xref: /petsc/src/mat/impls/sbaij/seq/sbaij.c (revision ed1caa07d30b4bf52fa72dab33adf574bdef2af0)
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;
1884d9d31abSKris Buschelman   default:
18929bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"unknown option");
19049b5e25fSSatish Balay   }
19149b5e25fSSatish Balay   PetscFunctionReturn(0);
19249b5e25fSSatish Balay }
19349b5e25fSSatish Balay 
1944a2ae208SSatish Balay #undef __FUNCT__
1954a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqSBAIJ"
19613f74950SBarry Smith PetscErrorCode MatGetRow_SeqSBAIJ(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **v)
19749b5e25fSSatish Balay {
19849b5e25fSSatish Balay   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
1996849ba73SBarry Smith   PetscErrorCode ierr;
20013f74950SBarry Smith   PetscInt       itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*cols_i,bs2;
20149b5e25fSSatish Balay   MatScalar      *aa,*aa_i;
20287828ca2SBarry Smith   PetscScalar    *v_i;
20349b5e25fSSatish Balay 
20449b5e25fSSatish Balay   PetscFunctionBegin;
205521d7252SBarry Smith   bs  = A->bs;
20649b5e25fSSatish Balay   ai  = a->i;
20749b5e25fSSatish Balay   aj  = a->j;
20849b5e25fSSatish Balay   aa  = a->a;
20949b5e25fSSatish Balay   bs2 = a->bs2;
21049b5e25fSSatish Balay 
21177431f27SBarry Smith   if (row < 0 || row >= A->m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE, "Row %D out of range", row);
21249b5e25fSSatish Balay 
21349b5e25fSSatish Balay   bn  = row/bs;   /* Block number */
21449b5e25fSSatish Balay   bp  = row % bs; /* Block position */
21549b5e25fSSatish Balay   M   = ai[bn+1] - ai[bn];
21649b5e25fSSatish Balay   *ncols = bs*M;
21749b5e25fSSatish Balay 
21849b5e25fSSatish Balay   if (v) {
21949b5e25fSSatish Balay     *v = 0;
22049b5e25fSSatish Balay     if (*ncols) {
22187828ca2SBarry Smith       ierr = PetscMalloc((*ncols+row)*sizeof(PetscScalar),v);CHKERRQ(ierr);
22249b5e25fSSatish Balay       for (i=0; i<M; i++) { /* for each block in the block row */
22349b5e25fSSatish Balay         v_i  = *v + i*bs;
22449b5e25fSSatish Balay         aa_i = aa + bs2*(ai[bn] + i);
22549b5e25fSSatish Balay         for (j=bp,k=0; j<bs2; j+=bs,k++) {v_i[k] = aa_i[j];}
22649b5e25fSSatish Balay       }
22749b5e25fSSatish Balay     }
22849b5e25fSSatish Balay   }
22949b5e25fSSatish Balay 
23049b5e25fSSatish Balay   if (cols) {
23149b5e25fSSatish Balay     *cols = 0;
23249b5e25fSSatish Balay     if (*ncols) {
23313f74950SBarry Smith       ierr = PetscMalloc((*ncols+row)*sizeof(PetscInt),cols);CHKERRQ(ierr);
23449b5e25fSSatish Balay       for (i=0; i<M; i++) { /* for each block in the block row */
23549b5e25fSSatish Balay         cols_i = *cols + i*bs;
23649b5e25fSSatish Balay         itmp  = bs*aj[ai[bn] + i];
23749b5e25fSSatish Balay         for (j=0; j<bs; j++) {cols_i[j] = itmp++;}
23849b5e25fSSatish Balay       }
23949b5e25fSSatish Balay     }
24049b5e25fSSatish Balay   }
24149b5e25fSSatish Balay 
24249b5e25fSSatish Balay   /*search column A(0:row-1,row) (=A(row,0:row-1)). Could be expensive! */
2435ddb2528SHong Zhang   /* this segment is currently removed, so only entries in the upper triangle are obtained */
2445ddb2528SHong Zhang #ifdef column_search
24549b5e25fSSatish Balay   v_i    = *v    + M*bs;
24649b5e25fSSatish Balay   cols_i = *cols + M*bs;
24749b5e25fSSatish Balay   for (i=0; i<bn; i++){ /* for each block row */
24849b5e25fSSatish Balay     M = ai[i+1] - ai[i];
24949b5e25fSSatish Balay     for (j=0; j<M; j++){
25049b5e25fSSatish Balay       itmp = aj[ai[i] + j];    /* block column value */
25149b5e25fSSatish Balay       if (itmp == bn){
25249b5e25fSSatish Balay         aa_i   = aa    + bs2*(ai[i] + j) + bs*bp;
25349b5e25fSSatish Balay         for (k=0; k<bs; k++) {
25449b5e25fSSatish Balay           *cols_i++ = i*bs+k;
25549b5e25fSSatish Balay           *v_i++    = aa_i[k];
25649b5e25fSSatish Balay         }
25749b5e25fSSatish Balay         *ncols += bs;
25849b5e25fSSatish Balay         break;
25949b5e25fSSatish Balay       }
26049b5e25fSSatish Balay     }
26149b5e25fSSatish Balay   }
2625ddb2528SHong Zhang #endif
26349b5e25fSSatish Balay   PetscFunctionReturn(0);
26449b5e25fSSatish Balay }
26549b5e25fSSatish Balay 
2664a2ae208SSatish Balay #undef __FUNCT__
2674a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqSBAIJ"
26813f74950SBarry Smith PetscErrorCode MatRestoreRow_SeqSBAIJ(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
26949b5e25fSSatish Balay {
270dfbe8321SBarry Smith   PetscErrorCode ierr;
27149b5e25fSSatish Balay 
27249b5e25fSSatish Balay   PetscFunctionBegin;
27349b5e25fSSatish Balay   if (idx) {if (*idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}}
27449b5e25fSSatish Balay   if (v)   {if (*v)   {ierr = PetscFree(*v);CHKERRQ(ierr);}}
27549b5e25fSSatish Balay   PetscFunctionReturn(0);
27649b5e25fSSatish Balay }
27749b5e25fSSatish Balay 
2784a2ae208SSatish Balay #undef __FUNCT__
2794a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqSBAIJ"
280dfbe8321SBarry Smith PetscErrorCode MatTranspose_SeqSBAIJ(Mat A,Mat *B)
28149b5e25fSSatish Balay {
282dfbe8321SBarry Smith   PetscErrorCode ierr;
28349b5e25fSSatish Balay   PetscFunctionBegin;
284999d9058SBarry Smith   ierr = MatDuplicate(A,MAT_COPY_VALUES,B);CHKERRQ(ierr);
2858115998fSBarry Smith   PetscFunctionReturn(0);
28649b5e25fSSatish Balay }
28749b5e25fSSatish Balay 
2884a2ae208SSatish Balay #undef __FUNCT__
2894a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_ASCII"
2906849ba73SBarry Smith static PetscErrorCode MatView_SeqSBAIJ_ASCII(Mat A,PetscViewer viewer)
29149b5e25fSSatish Balay {
29249b5e25fSSatish Balay   Mat_SeqSBAIJ      *a = (Mat_SeqSBAIJ*)A->data;
293dfbe8321SBarry Smith   PetscErrorCode    ierr;
294521d7252SBarry Smith   PetscInt          i,j,bs = A->bs,k,l,bs2=a->bs2;
295fb9695e5SSatish Balay   char              *name;
296f3ef73ceSBarry Smith   PetscViewerFormat format;
29749b5e25fSSatish Balay 
29849b5e25fSSatish Balay   PetscFunctionBegin;
29980fe4e49SBarry Smith   ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
300b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
301456192e2SBarry Smith   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
30277431f27SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  block size is %D\n",bs);CHKERRQ(ierr);
303fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_MATLAB) {
30429bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Matlab format not supported");
305fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_COMMON) {
306b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
30749b5e25fSSatish Balay     for (i=0; i<a->mbs; i++) {
30849b5e25fSSatish Balay       for (j=0; j<bs; j++) {
30977431f27SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i*bs+j);CHKERRQ(ierr);
31049b5e25fSSatish Balay         for (k=a->i[i]; k<a->i[i+1]; k++) {
31149b5e25fSSatish Balay           for (l=0; l<bs; l++) {
31249b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX)
31349b5e25fSSatish Balay             if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
31477431f27SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%D, %g + %g i) ",bs*a->j[k]+l,
31549b5e25fSSatish Balay                                             PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
31649b5e25fSSatish Balay             } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
31777431f27SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%D, %g - %g i) ",bs*a->j[k]+l,
31849b5e25fSSatish Balay                                             PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
31949b5e25fSSatish Balay             } else if (PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
32077431f27SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%D, %g) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
32149b5e25fSSatish Balay             }
32249b5e25fSSatish Balay #else
32349b5e25fSSatish Balay             if (a->a[bs2*k + l*bs + j] != 0.0) {
32477431f27SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%D, %g) ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr);
32549b5e25fSSatish Balay             }
32649b5e25fSSatish Balay #endif
32749b5e25fSSatish Balay           }
32849b5e25fSSatish Balay         }
329b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
33049b5e25fSSatish Balay       }
33149b5e25fSSatish Balay     }
332b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
33349b5e25fSSatish Balay   } else {
334b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
33549b5e25fSSatish Balay     for (i=0; i<a->mbs; i++) {
33649b5e25fSSatish Balay       for (j=0; j<bs; j++) {
33777431f27SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i*bs+j);CHKERRQ(ierr);
33849b5e25fSSatish Balay         for (k=a->i[i]; k<a->i[i+1]; k++) {
33949b5e25fSSatish Balay           for (l=0; l<bs; l++) {
34049b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX)
34149b5e25fSSatish Balay             if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0) {
34277431f27SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%D, %g + %g i) ",bs*a->j[k]+l,
34349b5e25fSSatish Balay                                             PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
34449b5e25fSSatish Balay             } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0) {
34577431f27SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%D, %g - %g i) ",bs*a->j[k]+l,
34649b5e25fSSatish Balay                                             PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
34749b5e25fSSatish Balay             } else {
34877431f27SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%D, %g) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
34949b5e25fSSatish Balay             }
35049b5e25fSSatish Balay #else
35177431f27SBarry Smith             ierr = PetscViewerASCIIPrintf(viewer," %D %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr);
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);
35949b5e25fSSatish Balay   }
360b0a32e0cSBarry Smith   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
36149b5e25fSSatish Balay   PetscFunctionReturn(0);
36249b5e25fSSatish Balay }
36349b5e25fSSatish Balay 
3644a2ae208SSatish Balay #undef __FUNCT__
3654a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Draw_Zoom"
3666849ba73SBarry Smith static PetscErrorCode MatView_SeqSBAIJ_Draw_Zoom(PetscDraw draw,void *Aa)
36749b5e25fSSatish Balay {
36849b5e25fSSatish Balay   Mat            A = (Mat) Aa;
36949b5e25fSSatish Balay   Mat_SeqSBAIJ   *a=(Mat_SeqSBAIJ*)A->data;
3706849ba73SBarry Smith   PetscErrorCode ierr;
37113f74950SBarry Smith   PetscInt       row,i,j,k,l,mbs=a->mbs,color,bs=A->bs,bs2=a->bs2;
37213f74950SBarry Smith   PetscMPIInt    rank;
37349b5e25fSSatish Balay   PetscReal      xl,yl,xr,yr,x_l,x_r,y_l,y_r;
37449b5e25fSSatish Balay   MatScalar      *aa;
37549b5e25fSSatish Balay   MPI_Comm       comm;
376b0a32e0cSBarry Smith   PetscViewer    viewer;
37749b5e25fSSatish Balay 
37849b5e25fSSatish Balay   PetscFunctionBegin;
37949b5e25fSSatish Balay   /*
38049b5e25fSSatish Balay     This is nasty. If this is called from an originally parallel matrix
38149b5e25fSSatish Balay     then all processes call this,but only the first has the matrix so the
38249b5e25fSSatish Balay     rest should return immediately.
38349b5e25fSSatish Balay   */
38449b5e25fSSatish Balay   ierr = PetscObjectGetComm((PetscObject)draw,&comm);CHKERRQ(ierr);
38549b5e25fSSatish Balay   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
38649b5e25fSSatish Balay   if (rank) PetscFunctionReturn(0);
38749b5e25fSSatish Balay 
38849b5e25fSSatish Balay   ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr);
38949b5e25fSSatish Balay 
390b0a32e0cSBarry Smith   ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr);
391b0a32e0cSBarry Smith   PetscDrawString(draw, .3*(xl+xr), .3*(yl+yr), PETSC_DRAW_BLACK, "symmetric");
39249b5e25fSSatish Balay 
39349b5e25fSSatish Balay   /* loop over matrix elements drawing boxes */
394b0a32e0cSBarry Smith   color = PETSC_DRAW_BLUE;
39549b5e25fSSatish Balay   for (i=0,row=0; i<mbs; i++,row+=bs) {
39649b5e25fSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
397c464158bSHong Zhang       y_l = A->m - row - 1.0; y_r = y_l + 1.0;
39849b5e25fSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
39949b5e25fSSatish Balay       aa = a->a + j*bs2;
40049b5e25fSSatish Balay       for (k=0; k<bs; k++) {
40149b5e25fSSatish Balay         for (l=0; l<bs; l++) {
40249b5e25fSSatish Balay           if (PetscRealPart(*aa++) >=  0.) continue;
403b0a32e0cSBarry Smith           ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
40449b5e25fSSatish Balay         }
40549b5e25fSSatish Balay       }
40649b5e25fSSatish Balay     }
40749b5e25fSSatish Balay   }
408b0a32e0cSBarry Smith   color = PETSC_DRAW_CYAN;
40949b5e25fSSatish Balay   for (i=0,row=0; i<mbs; i++,row+=bs) {
41049b5e25fSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
411c464158bSHong Zhang       y_l = A->m - row - 1.0; y_r = y_l + 1.0;
41249b5e25fSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
41349b5e25fSSatish Balay       aa = a->a + j*bs2;
41449b5e25fSSatish Balay       for (k=0; k<bs; k++) {
41549b5e25fSSatish Balay         for (l=0; l<bs; l++) {
41649b5e25fSSatish Balay           if (PetscRealPart(*aa++) != 0.) continue;
417b0a32e0cSBarry Smith           ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
41849b5e25fSSatish Balay         }
41949b5e25fSSatish Balay       }
42049b5e25fSSatish Balay     }
42149b5e25fSSatish Balay   }
42249b5e25fSSatish Balay 
423b0a32e0cSBarry Smith   color = PETSC_DRAW_RED;
42449b5e25fSSatish Balay   for (i=0,row=0; i<mbs; i++,row+=bs) {
42549b5e25fSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
426c464158bSHong Zhang       y_l = A->m - row - 1.0; y_r = y_l + 1.0;
42749b5e25fSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
42849b5e25fSSatish Balay       aa = a->a + j*bs2;
42949b5e25fSSatish Balay       for (k=0; k<bs; k++) {
43049b5e25fSSatish Balay         for (l=0; l<bs; l++) {
43149b5e25fSSatish Balay           if (PetscRealPart(*aa++) <= 0.) continue;
432b0a32e0cSBarry Smith           ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
43349b5e25fSSatish Balay         }
43449b5e25fSSatish Balay       }
43549b5e25fSSatish Balay     }
43649b5e25fSSatish Balay   }
43749b5e25fSSatish Balay   PetscFunctionReturn(0);
43849b5e25fSSatish Balay }
43949b5e25fSSatish Balay 
4404a2ae208SSatish Balay #undef __FUNCT__
4414a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Draw"
4426849ba73SBarry Smith static PetscErrorCode MatView_SeqSBAIJ_Draw(Mat A,PetscViewer viewer)
44349b5e25fSSatish Balay {
444dfbe8321SBarry Smith   PetscErrorCode ierr;
44549b5e25fSSatish Balay   PetscReal      xl,yl,xr,yr,w,h;
446b0a32e0cSBarry Smith   PetscDraw      draw;
44749b5e25fSSatish Balay   PetscTruth     isnull;
44849b5e25fSSatish Balay 
44949b5e25fSSatish Balay   PetscFunctionBegin;
45049b5e25fSSatish Balay 
451b0a32e0cSBarry Smith   ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
452b0a32e0cSBarry Smith   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
45349b5e25fSSatish Balay 
45449b5e25fSSatish Balay   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
455c464158bSHong Zhang   xr  = A->m; yr = A->m; h = yr/10.0; w = xr/10.0;
45649b5e25fSSatish Balay   xr += w;    yr += h;  xl = -w;     yl = -h;
457b0a32e0cSBarry Smith   ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
458b0a32e0cSBarry Smith   ierr = PetscDrawZoom(draw,MatView_SeqSBAIJ_Draw_Zoom,A);CHKERRQ(ierr);
45949b5e25fSSatish Balay   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr);
46049b5e25fSSatish Balay   PetscFunctionReturn(0);
46149b5e25fSSatish Balay }
46249b5e25fSSatish Balay 
4634a2ae208SSatish Balay #undef __FUNCT__
4644a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ"
465dfbe8321SBarry Smith PetscErrorCode MatView_SeqSBAIJ(Mat A,PetscViewer viewer)
46649b5e25fSSatish Balay {
467dfbe8321SBarry Smith   PetscErrorCode ierr;
46832077d6dSBarry Smith   PetscTruth     iascii,isdraw;
46949b5e25fSSatish Balay 
47049b5e25fSSatish Balay   PetscFunctionBegin;
47132077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
472fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
47332077d6dSBarry Smith   if (iascii){
47449b5e25fSSatish Balay     ierr = MatView_SeqSBAIJ_ASCII(A,viewer);CHKERRQ(ierr);
47549b5e25fSSatish Balay   } else if (isdraw) {
47649b5e25fSSatish Balay     ierr = MatView_SeqSBAIJ_Draw(A,viewer);CHKERRQ(ierr);
47749b5e25fSSatish Balay   } else {
478a5e6ed63SBarry Smith     Mat B;
479ceb03754SKris Buschelman     ierr = MatConvert(A,MATSEQAIJ,MAT_INITIAL_MATRIX,&B);CHKERRQ(ierr);
480a5e6ed63SBarry Smith     ierr = MatView(B,viewer);CHKERRQ(ierr);
481a5e6ed63SBarry Smith     ierr = MatDestroy(B);CHKERRQ(ierr);
48249b5e25fSSatish Balay   }
48349b5e25fSSatish Balay   PetscFunctionReturn(0);
48449b5e25fSSatish Balay }
48549b5e25fSSatish Balay 
48649b5e25fSSatish Balay 
4874a2ae208SSatish Balay #undef __FUNCT__
4884a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqSBAIJ"
48913f74950SBarry Smith PetscErrorCode MatGetValues_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],PetscScalar v[])
49049b5e25fSSatish Balay {
491045c9aa0SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
49213f74950SBarry Smith   PetscInt     *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j;
49313f74950SBarry Smith   PetscInt     *ai = a->i,*ailen = a->ilen;
49413f74950SBarry Smith   PetscInt     brow,bcol,ridx,cidx,bs=A->bs,bs2=a->bs2;
49549b5e25fSSatish Balay   MatScalar    *ap,*aa = a->a,zero = 0.0;
49649b5e25fSSatish Balay 
49749b5e25fSSatish Balay   PetscFunctionBegin;
49849b5e25fSSatish Balay   for (k=0; k<m; k++) { /* loop over rows */
49949b5e25fSSatish Balay     row  = im[k]; brow = row/bs;
50077431f27SBarry Smith     if (row < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative row: %D",row);
50177431f27SBarry Smith     if (row >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->m-1);
50249b5e25fSSatish Balay     rp   = aj + ai[brow] ; ap = aa + bs2*ai[brow] ;
50349b5e25fSSatish Balay     nrow = ailen[brow];
50449b5e25fSSatish Balay     for (l=0; l<n; l++) { /* loop over columns */
50577431f27SBarry Smith       if (in[l] < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative column: %D",in[l]);
50677431f27SBarry Smith       if (in[l] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",in[l],A->n-1);
50749b5e25fSSatish Balay       col  = in[l] ;
50849b5e25fSSatish Balay       bcol = col/bs;
50949b5e25fSSatish Balay       cidx = col%bs;
51049b5e25fSSatish Balay       ridx = row%bs;
51149b5e25fSSatish Balay       high = nrow;
51249b5e25fSSatish Balay       low  = 0; /* assume unsorted */
51349b5e25fSSatish Balay       while (high-low > 5) {
51449b5e25fSSatish Balay         t = (low+high)/2;
51549b5e25fSSatish Balay         if (rp[t] > bcol) high = t;
51649b5e25fSSatish Balay         else             low  = t;
51749b5e25fSSatish Balay       }
51849b5e25fSSatish Balay       for (i=low; i<high; i++) {
51949b5e25fSSatish Balay         if (rp[i] > bcol) break;
52049b5e25fSSatish Balay         if (rp[i] == bcol) {
52149b5e25fSSatish Balay           *v++ = ap[bs2*i+bs*cidx+ridx];
52249b5e25fSSatish Balay           goto finished;
52349b5e25fSSatish Balay         }
52449b5e25fSSatish Balay       }
52549b5e25fSSatish Balay       *v++ = zero;
52649b5e25fSSatish Balay        finished:;
52749b5e25fSSatish Balay     }
52849b5e25fSSatish Balay   }
52949b5e25fSSatish Balay   PetscFunctionReturn(0);
53049b5e25fSSatish Balay }
53149b5e25fSSatish Balay 
53249b5e25fSSatish Balay 
5334a2ae208SSatish Balay #undef __FUNCT__
5344a2ae208SSatish Balay #define __FUNCT__ "MatSetValuesBlocked_SeqSBAIJ"
53513f74950SBarry Smith PetscErrorCode MatSetValuesBlocked_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is)
53649b5e25fSSatish Balay {
5370880e062SHong Zhang   Mat_SeqSBAIJ    *a = (Mat_SeqSBAIJ*)A->data;
5386849ba73SBarry Smith   PetscErrorCode  ierr;
539e2ee6c50SBarry Smith   PetscInt        *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,lastcol = -1;
54013f74950SBarry Smith   PetscInt        *imax=a->imax,*ai=a->i,*ailen=a->ilen;
54113f74950SBarry Smith   PetscInt        *aj=a->j,nonew=a->nonew,bs2=a->bs2,bs=A->bs,stepval;
5420880e062SHong Zhang   PetscTruth      roworiented=a->roworiented;
543f15d580aSBarry Smith   const MatScalar *value = v;
544f15d580aSBarry Smith   MatScalar       *ap,*aa = a->a,*bap;
5450880e062SHong Zhang 
54649b5e25fSSatish Balay   PetscFunctionBegin;
5470880e062SHong Zhang   if (roworiented) {
5480880e062SHong Zhang     stepval = (n-1)*bs;
5490880e062SHong Zhang   } else {
5500880e062SHong Zhang     stepval = (m-1)*bs;
5510880e062SHong Zhang   }
5520880e062SHong Zhang   for (k=0; k<m; k++) { /* loop over added rows */
5530880e062SHong Zhang     row  = im[k];
5540880e062SHong Zhang     if (row < 0) continue;
5552515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
55677431f27SBarry Smith     if (row >= a->mbs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,a->mbs-1);
5570880e062SHong Zhang #endif
5580880e062SHong Zhang     rp   = aj + ai[row];
5590880e062SHong Zhang     ap   = aa + bs2*ai[row];
5600880e062SHong Zhang     rmax = imax[row];
5610880e062SHong Zhang     nrow = ailen[row];
5620880e062SHong Zhang     low  = 0;
563818f2c47SBarry Smith     high = nrow;
5640880e062SHong Zhang     for (l=0; l<n; l++) { /* loop over added columns */
5650880e062SHong Zhang       if (in[l] < 0) continue;
5660880e062SHong Zhang       col = in[l];
5672515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
56877431f27SBarry Smith       if (col >= a->nbs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",col,a->nbs-1);
569b1823623SSatish Balay #endif
5700880e062SHong Zhang       if (col < row) continue; /* ignore lower triangular block */
5710880e062SHong Zhang       if (roworiented) {
5720880e062SHong Zhang         value = v + k*(stepval+bs)*bs + l*bs;
5730880e062SHong Zhang       } else {
5740880e062SHong Zhang         value = v + l*(stepval+bs)*bs + k*bs;
5750880e062SHong Zhang       }
576818f2c47SBarry Smith       if (col < lastcol) low = 0; else high = nrow;
577e2ee6c50SBarry Smith       lastcol = col;
5780880e062SHong Zhang       while (high-low > 7) {
5790880e062SHong Zhang         t = (low+high)/2;
5800880e062SHong Zhang         if (rp[t] > col) high = t;
5810880e062SHong Zhang         else             low  = t;
5820880e062SHong Zhang       }
5830880e062SHong Zhang       for (i=low; i<high; i++) {
5840880e062SHong Zhang         if (rp[i] > col) break;
5850880e062SHong Zhang         if (rp[i] == col) {
5860880e062SHong Zhang           bap  = ap +  bs2*i;
5870880e062SHong Zhang           if (roworiented) {
5880880e062SHong Zhang             if (is == ADD_VALUES) {
5890880e062SHong Zhang               for (ii=0; ii<bs; ii++,value+=stepval) {
5900880e062SHong Zhang                 for (jj=ii; jj<bs2; jj+=bs) {
5910880e062SHong Zhang                   bap[jj] += *value++;
5920880e062SHong Zhang                 }
5930880e062SHong Zhang               }
5940880e062SHong Zhang             } else {
5950880e062SHong Zhang               for (ii=0; ii<bs; ii++,value+=stepval) {
5960880e062SHong Zhang                 for (jj=ii; jj<bs2; jj+=bs) {
5970880e062SHong Zhang                   bap[jj] = *value++;
5980880e062SHong Zhang                 }
5990880e062SHong Zhang                }
6000880e062SHong Zhang             }
6010880e062SHong Zhang           } else {
6020880e062SHong Zhang             if (is == ADD_VALUES) {
6030880e062SHong Zhang               for (ii=0; ii<bs; ii++,value+=stepval) {
6040880e062SHong Zhang                 for (jj=0; jj<bs; jj++) {
6050880e062SHong Zhang                   *bap++ += *value++;
6060880e062SHong Zhang                 }
6070880e062SHong Zhang               }
6080880e062SHong Zhang             } else {
6090880e062SHong Zhang               for (ii=0; ii<bs; ii++,value+=stepval) {
6100880e062SHong Zhang                 for (jj=0; jj<bs; jj++) {
6110880e062SHong Zhang                   *bap++  = *value++;
6120880e062SHong Zhang                 }
6130880e062SHong Zhang               }
6140880e062SHong Zhang             }
6150880e062SHong Zhang           }
6160880e062SHong Zhang           goto noinsert2;
6170880e062SHong Zhang         }
6180880e062SHong Zhang       }
6190880e062SHong Zhang       if (nonew == 1) goto noinsert2;
620085a36d4SBarry Smith       if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col);
621*ed1caa07SMatthew Knepley       MatSeqXAIJReallocateAIJ(a,bs2,nrow,row,col,rmax,aa,ai,aj,a->mbs,rp,ap,imax,nonew);
6220880e062SHong Zhang       N = nrow++ - 1;
6230880e062SHong Zhang       /* shift up all the later entries in this row */
6240880e062SHong Zhang       for (ii=N; ii>=i; ii--) {
6250880e062SHong Zhang         rp[ii+1] = rp[ii];
6260880e062SHong Zhang         ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr);
6270880e062SHong Zhang       }
6280880e062SHong Zhang       if (N >= i) {
6290880e062SHong Zhang         ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr);
6300880e062SHong Zhang       }
6310880e062SHong Zhang       rp[i] = col;
6320880e062SHong Zhang       bap   = ap +  bs2*i;
6330880e062SHong Zhang       if (roworiented) {
6340880e062SHong Zhang         for (ii=0; ii<bs; ii++,value+=stepval) {
6350880e062SHong Zhang           for (jj=ii; jj<bs2; jj+=bs) {
6360880e062SHong Zhang             bap[jj] = *value++;
6370880e062SHong Zhang           }
6380880e062SHong Zhang         }
6390880e062SHong Zhang       } else {
6400880e062SHong Zhang         for (ii=0; ii<bs; ii++,value+=stepval) {
6410880e062SHong Zhang           for (jj=0; jj<bs; jj++) {
6420880e062SHong Zhang             *bap++  = *value++;
6430880e062SHong Zhang           }
6440880e062SHong Zhang         }
6450880e062SHong Zhang        }
6460880e062SHong Zhang     noinsert2:;
6470880e062SHong Zhang       low = i;
6480880e062SHong Zhang     }
6490880e062SHong Zhang     ailen[row] = nrow;
6500880e062SHong Zhang   }
6510880e062SHong Zhang    PetscFunctionReturn(0);
65249b5e25fSSatish Balay }
65349b5e25fSSatish Balay 
6544a2ae208SSatish Balay #undef __FUNCT__
6554a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_SeqSBAIJ"
656dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_SeqSBAIJ(Mat A,MatAssemblyType mode)
65749b5e25fSSatish Balay {
65849b5e25fSSatish Balay   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
6596849ba73SBarry Smith   PetscErrorCode ierr;
66013f74950SBarry Smith   PetscInt       fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax;
66113f74950SBarry Smith   PetscInt       m = A->m,*ip,N,*ailen = a->ilen;
66213f74950SBarry Smith   PetscInt       mbs = a->mbs,bs2 = a->bs2,rmax = 0;
66349b5e25fSSatish Balay   MatScalar      *aa = a->a,*ap;
66449b5e25fSSatish Balay 
66549b5e25fSSatish Balay   PetscFunctionBegin;
66649b5e25fSSatish Balay   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0);
66749b5e25fSSatish Balay 
66849b5e25fSSatish Balay   if (m) rmax = ailen[0];
66949b5e25fSSatish Balay   for (i=1; i<mbs; i++) {
67049b5e25fSSatish Balay     /* move each row back by the amount of empty slots (fshift) before it*/
67149b5e25fSSatish Balay     fshift += imax[i-1] - ailen[i-1];
67249b5e25fSSatish Balay      rmax   = PetscMax(rmax,ailen[i]);
67349b5e25fSSatish Balay      if (fshift) {
67449b5e25fSSatish Balay        ip = aj + ai[i]; ap = aa + bs2*ai[i];
67549b5e25fSSatish Balay        N = ailen[i];
67649b5e25fSSatish Balay        for (j=0; j<N; j++) {
67749b5e25fSSatish Balay          ip[j-fshift] = ip[j];
67849b5e25fSSatish Balay          ierr = PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
67949b5e25fSSatish Balay        }
68049b5e25fSSatish Balay      }
68149b5e25fSSatish Balay      ai[i] = ai[i-1] + ailen[i-1];
68249b5e25fSSatish Balay   }
68349b5e25fSSatish Balay   if (mbs) {
68449b5e25fSSatish Balay     fshift += imax[mbs-1] - ailen[mbs-1];
68549b5e25fSSatish Balay      ai[mbs] = ai[mbs-1] + ailen[mbs-1];
68649b5e25fSSatish Balay   }
68749b5e25fSSatish Balay   /* reset ilen and imax for each row */
68849b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
68949b5e25fSSatish Balay     ailen[i] = imax[i] = ai[i+1] - ai[i];
69049b5e25fSSatish Balay   }
6916c6c5352SBarry Smith   a->nz = ai[mbs];
69249b5e25fSSatish Balay 
693b424e231SHong Zhang   /* diagonals may have moved, reset it */
694b424e231SHong Zhang   if (a->diag) {
69513f74950SBarry Smith     ierr = PetscMemcpy(a->diag,ai,(mbs+1)*sizeof(PetscInt));CHKERRQ(ierr);
69649b5e25fSSatish Balay   }
69763ba0a88SBarry Smith   ierr = PetscLogInfo((A,"MatAssemblyEnd_SeqSBAIJ:Matrix size: %D X %D, block size %D; storage space: %D unneeded, %D used\n",
69863ba0a88SBarry Smith                        m,A->m,A->bs,fshift*bs2,a->nz*bs2));CHKERRQ(ierr);
69963ba0a88SBarry Smith   ierr = PetscLogInfo((A,"MatAssemblyEnd_SeqSBAIJ:Number of mallocs during MatSetValues is %D\n",a->reallocs));CHKERRQ(ierr);
70063ba0a88SBarry Smith   ierr = PetscLogInfo((A,"MatAssemblyEnd_SeqSBAIJ:Most nonzeros blocks in any row is %D\n",rmax));CHKERRQ(ierr);
70149b5e25fSSatish Balay   a->reallocs          = 0;
70249b5e25fSSatish Balay   A->info.nz_unneeded  = (PetscReal)fshift*bs2;
70349b5e25fSSatish Balay   PetscFunctionReturn(0);
70449b5e25fSSatish Balay }
70549b5e25fSSatish Balay 
70649b5e25fSSatish Balay /*
70749b5e25fSSatish Balay    This function returns an array of flags which indicate the locations of contiguous
70849b5e25fSSatish Balay    blocks that should be zeroed. for eg: if bs = 3  and is = [0,1,2,3,5,6,7,8,9]
70949b5e25fSSatish Balay    then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)]
71049b5e25fSSatish Balay    Assume: sizes should be long enough to hold all the values.
71149b5e25fSSatish Balay */
7124a2ae208SSatish Balay #undef __FUNCT__
7134a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqSBAIJ_Check_Blocks"
71413f74950SBarry Smith PetscErrorCode MatZeroRows_SeqSBAIJ_Check_Blocks(PetscInt idx[],PetscInt n,PetscInt bs,PetscInt sizes[], PetscInt *bs_max)
71549b5e25fSSatish Balay {
71613f74950SBarry Smith   PetscInt   i,j,k,row;
71749b5e25fSSatish Balay   PetscTruth flg;
71849b5e25fSSatish Balay 
71949b5e25fSSatish Balay   PetscFunctionBegin;
72049b5e25fSSatish Balay    for (i=0,j=0; i<n; j++) {
72149b5e25fSSatish Balay      row = idx[i];
72249b5e25fSSatish Balay      if (row%bs!=0) { /* Not the begining of a block */
72349b5e25fSSatish Balay        sizes[j] = 1;
72449b5e25fSSatish Balay        i++;
72549b5e25fSSatish Balay      } else if (i+bs > n) { /* Beginning of a block, but complete block doesn't exist (at idx end) */
72649b5e25fSSatish Balay        sizes[j] = 1;         /* Also makes sure atleast 'bs' values exist for next else */
72749b5e25fSSatish Balay        i++;
72849b5e25fSSatish Balay      } else { /* Begining of the block, so check if the complete block exists */
72949b5e25fSSatish Balay        flg = PETSC_TRUE;
73049b5e25fSSatish Balay        for (k=1; k<bs; k++) {
73149b5e25fSSatish Balay          if (row+k != idx[i+k]) { /* break in the block */
73249b5e25fSSatish Balay            flg = PETSC_FALSE;
73349b5e25fSSatish Balay            break;
73449b5e25fSSatish Balay          }
73549b5e25fSSatish Balay        }
736abc0a331SBarry Smith        if (flg) { /* No break in the bs */
73749b5e25fSSatish Balay          sizes[j] = bs;
73849b5e25fSSatish Balay          i+= bs;
73949b5e25fSSatish Balay        } else {
74049b5e25fSSatish Balay          sizes[j] = 1;
74149b5e25fSSatish Balay          i++;
74249b5e25fSSatish Balay        }
74349b5e25fSSatish Balay      }
74449b5e25fSSatish Balay    }
74549b5e25fSSatish Balay    *bs_max = j;
74649b5e25fSSatish Balay    PetscFunctionReturn(0);
74749b5e25fSSatish Balay }
74849b5e25fSSatish Balay 
74949b5e25fSSatish Balay 
75049b5e25fSSatish Balay /* Only add/insert a(i,j) with i<=j (blocks).
75149b5e25fSSatish Balay    Any a(i,j) with i>j input by user is ingored.
75249b5e25fSSatish Balay */
75349b5e25fSSatish Balay 
7544a2ae208SSatish Balay #undef __FUNCT__
7554a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqSBAIJ"
75613f74950SBarry Smith PetscErrorCode MatSetValues_SeqSBAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is)
75749b5e25fSSatish Balay {
75849b5e25fSSatish Balay   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
7596849ba73SBarry Smith   PetscErrorCode ierr;
760e2ee6c50SBarry Smith   PetscInt       *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,lastcol = -1;
76113f74950SBarry Smith   PetscInt       *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented;
76213f74950SBarry Smith   PetscInt       *aj=a->j,nonew=a->nonew,bs=A->bs,brow,bcol;
76313f74950SBarry Smith   PetscInt       ridx,cidx,bs2=a->bs2;
76449b5e25fSSatish Balay   MatScalar      *ap,value,*aa=a->a,*bap;
76549b5e25fSSatish Balay 
76649b5e25fSSatish Balay   PetscFunctionBegin;
76749b5e25fSSatish Balay   for (k=0; k<m; k++) { /* loop over added rows */
76849b5e25fSSatish Balay     row  = im[k];       /* row number */
76949b5e25fSSatish Balay     brow = row/bs;      /* block row number */
77049b5e25fSSatish Balay     if (row < 0) continue;
7712515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
77277431f27SBarry Smith     if (row >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->m-1);
77349b5e25fSSatish Balay #endif
77449b5e25fSSatish Balay     rp   = aj + ai[brow]; /*ptr to beginning of column value of the row block*/
77549b5e25fSSatish Balay     ap   = aa + bs2*ai[brow]; /*ptr to beginning of element value of the row block*/
77649b5e25fSSatish Balay     rmax = imax[brow];  /* maximum space allocated for this row */
77749b5e25fSSatish Balay     nrow = ailen[brow]; /* actual length of this row */
77849b5e25fSSatish Balay     low  = 0;
77949b5e25fSSatish Balay 
78049b5e25fSSatish Balay     for (l=0; l<n; l++) { /* loop over added columns */
78149b5e25fSSatish Balay       if (in[l] < 0) continue;
7822515c552SBarry Smith #if defined(PETSC_USE_DEBUG)
78377431f27SBarry Smith       if (in[l] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",in[l],A->m-1);
78449b5e25fSSatish Balay #endif
78549b5e25fSSatish Balay       col = in[l];
78649b5e25fSSatish Balay       bcol = col/bs;              /* block col number */
78749b5e25fSSatish Balay 
788941593c8SHong Zhang       if (brow > bcol) {
789941593c8SHong Zhang         if (a->ignore_ltriangular){
790941593c8SHong Zhang           continue; /* ignore lower triangular values */
791941593c8SHong Zhang         } else {
792941593c8SHong 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)");
793941593c8SHong Zhang         }
794941593c8SHong Zhang       }
795f4989cb3SHong Zhang 
79649b5e25fSSatish Balay       ridx = row % bs; cidx = col % bs; /*row and col index inside the block */
7978549e402SHong Zhang       if ((brow==bcol && ridx<=cidx) || (brow<bcol)){
79849b5e25fSSatish Balay         /* element value a(k,l) */
79949b5e25fSSatish Balay         if (roworiented) {
80049b5e25fSSatish Balay           value = v[l + k*n];
80149b5e25fSSatish Balay         } else {
80249b5e25fSSatish Balay           value = v[k + l*m];
80349b5e25fSSatish Balay         }
80449b5e25fSSatish Balay 
80549b5e25fSSatish Balay         /* move pointer bap to a(k,l) quickly and add/insert value */
806e2ee6c50SBarry Smith         if (col < lastcol) low = 0; high = nrow;
807e2ee6c50SBarry Smith         lastcol = col;
80849b5e25fSSatish Balay         while (high-low > 7) {
80949b5e25fSSatish Balay           t = (low+high)/2;
81049b5e25fSSatish Balay           if (rp[t] > bcol) high = t;
81149b5e25fSSatish Balay           else              low  = t;
81249b5e25fSSatish Balay         }
81349b5e25fSSatish Balay         for (i=low; i<high; i++) {
81449b5e25fSSatish Balay           if (rp[i] > bcol) break;
81549b5e25fSSatish Balay           if (rp[i] == bcol) {
81649b5e25fSSatish Balay             bap  = ap +  bs2*i + bs*cidx + ridx;
81749b5e25fSSatish Balay             if (is == ADD_VALUES) *bap += value;
81849b5e25fSSatish Balay             else                  *bap  = value;
8198549e402SHong Zhang             /* for diag block, add/insert its symmetric element a(cidx,ridx) */
8208549e402SHong Zhang             if (brow == bcol && ridx < cidx){
8218549e402SHong Zhang               bap  = ap +  bs2*i + bs*ridx + cidx;
8228549e402SHong Zhang               if (is == ADD_VALUES) *bap += value;
8238549e402SHong Zhang               else                  *bap  = value;
8248549e402SHong Zhang             }
82549b5e25fSSatish Balay             goto noinsert1;
82649b5e25fSSatish Balay           }
82749b5e25fSSatish Balay         }
82849b5e25fSSatish Balay 
82949b5e25fSSatish Balay         if (nonew == 1) goto noinsert1;
830085a36d4SBarry Smith         if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero (%D, %D) in the matrix", row, col);
831*ed1caa07SMatthew Knepley         MatSeqXAIJReallocateAIJ(a,bs2,nrow,brow,bcol,rmax,aa,ai,aj,a->mbs,rp,ap,imax,nonew);
83249b5e25fSSatish Balay 
83349b5e25fSSatish Balay         N = nrow++ - 1;
83449b5e25fSSatish Balay         /* shift up all the later entries in this row */
83549b5e25fSSatish Balay         for (ii=N; ii>=i; ii--) {
83649b5e25fSSatish Balay           rp[ii+1] = rp[ii];
83749b5e25fSSatish Balay           ierr     = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr);
83849b5e25fSSatish Balay         }
83949b5e25fSSatish Balay         if (N>=i) {
84049b5e25fSSatish Balay           ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr);
84149b5e25fSSatish Balay         }
84249b5e25fSSatish Balay         rp[i]                      = bcol;
84349b5e25fSSatish Balay         ap[bs2*i + bs*cidx + ridx] = value;
84449b5e25fSSatish Balay       noinsert1:;
84549b5e25fSSatish Balay         low = i;
8468549e402SHong Zhang       }
84749b5e25fSSatish Balay     }   /* end of loop over added columns */
84849b5e25fSSatish Balay     ailen[brow] = nrow;
84949b5e25fSSatish Balay   }   /* end of loop over added rows */
85049b5e25fSSatish Balay   PetscFunctionReturn(0);
85149b5e25fSSatish Balay }
85249b5e25fSSatish Balay 
8536849ba73SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ(Mat,IS,MatFactorInfo*,Mat*);
8546849ba73SBarry Smith EXTERN PetscErrorCode MatCholeskyFactor_SeqSBAIJ(Mat,IS,MatFactorInfo*);
85513f74950SBarry Smith EXTERN PetscErrorCode MatIncreaseOverlap_SeqSBAIJ(Mat,PetscInt,IS[],PetscInt);
85613f74950SBarry Smith EXTERN PetscErrorCode MatGetSubMatrix_SeqSBAIJ(Mat,IS,IS,PetscInt,MatReuse,Mat*);
85713f74950SBarry Smith EXTERN PetscErrorCode MatGetSubMatrices_SeqSBAIJ(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat*[]);
8586849ba73SBarry Smith EXTERN PetscErrorCode MatScale_SeqSBAIJ(const PetscScalar*,Mat);
8596849ba73SBarry Smith EXTERN PetscErrorCode MatNorm_SeqSBAIJ(Mat,NormType,PetscReal *);
8606849ba73SBarry Smith EXTERN PetscErrorCode MatEqual_SeqSBAIJ(Mat,Mat,PetscTruth*);
8616849ba73SBarry Smith EXTERN PetscErrorCode MatGetDiagonal_SeqSBAIJ(Mat,Vec);
8626849ba73SBarry Smith EXTERN PetscErrorCode MatDiagonalScale_SeqSBAIJ(Mat,Vec,Vec);
8636849ba73SBarry Smith EXTERN PetscErrorCode MatGetInfo_SeqSBAIJ(Mat,MatInfoType,MatInfo *);
8646849ba73SBarry Smith EXTERN PetscErrorCode MatZeroEntries_SeqSBAIJ(Mat);
8656849ba73SBarry Smith EXTERN PetscErrorCode MatGetRowMax_SeqSBAIJ(Mat,Vec);
86649b5e25fSSatish Balay 
8676849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_N(Mat,Vec,Vec);
8686849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_1(Mat,Vec,Vec);
8696849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_2(Mat,Vec,Vec);
8706849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_3(Mat,Vec,Vec);
8716849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_4(Mat,Vec,Vec);
8726849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_5(Mat,Vec,Vec);
8736849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_6(Mat,Vec,Vec);
8746849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_7(Mat,Vec,Vec);
87549b5e25fSSatish Balay 
8766849ba73SBarry Smith EXTERN PetscErrorCode MatSolves_SeqSBAIJ_1(Mat,Vecs,Vecs);
877d59c15a7SBarry Smith 
8786849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_1_NaturalOrdering(Mat,Vec,Vec);
8796849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_2_NaturalOrdering(Mat,Vec,Vec);
8806849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_3_NaturalOrdering(Mat,Vec,Vec);
8816849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_4_NaturalOrdering(Mat,Vec,Vec);
8826849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_5_NaturalOrdering(Mat,Vec,Vec);
8836849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_6_NaturalOrdering(Mat,Vec,Vec);
8846849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_7_NaturalOrdering(Mat,Vec,Vec);
8856849ba73SBarry Smith EXTERN PetscErrorCode MatSolve_SeqSBAIJ_N_NaturalOrdering(Mat,Vec,Vec);
88607f98182SSatish Balay 
887af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat,MatFactorInfo*,Mat*);
888af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1(Mat,MatFactorInfo*,Mat*);
889af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_2(Mat,MatFactorInfo*,Mat*);
890af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_3(Mat,MatFactorInfo*,Mat*);
891af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_4(Mat,MatFactorInfo*,Mat*);
892af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_5(Mat,MatFactorInfo*,Mat*);
893af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_6(Mat,MatFactorInfo*,Mat*);
894af281ebdSHong Zhang EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_7(Mat,MatFactorInfo*,Mat*);
89513f74950SBarry Smith EXTERN PetscErrorCode MatGetInertia_SeqSBAIJ(Mat,PetscInt*,PetscInt*,PetscInt*);
89649b5e25fSSatish Balay 
8976849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_1(Mat,Vec,Vec);
8986849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_2(Mat,Vec,Vec);
8996849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_3(Mat,Vec,Vec);
9006849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_4(Mat,Vec,Vec);
9016849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_5(Mat,Vec,Vec);
9026849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_6(Mat,Vec,Vec);
9036849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_7(Mat,Vec,Vec);
9046849ba73SBarry Smith EXTERN PetscErrorCode MatMult_SeqSBAIJ_N(Mat,Vec,Vec);
90549b5e25fSSatish Balay 
9066849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_1(Mat,Vec,Vec,Vec);
9076849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_2(Mat,Vec,Vec,Vec);
9086849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_3(Mat,Vec,Vec,Vec);
9096849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_4(Mat,Vec,Vec,Vec);
9106849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_5(Mat,Vec,Vec,Vec);
9116849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_6(Mat,Vec,Vec,Vec);
9126849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_7(Mat,Vec,Vec,Vec);
9136849ba73SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqSBAIJ_N(Mat,Vec,Vec,Vec);
91449b5e25fSSatish Balay 
9154d101231SSatish Balay #ifdef HAVE_MatICCFactor
9166849ba73SBarry Smith /* modified from MatILUFactor_SeqSBAIJ, needs further work!  */
9174a2ae208SSatish Balay #undef __FUNCT__
9184d101231SSatish Balay #define __FUNCT__ "MatICCFactor_SeqSBAIJ"
919dfbe8321SBarry Smith PetscErrorCode MatICCFactor_SeqSBAIJ(Mat inA,IS row,MatFactorInfo *info)
92049b5e25fSSatish Balay {
9214ccecd49SHong Zhang   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)inA->data;
92249b5e25fSSatish Balay   Mat            outA;
923dfbe8321SBarry Smith   PetscErrorCode ierr;
92449b5e25fSSatish Balay   PetscTruth     row_identity,col_identity;
92549b5e25fSSatish Balay 
92649b5e25fSSatish Balay   PetscFunctionBegin;
92749b5e25fSSatish Balay   outA          = inA;
9281260e1f8SHong Zhang   inA->factor   = FACTOR_CHOLESKY;
92949b5e25fSSatish Balay 
93049b5e25fSSatish Balay   if (!a->diag) {
9311a3463dfSHong Zhang     ierr = MatMarkDiagonal_SeqSBAIJ(inA);CHKERRQ(ierr);
93249b5e25fSSatish Balay   }
93349b5e25fSSatish Balay   /*
93449b5e25fSSatish Balay     Blocksize 2, 3, 4, 5, 6 and 7 have a special faster factorization/solver
93549b5e25fSSatish Balay     for ILU(0) factorization with natural ordering
93649b5e25fSSatish Balay   */
93749b5e25fSSatish Balay   switch (a->bs) {
93849b5e25fSSatish Balay   case 1:
9390fe9e5beSBarry Smith     inA->ops->solve            = MatSolve_SeqSBAIJ_1_NaturalOrdering;
94007f98182SSatish Balay     inA->ops->solvetranspose   = MatSolve_SeqSBAIJ_1_NaturalOrdering;
941d59c15a7SBarry Smith     inA->ops->solves           = MatSolves_SeqSBAIJ_1;
94263ba0a88SBarry Smith     ierr = PetscLoginfo((inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering solvetrans BS=1\n"));CHKERRQ(ierr);
94349b5e25fSSatish Balay   case 2:
9441a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering;
9451a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_2_NaturalOrdering;
9461a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_2_NaturalOrdering;
94763ba0a88SBarry Smith     ierr = PetscLogInfo((inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=2\n"));CHKERRQ(ierr);
94849b5e25fSSatish Balay     break;
94949b5e25fSSatish Balay   case 3:
9501a3463dfSHong Zhang      inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering;
9511a3463dfSHong Zhang      inA->ops->solve           = MatSolve_SeqSBAIJ_3_NaturalOrdering;
9521a3463dfSHong Zhang      inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_3_NaturalOrdering;
95363ba0a88SBarry Smith      ierr = PetscLogInfo((inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=3\n"));CHKERRQ(ierr);
95449b5e25fSSatish Balay      break;
95549b5e25fSSatish Balay   case 4:
9561a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering;
9571a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_4_NaturalOrdering;
9581a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_4_NaturalOrdering;
95963ba0a88SBarry Smith     ierr = PetscLogInfo((inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=4\n"));CHKERRQ(ierr);
96049b5e25fSSatish Balay     break;
96149b5e25fSSatish Balay   case 5:
9621a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering;
9631a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_5_NaturalOrdering;
9641a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_5_NaturalOrdering;
96563ba0a88SBarry Smith     ierr = PetscLogInfo((inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=5\n"));CHKERRQ(ierr);
96649b5e25fSSatish Balay     break;
96749b5e25fSSatish Balay   case 6:
9681a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering;
9691a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_6_NaturalOrdering;
9701a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_6_NaturalOrdering;
97163ba0a88SBarry Smith     ierr = PetscLogInfo((inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=6\n"));CHKERRQ(ierr);
97249b5e25fSSatish Balay     break;
97349b5e25fSSatish Balay   case 7:
9741a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering;
9751a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_7_NaturalOrdering;
9761a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_7_NaturalOrdering;
97763ba0a88SBarry Smith     ierr = PetscLogInfo((inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=7\n"));CHKERRQ(ierr);
97849b5e25fSSatish Balay     break;
97949b5e25fSSatish Balay   default:
98049b5e25fSSatish Balay     a->row        = row;
9811a3463dfSHong Zhang     a->icol       = col;
98249b5e25fSSatish Balay     ierr          = PetscObjectReference((PetscObject)row);CHKERRQ(ierr);
98349b5e25fSSatish Balay     ierr          = PetscObjectReference((PetscObject)col);CHKERRQ(ierr);
98449b5e25fSSatish Balay 
98549b5e25fSSatish Balay     /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */
98649b5e25fSSatish Balay     ierr = ISInvertPermutation(col,PETSC_DECIDE, &(a->icol));CHKERRQ(ierr);
98752e6d16bSBarry Smith     ierr = PetscLogObjectParent(inA,a->icol);CHKERRQ(ierr);
98849b5e25fSSatish Balay 
98949b5e25fSSatish Balay     if (!a->solve_work) {
99087828ca2SBarry Smith       ierr = PetscMalloc((A->m+a->bs)*sizeof(PetscScalar),&a->solve_work);CHKERRQ(ierr);
99152e6d16bSBarry Smith       ierr = PetscLogObjectMemory(inA,(A->m+a->bs)*sizeof(PetscScalar));CHKERRQ(ierr);
99249b5e25fSSatish Balay     }
99349b5e25fSSatish Balay   }
99449b5e25fSSatish Balay 
995af281ebdSHong Zhang   ierr = MatCholeskyFactorNumeric(inA,info,&outA);CHKERRQ(ierr);
99649b5e25fSSatish Balay   PetscFunctionReturn(0);
99749b5e25fSSatish Balay }
998950f1e5bSHong Zhang #endif
999950f1e5bSHong Zhang 
10004a2ae208SSatish Balay #undef __FUNCT__
10014a2ae208SSatish Balay #define __FUNCT__ "MatPrintHelp_SeqSBAIJ"
1002dfbe8321SBarry Smith PetscErrorCode MatPrintHelp_SeqSBAIJ(Mat A)
100349b5e25fSSatish Balay {
100449b5e25fSSatish Balay   static PetscTruth called = PETSC_FALSE;
100549b5e25fSSatish Balay   MPI_Comm          comm = A->comm;
1006dfbe8321SBarry Smith   PetscErrorCode    ierr;
100749b5e25fSSatish Balay 
100849b5e25fSSatish Balay   PetscFunctionBegin;
100949b5e25fSSatish Balay   if (called) {PetscFunctionReturn(0);} else called = PETSC_TRUE;
101049b5e25fSSatish Balay   ierr = (*PetscHelpPrintf)(comm," Options for MATSEQSBAIJ and MATMPISBAIJ matrix formats (the defaults):\n");CHKERRQ(ierr);
101149b5e25fSSatish Balay   ierr = (*PetscHelpPrintf)(comm,"  -mat_block_size <block_size>\n");CHKERRQ(ierr);
1012941593c8SHong Zhang   ierr = (*PetscHelpPrintf)(comm,"  -mat_ignore_ltriangular: Ignore lower triangular values set by user\n");CHKERRQ(ierr);
101349b5e25fSSatish Balay   PetscFunctionReturn(0);
101449b5e25fSSatish Balay }
101549b5e25fSSatish Balay 
101649b5e25fSSatish Balay EXTERN_C_BEGIN
10174a2ae208SSatish Balay #undef __FUNCT__
10184a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices_SeqSBAIJ"
1019be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetColumnIndices_SeqSBAIJ(Mat mat,PetscInt *indices)
102049b5e25fSSatish Balay {
1021045c9aa0SHong Zhang   Mat_SeqSBAIJ *baij = (Mat_SeqSBAIJ *)mat->data;
102213f74950SBarry Smith   PetscInt     i,nz,n;
102349b5e25fSSatish Balay 
102449b5e25fSSatish Balay   PetscFunctionBegin;
10256c6c5352SBarry Smith   nz = baij->maxnz;
1026c464158bSHong Zhang   n  = mat->n;
102749b5e25fSSatish Balay   for (i=0; i<nz; i++) {
102849b5e25fSSatish Balay     baij->j[i] = indices[i];
102949b5e25fSSatish Balay   }
10306c6c5352SBarry Smith    baij->nz = nz;
103149b5e25fSSatish Balay    for (i=0; i<n; i++) {
103249b5e25fSSatish Balay      baij->ilen[i] = baij->imax[i];
103349b5e25fSSatish Balay    }
103449b5e25fSSatish Balay    PetscFunctionReturn(0);
103549b5e25fSSatish Balay }
103649b5e25fSSatish Balay EXTERN_C_END
103749b5e25fSSatish Balay 
10384a2ae208SSatish Balay #undef __FUNCT__
10394a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices"
104049b5e25fSSatish Balay /*@
104119585528SSatish Balay   MatSeqSBAIJSetColumnIndices - Set the column indices for all the rows
104249b5e25fSSatish Balay   in the matrix.
104349b5e25fSSatish Balay 
104449b5e25fSSatish Balay   Input Parameters:
104519585528SSatish Balay   +  mat     - the SeqSBAIJ matrix
104649b5e25fSSatish Balay   -  indices - the column indices
104749b5e25fSSatish Balay 
104849b5e25fSSatish Balay   Level: advanced
104949b5e25fSSatish Balay 
105049b5e25fSSatish Balay   Notes:
105149b5e25fSSatish Balay   This can be called if you have precomputed the nonzero structure of the
105249b5e25fSSatish Balay   matrix and want to provide it to the matrix object to improve the performance
105349b5e25fSSatish Balay   of the MatSetValues() operation.
105449b5e25fSSatish Balay 
105549b5e25fSSatish Balay   You MUST have set the correct numbers of nonzeros per row in the call to
1056d1be2dadSMatthew Knepley   MatCreateSeqSBAIJ(), and the columns indices MUST be sorted.
105749b5e25fSSatish Balay 
1058ab9f2c04SSatish Balay   MUST be called before any calls to MatSetValues()
105949b5e25fSSatish Balay 
1060ab9f2c04SSatish Balay   .seealso: MatCreateSeqSBAIJ
106149b5e25fSSatish Balay @*/
1062be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetColumnIndices(Mat mat,PetscInt *indices)
106349b5e25fSSatish Balay {
106413f74950SBarry Smith   PetscErrorCode ierr,(*f)(Mat,PetscInt *);
106549b5e25fSSatish Balay 
106649b5e25fSSatish Balay   PetscFunctionBegin;
10674482741eSBarry Smith   PetscValidHeaderSpecific(mat,MAT_COOKIE,1);
10684482741eSBarry Smith   PetscValidPointer(indices,2);
1069c134de8dSSatish Balay   ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqSBAIJSetColumnIndices_C",(void (**)(void))&f);CHKERRQ(ierr);
107049b5e25fSSatish Balay   if (f) {
107149b5e25fSSatish Balay     ierr = (*f)(mat,indices);CHKERRQ(ierr);
107249b5e25fSSatish Balay   } else {
1073e005ede5SBarry Smith     SETERRQ(PETSC_ERR_SUP,"Wrong type of matrix to set column indices");
107449b5e25fSSatish Balay   }
107549b5e25fSSatish Balay   PetscFunctionReturn(0);
107649b5e25fSSatish Balay }
107749b5e25fSSatish Balay 
10784a2ae208SSatish Balay #undef __FUNCT__
10794a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqSBAIJ"
1080dfbe8321SBarry Smith PetscErrorCode MatSetUpPreallocation_SeqSBAIJ(Mat A)
1081273d9f13SBarry Smith {
1082dfbe8321SBarry Smith   PetscErrorCode ierr;
1083273d9f13SBarry Smith 
1084273d9f13SBarry Smith   PetscFunctionBegin;
1085ab93d7beSBarry Smith   ierr =  MatSeqSBAIJSetPreallocation_SeqSBAIJ(A,1,PETSC_DEFAULT,0);CHKERRQ(ierr);
1086273d9f13SBarry Smith   PetscFunctionReturn(0);
1087273d9f13SBarry Smith }
1088273d9f13SBarry Smith 
1089a6ece127SHong Zhang #undef __FUNCT__
1090a6ece127SHong Zhang #define __FUNCT__ "MatGetArray_SeqSBAIJ"
1091dfbe8321SBarry Smith PetscErrorCode MatGetArray_SeqSBAIJ(Mat A,PetscScalar *array[])
1092a6ece127SHong Zhang {
1093a6ece127SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
1094a6ece127SHong Zhang   PetscFunctionBegin;
1095a6ece127SHong Zhang   *array = a->a;
1096a6ece127SHong Zhang   PetscFunctionReturn(0);
1097a6ece127SHong Zhang }
1098a6ece127SHong Zhang 
1099a6ece127SHong Zhang #undef __FUNCT__
1100a6ece127SHong Zhang #define __FUNCT__ "MatRestoreArray_SeqSBAIJ"
1101dfbe8321SBarry Smith PetscErrorCode MatRestoreArray_SeqSBAIJ(Mat A,PetscScalar *array[])
1102a6ece127SHong Zhang {
1103a6ece127SHong Zhang   PetscFunctionBegin;
1104a6ece127SHong Zhang   PetscFunctionReturn(0);
1105a6ece127SHong Zhang  }
1106a6ece127SHong Zhang 
110742ee4b1aSHong Zhang #include "petscblaslapack.h"
110842ee4b1aSHong Zhang #undef __FUNCT__
110942ee4b1aSHong Zhang #define __FUNCT__ "MatAXPY_SeqSBAIJ"
1110dfbe8321SBarry Smith PetscErrorCode MatAXPY_SeqSBAIJ(const PetscScalar *a,Mat X,Mat Y,MatStructure str)
111142ee4b1aSHong Zhang {
111242ee4b1aSHong Zhang   Mat_SeqSBAIJ   *x=(Mat_SeqSBAIJ *)X->data, *y=(Mat_SeqSBAIJ *)Y->data;
1113dfbe8321SBarry Smith   PetscErrorCode ierr;
1114521d7252SBarry Smith   PetscInt       i,bs=Y->bs,bs2,j;
11154ce68768SBarry Smith   PetscBLASInt   bnz = (PetscBLASInt)x->nz,one = 1;
111642ee4b1aSHong Zhang 
111742ee4b1aSHong Zhang   PetscFunctionBegin;
111842ee4b1aSHong Zhang   if (str == SAME_NONZERO_PATTERN) {
111971044d3cSBarry Smith     BLASaxpy_(&bnz,(PetscScalar*)a,x->a,&one,y->a,&one);
1120c537a176SHong Zhang   } else if (str == SUBSET_NONZERO_PATTERN) { /* nonzeros of X is a subset of Y's */
1121c4319e64SHong Zhang     if (y->xtoy && y->XtoY != X) {
1122c4319e64SHong Zhang       ierr = PetscFree(y->xtoy);CHKERRQ(ierr);
1123c4319e64SHong Zhang       ierr = MatDestroy(y->XtoY);CHKERRQ(ierr);
1124c537a176SHong Zhang     }
1125c4319e64SHong Zhang     if (!y->xtoy) { /* get xtoy */
1126c4319e64SHong Zhang       ierr = MatAXPYGetxtoy_Private(x->mbs,x->i,x->j,PETSC_NULL, y->i,y->j,PETSC_NULL, &y->xtoy);CHKERRQ(ierr);
1127c4319e64SHong Zhang       y->XtoY = X;
1128c537a176SHong Zhang     }
1129c4319e64SHong Zhang     bs2 = bs*bs;
11306c6c5352SBarry Smith     for (i=0; i<x->nz; i++) {
1131c4319e64SHong Zhang       j = 0;
1132c4319e64SHong Zhang       while (j < bs2){
11336550863bSHong Zhang         y->a[bs2*y->xtoy[i]+j] += (*a)*(x->a[bs2*i+j]);
1134c4319e64SHong Zhang         j++;
1135c537a176SHong Zhang       }
1136c4319e64SHong Zhang     }
113763ba0a88SBarry 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);
113842ee4b1aSHong Zhang   } else {
113942ee4b1aSHong Zhang     ierr = MatAXPY_Basic(a,X,Y,str);CHKERRQ(ierr);
114042ee4b1aSHong Zhang   }
114142ee4b1aSHong Zhang   PetscFunctionReturn(0);
114242ee4b1aSHong Zhang }
114342ee4b1aSHong Zhang 
1144efcf0fc3SBarry Smith #undef __FUNCT__
1145efcf0fc3SBarry Smith #define __FUNCT__ "MatIsSymmetric_SeqSBAIJ"
1146dfbe8321SBarry Smith PetscErrorCode MatIsSymmetric_SeqSBAIJ(Mat A,PetscReal tol,PetscTruth *flg)
1147efcf0fc3SBarry Smith {
1148efcf0fc3SBarry Smith   PetscFunctionBegin;
1149efcf0fc3SBarry Smith   *flg = PETSC_TRUE;
1150efcf0fc3SBarry Smith   PetscFunctionReturn(0);
1151efcf0fc3SBarry Smith }
1152efcf0fc3SBarry Smith 
1153efcf0fc3SBarry Smith #undef __FUNCT__
1154efcf0fc3SBarry Smith #define __FUNCT__ "MatIsStructurallySymmetric_SeqSBAIJ"
1155dfbe8321SBarry Smith PetscErrorCode MatIsStructurallySymmetric_SeqSBAIJ(Mat A,PetscTruth *flg)
1156efcf0fc3SBarry Smith {
1157efcf0fc3SBarry Smith    PetscFunctionBegin;
1158efcf0fc3SBarry Smith    *flg = PETSC_TRUE;
1159efcf0fc3SBarry Smith    PetscFunctionReturn(0);
1160efcf0fc3SBarry Smith }
1161efcf0fc3SBarry Smith 
1162efcf0fc3SBarry Smith #undef __FUNCT__
1163efcf0fc3SBarry Smith #define __FUNCT__ "MatIsHermitian_SeqSBAIJ"
1164dfbe8321SBarry Smith PetscErrorCode MatIsHermitian_SeqSBAIJ(Mat A,PetscTruth *flg)
1165efcf0fc3SBarry Smith  {
1166efcf0fc3SBarry Smith    PetscFunctionBegin;
1167efcf0fc3SBarry Smith    *flg = PETSC_FALSE;
1168efcf0fc3SBarry Smith    PetscFunctionReturn(0);
1169efcf0fc3SBarry Smith  }
1170efcf0fc3SBarry Smith 
117149b5e25fSSatish Balay /* -------------------------------------------------------------------*/
117249b5e25fSSatish Balay static struct _MatOps MatOps_Values = {MatSetValues_SeqSBAIJ,
117349b5e25fSSatish Balay        MatGetRow_SeqSBAIJ,
117449b5e25fSSatish Balay        MatRestoreRow_SeqSBAIJ,
117549b5e25fSSatish Balay        MatMult_SeqSBAIJ_N,
117697304618SKris Buschelman /* 4*/ MatMultAdd_SeqSBAIJ_N,
1177e005ede5SBarry Smith        MatMult_SeqSBAIJ_N,
1178e005ede5SBarry Smith        MatMultAdd_SeqSBAIJ_N,
117949b5e25fSSatish Balay        MatSolve_SeqSBAIJ_N,
118049b5e25fSSatish Balay        0,
118149b5e25fSSatish Balay        0,
118297304618SKris Buschelman /*10*/ 0,
118349b5e25fSSatish Balay        0,
11845f4ef2deSHong Zhang        MatCholeskyFactor_SeqSBAIJ,
1185d06b337dSHong Zhang        MatRelax_SeqSBAIJ,
118649b5e25fSSatish Balay        MatTranspose_SeqSBAIJ,
118797304618SKris Buschelman /*15*/ MatGetInfo_SeqSBAIJ,
118849b5e25fSSatish Balay        MatEqual_SeqSBAIJ,
118949b5e25fSSatish Balay        MatGetDiagonal_SeqSBAIJ,
119049b5e25fSSatish Balay        MatDiagonalScale_SeqSBAIJ,
119149b5e25fSSatish Balay        MatNorm_SeqSBAIJ,
119297304618SKris Buschelman /*20*/ 0,
119349b5e25fSSatish Balay        MatAssemblyEnd_SeqSBAIJ,
119449b5e25fSSatish Balay        0,
119549b5e25fSSatish Balay        MatSetOption_SeqSBAIJ,
119649b5e25fSSatish Balay        MatZeroEntries_SeqSBAIJ,
1197dcf5cc72SBarry Smith /*25*/ 0,
119849b5e25fSSatish Balay        0,
119949b5e25fSSatish Balay        0,
12005f4ef2deSHong Zhang        MatCholeskyFactorSymbolic_SeqSBAIJ,
12015f4ef2deSHong Zhang        MatCholeskyFactorNumeric_SeqSBAIJ_N,
120297304618SKris Buschelman /*30*/ MatSetUpPreallocation_SeqSBAIJ,
1203c464158bSHong Zhang        0,
12044d101231SSatish Balay        MatICCFactorSymbolic_SeqSBAIJ,
1205a6ece127SHong Zhang        MatGetArray_SeqSBAIJ,
1206a6ece127SHong Zhang        MatRestoreArray_SeqSBAIJ,
120797304618SKris Buschelman /*35*/ MatDuplicate_SeqSBAIJ,
120849b5e25fSSatish Balay        0,
120949b5e25fSSatish Balay        0,
121049b5e25fSSatish Balay        0,
1211950f1e5bSHong Zhang        0,
121297304618SKris Buschelman /*40*/ MatAXPY_SeqSBAIJ,
121349b5e25fSSatish Balay        MatGetSubMatrices_SeqSBAIJ,
121449b5e25fSSatish Balay        MatIncreaseOverlap_SeqSBAIJ,
121549b5e25fSSatish Balay        MatGetValues_SeqSBAIJ,
121649b5e25fSSatish Balay        0,
121797304618SKris Buschelman /*45*/ MatPrintHelp_SeqSBAIJ,
121849b5e25fSSatish Balay        MatScale_SeqSBAIJ,
121949b5e25fSSatish Balay        0,
122049b5e25fSSatish Balay        0,
122149b5e25fSSatish Balay        0,
1222521d7252SBarry Smith /*50*/ 0,
122349b5e25fSSatish Balay        MatGetRowIJ_SeqSBAIJ,
122449b5e25fSSatish Balay        MatRestoreRowIJ_SeqSBAIJ,
122549b5e25fSSatish Balay        0,
122649b5e25fSSatish Balay        0,
122797304618SKris Buschelman /*55*/ 0,
122849b5e25fSSatish Balay        0,
122949b5e25fSSatish Balay        0,
123049b5e25fSSatish Balay        0,
123149b5e25fSSatish Balay        MatSetValuesBlocked_SeqSBAIJ,
123297304618SKris Buschelman /*60*/ MatGetSubMatrix_SeqSBAIJ,
123349b5e25fSSatish Balay        0,
123449b5e25fSSatish Balay        0,
12358a124369SBarry Smith        MatGetPetscMaps_Petsc,
1236d959ec07SHong Zhang        0,
123797304618SKris Buschelman /*65*/ 0,
1238d959ec07SHong Zhang        0,
1239d959ec07SHong Zhang        0,
1240d959ec07SHong Zhang        0,
1241d959ec07SHong Zhang        0,
124297304618SKris Buschelman /*70*/ MatGetRowMax_SeqSBAIJ,
12433e0d88b5SBarry Smith        0,
12443e0d88b5SBarry Smith        0,
12453e0d88b5SBarry Smith        0,
12463e0d88b5SBarry Smith        0,
124797304618SKris Buschelman /*75*/ 0,
12483e0d88b5SBarry Smith        0,
12493e0d88b5SBarry Smith        0,
12503e0d88b5SBarry Smith        0,
12513e0d88b5SBarry Smith        0,
125297304618SKris Buschelman /*80*/ 0,
12533e0d88b5SBarry Smith        0,
12543e0d88b5SBarry Smith        0,
12553e0d88b5SBarry Smith #if !defined(PETSC_USE_COMPLEX)
125697304618SKris Buschelman        MatGetInertia_SeqSBAIJ,
12573e0d88b5SBarry Smith #else
125897304618SKris Buschelman        0,
12593e0d88b5SBarry Smith #endif
1260865e5f61SKris Buschelman        MatLoad_SeqSBAIJ,
1261865e5f61SKris Buschelman /*85*/ MatIsSymmetric_SeqSBAIJ,
1262865e5f61SKris Buschelman        MatIsHermitian_SeqSBAIJ,
1263efcf0fc3SBarry Smith        MatIsStructurallySymmetric_SeqSBAIJ,
1264865e5f61SKris Buschelman        0,
1265865e5f61SKris Buschelman        0,
1266865e5f61SKris Buschelman /*90*/ 0,
1267865e5f61SKris Buschelman        0,
1268865e5f61SKris Buschelman        0,
1269865e5f61SKris Buschelman        0,
1270865e5f61SKris Buschelman        0,
1271865e5f61SKris Buschelman /*95*/ 0,
1272865e5f61SKris Buschelman        0,
1273865e5f61SKris Buschelman        0,
1274865e5f61SKris Buschelman        0};
1275be1d678aSKris Buschelman 
127649b5e25fSSatish Balay EXTERN_C_BEGIN
12774a2ae208SSatish Balay #undef __FUNCT__
12784a2ae208SSatish Balay #define __FUNCT__ "MatStoreValues_SeqSBAIJ"
1279be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatStoreValues_SeqSBAIJ(Mat mat)
128049b5e25fSSatish Balay {
12814afc71dfSHong Zhang   Mat_SeqSBAIJ   *aij = (Mat_SeqSBAIJ *)mat->data;
1282521d7252SBarry Smith   PetscInt       nz = aij->i[mat->m]*mat->bs*aij->bs2;
1283dfbe8321SBarry Smith   PetscErrorCode ierr;
128449b5e25fSSatish Balay 
128549b5e25fSSatish Balay   PetscFunctionBegin;
128649b5e25fSSatish Balay   if (aij->nonew != 1) {
1287e005ede5SBarry Smith     SETERRQ(PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
128849b5e25fSSatish Balay   }
128949b5e25fSSatish Balay 
129049b5e25fSSatish Balay   /* allocate space for values if not already there */
129149b5e25fSSatish Balay   if (!aij->saved_values) {
129287828ca2SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&aij->saved_values);CHKERRQ(ierr);
129349b5e25fSSatish Balay   }
129449b5e25fSSatish Balay 
129549b5e25fSSatish Balay   /* copy values over */
129687828ca2SBarry Smith   ierr = PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(PetscScalar));CHKERRQ(ierr);
129749b5e25fSSatish Balay   PetscFunctionReturn(0);
129849b5e25fSSatish Balay }
129949b5e25fSSatish Balay EXTERN_C_END
130049b5e25fSSatish Balay 
130149b5e25fSSatish Balay EXTERN_C_BEGIN
13024a2ae208SSatish Balay #undef __FUNCT__
13034a2ae208SSatish Balay #define __FUNCT__ "MatRetrieveValues_SeqSBAIJ"
1304be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatRetrieveValues_SeqSBAIJ(Mat mat)
130549b5e25fSSatish Balay {
13064afc71dfSHong Zhang   Mat_SeqSBAIJ   *aij = (Mat_SeqSBAIJ *)mat->data;
13076849ba73SBarry Smith   PetscErrorCode ierr;
1308521d7252SBarry Smith   PetscInt       nz = aij->i[mat->m]*mat->bs*aij->bs2;
130949b5e25fSSatish Balay 
131049b5e25fSSatish Balay   PetscFunctionBegin;
131149b5e25fSSatish Balay   if (aij->nonew != 1) {
1312e005ede5SBarry Smith     SETERRQ(PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
131349b5e25fSSatish Balay   }
131449b5e25fSSatish Balay   if (!aij->saved_values) {
1315e005ede5SBarry Smith     SETERRQ(PETSC_ERR_ORDER,"Must call MatStoreValues(A);first");
131649b5e25fSSatish Balay   }
131749b5e25fSSatish Balay 
131849b5e25fSSatish Balay   /* copy values over */
131987828ca2SBarry Smith   ierr = PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(PetscScalar));CHKERRQ(ierr);
132049b5e25fSSatish Balay   PetscFunctionReturn(0);
132149b5e25fSSatish Balay }
132249b5e25fSSatish Balay EXTERN_C_END
132349b5e25fSSatish Balay 
13248549e402SHong Zhang EXTERN_C_BEGIN
13254a2ae208SSatish Balay #undef __FUNCT__
1326a23d5eceSKris Buschelman #define __FUNCT__ "MatSeqSBAIJSetPreallocation_SeqSBAIJ"
1327be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetPreallocation_SeqSBAIJ(Mat B,PetscInt bs,PetscInt nz,PetscInt *nnz)
132849b5e25fSSatish Balay {
1329c464158bSHong Zhang   Mat_SeqSBAIJ   *b = (Mat_SeqSBAIJ*)B->data;
13306849ba73SBarry Smith   PetscErrorCode ierr;
1331ab93d7beSBarry Smith   PetscInt       i,mbs,bs2;
1332ab93d7beSBarry Smith   PetscTruth     skipallocation = PETSC_FALSE,flg;
133349b5e25fSSatish Balay 
133449b5e25fSSatish Balay   PetscFunctionBegin;
1335273d9f13SBarry Smith   B->preallocated = PETSC_TRUE;
1336e82a3eeeSBarry Smith   ierr = PetscOptionsGetInt(B->prefix,"-mat_block_size",&bs,PETSC_NULL);CHKERRQ(ierr);
1337c464158bSHong Zhang   mbs  = B->m/bs;
133849b5e25fSSatish Balay   bs2  = bs*bs;
133949b5e25fSSatish Balay 
1340c464158bSHong Zhang   if (mbs*bs != B->m) {
134129bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_SIZ,"Number rows, cols must be divisible by blocksize");
134249b5e25fSSatish Balay   }
134349b5e25fSSatish Balay 
1344ab93d7beSBarry Smith   if (nz == MAT_SKIP_ALLOCATION) {
1345ab93d7beSBarry Smith     skipallocation = PETSC_TRUE;
1346ab93d7beSBarry Smith     nz             = 0;
1347ab93d7beSBarry Smith   }
1348ab93d7beSBarry Smith 
1349435da068SBarry Smith   if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 3;
135077431f27SBarry Smith   if (nz < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %D",nz);
135149b5e25fSSatish Balay   if (nnz) {
135249b5e25fSSatish Balay     for (i=0; i<mbs; i++) {
135377431f27SBarry Smith       if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %D value %D",i,nnz[i]);
135477431f27SBarry 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);
135549b5e25fSSatish Balay     }
135649b5e25fSSatish Balay   }
135749b5e25fSSatish Balay 
1358e82a3eeeSBarry Smith   ierr    = PetscOptionsHasName(B->prefix,"-mat_no_unroll",&flg);CHKERRQ(ierr);
135949b5e25fSSatish Balay   if (!flg) {
136049b5e25fSSatish Balay     switch (bs) {
136149b5e25fSSatish Balay     case 1:
13624ccecd49SHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1;
136349b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_1;
1364d59c15a7SBarry Smith       B->ops->solves          = MatSolves_SeqSBAIJ_1;
1365e005ede5SBarry Smith       B->ops->solvetranspose  = MatSolve_SeqSBAIJ_1;
136649b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_1;
136749b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_1;
136849b5e25fSSatish Balay       break;
136949b5e25fSSatish Balay     case 2:
13704ccecd49SHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2;
137149b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_2;
1372e005ede5SBarry Smith       B->ops->solvetranspose  = MatSolve_SeqSBAIJ_2;
137349b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_2;
137449b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_2;
137549b5e25fSSatish Balay       break;
137649b5e25fSSatish Balay     case 3:
13775f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3;
137849b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_3;
1379e005ede5SBarry Smith       B->ops->solvetranspose  = MatSolve_SeqSBAIJ_3;
138049b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_3;
138149b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_3;
138249b5e25fSSatish Balay       break;
138349b5e25fSSatish Balay     case 4:
13845f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4;
138549b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_4;
1386e005ede5SBarry Smith       B->ops->solvetranspose  = MatSolve_SeqSBAIJ_4;
138749b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_4;
138849b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_4;
138949b5e25fSSatish Balay       break;
139049b5e25fSSatish Balay     case 5:
13915f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5;
139249b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_5;
1393e005ede5SBarry Smith       B->ops->solvetranspose  = MatSolve_SeqSBAIJ_5;
139449b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_5;
139549b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_5;
139649b5e25fSSatish Balay       break;
139749b5e25fSSatish Balay     case 6:
13985f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6;
139949b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_6;
1400e005ede5SBarry Smith       B->ops->solvetranspose  = MatSolve_SeqSBAIJ_6;
140149b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_6;
140249b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_6;
140349b5e25fSSatish Balay       break;
140449b5e25fSSatish Balay     case 7:
1405de53e5efSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7;
140649b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_7;
1407e005ede5SBarry Smith       B->ops->solvetranspose  = MatSolve_SeqSBAIJ_7;
1408de53e5efSHong Zhang       B->ops->mult            = MatMult_SeqSBAIJ_7;
140949b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_7;
141049b5e25fSSatish Balay       break;
141149b5e25fSSatish Balay     }
141249b5e25fSSatish Balay   }
141349b5e25fSSatish Balay 
141449b5e25fSSatish Balay   b->mbs = mbs;
14154afc71dfSHong Zhang   b->nbs = mbs;
1416ab93d7beSBarry Smith   if (!skipallocation) {
1417ab93d7beSBarry Smith     /* b->ilen will count nonzeros in each block row so far. */
1418ab93d7beSBarry Smith     ierr   = PetscMalloc2(mbs,PetscInt,&b->imax,mbs,PetscInt,&b->ilen);CHKERRQ(ierr);
1419ab93d7beSBarry Smith     for (i=0; i<mbs; i++) { b->ilen[i] = 0;}
142049b5e25fSSatish Balay     if (!nnz) {
1421435da068SBarry Smith       if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5;
142249b5e25fSSatish Balay       else if (nz <= 0)        nz = 1;
142349b5e25fSSatish Balay       for (i=0; i<mbs; i++) {
14248cef66ccSHong Zhang         b->imax[i] = nz;
142549b5e25fSSatish Balay       }
1426153ea458SHong Zhang       nz = nz*mbs; /* total nz */
142749b5e25fSSatish Balay     } else {
142849b5e25fSSatish Balay       nz = 0;
14298cef66ccSHong Zhang       for (i=0; i<mbs; i++) {b->imax[i] = nnz[i]; nz += nnz[i];}
143049b5e25fSSatish Balay     }
14316c6c5352SBarry Smith     /* nz=(nz+mbs)/2; */ /* total diagonal and superdiagonal nonzero blocks */
143249b5e25fSSatish Balay 
143349b5e25fSSatish Balay     /* allocate the matrix space */
1434a96a251dSBarry Smith     ierr = PetscMalloc3(bs2*nz,PetscScalar,&b->a,nz,PetscInt,&b->j,B->m+1,PetscInt,&b->i);CHKERRQ(ierr);
14356c6c5352SBarry Smith     ierr = PetscMemzero(b->a,nz*bs2*sizeof(MatScalar));CHKERRQ(ierr);
143613f74950SBarry Smith     ierr = PetscMemzero(b->j,nz*sizeof(PetscInt));CHKERRQ(ierr);
143749b5e25fSSatish Balay     b->singlemalloc = PETSC_TRUE;
143849b5e25fSSatish Balay 
143949b5e25fSSatish Balay     /* pointer to beginning of each row */
144049b5e25fSSatish Balay     b->i[0] = 0;
144149b5e25fSSatish Balay     for (i=1; i<mbs+1; i++) {
144249b5e25fSSatish Balay       b->i[i] = b->i[i-1] + b->imax[i-1];
144349b5e25fSSatish Balay     }
1444ab93d7beSBarry Smith   }
144549b5e25fSSatish Balay 
1446521d7252SBarry Smith   B->bs               = bs;
144749b5e25fSSatish Balay   b->bs2              = bs2;
14486c6c5352SBarry Smith   b->nz             = 0;
14496c6c5352SBarry Smith   b->maxnz          = nz*bs2;
1450153ea458SHong Zhang 
145116cdd363SHong Zhang   b->inew             = 0;
145216cdd363SHong Zhang   b->jnew             = 0;
145316cdd363SHong Zhang   b->anew             = 0;
145416cdd363SHong Zhang   b->a2anew           = 0;
14551a3463dfSHong Zhang   b->permute          = PETSC_FALSE;
1456c464158bSHong Zhang   PetscFunctionReturn(0);
1457c464158bSHong Zhang }
1458a23d5eceSKris Buschelman EXTERN_C_END
1459153ea458SHong Zhang 
1460d769727bSBarry Smith EXTERN_C_BEGIN
1461be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqSBAIJ_SeqAIJ(Mat,const MatType,MatReuse,Mat*);
1462be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqSBAIJ_SeqBAIJ(Mat,const MatType,MatReuse,Mat*);
1463d769727bSBarry Smith EXTERN_C_END
1464d769727bSBarry Smith 
14650bad9183SKris Buschelman /*MC
1466fafad747SKris Buschelman   MATSEQSBAIJ - MATSEQSBAIJ = "seqsbaij" - A matrix type to be used for sequential symmetric block sparse matrices,
14670bad9183SKris Buschelman   based on block compressed sparse row format.  Only the upper triangular portion of the matrix is stored.
14680bad9183SKris Buschelman 
14690bad9183SKris Buschelman   Options Database Keys:
14700bad9183SKris Buschelman   . -mat_type seqsbaij - sets the matrix type to "seqsbaij" during a call to MatSetFromOptions()
14710bad9183SKris Buschelman 
14720bad9183SKris Buschelman   Level: beginner
14730bad9183SKris Buschelman 
14740bad9183SKris Buschelman   .seealso: MatCreateSeqSBAIJ
14750bad9183SKris Buschelman M*/
14760bad9183SKris Buschelman 
1477a23d5eceSKris Buschelman EXTERN_C_BEGIN
1478a23d5eceSKris Buschelman #undef __FUNCT__
1479a23d5eceSKris Buschelman #define __FUNCT__ "MatCreate_SeqSBAIJ"
1480be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_SeqSBAIJ(Mat B)
1481a23d5eceSKris Buschelman {
1482a23d5eceSKris Buschelman   Mat_SeqSBAIJ   *b;
1483dfbe8321SBarry Smith   PetscErrorCode ierr;
148413f74950SBarry Smith   PetscMPIInt    size;
1485941593c8SHong Zhang   PetscTruth     flg;
1486a23d5eceSKris Buschelman 
1487a23d5eceSKris Buschelman   PetscFunctionBegin;
1488a23d5eceSKris Buschelman   ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr);
1489a23d5eceSKris Buschelman   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1");
1490a23d5eceSKris Buschelman   B->m = B->M = PetscMax(B->m,B->M);
1491a23d5eceSKris Buschelman   B->n = B->N = PetscMax(B->n,B->N);
1492a23d5eceSKris Buschelman 
1493a23d5eceSKris Buschelman   ierr    = PetscNew(Mat_SeqSBAIJ,&b);CHKERRQ(ierr);
1494a23d5eceSKris Buschelman   B->data = (void*)b;
1495a23d5eceSKris Buschelman   ierr    = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
1496a23d5eceSKris Buschelman   B->ops->destroy     = MatDestroy_SeqSBAIJ;
1497a23d5eceSKris Buschelman   B->ops->view        = MatView_SeqSBAIJ;
1498a23d5eceSKris Buschelman   B->factor           = 0;
1499a23d5eceSKris Buschelman   B->mapping          = 0;
1500a23d5eceSKris Buschelman   b->row              = 0;
1501a23d5eceSKris Buschelman   b->icol             = 0;
1502a23d5eceSKris Buschelman   b->reallocs         = 0;
1503a23d5eceSKris Buschelman   b->saved_values     = 0;
1504a23d5eceSKris Buschelman 
1505a23d5eceSKris Buschelman   ierr = PetscMapCreateMPI(B->comm,B->m,B->M,&B->rmap);CHKERRQ(ierr);
1506a23d5eceSKris Buschelman   ierr = PetscMapCreateMPI(B->comm,B->n,B->N,&B->cmap);CHKERRQ(ierr);
1507a23d5eceSKris Buschelman 
1508a23d5eceSKris Buschelman   b->sorted           = PETSC_FALSE;
1509a23d5eceSKris Buschelman   b->roworiented      = PETSC_TRUE;
1510a23d5eceSKris Buschelman   b->nonew            = 0;
1511a23d5eceSKris Buschelman   b->diag             = 0;
1512a23d5eceSKris Buschelman   b->solve_work       = 0;
1513a23d5eceSKris Buschelman   b->mult_work        = 0;
1514a23d5eceSKris Buschelman   B->spptr            = 0;
1515a23d5eceSKris Buschelman   b->keepzeroedrows   = PETSC_FALSE;
1516a23d5eceSKris Buschelman   b->xtoy             = 0;
1517a23d5eceSKris Buschelman   b->XtoY             = 0;
1518a23d5eceSKris Buschelman 
1519a23d5eceSKris Buschelman   b->inew             = 0;
1520a23d5eceSKris Buschelman   b->jnew             = 0;
1521a23d5eceSKris Buschelman   b->anew             = 0;
1522a23d5eceSKris Buschelman   b->a2anew           = 0;
1523a23d5eceSKris Buschelman   b->permute          = PETSC_FALSE;
1524a23d5eceSKris Buschelman 
1525941593c8SHong Zhang   b->ignore_ltriangular = PETSC_FALSE;
1526941593c8SHong Zhang   ierr = PetscOptionsHasName(PETSC_NULL,"-mat_ignore_lower_triangular",&flg);CHKERRQ(ierr);
1527941593c8SHong Zhang   if (flg) b->ignore_ltriangular = PETSC_TRUE;
1528941593c8SHong Zhang 
1529a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C",
1530a23d5eceSKris Buschelman                                      "MatStoreValues_SeqSBAIJ",
1531a23d5eceSKris Buschelman                                      MatStoreValues_SeqSBAIJ);CHKERRQ(ierr);
1532a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C",
1533a23d5eceSKris Buschelman                                      "MatRetrieveValues_SeqSBAIJ",
1534a23d5eceSKris Buschelman                                      (void*)MatRetrieveValues_SeqSBAIJ);CHKERRQ(ierr);
1535a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqSBAIJSetColumnIndices_C",
1536a23d5eceSKris Buschelman                                      "MatSeqSBAIJSetColumnIndices_SeqSBAIJ",
1537a23d5eceSKris Buschelman                                      MatSeqSBAIJSetColumnIndices_SeqSBAIJ);CHKERRQ(ierr);
15384e5e7fe4SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqsbaij_seqaij_C",
15394e5e7fe4SHong Zhang                                      "MatConvert_SeqSBAIJ_SeqAIJ",
15404e5e7fe4SHong Zhang                                       MatConvert_SeqSBAIJ_SeqAIJ);CHKERRQ(ierr);
1541a0e1a404SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqsbaij_seqbaij_C",
1542a0e1a404SHong Zhang                                      "MatConvert_SeqSBAIJ_SeqBAIJ",
1543a0e1a404SHong Zhang                                       MatConvert_SeqSBAIJ_SeqBAIJ);CHKERRQ(ierr);
1544a23d5eceSKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqSBAIJSetPreallocation_C",
1545a23d5eceSKris Buschelman                                      "MatSeqSBAIJSetPreallocation_SeqSBAIJ",
1546a23d5eceSKris Buschelman                                      MatSeqSBAIJSetPreallocation_SeqSBAIJ);CHKERRQ(ierr);
154723ce1328SBarry Smith 
154823ce1328SBarry Smith   B->symmetric                  = PETSC_TRUE;
154923ce1328SBarry Smith   B->structurally_symmetric     = PETSC_TRUE;
155023ce1328SBarry Smith   B->symmetric_set              = PETSC_TRUE;
155123ce1328SBarry Smith   B->structurally_symmetric_set = PETSC_TRUE;
1552a23d5eceSKris Buschelman   PetscFunctionReturn(0);
1553a23d5eceSKris Buschelman }
1554a23d5eceSKris Buschelman EXTERN_C_END
1555a23d5eceSKris Buschelman 
1556a23d5eceSKris Buschelman #undef __FUNCT__
1557a23d5eceSKris Buschelman #define __FUNCT__ "MatSeqSBAIJSetPreallocation"
1558a23d5eceSKris Buschelman /*@C
1559a23d5eceSKris Buschelman    MatSeqSBAIJSetPreallocation - Creates a sparse symmetric matrix in block AIJ (block
1560a23d5eceSKris Buschelman    compressed row) format.  For good matrix assembly performance the
1561a23d5eceSKris Buschelman    user should preallocate the matrix storage by setting the parameter nz
1562a23d5eceSKris Buschelman    (or the array nnz).  By setting these parameters accurately, performance
1563a23d5eceSKris Buschelman    during matrix assembly can be increased by more than a factor of 50.
1564a23d5eceSKris Buschelman 
1565a23d5eceSKris Buschelman    Collective on Mat
1566a23d5eceSKris Buschelman 
1567a23d5eceSKris Buschelman    Input Parameters:
1568a23d5eceSKris Buschelman +  A - the symmetric matrix
1569a23d5eceSKris Buschelman .  bs - size of block
1570a23d5eceSKris Buschelman .  nz - number of block nonzeros per block row (same for all rows)
1571a23d5eceSKris Buschelman -  nnz - array containing the number of block nonzeros in the upper triangular plus
1572a23d5eceSKris Buschelman          diagonal portion of each block (possibly different for each block row) or PETSC_NULL
1573a23d5eceSKris Buschelman 
1574a23d5eceSKris Buschelman    Options Database Keys:
1575a23d5eceSKris Buschelman .   -mat_no_unroll - uses code that does not unroll the loops in the
1576a23d5eceSKris Buschelman                      block calculations (much slower)
1577a23d5eceSKris Buschelman .    -mat_block_size - size of the blocks to use
1578a23d5eceSKris Buschelman 
1579a23d5eceSKris Buschelman    Level: intermediate
1580a23d5eceSKris Buschelman 
1581a23d5eceSKris Buschelman    Notes:
1582a23d5eceSKris Buschelman    Specify the preallocated storage with either nz or nnz (not both).
1583a23d5eceSKris Buschelman    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
1584a23d5eceSKris Buschelman    allocation.  For additional details, see the users manual chapter on
1585a23d5eceSKris Buschelman    matrices.
1586a23d5eceSKris Buschelman 
158749a6f317SBarry Smith    If the nnz parameter is given then the nz parameter is ignored
158849a6f317SBarry Smith 
158949a6f317SBarry Smith 
1590a23d5eceSKris Buschelman .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ()
1591a23d5eceSKris Buschelman @*/
1592be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetPreallocation(Mat B,PetscInt bs,PetscInt nz,const PetscInt nnz[])
159313f74950SBarry Smith {
159413f74950SBarry Smith   PetscErrorCode ierr,(*f)(Mat,PetscInt,PetscInt,const PetscInt[]);
1595a23d5eceSKris Buschelman 
1596a23d5eceSKris Buschelman   PetscFunctionBegin;
1597a23d5eceSKris Buschelman   ierr = PetscObjectQueryFunction((PetscObject)B,"MatSeqSBAIJSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr);
1598a23d5eceSKris Buschelman   if (f) {
1599a23d5eceSKris Buschelman     ierr = (*f)(B,bs,nz,nnz);CHKERRQ(ierr);
1600a23d5eceSKris Buschelman   }
1601a23d5eceSKris Buschelman   PetscFunctionReturn(0);
1602a23d5eceSKris Buschelman }
160349b5e25fSSatish Balay 
16044a2ae208SSatish Balay #undef __FUNCT__
16054a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqSBAIJ"
1606c464158bSHong Zhang /*@C
1607c464158bSHong Zhang    MatCreateSeqSBAIJ - Creates a sparse symmetric matrix in block AIJ (block
1608c464158bSHong Zhang    compressed row) format.  For good matrix assembly performance the
1609c464158bSHong Zhang    user should preallocate the matrix storage by setting the parameter nz
1610c464158bSHong Zhang    (or the array nnz).  By setting these parameters accurately, performance
1611c464158bSHong Zhang    during matrix assembly can be increased by more than a factor of 50.
161249b5e25fSSatish Balay 
1613c464158bSHong Zhang    Collective on MPI_Comm
1614c464158bSHong Zhang 
1615c464158bSHong Zhang    Input Parameters:
1616c464158bSHong Zhang +  comm - MPI communicator, set to PETSC_COMM_SELF
1617c464158bSHong Zhang .  bs - size of block
1618c464158bSHong Zhang .  m - number of rows, or number of columns
1619c464158bSHong Zhang .  nz - number of block nonzeros per block row (same for all rows)
1620744e8345SSatish Balay -  nnz - array containing the number of block nonzeros in the upper triangular plus
1621744e8345SSatish Balay          diagonal portion of each block (possibly different for each block row) or PETSC_NULL
1622c464158bSHong Zhang 
1623c464158bSHong Zhang    Output Parameter:
1624c464158bSHong Zhang .  A - the symmetric matrix
1625c464158bSHong Zhang 
1626c464158bSHong Zhang    Options Database Keys:
1627c464158bSHong Zhang .   -mat_no_unroll - uses code that does not unroll the loops in the
1628c464158bSHong Zhang                      block calculations (much slower)
1629c464158bSHong Zhang .    -mat_block_size - size of the blocks to use
1630c464158bSHong Zhang 
1631c464158bSHong Zhang    Level: intermediate
1632c464158bSHong Zhang 
1633c464158bSHong Zhang    Notes:
1634c464158bSHong Zhang 
1635c464158bSHong Zhang    Specify the preallocated storage with either nz or nnz (not both).
1636c464158bSHong Zhang    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
1637c464158bSHong Zhang    allocation.  For additional details, see the users manual chapter on
1638c464158bSHong Zhang    matrices.
1639c464158bSHong Zhang 
164049a6f317SBarry Smith    If the nnz parameter is given then the nz parameter is ignored
164149a6f317SBarry Smith 
1642c464158bSHong Zhang .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ()
1643c464158bSHong Zhang @*/
1644be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJ(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],Mat *A)
1645c464158bSHong Zhang {
1646dfbe8321SBarry Smith   PetscErrorCode ierr;
1647c464158bSHong Zhang 
1648c464158bSHong Zhang   PetscFunctionBegin;
1649c464158bSHong Zhang   ierr = MatCreate(comm,m,n,m,n,A);CHKERRQ(ierr);
1650c464158bSHong Zhang   ierr = MatSetType(*A,MATSEQSBAIJ);CHKERRQ(ierr);
1651ab93d7beSBarry Smith   ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(*A,bs,nz,(PetscInt*)nnz);CHKERRQ(ierr);
165249b5e25fSSatish Balay   PetscFunctionReturn(0);
165349b5e25fSSatish Balay }
165449b5e25fSSatish Balay 
16554a2ae208SSatish Balay #undef __FUNCT__
16564a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqSBAIJ"
1657dfbe8321SBarry Smith PetscErrorCode MatDuplicate_SeqSBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B)
165849b5e25fSSatish Balay {
165949b5e25fSSatish Balay   Mat            C;
166049b5e25fSSatish Balay   Mat_SeqSBAIJ   *c,*a = (Mat_SeqSBAIJ*)A->data;
16616849ba73SBarry Smith   PetscErrorCode ierr;
1662b40805acSSatish Balay   PetscInt       i,mbs = a->mbs,nz = a->nz,bs2 =a->bs2;
166349b5e25fSSatish Balay 
166449b5e25fSSatish Balay   PetscFunctionBegin;
166529bbc08cSBarry Smith   if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,"Corrupt matrix");
166649b5e25fSSatish Balay 
166749b5e25fSSatish Balay   *B = 0;
1668692f9cbeSHong Zhang   ierr = MatCreate(A->comm,A->m,A->n,A->m,A->n,&C);CHKERRQ(ierr);
1669be5d1d56SKris Buschelman   ierr = MatSetType(C,A->type_name);CHKERRQ(ierr);
16701d5dac46SHong Zhang   ierr = PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
1671692f9cbeSHong Zhang   c    = (Mat_SeqSBAIJ*)C->data;
1672692f9cbeSHong Zhang 
1673273d9f13SBarry Smith   C->preallocated   = PETSC_TRUE;
167449b5e25fSSatish Balay   C->factor         = A->factor;
167549b5e25fSSatish Balay   c->row            = 0;
167649b5e25fSSatish Balay   c->icol           = 0;
167749b5e25fSSatish Balay   c->saved_values   = 0;
167849b5e25fSSatish Balay   c->keepzeroedrows = a->keepzeroedrows;
167949b5e25fSSatish Balay   C->assembled      = PETSC_TRUE;
168049b5e25fSSatish Balay 
168182327fa8SHong Zhang   C->M    = A->M;
168282327fa8SHong Zhang   C->N    = A->N;
1683521d7252SBarry Smith   C->bs   = A->bs;
168449b5e25fSSatish Balay   c->bs2  = a->bs2;
168549b5e25fSSatish Balay   c->mbs  = a->mbs;
168649b5e25fSSatish Balay   c->nbs  = a->nbs;
168749b5e25fSSatish Balay 
168813f74950SBarry Smith   ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&c->imax);CHKERRQ(ierr);
168913f74950SBarry Smith   ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&c->ilen);CHKERRQ(ierr);
169049b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
169149b5e25fSSatish Balay     c->imax[i] = a->imax[i];
169249b5e25fSSatish Balay     c->ilen[i] = a->ilen[i];
169349b5e25fSSatish Balay   }
1694b40805acSSatish Balay   ierr = PetscLogObjectMemory(C,2*(mbs+1)*sizeof(PetscInt)+sizeof(struct _p_Mat)+sizeof(Mat_SeqSBAIJ));CHKERRQ(ierr);
169549b5e25fSSatish Balay 
169649b5e25fSSatish Balay   /* allocate the matrix space */
1697b40805acSSatish Balay   ierr = PetscMalloc3(bs2*nz,MatScalar,&c->a,nz,PetscInt,&c->j,mbs+1,PetscInt,&c->i);CHKERRQ(ierr);
169849b5e25fSSatish Balay   c->singlemalloc = PETSC_TRUE;
169913f74950SBarry Smith   ierr = PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(PetscInt));CHKERRQ(ierr);
1700b40805acSSatish Balay   ierr = PetscLogObjectMemory(C,(mbs+1)*sizeof(PetscInt) + nz*(bs2*sizeof(MatScalar) + sizeof(PetscInt)));CHKERRQ(ierr);
170149b5e25fSSatish Balay   if (mbs > 0) {
170213f74950SBarry Smith     ierr = PetscMemcpy(c->j,a->j,nz*sizeof(PetscInt));CHKERRQ(ierr);
170349b5e25fSSatish Balay     if (cpvalues == MAT_COPY_VALUES) {
170449b5e25fSSatish Balay       ierr = PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr);
170549b5e25fSSatish Balay     } else {
170649b5e25fSSatish Balay       ierr = PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr);
170749b5e25fSSatish Balay     }
170849b5e25fSSatish Balay   }
170949b5e25fSSatish Balay 
171049b5e25fSSatish Balay   c->sorted      = a->sorted;
171149b5e25fSSatish Balay   c->roworiented = a->roworiented;
171249b5e25fSSatish Balay   c->nonew       = a->nonew;
171349b5e25fSSatish Balay 
171449b5e25fSSatish Balay   if (a->diag) {
171513f74950SBarry Smith     ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&c->diag);CHKERRQ(ierr);
171652e6d16bSBarry Smith     ierr = PetscLogObjectMemory(C,(mbs+1)*sizeof(PetscInt));CHKERRQ(ierr);
171749b5e25fSSatish Balay     for (i=0; i<mbs; i++) {
171849b5e25fSSatish Balay       c->diag[i] = a->diag[i];
171949b5e25fSSatish Balay     }
172049b5e25fSSatish Balay   } else c->diag        = 0;
17216c6c5352SBarry Smith   c->nz               = a->nz;
17226c6c5352SBarry Smith   c->maxnz            = a->maxnz;
172349b5e25fSSatish Balay   c->solve_work         = 0;
172449b5e25fSSatish Balay   c->mult_work          = 0;
172549b5e25fSSatish Balay   *B = C;
1726b0a32e0cSBarry Smith   ierr = PetscFListDuplicate(A->qlist,&C->qlist);CHKERRQ(ierr);
172749b5e25fSSatish Balay   PetscFunctionReturn(0);
172849b5e25fSSatish Balay }
172949b5e25fSSatish Balay 
17304a2ae208SSatish Balay #undef __FUNCT__
17314a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqSBAIJ"
1732dfbe8321SBarry Smith PetscErrorCode MatLoad_SeqSBAIJ(PetscViewer viewer,const MatType type,Mat *A)
173349b5e25fSSatish Balay {
173449b5e25fSSatish Balay   Mat_SeqSBAIJ   *a;
173549b5e25fSSatish Balay   Mat            B;
17366849ba73SBarry Smith   PetscErrorCode ierr;
173713f74950SBarry Smith   int            fd;
173813f74950SBarry Smith   PetscMPIInt    size;
173913f74950SBarry Smith   PetscInt       i,nz,header[4],*rowlengths=0,M,N,bs=1;
174013f74950SBarry Smith   PetscInt       *mask,mbs,*jj,j,rowcount,nzcount,k,*s_browlengths,maskcount;
174113f74950SBarry Smith   PetscInt       kmax,jcount,block,idx,point,nzcountb,extra_rows;
174213f74950SBarry Smith   PetscInt       *masked,nmask,tmp,bs2,ishift;
174387828ca2SBarry Smith   PetscScalar    *aa;
174449b5e25fSSatish Balay   MPI_Comm       comm = ((PetscObject)viewer)->comm;
174549b5e25fSSatish Balay 
174649b5e25fSSatish Balay   PetscFunctionBegin;
1747b0a32e0cSBarry Smith   ierr = PetscOptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,PETSC_NULL);CHKERRQ(ierr);
174849b5e25fSSatish Balay   bs2  = bs*bs;
174949b5e25fSSatish Balay 
175049b5e25fSSatish Balay   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
175129bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor");
1752b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
175349b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr);
1754552e946dSBarry Smith   if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not Mat object");
175549b5e25fSSatish Balay   M = header[1]; N = header[2]; nz = header[3];
175649b5e25fSSatish Balay 
175749b5e25fSSatish Balay   if (header[3] < 0) {
175829bbc08cSBarry Smith     SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as SeqSBAIJ");
175949b5e25fSSatish Balay   }
176049b5e25fSSatish Balay 
176129bbc08cSBarry Smith   if (M != N) SETERRQ(PETSC_ERR_SUP,"Can only do square matrices");
176249b5e25fSSatish Balay 
176349b5e25fSSatish Balay   /*
176449b5e25fSSatish Balay      This code adds extra rows to make sure the number of rows is
176549b5e25fSSatish Balay     divisible by the blocksize
176649b5e25fSSatish Balay   */
176749b5e25fSSatish Balay   mbs        = M/bs;
176849b5e25fSSatish Balay   extra_rows = bs - M + bs*(mbs);
176949b5e25fSSatish Balay   if (extra_rows == bs) extra_rows = 0;
177049b5e25fSSatish Balay   else                  mbs++;
177149b5e25fSSatish Balay   if (extra_rows) {
177263ba0a88SBarry Smith     ierr = PetscLogInfo((0,"MatLoad_SeqSBAIJ:Padding loaded matrix to match blocksize\n"));CHKERRQ(ierr);
177349b5e25fSSatish Balay   }
177449b5e25fSSatish Balay 
177549b5e25fSSatish Balay   /* read in row lengths */
177613f74950SBarry Smith   ierr = PetscMalloc((M+extra_rows)*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr);
177749b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
177849b5e25fSSatish Balay   for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1;
177949b5e25fSSatish Balay 
178049b5e25fSSatish Balay   /* read in column indices */
178113f74950SBarry Smith   ierr = PetscMalloc((nz+extra_rows)*sizeof(PetscInt),&jj);CHKERRQ(ierr);
178249b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT);CHKERRQ(ierr);
178349b5e25fSSatish Balay   for (i=0; i<extra_rows; i++) jj[nz+i] = M+i;
178449b5e25fSSatish Balay 
178549b5e25fSSatish Balay   /* loop over row lengths determining block row lengths */
178613f74950SBarry Smith   ierr     = PetscMalloc(mbs*sizeof(PetscInt),&s_browlengths);CHKERRQ(ierr);
178713f74950SBarry Smith   ierr     = PetscMemzero(s_browlengths,mbs*sizeof(PetscInt));CHKERRQ(ierr);
178813f74950SBarry Smith   ierr     = PetscMalloc(2*mbs*sizeof(PetscInt),&mask);CHKERRQ(ierr);
178913f74950SBarry Smith   ierr     = PetscMemzero(mask,mbs*sizeof(PetscInt));CHKERRQ(ierr);
179049b5e25fSSatish Balay   masked   = mask + mbs;
179149b5e25fSSatish Balay   rowcount = 0; nzcount = 0;
179249b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
179349b5e25fSSatish Balay     nmask = 0;
179449b5e25fSSatish Balay     for (j=0; j<bs; j++) {
179549b5e25fSSatish Balay       kmax = rowlengths[rowcount];
179649b5e25fSSatish Balay       for (k=0; k<kmax; k++) {
17972d703238SHong Zhang         tmp = jj[nzcount++]/bs;   /* block col. index */
179803630b6eSHong Zhang         if (!mask[tmp] && tmp >= i) {masked[nmask++] = tmp; mask[tmp] = 1;}
179949b5e25fSSatish Balay       }
180049b5e25fSSatish Balay       rowcount++;
180149b5e25fSSatish Balay     }
1802574b2666SHong Zhang     s_browlengths[i] += nmask;
1803574b2666SHong Zhang 
180449b5e25fSSatish Balay     /* zero out the mask elements we set */
180549b5e25fSSatish Balay     for (j=0; j<nmask; j++) mask[masked[j]] = 0;
180649b5e25fSSatish Balay   }
180749b5e25fSSatish Balay 
180849b5e25fSSatish Balay   /* create our matrix */
18099abb65ffSKris Buschelman   ierr = MatCreate(comm,M+extra_rows,N+extra_rows,M+extra_rows,N+extra_rows,&B);CHKERRQ(ierr);
18109abb65ffSKris Buschelman   ierr = MatSetType(B,type);CHKERRQ(ierr);
1811ab93d7beSBarry Smith   ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(B,bs,0,s_browlengths);CHKERRQ(ierr);
181249b5e25fSSatish Balay   a = (Mat_SeqSBAIJ*)B->data;
181349b5e25fSSatish Balay 
181449b5e25fSSatish Balay   /* set matrix "i" values */
181549b5e25fSSatish Balay   a->i[0] = 0;
181649b5e25fSSatish Balay   for (i=1; i<= mbs; i++) {
1817574b2666SHong Zhang     a->i[i]      = a->i[i-1] + s_browlengths[i-1];
1818574b2666SHong Zhang     a->ilen[i-1] = s_browlengths[i-1];
181949b5e25fSSatish Balay   }
18206c6c5352SBarry Smith   a->nz = a->i[mbs];
182149b5e25fSSatish Balay 
182249b5e25fSSatish Balay   /* read in nonzero values */
182387828ca2SBarry Smith   ierr = PetscMalloc((nz+extra_rows)*sizeof(PetscScalar),&aa);CHKERRQ(ierr);
182449b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);CHKERRQ(ierr);
182549b5e25fSSatish Balay   for (i=0; i<extra_rows; i++) aa[nz+i] = 1.0;
182649b5e25fSSatish Balay 
182749b5e25fSSatish Balay   /* set "a" and "j" values into matrix */
182849b5e25fSSatish Balay   nzcount = 0; jcount = 0;
182949b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
183049b5e25fSSatish Balay     nzcountb = nzcount;
183149b5e25fSSatish Balay     nmask    = 0;
183249b5e25fSSatish Balay     for (j=0; j<bs; j++) {
183349b5e25fSSatish Balay       kmax = rowlengths[i*bs+j];
183449b5e25fSSatish Balay       for (k=0; k<kmax; k++) {
18352d703238SHong Zhang         tmp = jj[nzcount++]/bs; /* block col. index */
183603630b6eSHong Zhang         if (!mask[tmp] && tmp >= i) { masked[nmask++] = tmp; mask[tmp] = 1;}
18372d703238SHong Zhang       }
18382d703238SHong Zhang     }
18392d703238SHong Zhang     /* sort the masked values */
18402d703238SHong Zhang     ierr = PetscSortInt(nmask,masked);CHKERRQ(ierr);
18412d703238SHong Zhang 
18422d703238SHong Zhang     /* set "j" values into matrix */
18432d703238SHong Zhang     maskcount = 1;
18442d703238SHong Zhang     for (j=0; j<nmask; j++) {
184549b5e25fSSatish Balay       a->j[jcount++]  = masked[j];
184649b5e25fSSatish Balay       mask[masked[j]] = maskcount++;
184749b5e25fSSatish Balay     }
1848574b2666SHong Zhang 
184949b5e25fSSatish Balay     /* set "a" values into matrix */
185049b5e25fSSatish Balay     ishift = bs2*a->i[i];
185149b5e25fSSatish Balay     for (j=0; j<bs; j++) {
185249b5e25fSSatish Balay       kmax = rowlengths[i*bs+j];
185349b5e25fSSatish Balay       for (k=0; k<kmax; k++) {
1854574b2666SHong Zhang         tmp       = jj[nzcountb]/bs ; /* block col. index */
1855574b2666SHong Zhang         if (tmp >= i){
185649b5e25fSSatish Balay           block     = mask[tmp] - 1;
185749b5e25fSSatish Balay           point     = jj[nzcountb] - bs*tmp;
185849b5e25fSSatish Balay           idx       = ishift + bs2*block + j + bs*point;
1859574b2666SHong Zhang           a->a[idx] = aa[nzcountb];
1860574b2666SHong Zhang         }
1861574b2666SHong Zhang         nzcountb++;
186249b5e25fSSatish Balay       }
186349b5e25fSSatish Balay     }
186449b5e25fSSatish Balay     /* zero out the mask elements we set */
186549b5e25fSSatish Balay     for (j=0; j<nmask; j++) mask[masked[j]] = 0;
186649b5e25fSSatish Balay   }
18676c6c5352SBarry Smith   if (jcount != a->nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Bad binary matrix");
186849b5e25fSSatish Balay 
186949b5e25fSSatish Balay   ierr = PetscFree(rowlengths);CHKERRQ(ierr);
1870574b2666SHong Zhang   ierr = PetscFree(s_browlengths);CHKERRQ(ierr);
187149b5e25fSSatish Balay   ierr = PetscFree(aa);CHKERRQ(ierr);
187249b5e25fSSatish Balay   ierr = PetscFree(jj);CHKERRQ(ierr);
187349b5e25fSSatish Balay   ierr = PetscFree(mask);CHKERRQ(ierr);
187449b5e25fSSatish Balay 
18759abb65ffSKris Buschelman   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
18769abb65ffSKris Buschelman   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
187749b5e25fSSatish Balay   ierr = MatView_Private(B);CHKERRQ(ierr);
18789abb65ffSKris Buschelman   *A = B;
187949b5e25fSSatish Balay   PetscFunctionReturn(0);
188049b5e25fSSatish Balay }
1881574b2666SHong Zhang 
1882d06b337dSHong Zhang #undef __FUNCT__
1883d06b337dSHong Zhang #define __FUNCT__ "MatRelax_SeqSBAIJ"
188413f74950SBarry Smith PetscErrorCode MatRelax_SeqSBAIJ(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx)
1885d06b337dSHong Zhang {
1886d06b337dSHong Zhang   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
1887d06b337dSHong Zhang   MatScalar      *aa=a->a,*v,*v1;
1888d06b337dSHong Zhang   PetscScalar    *x,*b,*t,sum,d;
18896849ba73SBarry Smith   PetscErrorCode ierr;
1890521d7252SBarry Smith   PetscInt       m=a->mbs,bs=A->bs,*ai=a->i,*aj=a->j;
1891521d7252SBarry Smith   PetscInt       nz,nz1,*vj,*vj1,i;
1892d06b337dSHong Zhang 
1893d06b337dSHong Zhang   PetscFunctionBegin;
1894b965ef7fSBarry Smith   its = its*lits;
189577431f27SBarry Smith   if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits);
1896d06b337dSHong Zhang 
1897d06b337dSHong Zhang   if (bs > 1)
1898d06b337dSHong Zhang     SETERRQ(PETSC_ERR_SUP,"SSOR for block size > 1 is not yet implemented");
1899d06b337dSHong Zhang 
19001ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1901d06b337dSHong Zhang   if (xx != bb) {
19021ebc52fbSHong Zhang     ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
1903d06b337dSHong Zhang   } else {
1904d06b337dSHong Zhang     b = x;
1905d06b337dSHong Zhang   }
1906d06b337dSHong Zhang 
1907d06b337dSHong Zhang   ierr = PetscMalloc(m*sizeof(PetscScalar),&t);CHKERRQ(ierr);
1908d06b337dSHong Zhang 
1909d06b337dSHong Zhang   if (flag & SOR_ZERO_INITIAL_GUESS) {
19102798e883SHong Zhang     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
1911d06b337dSHong Zhang       for (i=0; i<m; i++)
1912d06b337dSHong Zhang         t[i] = b[i];
1913d06b337dSHong Zhang 
1914d06b337dSHong Zhang       for (i=0; i<m; i++){
191544706875SHong Zhang         d  = *(aa + ai[i]);  /* diag[i] */
1916d06b337dSHong Zhang         v  = aa + ai[i] + 1;
1917d06b337dSHong Zhang         vj = aj + ai[i] + 1;
1918d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1;
1919d06b337dSHong Zhang         x[i] = omega*t[i]/d;
1920d06b337dSHong Zhang         while (nz--) t[*vj++] -= x[i]*(*v++); /* update rhs */
1921efee365bSSatish Balay         ierr = PetscLogFlops(2*nz-1);CHKERRQ(ierr);
1922d06b337dSHong Zhang       }
1923d06b337dSHong Zhang     }
1924d06b337dSHong Zhang 
19252798e883SHong Zhang     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
1926d06b337dSHong Zhang       for (i=0; i<m; i++)
1927d06b337dSHong Zhang         t[i] = b[i];
1928d06b337dSHong Zhang 
1929d06b337dSHong Zhang       for (i=0; i<m-1; i++){  /* update rhs */
1930d06b337dSHong Zhang         v  = aa + ai[i] + 1;
1931d06b337dSHong Zhang         vj = aj + ai[i] + 1;
1932d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1;
1933d06b337dSHong Zhang         while (nz--) t[*vj++] -= x[i]*(*v++);
1934efee365bSSatish Balay         ierr = PetscLogFlops(2*nz-1);CHKERRQ(ierr);
1935d06b337dSHong Zhang       }
1936d06b337dSHong Zhang       for (i=m-1; i>=0; i--){
1937d06b337dSHong Zhang         d  = *(aa + ai[i]);
1938d06b337dSHong Zhang         v  = aa + ai[i] + 1;
1939d06b337dSHong Zhang         vj = aj + ai[i] + 1;
1940d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1;
1941d06b337dSHong Zhang         sum = t[i];
1942d06b337dSHong Zhang         while (nz--) sum -= x[*vj++]*(*v++);
1943efee365bSSatish Balay         ierr = PetscLogFlops(2*nz-1);CHKERRQ(ierr);
1944d06b337dSHong Zhang         x[i] =   (1-omega)*x[i] + omega*sum/d;
1945d06b337dSHong Zhang       }
1946d06b337dSHong Zhang     }
1947d06b337dSHong Zhang     its--;
1948d06b337dSHong Zhang   }
1949d06b337dSHong Zhang 
1950d06b337dSHong Zhang   while (its--) {
195144706875SHong Zhang     /*
195244706875SHong Zhang        forward sweep:
195344706875SHong Zhang        for i=0,...,m-1:
195444706875SHong Zhang          sum[i] = (b[i] - U(i,:)x )/d[i];
195544706875SHong Zhang          x[i]   = (1-omega)x[i] + omega*sum[i];
195644706875SHong Zhang          b      = b - x[i]*U^T(i,:);
1957d06b337dSHong Zhang 
195844706875SHong Zhang     */
19592798e883SHong Zhang     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
1960d06b337dSHong Zhang       for (i=0; i<m; i++)
1961d06b337dSHong Zhang         t[i] = b[i];
1962d06b337dSHong Zhang 
1963d06b337dSHong Zhang       for (i=0; i<m; i++){
196444706875SHong Zhang         d  = *(aa + ai[i]);  /* diag[i] */
1965d06b337dSHong Zhang         v  = aa + ai[i] + 1; v1=v;
1966d06b337dSHong Zhang         vj = aj + ai[i] + 1; vj1=vj;
1967d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1; nz1=nz;
1968d06b337dSHong Zhang         sum = t[i];
1969d06b337dSHong Zhang         while (nz1--) sum -= (*v1++)*x[*vj1++];
1970d06b337dSHong Zhang         x[i] = (1-omega)*x[i] + omega*sum/d;
1971d06b337dSHong Zhang         while (nz--) t[*vj++] -= x[i]*(*v++);
1972efee365bSSatish Balay         ierr = PetscLogFlops(4*nz-2);CHKERRQ(ierr);
1973d06b337dSHong Zhang       }
1974d06b337dSHong Zhang     }
1975d06b337dSHong Zhang 
19762798e883SHong Zhang   if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
197744706875SHong Zhang       /*
197844706875SHong Zhang        backward sweep:
197944706875SHong Zhang        b = b - x[i]*U^T(i,:), i=0,...,n-2
198044706875SHong Zhang        for i=m-1,...,0:
198144706875SHong Zhang          sum[i] = (b[i] - U(i,:)x )/d[i];
198244706875SHong Zhang          x[i]   = (1-omega)x[i] + omega*sum[i];
198344706875SHong Zhang       */
1984d06b337dSHong Zhang       for (i=0; i<m; i++)
1985d06b337dSHong Zhang         t[i] = b[i];
1986d06b337dSHong Zhang 
1987d06b337dSHong Zhang       for (i=0; i<m-1; i++){  /* update rhs */
1988d06b337dSHong Zhang         v  = aa + ai[i] + 1;
1989d06b337dSHong Zhang         vj = aj + ai[i] + 1;
1990d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1;
1991d06b337dSHong Zhang         while (nz--) t[*vj++] -= x[i]*(*v++);
1992efee365bSSatish Balay         ierr = PetscLogFlops(2*nz-1);CHKERRQ(ierr);
1993d06b337dSHong Zhang       }
1994d06b337dSHong Zhang       for (i=m-1; i>=0; i--){
1995d06b337dSHong Zhang         d  = *(aa + ai[i]);
1996d06b337dSHong Zhang         v  = aa + ai[i] + 1;
1997d06b337dSHong Zhang         vj = aj + ai[i] + 1;
1998d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1;
1999d06b337dSHong Zhang         sum = t[i];
2000d06b337dSHong Zhang         while (nz--) sum -= x[*vj++]*(*v++);
2001efee365bSSatish Balay         ierr = PetscLogFlops(2*nz-1);CHKERRQ(ierr);
2002d06b337dSHong Zhang         x[i] =   (1-omega)*x[i] + omega*sum/d;
2003d06b337dSHong Zhang       }
2004d06b337dSHong Zhang     }
2005d06b337dSHong Zhang   }
2006d06b337dSHong Zhang 
2007d06b337dSHong Zhang   ierr = PetscFree(t);CHKERRQ(ierr);
20081ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2009d06b337dSHong Zhang   if (bb != xx) {
20101ebc52fbSHong Zhang     ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
2011d06b337dSHong Zhang   }
2012d06b337dSHong Zhang   PetscFunctionReturn(0);
2013d06b337dSHong Zhang }
2014d06b337dSHong Zhang 
2015d06b337dSHong Zhang 
2016d06b337dSHong Zhang 
2017d06b337dSHong Zhang 
201849b5e25fSSatish Balay 
201949b5e25fSSatish Balay 
2020