xref: /petsc/src/mat/impls/sbaij/seq/sbaij.c (revision e11b0e149a41efd678d96ff732835a7007915d7d)
173f4d377SMatthew Knepley /*$Id: sbaij.c,v 1.62 2001/08/07 03:03:01 balay Exp $*/
249b5e25fSSatish Balay 
349b5e25fSSatish Balay /*
449b5e25fSSatish Balay     Defines the basic matrix operations for the BAIJ (compressed row)
549b5e25fSSatish Balay   matrix storage format.
649b5e25fSSatish Balay */
79e070d67SMatthew Knepley #include "src/mat/impls/baij/seq/baij.h"         /*I "petscmat.h" I*/
849b5e25fSSatish Balay #include "src/vec/vecimpl.h"
949b5e25fSSatish Balay #include "src/inline/spops.h"
10aec82049SSatish Balay #include "src/mat/impls/sbaij/seq/sbaij.h"
1149b5e25fSSatish Balay 
1249b5e25fSSatish Balay #define CHUNKSIZE  10
1349b5e25fSSatish Balay 
1449b5e25fSSatish Balay /*
1549b5e25fSSatish Balay      Checks for missing diagonals
1649b5e25fSSatish Balay */
174a2ae208SSatish Balay #undef __FUNCT__
184a2ae208SSatish Balay #define __FUNCT__ "MatMissingDiagonal_SeqSBAIJ"
1949b5e25fSSatish Balay int MatMissingDiagonal_SeqSBAIJ(Mat A)
2049b5e25fSSatish Balay {
21045c9aa0SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
2249b5e25fSSatish Balay   int         *diag,*jj = a->j,i,ierr;
2349b5e25fSSatish Balay 
2449b5e25fSSatish Balay   PetscFunctionBegin;
25045c9aa0SHong Zhang   ierr = MatMarkDiagonal_SeqSBAIJ(A);CHKERRQ(ierr);
2649b5e25fSSatish Balay   diag = a->diag;
2749b5e25fSSatish Balay   for (i=0; i<a->mbs; i++) {
28*e11b0e14SHong Zhang     if (jj[diag[i]] != i) SETERRQ1(1,"Matrix is missing diagonal number %d",i);
2949b5e25fSSatish Balay   }
3049b5e25fSSatish Balay   PetscFunctionReturn(0);
3149b5e25fSSatish Balay }
3249b5e25fSSatish Balay 
334a2ae208SSatish Balay #undef __FUNCT__
344a2ae208SSatish Balay #define __FUNCT__ "MatMarkDiagonal_SeqSBAIJ"
3549b5e25fSSatish Balay int MatMarkDiagonal_SeqSBAIJ(Mat A)
3649b5e25fSSatish Balay {
37045c9aa0SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
38*e11b0e14SHong Zhang   int          i,m = a->mbs,ierr;
3949b5e25fSSatish Balay 
4049b5e25fSSatish Balay   PetscFunctionBegin;
4149b5e25fSSatish Balay   if (a->diag) PetscFunctionReturn(0);
4249b5e25fSSatish Balay 
43*e11b0e14SHong Zhang   ierr = PetscMalloc((m+1)*sizeof(int),&a->diag);CHKERRQ(ierr);
44b0a32e0cSBarry Smith   PetscLogObjectMemory(A,(m+1)*sizeof(int));
45*e11b0e14SHong Zhang   for (i=0; i<m; i++) a->diag[i] = a->i[i];
4649b5e25fSSatish Balay   PetscFunctionReturn(0);
4749b5e25fSSatish Balay }
4849b5e25fSSatish Balay 
4949b5e25fSSatish Balay 
5049b5e25fSSatish Balay extern int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int,int**,int**);
5149b5e25fSSatish Balay 
524a2ae208SSatish Balay #undef __FUNCT__
534a2ae208SSatish Balay #define __FUNCT__ "MatGetRowIJ_SeqSBAIJ"
5449b5e25fSSatish Balay static int MatGetRowIJ_SeqSBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja,PetscTruth *done)
5549b5e25fSSatish Balay {
5649b5e25fSSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
5749b5e25fSSatish Balay 
5849b5e25fSSatish Balay   PetscFunctionBegin;
59045c9aa0SHong Zhang   if (ia) {
6029bbc08cSBarry Smith     SETERRQ(1,"Function not yet written for SBAIJ format, only supports natural ordering");
6149b5e25fSSatish Balay   }
62045c9aa0SHong Zhang   *nn = a->mbs;
6349b5e25fSSatish Balay   PetscFunctionReturn(0);
6449b5e25fSSatish Balay }
6549b5e25fSSatish Balay 
664a2ae208SSatish Balay #undef __FUNCT__
674a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRowIJ_SeqSBAIJ"
68045c9aa0SHong Zhang static int MatRestoreRowIJ_SeqSBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja,PetscTruth *done)
6949b5e25fSSatish Balay {
7049b5e25fSSatish Balay   PetscFunctionBegin;
7149b5e25fSSatish Balay   if (!ia) PetscFunctionReturn(0);
7229bbc08cSBarry Smith   SETERRQ(1,"Function not yet written for SBAIJ format, only supports natural ordering");
7316cdd363SHong Zhang   /* PetscFunctionReturn(0); */
7449b5e25fSSatish Balay }
7549b5e25fSSatish Balay 
764a2ae208SSatish Balay #undef __FUNCT__
774a2ae208SSatish Balay #define __FUNCT__ "MatGetBlockSize_SeqSBAIJ"
7849b5e25fSSatish Balay int MatGetBlockSize_SeqSBAIJ(Mat mat,int *bs)
7949b5e25fSSatish Balay {
8049b5e25fSSatish Balay   Mat_SeqSBAIJ *sbaij = (Mat_SeqSBAIJ*)mat->data;
8149b5e25fSSatish Balay 
8249b5e25fSSatish Balay   PetscFunctionBegin;
8349b5e25fSSatish Balay   *bs = sbaij->bs;
8449b5e25fSSatish Balay   PetscFunctionReturn(0);
8549b5e25fSSatish Balay }
8649b5e25fSSatish Balay 
874a2ae208SSatish Balay #undef __FUNCT__
884a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqSBAIJ"
8949b5e25fSSatish Balay int MatDestroy_SeqSBAIJ(Mat A)
9049b5e25fSSatish Balay {
9149b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
9249b5e25fSSatish Balay   int         ierr;
9349b5e25fSSatish Balay 
9449b5e25fSSatish Balay   PetscFunctionBegin;
9549b5e25fSSatish Balay #if defined(PETSC_USE_LOG)
96b0a32e0cSBarry Smith   PetscLogObjectState((PetscObject)A,"Rows=%d, s_NZ=%d",A->m,a->s_nz);
9749b5e25fSSatish Balay #endif
9849b5e25fSSatish Balay   ierr = PetscFree(a->a);CHKERRQ(ierr);
9949b5e25fSSatish Balay   if (!a->singlemalloc) {
10049b5e25fSSatish Balay     ierr = PetscFree(a->i);CHKERRQ(ierr);
10163c8bf9fSHong Zhang     ierr = PetscFree(a->j);CHKERRQ(ierr);
10249b5e25fSSatish Balay   }
10349b5e25fSSatish Balay   if (a->row) {
10449b5e25fSSatish Balay     ierr = ISDestroy(a->row);CHKERRQ(ierr);
10549b5e25fSSatish Balay   }
10649b5e25fSSatish Balay   if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);}
10749b5e25fSSatish Balay   if (a->ilen) {ierr = PetscFree(a->ilen);CHKERRQ(ierr);}
10849b5e25fSSatish Balay   if (a->imax) {ierr = PetscFree(a->imax);CHKERRQ(ierr);}
10949b5e25fSSatish Balay   if (a->solve_work) {ierr = PetscFree(a->solve_work);CHKERRQ(ierr);}
11049b5e25fSSatish Balay   if (a->mult_work) {ierr = PetscFree(a->mult_work);CHKERRQ(ierr);}
11149b5e25fSSatish Balay   if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);}
11249b5e25fSSatish Balay   if (a->saved_values) {ierr = PetscFree(a->saved_values);CHKERRQ(ierr);}
1131a3463dfSHong Zhang 
1141a3463dfSHong Zhang   if (a->inew){
1151a3463dfSHong Zhang     ierr = PetscFree(a->inew);CHKERRQ(ierr);
1161a3463dfSHong Zhang     a->inew = 0;
1171a3463dfSHong Zhang   }
11849b5e25fSSatish Balay   ierr = PetscFree(a);CHKERRQ(ierr);
11949b5e25fSSatish Balay   PetscFunctionReturn(0);
12049b5e25fSSatish Balay }
12149b5e25fSSatish Balay 
1224a2ae208SSatish Balay #undef __FUNCT__
1234a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqSBAIJ"
12449b5e25fSSatish Balay int MatSetOption_SeqSBAIJ(Mat A,MatOption op)
12549b5e25fSSatish Balay {
126045c9aa0SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
12749b5e25fSSatish Balay 
12849b5e25fSSatish Balay   PetscFunctionBegin;
1294d9d31abSKris Buschelman   switch (op) {
1304d9d31abSKris Buschelman   case MAT_ROW_ORIENTED:
1314d9d31abSKris Buschelman     a->roworiented    = PETSC_TRUE;
1324d9d31abSKris Buschelman     break;
1334d9d31abSKris Buschelman   case MAT_COLUMN_ORIENTED:
1344d9d31abSKris Buschelman     a->roworiented    = PETSC_FALSE;
1354d9d31abSKris Buschelman     break;
1364d9d31abSKris Buschelman   case MAT_COLUMNS_SORTED:
1374d9d31abSKris Buschelman     a->sorted         = PETSC_TRUE;
1384d9d31abSKris Buschelman     break;
1394d9d31abSKris Buschelman   case MAT_COLUMNS_UNSORTED:
1404d9d31abSKris Buschelman     a->sorted         = PETSC_FALSE;
1414d9d31abSKris Buschelman     break;
1424d9d31abSKris Buschelman   case MAT_KEEP_ZEROED_ROWS:
1434d9d31abSKris Buschelman     a->keepzeroedrows = PETSC_TRUE;
1444d9d31abSKris Buschelman     break;
1454d9d31abSKris Buschelman   case MAT_NO_NEW_NONZERO_LOCATIONS:
1464d9d31abSKris Buschelman     a->nonew          = 1;
1474d9d31abSKris Buschelman     break;
1484d9d31abSKris Buschelman   case MAT_NEW_NONZERO_LOCATION_ERR:
1494d9d31abSKris Buschelman     a->nonew          = -1;
1504d9d31abSKris Buschelman     break;
1514d9d31abSKris Buschelman   case MAT_NEW_NONZERO_ALLOCATION_ERR:
1524d9d31abSKris Buschelman     a->nonew          = -2;
1534d9d31abSKris Buschelman     break;
1544d9d31abSKris Buschelman   case MAT_YES_NEW_NONZERO_LOCATIONS:
1554d9d31abSKris Buschelman     a->nonew          = 0;
1564d9d31abSKris Buschelman     break;
1574d9d31abSKris Buschelman   case MAT_ROWS_SORTED:
1584d9d31abSKris Buschelman   case MAT_ROWS_UNSORTED:
1594d9d31abSKris Buschelman   case MAT_YES_NEW_DIAGONALS:
1604d9d31abSKris Buschelman   case MAT_IGNORE_OFF_PROC_ENTRIES:
1614d9d31abSKris Buschelman   case MAT_USE_HASH_TABLE:
162d03495bdSKris Buschelman   case MAT_USE_SINGLE_PRECISION_SOLVES:
163b0a32e0cSBarry Smith     PetscLogInfo(A,"MatSetOption_SeqSBAIJ:Option ignored\n");
1644d9d31abSKris Buschelman     break;
1654d9d31abSKris Buschelman   case MAT_NO_NEW_DIAGONALS:
16629bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS");
1674d9d31abSKris Buschelman   default:
16829bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"unknown option");
16949b5e25fSSatish Balay   }
17049b5e25fSSatish Balay   PetscFunctionReturn(0);
17149b5e25fSSatish Balay }
17249b5e25fSSatish Balay 
1734a2ae208SSatish Balay #undef __FUNCT__
1744a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqSBAIJ"
17587828ca2SBarry Smith int MatGetRow_SeqSBAIJ(Mat A,int row,int *ncols,int **cols,PetscScalar **v)
17649b5e25fSSatish Balay {
17749b5e25fSSatish Balay   Mat_SeqSBAIJ  *a = (Mat_SeqSBAIJ*)A->data;
17882502324SSatish Balay   int          itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*cols_i,bs2,ierr;
17949b5e25fSSatish Balay   MatScalar    *aa,*aa_i;
18087828ca2SBarry Smith   PetscScalar  *v_i;
18149b5e25fSSatish Balay 
18249b5e25fSSatish Balay   PetscFunctionBegin;
18349b5e25fSSatish Balay   bs  = a->bs;
18449b5e25fSSatish Balay   ai  = a->i;
18549b5e25fSSatish Balay   aj  = a->j;
18649b5e25fSSatish Balay   aa  = a->a;
18749b5e25fSSatish Balay   bs2 = a->bs2;
18849b5e25fSSatish Balay 
189c464158bSHong Zhang   if (row < 0 || row >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row out of range");
19049b5e25fSSatish Balay 
19149b5e25fSSatish Balay   bn  = row/bs;   /* Block number */
19249b5e25fSSatish Balay   bp  = row % bs; /* Block position */
19349b5e25fSSatish Balay   M   = ai[bn+1] - ai[bn];
19449b5e25fSSatish Balay   *ncols = bs*M;
19549b5e25fSSatish Balay 
19649b5e25fSSatish Balay   if (v) {
19749b5e25fSSatish Balay     *v = 0;
19849b5e25fSSatish Balay     if (*ncols) {
19987828ca2SBarry Smith       ierr = PetscMalloc((*ncols+row)*sizeof(PetscScalar),v);CHKERRQ(ierr);
20049b5e25fSSatish Balay       for (i=0; i<M; i++) { /* for each block in the block row */
20149b5e25fSSatish Balay         v_i  = *v + i*bs;
20249b5e25fSSatish Balay         aa_i = aa + bs2*(ai[bn] + i);
20349b5e25fSSatish Balay         for (j=bp,k=0; j<bs2; j+=bs,k++) {v_i[k] = aa_i[j];}
20449b5e25fSSatish Balay       }
20549b5e25fSSatish Balay     }
20649b5e25fSSatish Balay   }
20749b5e25fSSatish Balay 
20849b5e25fSSatish Balay   if (cols) {
20949b5e25fSSatish Balay     *cols = 0;
21049b5e25fSSatish Balay     if (*ncols) {
21182502324SSatish Balay       ierr = PetscMalloc((*ncols+row)*sizeof(int),cols);CHKERRQ(ierr);
21249b5e25fSSatish Balay       for (i=0; i<M; i++) { /* for each block in the block row */
21349b5e25fSSatish Balay         cols_i = *cols + i*bs;
21449b5e25fSSatish Balay         itmp  = bs*aj[ai[bn] + i];
21549b5e25fSSatish Balay         for (j=0; j<bs; j++) {cols_i[j] = itmp++;}
21649b5e25fSSatish Balay       }
21749b5e25fSSatish Balay     }
21849b5e25fSSatish Balay   }
21949b5e25fSSatish Balay 
22049b5e25fSSatish Balay   /*search column A(0:row-1,row) (=A(row,0:row-1)). Could be expensive! */
2215ddb2528SHong Zhang   /* this segment is currently removed, so only entries in the upper triangle are obtained */
2225ddb2528SHong Zhang #ifdef column_search
22349b5e25fSSatish Balay   v_i    = *v    + M*bs;
22449b5e25fSSatish Balay   cols_i = *cols + M*bs;
22549b5e25fSSatish Balay   for (i=0; i<bn; i++){ /* for each block row */
22649b5e25fSSatish Balay     M = ai[i+1] - ai[i];
22749b5e25fSSatish Balay     for (j=0; j<M; j++){
22849b5e25fSSatish Balay       itmp = aj[ai[i] + j];    /* block column value */
22949b5e25fSSatish Balay       if (itmp == bn){
23049b5e25fSSatish Balay         aa_i   = aa    + bs2*(ai[i] + j) + bs*bp;
23149b5e25fSSatish Balay         for (k=0; k<bs; k++) {
23249b5e25fSSatish Balay           *cols_i++ = i*bs+k;
23349b5e25fSSatish Balay           *v_i++    = aa_i[k];
23449b5e25fSSatish Balay         }
23549b5e25fSSatish Balay         *ncols += bs;
23649b5e25fSSatish Balay         break;
23749b5e25fSSatish Balay       }
23849b5e25fSSatish Balay     }
23949b5e25fSSatish Balay   }
2405ddb2528SHong Zhang #endif
24149b5e25fSSatish Balay 
24249b5e25fSSatish Balay   PetscFunctionReturn(0);
24349b5e25fSSatish Balay }
24449b5e25fSSatish Balay 
2454a2ae208SSatish Balay #undef __FUNCT__
2464a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqSBAIJ"
24787828ca2SBarry Smith int MatRestoreRow_SeqSBAIJ(Mat A,int row,int *nz,int **idx,PetscScalar **v)
24849b5e25fSSatish Balay {
24949b5e25fSSatish Balay   int ierr;
25049b5e25fSSatish Balay 
25149b5e25fSSatish Balay   PetscFunctionBegin;
25249b5e25fSSatish Balay   if (idx) {if (*idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}}
25349b5e25fSSatish Balay   if (v)   {if (*v)   {ierr = PetscFree(*v);CHKERRQ(ierr);}}
25449b5e25fSSatish Balay   PetscFunctionReturn(0);
25549b5e25fSSatish Balay }
25649b5e25fSSatish Balay 
2574a2ae208SSatish Balay #undef __FUNCT__
2584a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqSBAIJ"
25949b5e25fSSatish Balay int MatTranspose_SeqSBAIJ(Mat A,Mat *B)
26049b5e25fSSatish Balay {
26149b5e25fSSatish Balay   PetscFunctionBegin;
26229bbc08cSBarry Smith   SETERRQ(1,"Matrix is symmetric. MatTranspose() should not be called");
26316cdd363SHong Zhang   /* PetscFunctionReturn(0); */
26449b5e25fSSatish Balay }
26549b5e25fSSatish Balay 
2664a2ae208SSatish Balay #undef __FUNCT__
2674a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Binary"
268b0a32e0cSBarry Smith static int MatView_SeqSBAIJ_Binary(Mat A,PetscViewer viewer)
26949b5e25fSSatish Balay {
27049b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
27149b5e25fSSatish Balay   int          i,fd,*col_lens,ierr,bs = a->bs,count,*jj,j,k,l,bs2=a->bs2;
27287828ca2SBarry Smith   PetscScalar  *aa;
27349b5e25fSSatish Balay   FILE         *file;
27449b5e25fSSatish Balay 
27549b5e25fSSatish Balay   PetscFunctionBegin;
276b0a32e0cSBarry Smith   ierr        = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
27782502324SSatish Balay   ierr        = PetscMalloc((4+A->m)*sizeof(int),&col_lens);CHKERRQ(ierr);
278552e946dSBarry Smith   col_lens[0] = MAT_FILE_COOKIE;
27949b5e25fSSatish Balay 
280c464158bSHong Zhang   col_lens[1] = A->m;
281c464158bSHong Zhang   col_lens[2] = A->m;
28249b5e25fSSatish Balay   col_lens[3] = a->s_nz*bs2;
28349b5e25fSSatish Balay 
28449b5e25fSSatish Balay   /* store lengths of each row and write (including header) to file */
28549b5e25fSSatish Balay   count = 0;
28649b5e25fSSatish Balay   for (i=0; i<a->mbs; i++) {
28749b5e25fSSatish Balay     for (j=0; j<bs; j++) {
28849b5e25fSSatish Balay       col_lens[4+count++] = bs*(a->i[i+1] - a->i[i]);
28949b5e25fSSatish Balay     }
29049b5e25fSSatish Balay   }
291c464158bSHong Zhang   ierr = PetscBinaryWrite(fd,col_lens,4+A->m,PETSC_INT,1);CHKERRQ(ierr);
29249b5e25fSSatish Balay   ierr = PetscFree(col_lens);CHKERRQ(ierr);
29349b5e25fSSatish Balay 
29449b5e25fSSatish Balay   /* store column indices (zero start index) */
29582502324SSatish Balay   ierr  = PetscMalloc((a->s_nz+1)*bs2*sizeof(int),&jj);CHKERRQ(ierr);
29649b5e25fSSatish Balay   count = 0;
29749b5e25fSSatish Balay   for (i=0; i<a->mbs; i++) {
29849b5e25fSSatish Balay     for (j=0; j<bs; j++) {
29949b5e25fSSatish Balay       for (k=a->i[i]; k<a->i[i+1]; k++) {
30049b5e25fSSatish Balay         for (l=0; l<bs; l++) {
30149b5e25fSSatish Balay           jj[count++] = bs*a->j[k] + l;
30249b5e25fSSatish Balay         }
30349b5e25fSSatish Balay       }
30449b5e25fSSatish Balay     }
30549b5e25fSSatish Balay   }
30649b5e25fSSatish Balay   ierr = PetscBinaryWrite(fd,jj,bs2*a->s_nz,PETSC_INT,0);CHKERRQ(ierr);
30749b5e25fSSatish Balay   ierr = PetscFree(jj);CHKERRQ(ierr);
30849b5e25fSSatish Balay 
30949b5e25fSSatish Balay   /* store nonzero values */
31087828ca2SBarry Smith   ierr  = PetscMalloc((a->s_nz+1)*bs2*sizeof(PetscScalar),&aa);CHKERRQ(ierr);
31149b5e25fSSatish Balay   count = 0;
31249b5e25fSSatish Balay   for (i=0; i<a->mbs; i++) {
31349b5e25fSSatish Balay     for (j=0; j<bs; j++) {
31449b5e25fSSatish Balay       for (k=a->i[i]; k<a->i[i+1]; k++) {
31549b5e25fSSatish Balay         for (l=0; l<bs; l++) {
31649b5e25fSSatish Balay           aa[count++] = a->a[bs2*k + l*bs + j];
31749b5e25fSSatish Balay         }
31849b5e25fSSatish Balay       }
31949b5e25fSSatish Balay     }
32049b5e25fSSatish Balay   }
32149b5e25fSSatish Balay   ierr = PetscBinaryWrite(fd,aa,bs2*a->s_nz,PETSC_SCALAR,0);CHKERRQ(ierr);
32249b5e25fSSatish Balay   ierr = PetscFree(aa);CHKERRQ(ierr);
32349b5e25fSSatish Balay 
324b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetInfoPointer(viewer,&file);CHKERRQ(ierr);
32549b5e25fSSatish Balay   if (file) {
32649b5e25fSSatish Balay     fprintf(file,"-matload_block_size %d\n",a->bs);
32749b5e25fSSatish Balay   }
32849b5e25fSSatish Balay   PetscFunctionReturn(0);
32949b5e25fSSatish Balay }
33049b5e25fSSatish Balay 
3314a2ae208SSatish Balay #undef __FUNCT__
3324a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_ASCII"
333b0a32e0cSBarry Smith static int MatView_SeqSBAIJ_ASCII(Mat A,PetscViewer viewer)
33449b5e25fSSatish Balay {
33549b5e25fSSatish Balay   Mat_SeqSBAIJ      *a = (Mat_SeqSBAIJ*)A->data;
336fb9695e5SSatish Balay   int               ierr,i,j,bs = a->bs,k,l,bs2=a->bs2;
337fb9695e5SSatish Balay   char              *name;
338f3ef73ceSBarry Smith   PetscViewerFormat format;
33949b5e25fSSatish Balay 
34049b5e25fSSatish Balay   PetscFunctionBegin;
34180fe4e49SBarry Smith   ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
342b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
343fb9695e5SSatish Balay   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_LONG) {
344b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  block size is %d\n",bs);CHKERRQ(ierr);
345fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_MATLAB) {
34629bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Matlab format not supported");
347fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_COMMON) {
348b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
34949b5e25fSSatish Balay     for (i=0; i<a->mbs; i++) {
35049b5e25fSSatish Balay       for (j=0; j<bs; j++) {
351b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr);
35249b5e25fSSatish Balay         for (k=a->i[i]; k<a->i[i+1]; k++) {
35349b5e25fSSatish Balay           for (l=0; l<bs; l++) {
35449b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX)
35549b5e25fSSatish Balay             if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
356b0a32e0cSBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," %d %g + %g i",bs*a->j[k]+l,
35749b5e25fSSatish Balay                       PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
35849b5e25fSSatish Balay             } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
359b0a32e0cSBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," %d %g - %g i",bs*a->j[k]+l,
36049b5e25fSSatish Balay                       PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
36149b5e25fSSatish Balay             } else if (PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
362b0a32e0cSBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
36349b5e25fSSatish Balay             }
36449b5e25fSSatish Balay #else
36549b5e25fSSatish Balay             if (a->a[bs2*k + l*bs + j] != 0.0) {
366b0a32e0cSBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr);
36749b5e25fSSatish Balay             }
36849b5e25fSSatish Balay #endif
36949b5e25fSSatish Balay           }
37049b5e25fSSatish Balay         }
371b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
37249b5e25fSSatish Balay       }
37349b5e25fSSatish Balay     }
374b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
37549b5e25fSSatish Balay   } else {
376b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
37749b5e25fSSatish Balay     for (i=0; i<a->mbs; i++) {
37849b5e25fSSatish Balay       for (j=0; j<bs; j++) {
379b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr);
38049b5e25fSSatish Balay         for (k=a->i[i]; k<a->i[i+1]; k++) {
38149b5e25fSSatish Balay           for (l=0; l<bs; l++) {
38249b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX)
38349b5e25fSSatish Balay             if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0) {
384b0a32e0cSBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," %d %g + %g i",bs*a->j[k]+l,
38549b5e25fSSatish Balay                 PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
38649b5e25fSSatish Balay             } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0) {
387b0a32e0cSBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," %d %g - %g i",bs*a->j[k]+l,
38849b5e25fSSatish Balay                 PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
38949b5e25fSSatish Balay             } else {
390b0a32e0cSBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
39149b5e25fSSatish Balay             }
39249b5e25fSSatish Balay #else
393b0a32e0cSBarry Smith             ierr = PetscViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr);
39449b5e25fSSatish Balay #endif
39549b5e25fSSatish Balay           }
39649b5e25fSSatish Balay         }
397b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
39849b5e25fSSatish Balay       }
39949b5e25fSSatish Balay     }
400b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
40149b5e25fSSatish Balay   }
402b0a32e0cSBarry Smith   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
40349b5e25fSSatish Balay   PetscFunctionReturn(0);
40449b5e25fSSatish Balay }
40549b5e25fSSatish Balay 
4064a2ae208SSatish Balay #undef __FUNCT__
4074a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Draw_Zoom"
408b0a32e0cSBarry Smith static int MatView_SeqSBAIJ_Draw_Zoom(PetscDraw draw,void *Aa)
40949b5e25fSSatish Balay {
41049b5e25fSSatish Balay   Mat          A = (Mat) Aa;
41149b5e25fSSatish Balay   Mat_SeqSBAIJ *a=(Mat_SeqSBAIJ*)A->data;
41249b5e25fSSatish Balay   int          row,ierr,i,j,k,l,mbs=a->mbs,color,bs=a->bs,bs2=a->bs2,rank;
41349b5e25fSSatish Balay   PetscReal    xl,yl,xr,yr,x_l,x_r,y_l,y_r;
41449b5e25fSSatish Balay   MatScalar    *aa;
41549b5e25fSSatish Balay   MPI_Comm     comm;
416b0a32e0cSBarry Smith   PetscViewer  viewer;
41749b5e25fSSatish Balay 
41849b5e25fSSatish Balay   PetscFunctionBegin;
41949b5e25fSSatish Balay   /*
42049b5e25fSSatish Balay       This is nasty. If this is called from an originally parallel matrix
42149b5e25fSSatish Balay    then all processes call this,but only the first has the matrix so the
42249b5e25fSSatish Balay    rest should return immediately.
42349b5e25fSSatish Balay   */
42449b5e25fSSatish Balay   ierr = PetscObjectGetComm((PetscObject)draw,&comm);CHKERRQ(ierr);
42549b5e25fSSatish Balay   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
42649b5e25fSSatish Balay   if (rank) PetscFunctionReturn(0);
42749b5e25fSSatish Balay 
42849b5e25fSSatish Balay   ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr);
42949b5e25fSSatish Balay 
430b0a32e0cSBarry Smith   ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr);
431b0a32e0cSBarry Smith   PetscDrawString(draw, .3*(xl+xr), .3*(yl+yr), PETSC_DRAW_BLACK, "symmetric");
43249b5e25fSSatish Balay 
43349b5e25fSSatish Balay   /* loop over matrix elements drawing boxes */
434b0a32e0cSBarry Smith   color = PETSC_DRAW_BLUE;
43549b5e25fSSatish Balay   for (i=0,row=0; i<mbs; i++,row+=bs) {
43649b5e25fSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
437c464158bSHong Zhang       y_l = A->m - row - 1.0; y_r = y_l + 1.0;
43849b5e25fSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
43949b5e25fSSatish Balay       aa = a->a + j*bs2;
44049b5e25fSSatish Balay       for (k=0; k<bs; k++) {
44149b5e25fSSatish Balay         for (l=0; l<bs; l++) {
44249b5e25fSSatish Balay           if (PetscRealPart(*aa++) >=  0.) continue;
443b0a32e0cSBarry Smith           ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
44449b5e25fSSatish Balay         }
44549b5e25fSSatish Balay       }
44649b5e25fSSatish Balay     }
44749b5e25fSSatish Balay   }
448b0a32e0cSBarry Smith   color = PETSC_DRAW_CYAN;
44949b5e25fSSatish Balay   for (i=0,row=0; i<mbs; i++,row+=bs) {
45049b5e25fSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
451c464158bSHong Zhang       y_l = A->m - row - 1.0; y_r = y_l + 1.0;
45249b5e25fSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
45349b5e25fSSatish Balay       aa = a->a + j*bs2;
45449b5e25fSSatish Balay       for (k=0; k<bs; k++) {
45549b5e25fSSatish Balay         for (l=0; l<bs; l++) {
45649b5e25fSSatish Balay           if (PetscRealPart(*aa++) != 0.) continue;
457b0a32e0cSBarry Smith           ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
45849b5e25fSSatish Balay         }
45949b5e25fSSatish Balay       }
46049b5e25fSSatish Balay     }
46149b5e25fSSatish Balay   }
46249b5e25fSSatish Balay 
463b0a32e0cSBarry Smith   color = PETSC_DRAW_RED;
46449b5e25fSSatish Balay   for (i=0,row=0; i<mbs; i++,row+=bs) {
46549b5e25fSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
466c464158bSHong Zhang       y_l = A->m - row - 1.0; y_r = y_l + 1.0;
46749b5e25fSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
46849b5e25fSSatish Balay       aa = a->a + j*bs2;
46949b5e25fSSatish Balay       for (k=0; k<bs; k++) {
47049b5e25fSSatish Balay         for (l=0; l<bs; l++) {
47149b5e25fSSatish Balay           if (PetscRealPart(*aa++) <= 0.) continue;
472b0a32e0cSBarry Smith           ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
47349b5e25fSSatish Balay         }
47449b5e25fSSatish Balay       }
47549b5e25fSSatish Balay     }
47649b5e25fSSatish Balay   }
47749b5e25fSSatish Balay   PetscFunctionReturn(0);
47849b5e25fSSatish Balay }
47949b5e25fSSatish Balay 
4804a2ae208SSatish Balay #undef __FUNCT__
4814a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Draw"
482b0a32e0cSBarry Smith static int MatView_SeqSBAIJ_Draw(Mat A,PetscViewer viewer)
48349b5e25fSSatish Balay {
48449b5e25fSSatish Balay   int          ierr;
48549b5e25fSSatish Balay   PetscReal    xl,yl,xr,yr,w,h;
486b0a32e0cSBarry Smith   PetscDraw    draw;
48749b5e25fSSatish Balay   PetscTruth   isnull;
48849b5e25fSSatish Balay 
48949b5e25fSSatish Balay   PetscFunctionBegin;
49049b5e25fSSatish Balay 
491b0a32e0cSBarry Smith   ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
492b0a32e0cSBarry Smith   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
49349b5e25fSSatish Balay 
49449b5e25fSSatish Balay   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
495c464158bSHong Zhang   xr  = A->m; yr = A->m; h = yr/10.0; w = xr/10.0;
49649b5e25fSSatish Balay   xr += w;    yr += h;  xl = -w;     yl = -h;
497b0a32e0cSBarry Smith   ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
498b0a32e0cSBarry Smith   ierr = PetscDrawZoom(draw,MatView_SeqSBAIJ_Draw_Zoom,A);CHKERRQ(ierr);
49949b5e25fSSatish Balay   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr);
50049b5e25fSSatish Balay   PetscFunctionReturn(0);
50149b5e25fSSatish Balay }
50249b5e25fSSatish Balay 
5034a2ae208SSatish Balay #undef __FUNCT__
5044a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ"
505b0a32e0cSBarry Smith int MatView_SeqSBAIJ(Mat A,PetscViewer viewer)
50649b5e25fSSatish Balay {
50749b5e25fSSatish Balay   int        ierr;
50849b5e25fSSatish Balay   PetscTruth issocket,isascii,isbinary,isdraw;
50949b5e25fSSatish Balay 
51049b5e25fSSatish Balay   PetscFunctionBegin;
511b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);CHKERRQ(ierr);
512b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
513fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr);
514fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
51549b5e25fSSatish Balay   if (issocket) {
51629bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Socket viewer not supported");
51749b5e25fSSatish Balay   } else if (isascii){
51849b5e25fSSatish Balay     ierr = MatView_SeqSBAIJ_ASCII(A,viewer);CHKERRQ(ierr);
51949b5e25fSSatish Balay   } else if (isbinary) {
52049b5e25fSSatish Balay     ierr = MatView_SeqSBAIJ_Binary(A,viewer);CHKERRQ(ierr);
52149b5e25fSSatish Balay   } else if (isdraw) {
52249b5e25fSSatish Balay     ierr = MatView_SeqSBAIJ_Draw(A,viewer);CHKERRQ(ierr);
52349b5e25fSSatish Balay   } else {
52429bbc08cSBarry Smith     SETERRQ1(1,"Viewer type %s not supported by SeqBAIJ matrices",((PetscObject)viewer)->type_name);
52549b5e25fSSatish Balay   }
52649b5e25fSSatish Balay   PetscFunctionReturn(0);
52749b5e25fSSatish Balay }
52849b5e25fSSatish Balay 
52949b5e25fSSatish Balay 
5304a2ae208SSatish Balay #undef __FUNCT__
5314a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqSBAIJ"
53287828ca2SBarry Smith int MatGetValues_SeqSBAIJ(Mat A,int m,int *im,int n,int *in,PetscScalar *v)
53349b5e25fSSatish Balay {
534045c9aa0SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
53549b5e25fSSatish Balay   int        *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j;
53649b5e25fSSatish Balay   int        *ai = a->i,*ailen = a->ilen;
53749b5e25fSSatish Balay   int        brow,bcol,ridx,cidx,bs=a->bs,bs2=a->bs2;
53849b5e25fSSatish Balay   MatScalar  *ap,*aa = a->a,zero = 0.0;
53949b5e25fSSatish Balay 
54049b5e25fSSatish Balay   PetscFunctionBegin;
54149b5e25fSSatish Balay   for (k=0; k<m; k++) { /* loop over rows */
54249b5e25fSSatish Balay     row  = im[k]; brow = row/bs;
54329bbc08cSBarry Smith     if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative row");
544c464158bSHong Zhang     if (row >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
54549b5e25fSSatish Balay     rp   = aj + ai[brow] ; ap = aa + bs2*ai[brow] ;
54649b5e25fSSatish Balay     nrow = ailen[brow];
54749b5e25fSSatish Balay     for (l=0; l<n; l++) { /* loop over columns */
54829bbc08cSBarry Smith       if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative column");
549c464158bSHong Zhang       if (in[l] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
55049b5e25fSSatish Balay       col  = in[l] ;
55149b5e25fSSatish Balay       bcol = col/bs;
55249b5e25fSSatish Balay       cidx = col%bs;
55349b5e25fSSatish Balay       ridx = row%bs;
55449b5e25fSSatish Balay       high = nrow;
55549b5e25fSSatish Balay       low  = 0; /* assume unsorted */
55649b5e25fSSatish Balay       while (high-low > 5) {
55749b5e25fSSatish Balay         t = (low+high)/2;
55849b5e25fSSatish Balay         if (rp[t] > bcol) high = t;
55949b5e25fSSatish Balay         else             low  = t;
56049b5e25fSSatish Balay       }
56149b5e25fSSatish Balay       for (i=low; i<high; i++) {
56249b5e25fSSatish Balay         if (rp[i] > bcol) break;
56349b5e25fSSatish Balay         if (rp[i] == bcol) {
56449b5e25fSSatish Balay           *v++ = ap[bs2*i+bs*cidx+ridx];
56549b5e25fSSatish Balay           goto finished;
56649b5e25fSSatish Balay         }
56749b5e25fSSatish Balay       }
56849b5e25fSSatish Balay       *v++ = zero;
56949b5e25fSSatish Balay       finished:;
57049b5e25fSSatish Balay     }
57149b5e25fSSatish Balay   }
57249b5e25fSSatish Balay   PetscFunctionReturn(0);
57349b5e25fSSatish Balay }
57449b5e25fSSatish Balay 
57549b5e25fSSatish Balay 
5764a2ae208SSatish Balay #undef __FUNCT__
5774a2ae208SSatish Balay #define __FUNCT__ "MatSetValuesBlocked_SeqSBAIJ"
57887828ca2SBarry Smith int MatSetValuesBlocked_SeqSBAIJ(Mat A,int m,int *im,int n,int *in,PetscScalar *v,InsertMode is)
57949b5e25fSSatish Balay {
58049b5e25fSSatish Balay   PetscFunctionBegin;
58129bbc08cSBarry Smith   SETERRQ(1,"Function not yet written for SBAIJ format");
58216cdd363SHong Zhang   /* PetscFunctionReturn(0); */
58349b5e25fSSatish Balay }
58449b5e25fSSatish Balay 
5854a2ae208SSatish Balay #undef __FUNCT__
5864a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_SeqSBAIJ"
58749b5e25fSSatish Balay int MatAssemblyEnd_SeqSBAIJ(Mat A,MatAssemblyType mode)
58849b5e25fSSatish Balay {
58949b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
59049b5e25fSSatish Balay   int        fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax;
591c464158bSHong Zhang   int        m = A->m,*ip,N,*ailen = a->ilen;
59249b5e25fSSatish Balay   int        mbs = a->mbs,bs2 = a->bs2,rmax = 0,ierr;
59349b5e25fSSatish Balay   MatScalar  *aa = a->a,*ap;
59449b5e25fSSatish Balay 
59549b5e25fSSatish Balay   PetscFunctionBegin;
59649b5e25fSSatish Balay   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0);
59749b5e25fSSatish Balay 
59849b5e25fSSatish Balay   if (m) rmax = ailen[0];
59949b5e25fSSatish Balay   for (i=1; i<mbs; i++) {
60049b5e25fSSatish Balay     /* move each row back by the amount of empty slots (fshift) before it*/
60149b5e25fSSatish Balay     fshift += imax[i-1] - ailen[i-1];
60249b5e25fSSatish Balay     rmax   = PetscMax(rmax,ailen[i]);
60349b5e25fSSatish Balay     if (fshift) {
60449b5e25fSSatish Balay       ip = aj + ai[i]; ap = aa + bs2*ai[i];
60549b5e25fSSatish Balay       N = ailen[i];
60649b5e25fSSatish Balay       for (j=0; j<N; j++) {
60749b5e25fSSatish Balay         ip[j-fshift] = ip[j];
60849b5e25fSSatish Balay         ierr = PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
60949b5e25fSSatish Balay       }
61049b5e25fSSatish Balay     }
61149b5e25fSSatish Balay     ai[i] = ai[i-1] + ailen[i-1];
61249b5e25fSSatish Balay   }
61349b5e25fSSatish Balay   if (mbs) {
61449b5e25fSSatish Balay     fshift += imax[mbs-1] - ailen[mbs-1];
61549b5e25fSSatish Balay     ai[mbs] = ai[mbs-1] + ailen[mbs-1];
61649b5e25fSSatish Balay   }
61749b5e25fSSatish Balay   /* reset ilen and imax for each row */
61849b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
61949b5e25fSSatish Balay     ailen[i] = imax[i] = ai[i+1] - ai[i];
62049b5e25fSSatish Balay   }
62149b5e25fSSatish Balay   a->s_nz = ai[mbs];
62249b5e25fSSatish Balay 
62349b5e25fSSatish Balay   /* diagonals may have moved, so kill the diagonal pointers */
62449b5e25fSSatish Balay   if (fshift && a->diag) {
62549b5e25fSSatish Balay     ierr = PetscFree(a->diag);CHKERRQ(ierr);
626b0a32e0cSBarry Smith     PetscLogObjectMemory(A,-(m+1)*sizeof(int));
62749b5e25fSSatish Balay     a->diag = 0;
62849b5e25fSSatish Balay   }
629b0a32e0cSBarry Smith   PetscLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Matrix size: %d X %d, block size %d; storage space: %d unneeded, %d used\n",
630c464158bSHong Zhang            m,A->m,a->bs,fshift*bs2,a->s_nz*bs2);
631b0a32e0cSBarry Smith   PetscLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Number of mallocs during MatSetValues is %d\n",
63249b5e25fSSatish Balay            a->reallocs);
633b0a32e0cSBarry Smith   PetscLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Most nonzeros blocks in any row is %d\n",rmax);
63449b5e25fSSatish Balay   a->reallocs          = 0;
63549b5e25fSSatish Balay   A->info.nz_unneeded  = (PetscReal)fshift*bs2;
63649b5e25fSSatish Balay 
63749b5e25fSSatish Balay   PetscFunctionReturn(0);
63849b5e25fSSatish Balay }
63949b5e25fSSatish Balay 
64049b5e25fSSatish Balay /*
64149b5e25fSSatish Balay    This function returns an array of flags which indicate the locations of contiguous
64249b5e25fSSatish Balay    blocks that should be zeroed. for eg: if bs = 3  and is = [0,1,2,3,5,6,7,8,9]
64349b5e25fSSatish Balay    then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)]
64449b5e25fSSatish Balay    Assume: sizes should be long enough to hold all the values.
64549b5e25fSSatish Balay */
6464a2ae208SSatish Balay #undef __FUNCT__
6474a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqSBAIJ_Check_Blocks"
64849b5e25fSSatish Balay static int MatZeroRows_SeqSBAIJ_Check_Blocks(int idx[],int n,int bs,int sizes[], int *bs_max)
64949b5e25fSSatish Balay {
65049b5e25fSSatish Balay   int        i,j,k,row;
65149b5e25fSSatish Balay   PetscTruth flg;
65249b5e25fSSatish Balay 
65349b5e25fSSatish Balay   PetscFunctionBegin;
65449b5e25fSSatish Balay   for (i=0,j=0; i<n; j++) {
65549b5e25fSSatish Balay     row = idx[i];
65649b5e25fSSatish Balay     if (row%bs!=0) { /* Not the begining of a block */
65749b5e25fSSatish Balay       sizes[j] = 1;
65849b5e25fSSatish Balay       i++;
65949b5e25fSSatish Balay     } else if (i+bs > n) { /* Beginning of a block, but complete block doesn't exist (at idx end) */
66049b5e25fSSatish Balay       sizes[j] = 1;         /* Also makes sure atleast 'bs' values exist for next else */
66149b5e25fSSatish Balay       i++;
66249b5e25fSSatish Balay     } else { /* Begining of the block, so check if the complete block exists */
66349b5e25fSSatish Balay       flg = PETSC_TRUE;
66449b5e25fSSatish Balay       for (k=1; k<bs; k++) {
66549b5e25fSSatish Balay         if (row+k != idx[i+k]) { /* break in the block */
66649b5e25fSSatish Balay           flg = PETSC_FALSE;
66749b5e25fSSatish Balay           break;
66849b5e25fSSatish Balay         }
66949b5e25fSSatish Balay       }
67049b5e25fSSatish Balay       if (flg == PETSC_TRUE) { /* No break in the bs */
67149b5e25fSSatish Balay         sizes[j] = bs;
67249b5e25fSSatish Balay         i+= bs;
67349b5e25fSSatish Balay       } else {
67449b5e25fSSatish Balay         sizes[j] = 1;
67549b5e25fSSatish Balay         i++;
67649b5e25fSSatish Balay       }
67749b5e25fSSatish Balay     }
67849b5e25fSSatish Balay   }
67949b5e25fSSatish Balay   *bs_max = j;
68049b5e25fSSatish Balay   PetscFunctionReturn(0);
68149b5e25fSSatish Balay }
68249b5e25fSSatish Balay 
6834a2ae208SSatish Balay #undef __FUNCT__
6844a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqSBAIJ"
68587828ca2SBarry Smith int MatZeroRows_SeqSBAIJ(Mat A,IS is,PetscScalar *diag)
68649b5e25fSSatish Balay {
68749b5e25fSSatish Balay   Mat_SeqSBAIJ  *sbaij=(Mat_SeqSBAIJ*)A->data;
68849b5e25fSSatish Balay   int           ierr,i,j,k,count,is_n,*is_idx,*rows;
68949b5e25fSSatish Balay   int           bs=sbaij->bs,bs2=sbaij->bs2,*sizes,row,bs_max;
69087828ca2SBarry Smith   PetscScalar   zero = 0.0;
69149b5e25fSSatish Balay   MatScalar     *aa;
69249b5e25fSSatish Balay 
69349b5e25fSSatish Balay   PetscFunctionBegin;
69449b5e25fSSatish Balay   /* Make a copy of the IS and  sort it */
69549b5e25fSSatish Balay   ierr = ISGetSize(is,&is_n);CHKERRQ(ierr);
69649b5e25fSSatish Balay   ierr = ISGetIndices(is,&is_idx);CHKERRQ(ierr);
69749b5e25fSSatish Balay 
69849b5e25fSSatish Balay   /* allocate memory for rows,sizes */
699b0a32e0cSBarry Smith   ierr = PetscMalloc((3*is_n+1)*sizeof(int),&rows);CHKERRQ(ierr);
70049b5e25fSSatish Balay   sizes = rows + is_n;
70149b5e25fSSatish Balay 
70249b5e25fSSatish Balay   /* initialize copy IS values to rows, and sort them */
70349b5e25fSSatish Balay   for (i=0; i<is_n; i++) { rows[i] = is_idx[i]; }
70449b5e25fSSatish Balay   ierr = PetscSortInt(is_n,rows);CHKERRQ(ierr);
70549b5e25fSSatish Balay   if (sbaij->keepzeroedrows) { /* do not change nonzero structure */
70649b5e25fSSatish Balay     for (i=0; i<is_n; i++) { sizes[i] = 1; } /* sizes: size of blocks, = 1 or bs */
70749b5e25fSSatish Balay     bs_max = is_n;            /* bs_max: num. of contiguous block row in the row */
70849b5e25fSSatish Balay   } else {
70949b5e25fSSatish Balay     ierr = MatZeroRows_SeqSBAIJ_Check_Blocks(rows,is_n,bs,sizes,&bs_max);CHKERRQ(ierr);
71049b5e25fSSatish Balay   }
71149b5e25fSSatish Balay   ierr = ISRestoreIndices(is,&is_idx);CHKERRQ(ierr);
71249b5e25fSSatish Balay 
71349b5e25fSSatish Balay   for (i=0,j=0; i<bs_max; j+=sizes[i],i++) {
71449b5e25fSSatish Balay     row   = rows[j];                  /* row to be zeroed */
715c464158bSHong Zhang     if (row < 0 || row > A->m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"row %d out of range",row);
71649b5e25fSSatish Balay     count = (sbaij->i[row/bs +1] - sbaij->i[row/bs])*bs; /* num. of elements in the row */
71749b5e25fSSatish Balay     aa    = sbaij->a + sbaij->i[row/bs]*bs2 + (row%bs);
71849b5e25fSSatish Balay     if (sizes[i] == bs && !sbaij->keepzeroedrows) {
71949b5e25fSSatish Balay       if (diag) {
72049b5e25fSSatish Balay         if (sbaij->ilen[row/bs] > 0) {
72149b5e25fSSatish Balay           sbaij->ilen[row/bs] = 1;
72249b5e25fSSatish Balay           sbaij->j[sbaij->i[row/bs]] = row/bs;
72349b5e25fSSatish Balay           ierr = PetscMemzero(aa,count*bs*sizeof(MatScalar));CHKERRQ(ierr);
72449b5e25fSSatish Balay         }
72549b5e25fSSatish Balay         /* Now insert all the diagoanl values for this bs */
72649b5e25fSSatish Balay         for (k=0; k<bs; k++) {
72749b5e25fSSatish Balay           ierr = (*A->ops->setvalues)(A,1,rows+j+k,1,rows+j+k,diag,INSERT_VALUES);CHKERRQ(ierr);
72849b5e25fSSatish Balay         }
72949b5e25fSSatish Balay       } else { /* (!diag) */
73049b5e25fSSatish Balay         sbaij->ilen[row/bs] = 0;
73149b5e25fSSatish Balay       } /* end (!diag) */
73249b5e25fSSatish Balay     } else { /* (sizes[i] != bs), broken block */
73349b5e25fSSatish Balay #if defined (PETSC_USE_DEBUG)
73429bbc08cSBarry Smith       if (sizes[i] != 1) SETERRQ(1,"Internal Error. Value should be 1");
73549b5e25fSSatish Balay #endif
73649b5e25fSSatish Balay       for (k=0; k<count; k++) {
73749b5e25fSSatish Balay         aa[0] = zero;
73849b5e25fSSatish Balay         aa+=bs;
73949b5e25fSSatish Balay       }
74049b5e25fSSatish Balay       if (diag) {
74149b5e25fSSatish Balay         ierr = (*A->ops->setvalues)(A,1,rows+j,1,rows+j,diag,INSERT_VALUES);CHKERRQ(ierr);
74249b5e25fSSatish Balay       }
74349b5e25fSSatish Balay     }
74449b5e25fSSatish Balay   }
74549b5e25fSSatish Balay 
74649b5e25fSSatish Balay   ierr = PetscFree(rows);CHKERRQ(ierr);
74749b5e25fSSatish Balay   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
74849b5e25fSSatish Balay   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
74949b5e25fSSatish Balay   PetscFunctionReturn(0);
75049b5e25fSSatish Balay }
75149b5e25fSSatish Balay 
75249b5e25fSSatish Balay /* Only add/insert a(i,j) with i<=j (blocks).
75349b5e25fSSatish Balay    Any a(i,j) with i>j input by user is ingored.
75449b5e25fSSatish Balay */
75549b5e25fSSatish Balay 
7564a2ae208SSatish Balay #undef __FUNCT__
7574a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqSBAIJ"
75887828ca2SBarry Smith int MatSetValues_SeqSBAIJ(Mat A,int m,int *im,int n,int *in,PetscScalar *v,InsertMode is)
75949b5e25fSSatish Balay {
76049b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
76149b5e25fSSatish Balay   int         *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted=a->sorted;
76249b5e25fSSatish Balay   int         *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented;
76349b5e25fSSatish Balay   int         *aj=a->j,nonew=a->nonew,bs=a->bs,brow,bcol;
76449b5e25fSSatish Balay   int         ridx,cidx,bs2=a->bs2,ierr;
76549b5e25fSSatish Balay   MatScalar   *ap,value,*aa=a->a,*bap;
76649b5e25fSSatish Balay 
76749b5e25fSSatish Balay   PetscFunctionBegin;
76849b5e25fSSatish Balay 
76949b5e25fSSatish Balay   for (k=0; k<m; k++) { /* loop over added rows */
77049b5e25fSSatish Balay     row  = im[k];       /* row number */
77149b5e25fSSatish Balay     brow = row/bs;      /* block row number */
77249b5e25fSSatish Balay     if (row < 0) continue;
77349b5e25fSSatish Balay #if defined(PETSC_USE_BOPT_g)
774c464158bSHong Zhang     if (row >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",row,A->m);
77549b5e25fSSatish Balay #endif
77649b5e25fSSatish Balay     rp   = aj + ai[brow]; /*ptr to beginning of column value of the row block*/
77749b5e25fSSatish Balay     ap   = aa + bs2*ai[brow]; /*ptr to beginning of element value of the row block*/
77849b5e25fSSatish Balay     rmax = imax[brow];  /* maximum space allocated for this row */
77949b5e25fSSatish Balay     nrow = ailen[brow]; /* actual length of this row */
78049b5e25fSSatish Balay     low  = 0;
78149b5e25fSSatish Balay 
78249b5e25fSSatish Balay     for (l=0; l<n; l++) { /* loop over added columns */
78349b5e25fSSatish Balay       if (in[l] < 0) continue;
78449b5e25fSSatish Balay #if defined(PETSC_USE_BOPT_g)
785c464158bSHong Zhang       if (in[l] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",in[l],A->m);
78649b5e25fSSatish Balay #endif
78749b5e25fSSatish Balay       col = in[l];
78849b5e25fSSatish Balay       bcol = col/bs;              /* block col number */
78949b5e25fSSatish Balay 
79049b5e25fSSatish Balay       if (brow <= bcol){
79149b5e25fSSatish Balay         ridx = row % bs; cidx = col % bs; /*row and col index inside the block */
7928549e402SHong Zhang         if ((brow==bcol && ridx<=cidx) || (brow<bcol)){
79349b5e25fSSatish Balay           /* element value a(k,l) */
79449b5e25fSSatish Balay           if (roworiented) {
79549b5e25fSSatish Balay             value = v[l + k*n];
79649b5e25fSSatish Balay           } else {
79749b5e25fSSatish Balay             value = v[k + l*m];
79849b5e25fSSatish Balay           }
79949b5e25fSSatish Balay 
80049b5e25fSSatish Balay           /* move pointer bap to a(k,l) quickly and add/insert value */
80149b5e25fSSatish Balay           if (!sorted) low = 0; high = nrow;
80249b5e25fSSatish Balay           while (high-low > 7) {
80349b5e25fSSatish Balay             t = (low+high)/2;
80449b5e25fSSatish Balay             if (rp[t] > bcol) high = t;
80549b5e25fSSatish Balay             else              low  = t;
80649b5e25fSSatish Balay           }
80749b5e25fSSatish Balay           for (i=low; i<high; i++) {
80849b5e25fSSatish Balay             /* printf("The loop of i=low.., rp[%d]=%d\n",i,rp[i]); */
80949b5e25fSSatish Balay             if (rp[i] > bcol) break;
81049b5e25fSSatish Balay             if (rp[i] == bcol) {
81149b5e25fSSatish Balay               bap  = ap +  bs2*i + bs*cidx + ridx;
81249b5e25fSSatish Balay               if (is == ADD_VALUES) *bap += value;
81349b5e25fSSatish Balay               else                  *bap  = value;
8148549e402SHong Zhang               /* for diag block, add/insert its symmetric element a(cidx,ridx) */
8158549e402SHong Zhang               if (brow == bcol && ridx < cidx){
8168549e402SHong Zhang                 bap  = ap +  bs2*i + bs*ridx + cidx;
8178549e402SHong Zhang                 if (is == ADD_VALUES) *bap += value;
8188549e402SHong Zhang                 else                  *bap  = value;
8198549e402SHong Zhang               }
82049b5e25fSSatish Balay               goto noinsert1;
82149b5e25fSSatish Balay             }
82249b5e25fSSatish Balay           }
82349b5e25fSSatish Balay 
82449b5e25fSSatish Balay       if (nonew == 1) goto noinsert1;
82529bbc08cSBarry Smith       else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix");
82649b5e25fSSatish Balay       if (nrow >= rmax) {
82749b5e25fSSatish Balay         /* there is no extra room in row, therefore enlarge */
82849b5e25fSSatish Balay         int       new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j;
82949b5e25fSSatish Balay         MatScalar *new_a;
83049b5e25fSSatish Balay 
83129bbc08cSBarry Smith         if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix");
83249b5e25fSSatish Balay 
83349b5e25fSSatish Balay         /* Malloc new storage space */
83449b5e25fSSatish Balay         len   = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int);
83582502324SSatish Balay         ierr  = PetscMalloc(len,&new_a);CHKERRQ(ierr);
83649b5e25fSSatish Balay         new_j = (int*)(new_a + bs2*new_nz);
83749b5e25fSSatish Balay         new_i = new_j + new_nz;
83849b5e25fSSatish Balay 
83949b5e25fSSatish Balay         /* copy over old data into new slots */
84049b5e25fSSatish Balay         for (ii=0; ii<brow+1; ii++) {new_i[ii] = ai[ii];}
84149b5e25fSSatish Balay         for (ii=brow+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;}
84249b5e25fSSatish Balay         ierr = PetscMemcpy(new_j,aj,(ai[brow]+nrow)*sizeof(int));CHKERRQ(ierr);
84349b5e25fSSatish Balay         len  = (new_nz - CHUNKSIZE - ai[brow] - nrow);
84449b5e25fSSatish Balay         ierr = PetscMemcpy(new_j+ai[brow]+nrow+CHUNKSIZE,aj+ai[brow]+nrow,len*sizeof(int));CHKERRQ(ierr);
84549b5e25fSSatish Balay         ierr = PetscMemcpy(new_a,aa,(ai[brow]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr);
84649b5e25fSSatish Balay         ierr = PetscMemzero(new_a+bs2*(ai[brow]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr);
84749b5e25fSSatish Balay         ierr = PetscMemcpy(new_a+bs2*(ai[brow]+nrow+CHUNKSIZE),aa+bs2*(ai[brow]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr);
84849b5e25fSSatish Balay         /* free up old matrix storage */
84949b5e25fSSatish Balay         ierr = PetscFree(a->a);CHKERRQ(ierr);
85049b5e25fSSatish Balay         if (!a->singlemalloc) {
85149b5e25fSSatish Balay           ierr = PetscFree(a->i);CHKERRQ(ierr);
85249b5e25fSSatish Balay           ierr = PetscFree(a->j);CHKERRQ(ierr);
85349b5e25fSSatish Balay         }
85449b5e25fSSatish Balay         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
85549b5e25fSSatish Balay         a->singlemalloc = PETSC_TRUE;
85649b5e25fSSatish Balay 
85749b5e25fSSatish Balay         rp   = aj + ai[brow]; ap = aa + bs2*ai[brow];
85849b5e25fSSatish Balay         rmax = imax[brow] = imax[brow] + CHUNKSIZE;
859b0a32e0cSBarry Smith         PetscLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar)));
86049b5e25fSSatish Balay         a->s_maxnz += bs2*CHUNKSIZE;
86149b5e25fSSatish Balay         a->reallocs++;
86249b5e25fSSatish Balay         a->s_nz++;
86349b5e25fSSatish Balay       }
86449b5e25fSSatish Balay 
86549b5e25fSSatish Balay       N = nrow++ - 1;
86649b5e25fSSatish Balay       /* shift up all the later entries in this row */
86749b5e25fSSatish Balay       for (ii=N; ii>=i; ii--) {
86849b5e25fSSatish Balay         rp[ii+1] = rp[ii];
86949b5e25fSSatish Balay         ierr     = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr);
87049b5e25fSSatish Balay       }
87149b5e25fSSatish Balay       if (N>=i) {
87249b5e25fSSatish Balay         ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr);
87349b5e25fSSatish Balay       }
87449b5e25fSSatish Balay       rp[i]                      = bcol;
87549b5e25fSSatish Balay       ap[bs2*i + bs*cidx + ridx] = value;
87649b5e25fSSatish Balay       noinsert1:;
87749b5e25fSSatish Balay       low = i;
87849b5e25fSSatish Balay       /* } */
8798549e402SHong Zhang         }
88049b5e25fSSatish Balay       } /* end of if .. if..  */
88149b5e25fSSatish Balay     }                     /* end of loop over added columns */
88249b5e25fSSatish Balay     ailen[brow] = nrow;
88349b5e25fSSatish Balay   }                       /* end of loop over added rows */
88449b5e25fSSatish Balay 
88549b5e25fSSatish Balay   PetscFunctionReturn(0);
88649b5e25fSSatish Balay }
88749b5e25fSSatish Balay 
888915354c9SHong Zhang extern int MatCholeskyFactorSymbolic_SeqSBAIJ(Mat,IS,PetscReal,Mat*);
889915354c9SHong Zhang extern int MatCholeskyFactor_SeqSBAIJ(Mat,IS,PetscReal);
89049b5e25fSSatish Balay extern int MatIncreaseOverlap_SeqSBAIJ(Mat,int,IS*,int);
89149b5e25fSSatish Balay extern int MatGetSubMatrix_SeqSBAIJ(Mat,IS,IS,int,MatReuse,Mat*);
89249b5e25fSSatish Balay extern int MatGetSubMatrices_SeqSBAIJ(Mat,int,IS*,IS*,MatReuse,Mat**);
89349b5e25fSSatish Balay extern int MatMultTranspose_SeqSBAIJ(Mat,Vec,Vec);
89449b5e25fSSatish Balay extern int MatMultTransposeAdd_SeqSBAIJ(Mat,Vec,Vec,Vec);
89587828ca2SBarry Smith extern int MatScale_SeqSBAIJ(PetscScalar*,Mat);
89649b5e25fSSatish Balay extern int MatNorm_SeqSBAIJ(Mat,NormType,PetscReal *);
89749b5e25fSSatish Balay extern int MatEqual_SeqSBAIJ(Mat,Mat,PetscTruth*);
89849b5e25fSSatish Balay extern int MatGetDiagonal_SeqSBAIJ(Mat,Vec);
89949b5e25fSSatish Balay extern int MatDiagonalScale_SeqSBAIJ(Mat,Vec,Vec);
90049b5e25fSSatish Balay extern int MatGetInfo_SeqSBAIJ(Mat,MatInfoType,MatInfo *);
90149b5e25fSSatish Balay extern int MatZeroEntries_SeqSBAIJ(Mat);
90213a02c98SHong Zhang extern int MatGetRowMax_SeqSBAIJ(Mat,Vec);
90349b5e25fSSatish Balay 
90449b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_N(Mat,Vec,Vec);
90549b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_1(Mat,Vec,Vec);
90649b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_2(Mat,Vec,Vec);
90749b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_3(Mat,Vec,Vec);
90849b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_4(Mat,Vec,Vec);
90949b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_5(Mat,Vec,Vec);
91049b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_6(Mat,Vec,Vec);
91149b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_7(Mat,Vec,Vec);
91249b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_7(Mat,Vec,Vec);
91349b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_6(Mat,Vec,Vec);
91449b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_5(Mat,Vec,Vec);
91549b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_4(Mat,Vec,Vec);
91649b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_3(Mat,Vec,Vec);
91749b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_2(Mat,Vec,Vec);
91849b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_1(Mat,Vec,Vec);
91949b5e25fSSatish Balay 
92007f98182SSatish Balay extern int MatSolve_SeqSBAIJ_1_NaturalOrdering(Mat,Vec,Vec);
92107f98182SSatish Balay extern int MatSolve_SeqSBAIJ_2_NaturalOrdering(Mat,Vec,Vec);
92207f98182SSatish Balay extern int MatSolve_SeqSBAIJ_3_NaturalOrdering(Mat,Vec,Vec);
92307f98182SSatish Balay extern int MatSolve_SeqSBAIJ_4_NaturalOrdering(Mat,Vec,Vec);
92407f98182SSatish Balay extern int MatSolve_SeqSBAIJ_5_NaturalOrdering(Mat,Vec,Vec);
92507f98182SSatish Balay extern int MatSolve_SeqSBAIJ_6_NaturalOrdering(Mat,Vec,Vec);
92607f98182SSatish Balay extern int MatSolve_SeqSBAIJ_7_NaturalOrdering(Mat,Vec,Vec);
92707f98182SSatish Balay extern int MatSolve_SeqSBAIJ_N_NaturalOrdering(Mat,Vec,Vec);
92807f98182SSatish Balay 
9295f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat,Mat*);
9305f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_1(Mat,Mat*);
9315f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_2(Mat,Mat*);
9325f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_3(Mat,Mat*);
9335f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_4(Mat,Mat*);
9345f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_5(Mat,Mat*);
9355f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_6(Mat,Mat*);
93657d5e339SHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_7(Mat,Mat*);
93749b5e25fSSatish Balay 
93849b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_1(Mat,Vec,Vec);
93949b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_2(Mat,Vec,Vec);
94049b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_3(Mat,Vec,Vec);
94149b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_4(Mat,Vec,Vec);
94249b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_5(Mat,Vec,Vec);
94349b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_6(Mat,Vec,Vec);
94449b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_7(Mat,Vec,Vec);
94549b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_N(Mat,Vec,Vec);
94649b5e25fSSatish Balay 
94749b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_1(Mat,Vec,Vec,Vec);
94849b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_2(Mat,Vec,Vec,Vec);
94949b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_3(Mat,Vec,Vec,Vec);
95049b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_4(Mat,Vec,Vec,Vec);
95149b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_5(Mat,Vec,Vec,Vec);
95249b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_6(Mat,Vec,Vec,Vec);
95349b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_7(Mat,Vec,Vec,Vec);
95449b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_N(Mat,Vec,Vec,Vec);
95549b5e25fSSatish Balay 
9564d101231SSatish Balay #ifdef HAVE_MatICCFactor
9571a3463dfSHong Zhang /* modefied from MatILUFactor_SeqSBAIJ, needs further work!  */
9584a2ae208SSatish Balay #undef __FUNCT__
9594d101231SSatish Balay #define __FUNCT__ "MatICCFactor_SeqSBAIJ"
9604d101231SSatish Balay int MatICCFactor_SeqSBAIJ(Mat inA,IS row,PetscReal fill,int level)
96149b5e25fSSatish Balay {
9624ccecd49SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)inA->data;
96349b5e25fSSatish Balay   Mat         outA;
96449b5e25fSSatish Balay   int         ierr;
96549b5e25fSSatish Balay   PetscTruth  row_identity,col_identity;
96649b5e25fSSatish Balay 
96749b5e25fSSatish Balay   PetscFunctionBegin;
9681a3463dfSHong Zhang   /*
96929bbc08cSBarry Smith   if (level != 0) SETERRQ(PETSC_ERR_SUP,"Only levels = 0 supported for in-place ILU");
97049b5e25fSSatish Balay   ierr = ISIdentity(row,&row_identity);CHKERRQ(ierr);
97149b5e25fSSatish Balay   ierr = ISIdentity(col,&col_identity);CHKERRQ(ierr);
97249b5e25fSSatish Balay   if (!row_identity || !col_identity) {
97329bbc08cSBarry Smith     SETERRQ(1,"Row and column permutations must be identity for in-place ICC");
97449b5e25fSSatish Balay   }
9751a3463dfSHong Zhang   */
97649b5e25fSSatish Balay 
97749b5e25fSSatish Balay   outA          = inA;
9781260e1f8SHong Zhang   inA->factor   = FACTOR_CHOLESKY;
97949b5e25fSSatish Balay 
98049b5e25fSSatish Balay   if (!a->diag) {
9811a3463dfSHong Zhang     ierr = MatMarkDiagonal_SeqSBAIJ(inA);CHKERRQ(ierr);
98249b5e25fSSatish Balay   }
98349b5e25fSSatish Balay   /*
98449b5e25fSSatish Balay       Blocksize 2, 3, 4, 5, 6 and 7 have a special faster factorization/solver
98549b5e25fSSatish Balay       for ILU(0) factorization with natural ordering
98649b5e25fSSatish Balay   */
98749b5e25fSSatish Balay   switch (a->bs) {
98849b5e25fSSatish Balay   case 1:
9890fe9e5beSBarry Smith     inA->ops->solve            = MatSolve_SeqSBAIJ_1_NaturalOrdering;
99007f98182SSatish Balay     inA->ops->solvetranspose   = MatSolve_SeqSBAIJ_1_NaturalOrdering;
9914d101231SSatish Balay     PetscLoginfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering solvetrans BS=1\n");
99249b5e25fSSatish Balay   case 2:
9931a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering;
9941a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_2_NaturalOrdering;
9951a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_2_NaturalOrdering;
9964d101231SSatish Balay     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=2\n");
99749b5e25fSSatish Balay     break;
99849b5e25fSSatish Balay   case 3:
9991a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering;
10001a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_3_NaturalOrdering;
10011a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_3_NaturalOrdering;
10024d101231SSatish Balay     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=3\n");
100349b5e25fSSatish Balay     break;
100449b5e25fSSatish Balay   case 4:
10051a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering;
10061a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_4_NaturalOrdering;
10071a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_4_NaturalOrdering;
10084d101231SSatish Balay     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=4\n");
100949b5e25fSSatish Balay     break;
101049b5e25fSSatish Balay   case 5:
10111a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering;
10121a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_5_NaturalOrdering;
10131a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_5_NaturalOrdering;
10144d101231SSatish Balay     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=5\n");
101549b5e25fSSatish Balay     break;
101649b5e25fSSatish Balay   case 6:
10171a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering;
10181a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_6_NaturalOrdering;
10191a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_6_NaturalOrdering;
10204d101231SSatish Balay     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=6\n");
102149b5e25fSSatish Balay     break;
102249b5e25fSSatish Balay   case 7:
10231a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering;
10241a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_7_NaturalOrdering;
10251a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_7_NaturalOrdering;
10264d101231SSatish Balay     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=7\n");
102749b5e25fSSatish Balay     break;
102849b5e25fSSatish Balay   default:
102949b5e25fSSatish Balay     a->row        = row;
10301a3463dfSHong Zhang     a->icol       = col;
103149b5e25fSSatish Balay     ierr          = PetscObjectReference((PetscObject)row);CHKERRQ(ierr);
103249b5e25fSSatish Balay     ierr          = PetscObjectReference((PetscObject)col);CHKERRQ(ierr);
103349b5e25fSSatish Balay 
103449b5e25fSSatish Balay     /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */
103549b5e25fSSatish Balay     ierr = ISInvertPermutation(col,PETSC_DECIDE, &(a->icol));CHKERRQ(ierr);
1036b0a32e0cSBarry Smith     PetscLogObjectParent(inA,a->icol);
103749b5e25fSSatish Balay 
103849b5e25fSSatish Balay     if (!a->solve_work) {
103987828ca2SBarry Smith       ierr = PetscMalloc((A->m+a->bs)*sizeof(PetscScalar),&a->solve_work);CHKERRQ(ierr);
104087828ca2SBarry Smith       PetscLogObjectMemory(inA,(A->m+a->bs)*sizeof(PetscScalar));
104149b5e25fSSatish Balay     }
104249b5e25fSSatish Balay   }
104349b5e25fSSatish Balay 
10441a3463dfSHong Zhang   ierr = MatCholeskyFactorNumeric(inA,&outA);CHKERRQ(ierr);
104549b5e25fSSatish Balay 
104649b5e25fSSatish Balay   PetscFunctionReturn(0);
104749b5e25fSSatish Balay }
1048950f1e5bSHong Zhang #endif
1049950f1e5bSHong Zhang 
10504a2ae208SSatish Balay #undef __FUNCT__
10514a2ae208SSatish Balay #define __FUNCT__ "MatPrintHelp_SeqSBAIJ"
105249b5e25fSSatish Balay int MatPrintHelp_SeqSBAIJ(Mat A)
105349b5e25fSSatish Balay {
105449b5e25fSSatish Balay   static PetscTruth called = PETSC_FALSE;
105549b5e25fSSatish Balay   MPI_Comm          comm = A->comm;
105649b5e25fSSatish Balay   int               ierr;
105749b5e25fSSatish Balay 
105849b5e25fSSatish Balay   PetscFunctionBegin;
105949b5e25fSSatish Balay   if (called) {PetscFunctionReturn(0);} else called = PETSC_TRUE;
106049b5e25fSSatish Balay   ierr = (*PetscHelpPrintf)(comm," Options for MATSEQSBAIJ and MATMPISBAIJ matrix formats (the defaults):\n");CHKERRQ(ierr);
106149b5e25fSSatish Balay   ierr = (*PetscHelpPrintf)(comm,"  -mat_block_size <block_size>\n");CHKERRQ(ierr);
106249b5e25fSSatish Balay   PetscFunctionReturn(0);
106349b5e25fSSatish Balay }
106449b5e25fSSatish Balay 
106549b5e25fSSatish Balay EXTERN_C_BEGIN
10664a2ae208SSatish Balay #undef __FUNCT__
10674a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices_SeqSBAIJ"
106849b5e25fSSatish Balay int MatSeqSBAIJSetColumnIndices_SeqSBAIJ(Mat mat,int *indices)
106949b5e25fSSatish Balay {
1070045c9aa0SHong Zhang   Mat_SeqSBAIJ *baij = (Mat_SeqSBAIJ *)mat->data;
107149b5e25fSSatish Balay   int         i,nz,n;
107249b5e25fSSatish Balay 
107349b5e25fSSatish Balay   PetscFunctionBegin;
1074045c9aa0SHong Zhang   nz = baij->s_maxnz;
1075c464158bSHong Zhang   n  = mat->n;
107649b5e25fSSatish Balay   for (i=0; i<nz; i++) {
107749b5e25fSSatish Balay     baij->j[i] = indices[i];
107849b5e25fSSatish Balay   }
1079045c9aa0SHong Zhang   baij->s_nz = nz;
108049b5e25fSSatish Balay   for (i=0; i<n; i++) {
108149b5e25fSSatish Balay     baij->ilen[i] = baij->imax[i];
108249b5e25fSSatish Balay   }
108349b5e25fSSatish Balay 
108449b5e25fSSatish Balay   PetscFunctionReturn(0);
108549b5e25fSSatish Balay }
108649b5e25fSSatish Balay EXTERN_C_END
108749b5e25fSSatish Balay 
10884a2ae208SSatish Balay #undef __FUNCT__
10894a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices"
109049b5e25fSSatish Balay /*@
109119585528SSatish Balay     MatSeqSBAIJSetColumnIndices - Set the column indices for all the rows
109249b5e25fSSatish Balay        in the matrix.
109349b5e25fSSatish Balay 
109449b5e25fSSatish Balay   Input Parameters:
109519585528SSatish Balay +  mat     - the SeqSBAIJ matrix
109649b5e25fSSatish Balay -  indices - the column indices
109749b5e25fSSatish Balay 
109849b5e25fSSatish Balay   Level: advanced
109949b5e25fSSatish Balay 
110049b5e25fSSatish Balay   Notes:
110149b5e25fSSatish Balay     This can be called if you have precomputed the nonzero structure of the
110249b5e25fSSatish Balay   matrix and want to provide it to the matrix object to improve the performance
110349b5e25fSSatish Balay   of the MatSetValues() operation.
110449b5e25fSSatish Balay 
110549b5e25fSSatish Balay     You MUST have set the correct numbers of nonzeros per row in the call to
110619585528SSatish Balay   MatCreateSeqSBAIJ().
110749b5e25fSSatish Balay 
1108ab9f2c04SSatish Balay     MUST be called before any calls to MatSetValues()
110949b5e25fSSatish Balay 
1110ab9f2c04SSatish Balay .seealso: MatCreateSeqSBAIJ
111149b5e25fSSatish Balay @*/
111249b5e25fSSatish Balay int MatSeqSBAIJSetColumnIndices(Mat mat,int *indices)
111349b5e25fSSatish Balay {
111449b5e25fSSatish Balay   int ierr,(*f)(Mat,int *);
111549b5e25fSSatish Balay 
111649b5e25fSSatish Balay   PetscFunctionBegin;
111749b5e25fSSatish Balay   PetscValidHeaderSpecific(mat,MAT_COOKIE);
1118c134de8dSSatish Balay   ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqSBAIJSetColumnIndices_C",(void (**)(void))&f);CHKERRQ(ierr);
111949b5e25fSSatish Balay   if (f) {
112049b5e25fSSatish Balay     ierr = (*f)(mat,indices);CHKERRQ(ierr);
112149b5e25fSSatish Balay   } else {
112229bbc08cSBarry Smith     SETERRQ(1,"Wrong type of matrix to set column indices");
112349b5e25fSSatish Balay   }
112449b5e25fSSatish Balay   PetscFunctionReturn(0);
112549b5e25fSSatish Balay }
112649b5e25fSSatish Balay 
11274a2ae208SSatish Balay #undef __FUNCT__
11284a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqSBAIJ"
1129273d9f13SBarry Smith int MatSetUpPreallocation_SeqSBAIJ(Mat A)
1130273d9f13SBarry Smith {
1131273d9f13SBarry Smith   int        ierr;
1132273d9f13SBarry Smith 
1133273d9f13SBarry Smith   PetscFunctionBegin;
1134273d9f13SBarry Smith   ierr =  MatSeqSBAIJSetPreallocation(A,1,PETSC_DEFAULT,0);CHKERRQ(ierr);
1135273d9f13SBarry Smith   PetscFunctionReturn(0);
1136273d9f13SBarry Smith }
1137273d9f13SBarry Smith 
113849b5e25fSSatish Balay /* -------------------------------------------------------------------*/
113949b5e25fSSatish Balay static struct _MatOps MatOps_Values = {MatSetValues_SeqSBAIJ,
114049b5e25fSSatish Balay        MatGetRow_SeqSBAIJ,
114149b5e25fSSatish Balay        MatRestoreRow_SeqSBAIJ,
114249b5e25fSSatish Balay        MatMult_SeqSBAIJ_N,
114349b5e25fSSatish Balay        MatMultAdd_SeqSBAIJ_N,
114449b5e25fSSatish Balay        MatMultTranspose_SeqSBAIJ,
114549b5e25fSSatish Balay        MatMultTransposeAdd_SeqSBAIJ,
114649b5e25fSSatish Balay        MatSolve_SeqSBAIJ_N,
114749b5e25fSSatish Balay        0,
114849b5e25fSSatish Balay        0,
114949b5e25fSSatish Balay        0,
115049b5e25fSSatish Balay        0,
11515f4ef2deSHong Zhang        MatCholeskyFactor_SeqSBAIJ,
1152d06b337dSHong Zhang        MatRelax_SeqSBAIJ,
115349b5e25fSSatish Balay        MatTranspose_SeqSBAIJ,
115449b5e25fSSatish Balay        MatGetInfo_SeqSBAIJ,
115549b5e25fSSatish Balay        MatEqual_SeqSBAIJ,
115649b5e25fSSatish Balay        MatGetDiagonal_SeqSBAIJ,
115749b5e25fSSatish Balay        MatDiagonalScale_SeqSBAIJ,
115849b5e25fSSatish Balay        MatNorm_SeqSBAIJ,
115949b5e25fSSatish Balay        0,
116049b5e25fSSatish Balay        MatAssemblyEnd_SeqSBAIJ,
116149b5e25fSSatish Balay        0,
116249b5e25fSSatish Balay        MatSetOption_SeqSBAIJ,
116349b5e25fSSatish Balay        MatZeroEntries_SeqSBAIJ,
116449b5e25fSSatish Balay        MatZeroRows_SeqSBAIJ,
116549b5e25fSSatish Balay        0,
116649b5e25fSSatish Balay        0,
11675f4ef2deSHong Zhang        MatCholeskyFactorSymbolic_SeqSBAIJ,
11685f4ef2deSHong Zhang        MatCholeskyFactorNumeric_SeqSBAIJ_N,
1169273d9f13SBarry Smith        MatSetUpPreallocation_SeqSBAIJ,
1170c464158bSHong Zhang        0,
11714d101231SSatish Balay        MatICCFactorSymbolic_SeqSBAIJ,
117249b5e25fSSatish Balay        0,
117349b5e25fSSatish Balay        0,
117449b5e25fSSatish Balay        MatDuplicate_SeqSBAIJ,
117549b5e25fSSatish Balay        0,
117649b5e25fSSatish Balay        0,
117749b5e25fSSatish Balay        0,
1178950f1e5bSHong Zhang        0,
117949b5e25fSSatish Balay        0,
118049b5e25fSSatish Balay        MatGetSubMatrices_SeqSBAIJ,
118149b5e25fSSatish Balay        MatIncreaseOverlap_SeqSBAIJ,
118249b5e25fSSatish Balay        MatGetValues_SeqSBAIJ,
118349b5e25fSSatish Balay        0,
118449b5e25fSSatish Balay        MatPrintHelp_SeqSBAIJ,
118549b5e25fSSatish Balay        MatScale_SeqSBAIJ,
118649b5e25fSSatish Balay        0,
118749b5e25fSSatish Balay        0,
118849b5e25fSSatish Balay        0,
118949b5e25fSSatish Balay        MatGetBlockSize_SeqSBAIJ,
119049b5e25fSSatish Balay        MatGetRowIJ_SeqSBAIJ,
119149b5e25fSSatish Balay        MatRestoreRowIJ_SeqSBAIJ,
119249b5e25fSSatish Balay        0,
119349b5e25fSSatish Balay        0,
119449b5e25fSSatish Balay        0,
119549b5e25fSSatish Balay        0,
119649b5e25fSSatish Balay        0,
119749b5e25fSSatish Balay        0,
119849b5e25fSSatish Balay        MatSetValuesBlocked_SeqSBAIJ,
119949b5e25fSSatish Balay        MatGetSubMatrix_SeqSBAIJ,
120049b5e25fSSatish Balay        0,
120149b5e25fSSatish Balay        0,
12028a124369SBarry Smith        MatGetPetscMaps_Petsc,
1203d959ec07SHong Zhang        0,
1204d959ec07SHong Zhang        0,
1205d959ec07SHong Zhang        0,
1206d959ec07SHong Zhang        0,
1207d959ec07SHong Zhang        0,
1208d959ec07SHong Zhang        0,
120924d5174aSHong Zhang        MatGetRowMax_SeqSBAIJ};
121049b5e25fSSatish Balay 
121149b5e25fSSatish Balay EXTERN_C_BEGIN
12124a2ae208SSatish Balay #undef __FUNCT__
12134a2ae208SSatish Balay #define __FUNCT__ "MatStoreValues_SeqSBAIJ"
121449b5e25fSSatish Balay int MatStoreValues_SeqSBAIJ(Mat mat)
121549b5e25fSSatish Balay {
12164afc71dfSHong Zhang   Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data;
1217c464158bSHong Zhang   int         nz = aij->i[mat->m]*aij->bs*aij->bs2;
121849b5e25fSSatish Balay   int         ierr;
121949b5e25fSSatish Balay 
122049b5e25fSSatish Balay   PetscFunctionBegin;
122149b5e25fSSatish Balay   if (aij->nonew != 1) {
122229bbc08cSBarry Smith     SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
122349b5e25fSSatish Balay   }
122449b5e25fSSatish Balay 
122549b5e25fSSatish Balay   /* allocate space for values if not already there */
122649b5e25fSSatish Balay   if (!aij->saved_values) {
122787828ca2SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&aij->saved_values);CHKERRQ(ierr);
122849b5e25fSSatish Balay   }
122949b5e25fSSatish Balay 
123049b5e25fSSatish Balay   /* copy values over */
123187828ca2SBarry Smith   ierr = PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(PetscScalar));CHKERRQ(ierr);
123249b5e25fSSatish Balay   PetscFunctionReturn(0);
123349b5e25fSSatish Balay }
123449b5e25fSSatish Balay EXTERN_C_END
123549b5e25fSSatish Balay 
123649b5e25fSSatish Balay EXTERN_C_BEGIN
12374a2ae208SSatish Balay #undef __FUNCT__
12384a2ae208SSatish Balay #define __FUNCT__ "MatRetrieveValues_SeqSBAIJ"
123949b5e25fSSatish Balay int MatRetrieveValues_SeqSBAIJ(Mat mat)
124049b5e25fSSatish Balay {
12414afc71dfSHong Zhang   Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data;
1242c464158bSHong Zhang   int         nz = aij->i[mat->m]*aij->bs*aij->bs2,ierr;
124349b5e25fSSatish Balay 
124449b5e25fSSatish Balay   PetscFunctionBegin;
124549b5e25fSSatish Balay   if (aij->nonew != 1) {
124629bbc08cSBarry Smith     SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
124749b5e25fSSatish Balay   }
124849b5e25fSSatish Balay   if (!aij->saved_values) {
124929bbc08cSBarry Smith     SETERRQ(1,"Must call MatStoreValues(A);first");
125049b5e25fSSatish Balay   }
125149b5e25fSSatish Balay 
125249b5e25fSSatish Balay   /* copy values over */
125387828ca2SBarry Smith   ierr = PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(PetscScalar));CHKERRQ(ierr);
125449b5e25fSSatish Balay   PetscFunctionReturn(0);
125549b5e25fSSatish Balay }
125649b5e25fSSatish Balay EXTERN_C_END
125749b5e25fSSatish Balay 
12588549e402SHong Zhang EXTERN_C_BEGIN
12594a2ae208SSatish Balay #undef __FUNCT__
12604a2ae208SSatish Balay #define __FUNCT__ "MatCreate_SeqSBAIJ"
1261c464158bSHong Zhang int MatCreate_SeqSBAIJ(Mat B)
1262c464158bSHong Zhang {
1263c464158bSHong Zhang   Mat_SeqSBAIJ *b;
12640c955e93SHong Zhang   int          ierr,size;
1265c464158bSHong Zhang 
1266c464158bSHong Zhang   PetscFunctionBegin;
1267c464158bSHong Zhang   ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr);
1268c464158bSHong Zhang   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1");
1269273d9f13SBarry Smith   B->m = B->M = PetscMax(B->m,B->M);
1270273d9f13SBarry Smith   B->n = B->N = PetscMax(B->n,B->N);
1271c464158bSHong Zhang 
1272b0a32e0cSBarry Smith   ierr    = PetscNew(Mat_SeqSBAIJ,&b);CHKERRQ(ierr);
1273b0a32e0cSBarry Smith   B->data = (void*)b;
1274c464158bSHong Zhang   ierr    = PetscMemzero(b,sizeof(Mat_SeqSBAIJ));CHKERRQ(ierr);
1275c464158bSHong Zhang   ierr    = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
1276c464158bSHong Zhang   B->ops->destroy     = MatDestroy_SeqSBAIJ;
1277c464158bSHong Zhang   B->ops->view        = MatView_SeqSBAIJ;
1278c464158bSHong Zhang   B->factor           = 0;
1279c464158bSHong Zhang   B->lupivotthreshold = 1.0;
1280c464158bSHong Zhang   B->mapping          = 0;
1281c464158bSHong Zhang   b->row              = 0;
1282c464158bSHong Zhang   b->icol             = 0;
1283c464158bSHong Zhang   b->reallocs         = 0;
1284c464158bSHong Zhang   b->saved_values     = 0;
1285c464158bSHong Zhang 
12868a124369SBarry Smith   ierr = PetscMapCreateMPI(B->comm,B->m,B->M,&B->rmap);CHKERRQ(ierr);
12878a124369SBarry Smith   ierr = PetscMapCreateMPI(B->comm,B->n,B->N,&B->cmap);CHKERRQ(ierr);
1288c464158bSHong Zhang 
1289c464158bSHong Zhang   b->sorted           = PETSC_FALSE;
1290c464158bSHong Zhang   b->roworiented      = PETSC_TRUE;
1291c464158bSHong Zhang   b->nonew            = 0;
1292c464158bSHong Zhang   b->diag             = 0;
1293c464158bSHong Zhang   b->solve_work       = 0;
1294c464158bSHong Zhang   b->mult_work        = 0;
1295c464158bSHong Zhang   b->spptr            = 0;
1296c464158bSHong Zhang   b->keepzeroedrows   = PETSC_FALSE;
1297c464158bSHong Zhang 
1298c464158bSHong Zhang   b->inew             = 0;
1299c464158bSHong Zhang   b->jnew             = 0;
1300c464158bSHong Zhang   b->anew             = 0;
1301c464158bSHong Zhang   b->a2anew           = 0;
1302c464158bSHong Zhang   b->permute          = PETSC_FALSE;
1303c464158bSHong Zhang 
1304c464158bSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C",
1305c464158bSHong Zhang                                      "MatStoreValues_SeqSBAIJ",
1306c464158bSHong Zhang                                      (void*)MatStoreValues_SeqSBAIJ);CHKERRQ(ierr);
1307c464158bSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C",
1308c464158bSHong Zhang                                      "MatRetrieveValues_SeqSBAIJ",
1309c464158bSHong Zhang                                      (void*)MatRetrieveValues_SeqSBAIJ);CHKERRQ(ierr);
1310c464158bSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqSBAIJSetColumnIndices_C",
1311c464158bSHong Zhang                                      "MatSeqSBAIJSetColumnIndices_SeqSBAIJ",
1312c464158bSHong Zhang                                      (void*)MatSeqSBAIJSetColumnIndices_SeqSBAIJ);CHKERRQ(ierr);
1313c464158bSHong Zhang   PetscFunctionReturn(0);
1314c464158bSHong Zhang }
13158549e402SHong Zhang EXTERN_C_END
1316c464158bSHong Zhang 
13174a2ae208SSatish Balay #undef __FUNCT__
13184a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetPreallocation"
131949b5e25fSSatish Balay /*@C
1320273d9f13SBarry Smith    MatSeqSBAIJSetPreallocation - Creates a sparse symmetric matrix in block AIJ (block
132149b5e25fSSatish Balay    compressed row) format.  For good matrix assembly performance the
132249b5e25fSSatish Balay    user should preallocate the matrix storage by setting the parameter nz
132349b5e25fSSatish Balay    (or the array nnz).  By setting these parameters accurately, performance
132449b5e25fSSatish Balay    during matrix assembly can be increased by more than a factor of 50.
132549b5e25fSSatish Balay 
1326c464158bSHong Zhang    Collective on Mat
132749b5e25fSSatish Balay 
132849b5e25fSSatish Balay    Input Parameters:
1329c464158bSHong Zhang +  A - the symmetric matrix
133049b5e25fSSatish Balay .  bs - size of block
133149b5e25fSSatish Balay .  nz - number of block nonzeros per block row (same for all rows)
133249b5e25fSSatish Balay -  nnz - array containing the number of block nonzeros in the various block rows
133349b5e25fSSatish Balay          (possibly different for each block row) or PETSC_NULL
133449b5e25fSSatish Balay 
133549b5e25fSSatish Balay    Options Database Keys:
133649b5e25fSSatish Balay .   -mat_no_unroll - uses code that does not unroll the loops in the
133749b5e25fSSatish Balay                      block calculations (much slower)
133849b5e25fSSatish Balay .    -mat_block_size - size of the blocks to use
133949b5e25fSSatish Balay 
134049b5e25fSSatish Balay    Level: intermediate
134149b5e25fSSatish Balay 
134249b5e25fSSatish Balay    Notes:
1343c464158bSHong Zhang    The block AIJ format is compatible with standard Fortran 77
134449b5e25fSSatish Balay    storage.  That is, the stored row and column indices can begin at
134549b5e25fSSatish Balay    either one (as in Fortran) or zero.  See the users' manual for details.
134649b5e25fSSatish Balay 
134749b5e25fSSatish Balay    Specify the preallocated storage with either nz or nnz (not both).
134849b5e25fSSatish Balay    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
134949b5e25fSSatish Balay    allocation.  For additional details, see the users manual chapter on
135049b5e25fSSatish Balay    matrices.
135149b5e25fSSatish Balay 
1352a209d233SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ()
135349b5e25fSSatish Balay @*/
1354273d9f13SBarry Smith int MatSeqSBAIJSetPreallocation(Mat B,int bs,int nz,int *nnz)
135549b5e25fSSatish Balay {
1356c464158bSHong Zhang   Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ*)B->data;
13570c955e93SHong Zhang   int          i,len,ierr,mbs,bs2;
135849b5e25fSSatish Balay   PetscTruth   flg;
13594afc71dfSHong Zhang   int          s_nz;
136049b5e25fSSatish Balay 
136149b5e25fSSatish Balay   PetscFunctionBegin;
1362273d9f13SBarry Smith   B->preallocated = PETSC_TRUE;
1363b0a32e0cSBarry Smith   ierr = PetscOptionsGetInt(PETSC_NULL,"-mat_block_size",&bs,PETSC_NULL);CHKERRQ(ierr);
1364c464158bSHong Zhang   mbs  = B->m/bs;
136549b5e25fSSatish Balay   bs2  = bs*bs;
136649b5e25fSSatish Balay 
1367c464158bSHong Zhang   if (mbs*bs != B->m) {
136829bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_SIZ,"Number rows, cols must be divisible by blocksize");
136949b5e25fSSatish Balay   }
137049b5e25fSSatish Balay 
1371435da068SBarry Smith   if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 3;
1372435da068SBarry Smith   if (nz < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %d",nz);
137349b5e25fSSatish Balay   if (nnz) {
137449b5e25fSSatish Balay     for (i=0; i<mbs; i++) {
137529bbc08cSBarry Smith       if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %d value %d",i,nnz[i]);
137680fe4e49SBarry 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);
137749b5e25fSSatish Balay     }
137849b5e25fSSatish Balay   }
137949b5e25fSSatish Balay 
1380b0a32e0cSBarry Smith   ierr    = PetscOptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg);CHKERRQ(ierr);
138149b5e25fSSatish Balay   if (!flg) {
138249b5e25fSSatish Balay     switch (bs) {
138349b5e25fSSatish Balay     case 1:
13844ccecd49SHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1;
138549b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_1;
138649b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_1;
138749b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_1;
138849b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_1;
138949b5e25fSSatish Balay       break;
139049b5e25fSSatish Balay     case 2:
13914ccecd49SHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2;
139249b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_2;
139349b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_2;
139449b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_2;
139549b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_2;
139649b5e25fSSatish Balay       break;
139749b5e25fSSatish Balay     case 3:
13985f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3;
139949b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_3;
140049b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_3;
140149b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_3;
140249b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_3;
140349b5e25fSSatish Balay       break;
140449b5e25fSSatish Balay     case 4:
14055f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4;
140649b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_4;
140749b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_4;
140849b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_4;
140949b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_4;
141049b5e25fSSatish Balay       break;
141149b5e25fSSatish Balay     case 5:
14125f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5;
141349b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_5;
141449b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_5;
141549b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_5;
141649b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_5;
141749b5e25fSSatish Balay       break;
141849b5e25fSSatish Balay     case 6:
14195f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6;
142049b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_6;
142149b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_6;
142249b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_6;
142349b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_6;
142449b5e25fSSatish Balay       break;
142549b5e25fSSatish Balay     case 7:
1426de53e5efSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7;
142749b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_7;
142849b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_7;
1429de53e5efSHong Zhang       B->ops->mult            = MatMult_SeqSBAIJ_7;
143049b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_7;
143149b5e25fSSatish Balay       break;
143249b5e25fSSatish Balay     }
143349b5e25fSSatish Balay   }
143449b5e25fSSatish Balay 
143549b5e25fSSatish Balay   b->mbs = mbs;
14364afc71dfSHong Zhang   b->nbs = mbs;
1437b0a32e0cSBarry Smith   ierr   = PetscMalloc((mbs+1)*sizeof(int),&b->imax);CHKERRQ(ierr);
143849b5e25fSSatish Balay   if (!nnz) {
1439435da068SBarry Smith     if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5;
144049b5e25fSSatish Balay     else if (nz <= 0)        nz = 1;
144149b5e25fSSatish Balay     for (i=0; i<mbs; i++) {
14428cef66ccSHong Zhang       b->imax[i] = nz;
144349b5e25fSSatish Balay     }
1444153ea458SHong Zhang     nz = nz*mbs; /* total nz */
144549b5e25fSSatish Balay   } else {
144649b5e25fSSatish Balay     nz = 0;
14478cef66ccSHong Zhang     for (i=0; i<mbs; i++) {b->imax[i] = nnz[i]; nz += nnz[i];}
144849b5e25fSSatish Balay   }
14498cef66ccSHong Zhang   /* s_nz=(nz+mbs)/2; */ /* total diagonal and superdiagonal nonzero blocks */
14508cef66ccSHong Zhang   s_nz = nz;
145149b5e25fSSatish Balay 
145249b5e25fSSatish Balay   /* allocate the matrix space */
1453c464158bSHong Zhang   len  = s_nz*sizeof(int) + s_nz*bs2*sizeof(MatScalar) + (B->m+1)*sizeof(int);
145482502324SSatish Balay   ierr = PetscMalloc(len,&b->a);CHKERRQ(ierr);
145549b5e25fSSatish Balay   ierr = PetscMemzero(b->a,s_nz*bs2*sizeof(MatScalar));CHKERRQ(ierr);
145649b5e25fSSatish Balay   b->j = (int*)(b->a + s_nz*bs2);
145749b5e25fSSatish Balay   ierr = PetscMemzero(b->j,s_nz*sizeof(int));CHKERRQ(ierr);
145849b5e25fSSatish Balay   b->i = b->j + s_nz;
145949b5e25fSSatish Balay   b->singlemalloc = PETSC_TRUE;
146049b5e25fSSatish Balay 
146149b5e25fSSatish Balay   /* pointer to beginning of each row */
146249b5e25fSSatish Balay   b->i[0] = 0;
146349b5e25fSSatish Balay   for (i=1; i<mbs+1; i++) {
146449b5e25fSSatish Balay     b->i[i] = b->i[i-1] + b->imax[i-1];
146549b5e25fSSatish Balay   }
146649b5e25fSSatish Balay 
146749b5e25fSSatish Balay   /* b->ilen will count nonzeros in each block row so far. */
14685033bc48SSatish Balay   ierr = PetscMalloc((mbs+1)*sizeof(int),&b->ilen);CHKERRQ(ierr);
1469b0a32e0cSBarry Smith   PetscLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqSBAIJ));
147049b5e25fSSatish Balay   for (i=0; i<mbs; i++) { b->ilen[i] = 0;}
147149b5e25fSSatish Balay 
147249b5e25fSSatish Balay   b->bs               = bs;
147349b5e25fSSatish Balay   b->bs2              = bs2;
147449b5e25fSSatish Balay   b->s_nz             = 0;
147549b5e25fSSatish Balay   b->s_maxnz          = s_nz*bs2;
1476153ea458SHong Zhang 
147716cdd363SHong Zhang   b->inew             = 0;
147816cdd363SHong Zhang   b->jnew             = 0;
147916cdd363SHong Zhang   b->anew             = 0;
148016cdd363SHong Zhang   b->a2anew           = 0;
14811a3463dfSHong Zhang   b->permute          = PETSC_FALSE;
1482c464158bSHong Zhang   PetscFunctionReturn(0);
1483c464158bSHong Zhang }
1484153ea458SHong Zhang 
148549b5e25fSSatish Balay 
14864a2ae208SSatish Balay #undef __FUNCT__
14874a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqSBAIJ"
1488c464158bSHong Zhang /*@C
1489c464158bSHong Zhang    MatCreateSeqSBAIJ - Creates a sparse symmetric matrix in block AIJ (block
1490c464158bSHong Zhang    compressed row) format.  For good matrix assembly performance the
1491c464158bSHong Zhang    user should preallocate the matrix storage by setting the parameter nz
1492c464158bSHong Zhang    (or the array nnz).  By setting these parameters accurately, performance
1493c464158bSHong Zhang    during matrix assembly can be increased by more than a factor of 50.
149449b5e25fSSatish Balay 
1495c464158bSHong Zhang    Collective on MPI_Comm
1496c464158bSHong Zhang 
1497c464158bSHong Zhang    Input Parameters:
1498c464158bSHong Zhang +  comm - MPI communicator, set to PETSC_COMM_SELF
1499c464158bSHong Zhang .  bs - size of block
1500c464158bSHong Zhang .  m - number of rows, or number of columns
1501c464158bSHong Zhang .  nz - number of block nonzeros per block row (same for all rows)
1502c464158bSHong Zhang -  nnz - array containing the number of block nonzeros in the various block rows
1503c464158bSHong Zhang          (possibly different for each block row) or PETSC_NULL
1504c464158bSHong Zhang 
1505c464158bSHong Zhang    Output Parameter:
1506c464158bSHong Zhang .  A - the symmetric matrix
1507c464158bSHong Zhang 
1508c464158bSHong Zhang    Options Database Keys:
1509c464158bSHong Zhang .   -mat_no_unroll - uses code that does not unroll the loops in the
1510c464158bSHong Zhang                      block calculations (much slower)
1511c464158bSHong Zhang .    -mat_block_size - size of the blocks to use
1512c464158bSHong Zhang 
1513c464158bSHong Zhang    Level: intermediate
1514c464158bSHong Zhang 
1515c464158bSHong Zhang    Notes:
1516c464158bSHong Zhang    The block AIJ format is compatible with standard Fortran 77
1517c464158bSHong Zhang    storage.  That is, the stored row and column indices can begin at
1518c464158bSHong Zhang    either one (as in Fortran) or zero.  See the users' manual for details.
1519c464158bSHong Zhang 
1520c464158bSHong Zhang    Specify the preallocated storage with either nz or nnz (not both).
1521c464158bSHong Zhang    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
1522c464158bSHong Zhang    allocation.  For additional details, see the users manual chapter on
1523c464158bSHong Zhang    matrices.
1524c464158bSHong Zhang 
1525c464158bSHong Zhang .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ()
1526c464158bSHong Zhang @*/
1527c464158bSHong Zhang int MatCreateSeqSBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,int *nnz,Mat *A)
1528c464158bSHong Zhang {
1529c464158bSHong Zhang   int ierr;
1530c464158bSHong Zhang 
1531c464158bSHong Zhang   PetscFunctionBegin;
1532c464158bSHong Zhang   ierr = MatCreate(comm,m,n,m,n,A);CHKERRQ(ierr);
1533c464158bSHong Zhang   ierr = MatSetType(*A,MATSEQSBAIJ);CHKERRQ(ierr);
1534273d9f13SBarry Smith   ierr = MatSeqSBAIJSetPreallocation(*A,bs,nz,nnz);CHKERRQ(ierr);
153549b5e25fSSatish Balay   PetscFunctionReturn(0);
153649b5e25fSSatish Balay }
153749b5e25fSSatish Balay 
15384a2ae208SSatish Balay #undef __FUNCT__
15394a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqSBAIJ"
154049b5e25fSSatish Balay int MatDuplicate_SeqSBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B)
154149b5e25fSSatish Balay {
154249b5e25fSSatish Balay   Mat         C;
154349b5e25fSSatish Balay   Mat_SeqSBAIJ *c,*a = (Mat_SeqSBAIJ*)A->data;
154449b5e25fSSatish Balay   int         i,len,mbs = a->mbs,nz = a->s_nz,bs2 =a->bs2,ierr;
154549b5e25fSSatish Balay 
154649b5e25fSSatish Balay   PetscFunctionBegin;
154729bbc08cSBarry Smith   if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,"Corrupt matrix");
154849b5e25fSSatish Balay 
154949b5e25fSSatish Balay   *B = 0;
1550692f9cbeSHong Zhang   ierr = MatCreate(A->comm,A->m,A->n,A->m,A->n,&C);CHKERRQ(ierr);
1551692f9cbeSHong Zhang   ierr = MatSetType(C,MATSEQSBAIJ);CHKERRQ(ierr);
1552692f9cbeSHong Zhang   c    = (Mat_SeqSBAIJ*)C->data;
1553692f9cbeSHong Zhang 
155449b5e25fSSatish Balay   ierr              = PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
1555273d9f13SBarry Smith   C->preallocated   = PETSC_TRUE;
155649b5e25fSSatish Balay   C->factor         = A->factor;
155749b5e25fSSatish Balay   c->row            = 0;
155849b5e25fSSatish Balay   c->icol           = 0;
155949b5e25fSSatish Balay   c->saved_values   = 0;
156049b5e25fSSatish Balay   c->keepzeroedrows = a->keepzeroedrows;
156149b5e25fSSatish Balay   C->assembled      = PETSC_TRUE;
156249b5e25fSSatish Balay 
156349b5e25fSSatish Balay   c->bs         = a->bs;
156449b5e25fSSatish Balay   c->bs2        = a->bs2;
156549b5e25fSSatish Balay   c->mbs        = a->mbs;
156649b5e25fSSatish Balay   c->nbs        = a->nbs;
156749b5e25fSSatish Balay 
1568b0a32e0cSBarry Smith   ierr = PetscMalloc((mbs+1)*sizeof(int),&c->imax);CHKERRQ(ierr);
1569b0a32e0cSBarry Smith   ierr = PetscMalloc((mbs+1)*sizeof(int),&c->ilen);CHKERRQ(ierr);
157049b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
157149b5e25fSSatish Balay     c->imax[i] = a->imax[i];
157249b5e25fSSatish Balay     c->ilen[i] = a->ilen[i];
157349b5e25fSSatish Balay   }
157449b5e25fSSatish Balay 
157549b5e25fSSatish Balay   /* allocate the matrix space */
157649b5e25fSSatish Balay   c->singlemalloc = PETSC_TRUE;
157749b5e25fSSatish Balay   len  = (mbs+1)*sizeof(int) + nz*(bs2*sizeof(MatScalar) + sizeof(int));
157882502324SSatish Balay   ierr = PetscMalloc(len,&c->a);CHKERRQ(ierr);
157949b5e25fSSatish Balay   c->j = (int*)(c->a + nz*bs2);
158049b5e25fSSatish Balay   c->i = c->j + nz;
158149b5e25fSSatish Balay   ierr = PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int));CHKERRQ(ierr);
158249b5e25fSSatish Balay   if (mbs > 0) {
158349b5e25fSSatish Balay     ierr = PetscMemcpy(c->j,a->j,nz*sizeof(int));CHKERRQ(ierr);
158449b5e25fSSatish Balay     if (cpvalues == MAT_COPY_VALUES) {
158549b5e25fSSatish Balay       ierr = PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr);
158649b5e25fSSatish Balay     } else {
158749b5e25fSSatish Balay       ierr = PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr);
158849b5e25fSSatish Balay     }
158949b5e25fSSatish Balay   }
159049b5e25fSSatish Balay 
1591b0a32e0cSBarry Smith   PetscLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqSBAIJ));
159249b5e25fSSatish Balay   c->sorted      = a->sorted;
159349b5e25fSSatish Balay   c->roworiented = a->roworiented;
159449b5e25fSSatish Balay   c->nonew       = a->nonew;
159549b5e25fSSatish Balay 
159649b5e25fSSatish Balay   if (a->diag) {
1597b0a32e0cSBarry Smith     ierr = PetscMalloc((mbs+1)*sizeof(int),&c->diag);CHKERRQ(ierr);
1598b0a32e0cSBarry Smith     PetscLogObjectMemory(C,(mbs+1)*sizeof(int));
159949b5e25fSSatish Balay     for (i=0; i<mbs; i++) {
160049b5e25fSSatish Balay       c->diag[i] = a->diag[i];
160149b5e25fSSatish Balay     }
160249b5e25fSSatish Balay   } else c->diag        = 0;
160349b5e25fSSatish Balay   c->s_nz               = a->s_nz;
160449b5e25fSSatish Balay   c->s_maxnz            = a->s_maxnz;
160549b5e25fSSatish Balay   c->solve_work         = 0;
160649b5e25fSSatish Balay   c->spptr              = 0;      /* Dangerous -I'm throwing away a->spptr */
160749b5e25fSSatish Balay   c->mult_work          = 0;
160849b5e25fSSatish Balay   *B = C;
1609b0a32e0cSBarry Smith   ierr = PetscFListDuplicate(A->qlist,&C->qlist);CHKERRQ(ierr);
161049b5e25fSSatish Balay   PetscFunctionReturn(0);
161149b5e25fSSatish Balay }
161249b5e25fSSatish Balay 
16138549e402SHong Zhang EXTERN_C_BEGIN
16144a2ae208SSatish Balay #undef __FUNCT__
16154a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqSBAIJ"
1616b0a32e0cSBarry Smith int MatLoad_SeqSBAIJ(PetscViewer viewer,MatType type,Mat *A)
161749b5e25fSSatish Balay {
161849b5e25fSSatish Balay   Mat_SeqSBAIJ *a;
161949b5e25fSSatish Balay   Mat          B;
162049b5e25fSSatish Balay   int          i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1;
1621574b2666SHong Zhang   int          *mask,mbs,*jj,j,rowcount,nzcount,k,*browlengths,*s_browlengths,maskcount;
162249b5e25fSSatish Balay   int          kmax,jcount,block,idx,point,nzcountb,extra_rows;
162349b5e25fSSatish Balay   int          *masked,nmask,tmp,bs2,ishift;
162487828ca2SBarry Smith   PetscScalar  *aa;
162549b5e25fSSatish Balay   MPI_Comm     comm = ((PetscObject)viewer)->comm;
162649b5e25fSSatish Balay 
162749b5e25fSSatish Balay   PetscFunctionBegin;
1628b0a32e0cSBarry Smith   ierr = PetscOptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,PETSC_NULL);CHKERRQ(ierr);
162949b5e25fSSatish Balay   bs2  = bs*bs;
163049b5e25fSSatish Balay 
163149b5e25fSSatish Balay   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
163229bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor");
1633b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
163449b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr);
1635552e946dSBarry Smith   if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not Mat object");
163649b5e25fSSatish Balay   M = header[1]; N = header[2]; nz = header[3];
163749b5e25fSSatish Balay 
163849b5e25fSSatish Balay   if (header[3] < 0) {
163929bbc08cSBarry Smith     SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as SeqSBAIJ");
164049b5e25fSSatish Balay   }
164149b5e25fSSatish Balay 
164229bbc08cSBarry Smith   if (M != N) SETERRQ(PETSC_ERR_SUP,"Can only do square matrices");
164349b5e25fSSatish Balay 
164449b5e25fSSatish Balay   /*
164549b5e25fSSatish Balay      This code adds extra rows to make sure the number of rows is
164649b5e25fSSatish Balay     divisible by the blocksize
164749b5e25fSSatish Balay   */
164849b5e25fSSatish Balay   mbs        = M/bs;
164949b5e25fSSatish Balay   extra_rows = bs - M + bs*(mbs);
165049b5e25fSSatish Balay   if (extra_rows == bs) extra_rows = 0;
165149b5e25fSSatish Balay   else                  mbs++;
165249b5e25fSSatish Balay   if (extra_rows) {
1653b0a32e0cSBarry Smith     PetscLogInfo(0,"MatLoad_SeqSBAIJ:Padding loaded matrix to match blocksize\n");
165449b5e25fSSatish Balay   }
165549b5e25fSSatish Balay 
165649b5e25fSSatish Balay   /* read in row lengths */
1657b0a32e0cSBarry Smith   ierr = PetscMalloc((M+extra_rows)*sizeof(int),&rowlengths);CHKERRQ(ierr);
165849b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
165949b5e25fSSatish Balay   for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1;
166049b5e25fSSatish Balay 
166149b5e25fSSatish Balay   /* read in column indices */
1662b0a32e0cSBarry Smith   ierr = PetscMalloc((nz+extra_rows)*sizeof(int),&jj);CHKERRQ(ierr);
166349b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT);CHKERRQ(ierr);
166449b5e25fSSatish Balay   for (i=0; i<extra_rows; i++) jj[nz+i] = M+i;
166549b5e25fSSatish Balay 
166649b5e25fSSatish Balay   /* loop over row lengths determining block row lengths */
166782502324SSatish Balay   ierr     = PetscMalloc(mbs*sizeof(int),&browlengths);CHKERRQ(ierr);
166882502324SSatish Balay   ierr     = PetscMalloc(mbs*sizeof(int),&s_browlengths);CHKERRQ(ierr);
1669574b2666SHong Zhang   ierr     = PetscMemzero(s_browlengths,mbs*sizeof(int));CHKERRQ(ierr);
167082502324SSatish Balay   ierr     = PetscMalloc(2*mbs*sizeof(int),&mask);CHKERRQ(ierr);
167149b5e25fSSatish Balay   ierr     = PetscMemzero(mask,mbs*sizeof(int));CHKERRQ(ierr);
167249b5e25fSSatish Balay   masked   = mask + mbs;
167349b5e25fSSatish Balay   rowcount = 0; nzcount = 0;
167449b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
167549b5e25fSSatish Balay     nmask = 0;
167649b5e25fSSatish Balay     for (j=0; j<bs; j++) {
167749b5e25fSSatish Balay       kmax = rowlengths[rowcount];
167849b5e25fSSatish Balay       for (k=0; k<kmax; k++) {
16792d703238SHong Zhang         tmp = jj[nzcount++]/bs;   /* block col. index */
168003630b6eSHong Zhang         if (!mask[tmp] && tmp >= i) {masked[nmask++] = tmp; mask[tmp] = 1;}
168149b5e25fSSatish Balay       }
168249b5e25fSSatish Balay       rowcount++;
168349b5e25fSSatish Balay     }
1684574b2666SHong Zhang     s_browlengths[i] += nmask;
1685574b2666SHong Zhang     browlengths[i]    = 2*s_browlengths[i];
1686574b2666SHong Zhang 
168749b5e25fSSatish Balay     /* zero out the mask elements we set */
168849b5e25fSSatish Balay     for (j=0; j<nmask; j++) mask[masked[j]] = 0;
168949b5e25fSSatish Balay   }
169049b5e25fSSatish Balay 
169149b5e25fSSatish Balay   /* create our matrix */
169249b5e25fSSatish Balay   ierr = MatCreateSeqSBAIJ(comm,bs,M+extra_rows,N+extra_rows,0,browlengths,A);CHKERRQ(ierr);
169349b5e25fSSatish Balay   B = *A;
169449b5e25fSSatish Balay   a = (Mat_SeqSBAIJ*)B->data;
169549b5e25fSSatish Balay 
169649b5e25fSSatish Balay   /* set matrix "i" values */
169749b5e25fSSatish Balay   a->i[0] = 0;
169849b5e25fSSatish Balay   for (i=1; i<= mbs; i++) {
1699574b2666SHong Zhang     a->i[i]      = a->i[i-1] + s_browlengths[i-1];
1700574b2666SHong Zhang     a->ilen[i-1] = s_browlengths[i-1];
170149b5e25fSSatish Balay   }
170249b5e25fSSatish Balay   a->s_nz         = 0;
1703574b2666SHong Zhang   for (i=0; i<mbs; i++) a->s_nz += s_browlengths[i];
170449b5e25fSSatish Balay 
170549b5e25fSSatish Balay   /* read in nonzero values */
170687828ca2SBarry Smith   ierr = PetscMalloc((nz+extra_rows)*sizeof(PetscScalar),&aa);CHKERRQ(ierr);
170749b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);CHKERRQ(ierr);
170849b5e25fSSatish Balay   for (i=0; i<extra_rows; i++) aa[nz+i] = 1.0;
170949b5e25fSSatish Balay 
171049b5e25fSSatish Balay   /* set "a" and "j" values into matrix */
171149b5e25fSSatish Balay   nzcount = 0; jcount = 0;
171249b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
171349b5e25fSSatish Balay     nzcountb = nzcount;
171449b5e25fSSatish Balay     nmask    = 0;
171549b5e25fSSatish Balay     for (j=0; j<bs; j++) {
171649b5e25fSSatish Balay       kmax = rowlengths[i*bs+j];
171749b5e25fSSatish Balay       for (k=0; k<kmax; k++) {
17182d703238SHong Zhang         tmp = jj[nzcount++]/bs; /* block col. index */
171903630b6eSHong Zhang         if (!mask[tmp] && tmp >= i) { masked[nmask++] = tmp; mask[tmp] = 1;}
17202d703238SHong Zhang       }
17212d703238SHong Zhang     }
17222d703238SHong Zhang     /* sort the masked values */
17232d703238SHong Zhang     ierr = PetscSortInt(nmask,masked);CHKERRQ(ierr);
17242d703238SHong Zhang 
17252d703238SHong Zhang     /* set "j" values into matrix */
17262d703238SHong Zhang     maskcount = 1;
17272d703238SHong Zhang     for (j=0; j<nmask; j++) {
172849b5e25fSSatish Balay       a->j[jcount++]  = masked[j];
172949b5e25fSSatish Balay       mask[masked[j]] = maskcount++;
173049b5e25fSSatish Balay     }
1731574b2666SHong Zhang 
173249b5e25fSSatish Balay     /* set "a" values into matrix */
173349b5e25fSSatish Balay     ishift = bs2*a->i[i];
173449b5e25fSSatish Balay     for (j=0; j<bs; j++) {
173549b5e25fSSatish Balay       kmax = rowlengths[i*bs+j];
173649b5e25fSSatish Balay       for (k=0; k<kmax; k++) {
1737574b2666SHong Zhang         tmp       = jj[nzcountb]/bs ; /* block col. index */
1738574b2666SHong Zhang         if (tmp >= i){
173949b5e25fSSatish Balay           block     = mask[tmp] - 1;
174049b5e25fSSatish Balay           point     = jj[nzcountb] - bs*tmp;
174149b5e25fSSatish Balay           idx       = ishift + bs2*block + j + bs*point;
1742574b2666SHong Zhang           a->a[idx] = aa[nzcountb];
1743574b2666SHong Zhang         }
1744574b2666SHong Zhang         nzcountb++;
174549b5e25fSSatish Balay       }
174649b5e25fSSatish Balay     }
174749b5e25fSSatish Balay     /* zero out the mask elements we set */
174849b5e25fSSatish Balay     for (j=0; j<nmask; j++) mask[masked[j]] = 0;
174949b5e25fSSatish Balay   }
175029bbc08cSBarry Smith   if (jcount != a->s_nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Bad binary matrix");
175149b5e25fSSatish Balay 
175249b5e25fSSatish Balay   ierr = PetscFree(rowlengths);CHKERRQ(ierr);
175349b5e25fSSatish Balay   ierr = PetscFree(browlengths);CHKERRQ(ierr);
1754574b2666SHong Zhang   ierr = PetscFree(s_browlengths);CHKERRQ(ierr);
175549b5e25fSSatish Balay   ierr = PetscFree(aa);CHKERRQ(ierr);
175649b5e25fSSatish Balay   ierr = PetscFree(jj);CHKERRQ(ierr);
175749b5e25fSSatish Balay   ierr = PetscFree(mask);CHKERRQ(ierr);
175849b5e25fSSatish Balay 
175949b5e25fSSatish Balay   B->assembled = PETSC_TRUE;
176049b5e25fSSatish Balay   ierr = MatView_Private(B);CHKERRQ(ierr);
176149b5e25fSSatish Balay   PetscFunctionReturn(0);
176249b5e25fSSatish Balay }
17638549e402SHong Zhang EXTERN_C_END
1764574b2666SHong Zhang 
1765d06b337dSHong Zhang #undef __FUNCT__
1766d06b337dSHong Zhang #define __FUNCT__ "MatRelax_SeqSBAIJ"
1767c14dc6b6SHong Zhang int MatRelax_SeqSBAIJ(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,int its,int lits,Vec xx)
1768d06b337dSHong Zhang {
1769d06b337dSHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
1770d06b337dSHong Zhang   MatScalar    *aa=a->a,*v,*v1;
1771d06b337dSHong Zhang   PetscScalar  *x,*b,*t,sum,d;
1772d06b337dSHong Zhang   int          m=a->mbs,bs=a->bs,*ai=a->i,*aj=a->j,ierr;
1773d06b337dSHong Zhang   int          nz,nz1,*vj,*vj1,i;
1774d06b337dSHong Zhang 
1775d06b337dSHong Zhang   PetscFunctionBegin;
1776b965ef7fSBarry Smith   its = its*lits;
177791723122SBarry Smith   if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %d and local its %d both positive",its,lits);
1778d06b337dSHong Zhang 
1779d06b337dSHong Zhang   if (bs > 1)
1780d06b337dSHong Zhang     SETERRQ(PETSC_ERR_SUP,"SSOR for block size > 1 is not yet implemented");
1781d06b337dSHong Zhang 
1782d06b337dSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1783d06b337dSHong Zhang   if (xx != bb) {
1784d06b337dSHong Zhang     ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
1785d06b337dSHong Zhang   } else {
1786d06b337dSHong Zhang     b = x;
1787d06b337dSHong Zhang   }
1788d06b337dSHong Zhang 
1789d06b337dSHong Zhang   ierr = PetscMalloc(m*sizeof(PetscScalar),&t);CHKERRQ(ierr);
1790d06b337dSHong Zhang 
1791d06b337dSHong Zhang   if (flag & SOR_ZERO_INITIAL_GUESS) {
17922798e883SHong Zhang     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
1793d06b337dSHong Zhang       for (i=0; i<m; i++)
1794d06b337dSHong Zhang         t[i] = b[i];
1795d06b337dSHong Zhang 
1796d06b337dSHong Zhang       for (i=0; i<m; i++){
179744706875SHong Zhang         d  = *(aa + ai[i]);  /* diag[i] */
1798d06b337dSHong Zhang         v  = aa + ai[i] + 1;
1799d06b337dSHong Zhang         vj = aj + ai[i] + 1;
1800d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1;
1801d06b337dSHong Zhang         x[i] = omega*t[i]/d;
1802d06b337dSHong Zhang         while (nz--) t[*vj++] -= x[i]*(*v++); /* update rhs */
180344706875SHong Zhang         PetscLogFlops(2*nz-1);
1804d06b337dSHong Zhang       }
1805d06b337dSHong Zhang     }
1806d06b337dSHong Zhang 
18072798e883SHong Zhang     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
1808d06b337dSHong Zhang       for (i=0; i<m; i++)
1809d06b337dSHong Zhang         t[i] = b[i];
1810d06b337dSHong Zhang 
1811d06b337dSHong Zhang       for (i=0; i<m-1; i++){  /* update rhs */
1812d06b337dSHong Zhang         v  = aa + ai[i] + 1;
1813d06b337dSHong Zhang         vj = aj + ai[i] + 1;
1814d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1;
1815d06b337dSHong Zhang         while (nz--) t[*vj++] -= x[i]*(*v++);
181644706875SHong Zhang         PetscLogFlops(2*nz-1);
1817d06b337dSHong Zhang       }
1818d06b337dSHong Zhang       for (i=m-1; i>=0; i--){
1819d06b337dSHong Zhang         d  = *(aa + ai[i]);
1820d06b337dSHong Zhang         v  = aa + ai[i] + 1;
1821d06b337dSHong Zhang         vj = aj + ai[i] + 1;
1822d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1;
1823d06b337dSHong Zhang         sum = t[i];
1824d06b337dSHong Zhang         while (nz--) sum -= x[*vj++]*(*v++);
182544706875SHong Zhang         PetscLogFlops(2*nz-1);
1826d06b337dSHong Zhang         x[i] =   (1-omega)*x[i] + omega*sum/d;
1827d06b337dSHong Zhang       }
1828d06b337dSHong Zhang     }
1829d06b337dSHong Zhang     its--;
1830d06b337dSHong Zhang   }
1831d06b337dSHong Zhang 
1832d06b337dSHong Zhang   while (its--) {
183344706875SHong Zhang     /*
183444706875SHong Zhang        forward sweep:
183544706875SHong Zhang        for i=0,...,m-1:
183644706875SHong Zhang          sum[i] = (b[i] - U(i,:)x )/d[i];
183744706875SHong Zhang          x[i]   = (1-omega)x[i] + omega*sum[i];
183844706875SHong Zhang          b      = b - x[i]*U^T(i,:);
1839d06b337dSHong Zhang 
184044706875SHong Zhang     */
18412798e883SHong Zhang     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
1842d06b337dSHong Zhang       for (i=0; i<m; i++)
1843d06b337dSHong Zhang         t[i] = b[i];
1844d06b337dSHong Zhang 
1845d06b337dSHong Zhang       for (i=0; i<m; i++){
184644706875SHong Zhang         d  = *(aa + ai[i]);  /* diag[i] */
1847d06b337dSHong Zhang         v  = aa + ai[i] + 1; v1=v;
1848d06b337dSHong Zhang         vj = aj + ai[i] + 1; vj1=vj;
1849d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1; nz1=nz;
1850d06b337dSHong Zhang         sum = t[i];
1851d06b337dSHong Zhang         while (nz1--) sum -= (*v1++)*x[*vj1++];
1852d06b337dSHong Zhang         x[i] = (1-omega)*x[i] + omega*sum/d;
1853d06b337dSHong Zhang         while (nz--) t[*vj++] -= x[i]*(*v++);
185444706875SHong Zhang         PetscLogFlops(4*nz-2);
1855d06b337dSHong Zhang       }
1856d06b337dSHong Zhang     }
1857d06b337dSHong Zhang 
18582798e883SHong Zhang   if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
185944706875SHong Zhang       /*
186044706875SHong Zhang        backward sweep:
186144706875SHong Zhang        b = b - x[i]*U^T(i,:), i=0,...,n-2
186244706875SHong Zhang        for i=m-1,...,0:
186344706875SHong Zhang          sum[i] = (b[i] - U(i,:)x )/d[i];
186444706875SHong Zhang          x[i]   = (1-omega)x[i] + omega*sum[i];
186544706875SHong Zhang       */
1866d06b337dSHong Zhang       for (i=0; i<m; i++)
1867d06b337dSHong Zhang         t[i] = b[i];
1868d06b337dSHong Zhang 
1869d06b337dSHong Zhang       for (i=0; i<m-1; i++){  /* update rhs */
1870d06b337dSHong Zhang         v  = aa + ai[i] + 1;
1871d06b337dSHong Zhang         vj = aj + ai[i] + 1;
1872d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1;
1873d06b337dSHong Zhang         while (nz--) t[*vj++] -= x[i]*(*v++);
187444706875SHong Zhang         PetscLogFlops(2*nz-1);
1875d06b337dSHong Zhang       }
1876d06b337dSHong Zhang       for (i=m-1; i>=0; i--){
1877d06b337dSHong Zhang         d  = *(aa + ai[i]);
1878d06b337dSHong Zhang         v  = aa + ai[i] + 1;
1879d06b337dSHong Zhang         vj = aj + ai[i] + 1;
1880d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1;
1881d06b337dSHong Zhang         sum = t[i];
1882d06b337dSHong Zhang         while (nz--) sum -= x[*vj++]*(*v++);
188344706875SHong Zhang         PetscLogFlops(2*nz-1);
1884d06b337dSHong Zhang         x[i] =   (1-omega)*x[i] + omega*sum/d;
1885d06b337dSHong Zhang       }
1886d06b337dSHong Zhang     }
1887d06b337dSHong Zhang   }
1888d06b337dSHong Zhang 
1889d06b337dSHong Zhang   ierr = PetscFree(t); CHKERRQ(ierr);
1890d06b337dSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
1891d06b337dSHong Zhang   if (bb != xx) {
1892d06b337dSHong Zhang     ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
1893d06b337dSHong Zhang   }
18942798e883SHong Zhang 
1895d06b337dSHong Zhang   PetscFunctionReturn(0);
1896d06b337dSHong Zhang }
1897d06b337dSHong Zhang 
1898d06b337dSHong Zhang 
1899d06b337dSHong Zhang 
1900d06b337dSHong Zhang 
190149b5e25fSSatish Balay 
190249b5e25fSSatish Balay 
1903