xref: /petsc/src/mat/impls/sbaij/seq/sbaij.c (revision f5edf698d856978453db75eb7a802a869d6e16dd)
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;
3913f74950SBarry Smith   PetscInt       i,mbs = a->mbs;
4049b5e25fSSatish Balay 
4149b5e25fSSatish Balay   PetscFunctionBegin;
4249b5e25fSSatish Balay   if (a->diag) PetscFunctionReturn(0);
4349b5e25fSSatish Balay 
4413f74950SBarry Smith   ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&a->diag);CHKERRQ(ierr);
4552e6d16bSBarry Smith   ierr = PetscLogObjectMemory(A,(mbs+1)*sizeof(PetscInt));CHKERRQ(ierr);
46b424e231SHong Zhang   for (i=0; i<mbs; i++) a->diag[i] = a->i[i];
4749b5e25fSSatish Balay   PetscFunctionReturn(0);
4849b5e25fSSatish Balay }
4949b5e25fSSatish Balay 
504a2ae208SSatish Balay #undef __FUNCT__
514a2ae208SSatish Balay #define __FUNCT__ "MatGetRowIJ_SeqSBAIJ"
5213f74950SBarry Smith static PetscErrorCode MatGetRowIJ_SeqSBAIJ(Mat A,PetscInt oshift,PetscTruth symmetric,PetscInt *nn,PetscInt *ia[],PetscInt *ja[],PetscTruth *done)
5349b5e25fSSatish Balay {
54a6ece127SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
5513f74950SBarry Smith   PetscInt     n = a->mbs,i;
5649b5e25fSSatish Balay 
5749b5e25fSSatish Balay   PetscFunctionBegin;
58d3e5a4abSHong Zhang   *nn = n;
59a1373b80SHong Zhang   if (!ia) PetscFunctionReturn(0);
60a1373b80SHong Zhang 
61a6ece127SHong Zhang   if (oshift == 1) {
62a6ece127SHong Zhang     /* temporarily add 1 to i and j indices */
6313f74950SBarry Smith     PetscInt nz = a->i[n];
646c6c5352SBarry Smith     for (i=0; i<nz; i++) a->j[i]++;
65a1373b80SHong Zhang     for (i=0; i<n+1; i++) a->i[i]++;
66a1373b80SHong Zhang     *ia = a->i; *ja = a->j;
67a1373b80SHong Zhang   } else {
68a1373b80SHong Zhang     *ia = a->i; *ja = a->j;
69a6ece127SHong Zhang   }
7049b5e25fSSatish Balay   PetscFunctionReturn(0);
7149b5e25fSSatish Balay }
7249b5e25fSSatish Balay 
734a2ae208SSatish Balay #undef __FUNCT__
744a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRowIJ_SeqSBAIJ"
7513f74950SBarry Smith static PetscErrorCode MatRestoreRowIJ_SeqSBAIJ(Mat A,PetscInt oshift,PetscTruth symmetric,PetscInt *nn,PetscInt *ia[],PetscInt *ja[],PetscTruth *done)
7649b5e25fSSatish Balay {
77b7aaefc3SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
7813f74950SBarry Smith   PetscInt     i,n = a->mbs;
79a6ece127SHong Zhang 
8049b5e25fSSatish Balay   PetscFunctionBegin;
8149b5e25fSSatish Balay   if (!ia) PetscFunctionReturn(0);
82a6ece127SHong Zhang 
83a6ece127SHong Zhang   if (oshift == 1) {
8413f74950SBarry Smith     PetscInt nz = a->i[n]-1;
856c6c5352SBarry Smith     for (i=0; i<nz; i++) a->j[i]--;
86a6ece127SHong Zhang     for (i=0; i<n+1; i++) a->i[i]--;
87a6ece127SHong Zhang   }
88a6ece127SHong Zhang   PetscFunctionReturn(0);
8949b5e25fSSatish Balay }
9049b5e25fSSatish Balay 
914a2ae208SSatish Balay #undef __FUNCT__
924a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqSBAIJ"
93dfbe8321SBarry Smith PetscErrorCode MatDestroy_SeqSBAIJ(Mat A)
9449b5e25fSSatish Balay {
9549b5e25fSSatish Balay   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
96dfbe8321SBarry Smith   PetscErrorCode ierr;
9749b5e25fSSatish Balay 
9849b5e25fSSatish Balay   PetscFunctionBegin;
9977431f27SBarry Smith   PetscLogObjectState((PetscObject)A,"Rows=%D, NZ=%D",A->m,a->nz);
100085a36d4SBarry Smith   ierr = MatSeqXAIJFreeAIJ(a->singlemalloc,&a->a,&a->j,&a->i);CHKERRQ(ierr);
10149b5e25fSSatish Balay   if (a->row) {
10249b5e25fSSatish Balay     ierr = ISDestroy(a->row);CHKERRQ(ierr);
10349b5e25fSSatish Balay   }
10449b5e25fSSatish Balay   if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);}
105ab93d7beSBarry Smith   if (a->imax) {ierr = PetscFree2(a->imax,a->ilen);CHKERRQ(ierr);}
10649b5e25fSSatish Balay   if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);}
107d59c15a7SBarry Smith   if (a->solve_work)  {ierr = PetscFree(a->solve_work);CHKERRQ(ierr);}
108d59c15a7SBarry Smith   if (a->solves_work) {ierr = PetscFree(a->solves_work);CHKERRQ(ierr);}
109d59c15a7SBarry Smith   if (a->mult_work)   {ierr = PetscFree(a->mult_work);CHKERRQ(ierr);}
11049b5e25fSSatish Balay   if (a->saved_values) {ierr = PetscFree(a->saved_values);CHKERRQ(ierr);}
111c4319e64SHong Zhang   if (a->xtoy) {ierr = PetscFree(a->xtoy);CHKERRQ(ierr);}
1121a3463dfSHong Zhang 
1131a3463dfSHong Zhang   if (a->inew){
1141a3463dfSHong Zhang     ierr = PetscFree(a->inew);CHKERRQ(ierr);
1151a3463dfSHong Zhang     a->inew = 0;
1161a3463dfSHong Zhang   }
11749b5e25fSSatish Balay   ierr = PetscFree(a);CHKERRQ(ierr);
118901853e0SKris Buschelman 
119901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)A,"MatStoreValues_C","",PETSC_NULL);CHKERRQ(ierr);
120901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)A,"MatRetrieveValues_C","",PETSC_NULL);CHKERRQ(ierr);
121901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)A,"MatSeqSBAIJSetColumnIndices_C","",PETSC_NULL);CHKERRQ(ierr);
122901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqsbaij_seqaij_C","",PETSC_NULL);CHKERRQ(ierr);
123901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_seqsbaij_seqbaij_C","",PETSC_NULL);CHKERRQ(ierr);
124901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)A,"MatSeqSBAIJSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr);
12549b5e25fSSatish Balay   PetscFunctionReturn(0);
12649b5e25fSSatish Balay }
12749b5e25fSSatish Balay 
1284a2ae208SSatish Balay #undef __FUNCT__
1294a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqSBAIJ"
130dfbe8321SBarry Smith PetscErrorCode MatSetOption_SeqSBAIJ(Mat A,MatOption op)
13149b5e25fSSatish Balay {
132045c9aa0SHong Zhang   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
13363ba0a88SBarry Smith   PetscErrorCode ierr;
13449b5e25fSSatish Balay 
13549b5e25fSSatish Balay   PetscFunctionBegin;
1364d9d31abSKris Buschelman   switch (op) {
1374d9d31abSKris Buschelman   case MAT_ROW_ORIENTED:
1384d9d31abSKris Buschelman     a->roworiented = PETSC_TRUE;
1394d9d31abSKris Buschelman     break;
1404d9d31abSKris Buschelman   case MAT_COLUMN_ORIENTED:
1414d9d31abSKris Buschelman     a->roworiented = PETSC_FALSE;
1424d9d31abSKris Buschelman     break;
1434d9d31abSKris Buschelman   case MAT_COLUMNS_SORTED:
1444d9d31abSKris Buschelman     a->sorted = PETSC_TRUE;
1454d9d31abSKris Buschelman     break;
1464d9d31abSKris Buschelman   case MAT_COLUMNS_UNSORTED:
1474d9d31abSKris Buschelman     a->sorted = PETSC_FALSE;
1484d9d31abSKris Buschelman     break;
1494d9d31abSKris Buschelman   case MAT_KEEP_ZEROED_ROWS:
1504d9d31abSKris Buschelman     a->keepzeroedrows = PETSC_TRUE;
1514d9d31abSKris Buschelman     break;
1524d9d31abSKris Buschelman   case MAT_NO_NEW_NONZERO_LOCATIONS:
1534d9d31abSKris Buschelman     a->nonew = 1;
1544d9d31abSKris Buschelman     break;
1554d9d31abSKris Buschelman   case MAT_NEW_NONZERO_LOCATION_ERR:
1564d9d31abSKris Buschelman     a->nonew = -1;
1574d9d31abSKris Buschelman     break;
1584d9d31abSKris Buschelman   case MAT_NEW_NONZERO_ALLOCATION_ERR:
1594d9d31abSKris Buschelman     a->nonew = -2;
1604d9d31abSKris Buschelman     break;
1614d9d31abSKris Buschelman   case MAT_YES_NEW_NONZERO_LOCATIONS:
1624d9d31abSKris Buschelman     a->nonew = 0;
1634d9d31abSKris Buschelman     break;
1644d9d31abSKris Buschelman   case MAT_ROWS_SORTED:
1654d9d31abSKris Buschelman   case MAT_ROWS_UNSORTED:
1664d9d31abSKris Buschelman   case MAT_YES_NEW_DIAGONALS:
1674d9d31abSKris Buschelman   case MAT_IGNORE_OFF_PROC_ENTRIES:
1684d9d31abSKris Buschelman   case MAT_USE_HASH_TABLE:
16963ba0a88SBarry Smith     ierr = PetscLogInfo((A,"MatSetOption_SeqSBAIJ:Option ignored\n"));CHKERRQ(ierr);
1704d9d31abSKris Buschelman     break;
1714d9d31abSKris Buschelman   case MAT_NO_NEW_DIAGONALS:
17229bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS");
1739a4540c5SBarry Smith   case MAT_NOT_SYMMETRIC:
1749a4540c5SBarry Smith   case MAT_NOT_STRUCTURALLY_SYMMETRIC:
1759a4540c5SBarry Smith   case MAT_HERMITIAN:
1769a4540c5SBarry Smith     SETERRQ(PETSC_ERR_SUP,"Matrix must be symmetric");
17777e54ba9SKris Buschelman   case MAT_SYMMETRIC:
17877e54ba9SKris Buschelman   case MAT_STRUCTURALLY_SYMMETRIC:
1799a4540c5SBarry Smith   case MAT_NOT_HERMITIAN:
1809a4540c5SBarry Smith   case MAT_SYMMETRY_ETERNAL:
1819a4540c5SBarry Smith   case MAT_NOT_SYMMETRY_ETERNAL:
182941593c8SHong Zhang   case MAT_IGNORE_LOWER_TRIANGULAR:
183941593c8SHong Zhang     a->ignore_ltriangular = PETSC_TRUE;
184941593c8SHong Zhang     break;
185941593c8SHong Zhang   case MAT_ERROR_LOWER_TRIANGULAR:
186941593c8SHong Zhang     a->ignore_ltriangular = PETSC_FALSE;
18777e54ba9SKris Buschelman     break;
188*f5edf698SHong Zhang   case MAT_GETROW_UPPERTRIANGULAR:
189*f5edf698SHong Zhang     a->getrow_utriangular = PETSC_TRUE;
190*f5edf698SHong Zhang     break;
1914d9d31abSKris Buschelman   default:
19229bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"unknown option");
19349b5e25fSSatish Balay   }
19449b5e25fSSatish Balay   PetscFunctionReturn(0);
19549b5e25fSSatish Balay }
19649b5e25fSSatish Balay 
1974a2ae208SSatish Balay #undef __FUNCT__
1984a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqSBAIJ"
19913f74950SBarry Smith PetscErrorCode MatGetRow_SeqSBAIJ(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **v)
20049b5e25fSSatish Balay {
20149b5e25fSSatish Balay   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
2026849ba73SBarry Smith   PetscErrorCode ierr;
20313f74950SBarry Smith   PetscInt       itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*cols_i,bs2;
20449b5e25fSSatish Balay   MatScalar      *aa,*aa_i;
20587828ca2SBarry Smith   PetscScalar    *v_i;
20649b5e25fSSatish Balay 
20749b5e25fSSatish Balay   PetscFunctionBegin;
208*f5edf698SHong Zhang   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) or MatGetRowUpperTriangular()");
209*f5edf698SHong Zhang   /* Get the upper triangular part of the row */
210521d7252SBarry Smith   bs  = A->bs;
21149b5e25fSSatish Balay   ai  = a->i;
21249b5e25fSSatish Balay   aj  = a->j;
21349b5e25fSSatish Balay   aa  = a->a;
21449b5e25fSSatish Balay   bs2 = a->bs2;
21549b5e25fSSatish Balay 
21677431f27SBarry Smith   if (row < 0 || row >= A->m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE, "Row %D out of range", row);
21749b5e25fSSatish Balay 
21849b5e25fSSatish Balay   bn  = row/bs;   /* Block number */
21949b5e25fSSatish Balay   bp  = row % bs; /* Block position */
22049b5e25fSSatish Balay   M   = ai[bn+1] - ai[bn];
22149b5e25fSSatish Balay   *ncols = bs*M;
22249b5e25fSSatish Balay 
22349b5e25fSSatish Balay   if (v) {
22449b5e25fSSatish Balay     *v = 0;
22549b5e25fSSatish Balay     if (*ncols) {
22687828ca2SBarry Smith       ierr = PetscMalloc((*ncols+row)*sizeof(PetscScalar),v);CHKERRQ(ierr);
22749b5e25fSSatish Balay       for (i=0; i<M; i++) { /* for each block in the block row */
22849b5e25fSSatish Balay         v_i  = *v + i*bs;
22949b5e25fSSatish Balay         aa_i = aa + bs2*(ai[bn] + i);
23049b5e25fSSatish Balay         for (j=bp,k=0; j<bs2; j+=bs,k++) {v_i[k] = aa_i[j];}
23149b5e25fSSatish Balay       }
23249b5e25fSSatish Balay     }
23349b5e25fSSatish Balay   }
23449b5e25fSSatish Balay 
23549b5e25fSSatish Balay   if (cols) {
23649b5e25fSSatish Balay     *cols = 0;
23749b5e25fSSatish Balay     if (*ncols) {
23813f74950SBarry Smith       ierr = PetscMalloc((*ncols+row)*sizeof(PetscInt),cols);CHKERRQ(ierr);
23949b5e25fSSatish Balay       for (i=0; i<M; i++) { /* for each block in the block row */
24049b5e25fSSatish Balay         cols_i = *cols + i*bs;
24149b5e25fSSatish Balay         itmp  = bs*aj[ai[bn] + i];
24249b5e25fSSatish Balay         for (j=0; j<bs; j++) {cols_i[j] = itmp++;}
24349b5e25fSSatish Balay       }
24449b5e25fSSatish Balay     }
24549b5e25fSSatish Balay   }
24649b5e25fSSatish Balay 
24749b5e25fSSatish Balay   /*search column A(0:row-1,row) (=A(row,0:row-1)). Could be expensive! */
2485ddb2528SHong Zhang   /* this segment is currently removed, so only entries in the upper triangle are obtained */
2495ddb2528SHong Zhang #ifdef column_search
25049b5e25fSSatish Balay   v_i    = *v    + M*bs;
25149b5e25fSSatish Balay   cols_i = *cols + M*bs;
25249b5e25fSSatish Balay   for (i=0; i<bn; i++){ /* for each block row */
25349b5e25fSSatish Balay     M = ai[i+1] - ai[i];
25449b5e25fSSatish Balay     for (j=0; j<M; j++){
25549b5e25fSSatish Balay       itmp = aj[ai[i] + j];    /* block column value */
25649b5e25fSSatish Balay       if (itmp == bn){
25749b5e25fSSatish Balay         aa_i   = aa    + bs2*(ai[i] + j) + bs*bp;
25849b5e25fSSatish Balay         for (k=0; k<bs; k++) {
25949b5e25fSSatish Balay           *cols_i++ = i*bs+k;
26049b5e25fSSatish Balay           *v_i++    = aa_i[k];
26149b5e25fSSatish Balay         }
26249b5e25fSSatish Balay         *ncols += bs;
26349b5e25fSSatish Balay         break;
26449b5e25fSSatish Balay       }
26549b5e25fSSatish Balay     }
26649b5e25fSSatish Balay   }
2675ddb2528SHong Zhang #endif
26849b5e25fSSatish Balay   PetscFunctionReturn(0);
26949b5e25fSSatish Balay }
27049b5e25fSSatish Balay 
2714a2ae208SSatish Balay #undef __FUNCT__
2724a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqSBAIJ"
27313f74950SBarry Smith PetscErrorCode MatRestoreRow_SeqSBAIJ(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
27449b5e25fSSatish Balay {
275dfbe8321SBarry Smith   PetscErrorCode ierr;
27649b5e25fSSatish Balay 
27749b5e25fSSatish Balay   PetscFunctionBegin;
27849b5e25fSSatish Balay   if (idx) {if (*idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}}
27949b5e25fSSatish Balay   if (v)   {if (*v)   {ierr = PetscFree(*v);CHKERRQ(ierr);}}
28049b5e25fSSatish Balay   PetscFunctionReturn(0);
28149b5e25fSSatish Balay }
28249b5e25fSSatish Balay 
2834a2ae208SSatish Balay #undef __FUNCT__
284*f5edf698SHong Zhang #define __FUNCT__ "MatGetRowUpperTriangular_SeqSBAIJ"
285*f5edf698SHong Zhang PetscErrorCode MatGetRowUpperTriangular_SeqSBAIJ(Mat A)
286*f5edf698SHong Zhang {
287*f5edf698SHong Zhang   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
288*f5edf698SHong Zhang 
289*f5edf698SHong Zhang   PetscFunctionBegin;
290*f5edf698SHong Zhang   a->getrow_utriangular = PETSC_TRUE;
291*f5edf698SHong Zhang   PetscFunctionReturn(0);
292*f5edf698SHong Zhang }
293*f5edf698SHong Zhang #undef __FUNCT__
294*f5edf698SHong Zhang #define __FUNCT__ "MatRestoreRowUpperTriangular_SeqSBAIJ"
295*f5edf698SHong Zhang PetscErrorCode MatRestoreRowUpperTriangular_SeqSBAIJ(Mat A)
296*f5edf698SHong Zhang {
297*f5edf698SHong Zhang   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
298*f5edf698SHong Zhang 
299*f5edf698SHong Zhang   PetscFunctionBegin;
300*f5edf698SHong Zhang   a->getrow_utriangular = PETSC_FALSE;
301*f5edf698SHong Zhang   PetscFunctionReturn(0);
302*f5edf698SHong Zhang }
303*f5edf698SHong Zhang 
304*f5edf698SHong Zhang #undef __FUNCT__
3054a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqSBAIJ"
306dfbe8321SBarry Smith PetscErrorCode MatTranspose_SeqSBAIJ(Mat A,Mat *B)
30749b5e25fSSatish Balay {
308dfbe8321SBarry Smith   PetscErrorCode ierr;
30949b5e25fSSatish Balay   PetscFunctionBegin;
310999d9058SBarry Smith   ierr = MatDuplicate(A,MAT_COPY_VALUES,B);CHKERRQ(ierr);
3118115998fSBarry Smith   PetscFunctionReturn(0);
31249b5e25fSSatish Balay }
31349b5e25fSSatish Balay 
3144a2ae208SSatish Balay #undef __FUNCT__
3154a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_ASCII"
3166849ba73SBarry Smith static PetscErrorCode MatView_SeqSBAIJ_ASCII(Mat A,PetscViewer viewer)
31749b5e25fSSatish Balay {
31849b5e25fSSatish Balay   Mat_SeqSBAIJ      *a = (Mat_SeqSBAIJ*)A->data;
319dfbe8321SBarry Smith   PetscErrorCode    ierr;
320521d7252SBarry Smith   PetscInt          i,j,bs = A->bs,k,l,bs2=a->bs2;
3212dcb1b2aSMatthew Knepley   const char        *name;
322f3ef73ceSBarry Smith   PetscViewerFormat format;
32349b5e25fSSatish Balay 
32449b5e25fSSatish Balay   PetscFunctionBegin;
32580fe4e49SBarry Smith   ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
326b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
327456192e2SBarry Smith   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
32877431f27SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  block size is %D\n",bs);CHKERRQ(ierr);
329fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_MATLAB) {
33029bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Matlab format not supported");
331fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_COMMON) {
332b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
33349b5e25fSSatish Balay     for (i=0; i<a->mbs; i++) {
33449b5e25fSSatish Balay       for (j=0; j<bs; j++) {
33577431f27SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i*bs+j);CHKERRQ(ierr);
33649b5e25fSSatish Balay         for (k=a->i[i]; k<a->i[i+1]; k++) {
33749b5e25fSSatish Balay           for (l=0; l<bs; l++) {
33849b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX)
33949b5e25fSSatish Balay             if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
34077431f27SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%D, %g + %g i) ",bs*a->j[k]+l,
34149b5e25fSSatish Balay                                             PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
34249b5e25fSSatish Balay             } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
34377431f27SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%D, %g - %g i) ",bs*a->j[k]+l,
34449b5e25fSSatish Balay                                             PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
34549b5e25fSSatish Balay             } else if (PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
34677431f27SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%D, %g) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
34749b5e25fSSatish Balay             }
34849b5e25fSSatish Balay #else
34949b5e25fSSatish Balay             if (a->a[bs2*k + l*bs + j] != 0.0) {
35077431f27SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%D, %g) ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr);
35149b5e25fSSatish Balay             }
35249b5e25fSSatish Balay #endif
35349b5e25fSSatish Balay           }
35449b5e25fSSatish Balay         }
355b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
35649b5e25fSSatish Balay       }
35749b5e25fSSatish Balay     }
358b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
359c1490034SHong Zhang   } else if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) {
360c1490034SHong Zhang      PetscFunctionReturn(0);
36149b5e25fSSatish Balay   } else {
362b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
36349b5e25fSSatish Balay     for (i=0; i<a->mbs; i++) {
36449b5e25fSSatish Balay       for (j=0; j<bs; j++) {
36577431f27SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i*bs+j);CHKERRQ(ierr);
36649b5e25fSSatish Balay         for (k=a->i[i]; k<a->i[i+1]; k++) {
36749b5e25fSSatish Balay           for (l=0; l<bs; l++) {
36849b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX)
36949b5e25fSSatish Balay             if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0) {
37077431f27SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%D, %g + %g i) ",bs*a->j[k]+l,
37149b5e25fSSatish Balay                                             PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
37249b5e25fSSatish Balay             } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0) {
37377431f27SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%D, %g - %g i) ",bs*a->j[k]+l,
37449b5e25fSSatish Balay                                             PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
37549b5e25fSSatish Balay             } else {
37677431f27SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%D, %g) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
37749b5e25fSSatish Balay             }
37849b5e25fSSatish Balay #else
37977431f27SBarry Smith             ierr = PetscViewerASCIIPrintf(viewer," %D %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr);
38049b5e25fSSatish Balay #endif
38149b5e25fSSatish Balay           }
38249b5e25fSSatish Balay         }
383b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
38449b5e25fSSatish Balay       }
38549b5e25fSSatish Balay     }
386b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
38749b5e25fSSatish Balay   }
388b0a32e0cSBarry Smith   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
38949b5e25fSSatish Balay   PetscFunctionReturn(0);
39049b5e25fSSatish Balay }
39149b5e25fSSatish Balay 
3924a2ae208SSatish Balay #undef __FUNCT__
3934a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Draw_Zoom"
3946849ba73SBarry Smith static PetscErrorCode MatView_SeqSBAIJ_Draw_Zoom(PetscDraw draw,void *Aa)
39549b5e25fSSatish Balay {
39649b5e25fSSatish Balay   Mat            A = (Mat) Aa;
39749b5e25fSSatish Balay   Mat_SeqSBAIJ   *a=(Mat_SeqSBAIJ*)A->data;
3986849ba73SBarry Smith   PetscErrorCode ierr;
39913f74950SBarry Smith   PetscInt       row,i,j,k,l,mbs=a->mbs,color,bs=A->bs,bs2=a->bs2;
40013f74950SBarry Smith   PetscMPIInt    rank;
40149b5e25fSSatish Balay   PetscReal      xl,yl,xr,yr,x_l,x_r,y_l,y_r;
40249b5e25fSSatish Balay   MatScalar      *aa;
40349b5e25fSSatish Balay   MPI_Comm       comm;
404b0a32e0cSBarry Smith   PetscViewer    viewer;
40549b5e25fSSatish Balay 
40649b5e25fSSatish Balay   PetscFunctionBegin;
40749b5e25fSSatish Balay   /*
40849b5e25fSSatish Balay     This is nasty. If this is called from an originally parallel matrix
40949b5e25fSSatish Balay     then all processes call this,but only the first has the matrix so the
41049b5e25fSSatish Balay     rest should return immediately.
41149b5e25fSSatish Balay   */
41249b5e25fSSatish Balay   ierr = PetscObjectGetComm((PetscObject)draw,&comm);CHKERRQ(ierr);
41349b5e25fSSatish Balay   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
41449b5e25fSSatish Balay   if (rank) PetscFunctionReturn(0);
41549b5e25fSSatish Balay 
41649b5e25fSSatish Balay   ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr);
41749b5e25fSSatish Balay 
418b0a32e0cSBarry Smith   ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr);
419b0a32e0cSBarry Smith   PetscDrawString(draw, .3*(xl+xr), .3*(yl+yr), PETSC_DRAW_BLACK, "symmetric");
42049b5e25fSSatish Balay 
42149b5e25fSSatish Balay   /* loop over matrix elements drawing boxes */
422b0a32e0cSBarry Smith   color = PETSC_DRAW_BLUE;
42349b5e25fSSatish Balay   for (i=0,row=0; i<mbs; i++,row+=bs) {
42449b5e25fSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
425c464158bSHong Zhang       y_l = A->m - row - 1.0; y_r = y_l + 1.0;
42649b5e25fSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
42749b5e25fSSatish Balay       aa = a->a + j*bs2;
42849b5e25fSSatish Balay       for (k=0; k<bs; k++) {
42949b5e25fSSatish Balay         for (l=0; l<bs; l++) {
43049b5e25fSSatish Balay           if (PetscRealPart(*aa++) >=  0.) continue;
431b0a32e0cSBarry Smith           ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
43249b5e25fSSatish Balay         }
43349b5e25fSSatish Balay       }
43449b5e25fSSatish Balay     }
43549b5e25fSSatish Balay   }
436b0a32e0cSBarry Smith   color = PETSC_DRAW_CYAN;
43749b5e25fSSatish Balay   for (i=0,row=0; i<mbs; i++,row+=bs) {
43849b5e25fSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
439c464158bSHong Zhang       y_l = A->m - row - 1.0; y_r = y_l + 1.0;
44049b5e25fSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
44149b5e25fSSatish Balay       aa = a->a + j*bs2;
44249b5e25fSSatish Balay       for (k=0; k<bs; k++) {
44349b5e25fSSatish Balay         for (l=0; l<bs; l++) {
44449b5e25fSSatish Balay           if (PetscRealPart(*aa++) != 0.) continue;
445b0a32e0cSBarry Smith           ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
44649b5e25fSSatish Balay         }
44749b5e25fSSatish Balay       }
44849b5e25fSSatish Balay     }
44949b5e25fSSatish Balay   }
45049b5e25fSSatish Balay 
451b0a32e0cSBarry Smith   color = PETSC_DRAW_RED;
45249b5e25fSSatish Balay   for (i=0,row=0; i<mbs; i++,row+=bs) {
45349b5e25fSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
454c464158bSHong Zhang       y_l = A->m - row - 1.0; y_r = y_l + 1.0;
45549b5e25fSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
45649b5e25fSSatish Balay       aa = a->a + j*bs2;
45749b5e25fSSatish Balay       for (k=0; k<bs; k++) {
45849b5e25fSSatish Balay         for (l=0; l<bs; l++) {
45949b5e25fSSatish Balay           if (PetscRealPart(*aa++) <= 0.) continue;
460b0a32e0cSBarry Smith           ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
46149b5e25fSSatish Balay         }
46249b5e25fSSatish Balay       }
46349b5e25fSSatish Balay     }
46449b5e25fSSatish Balay   }
46549b5e25fSSatish Balay   PetscFunctionReturn(0);
46649b5e25fSSatish Balay }
46749b5e25fSSatish Balay 
4684a2ae208SSatish Balay #undef __FUNCT__
4694a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Draw"
4706849ba73SBarry Smith static PetscErrorCode MatView_SeqSBAIJ_Draw(Mat A,PetscViewer viewer)
47149b5e25fSSatish Balay {
472dfbe8321SBarry Smith   PetscErrorCode ierr;
47349b5e25fSSatish Balay   PetscReal      xl,yl,xr,yr,w,h;
474b0a32e0cSBarry Smith   PetscDraw      draw;
47549b5e25fSSatish Balay   PetscTruth     isnull;
47649b5e25fSSatish Balay 
47749b5e25fSSatish Balay   PetscFunctionBegin;
47849b5e25fSSatish Balay 
479b0a32e0cSBarry Smith   ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
480b0a32e0cSBarry Smith   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
48149b5e25fSSatish Balay 
48249b5e25fSSatish Balay   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
483c464158bSHong Zhang   xr  = A->m; yr = A->m; h = yr/10.0; w = xr/10.0;
48449b5e25fSSatish Balay   xr += w;    yr += h;  xl = -w;     yl = -h;
485b0a32e0cSBarry Smith   ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
486b0a32e0cSBarry Smith   ierr = PetscDrawZoom(draw,MatView_SeqSBAIJ_Draw_Zoom,A);CHKERRQ(ierr);
48749b5e25fSSatish Balay   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr);
48849b5e25fSSatish Balay   PetscFunctionReturn(0);
48949b5e25fSSatish Balay }
49049b5e25fSSatish Balay 
4914a2ae208SSatish Balay #undef __FUNCT__
4924a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ"
493dfbe8321SBarry Smith PetscErrorCode MatView_SeqSBAIJ(Mat A,PetscViewer viewer)
49449b5e25fSSatish Balay {
495dfbe8321SBarry Smith   PetscErrorCode ierr;
49632077d6dSBarry Smith   PetscTruth     iascii,isdraw;
49749b5e25fSSatish Balay 
49849b5e25fSSatish Balay   PetscFunctionBegin;
49932077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
500fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
50132077d6dSBarry Smith   if (iascii){
50249b5e25fSSatish Balay     ierr = MatView_SeqSBAIJ_ASCII(A,viewer);CHKERRQ(ierr);
50349b5e25fSSatish Balay   } else if (isdraw) {
50449b5e25fSSatish Balay     ierr = MatView_SeqSBAIJ_Draw(A,viewer);CHKERRQ(ierr);
50549b5e25fSSatish Balay   } else {
506a5e6ed63SBarry Smith     Mat B;
507ceb03754SKris Buschelman     ierr = MatConvert(A,MATSEQAIJ,MAT_INITIAL_MATRIX,&B);CHKERRQ(ierr);
508a5e6ed63SBarry Smith     ierr = MatView(B,viewer);CHKERRQ(ierr);
509a5e6ed63SBarry Smith     ierr = MatDestroy(B);CHKERRQ(ierr);
51049b5e25fSSatish Balay   }
51149b5e25fSSatish Balay   PetscFunctionReturn(0);
51249b5e25fSSatish Balay }
51349b5e25fSSatish Balay 
51449b5e25fSSatish Balay 
5154a2ae208SSatish Balay #undef __FUNCT__
5164a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqSBAIJ"
51713f74950SBarry Smith PetscErrorCode MatGetValues_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],PetscScalar v[])
51849b5e25fSSatish Balay {
519045c9aa0SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
52013f74950SBarry Smith   PetscInt     *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j;
52113f74950SBarry Smith   PetscInt     *ai = a->i,*ailen = a->ilen;
52213f74950SBarry Smith   PetscInt     brow,bcol,ridx,cidx,bs=A->bs,bs2=a->bs2;
52349b5e25fSSatish Balay   MatScalar    *ap,*aa = a->a,zero = 0.0;
52449b5e25fSSatish Balay 
52549b5e25fSSatish Balay   PetscFunctionBegin;
52649b5e25fSSatish Balay   for (k=0; k<m; k++) { /* loop over rows */
52749b5e25fSSatish Balay     row  = im[k]; brow = row/bs;
52877431f27SBarry Smith     if (row < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative row: %D",row);
52977431f27SBarry Smith     if (row >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->m-1);
53049b5e25fSSatish Balay     rp   = aj + ai[brow] ; ap = aa + bs2*ai[brow] ;
53149b5e25fSSatish Balay     nrow = ailen[brow];
53249b5e25fSSatish Balay     for (l=0; l<n; l++) { /* loop over columns */
53377431f27SBarry Smith       if (in[l] < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative column: %D",in[l]);
53477431f27SBarry Smith       if (in[l] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",in[l],A->n-1);
53549b5e25fSSatish Balay       col  = in[l] ;
53649b5e25fSSatish Balay       bcol = col/bs;
53749b5e25fSSatish Balay       cidx = col%bs;
53849b5e25fSSatish Balay       ridx = row%bs;
53949b5e25fSSatish Balay       high = nrow;
54049b5e25fSSatish Balay       low  = 0; /* assume unsorted */
54149b5e25fSSatish Balay       while (high-low > 5) {
54249b5e25fSSatish Balay         t = (low+high)/2;
54349b5e25fSSatish Balay         if (rp[t] > bcol) high = t;
54449b5e25fSSatish Balay         else             low  = t;
54549b5e25fSSatish Balay       }
54649b5e25fSSatish Balay       for (i=low; i<high; i++) {
54749b5e25fSSatish Balay         if (rp[i] > bcol) break;
54849b5e25fSSatish Balay         if (rp[i] == bcol) {
54949b5e25fSSatish Balay           *v++ = ap[bs2*i+bs*cidx+ridx];
55049b5e25fSSatish Balay           goto finished;
55149b5e25fSSatish Balay         }
55249b5e25fSSatish Balay       }
55349b5e25fSSatish Balay       *v++ = zero;
55449b5e25fSSatish Balay        finished:;
55549b5e25fSSatish Balay     }
55649b5e25fSSatish Balay   }
55749b5e25fSSatish Balay   PetscFunctionReturn(0);
55849b5e25fSSatish Balay }
55949b5e25fSSatish Balay 
56049b5e25fSSatish Balay 
5614a2ae208SSatish Balay #undef __FUNCT__
5624a2ae208SSatish Balay #define __FUNCT__ "MatSetValuesBlocked_SeqSBAIJ"
56313f74950SBarry Smith PetscErrorCode MatSetValuesBlocked_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is)
56449b5e25fSSatish Balay {
5650880e062SHong Zhang   Mat_SeqSBAIJ    *a = (Mat_SeqSBAIJ*)A->data;
5666849ba73SBarry Smith   PetscErrorCode  ierr;
567e2ee6c50SBarry Smith   PetscInt        *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,lastcol = -1;
56813f74950SBarry Smith   PetscInt        *imax=a->imax,*ai=a->i,*ailen=a->ilen;
56913f74950SBarry Smith   PetscInt        *aj=a->j,nonew=a->nonew,bs2=a->bs2,bs=A->bs,stepval;
5700880e062SHong Zhang   PetscTruth      roworiented=a->roworiented;
571f15d580aSBarry Smith   const MatScalar *value = v;
572f15d580aSBarry Smith   MatScalar       *ap,*aa = a->a,*bap;
5730880e062SHong Zhang 
57449b5e25fSSatish Balay   PetscFunctionBegin;
5750880e062SHong Zhang   if (roworiented) {
5760880e062SHong Zhang     stepval = (n-1)*bs;
5770880e062SHong Zhang   } else {
5780880e062SHong Zhang     stepval = (m-1)*bs;
5790880e062SHong Zhang   }
5800880e062SHong Zhang   for (k=0; k<m; k++) { /* loop over added rows */
5810880e062SHong Zhang     row  = im[k];
5820880e062SHong Zhang     if (row < 0) continue;
5832515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
58477431f27SBarry Smith     if (row >= a->mbs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,a->mbs-1);
5850880e062SHong Zhang #endif
5860880e062SHong Zhang     rp   = aj + ai[row];
5870880e062SHong Zhang     ap   = aa + bs2*ai[row];
5880880e062SHong Zhang     rmax = imax[row];
5890880e062SHong Zhang     nrow = ailen[row];
5900880e062SHong Zhang     low  = 0;
591818f2c47SBarry Smith     high = nrow;
5920880e062SHong Zhang     for (l=0; l<n; l++) { /* loop over added columns */
5930880e062SHong Zhang       if (in[l] < 0) continue;
5940880e062SHong Zhang       col = in[l];
5952515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
59677431f27SBarry Smith       if (col >= a->nbs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",col,a->nbs-1);
597b1823623SSatish Balay #endif
5980880e062SHong Zhang       if (col < row) continue; /* ignore lower triangular block */
5990880e062SHong Zhang       if (roworiented) {
6000880e062SHong Zhang         value = v + k*(stepval+bs)*bs + l*bs;
6010880e062SHong Zhang       } else {
6020880e062SHong Zhang         value = v + l*(stepval+bs)*bs + k*bs;
6030880e062SHong Zhang       }
6047cd84e04SBarry Smith       if (col <= lastcol) low = 0; else high = nrow;
605e2ee6c50SBarry Smith       lastcol = col;
6060880e062SHong Zhang       while (high-low > 7) {
6070880e062SHong Zhang         t = (low+high)/2;
6080880e062SHong Zhang         if (rp[t] > col) high = t;
6090880e062SHong Zhang         else             low  = t;
6100880e062SHong Zhang       }
6110880e062SHong Zhang       for (i=low; i<high; i++) {
6120880e062SHong Zhang         if (rp[i] > col) break;
6130880e062SHong Zhang         if (rp[i] == col) {
6140880e062SHong Zhang           bap  = ap +  bs2*i;
6150880e062SHong Zhang           if (roworiented) {
6160880e062SHong Zhang             if (is == ADD_VALUES) {
6170880e062SHong Zhang               for (ii=0; ii<bs; ii++,value+=stepval) {
6180880e062SHong Zhang                 for (jj=ii; jj<bs2; jj+=bs) {
6190880e062SHong Zhang                   bap[jj] += *value++;
6200880e062SHong Zhang                 }
6210880e062SHong Zhang               }
6220880e062SHong Zhang             } else {
6230880e062SHong Zhang               for (ii=0; ii<bs; ii++,value+=stepval) {
6240880e062SHong Zhang                 for (jj=ii; jj<bs2; jj+=bs) {
6250880e062SHong Zhang                   bap[jj] = *value++;
6260880e062SHong Zhang                 }
6270880e062SHong Zhang                }
6280880e062SHong Zhang             }
6290880e062SHong Zhang           } else {
6300880e062SHong Zhang             if (is == ADD_VALUES) {
6310880e062SHong Zhang               for (ii=0; ii<bs; ii++,value+=stepval) {
6320880e062SHong Zhang                 for (jj=0; jj<bs; jj++) {
6330880e062SHong Zhang                   *bap++ += *value++;
6340880e062SHong Zhang                 }
6350880e062SHong Zhang               }
6360880e062SHong Zhang             } else {
6370880e062SHong Zhang               for (ii=0; ii<bs; ii++,value+=stepval) {
6380880e062SHong Zhang                 for (jj=0; jj<bs; jj++) {
6390880e062SHong Zhang                   *bap++  = *value++;
6400880e062SHong Zhang                 }
6410880e062SHong Zhang               }
6420880e062SHong Zhang             }
6430880e062SHong Zhang           }
6440880e062SHong Zhang           goto noinsert2;
6450880e062SHong Zhang         }
6460880e062SHong Zhang       }
6470880e062SHong Zhang       if (nonew == 1) goto noinsert2;
648085a36d4SBarry Smith       if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col);
649ed1caa07SMatthew Knepley       MatSeqXAIJReallocateAIJ(a,bs2,nrow,row,col,rmax,aa,ai,aj,a->mbs,rp,ap,imax,nonew);
650c03d1d03SSatish Balay       N = nrow++ - 1; high++;
6510880e062SHong Zhang       /* shift up all the later entries in this row */
6520880e062SHong Zhang       for (ii=N; ii>=i; ii--) {
6530880e062SHong Zhang         rp[ii+1] = rp[ii];
6540880e062SHong Zhang         ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr);
6550880e062SHong Zhang       }
6560880e062SHong Zhang       if (N >= i) {
6570880e062SHong Zhang         ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr);
6580880e062SHong Zhang       }
6590880e062SHong Zhang       rp[i] = col;
6600880e062SHong Zhang       bap   = ap +  bs2*i;
6610880e062SHong Zhang       if (roworiented) {
6620880e062SHong Zhang         for (ii=0; ii<bs; ii++,value+=stepval) {
6630880e062SHong Zhang           for (jj=ii; jj<bs2; jj+=bs) {
6640880e062SHong Zhang             bap[jj] = *value++;
6650880e062SHong Zhang           }
6660880e062SHong Zhang         }
6670880e062SHong Zhang       } else {
6680880e062SHong Zhang         for (ii=0; ii<bs; ii++,value+=stepval) {
6690880e062SHong Zhang           for (jj=0; jj<bs; jj++) {
6700880e062SHong Zhang             *bap++  = *value++;
6710880e062SHong Zhang           }
6720880e062SHong Zhang         }
6730880e062SHong Zhang        }
6740880e062SHong Zhang     noinsert2:;
6750880e062SHong Zhang       low = i;
6760880e062SHong Zhang     }
6770880e062SHong Zhang     ailen[row] = nrow;
6780880e062SHong Zhang   }
6790880e062SHong Zhang    PetscFunctionReturn(0);
68049b5e25fSSatish Balay }
68149b5e25fSSatish Balay 
6824a2ae208SSatish Balay #undef __FUNCT__
6834a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_SeqSBAIJ"
684dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_SeqSBAIJ(Mat A,MatAssemblyType mode)
68549b5e25fSSatish Balay {
68649b5e25fSSatish Balay   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
6876849ba73SBarry Smith   PetscErrorCode ierr;
68813f74950SBarry Smith   PetscInt       fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax;
68913f74950SBarry Smith   PetscInt       m = A->m,*ip,N,*ailen = a->ilen;
69013f74950SBarry Smith   PetscInt       mbs = a->mbs,bs2 = a->bs2,rmax = 0;
69149b5e25fSSatish Balay   MatScalar      *aa = a->a,*ap;
69249b5e25fSSatish Balay 
69349b5e25fSSatish Balay   PetscFunctionBegin;
69449b5e25fSSatish Balay   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0);
69549b5e25fSSatish Balay 
69649b5e25fSSatish Balay   if (m) rmax = ailen[0];
69749b5e25fSSatish Balay   for (i=1; i<mbs; i++) {
69849b5e25fSSatish Balay     /* move each row back by the amount of empty slots (fshift) before it*/
69949b5e25fSSatish Balay     fshift += imax[i-1] - ailen[i-1];
70049b5e25fSSatish Balay      rmax   = PetscMax(rmax,ailen[i]);
70149b5e25fSSatish Balay      if (fshift) {
70249b5e25fSSatish Balay        ip = aj + ai[i]; ap = aa + bs2*ai[i];
70349b5e25fSSatish Balay        N = ailen[i];
70449b5e25fSSatish Balay        for (j=0; j<N; j++) {
70549b5e25fSSatish Balay          ip[j-fshift] = ip[j];
70649b5e25fSSatish Balay          ierr = PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
70749b5e25fSSatish Balay        }
70849b5e25fSSatish Balay      }
70949b5e25fSSatish Balay      ai[i] = ai[i-1] + ailen[i-1];
71049b5e25fSSatish Balay   }
71149b5e25fSSatish Balay   if (mbs) {
71249b5e25fSSatish Balay     fshift += imax[mbs-1] - ailen[mbs-1];
71349b5e25fSSatish Balay      ai[mbs] = ai[mbs-1] + ailen[mbs-1];
71449b5e25fSSatish Balay   }
71549b5e25fSSatish Balay   /* reset ilen and imax for each row */
71649b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
71749b5e25fSSatish Balay     ailen[i] = imax[i] = ai[i+1] - ai[i];
71849b5e25fSSatish Balay   }
7196c6c5352SBarry Smith   a->nz = ai[mbs];
72049b5e25fSSatish Balay 
721b424e231SHong Zhang   /* diagonals may have moved, reset it */
722b424e231SHong Zhang   if (a->diag) {
72313f74950SBarry Smith     ierr = PetscMemcpy(a->diag,ai,(mbs+1)*sizeof(PetscInt));CHKERRQ(ierr);
72449b5e25fSSatish Balay   }
72563ba0a88SBarry Smith   ierr = PetscLogInfo((A,"MatAssemblyEnd_SeqSBAIJ:Matrix size: %D X %D, block size %D; storage space: %D unneeded, %D used\n",
72663ba0a88SBarry Smith                        m,A->m,A->bs,fshift*bs2,a->nz*bs2));CHKERRQ(ierr);
72763ba0a88SBarry Smith   ierr = PetscLogInfo((A,"MatAssemblyEnd_SeqSBAIJ:Number of mallocs during MatSetValues is %D\n",a->reallocs));CHKERRQ(ierr);
72863ba0a88SBarry Smith   ierr = PetscLogInfo((A,"MatAssemblyEnd_SeqSBAIJ:Most nonzeros blocks in any row is %D\n",rmax));CHKERRQ(ierr);
72949b5e25fSSatish Balay   a->reallocs          = 0;
73049b5e25fSSatish Balay   A->info.nz_unneeded  = (PetscReal)fshift*bs2;
73149b5e25fSSatish Balay   PetscFunctionReturn(0);
73249b5e25fSSatish Balay }
73349b5e25fSSatish Balay 
73449b5e25fSSatish Balay /*
73549b5e25fSSatish Balay    This function returns an array of flags which indicate the locations of contiguous
73649b5e25fSSatish Balay    blocks that should be zeroed. for eg: if bs = 3  and is = [0,1,2,3,5,6,7,8,9]
73749b5e25fSSatish Balay    then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)]
73849b5e25fSSatish Balay    Assume: sizes should be long enough to hold all the values.
73949b5e25fSSatish Balay */
7404a2ae208SSatish Balay #undef __FUNCT__
7414a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqSBAIJ_Check_Blocks"
74213f74950SBarry Smith PetscErrorCode MatZeroRows_SeqSBAIJ_Check_Blocks(PetscInt idx[],PetscInt n,PetscInt bs,PetscInt sizes[], PetscInt *bs_max)
74349b5e25fSSatish Balay {
74413f74950SBarry Smith   PetscInt   i,j,k,row;
74549b5e25fSSatish Balay   PetscTruth flg;
74649b5e25fSSatish Balay 
74749b5e25fSSatish Balay   PetscFunctionBegin;
74849b5e25fSSatish Balay    for (i=0,j=0; i<n; j++) {
74949b5e25fSSatish Balay      row = idx[i];
75049b5e25fSSatish Balay      if (row%bs!=0) { /* Not the begining of a block */
75149b5e25fSSatish Balay        sizes[j] = 1;
75249b5e25fSSatish Balay        i++;
75349b5e25fSSatish Balay      } else if (i+bs > n) { /* Beginning of a block, but complete block doesn't exist (at idx end) */
75449b5e25fSSatish Balay        sizes[j] = 1;         /* Also makes sure atleast 'bs' values exist for next else */
75549b5e25fSSatish Balay        i++;
75649b5e25fSSatish Balay      } else { /* Begining of the block, so check if the complete block exists */
75749b5e25fSSatish Balay        flg = PETSC_TRUE;
75849b5e25fSSatish Balay        for (k=1; k<bs; k++) {
75949b5e25fSSatish Balay          if (row+k != idx[i+k]) { /* break in the block */
76049b5e25fSSatish Balay            flg = PETSC_FALSE;
76149b5e25fSSatish Balay            break;
76249b5e25fSSatish Balay          }
76349b5e25fSSatish Balay        }
764abc0a331SBarry Smith        if (flg) { /* No break in the bs */
76549b5e25fSSatish Balay          sizes[j] = bs;
76649b5e25fSSatish Balay          i+= bs;
76749b5e25fSSatish Balay        } else {
76849b5e25fSSatish Balay          sizes[j] = 1;
76949b5e25fSSatish Balay          i++;
77049b5e25fSSatish Balay        }
77149b5e25fSSatish Balay      }
77249b5e25fSSatish Balay    }
77349b5e25fSSatish Balay    *bs_max = j;
77449b5e25fSSatish Balay    PetscFunctionReturn(0);
77549b5e25fSSatish Balay }
77649b5e25fSSatish Balay 
77749b5e25fSSatish Balay 
77849b5e25fSSatish Balay /* Only add/insert a(i,j) with i<=j (blocks).
77949b5e25fSSatish Balay    Any a(i,j) with i>j input by user is ingored.
78049b5e25fSSatish Balay */
78149b5e25fSSatish Balay 
7824a2ae208SSatish Balay #undef __FUNCT__
7834a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqSBAIJ"
78413f74950SBarry Smith PetscErrorCode MatSetValues_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is)
78549b5e25fSSatish Balay {
78649b5e25fSSatish Balay   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
7876849ba73SBarry Smith   PetscErrorCode ierr;
788e2ee6c50SBarry Smith   PetscInt       *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,lastcol = -1;
78913f74950SBarry Smith   PetscInt       *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented;
79013f74950SBarry Smith   PetscInt       *aj=a->j,nonew=a->nonew,bs=A->bs,brow,bcol;
79113f74950SBarry Smith   PetscInt       ridx,cidx,bs2=a->bs2;
79249b5e25fSSatish Balay   MatScalar      *ap,value,*aa=a->a,*bap;
79349b5e25fSSatish Balay 
79449b5e25fSSatish Balay   PetscFunctionBegin;
79549b5e25fSSatish Balay   for (k=0; k<m; k++) { /* loop over added rows */
79649b5e25fSSatish Balay     row  = im[k];       /* row number */
79749b5e25fSSatish Balay     brow = row/bs;      /* block row number */
79849b5e25fSSatish Balay     if (row < 0) continue;
7992515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
80077431f27SBarry Smith     if (row >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->m-1);
80149b5e25fSSatish Balay #endif
80249b5e25fSSatish Balay     rp   = aj + ai[brow]; /*ptr to beginning of column value of the row block*/
80349b5e25fSSatish Balay     ap   = aa + bs2*ai[brow]; /*ptr to beginning of element value of the row block*/
80449b5e25fSSatish Balay     rmax = imax[brow];  /* maximum space allocated for this row */
80549b5e25fSSatish Balay     nrow = ailen[brow]; /* actual length of this row */
80649b5e25fSSatish Balay     low  = 0;
80749b5e25fSSatish Balay 
80849b5e25fSSatish Balay     for (l=0; l<n; l++) { /* loop over added columns */
80949b5e25fSSatish Balay       if (in[l] < 0) continue;
8102515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
81177431f27SBarry Smith       if (in[l] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",in[l],A->m-1);
81249b5e25fSSatish Balay #endif
81349b5e25fSSatish Balay       col = in[l];
81449b5e25fSSatish Balay       bcol = col/bs;              /* block col number */
81549b5e25fSSatish Balay 
816941593c8SHong Zhang       if (brow > bcol) {
817941593c8SHong Zhang         if (a->ignore_ltriangular){
818941593c8SHong Zhang           continue; /* ignore lower triangular values */
819941593c8SHong Zhang         } else {
820941593c8SHong Zhang           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)");
821941593c8SHong Zhang         }
822941593c8SHong Zhang       }
823f4989cb3SHong Zhang 
82449b5e25fSSatish Balay       ridx = row % bs; cidx = col % bs; /*row and col index inside the block */
8258549e402SHong Zhang       if ((brow==bcol && ridx<=cidx) || (brow<bcol)){
82649b5e25fSSatish Balay         /* element value a(k,l) */
82749b5e25fSSatish Balay         if (roworiented) {
82849b5e25fSSatish Balay           value = v[l + k*n];
82949b5e25fSSatish Balay         } else {
83049b5e25fSSatish Balay           value = v[k + l*m];
83149b5e25fSSatish Balay         }
83249b5e25fSSatish Balay 
83349b5e25fSSatish Balay         /* move pointer bap to a(k,l) quickly and add/insert value */
8347cd84e04SBarry Smith         if (col <= lastcol) low = 0; high = nrow;
835e2ee6c50SBarry Smith         lastcol = col;
83649b5e25fSSatish Balay         while (high-low > 7) {
83749b5e25fSSatish Balay           t = (low+high)/2;
83849b5e25fSSatish Balay           if (rp[t] > bcol) high = t;
83949b5e25fSSatish Balay           else              low  = t;
84049b5e25fSSatish Balay         }
84149b5e25fSSatish Balay         for (i=low; i<high; i++) {
84249b5e25fSSatish Balay           if (rp[i] > bcol) break;
84349b5e25fSSatish Balay           if (rp[i] == bcol) {
84449b5e25fSSatish Balay             bap  = ap +  bs2*i + bs*cidx + ridx;
84549b5e25fSSatish Balay             if (is == ADD_VALUES) *bap += value;
84649b5e25fSSatish Balay             else                  *bap  = value;
8478549e402SHong Zhang             /* for diag block, add/insert its symmetric element a(cidx,ridx) */
8488549e402SHong Zhang             if (brow == bcol && ridx < cidx){
8498549e402SHong Zhang               bap  = ap +  bs2*i + bs*ridx + cidx;
8508549e402SHong Zhang               if (is == ADD_VALUES) *bap += value;
8518549e402SHong Zhang               else                  *bap  = value;
8528549e402SHong Zhang             }
85349b5e25fSSatish Balay             goto noinsert1;
85449b5e25fSSatish Balay           }
85549b5e25fSSatish Balay         }
85649b5e25fSSatish Balay 
85749b5e25fSSatish Balay         if (nonew == 1) goto noinsert1;
858085a36d4SBarry Smith         if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col);
859ed1caa07SMatthew Knepley         MatSeqXAIJReallocateAIJ(a,bs2,nrow,brow,bcol,rmax,aa,ai,aj,a->mbs,rp,ap,imax,nonew);
86049b5e25fSSatish Balay 
861c03d1d03SSatish Balay         N = nrow++ - 1; high++;
86249b5e25fSSatish Balay         /* shift up all the later entries in this row */
86349b5e25fSSatish Balay         for (ii=N; ii>=i; ii--) {
86449b5e25fSSatish Balay           rp[ii+1] = rp[ii];
86549b5e25fSSatish Balay           ierr     = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr);
86649b5e25fSSatish Balay         }
86749b5e25fSSatish Balay         if (N>=i) {
86849b5e25fSSatish Balay           ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr);
86949b5e25fSSatish Balay         }
87049b5e25fSSatish Balay         rp[i]                      = bcol;
87149b5e25fSSatish Balay         ap[bs2*i + bs*cidx + ridx] = value;
87249b5e25fSSatish Balay       noinsert1:;
87349b5e25fSSatish Balay         low = i;
8748549e402SHong Zhang       }
87549b5e25fSSatish Balay     }   /* end of loop over added columns */
87649b5e25fSSatish Balay     ailen[brow] = nrow;
87749b5e25fSSatish Balay   }   /* end of loop over added rows */
87849b5e25fSSatish Balay   PetscFunctionReturn(0);
87949b5e25fSSatish Balay }
88049b5e25fSSatish Balay 
8816849ba73SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ(Mat,IS,MatFactorInfo*,Mat*);
8826849ba73SBarry Smith EXTERN PetscErrorCode MatCholeskyFactor_SeqSBAIJ(Mat,IS,MatFactorInfo*);
88313f74950SBarry Smith EXTERN PetscErrorCode MatIncreaseOverlap_SeqSBAIJ(Mat,PetscInt,IS[],PetscInt);
88413f74950SBarry Smith EXTERN PetscErrorCode MatGetSubMatrix_SeqSBAIJ(Mat,IS,IS,PetscInt,MatReuse,Mat*);
88513f74950SBarry Smith EXTERN PetscErrorCode MatGetSubMatrices_SeqSBAIJ(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat*[]);
886f4df32b1SMatthew Knepley EXTERN PetscErrorCode MatScale_SeqSBAIJ(Mat,PetscScalar);
8876849ba73SBarry Smith EXTERN PetscErrorCode MatNorm_SeqSBAIJ(Mat,NormType,PetscReal *);
8886849ba73SBarry Smith EXTERN PetscErrorCode MatEqual_SeqSBAIJ(Mat,Mat,PetscTruth*);
8896849ba73SBarry Smith EXTERN PetscErrorCode MatGetDiagonal_SeqSBAIJ(Mat,Vec);
8906849ba73SBarry Smith EXTERN PetscErrorCode MatDiagonalScale_SeqSBAIJ(Mat,Vec,Vec);
8916849ba73SBarry Smith EXTERN PetscErrorCode MatGetInfo_SeqSBAIJ(Mat,MatInfoType,MatInfo *);
8926849ba73SBarry Smith EXTERN PetscErrorCode MatZeroEntries_SeqSBAIJ(Mat);
8936849ba73SBarry Smith EXTERN PetscErrorCode MatGetRowMax_SeqSBAIJ(Mat,Vec);
89449b5e25fSSatish Balay 
8956849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_N(Mat,Vec,Vec);
8966849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_1(Mat,Vec,Vec);
8976849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_2(Mat,Vec,Vec);
8986849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_3(Mat,Vec,Vec);
8996849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_4(Mat,Vec,Vec);
9006849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_5(Mat,Vec,Vec);
9016849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_6(Mat,Vec,Vec);
9026849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_7(Mat,Vec,Vec);
90349b5e25fSSatish Balay 
9046849ba73SBarry Smith EXTERN PetscErrorCode MatSolves_SeqSBAIJ_1(Mat,Vecs,Vecs);
905d59c15a7SBarry Smith 
9066849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_1_NaturalOrdering(Mat,Vec,Vec);
9076849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_2_NaturalOrdering(Mat,Vec,Vec);
9086849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_3_NaturalOrdering(Mat,Vec,Vec);
9096849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_4_NaturalOrdering(Mat,Vec,Vec);
9106849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_5_NaturalOrdering(Mat,Vec,Vec);
9116849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_6_NaturalOrdering(Mat,Vec,Vec);
9126849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_7_NaturalOrdering(Mat,Vec,Vec);
9136849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_N_NaturalOrdering(Mat,Vec,Vec);
91407f98182SSatish Balay 
915af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat,MatFactorInfo*,Mat*);
916af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1(Mat,MatFactorInfo*,Mat*);
917af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_2(Mat,MatFactorInfo*,Mat*);
918af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_3(Mat,MatFactorInfo*,Mat*);
919af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_4(Mat,MatFactorInfo*,Mat*);
920af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_5(Mat,MatFactorInfo*,Mat*);
921af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_6(Mat,MatFactorInfo*,Mat*);
922af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_7(Mat,MatFactorInfo*,Mat*);
92313f74950SBarry Smith EXTERN PetscErrorCode MatGetInertia_SeqSBAIJ(Mat,PetscInt*,PetscInt*,PetscInt*);
92449b5e25fSSatish Balay 
9256849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_1(Mat,Vec,Vec);
9266849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_2(Mat,Vec,Vec);
9276849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_3(Mat,Vec,Vec);
9286849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_4(Mat,Vec,Vec);
9296849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_5(Mat,Vec,Vec);
9306849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_6(Mat,Vec,Vec);
9316849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_7(Mat,Vec,Vec);
9326849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_N(Mat,Vec,Vec);
93349b5e25fSSatish Balay 
9346849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_1(Mat,Vec,Vec,Vec);
9356849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_2(Mat,Vec,Vec,Vec);
9366849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_3(Mat,Vec,Vec,Vec);
9376849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_4(Mat,Vec,Vec,Vec);
9386849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_5(Mat,Vec,Vec,Vec);
9396849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_6(Mat,Vec,Vec,Vec);
9406849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_7(Mat,Vec,Vec,Vec);
9416849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_N(Mat,Vec,Vec,Vec);
94249b5e25fSSatish Balay 
9434d101231SSatish Balay #ifdef HAVE_MatICCFactor
9446849ba73SBarry Smith /* modified from MatILUFactor_SeqSBAIJ, needs further work!  */
9454a2ae208SSatish Balay #undef __FUNCT__
9464d101231SSatish Balay #define __FUNCT__ "MatICCFactor_SeqSBAIJ"
947dfbe8321SBarry Smith PetscErrorCode MatICCFactor_SeqSBAIJ(Mat inA,IS row,MatFactorInfo *info)
94849b5e25fSSatish Balay {
9494ccecd49SHong Zhang   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)inA->data;
95049b5e25fSSatish Balay   Mat            outA;
951dfbe8321SBarry Smith   PetscErrorCode ierr;
95249b5e25fSSatish Balay   PetscTruth     row_identity,col_identity;
95349b5e25fSSatish Balay 
95449b5e25fSSatish Balay   PetscFunctionBegin;
95549b5e25fSSatish Balay   outA          = inA;
9561260e1f8SHong Zhang   inA->factor   = FACTOR_CHOLESKY;
95749b5e25fSSatish Balay 
95849b5e25fSSatish Balay   if (!a->diag) {
9591a3463dfSHong Zhang     ierr = MatMarkDiagonal_SeqSBAIJ(inA);CHKERRQ(ierr);
96049b5e25fSSatish Balay   }
96149b5e25fSSatish Balay   /*
96249b5e25fSSatish Balay     Blocksize 2, 3, 4, 5, 6 and 7 have a special faster factorization/solver
96349b5e25fSSatish Balay     for ILU(0) factorization with natural ordering
96449b5e25fSSatish Balay   */
96549b5e25fSSatish Balay   switch (a->bs) {
96649b5e25fSSatish Balay   case 1:
9670fe9e5beSBarry Smith     inA->ops->solve            = MatSolve_SeqSBAIJ_1_NaturalOrdering;
96807f98182SSatish Balay     inA->ops->solvetranspose   = MatSolve_SeqSBAIJ_1_NaturalOrdering;
969d59c15a7SBarry Smith     inA->ops->solves           = MatSolves_SeqSBAIJ_1;
97063ba0a88SBarry Smith     ierr = PetscLoginfo((inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering solvetrans BS=1\n"));CHKERRQ(ierr);
97149b5e25fSSatish Balay   case 2:
9721a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering;
9731a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_2_NaturalOrdering;
9741a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_2_NaturalOrdering;
97563ba0a88SBarry Smith     ierr = PetscLogInfo((inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=2\n"));CHKERRQ(ierr);
97649b5e25fSSatish Balay     break;
97749b5e25fSSatish Balay   case 3:
9781a3463dfSHong Zhang      inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering;
9791a3463dfSHong Zhang      inA->ops->solve           = MatSolve_SeqSBAIJ_3_NaturalOrdering;
9801a3463dfSHong Zhang      inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_3_NaturalOrdering;
98163ba0a88SBarry Smith      ierr = PetscLogInfo((inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=3\n"));CHKERRQ(ierr);
98249b5e25fSSatish Balay      break;
98349b5e25fSSatish Balay   case 4:
9841a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering;
9851a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_4_NaturalOrdering;
9861a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_4_NaturalOrdering;
98763ba0a88SBarry Smith     ierr = PetscLogInfo((inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=4\n"));CHKERRQ(ierr);
98849b5e25fSSatish Balay     break;
98949b5e25fSSatish Balay   case 5:
9901a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering;
9911a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_5_NaturalOrdering;
9921a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_5_NaturalOrdering;
99363ba0a88SBarry Smith     ierr = PetscLogInfo((inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=5\n"));CHKERRQ(ierr);
99449b5e25fSSatish Balay     break;
99549b5e25fSSatish Balay   case 6:
9961a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering;
9971a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_6_NaturalOrdering;
9981a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_6_NaturalOrdering;
99963ba0a88SBarry Smith     ierr = PetscLogInfo((inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=6\n"));CHKERRQ(ierr);
100049b5e25fSSatish Balay     break;
100149b5e25fSSatish Balay   case 7:
10021a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering;
10031a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_7_NaturalOrdering;
10041a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_7_NaturalOrdering;
100563ba0a88SBarry Smith     ierr = PetscLogInfo((inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=7\n"));CHKERRQ(ierr);
100649b5e25fSSatish Balay     break;
100749b5e25fSSatish Balay   default:
100849b5e25fSSatish Balay     a->row        = row;
10091a3463dfSHong Zhang     a->icol       = col;
101049b5e25fSSatish Balay     ierr          = PetscObjectReference((PetscObject)row);CHKERRQ(ierr);
101149b5e25fSSatish Balay     ierr          = PetscObjectReference((PetscObject)col);CHKERRQ(ierr);
101249b5e25fSSatish Balay 
101349b5e25fSSatish Balay     /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */
101449b5e25fSSatish Balay     ierr = ISInvertPermutation(col,PETSC_DECIDE, &(a->icol));CHKERRQ(ierr);
101552e6d16bSBarry Smith     ierr = PetscLogObjectParent(inA,a->icol);CHKERRQ(ierr);
101649b5e25fSSatish Balay 
101749b5e25fSSatish Balay     if (!a->solve_work) {
101887828ca2SBarry Smith       ierr = PetscMalloc((A->m+a->bs)*sizeof(PetscScalar),&a->solve_work);CHKERRQ(ierr);
101952e6d16bSBarry Smith       ierr = PetscLogObjectMemory(inA,(A->m+a->bs)*sizeof(PetscScalar));CHKERRQ(ierr);
102049b5e25fSSatish Balay     }
102149b5e25fSSatish Balay   }
102249b5e25fSSatish Balay 
1023af281ebdSHong Zhang   ierr = MatCholeskyFactorNumeric(inA,info,&outA);CHKERRQ(ierr);
102449b5e25fSSatish Balay   PetscFunctionReturn(0);
102549b5e25fSSatish Balay }
1026950f1e5bSHong Zhang #endif
1027950f1e5bSHong Zhang 
10284a2ae208SSatish Balay #undef __FUNCT__
10294a2ae208SSatish Balay #define __FUNCT__ "MatPrintHelp_SeqSBAIJ"
1030dfbe8321SBarry Smith PetscErrorCode MatPrintHelp_SeqSBAIJ(Mat A)
103149b5e25fSSatish Balay {
103249b5e25fSSatish Balay   static PetscTruth called = PETSC_FALSE;
103349b5e25fSSatish Balay   MPI_Comm          comm = A->comm;
1034dfbe8321SBarry Smith   PetscErrorCode    ierr;
103549b5e25fSSatish Balay 
103649b5e25fSSatish Balay   PetscFunctionBegin;
103749b5e25fSSatish Balay   if (called) {PetscFunctionReturn(0);} else called = PETSC_TRUE;
103849b5e25fSSatish Balay   ierr = (*PetscHelpPrintf)(comm," Options for MATSEQSBAIJ and MATMPISBAIJ matrix formats (the defaults):\n");CHKERRQ(ierr);
103949b5e25fSSatish Balay   ierr = (*PetscHelpPrintf)(comm,"  -mat_block_size <block_size>\n");CHKERRQ(ierr);
1040*f5edf698SHong Zhang   ierr = (*PetscHelpPrintf)(comm,"  -mat_ignore_lower_triangular: Ignore lower triangular values set by user\n");CHKERRQ(ierr);
1041*f5edf698SHong Zhang   ierr = (*PetscHelpPrintf)(comm,"  -mat_getrow_uppertriangular: Enable MatGetRow() for retrieving the upper triangular part of the row\n");CHKERRQ(ierr);
104249b5e25fSSatish Balay   PetscFunctionReturn(0);
104349b5e25fSSatish Balay }
104449b5e25fSSatish Balay 
104549b5e25fSSatish Balay EXTERN_C_BEGIN
10464a2ae208SSatish Balay #undef __FUNCT__
10474a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices_SeqSBAIJ"
1048be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetColumnIndices_SeqSBAIJ(Mat mat,PetscInt *indices)
104949b5e25fSSatish Balay {
1050045c9aa0SHong Zhang   Mat_SeqSBAIJ *baij = (Mat_SeqSBAIJ *)mat->data;
105113f74950SBarry Smith   PetscInt     i,nz,n;
105249b5e25fSSatish Balay 
105349b5e25fSSatish Balay   PetscFunctionBegin;
10546c6c5352SBarry Smith   nz = baij->maxnz;
1055c464158bSHong Zhang   n  = mat->n;
105649b5e25fSSatish Balay   for (i=0; i<nz; i++) {
105749b5e25fSSatish Balay     baij->j[i] = indices[i];
105849b5e25fSSatish Balay   }
10596c6c5352SBarry Smith    baij->nz = nz;
106049b5e25fSSatish Balay    for (i=0; i<n; i++) {
106149b5e25fSSatish Balay      baij->ilen[i] = baij->imax[i];
106249b5e25fSSatish Balay    }
106349b5e25fSSatish Balay    PetscFunctionReturn(0);
106449b5e25fSSatish Balay }
106549b5e25fSSatish Balay EXTERN_C_END
106649b5e25fSSatish Balay 
10674a2ae208SSatish Balay #undef __FUNCT__
10684a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices"
106949b5e25fSSatish Balay /*@
107019585528SSatish Balay   MatSeqSBAIJSetColumnIndices - Set the column indices for all the rows
107149b5e25fSSatish Balay   in the matrix.
107249b5e25fSSatish Balay 
107349b5e25fSSatish Balay   Input Parameters:
107419585528SSatish Balay   +  mat     - the SeqSBAIJ matrix
107549b5e25fSSatish Balay   -  indices - the column indices
107649b5e25fSSatish Balay 
107749b5e25fSSatish Balay   Level: advanced
107849b5e25fSSatish Balay 
107949b5e25fSSatish Balay   Notes:
108049b5e25fSSatish Balay   This can be called if you have precomputed the nonzero structure of the
108149b5e25fSSatish Balay   matrix and want to provide it to the matrix object to improve the performance
108249b5e25fSSatish Balay   of the MatSetValues() operation.
108349b5e25fSSatish Balay 
108449b5e25fSSatish Balay   You MUST have set the correct numbers of nonzeros per row in the call to
1085d1be2dadSMatthew Knepley   MatCreateSeqSBAIJ(), and the columns indices MUST be sorted.
108649b5e25fSSatish Balay 
1087ab9f2c04SSatish Balay   MUST be called before any calls to MatSetValues()
108849b5e25fSSatish Balay 
1089ab9f2c04SSatish Balay   .seealso: MatCreateSeqSBAIJ
109049b5e25fSSatish Balay @*/
1091be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetColumnIndices(Mat mat,PetscInt *indices)
109249b5e25fSSatish Balay {
109313f74950SBarry Smith   PetscErrorCode ierr,(*f)(Mat,PetscInt *);
109449b5e25fSSatish Balay 
109549b5e25fSSatish Balay   PetscFunctionBegin;
10964482741eSBarry Smith   PetscValidHeaderSpecific(mat,MAT_COOKIE,1);
10974482741eSBarry Smith   PetscValidPointer(indices,2);
1098c134de8dSSatish Balay   ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqSBAIJSetColumnIndices_C",(void (**)(void))&f);CHKERRQ(ierr);
109949b5e25fSSatish Balay   if (f) {
110049b5e25fSSatish Balay     ierr = (*f)(mat,indices);CHKERRQ(ierr);
110149b5e25fSSatish Balay   } else {
1102e005ede5SBarry Smith     SETERRQ(PETSC_ERR_SUP,"Wrong type of matrix to set column indices");
110349b5e25fSSatish Balay   }
110449b5e25fSSatish Balay   PetscFunctionReturn(0);
110549b5e25fSSatish Balay }
110649b5e25fSSatish Balay 
11074a2ae208SSatish Balay #undef __FUNCT__
11083c896bc6SHong Zhang #define __FUNCT__ "MatCopy_SeqSBAIJ"
11093c896bc6SHong Zhang PetscErrorCode MatCopy_SeqSBAIJ(Mat A,Mat B,MatStructure str)
11103c896bc6SHong Zhang {
11113c896bc6SHong Zhang   PetscErrorCode ierr;
11123c896bc6SHong Zhang 
11133c896bc6SHong Zhang   PetscFunctionBegin;
11143c896bc6SHong Zhang   /* If the two matrices have the same copy implementation, use fast copy. */
11153c896bc6SHong Zhang   if (str == SAME_NONZERO_PATTERN && (A->ops->copy == B->ops->copy)) {
11163c896bc6SHong Zhang     Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
11173c896bc6SHong Zhang     Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ*)B->data;
11183c896bc6SHong Zhang 
11193c896bc6SHong Zhang     if (a->i[A->m] != b->i[B->m]) {
11203c896bc6SHong Zhang       SETERRQ(PETSC_ERR_ARG_INCOMP,"Number of nonzeros in two matrices are different");
11213c896bc6SHong Zhang     }
11223c896bc6SHong Zhang     ierr = PetscMemcpy(b->a,a->a,(a->i[A->m])*sizeof(PetscScalar));CHKERRQ(ierr);
11233c896bc6SHong Zhang   } else {
1124*f5edf698SHong Zhang     ierr = MatGetRowUpperTriangular(A);CHKERRQ(ierr);
11253c896bc6SHong Zhang     ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr);
1126*f5edf698SHong Zhang     ierr = MatRestoreRowUpperTriangular(A);CHKERRQ(ierr);
11273c896bc6SHong Zhang   }
11283c896bc6SHong Zhang   PetscFunctionReturn(0);
11293c896bc6SHong Zhang }
11303c896bc6SHong Zhang 
11313c896bc6SHong Zhang #undef __FUNCT__
11324a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqSBAIJ"
1133dfbe8321SBarry Smith PetscErrorCode MatSetUpPreallocation_SeqSBAIJ(Mat A)
1134273d9f13SBarry Smith {
1135dfbe8321SBarry Smith   PetscErrorCode ierr;
1136273d9f13SBarry Smith 
1137273d9f13SBarry Smith   PetscFunctionBegin;
1138ab93d7beSBarry Smith   ierr =  MatSeqSBAIJSetPreallocation_SeqSBAIJ(A,1,PETSC_DEFAULT,0);CHKERRQ(ierr);
1139273d9f13SBarry Smith   PetscFunctionReturn(0);
1140273d9f13SBarry Smith }
1141273d9f13SBarry Smith 
1142a6ece127SHong Zhang #undef __FUNCT__
1143a6ece127SHong Zhang #define __FUNCT__ "MatGetArray_SeqSBAIJ"
1144dfbe8321SBarry Smith PetscErrorCode MatGetArray_SeqSBAIJ(Mat A,PetscScalar *array[])
1145a6ece127SHong Zhang {
1146a6ece127SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
1147a6ece127SHong Zhang   PetscFunctionBegin;
1148a6ece127SHong Zhang   *array = a->a;
1149a6ece127SHong Zhang   PetscFunctionReturn(0);
1150a6ece127SHong Zhang }
1151a6ece127SHong Zhang 
1152a6ece127SHong Zhang #undef __FUNCT__
1153a6ece127SHong Zhang #define __FUNCT__ "MatRestoreArray_SeqSBAIJ"
1154dfbe8321SBarry Smith PetscErrorCode MatRestoreArray_SeqSBAIJ(Mat A,PetscScalar *array[])
1155a6ece127SHong Zhang {
1156a6ece127SHong Zhang   PetscFunctionBegin;
1157a6ece127SHong Zhang   PetscFunctionReturn(0);
1158a6ece127SHong Zhang  }
1159a6ece127SHong Zhang 
116042ee4b1aSHong Zhang #include "petscblaslapack.h"
116142ee4b1aSHong Zhang #undef __FUNCT__
116242ee4b1aSHong Zhang #define __FUNCT__ "MatAXPY_SeqSBAIJ"
1163f4df32b1SMatthew Knepley PetscErrorCode MatAXPY_SeqSBAIJ(Mat Y,PetscScalar a,Mat X,MatStructure str)
116442ee4b1aSHong Zhang {
116542ee4b1aSHong Zhang   Mat_SeqSBAIJ   *x=(Mat_SeqSBAIJ *)X->data, *y=(Mat_SeqSBAIJ *)Y->data;
1166dfbe8321SBarry Smith   PetscErrorCode ierr;
1167521d7252SBarry Smith   PetscInt       i,bs=Y->bs,bs2,j;
11684ce68768SBarry Smith   PetscBLASInt   bnz = (PetscBLASInt)x->nz,one = 1;
116942ee4b1aSHong Zhang 
117042ee4b1aSHong Zhang   PetscFunctionBegin;
117142ee4b1aSHong Zhang   if (str == SAME_NONZERO_PATTERN) {
1172f4df32b1SMatthew Knepley     PetscScalar alpha = a;
1173f4df32b1SMatthew Knepley     BLASaxpy_(&bnz,&alpha,x->a,&one,y->a,&one);
1174c537a176SHong Zhang   } else if (str == SUBSET_NONZERO_PATTERN) { /* nonzeros of X is a subset of Y's */
1175c4319e64SHong Zhang     if (y->xtoy && y->XtoY != X) {
1176c4319e64SHong Zhang       ierr = PetscFree(y->xtoy);CHKERRQ(ierr);
1177c4319e64SHong Zhang       ierr = MatDestroy(y->XtoY);CHKERRQ(ierr);
1178c537a176SHong Zhang     }
1179c4319e64SHong Zhang     if (!y->xtoy) { /* get xtoy */
1180c4319e64SHong Zhang       ierr = MatAXPYGetxtoy_Private(x->mbs,x->i,x->j,PETSC_NULL, y->i,y->j,PETSC_NULL, &y->xtoy);CHKERRQ(ierr);
1181c4319e64SHong Zhang       y->XtoY = X;
1182c537a176SHong Zhang     }
1183c4319e64SHong Zhang     bs2 = bs*bs;
11846c6c5352SBarry Smith     for (i=0; i<x->nz; i++) {
1185c4319e64SHong Zhang       j = 0;
1186c4319e64SHong Zhang       while (j < bs2){
1187f4df32b1SMatthew Knepley         y->a[bs2*y->xtoy[i]+j] += a*(x->a[bs2*i+j]);
1188c4319e64SHong Zhang         j++;
1189c537a176SHong Zhang       }
1190c4319e64SHong Zhang     }
119163ba0a88SBarry Smith     ierr = PetscLogInfo((0,"MatAXPY_SeqSBAIJ: 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);
119242ee4b1aSHong Zhang   } else {
1193*f5edf698SHong Zhang     ierr = MatGetRowUpperTriangular(X);CHKERRQ(ierr);
1194f4df32b1SMatthew Knepley     ierr = MatAXPY_Basic(Y,a,X,str);CHKERRQ(ierr);
1195*f5edf698SHong Zhang     ierr = MatRestoreRowUpperTriangular(X);CHKERRQ(ierr);
119642ee4b1aSHong Zhang   }
119742ee4b1aSHong Zhang   PetscFunctionReturn(0);
119842ee4b1aSHong Zhang }
119942ee4b1aSHong Zhang 
1200efcf0fc3SBarry Smith #undef __FUNCT__
1201efcf0fc3SBarry Smith #define __FUNCT__ "MatIsSymmetric_SeqSBAIJ"
1202dfbe8321SBarry Smith PetscErrorCode MatIsSymmetric_SeqSBAIJ(Mat A,PetscReal tol,PetscTruth *flg)
1203efcf0fc3SBarry Smith {
1204efcf0fc3SBarry Smith   PetscFunctionBegin;
1205efcf0fc3SBarry Smith   *flg = PETSC_TRUE;
1206efcf0fc3SBarry Smith   PetscFunctionReturn(0);
1207efcf0fc3SBarry Smith }
1208efcf0fc3SBarry Smith 
1209efcf0fc3SBarry Smith #undef __FUNCT__
1210efcf0fc3SBarry Smith #define __FUNCT__ "MatIsStructurallySymmetric_SeqSBAIJ"
1211dfbe8321SBarry Smith PetscErrorCode MatIsStructurallySymmetric_SeqSBAIJ(Mat A,PetscTruth *flg)
1212efcf0fc3SBarry Smith {
1213efcf0fc3SBarry Smith    PetscFunctionBegin;
1214efcf0fc3SBarry Smith    *flg = PETSC_TRUE;
1215efcf0fc3SBarry Smith    PetscFunctionReturn(0);
1216efcf0fc3SBarry Smith }
1217efcf0fc3SBarry Smith 
1218efcf0fc3SBarry Smith #undef __FUNCT__
1219efcf0fc3SBarry Smith #define __FUNCT__ "MatIsHermitian_SeqSBAIJ"
1220dfbe8321SBarry Smith PetscErrorCode MatIsHermitian_SeqSBAIJ(Mat A,PetscTruth *flg)
1221efcf0fc3SBarry Smith  {
1222efcf0fc3SBarry Smith    PetscFunctionBegin;
1223efcf0fc3SBarry Smith    *flg = PETSC_FALSE;
1224efcf0fc3SBarry Smith    PetscFunctionReturn(0);
1225efcf0fc3SBarry Smith  }
1226efcf0fc3SBarry Smith 
122799cafbc1SBarry Smith #undef __FUNCT__
122899cafbc1SBarry Smith #define __FUNCT__ "MatRealPart_SeqSBAIJ"
122999cafbc1SBarry Smith PetscErrorCode MatRealPart_SeqSBAIJ(Mat A)
123099cafbc1SBarry Smith {
123199cafbc1SBarry Smith   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
123299cafbc1SBarry Smith   PetscInt       i,nz = a->bs2*a->i[a->mbs];
123399cafbc1SBarry Smith   PetscScalar    *aa = a->a;
123499cafbc1SBarry Smith 
123599cafbc1SBarry Smith   PetscFunctionBegin;
123699cafbc1SBarry Smith   for (i=0; i<nz; i++) aa[i] = PetscRealPart(aa[i]);
123799cafbc1SBarry Smith   PetscFunctionReturn(0);
123899cafbc1SBarry Smith }
123999cafbc1SBarry Smith 
124099cafbc1SBarry Smith #undef __FUNCT__
124199cafbc1SBarry Smith #define __FUNCT__ "MatImaginaryPart_SeqSBAIJ"
124299cafbc1SBarry Smith PetscErrorCode MatImaginaryPart_SeqSBAIJ(Mat A)
124399cafbc1SBarry Smith {
124499cafbc1SBarry Smith   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
124599cafbc1SBarry Smith   PetscInt       i,nz = a->bs2*a->i[a->mbs];
124699cafbc1SBarry Smith   PetscScalar    *aa = a->a;
124799cafbc1SBarry Smith 
124899cafbc1SBarry Smith   PetscFunctionBegin;
124999cafbc1SBarry Smith   for (i=0; i<nz; i++) aa[i] = PetscImaginaryPart(aa[i]);
125099cafbc1SBarry Smith   PetscFunctionReturn(0);
125199cafbc1SBarry Smith }
125299cafbc1SBarry Smith 
125349b5e25fSSatish Balay /* -------------------------------------------------------------------*/
125449b5e25fSSatish Balay static struct _MatOps MatOps_Values = {MatSetValues_SeqSBAIJ,
125549b5e25fSSatish Balay        MatGetRow_SeqSBAIJ,
125649b5e25fSSatish Balay        MatRestoreRow_SeqSBAIJ,
125749b5e25fSSatish Balay        MatMult_SeqSBAIJ_N,
125897304618SKris Buschelman /* 4*/ MatMultAdd_SeqSBAIJ_N,
1259431c96f7SBarry Smith        MatMult_SeqSBAIJ_N,       /* transpose versions are same as non-transpose versions */
1260e005ede5SBarry Smith        MatMultAdd_SeqSBAIJ_N,
126149b5e25fSSatish Balay        MatSolve_SeqSBAIJ_N,
126249b5e25fSSatish Balay        0,
126349b5e25fSSatish Balay        0,
126497304618SKris Buschelman /*10*/ 0,
126549b5e25fSSatish Balay        0,
12665f4ef2deSHong Zhang        MatCholeskyFactor_SeqSBAIJ,
1267d06b337dSHong Zhang        MatRelax_SeqSBAIJ,
126849b5e25fSSatish Balay        MatTranspose_SeqSBAIJ,
126997304618SKris Buschelman /*15*/ MatGetInfo_SeqSBAIJ,
127049b5e25fSSatish Balay        MatEqual_SeqSBAIJ,
127149b5e25fSSatish Balay        MatGetDiagonal_SeqSBAIJ,
127249b5e25fSSatish Balay        MatDiagonalScale_SeqSBAIJ,
127349b5e25fSSatish Balay        MatNorm_SeqSBAIJ,
127497304618SKris Buschelman /*20*/ 0,
127549b5e25fSSatish Balay        MatAssemblyEnd_SeqSBAIJ,
127649b5e25fSSatish Balay        0,
127749b5e25fSSatish Balay        MatSetOption_SeqSBAIJ,
127849b5e25fSSatish Balay        MatZeroEntries_SeqSBAIJ,
1279dcf5cc72SBarry Smith /*25*/ 0,
128049b5e25fSSatish Balay        0,
128149b5e25fSSatish Balay        0,
12825f4ef2deSHong Zhang        MatCholeskyFactorSymbolic_SeqSBAIJ,
12835f4ef2deSHong Zhang        MatCholeskyFactorNumeric_SeqSBAIJ_N,
128497304618SKris Buschelman /*30*/ MatSetUpPreallocation_SeqSBAIJ,
1285c464158bSHong Zhang        0,
12864d101231SSatish Balay        MatICCFactorSymbolic_SeqSBAIJ,
1287a6ece127SHong Zhang        MatGetArray_SeqSBAIJ,
1288a6ece127SHong Zhang        MatRestoreArray_SeqSBAIJ,
128997304618SKris Buschelman /*35*/ MatDuplicate_SeqSBAIJ,
129049b5e25fSSatish Balay        0,
129149b5e25fSSatish Balay        0,
129249b5e25fSSatish Balay        0,
1293950f1e5bSHong Zhang        0,
129497304618SKris Buschelman /*40*/ MatAXPY_SeqSBAIJ,
129549b5e25fSSatish Balay        MatGetSubMatrices_SeqSBAIJ,
129649b5e25fSSatish Balay        MatIncreaseOverlap_SeqSBAIJ,
129749b5e25fSSatish Balay        MatGetValues_SeqSBAIJ,
12983c896bc6SHong Zhang        MatCopy_SeqSBAIJ,
129997304618SKris Buschelman /*45*/ MatPrintHelp_SeqSBAIJ,
130049b5e25fSSatish Balay        MatScale_SeqSBAIJ,
130149b5e25fSSatish Balay        0,
130249b5e25fSSatish Balay        0,
130349b5e25fSSatish Balay        0,
1304521d7252SBarry Smith /*50*/ 0,
130549b5e25fSSatish Balay        MatGetRowIJ_SeqSBAIJ,
130649b5e25fSSatish Balay        MatRestoreRowIJ_SeqSBAIJ,
130749b5e25fSSatish Balay        0,
130849b5e25fSSatish Balay        0,
130997304618SKris Buschelman /*55*/ 0,
131049b5e25fSSatish Balay        0,
131149b5e25fSSatish Balay        0,
131249b5e25fSSatish Balay        0,
131349b5e25fSSatish Balay        MatSetValuesBlocked_SeqSBAIJ,
131497304618SKris Buschelman /*60*/ MatGetSubMatrix_SeqSBAIJ,
131549b5e25fSSatish Balay        0,
131649b5e25fSSatish Balay        0,
13178a124369SBarry Smith        MatGetPetscMaps_Petsc,
1318d959ec07SHong Zhang        0,
131997304618SKris Buschelman /*65*/ 0,
1320d959ec07SHong Zhang        0,
1321d959ec07SHong Zhang        0,
1322d959ec07SHong Zhang        0,
1323d959ec07SHong Zhang        0,
132497304618SKris Buschelman /*70*/ MatGetRowMax_SeqSBAIJ,
13253e0d88b5SBarry Smith        0,
13263e0d88b5SBarry Smith        0,
13273e0d88b5SBarry Smith        0,
13283e0d88b5SBarry Smith        0,
132997304618SKris Buschelman /*75*/ 0,
13303e0d88b5SBarry Smith        0,
13313e0d88b5SBarry Smith        0,
13323e0d88b5SBarry Smith        0,
13333e0d88b5SBarry Smith        0,
133497304618SKris Buschelman /*80*/ 0,
13353e0d88b5SBarry Smith        0,
13363e0d88b5SBarry Smith        0,
13373e0d88b5SBarry Smith #if !defined(PETSC_USE_COMPLEX)
133897304618SKris Buschelman        MatGetInertia_SeqSBAIJ,
13393e0d88b5SBarry Smith #else
134097304618SKris Buschelman        0,
13413e0d88b5SBarry Smith #endif
1342865e5f61SKris Buschelman        MatLoad_SeqSBAIJ,
1343865e5f61SKris Buschelman /*85*/ MatIsSymmetric_SeqSBAIJ,
1344865e5f61SKris Buschelman        MatIsHermitian_SeqSBAIJ,
1345efcf0fc3SBarry Smith        MatIsStructurallySymmetric_SeqSBAIJ,
1346865e5f61SKris Buschelman        0,
1347865e5f61SKris Buschelman        0,
1348865e5f61SKris Buschelman /*90*/ 0,
1349865e5f61SKris Buschelman        0,
1350865e5f61SKris Buschelman        0,
1351865e5f61SKris Buschelman        0,
1352865e5f61SKris Buschelman        0,
1353865e5f61SKris Buschelman /*95*/ 0,
1354865e5f61SKris Buschelman        0,
1355865e5f61SKris Buschelman        0,
135699cafbc1SBarry Smith        0,
135799cafbc1SBarry Smith        0,
135899cafbc1SBarry Smith /*100*/0,
135999cafbc1SBarry Smith        0,
136099cafbc1SBarry Smith        0,
136199cafbc1SBarry Smith        0,
136299cafbc1SBarry Smith        0,
136399cafbc1SBarry Smith /*105*/0,
136499cafbc1SBarry Smith        MatRealPart_SeqSBAIJ,
1365*f5edf698SHong Zhang        MatImaginaryPart_SeqSBAIJ,
1366*f5edf698SHong Zhang        MatGetRowUpperTriangular_SeqSBAIJ,
1367*f5edf698SHong Zhang        MatRestoreRowUpperTriangular_SeqSBAIJ
136899cafbc1SBarry Smith };
1369be1d678aSKris Buschelman 
137049b5e25fSSatish Balay EXTERN_C_BEGIN
13714a2ae208SSatish Balay #undef __FUNCT__
13724a2ae208SSatish Balay #define __FUNCT__ "MatStoreValues_SeqSBAIJ"
1373be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatStoreValues_SeqSBAIJ(Mat mat)
137449b5e25fSSatish Balay {
13754afc71dfSHong Zhang   Mat_SeqSBAIJ   *aij = (Mat_SeqSBAIJ *)mat->data;
1376521d7252SBarry Smith   PetscInt       nz = aij->i[mat->m]*mat->bs*aij->bs2;
1377dfbe8321SBarry Smith   PetscErrorCode ierr;
137849b5e25fSSatish Balay 
137949b5e25fSSatish Balay   PetscFunctionBegin;
138049b5e25fSSatish Balay   if (aij->nonew != 1) {
1381e005ede5SBarry Smith     SETERRQ(PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
138249b5e25fSSatish Balay   }
138349b5e25fSSatish Balay 
138449b5e25fSSatish Balay   /* allocate space for values if not already there */
138549b5e25fSSatish Balay   if (!aij->saved_values) {
138687828ca2SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&aij->saved_values);CHKERRQ(ierr);
138749b5e25fSSatish Balay   }
138849b5e25fSSatish Balay 
138949b5e25fSSatish Balay   /* copy values over */
139087828ca2SBarry Smith   ierr = PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(PetscScalar));CHKERRQ(ierr);
139149b5e25fSSatish Balay   PetscFunctionReturn(0);
139249b5e25fSSatish Balay }
139349b5e25fSSatish Balay EXTERN_C_END
139449b5e25fSSatish Balay 
139549b5e25fSSatish Balay EXTERN_C_BEGIN
13964a2ae208SSatish Balay #undef __FUNCT__
13974a2ae208SSatish Balay #define __FUNCT__ "MatRetrieveValues_SeqSBAIJ"
1398be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatRetrieveValues_SeqSBAIJ(Mat mat)
139949b5e25fSSatish Balay {
14004afc71dfSHong Zhang   Mat_SeqSBAIJ   *aij = (Mat_SeqSBAIJ *)mat->data;
14016849ba73SBarry Smith   PetscErrorCode ierr;
1402521d7252SBarry Smith   PetscInt       nz = aij->i[mat->m]*mat->bs*aij->bs2;
140349b5e25fSSatish Balay 
140449b5e25fSSatish Balay   PetscFunctionBegin;
140549b5e25fSSatish Balay   if (aij->nonew != 1) {
1406e005ede5SBarry Smith     SETERRQ(PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
140749b5e25fSSatish Balay   }
140849b5e25fSSatish Balay   if (!aij->saved_values) {
1409e005ede5SBarry Smith     SETERRQ(PETSC_ERR_ORDER,"Must call MatStoreValues(A);first");
141049b5e25fSSatish Balay   }
141149b5e25fSSatish Balay 
141249b5e25fSSatish Balay   /* copy values over */
141387828ca2SBarry Smith   ierr = PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(PetscScalar));CHKERRQ(ierr);
141449b5e25fSSatish Balay   PetscFunctionReturn(0);
141549b5e25fSSatish Balay }
141649b5e25fSSatish Balay EXTERN_C_END
141749b5e25fSSatish Balay 
14188549e402SHong Zhang EXTERN_C_BEGIN
14194a2ae208SSatish Balay #undef __FUNCT__
1420a23d5eceSKris Buschelman #define __FUNCT__ "MatSeqSBAIJSetPreallocation_SeqSBAIJ"
1421be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetPreallocation_SeqSBAIJ(Mat B,PetscInt bs,PetscInt nz,PetscInt *nnz)
142249b5e25fSSatish Balay {
1423c464158bSHong Zhang   Mat_SeqSBAIJ   *b = (Mat_SeqSBAIJ*)B->data;
14246849ba73SBarry Smith   PetscErrorCode ierr;
1425ab93d7beSBarry Smith   PetscInt       i,mbs,bs2;
1426ab93d7beSBarry Smith   PetscTruth     skipallocation = PETSC_FALSE,flg;
142749b5e25fSSatish Balay 
142849b5e25fSSatish Balay   PetscFunctionBegin;
1429273d9f13SBarry Smith   B->preallocated = PETSC_TRUE;
1430e82a3eeeSBarry Smith   ierr = PetscOptionsGetInt(B->prefix,"-mat_block_size",&bs,PETSC_NULL);CHKERRQ(ierr);
1431c464158bSHong Zhang   mbs  = B->m/bs;
143249b5e25fSSatish Balay   bs2  = bs*bs;
143349b5e25fSSatish Balay 
1434c464158bSHong Zhang   if (mbs*bs != B->m) {
143529bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_SIZ,"Number rows, cols must be divisible by blocksize");
143649b5e25fSSatish Balay   }
143749b5e25fSSatish Balay 
1438ab93d7beSBarry Smith   if (nz == MAT_SKIP_ALLOCATION) {
1439ab93d7beSBarry Smith     skipallocation = PETSC_TRUE;
1440ab93d7beSBarry Smith     nz             = 0;
1441ab93d7beSBarry Smith   }
1442ab93d7beSBarry Smith 
1443435da068SBarry Smith   if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 3;
144477431f27SBarry Smith   if (nz < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %D",nz);
144549b5e25fSSatish Balay   if (nnz) {
144649b5e25fSSatish Balay     for (i=0; i<mbs; i++) {
144777431f27SBarry Smith       if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %D value %D",i,nnz[i]);
144877431f27SBarry 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);
144949b5e25fSSatish Balay     }
145049b5e25fSSatish Balay   }
145149b5e25fSSatish Balay 
1452e82a3eeeSBarry Smith   ierr    = PetscOptionsHasName(B->prefix,"-mat_no_unroll",&flg);CHKERRQ(ierr);
145349b5e25fSSatish Balay   if (!flg) {
145449b5e25fSSatish Balay     switch (bs) {
145549b5e25fSSatish Balay     case 1:
14564ccecd49SHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1;
145749b5e25fSSatish Balay       B->ops->solve            = MatSolve_SeqSBAIJ_1;
1458d59c15a7SBarry Smith       B->ops->solves           = MatSolves_SeqSBAIJ_1;
1459e005ede5SBarry Smith       B->ops->solvetranspose   = MatSolve_SeqSBAIJ_1;
146049b5e25fSSatish Balay       B->ops->mult             = MatMult_SeqSBAIJ_1;
146149b5e25fSSatish Balay       B->ops->multadd          = MatMultAdd_SeqSBAIJ_1;
1462431c96f7SBarry Smith       B->ops->multtranspose    = MatMult_SeqSBAIJ_1;
1463431c96f7SBarry Smith       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_1;
146449b5e25fSSatish Balay       break;
146549b5e25fSSatish Balay     case 2:
14664ccecd49SHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2;
146749b5e25fSSatish Balay       B->ops->solve            = MatSolve_SeqSBAIJ_2;
1468e005ede5SBarry Smith       B->ops->solvetranspose   = MatSolve_SeqSBAIJ_2;
146949b5e25fSSatish Balay       B->ops->mult             = MatMult_SeqSBAIJ_2;
147049b5e25fSSatish Balay       B->ops->multadd          = MatMultAdd_SeqSBAIJ_2;
1471431c96f7SBarry Smith       B->ops->multtranspose    = MatMult_SeqSBAIJ_2;
1472431c96f7SBarry Smith       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_2;
147349b5e25fSSatish Balay       break;
147449b5e25fSSatish Balay     case 3:
14755f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3;
147649b5e25fSSatish Balay       B->ops->solve            = MatSolve_SeqSBAIJ_3;
1477e005ede5SBarry Smith       B->ops->solvetranspose   = MatSolve_SeqSBAIJ_3;
147849b5e25fSSatish Balay       B->ops->mult             = MatMult_SeqSBAIJ_3;
147949b5e25fSSatish Balay       B->ops->multadd          = MatMultAdd_SeqSBAIJ_3;
1480431c96f7SBarry Smith       B->ops->multtranspose    = MatMult_SeqSBAIJ_3;
1481431c96f7SBarry Smith       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_3;
148249b5e25fSSatish Balay       break;
148349b5e25fSSatish Balay     case 4:
14845f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4;
148549b5e25fSSatish Balay       B->ops->solve            = MatSolve_SeqSBAIJ_4;
1486e005ede5SBarry Smith       B->ops->solvetranspose   = MatSolve_SeqSBAIJ_4;
148749b5e25fSSatish Balay       B->ops->mult             = MatMult_SeqSBAIJ_4;
148849b5e25fSSatish Balay       B->ops->multadd          = MatMultAdd_SeqSBAIJ_4;
1489431c96f7SBarry Smith       B->ops->multtranspose    = MatMult_SeqSBAIJ_4;
1490431c96f7SBarry Smith       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_4;
149149b5e25fSSatish Balay       break;
149249b5e25fSSatish Balay     case 5:
14935f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5;
149449b5e25fSSatish Balay       B->ops->solve            = MatSolve_SeqSBAIJ_5;
1495e005ede5SBarry Smith       B->ops->solvetranspose   = MatSolve_SeqSBAIJ_5;
149649b5e25fSSatish Balay       B->ops->mult             = MatMult_SeqSBAIJ_5;
149749b5e25fSSatish Balay       B->ops->multadd          = MatMultAdd_SeqSBAIJ_5;
1498431c96f7SBarry Smith       B->ops->multtranspose    = MatMult_SeqSBAIJ_5;
1499431c96f7SBarry Smith       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_5;
150049b5e25fSSatish Balay       break;
150149b5e25fSSatish Balay     case 6:
15025f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6;
150349b5e25fSSatish Balay       B->ops->solve            = MatSolve_SeqSBAIJ_6;
1504e005ede5SBarry Smith       B->ops->solvetranspose   = MatSolve_SeqSBAIJ_6;
150549b5e25fSSatish Balay       B->ops->mult             = MatMult_SeqSBAIJ_6;
150649b5e25fSSatish Balay       B->ops->multadd          = MatMultAdd_SeqSBAIJ_6;
1507431c96f7SBarry Smith       B->ops->multtranspose    = MatMult_SeqSBAIJ_6;
1508431c96f7SBarry Smith       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_6;
150949b5e25fSSatish Balay       break;
151049b5e25fSSatish Balay     case 7:
1511de53e5efSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7;
151249b5e25fSSatish Balay       B->ops->solve            = MatSolve_SeqSBAIJ_7;
1513e005ede5SBarry Smith       B->ops->solvetranspose   = MatSolve_SeqSBAIJ_7;
1514de53e5efSHong Zhang       B->ops->mult             = MatMult_SeqSBAIJ_7;
151549b5e25fSSatish Balay       B->ops->multadd          = MatMultAdd_SeqSBAIJ_7;
1516431c96f7SBarry Smith       B->ops->multtranspose    = MatMult_SeqSBAIJ_7;
1517431c96f7SBarry Smith       B->ops->multtransposeadd = MatMultAdd_SeqSBAIJ_7;
151849b5e25fSSatish Balay       break;
151949b5e25fSSatish Balay     }
152049b5e25fSSatish Balay   }
152149b5e25fSSatish Balay 
152249b5e25fSSatish Balay   b->mbs = mbs;
15234afc71dfSHong Zhang   b->nbs = mbs;
1524ab93d7beSBarry Smith   if (!skipallocation) {
1525ab93d7beSBarry Smith     /* b->ilen will count nonzeros in each block row so far. */
1526ab93d7beSBarry Smith     ierr   = PetscMalloc2(mbs,PetscInt,&b->imax,mbs,PetscInt,&b->ilen);CHKERRQ(ierr);
1527ab93d7beSBarry Smith     for (i=0; i<mbs; i++) { b->ilen[i] = 0;}
152849b5e25fSSatish Balay     if (!nnz) {
1529435da068SBarry Smith       if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5;
153049b5e25fSSatish Balay       else if (nz <= 0)        nz = 1;
153149b5e25fSSatish Balay       for (i=0; i<mbs; i++) {
15328cef66ccSHong Zhang         b->imax[i] = nz;
153349b5e25fSSatish Balay       }
1534153ea458SHong Zhang       nz = nz*mbs; /* total nz */
153549b5e25fSSatish Balay     } else {
153649b5e25fSSatish Balay       nz = 0;
15378cef66ccSHong Zhang       for (i=0; i<mbs; i++) {b->imax[i] = nnz[i]; nz += nnz[i];}
153849b5e25fSSatish Balay     }
15396c6c5352SBarry Smith     /* nz=(nz+mbs)/2; */ /* total diagonal and superdiagonal nonzero blocks */
154049b5e25fSSatish Balay 
154149b5e25fSSatish Balay     /* allocate the matrix space */
1542a96a251dSBarry Smith     ierr = PetscMalloc3(bs2*nz,PetscScalar,&b->a,nz,PetscInt,&b->j,B->m+1,PetscInt,&b->i);CHKERRQ(ierr);
15436c6c5352SBarry Smith     ierr = PetscMemzero(b->a,nz*bs2*sizeof(MatScalar));CHKERRQ(ierr);
154413f74950SBarry Smith     ierr = PetscMemzero(b->j,nz*sizeof(PetscInt));CHKERRQ(ierr);
154549b5e25fSSatish Balay     b->singlemalloc = PETSC_TRUE;
154649b5e25fSSatish Balay 
154749b5e25fSSatish Balay     /* pointer to beginning of each row */
154849b5e25fSSatish Balay     b->i[0] = 0;
154949b5e25fSSatish Balay     for (i=1; i<mbs+1; i++) {
155049b5e25fSSatish Balay       b->i[i] = b->i[i-1] + b->imax[i-1];
155149b5e25fSSatish Balay     }
1552ab93d7beSBarry Smith   }
155349b5e25fSSatish Balay 
1554521d7252SBarry Smith   B->bs               = bs;
155549b5e25fSSatish Balay   b->bs2              = bs2;
15566c6c5352SBarry Smith   b->nz             = 0;
15576c6c5352SBarry Smith   b->maxnz          = nz*bs2;
1558153ea458SHong Zhang 
155916cdd363SHong Zhang   b->inew             = 0;
156016cdd363SHong Zhang   b->jnew             = 0;
156116cdd363SHong Zhang   b->anew             = 0;
156216cdd363SHong Zhang   b->a2anew           = 0;
15631a3463dfSHong Zhang   b->permute          = PETSC_FALSE;
1564c464158bSHong Zhang   PetscFunctionReturn(0);
1565c464158bSHong Zhang }
1566a23d5eceSKris Buschelman EXTERN_C_END
1567153ea458SHong Zhang 
1568d769727bSBarry Smith EXTERN_C_BEGIN
1569f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqSBAIJ_SeqAIJ(Mat, MatType,MatReuse,Mat*);
1570f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqSBAIJ_SeqBAIJ(Mat, MatType,MatReuse,Mat*);
1571d769727bSBarry Smith EXTERN_C_END
1572d769727bSBarry Smith 
15730bad9183SKris Buschelman /*MC
1574fafad747SKris Buschelman   MATSEQSBAIJ - MATSEQSBAIJ = "seqsbaij" - A matrix type to be used for sequential symmetric block sparse matrices,
15750bad9183SKris Buschelman   based on block compressed sparse row format.  Only the upper triangular portion of the matrix is stored.
15760bad9183SKris Buschelman 
15770bad9183SKris Buschelman   Options Database Keys:
15780bad9183SKris Buschelman   . -mat_type seqsbaij - sets the matrix type to "seqsbaij" during a call to MatSetFromOptions()
15790bad9183SKris Buschelman 
15800bad9183SKris Buschelman   Level: beginner
15810bad9183SKris Buschelman 
15820bad9183SKris Buschelman   .seealso: MatCreateSeqSBAIJ
15830bad9183SKris Buschelman M*/
15840bad9183SKris Buschelman 
1585a23d5eceSKris Buschelman EXTERN_C_BEGIN
1586a23d5eceSKris Buschelman #undef __FUNCT__
1587a23d5eceSKris Buschelman #define __FUNCT__ "MatCreate_SeqSBAIJ"
1588be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_SeqSBAIJ(Mat B)
1589a23d5eceSKris Buschelman {
1590a23d5eceSKris Buschelman   Mat_SeqSBAIJ   *b;
1591dfbe8321SBarry Smith   PetscErrorCode ierr;
159213f74950SBarry Smith   PetscMPIInt    size;
1593941593c8SHong Zhang   PetscTruth     flg;
1594a23d5eceSKris Buschelman 
1595a23d5eceSKris Buschelman   PetscFunctionBegin;
1596a23d5eceSKris Buschelman   ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr);
1597a23d5eceSKris Buschelman   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1");
1598a23d5eceSKris Buschelman   B->m = B->M = PetscMax(B->m,B->M);
1599a23d5eceSKris Buschelman   B->n = B->N = PetscMax(B->n,B->N);
1600a23d5eceSKris Buschelman 
1601a23d5eceSKris Buschelman   ierr    = PetscNew(Mat_SeqSBAIJ,&b);CHKERRQ(ierr);
1602a23d5eceSKris Buschelman   B->data = (void*)b;
1603a23d5eceSKris Buschelman   ierr    = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
1604a23d5eceSKris Buschelman   B->ops->destroy     = MatDestroy_SeqSBAIJ;
1605a23d5eceSKris Buschelman   B->ops->view        = MatView_SeqSBAIJ;
1606a23d5eceSKris Buschelman   B->factor           = 0;
1607a23d5eceSKris Buschelman   B->mapping          = 0;
1608a23d5eceSKris Buschelman   b->row              = 0;
1609a23d5eceSKris Buschelman   b->icol             = 0;
1610a23d5eceSKris Buschelman   b->reallocs         = 0;
1611a23d5eceSKris Buschelman   b->saved_values     = 0;
1612a23d5eceSKris Buschelman 
1613a23d5eceSKris Buschelman   ierr = PetscMapCreateMPI(B->comm,B->m,B->M,&B->rmap);CHKERRQ(ierr);
1614a23d5eceSKris Buschelman   ierr = PetscMapCreateMPI(B->comm,B->n,B->N,&B->cmap);CHKERRQ(ierr);
1615a23d5eceSKris Buschelman 
1616a23d5eceSKris Buschelman   b->sorted           = PETSC_FALSE;
1617a23d5eceSKris Buschelman   b->roworiented      = PETSC_TRUE;
1618a23d5eceSKris Buschelman   b->nonew            = 0;
1619a23d5eceSKris Buschelman   b->diag             = 0;
1620a23d5eceSKris Buschelman   b->solve_work       = 0;
1621a23d5eceSKris Buschelman   b->mult_work        = 0;
1622a23d5eceSKris Buschelman   B->spptr            = 0;
1623a23d5eceSKris Buschelman   b->keepzeroedrows   = PETSC_FALSE;
1624a23d5eceSKris Buschelman   b->xtoy             = 0;
1625a23d5eceSKris Buschelman   b->XtoY             = 0;
1626a23d5eceSKris Buschelman 
1627a23d5eceSKris Buschelman   b->inew             = 0;
1628a23d5eceSKris Buschelman   b->jnew             = 0;
1629a23d5eceSKris Buschelman   b->anew             = 0;
1630a23d5eceSKris Buschelman   b->a2anew           = 0;
1631a23d5eceSKris Buschelman   b->permute          = PETSC_FALSE;
1632a23d5eceSKris Buschelman 
1633941593c8SHong Zhang   b->ignore_ltriangular = PETSC_FALSE;
1634941593c8SHong Zhang   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_ignore_lower_triangular",&flg);CHKERRQ(ierr);
1635941593c8SHong Zhang   if (flg) b->ignore_ltriangular = PETSC_TRUE;
1636941593c8SHong Zhang 
1637*f5edf698SHong Zhang   b->getrow_utriangular = PETSC_FALSE;
1638*f5edf698SHong Zhang   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_getrow_uppertriangular",&flg);CHKERRQ(ierr);
1639*f5edf698SHong Zhang   if (flg) b->getrow_utriangular = PETSC_TRUE;
1640*f5edf698SHong Zhang 
1641a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C",
1642a23d5eceSKris Buschelman                                      "MatStoreValues_SeqSBAIJ",
1643a23d5eceSKris Buschelman                                      MatStoreValues_SeqSBAIJ);CHKERRQ(ierr);
1644a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C",
1645a23d5eceSKris Buschelman                                      "MatRetrieveValues_SeqSBAIJ",
1646a23d5eceSKris Buschelman                                      (void*)MatRetrieveValues_SeqSBAIJ);CHKERRQ(ierr);
1647a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqSBAIJSetColumnIndices_C",
1648a23d5eceSKris Buschelman                                      "MatSeqSBAIJSetColumnIndices_SeqSBAIJ",
1649a23d5eceSKris Buschelman                                      MatSeqSBAIJSetColumnIndices_SeqSBAIJ);CHKERRQ(ierr);
16504e5e7fe4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqsbaij_seqaij_C",
16514e5e7fe4SHong Zhang                                      "MatConvert_SeqSBAIJ_SeqAIJ",
16524e5e7fe4SHong Zhang                                       MatConvert_SeqSBAIJ_SeqAIJ);CHKERRQ(ierr);
1653a0e1a404SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqsbaij_seqbaij_C",
1654a0e1a404SHong Zhang                                      "MatConvert_SeqSBAIJ_SeqBAIJ",
1655a0e1a404SHong Zhang                                       MatConvert_SeqSBAIJ_SeqBAIJ);CHKERRQ(ierr);
1656a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqSBAIJSetPreallocation_C",
1657a23d5eceSKris Buschelman                                      "MatSeqSBAIJSetPreallocation_SeqSBAIJ",
1658a23d5eceSKris Buschelman                                      MatSeqSBAIJSetPreallocation_SeqSBAIJ);CHKERRQ(ierr);
165923ce1328SBarry Smith 
166023ce1328SBarry Smith   B->symmetric                  = PETSC_TRUE;
166123ce1328SBarry Smith   B->structurally_symmetric     = PETSC_TRUE;
166223ce1328SBarry Smith   B->symmetric_set              = PETSC_TRUE;
166323ce1328SBarry Smith   B->structurally_symmetric_set = PETSC_TRUE;
1664a23d5eceSKris Buschelman   PetscFunctionReturn(0);
1665a23d5eceSKris Buschelman }
1666a23d5eceSKris Buschelman EXTERN_C_END
1667a23d5eceSKris Buschelman 
1668a23d5eceSKris Buschelman #undef __FUNCT__
1669a23d5eceSKris Buschelman #define __FUNCT__ "MatSeqSBAIJSetPreallocation"
1670a23d5eceSKris Buschelman /*@C
1671a23d5eceSKris Buschelman    MatSeqSBAIJSetPreallocation - Creates a sparse symmetric matrix in block AIJ (block
1672a23d5eceSKris Buschelman    compressed row) format.  For good matrix assembly performance the
1673a23d5eceSKris Buschelman    user should preallocate the matrix storage by setting the parameter nz
1674a23d5eceSKris Buschelman    (or the array nnz).  By setting these parameters accurately, performance
1675a23d5eceSKris Buschelman    during matrix assembly can be increased by more than a factor of 50.
1676a23d5eceSKris Buschelman 
1677a23d5eceSKris Buschelman    Collective on Mat
1678a23d5eceSKris Buschelman 
1679a23d5eceSKris Buschelman    Input Parameters:
1680a23d5eceSKris Buschelman +  A - the symmetric matrix
1681a23d5eceSKris Buschelman .  bs - size of block
1682a23d5eceSKris Buschelman .  nz - number of block nonzeros per block row (same for all rows)
1683a23d5eceSKris Buschelman -  nnz - array containing the number of block nonzeros in the upper triangular plus
1684a23d5eceSKris Buschelman          diagonal portion of each block (possibly different for each block row) or PETSC_NULL
1685a23d5eceSKris Buschelman 
1686a23d5eceSKris Buschelman    Options Database Keys:
1687a23d5eceSKris Buschelman .   -mat_no_unroll - uses code that does not unroll the loops in the
1688a23d5eceSKris Buschelman                      block calculations (much slower)
1689a23d5eceSKris Buschelman .    -mat_block_size - size of the blocks to use
1690a23d5eceSKris Buschelman 
1691a23d5eceSKris Buschelman    Level: intermediate
1692a23d5eceSKris Buschelman 
1693a23d5eceSKris Buschelman    Notes:
1694a23d5eceSKris Buschelman    Specify the preallocated storage with either nz or nnz (not both).
1695a23d5eceSKris Buschelman    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
1696a23d5eceSKris Buschelman    allocation.  For additional details, see the users manual chapter on
1697a23d5eceSKris Buschelman    matrices.
1698a23d5eceSKris Buschelman 
169949a6f317SBarry Smith    If the nnz parameter is given then the nz parameter is ignored
170049a6f317SBarry Smith 
170149a6f317SBarry Smith 
1702a23d5eceSKris Buschelman .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ()
1703a23d5eceSKris Buschelman @*/
1704be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetPreallocation(Mat B,PetscInt bs,PetscInt nz,const PetscInt nnz[])
170513f74950SBarry Smith {
170613f74950SBarry Smith   PetscErrorCode ierr,(*f)(Mat,PetscInt,PetscInt,const PetscInt[]);
1707a23d5eceSKris Buschelman 
1708a23d5eceSKris Buschelman   PetscFunctionBegin;
1709a23d5eceSKris Buschelman   ierr = PetscObjectQueryFunction((PetscObject)B,"MatSeqSBAIJSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr);
1710a23d5eceSKris Buschelman   if (f) {
1711a23d5eceSKris Buschelman     ierr = (*f)(B,bs,nz,nnz);CHKERRQ(ierr);
1712a23d5eceSKris Buschelman   }
1713a23d5eceSKris Buschelman   PetscFunctionReturn(0);
1714a23d5eceSKris Buschelman }
171549b5e25fSSatish Balay 
17164a2ae208SSatish Balay #undef __FUNCT__
17174a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqSBAIJ"
1718c464158bSHong Zhang /*@C
1719c464158bSHong Zhang    MatCreateSeqSBAIJ - Creates a sparse symmetric matrix in block AIJ (block
1720c464158bSHong Zhang    compressed row) format.  For good matrix assembly performance the
1721c464158bSHong Zhang    user should preallocate the matrix storage by setting the parameter nz
1722c464158bSHong Zhang    (or the array nnz).  By setting these parameters accurately, performance
1723c464158bSHong Zhang    during matrix assembly can be increased by more than a factor of 50.
172449b5e25fSSatish Balay 
1725c464158bSHong Zhang    Collective on MPI_Comm
1726c464158bSHong Zhang 
1727c464158bSHong Zhang    Input Parameters:
1728c464158bSHong Zhang +  comm - MPI communicator, set to PETSC_COMM_SELF
1729c464158bSHong Zhang .  bs - size of block
1730c464158bSHong Zhang .  m - number of rows, or number of columns
1731c464158bSHong Zhang .  nz - number of block nonzeros per block row (same for all rows)
1732744e8345SSatish Balay -  nnz - array containing the number of block nonzeros in the upper triangular plus
1733744e8345SSatish Balay          diagonal portion of each block (possibly different for each block row) or PETSC_NULL
1734c464158bSHong Zhang 
1735c464158bSHong Zhang    Output Parameter:
1736c464158bSHong Zhang .  A - the symmetric matrix
1737c464158bSHong Zhang 
1738c464158bSHong Zhang    Options Database Keys:
1739c464158bSHong Zhang .   -mat_no_unroll - uses code that does not unroll the loops in the
1740c464158bSHong Zhang                      block calculations (much slower)
1741c464158bSHong Zhang .    -mat_block_size - size of the blocks to use
1742c464158bSHong Zhang 
1743c464158bSHong Zhang    Level: intermediate
1744c464158bSHong Zhang 
1745c464158bSHong Zhang    Notes:
1746c464158bSHong Zhang 
1747c464158bSHong Zhang    Specify the preallocated storage with either nz or nnz (not both).
1748c464158bSHong Zhang    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
1749c464158bSHong Zhang    allocation.  For additional details, see the users manual chapter on
1750c464158bSHong Zhang    matrices.
1751c464158bSHong Zhang 
175249a6f317SBarry Smith    If the nnz parameter is given then the nz parameter is ignored
175349a6f317SBarry Smith 
1754c464158bSHong Zhang .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ()
1755c464158bSHong Zhang @*/
1756be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJ(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],Mat *A)
1757c464158bSHong Zhang {
1758dfbe8321SBarry Smith   PetscErrorCode ierr;
1759c464158bSHong Zhang 
1760c464158bSHong Zhang   PetscFunctionBegin;
1761f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,A);CHKERRQ(ierr);
1762f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr);
1763c464158bSHong Zhang   ierr = MatSetType(*A,MATSEQSBAIJ);CHKERRQ(ierr);
1764ab93d7beSBarry Smith   ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(*A,bs,nz,(PetscInt*)nnz);CHKERRQ(ierr);
176549b5e25fSSatish Balay   PetscFunctionReturn(0);
176649b5e25fSSatish Balay }
176749b5e25fSSatish Balay 
17684a2ae208SSatish Balay #undef __FUNCT__
17694a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqSBAIJ"
1770dfbe8321SBarry Smith PetscErrorCode MatDuplicate_SeqSBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B)
177149b5e25fSSatish Balay {
177249b5e25fSSatish Balay   Mat            C;
177349b5e25fSSatish Balay   Mat_SeqSBAIJ   *c,*a = (Mat_SeqSBAIJ*)A->data;
17746849ba73SBarry Smith   PetscErrorCode ierr;
1775b40805acSSatish Balay   PetscInt       i,mbs = a->mbs,nz = a->nz,bs2 =a->bs2;
177649b5e25fSSatish Balay 
177749b5e25fSSatish Balay   PetscFunctionBegin;
177829bbc08cSBarry Smith   if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,"Corrupt matrix");
177949b5e25fSSatish Balay 
178049b5e25fSSatish Balay   *B = 0;
1781f69a0ea3SMatthew Knepley   ierr = MatCreate(A->comm,&C);CHKERRQ(ierr);
1782f69a0ea3SMatthew Knepley   ierr = MatSetSizes(C,A->m,A->n,A->m,A->n);CHKERRQ(ierr);
1783be5d1d56SKris Buschelman   ierr = MatSetType(C,A->type_name);CHKERRQ(ierr);
17841d5dac46SHong Zhang   ierr = PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
1785692f9cbeSHong Zhang   c    = (Mat_SeqSBAIJ*)C->data;
1786692f9cbeSHong Zhang 
1787273d9f13SBarry Smith   C->preallocated   = PETSC_TRUE;
178849b5e25fSSatish Balay   C->factor         = A->factor;
178949b5e25fSSatish Balay   c->row            = 0;
179049b5e25fSSatish Balay   c->icol           = 0;
179149b5e25fSSatish Balay   c->saved_values   = 0;
179249b5e25fSSatish Balay   c->keepzeroedrows = a->keepzeroedrows;
179349b5e25fSSatish Balay   C->assembled      = PETSC_TRUE;
179449b5e25fSSatish Balay 
179582327fa8SHong Zhang   C->M    = A->M;
179682327fa8SHong Zhang   C->N    = A->N;
1797521d7252SBarry Smith   C->bs   = A->bs;
179849b5e25fSSatish Balay   c->bs2  = a->bs2;
179949b5e25fSSatish Balay   c->mbs  = a->mbs;
180049b5e25fSSatish Balay   c->nbs  = a->nbs;
180149b5e25fSSatish Balay 
180213f74950SBarry Smith   ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&c->imax);CHKERRQ(ierr);
180313f74950SBarry Smith   ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&c->ilen);CHKERRQ(ierr);
180449b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
180549b5e25fSSatish Balay     c->imax[i] = a->imax[i];
180649b5e25fSSatish Balay     c->ilen[i] = a->ilen[i];
180749b5e25fSSatish Balay   }
1808b40805acSSatish Balay   ierr = PetscLogObjectMemory(C,2*(mbs+1)*sizeof(PetscInt)+sizeof(struct _p_Mat)+sizeof(Mat_SeqSBAIJ));CHKERRQ(ierr);
180949b5e25fSSatish Balay 
181049b5e25fSSatish Balay   /* allocate the matrix space */
1811b40805acSSatish Balay   ierr = PetscMalloc3(bs2*nz,MatScalar,&c->a,nz,PetscInt,&c->j,mbs+1,PetscInt,&c->i);CHKERRQ(ierr);
181249b5e25fSSatish Balay   c->singlemalloc = PETSC_TRUE;
181313f74950SBarry Smith   ierr = PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(PetscInt));CHKERRQ(ierr);
1814b40805acSSatish Balay   ierr = PetscLogObjectMemory(C,(mbs+1)*sizeof(PetscInt) + nz*(bs2*sizeof(MatScalar) + sizeof(PetscInt)));CHKERRQ(ierr);
181549b5e25fSSatish Balay   if (mbs > 0) {
181613f74950SBarry Smith     ierr = PetscMemcpy(c->j,a->j,nz*sizeof(PetscInt));CHKERRQ(ierr);
181749b5e25fSSatish Balay     if (cpvalues == MAT_COPY_VALUES) {
181849b5e25fSSatish Balay       ierr = PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr);
181949b5e25fSSatish Balay     } else {
182049b5e25fSSatish Balay       ierr = PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr);
182149b5e25fSSatish Balay     }
182249b5e25fSSatish Balay   }
182349b5e25fSSatish Balay 
182449b5e25fSSatish Balay   c->sorted      = a->sorted;
182549b5e25fSSatish Balay   c->roworiented = a->roworiented;
182649b5e25fSSatish Balay   c->nonew       = a->nonew;
182749b5e25fSSatish Balay 
182849b5e25fSSatish Balay   if (a->diag) {
182913f74950SBarry Smith     ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&c->diag);CHKERRQ(ierr);
183052e6d16bSBarry Smith     ierr = PetscLogObjectMemory(C,(mbs+1)*sizeof(PetscInt));CHKERRQ(ierr);
183149b5e25fSSatish Balay     for (i=0; i<mbs; i++) {
183249b5e25fSSatish Balay       c->diag[i] = a->diag[i];
183349b5e25fSSatish Balay     }
183449b5e25fSSatish Balay   } else c->diag        = 0;
18356c6c5352SBarry Smith   c->nz               = a->nz;
18366c6c5352SBarry Smith   c->maxnz            = a->maxnz;
183749b5e25fSSatish Balay   c->solve_work         = 0;
183849b5e25fSSatish Balay   c->mult_work          = 0;
183949b5e25fSSatish Balay   *B = C;
1840b0a32e0cSBarry Smith   ierr = PetscFListDuplicate(A->qlist,&C->qlist);CHKERRQ(ierr);
184149b5e25fSSatish Balay   PetscFunctionReturn(0);
184249b5e25fSSatish Balay }
184349b5e25fSSatish Balay 
18444a2ae208SSatish Balay #undef __FUNCT__
18454a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqSBAIJ"
1846f69a0ea3SMatthew Knepley PetscErrorCode MatLoad_SeqSBAIJ(PetscViewer viewer, MatType type,Mat *A)
184749b5e25fSSatish Balay {
184849b5e25fSSatish Balay   Mat_SeqSBAIJ   *a;
184949b5e25fSSatish Balay   Mat            B;
18506849ba73SBarry Smith   PetscErrorCode ierr;
185113f74950SBarry Smith   int            fd;
185213f74950SBarry Smith   PetscMPIInt    size;
185313f74950SBarry Smith   PetscInt       i,nz,header[4],*rowlengths=0,M,N,bs=1;
185413f74950SBarry Smith   PetscInt       *mask,mbs,*jj,j,rowcount,nzcount,k,*s_browlengths,maskcount;
185513f74950SBarry Smith   PetscInt       kmax,jcount,block,idx,point,nzcountb,extra_rows;
185613f74950SBarry Smith   PetscInt       *masked,nmask,tmp,bs2,ishift;
185787828ca2SBarry Smith   PetscScalar    *aa;
185849b5e25fSSatish Balay   MPI_Comm       comm = ((PetscObject)viewer)->comm;
185949b5e25fSSatish Balay 
186049b5e25fSSatish Balay   PetscFunctionBegin;
1861b0a32e0cSBarry Smith   ierr = PetscOptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,PETSC_NULL);CHKERRQ(ierr);
186249b5e25fSSatish Balay   bs2  = bs*bs;
186349b5e25fSSatish Balay 
186449b5e25fSSatish Balay   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
186529bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor");
1866b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
186749b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr);
1868552e946dSBarry Smith   if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not Mat object");
186949b5e25fSSatish Balay   M = header[1]; N = header[2]; nz = header[3];
187049b5e25fSSatish Balay 
187149b5e25fSSatish Balay   if (header[3] < 0) {
187229bbc08cSBarry Smith     SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as SeqSBAIJ");
187349b5e25fSSatish Balay   }
187449b5e25fSSatish Balay 
187529bbc08cSBarry Smith   if (M != N) SETERRQ(PETSC_ERR_SUP,"Can only do square matrices");
187649b5e25fSSatish Balay 
187749b5e25fSSatish Balay   /*
187849b5e25fSSatish Balay      This code adds extra rows to make sure the number of rows is
187949b5e25fSSatish Balay     divisible by the blocksize
188049b5e25fSSatish Balay   */
188149b5e25fSSatish Balay   mbs        = M/bs;
188249b5e25fSSatish Balay   extra_rows = bs - M + bs*(mbs);
188349b5e25fSSatish Balay   if (extra_rows == bs) extra_rows = 0;
188449b5e25fSSatish Balay   else                  mbs++;
188549b5e25fSSatish Balay   if (extra_rows) {
188663ba0a88SBarry Smith     ierr = PetscLogInfo((0,"MatLoad_SeqSBAIJ:Padding loaded matrix to match blocksize\n"));CHKERRQ(ierr);
188749b5e25fSSatish Balay   }
188849b5e25fSSatish Balay 
188949b5e25fSSatish Balay   /* read in row lengths */
189013f74950SBarry Smith   ierr = PetscMalloc((M+extra_rows)*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr);
189149b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
189249b5e25fSSatish Balay   for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1;
189349b5e25fSSatish Balay 
189449b5e25fSSatish Balay   /* read in column indices */
189513f74950SBarry Smith   ierr = PetscMalloc((nz+extra_rows)*sizeof(PetscInt),&jj);CHKERRQ(ierr);
189649b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT);CHKERRQ(ierr);
189749b5e25fSSatish Balay   for (i=0; i<extra_rows; i++) jj[nz+i] = M+i;
189849b5e25fSSatish Balay 
189949b5e25fSSatish Balay   /* loop over row lengths determining block row lengths */
190013f74950SBarry Smith   ierr     = PetscMalloc(mbs*sizeof(PetscInt),&s_browlengths);CHKERRQ(ierr);
190113f74950SBarry Smith   ierr     = PetscMemzero(s_browlengths,mbs*sizeof(PetscInt));CHKERRQ(ierr);
190213f74950SBarry Smith   ierr     = PetscMalloc(2*mbs*sizeof(PetscInt),&mask);CHKERRQ(ierr);
190313f74950SBarry Smith   ierr     = PetscMemzero(mask,mbs*sizeof(PetscInt));CHKERRQ(ierr);
190449b5e25fSSatish Balay   masked   = mask + mbs;
190549b5e25fSSatish Balay   rowcount = 0; nzcount = 0;
190649b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
190749b5e25fSSatish Balay     nmask = 0;
190849b5e25fSSatish Balay     for (j=0; j<bs; j++) {
190949b5e25fSSatish Balay       kmax = rowlengths[rowcount];
191049b5e25fSSatish Balay       for (k=0; k<kmax; k++) {
19112d703238SHong Zhang         tmp = jj[nzcount++]/bs;   /* block col. index */
191203630b6eSHong Zhang         if (!mask[tmp] && tmp >= i) {masked[nmask++] = tmp; mask[tmp] = 1;}
191349b5e25fSSatish Balay       }
191449b5e25fSSatish Balay       rowcount++;
191549b5e25fSSatish Balay     }
1916574b2666SHong Zhang     s_browlengths[i] += nmask;
1917574b2666SHong Zhang 
191849b5e25fSSatish Balay     /* zero out the mask elements we set */
191949b5e25fSSatish Balay     for (j=0; j<nmask; j++) mask[masked[j]] = 0;
192049b5e25fSSatish Balay   }
192149b5e25fSSatish Balay 
192249b5e25fSSatish Balay   /* create our matrix */
1923f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,&B);CHKERRQ(ierr);
1924f69a0ea3SMatthew Knepley   ierr = MatSetSizes(B,M+extra_rows,N+extra_rows,M+extra_rows,N+extra_rows);CHKERRQ(ierr);
19259abb65ffSKris Buschelman   ierr = MatSetType(B,type);CHKERRQ(ierr);
1926ab93d7beSBarry Smith   ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(B,bs,0,s_browlengths);CHKERRQ(ierr);
192749b5e25fSSatish Balay   a = (Mat_SeqSBAIJ*)B->data;
192849b5e25fSSatish Balay 
192949b5e25fSSatish Balay   /* set matrix "i" values */
193049b5e25fSSatish Balay   a->i[0] = 0;
193149b5e25fSSatish Balay   for (i=1; i<= mbs; i++) {
1932574b2666SHong Zhang     a->i[i]      = a->i[i-1] + s_browlengths[i-1];
1933574b2666SHong Zhang     a->ilen[i-1] = s_browlengths[i-1];
193449b5e25fSSatish Balay   }
19356c6c5352SBarry Smith   a->nz = a->i[mbs];
193649b5e25fSSatish Balay 
193749b5e25fSSatish Balay   /* read in nonzero values */
193887828ca2SBarry Smith   ierr = PetscMalloc((nz+extra_rows)*sizeof(PetscScalar),&aa);CHKERRQ(ierr);
193949b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);CHKERRQ(ierr);
194049b5e25fSSatish Balay   for (i=0; i<extra_rows; i++) aa[nz+i] = 1.0;
194149b5e25fSSatish Balay 
194249b5e25fSSatish Balay   /* set "a" and "j" values into matrix */
194349b5e25fSSatish Balay   nzcount = 0; jcount = 0;
194449b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
194549b5e25fSSatish Balay     nzcountb = nzcount;
194649b5e25fSSatish Balay     nmask    = 0;
194749b5e25fSSatish Balay     for (j=0; j<bs; j++) {
194849b5e25fSSatish Balay       kmax = rowlengths[i*bs+j];
194949b5e25fSSatish Balay       for (k=0; k<kmax; k++) {
19502d703238SHong Zhang         tmp = jj[nzcount++]/bs; /* block col. index */
195103630b6eSHong Zhang         if (!mask[tmp] && tmp >= i) { masked[nmask++] = tmp; mask[tmp] = 1;}
19522d703238SHong Zhang       }
19532d703238SHong Zhang     }
19542d703238SHong Zhang     /* sort the masked values */
19552d703238SHong Zhang     ierr = PetscSortInt(nmask,masked);CHKERRQ(ierr);
19562d703238SHong Zhang 
19572d703238SHong Zhang     /* set "j" values into matrix */
19582d703238SHong Zhang     maskcount = 1;
19592d703238SHong Zhang     for (j=0; j<nmask; j++) {
196049b5e25fSSatish Balay       a->j[jcount++]  = masked[j];
196149b5e25fSSatish Balay       mask[masked[j]] = maskcount++;
196249b5e25fSSatish Balay     }
1963574b2666SHong Zhang 
196449b5e25fSSatish Balay     /* set "a" values into matrix */
196549b5e25fSSatish Balay     ishift = bs2*a->i[i];
196649b5e25fSSatish Balay     for (j=0; j<bs; j++) {
196749b5e25fSSatish Balay       kmax = rowlengths[i*bs+j];
196849b5e25fSSatish Balay       for (k=0; k<kmax; k++) {
1969574b2666SHong Zhang         tmp       = jj[nzcountb]/bs ; /* block col. index */
1970574b2666SHong Zhang         if (tmp >= i){
197149b5e25fSSatish Balay           block     = mask[tmp] - 1;
197249b5e25fSSatish Balay           point     = jj[nzcountb] - bs*tmp;
197349b5e25fSSatish Balay           idx       = ishift + bs2*block + j + bs*point;
1974574b2666SHong Zhang           a->a[idx] = aa[nzcountb];
1975574b2666SHong Zhang         }
1976574b2666SHong Zhang         nzcountb++;
197749b5e25fSSatish Balay       }
197849b5e25fSSatish Balay     }
197949b5e25fSSatish Balay     /* zero out the mask elements we set */
198049b5e25fSSatish Balay     for (j=0; j<nmask; j++) mask[masked[j]] = 0;
198149b5e25fSSatish Balay   }
19826c6c5352SBarry Smith   if (jcount != a->nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Bad binary matrix");
198349b5e25fSSatish Balay 
198449b5e25fSSatish Balay   ierr = PetscFree(rowlengths);CHKERRQ(ierr);
1985574b2666SHong Zhang   ierr = PetscFree(s_browlengths);CHKERRQ(ierr);
198649b5e25fSSatish Balay   ierr = PetscFree(aa);CHKERRQ(ierr);
198749b5e25fSSatish Balay   ierr = PetscFree(jj);CHKERRQ(ierr);
198849b5e25fSSatish Balay   ierr = PetscFree(mask);CHKERRQ(ierr);
198949b5e25fSSatish Balay 
19909abb65ffSKris Buschelman   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
19919abb65ffSKris Buschelman   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
199249b5e25fSSatish Balay   ierr = MatView_Private(B);CHKERRQ(ierr);
19939abb65ffSKris Buschelman   *A = B;
199449b5e25fSSatish Balay   PetscFunctionReturn(0);
199549b5e25fSSatish Balay }
1996574b2666SHong Zhang 
1997d06b337dSHong Zhang #undef __FUNCT__
1998d06b337dSHong Zhang #define __FUNCT__ "MatRelax_SeqSBAIJ"
199913f74950SBarry Smith PetscErrorCode MatRelax_SeqSBAIJ(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx)
2000d06b337dSHong Zhang {
2001d06b337dSHong Zhang   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
2002d06b337dSHong Zhang   MatScalar      *aa=a->a,*v,*v1;
2003d06b337dSHong Zhang   PetscScalar    *x,*b,*t,sum,d;
20046849ba73SBarry Smith   PetscErrorCode ierr;
2005521d7252SBarry Smith   PetscInt       m=a->mbs,bs=A->bs,*ai=a->i,*aj=a->j;
2006521d7252SBarry Smith   PetscInt       nz,nz1,*vj,*vj1,i;
2007d06b337dSHong Zhang 
2008d06b337dSHong Zhang   PetscFunctionBegin;
2009b965ef7fSBarry Smith   its = its*lits;
201077431f27SBarry Smith   if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits);
2011d06b337dSHong Zhang 
2012d06b337dSHong Zhang   if (bs > 1)
2013d06b337dSHong Zhang     SETERRQ(PETSC_ERR_SUP,"SSOR for block size > 1 is not yet implemented");
2014d06b337dSHong Zhang 
20151ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2016d06b337dSHong Zhang   if (xx != bb) {
20171ebc52fbSHong Zhang     ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
2018d06b337dSHong Zhang   } else {
2019d06b337dSHong Zhang     b = x;
2020d06b337dSHong Zhang   }
2021d06b337dSHong Zhang 
2022d06b337dSHong Zhang   ierr = PetscMalloc(m*sizeof(PetscScalar),&t);CHKERRQ(ierr);
2023d06b337dSHong Zhang 
2024d06b337dSHong Zhang   if (flag & SOR_ZERO_INITIAL_GUESS) {
20252798e883SHong Zhang     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
2026d06b337dSHong Zhang       for (i=0; i<m; i++)
2027d06b337dSHong Zhang         t[i] = b[i];
2028d06b337dSHong Zhang 
2029d06b337dSHong Zhang       for (i=0; i<m; i++){
203044706875SHong Zhang         d  = *(aa + ai[i]);  /* diag[i] */
2031d06b337dSHong Zhang         v  = aa + ai[i] + 1;
2032d06b337dSHong Zhang         vj = aj + ai[i] + 1;
2033d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1;
2034d06b337dSHong Zhang         x[i] = omega*t[i]/d;
2035d06b337dSHong Zhang         while (nz--) t[*vj++] -= x[i]*(*v++); /* update rhs */
2036efee365bSSatish Balay         ierr = PetscLogFlops(2*nz-1);CHKERRQ(ierr);
2037d06b337dSHong Zhang       }
2038d06b337dSHong Zhang     }
2039d06b337dSHong Zhang 
20402798e883SHong Zhang     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
2041d06b337dSHong Zhang       for (i=0; i<m; i++)
2042d06b337dSHong Zhang         t[i] = b[i];
2043d06b337dSHong Zhang 
2044d06b337dSHong Zhang       for (i=0; i<m-1; i++){  /* update rhs */
2045d06b337dSHong Zhang         v  = aa + ai[i] + 1;
2046d06b337dSHong Zhang         vj = aj + ai[i] + 1;
2047d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1;
2048d06b337dSHong Zhang         while (nz--) t[*vj++] -= x[i]*(*v++);
2049efee365bSSatish Balay         ierr = PetscLogFlops(2*nz-1);CHKERRQ(ierr);
2050d06b337dSHong Zhang       }
2051d06b337dSHong Zhang       for (i=m-1; i>=0; i--){
2052d06b337dSHong Zhang         d  = *(aa + ai[i]);
2053d06b337dSHong Zhang         v  = aa + ai[i] + 1;
2054d06b337dSHong Zhang         vj = aj + ai[i] + 1;
2055d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1;
2056d06b337dSHong Zhang         sum = t[i];
2057d06b337dSHong Zhang         while (nz--) sum -= x[*vj++]*(*v++);
2058efee365bSSatish Balay         ierr = PetscLogFlops(2*nz-1);CHKERRQ(ierr);
2059d06b337dSHong Zhang         x[i] =   (1-omega)*x[i] + omega*sum/d;
2060d06b337dSHong Zhang       }
2061d06b337dSHong Zhang     }
2062d06b337dSHong Zhang     its--;
2063d06b337dSHong Zhang   }
2064d06b337dSHong Zhang 
2065d06b337dSHong Zhang   while (its--) {
206644706875SHong Zhang     /*
206744706875SHong Zhang        forward sweep:
206844706875SHong Zhang        for i=0,...,m-1:
206944706875SHong Zhang          sum[i] = (b[i] - U(i,:)x )/d[i];
207044706875SHong Zhang          x[i]   = (1-omega)x[i] + omega*sum[i];
207144706875SHong Zhang          b      = b - x[i]*U^T(i,:);
2072d06b337dSHong Zhang 
207344706875SHong Zhang     */
20742798e883SHong Zhang     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
2075d06b337dSHong Zhang       for (i=0; i<m; i++)
2076d06b337dSHong Zhang         t[i] = b[i];
2077d06b337dSHong Zhang 
2078d06b337dSHong Zhang       for (i=0; i<m; i++){
207944706875SHong Zhang         d  = *(aa + ai[i]);  /* diag[i] */
2080d06b337dSHong Zhang         v  = aa + ai[i] + 1; v1=v;
2081d06b337dSHong Zhang         vj = aj + ai[i] + 1; vj1=vj;
2082d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1; nz1=nz;
2083d06b337dSHong Zhang         sum = t[i];
2084d06b337dSHong Zhang         while (nz1--) sum -= (*v1++)*x[*vj1++];
2085d06b337dSHong Zhang         x[i] = (1-omega)*x[i] + omega*sum/d;
2086d06b337dSHong Zhang         while (nz--) t[*vj++] -= x[i]*(*v++);
2087efee365bSSatish Balay         ierr = PetscLogFlops(4*nz-2);CHKERRQ(ierr);
2088d06b337dSHong Zhang       }
2089d06b337dSHong Zhang     }
2090d06b337dSHong Zhang 
20912798e883SHong Zhang   if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
209244706875SHong Zhang       /*
209344706875SHong Zhang        backward sweep:
209444706875SHong Zhang        b = b - x[i]*U^T(i,:), i=0,...,n-2
209544706875SHong Zhang        for i=m-1,...,0:
209644706875SHong Zhang          sum[i] = (b[i] - U(i,:)x )/d[i];
209744706875SHong Zhang          x[i]   = (1-omega)x[i] + omega*sum[i];
209844706875SHong Zhang       */
2099d06b337dSHong Zhang       for (i=0; i<m; i++)
2100d06b337dSHong Zhang         t[i] = b[i];
2101d06b337dSHong Zhang 
2102d06b337dSHong Zhang       for (i=0; i<m-1; i++){  /* update rhs */
2103d06b337dSHong Zhang         v  = aa + ai[i] + 1;
2104d06b337dSHong Zhang         vj = aj + ai[i] + 1;
2105d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1;
2106d06b337dSHong Zhang         while (nz--) t[*vj++] -= x[i]*(*v++);
2107efee365bSSatish Balay         ierr = PetscLogFlops(2*nz-1);CHKERRQ(ierr);
2108d06b337dSHong Zhang       }
2109d06b337dSHong Zhang       for (i=m-1; i>=0; i--){
2110d06b337dSHong Zhang         d  = *(aa + ai[i]);
2111d06b337dSHong Zhang         v  = aa + ai[i] + 1;
2112d06b337dSHong Zhang         vj = aj + ai[i] + 1;
2113d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1;
2114d06b337dSHong Zhang         sum = t[i];
2115d06b337dSHong Zhang         while (nz--) sum -= x[*vj++]*(*v++);
2116efee365bSSatish Balay         ierr = PetscLogFlops(2*nz-1);CHKERRQ(ierr);
2117d06b337dSHong Zhang         x[i] =   (1-omega)*x[i] + omega*sum/d;
2118d06b337dSHong Zhang       }
2119d06b337dSHong Zhang     }
2120d06b337dSHong Zhang   }
2121d06b337dSHong Zhang 
2122d06b337dSHong Zhang   ierr = PetscFree(t);CHKERRQ(ierr);
21231ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2124d06b337dSHong Zhang   if (bb != xx) {
21251ebc52fbSHong Zhang     ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
2126d06b337dSHong Zhang   }
2127d06b337dSHong Zhang   PetscFunctionReturn(0);
2128d06b337dSHong Zhang }
2129d06b337dSHong Zhang 
2130d06b337dSHong Zhang 
2131d06b337dSHong Zhang 
2132d06b337dSHong Zhang 
213349b5e25fSSatish Balay 
213449b5e25fSSatish Balay 
2135