xref: /petsc/src/mat/impls/sbaij/seq/sbaij.c (revision 0880e06233d9681dd0ed22d6e02e81e220065bc9)
173f4d377SMatthew Knepley /*$Id: sbaij.c,v 1.62 2001/08/07 03:03:01 balay Exp $*/
249b5e25fSSatish Balay 
349b5e25fSSatish Balay /*
4a1373b80SHong Zhang     Defines the basic matrix operations for the SBAIJ (compressed row)
549b5e25fSSatish Balay   matrix storage format.
649b5e25fSSatish Balay */
79e070d67SMatthew Knepley #include "src/mat/impls/baij/seq/baij.h"         /*I "petscmat.h" I*/
849b5e25fSSatish Balay #include "src/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++) {
28e11b0e14SHong 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;
38b424e231SHong Zhang   int          i,mbs = a->mbs,ierr;
3949b5e25fSSatish Balay 
4049b5e25fSSatish Balay   PetscFunctionBegin;
4149b5e25fSSatish Balay   if (a->diag) PetscFunctionReturn(0);
4249b5e25fSSatish Balay 
43b424e231SHong Zhang   ierr = PetscMalloc((mbs+1)*sizeof(int),&a->diag);CHKERRQ(ierr);
44b424e231SHong Zhang   PetscLogObjectMemory(A,(mbs+1)*sizeof(int));
45b424e231SHong Zhang   for (i=0; i<mbs; i++) a->diag[i] = a->i[i];
4649b5e25fSSatish Balay   PetscFunctionReturn(0);
4749b5e25fSSatish Balay }
4849b5e25fSSatish Balay 
49a1373b80SHong Zhang /* extern int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int,int**,int**); */
5049b5e25fSSatish Balay 
514a2ae208SSatish Balay #undef __FUNCT__
524a2ae208SSatish Balay #define __FUNCT__ "MatGetRowIJ_SeqSBAIJ"
5349b5e25fSSatish Balay static int MatGetRowIJ_SeqSBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja,PetscTruth *done)
5449b5e25fSSatish Balay {
5549b5e25fSSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
5649b5e25fSSatish Balay 
5749b5e25fSSatish Balay   PetscFunctionBegin;
58a1373b80SHong Zhang 
59045c9aa0SHong Zhang   if (ia) {
6029bbc08cSBarry Smith     SETERRQ(1,"Function not yet written for SBAIJ format, only supports natural ordering");
6149b5e25fSSatish Balay   }
62a1373b80SHong Zhang 
63a1373b80SHong Zhang #ifdef NEW
64a1373b80SHong Zhang   int         ierr;
65a1373b80SHong Zhang 
66a1373b80SHong Zhang   if (!ia) PetscFunctionReturn(0);
67a1373b80SHong Zhang 
68a1373b80SHong Zhang /*
69a1373b80SHong Zhang   if (symmetric) {
70a1373b80SHong Zhang     ierr = MatToSymmetricIJ_SeqAIJ(n,a->i,a->j,0,oshift,ia,ja);CHKERRQ(ierr);
71a1373b80SHong Zhang   } else if (oshift == 1) {
72a1373b80SHong Zhang     int nz = a->i[n];
73a1373b80SHong Zhang     for (i=0; i<nz; i++) a->j[i]++;
74a1373b80SHong Zhang     for (i=0; i<n+1; i++) a->i[i]++;
75a1373b80SHong Zhang     *ia = a->i; *ja = a->j;
76a1373b80SHong Zhang   } else {
77a1373b80SHong Zhang */
78a1373b80SHong Zhang     *ia = a->i; *ja = a->j;
79a1373b80SHong Zhang     /* } */
80a1373b80SHong Zhang #endif
81045c9aa0SHong Zhang   *nn = a->mbs;
8249b5e25fSSatish Balay   PetscFunctionReturn(0);
8349b5e25fSSatish Balay }
8449b5e25fSSatish Balay 
854a2ae208SSatish Balay #undef __FUNCT__
864a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRowIJ_SeqSBAIJ"
87045c9aa0SHong Zhang static int MatRestoreRowIJ_SeqSBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja,PetscTruth *done)
8849b5e25fSSatish Balay {
8949b5e25fSSatish Balay   PetscFunctionBegin;
9049b5e25fSSatish Balay   if (!ia) PetscFunctionReturn(0);
9129bbc08cSBarry Smith   SETERRQ(1,"Function not yet written for SBAIJ format, only supports natural ordering");
9216cdd363SHong Zhang   /* PetscFunctionReturn(0); */
9349b5e25fSSatish Balay }
9449b5e25fSSatish Balay 
954a2ae208SSatish Balay #undef __FUNCT__
964a2ae208SSatish Balay #define __FUNCT__ "MatGetBlockSize_SeqSBAIJ"
9749b5e25fSSatish Balay int MatGetBlockSize_SeqSBAIJ(Mat mat,int *bs)
9849b5e25fSSatish Balay {
9949b5e25fSSatish Balay   Mat_SeqSBAIJ *sbaij = (Mat_SeqSBAIJ*)mat->data;
10049b5e25fSSatish Balay 
10149b5e25fSSatish Balay   PetscFunctionBegin;
10249b5e25fSSatish Balay   *bs = sbaij->bs;
10349b5e25fSSatish Balay   PetscFunctionReturn(0);
10449b5e25fSSatish Balay }
10549b5e25fSSatish Balay 
1064a2ae208SSatish Balay #undef __FUNCT__
1074a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqSBAIJ"
10849b5e25fSSatish Balay int MatDestroy_SeqSBAIJ(Mat A)
10949b5e25fSSatish Balay {
11049b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
11149b5e25fSSatish Balay   int         ierr;
11249b5e25fSSatish Balay 
11349b5e25fSSatish Balay   PetscFunctionBegin;
11449b5e25fSSatish Balay #if defined(PETSC_USE_LOG)
115b0a32e0cSBarry Smith   PetscLogObjectState((PetscObject)A,"Rows=%d, s_NZ=%d",A->m,a->s_nz);
11649b5e25fSSatish Balay #endif
11749b5e25fSSatish Balay   ierr = PetscFree(a->a);CHKERRQ(ierr);
11849b5e25fSSatish Balay   if (!a->singlemalloc) {
11949b5e25fSSatish Balay     ierr = PetscFree(a->i);CHKERRQ(ierr);
12063c8bf9fSHong Zhang     ierr = PetscFree(a->j);CHKERRQ(ierr);
12149b5e25fSSatish Balay   }
12249b5e25fSSatish Balay   if (a->row) {
12349b5e25fSSatish Balay     ierr = ISDestroy(a->row);CHKERRQ(ierr);
12449b5e25fSSatish Balay   }
12549b5e25fSSatish Balay   if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);}
12649b5e25fSSatish Balay   if (a->ilen) {ierr = PetscFree(a->ilen);CHKERRQ(ierr);}
12749b5e25fSSatish Balay   if (a->imax) {ierr = PetscFree(a->imax);CHKERRQ(ierr);}
12849b5e25fSSatish Balay   if (a->solve_work) {ierr = PetscFree(a->solve_work);CHKERRQ(ierr);}
12949b5e25fSSatish Balay   if (a->mult_work) {ierr = PetscFree(a->mult_work);CHKERRQ(ierr);}
13049b5e25fSSatish Balay   if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);}
13149b5e25fSSatish Balay   if (a->saved_values) {ierr = PetscFree(a->saved_values);CHKERRQ(ierr);}
1321a3463dfSHong Zhang 
1331a3463dfSHong Zhang   if (a->inew){
1341a3463dfSHong Zhang     ierr = PetscFree(a->inew);CHKERRQ(ierr);
1351a3463dfSHong Zhang     a->inew = 0;
1361a3463dfSHong Zhang   }
13749b5e25fSSatish Balay   ierr = PetscFree(a);CHKERRQ(ierr);
13849b5e25fSSatish Balay   PetscFunctionReturn(0);
13949b5e25fSSatish Balay }
14049b5e25fSSatish Balay 
1414a2ae208SSatish Balay #undef __FUNCT__
1424a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqSBAIJ"
14349b5e25fSSatish Balay int MatSetOption_SeqSBAIJ(Mat A,MatOption op)
14449b5e25fSSatish Balay {
145045c9aa0SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
14649b5e25fSSatish Balay 
14749b5e25fSSatish Balay   PetscFunctionBegin;
1484d9d31abSKris Buschelman   switch (op) {
1494d9d31abSKris Buschelman   case MAT_ROW_ORIENTED:
1504d9d31abSKris Buschelman     a->roworiented    = PETSC_TRUE;
1514d9d31abSKris Buschelman     break;
1524d9d31abSKris Buschelman   case MAT_COLUMN_ORIENTED:
1534d9d31abSKris Buschelman     a->roworiented    = PETSC_FALSE;
1544d9d31abSKris Buschelman     break;
1554d9d31abSKris Buschelman   case MAT_COLUMNS_SORTED:
1564d9d31abSKris Buschelman     a->sorted         = PETSC_TRUE;
1574d9d31abSKris Buschelman     break;
1584d9d31abSKris Buschelman   case MAT_COLUMNS_UNSORTED:
1594d9d31abSKris Buschelman     a->sorted         = PETSC_FALSE;
1604d9d31abSKris Buschelman     break;
1614d9d31abSKris Buschelman   case MAT_KEEP_ZEROED_ROWS:
1624d9d31abSKris Buschelman     a->keepzeroedrows = PETSC_TRUE;
1634d9d31abSKris Buschelman     break;
1644d9d31abSKris Buschelman   case MAT_NO_NEW_NONZERO_LOCATIONS:
1654d9d31abSKris Buschelman     a->nonew          = 1;
1664d9d31abSKris Buschelman     break;
1674d9d31abSKris Buschelman   case MAT_NEW_NONZERO_LOCATION_ERR:
1684d9d31abSKris Buschelman     a->nonew          = -1;
1694d9d31abSKris Buschelman     break;
1704d9d31abSKris Buschelman   case MAT_NEW_NONZERO_ALLOCATION_ERR:
1714d9d31abSKris Buschelman     a->nonew          = -2;
1724d9d31abSKris Buschelman     break;
1734d9d31abSKris Buschelman   case MAT_YES_NEW_NONZERO_LOCATIONS:
1744d9d31abSKris Buschelman     a->nonew          = 0;
1754d9d31abSKris Buschelman     break;
1764d9d31abSKris Buschelman   case MAT_ROWS_SORTED:
1774d9d31abSKris Buschelman   case MAT_ROWS_UNSORTED:
1784d9d31abSKris Buschelman   case MAT_YES_NEW_DIAGONALS:
1794d9d31abSKris Buschelman   case MAT_IGNORE_OFF_PROC_ENTRIES:
1804d9d31abSKris Buschelman   case MAT_USE_HASH_TABLE:
181d03495bdSKris Buschelman   case MAT_USE_SINGLE_PRECISION_SOLVES:
182b0a32e0cSBarry Smith     PetscLogInfo(A,"MatSetOption_SeqSBAIJ:Option ignored\n");
1834d9d31abSKris Buschelman     break;
1844d9d31abSKris Buschelman   case MAT_NO_NEW_DIAGONALS:
18529bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS");
1864d9d31abSKris Buschelman   default:
18729bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"unknown option");
18849b5e25fSSatish Balay   }
18949b5e25fSSatish Balay   PetscFunctionReturn(0);
19049b5e25fSSatish Balay }
19149b5e25fSSatish Balay 
1924a2ae208SSatish Balay #undef __FUNCT__
1934a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqSBAIJ"
19487828ca2SBarry Smith int MatGetRow_SeqSBAIJ(Mat A,int row,int *ncols,int **cols,PetscScalar **v)
19549b5e25fSSatish Balay {
19649b5e25fSSatish Balay   Mat_SeqSBAIJ  *a = (Mat_SeqSBAIJ*)A->data;
19782502324SSatish Balay   int          itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*cols_i,bs2,ierr;
19849b5e25fSSatish Balay   MatScalar    *aa,*aa_i;
19987828ca2SBarry Smith   PetscScalar  *v_i;
20049b5e25fSSatish Balay 
20149b5e25fSSatish Balay   PetscFunctionBegin;
20249b5e25fSSatish Balay   bs  = a->bs;
20349b5e25fSSatish Balay   ai  = a->i;
20449b5e25fSSatish Balay   aj  = a->j;
20549b5e25fSSatish Balay   aa  = a->a;
20649b5e25fSSatish Balay   bs2 = a->bs2;
20749b5e25fSSatish Balay 
208c464158bSHong Zhang   if (row < 0 || row >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row out of range");
20949b5e25fSSatish Balay 
21049b5e25fSSatish Balay   bn  = row/bs;   /* Block number */
21149b5e25fSSatish Balay   bp  = row % bs; /* Block position */
21249b5e25fSSatish Balay   M   = ai[bn+1] - ai[bn];
21349b5e25fSSatish Balay   *ncols = bs*M;
21449b5e25fSSatish Balay 
21549b5e25fSSatish Balay   if (v) {
21649b5e25fSSatish Balay     *v = 0;
21749b5e25fSSatish Balay     if (*ncols) {
21887828ca2SBarry Smith       ierr = PetscMalloc((*ncols+row)*sizeof(PetscScalar),v);CHKERRQ(ierr);
21949b5e25fSSatish Balay       for (i=0; i<M; i++) { /* for each block in the block row */
22049b5e25fSSatish Balay         v_i  = *v + i*bs;
22149b5e25fSSatish Balay         aa_i = aa + bs2*(ai[bn] + i);
22249b5e25fSSatish Balay         for (j=bp,k=0; j<bs2; j+=bs,k++) {v_i[k] = aa_i[j];}
22349b5e25fSSatish Balay       }
22449b5e25fSSatish Balay     }
22549b5e25fSSatish Balay   }
22649b5e25fSSatish Balay 
22749b5e25fSSatish Balay   if (cols) {
22849b5e25fSSatish Balay     *cols = 0;
22949b5e25fSSatish Balay     if (*ncols) {
23082502324SSatish Balay       ierr = PetscMalloc((*ncols+row)*sizeof(int),cols);CHKERRQ(ierr);
23149b5e25fSSatish Balay       for (i=0; i<M; i++) { /* for each block in the block row */
23249b5e25fSSatish Balay         cols_i = *cols + i*bs;
23349b5e25fSSatish Balay         itmp  = bs*aj[ai[bn] + i];
23449b5e25fSSatish Balay         for (j=0; j<bs; j++) {cols_i[j] = itmp++;}
23549b5e25fSSatish Balay       }
23649b5e25fSSatish Balay     }
23749b5e25fSSatish Balay   }
23849b5e25fSSatish Balay 
23949b5e25fSSatish Balay   /*search column A(0:row-1,row) (=A(row,0:row-1)). Could be expensive! */
2405ddb2528SHong Zhang   /* this segment is currently removed, so only entries in the upper triangle are obtained */
2415ddb2528SHong Zhang #ifdef column_search
24249b5e25fSSatish Balay   v_i    = *v    + M*bs;
24349b5e25fSSatish Balay   cols_i = *cols + M*bs;
24449b5e25fSSatish Balay   for (i=0; i<bn; i++){ /* for each block row */
24549b5e25fSSatish Balay     M = ai[i+1] - ai[i];
24649b5e25fSSatish Balay     for (j=0; j<M; j++){
24749b5e25fSSatish Balay       itmp = aj[ai[i] + j];    /* block column value */
24849b5e25fSSatish Balay       if (itmp == bn){
24949b5e25fSSatish Balay         aa_i   = aa    + bs2*(ai[i] + j) + bs*bp;
25049b5e25fSSatish Balay         for (k=0; k<bs; k++) {
25149b5e25fSSatish Balay           *cols_i++ = i*bs+k;
25249b5e25fSSatish Balay           *v_i++    = aa_i[k];
25349b5e25fSSatish Balay         }
25449b5e25fSSatish Balay         *ncols += bs;
25549b5e25fSSatish Balay         break;
25649b5e25fSSatish Balay       }
25749b5e25fSSatish Balay     }
25849b5e25fSSatish Balay   }
2595ddb2528SHong Zhang #endif
26049b5e25fSSatish Balay 
26149b5e25fSSatish Balay   PetscFunctionReturn(0);
26249b5e25fSSatish Balay }
26349b5e25fSSatish Balay 
2644a2ae208SSatish Balay #undef __FUNCT__
2654a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqSBAIJ"
26687828ca2SBarry Smith int MatRestoreRow_SeqSBAIJ(Mat A,int row,int *nz,int **idx,PetscScalar **v)
26749b5e25fSSatish Balay {
26849b5e25fSSatish Balay   int ierr;
26949b5e25fSSatish Balay 
27049b5e25fSSatish Balay   PetscFunctionBegin;
27149b5e25fSSatish Balay   if (idx) {if (*idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}}
27249b5e25fSSatish Balay   if (v)   {if (*v)   {ierr = PetscFree(*v);CHKERRQ(ierr);}}
27349b5e25fSSatish Balay   PetscFunctionReturn(0);
27449b5e25fSSatish Balay }
27549b5e25fSSatish Balay 
2764a2ae208SSatish Balay #undef __FUNCT__
2774a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqSBAIJ"
27849b5e25fSSatish Balay int MatTranspose_SeqSBAIJ(Mat A,Mat *B)
27949b5e25fSSatish Balay {
2808115998fSBarry Smith   int ierr;
28149b5e25fSSatish Balay   PetscFunctionBegin;
282999d9058SBarry Smith   ierr = MatDuplicate(A,MAT_COPY_VALUES,B);CHKERRQ(ierr);
2838115998fSBarry Smith   PetscFunctionReturn(0);
28449b5e25fSSatish Balay }
28549b5e25fSSatish Balay 
2864a2ae208SSatish Balay #undef __FUNCT__
2874a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Binary"
288b0a32e0cSBarry Smith static int MatView_SeqSBAIJ_Binary(Mat A,PetscViewer viewer)
28949b5e25fSSatish Balay {
29049b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
29149b5e25fSSatish Balay   int          i,fd,*col_lens,ierr,bs = a->bs,count,*jj,j,k,l,bs2=a->bs2;
29287828ca2SBarry Smith   PetscScalar  *aa;
29349b5e25fSSatish Balay   FILE         *file;
29449b5e25fSSatish Balay 
29549b5e25fSSatish Balay   PetscFunctionBegin;
296b0a32e0cSBarry Smith   ierr        = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
29782502324SSatish Balay   ierr        = PetscMalloc((4+A->m)*sizeof(int),&col_lens);CHKERRQ(ierr);
298552e946dSBarry Smith   col_lens[0] = MAT_FILE_COOKIE;
29949b5e25fSSatish Balay 
300c464158bSHong Zhang   col_lens[1] = A->m;
301c464158bSHong Zhang   col_lens[2] = A->m;
30249b5e25fSSatish Balay   col_lens[3] = a->s_nz*bs2;
30349b5e25fSSatish Balay 
30449b5e25fSSatish Balay   /* store lengths of each row and write (including header) to file */
30549b5e25fSSatish Balay   count = 0;
30649b5e25fSSatish Balay   for (i=0; i<a->mbs; i++) {
30749b5e25fSSatish Balay     for (j=0; j<bs; j++) {
30849b5e25fSSatish Balay       col_lens[4+count++] = bs*(a->i[i+1] - a->i[i]);
30949b5e25fSSatish Balay     }
31049b5e25fSSatish Balay   }
311c464158bSHong Zhang   ierr = PetscBinaryWrite(fd,col_lens,4+A->m,PETSC_INT,1);CHKERRQ(ierr);
31249b5e25fSSatish Balay   ierr = PetscFree(col_lens);CHKERRQ(ierr);
31349b5e25fSSatish Balay 
31449b5e25fSSatish Balay   /* store column indices (zero start index) */
31582502324SSatish Balay   ierr  = PetscMalloc((a->s_nz+1)*bs2*sizeof(int),&jj);CHKERRQ(ierr);
31649b5e25fSSatish Balay   count = 0;
31749b5e25fSSatish Balay   for (i=0; i<a->mbs; i++) {
31849b5e25fSSatish Balay     for (j=0; j<bs; j++) {
31949b5e25fSSatish Balay       for (k=a->i[i]; k<a->i[i+1]; k++) {
32049b5e25fSSatish Balay         for (l=0; l<bs; l++) {
32149b5e25fSSatish Balay           jj[count++] = bs*a->j[k] + l;
32249b5e25fSSatish Balay         }
32349b5e25fSSatish Balay       }
32449b5e25fSSatish Balay     }
32549b5e25fSSatish Balay   }
32649b5e25fSSatish Balay   ierr = PetscBinaryWrite(fd,jj,bs2*a->s_nz,PETSC_INT,0);CHKERRQ(ierr);
32749b5e25fSSatish Balay   ierr = PetscFree(jj);CHKERRQ(ierr);
32849b5e25fSSatish Balay 
32949b5e25fSSatish Balay   /* store nonzero values */
33087828ca2SBarry Smith   ierr  = PetscMalloc((a->s_nz+1)*bs2*sizeof(PetscScalar),&aa);CHKERRQ(ierr);
33149b5e25fSSatish Balay   count = 0;
33249b5e25fSSatish Balay   for (i=0; i<a->mbs; i++) {
33349b5e25fSSatish Balay     for (j=0; j<bs; j++) {
33449b5e25fSSatish Balay       for (k=a->i[i]; k<a->i[i+1]; k++) {
33549b5e25fSSatish Balay         for (l=0; l<bs; l++) {
33649b5e25fSSatish Balay           aa[count++] = a->a[bs2*k + l*bs + j];
33749b5e25fSSatish Balay         }
33849b5e25fSSatish Balay       }
33949b5e25fSSatish Balay     }
34049b5e25fSSatish Balay   }
34149b5e25fSSatish Balay   ierr = PetscBinaryWrite(fd,aa,bs2*a->s_nz,PETSC_SCALAR,0);CHKERRQ(ierr);
34249b5e25fSSatish Balay   ierr = PetscFree(aa);CHKERRQ(ierr);
34349b5e25fSSatish Balay 
344b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetInfoPointer(viewer,&file);CHKERRQ(ierr);
34549b5e25fSSatish Balay   if (file) {
34649b5e25fSSatish Balay     fprintf(file,"-matload_block_size %d\n",a->bs);
34749b5e25fSSatish Balay   }
34849b5e25fSSatish Balay   PetscFunctionReturn(0);
34949b5e25fSSatish Balay }
35049b5e25fSSatish Balay 
3514a2ae208SSatish Balay #undef __FUNCT__
3524a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_ASCII"
353b0a32e0cSBarry Smith static int MatView_SeqSBAIJ_ASCII(Mat A,PetscViewer viewer)
35449b5e25fSSatish Balay {
35549b5e25fSSatish Balay   Mat_SeqSBAIJ      *a = (Mat_SeqSBAIJ*)A->data;
356fb9695e5SSatish Balay   int               ierr,i,j,bs = a->bs,k,l,bs2=a->bs2;
357fb9695e5SSatish Balay   char              *name;
358f3ef73ceSBarry Smith   PetscViewerFormat format;
35949b5e25fSSatish Balay 
36049b5e25fSSatish Balay   PetscFunctionBegin;
36180fe4e49SBarry Smith   ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr);
362b0a32e0cSBarry Smith   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
363456192e2SBarry Smith   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
364b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  block size is %d\n",bs);CHKERRQ(ierr);
365fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_MATLAB) {
36629bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Matlab format not supported");
367fb9695e5SSatish Balay   } else if (format == PETSC_VIEWER_ASCII_COMMON) {
368b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
36949b5e25fSSatish Balay     for (i=0; i<a->mbs; i++) {
37049b5e25fSSatish Balay       for (j=0; j<bs; j++) {
371b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr);
37249b5e25fSSatish Balay         for (k=a->i[i]; k<a->i[i+1]; k++) {
37349b5e25fSSatish Balay           for (l=0; l<bs; l++) {
37449b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX)
37549b5e25fSSatish Balay             if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
37662b941d6SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%d, %g + %g i) ",bs*a->j[k]+l,
37749b5e25fSSatish Balay                       PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
37849b5e25fSSatish Balay             } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
37962b941d6SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%d, %g - %g i) ",bs*a->j[k]+l,
38049b5e25fSSatish Balay                       PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
38149b5e25fSSatish Balay             } else if (PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
38262b941d6SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
38349b5e25fSSatish Balay             }
38449b5e25fSSatish Balay #else
38549b5e25fSSatish Balay             if (a->a[bs2*k + l*bs + j] != 0.0) {
38662b941d6SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr);
38749b5e25fSSatish Balay             }
38849b5e25fSSatish Balay #endif
38949b5e25fSSatish Balay           }
39049b5e25fSSatish Balay         }
391b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
39249b5e25fSSatish Balay       }
39349b5e25fSSatish Balay     }
394b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
39549b5e25fSSatish Balay   } else {
396b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
39749b5e25fSSatish Balay     for (i=0; i<a->mbs; i++) {
39849b5e25fSSatish Balay       for (j=0; j<bs; j++) {
399b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr);
40049b5e25fSSatish Balay         for (k=a->i[i]; k<a->i[i+1]; k++) {
40149b5e25fSSatish Balay           for (l=0; l<bs; l++) {
40249b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX)
40349b5e25fSSatish Balay             if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0) {
40462b941d6SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%d, %g + %g i) ",bs*a->j[k]+l,
40549b5e25fSSatish Balay                 PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
40649b5e25fSSatish Balay             } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0) {
40762b941d6SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%d, %g - %g i) ",bs*a->j[k]+l,
40849b5e25fSSatish Balay                 PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
40949b5e25fSSatish Balay             } else {
41062b941d6SBarry Smith               ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
41149b5e25fSSatish Balay             }
41249b5e25fSSatish Balay #else
413b0a32e0cSBarry Smith             ierr = PetscViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr);
41449b5e25fSSatish Balay #endif
41549b5e25fSSatish Balay           }
41649b5e25fSSatish Balay         }
417b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
41849b5e25fSSatish Balay       }
41949b5e25fSSatish Balay     }
420b0a32e0cSBarry Smith     ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
42149b5e25fSSatish Balay   }
422b0a32e0cSBarry Smith   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
42349b5e25fSSatish Balay   PetscFunctionReturn(0);
42449b5e25fSSatish Balay }
42549b5e25fSSatish Balay 
4264a2ae208SSatish Balay #undef __FUNCT__
4274a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Draw_Zoom"
428b0a32e0cSBarry Smith static int MatView_SeqSBAIJ_Draw_Zoom(PetscDraw draw,void *Aa)
42949b5e25fSSatish Balay {
43049b5e25fSSatish Balay   Mat          A = (Mat) Aa;
43149b5e25fSSatish Balay   Mat_SeqSBAIJ *a=(Mat_SeqSBAIJ*)A->data;
43249b5e25fSSatish Balay   int          row,ierr,i,j,k,l,mbs=a->mbs,color,bs=a->bs,bs2=a->bs2,rank;
43349b5e25fSSatish Balay   PetscReal    xl,yl,xr,yr,x_l,x_r,y_l,y_r;
43449b5e25fSSatish Balay   MatScalar    *aa;
43549b5e25fSSatish Balay   MPI_Comm     comm;
436b0a32e0cSBarry Smith   PetscViewer  viewer;
43749b5e25fSSatish Balay 
43849b5e25fSSatish Balay   PetscFunctionBegin;
43949b5e25fSSatish Balay   /*
44049b5e25fSSatish Balay       This is nasty. If this is called from an originally parallel matrix
44149b5e25fSSatish Balay    then all processes call this,but only the first has the matrix so the
44249b5e25fSSatish Balay    rest should return immediately.
44349b5e25fSSatish Balay   */
44449b5e25fSSatish Balay   ierr = PetscObjectGetComm((PetscObject)draw,&comm);CHKERRQ(ierr);
44549b5e25fSSatish Balay   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
44649b5e25fSSatish Balay   if (rank) PetscFunctionReturn(0);
44749b5e25fSSatish Balay 
44849b5e25fSSatish Balay   ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr);
44949b5e25fSSatish Balay 
450b0a32e0cSBarry Smith   ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr);
451b0a32e0cSBarry Smith   PetscDrawString(draw, .3*(xl+xr), .3*(yl+yr), PETSC_DRAW_BLACK, "symmetric");
45249b5e25fSSatish Balay 
45349b5e25fSSatish Balay   /* loop over matrix elements drawing boxes */
454b0a32e0cSBarry Smith   color = PETSC_DRAW_BLUE;
45549b5e25fSSatish Balay   for (i=0,row=0; i<mbs; i++,row+=bs) {
45649b5e25fSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
457c464158bSHong Zhang       y_l = A->m - row - 1.0; y_r = y_l + 1.0;
45849b5e25fSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
45949b5e25fSSatish Balay       aa = a->a + j*bs2;
46049b5e25fSSatish Balay       for (k=0; k<bs; k++) {
46149b5e25fSSatish Balay         for (l=0; l<bs; l++) {
46249b5e25fSSatish Balay           if (PetscRealPart(*aa++) >=  0.) continue;
463b0a32e0cSBarry Smith           ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
46449b5e25fSSatish Balay         }
46549b5e25fSSatish Balay       }
46649b5e25fSSatish Balay     }
46749b5e25fSSatish Balay   }
468b0a32e0cSBarry Smith   color = PETSC_DRAW_CYAN;
46949b5e25fSSatish Balay   for (i=0,row=0; i<mbs; i++,row+=bs) {
47049b5e25fSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
471c464158bSHong Zhang       y_l = A->m - row - 1.0; y_r = y_l + 1.0;
47249b5e25fSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
47349b5e25fSSatish Balay       aa = a->a + j*bs2;
47449b5e25fSSatish Balay       for (k=0; k<bs; k++) {
47549b5e25fSSatish Balay         for (l=0; l<bs; l++) {
47649b5e25fSSatish Balay           if (PetscRealPart(*aa++) != 0.) continue;
477b0a32e0cSBarry Smith           ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
47849b5e25fSSatish Balay         }
47949b5e25fSSatish Balay       }
48049b5e25fSSatish Balay     }
48149b5e25fSSatish Balay   }
48249b5e25fSSatish Balay 
483b0a32e0cSBarry Smith   color = PETSC_DRAW_RED;
48449b5e25fSSatish Balay   for (i=0,row=0; i<mbs; i++,row+=bs) {
48549b5e25fSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
486c464158bSHong Zhang       y_l = A->m - row - 1.0; y_r = y_l + 1.0;
48749b5e25fSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
48849b5e25fSSatish Balay       aa = a->a + j*bs2;
48949b5e25fSSatish Balay       for (k=0; k<bs; k++) {
49049b5e25fSSatish Balay         for (l=0; l<bs; l++) {
49149b5e25fSSatish Balay           if (PetscRealPart(*aa++) <= 0.) continue;
492b0a32e0cSBarry Smith           ierr = PetscDrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
49349b5e25fSSatish Balay         }
49449b5e25fSSatish Balay       }
49549b5e25fSSatish Balay     }
49649b5e25fSSatish Balay   }
49749b5e25fSSatish Balay   PetscFunctionReturn(0);
49849b5e25fSSatish Balay }
49949b5e25fSSatish Balay 
5004a2ae208SSatish Balay #undef __FUNCT__
5014a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ_Draw"
502b0a32e0cSBarry Smith static int MatView_SeqSBAIJ_Draw(Mat A,PetscViewer viewer)
50349b5e25fSSatish Balay {
50449b5e25fSSatish Balay   int          ierr;
50549b5e25fSSatish Balay   PetscReal    xl,yl,xr,yr,w,h;
506b0a32e0cSBarry Smith   PetscDraw    draw;
50749b5e25fSSatish Balay   PetscTruth   isnull;
50849b5e25fSSatish Balay 
50949b5e25fSSatish Balay   PetscFunctionBegin;
51049b5e25fSSatish Balay 
511b0a32e0cSBarry Smith   ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
512b0a32e0cSBarry Smith   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
51349b5e25fSSatish Balay 
51449b5e25fSSatish Balay   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
515c464158bSHong Zhang   xr  = A->m; yr = A->m; h = yr/10.0; w = xr/10.0;
51649b5e25fSSatish Balay   xr += w;    yr += h;  xl = -w;     yl = -h;
517b0a32e0cSBarry Smith   ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
518b0a32e0cSBarry Smith   ierr = PetscDrawZoom(draw,MatView_SeqSBAIJ_Draw_Zoom,A);CHKERRQ(ierr);
51949b5e25fSSatish Balay   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr);
52049b5e25fSSatish Balay   PetscFunctionReturn(0);
52149b5e25fSSatish Balay }
52249b5e25fSSatish Balay 
5234a2ae208SSatish Balay #undef __FUNCT__
5244a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqSBAIJ"
525b0a32e0cSBarry Smith int MatView_SeqSBAIJ(Mat A,PetscViewer viewer)
52649b5e25fSSatish Balay {
52749b5e25fSSatish Balay   int        ierr;
52849b5e25fSSatish Balay   PetscTruth issocket,isascii,isbinary,isdraw;
52949b5e25fSSatish Balay 
53049b5e25fSSatish Balay   PetscFunctionBegin;
531b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);CHKERRQ(ierr);
532b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
533fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr);
534fb9695e5SSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr);
53549b5e25fSSatish Balay   if (issocket) {
53629bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Socket viewer not supported");
53749b5e25fSSatish Balay   } else if (isascii){
53849b5e25fSSatish Balay     ierr = MatView_SeqSBAIJ_ASCII(A,viewer);CHKERRQ(ierr);
53949b5e25fSSatish Balay   } else if (isbinary) {
54049b5e25fSSatish Balay     ierr = MatView_SeqSBAIJ_Binary(A,viewer);CHKERRQ(ierr);
54149b5e25fSSatish Balay   } else if (isdraw) {
54249b5e25fSSatish Balay     ierr = MatView_SeqSBAIJ_Draw(A,viewer);CHKERRQ(ierr);
54349b5e25fSSatish Balay   } else {
54429bbc08cSBarry Smith     SETERRQ1(1,"Viewer type %s not supported by SeqBAIJ matrices",((PetscObject)viewer)->type_name);
54549b5e25fSSatish Balay   }
54649b5e25fSSatish Balay   PetscFunctionReturn(0);
54749b5e25fSSatish Balay }
54849b5e25fSSatish Balay 
54949b5e25fSSatish Balay 
5504a2ae208SSatish Balay #undef __FUNCT__
5514a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqSBAIJ"
55287828ca2SBarry Smith int MatGetValues_SeqSBAIJ(Mat A,int m,int *im,int n,int *in,PetscScalar *v)
55349b5e25fSSatish Balay {
554045c9aa0SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
55549b5e25fSSatish Balay   int        *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j;
55649b5e25fSSatish Balay   int        *ai = a->i,*ailen = a->ilen;
55749b5e25fSSatish Balay   int        brow,bcol,ridx,cidx,bs=a->bs,bs2=a->bs2;
55849b5e25fSSatish Balay   MatScalar  *ap,*aa = a->a,zero = 0.0;
55949b5e25fSSatish Balay 
56049b5e25fSSatish Balay   PetscFunctionBegin;
56149b5e25fSSatish Balay   for (k=0; k<m; k++) { /* loop over rows */
56249b5e25fSSatish Balay     row  = im[k]; brow = row/bs;
56329bbc08cSBarry Smith     if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative row");
564c464158bSHong Zhang     if (row >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
56549b5e25fSSatish Balay     rp   = aj + ai[brow] ; ap = aa + bs2*ai[brow] ;
56649b5e25fSSatish Balay     nrow = ailen[brow];
56749b5e25fSSatish Balay     for (l=0; l<n; l++) { /* loop over columns */
56829bbc08cSBarry Smith       if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative column");
569c464158bSHong Zhang       if (in[l] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
57049b5e25fSSatish Balay       col  = in[l] ;
57149b5e25fSSatish Balay       bcol = col/bs;
57249b5e25fSSatish Balay       cidx = col%bs;
57349b5e25fSSatish Balay       ridx = row%bs;
57449b5e25fSSatish Balay       high = nrow;
57549b5e25fSSatish Balay       low  = 0; /* assume unsorted */
57649b5e25fSSatish Balay       while (high-low > 5) {
57749b5e25fSSatish Balay         t = (low+high)/2;
57849b5e25fSSatish Balay         if (rp[t] > bcol) high = t;
57949b5e25fSSatish Balay         else             low  = t;
58049b5e25fSSatish Balay       }
58149b5e25fSSatish Balay       for (i=low; i<high; i++) {
58249b5e25fSSatish Balay         if (rp[i] > bcol) break;
58349b5e25fSSatish Balay         if (rp[i] == bcol) {
58449b5e25fSSatish Balay           *v++ = ap[bs2*i+bs*cidx+ridx];
58549b5e25fSSatish Balay           goto finished;
58649b5e25fSSatish Balay         }
58749b5e25fSSatish Balay       }
58849b5e25fSSatish Balay       *v++ = zero;
58949b5e25fSSatish Balay       finished:;
59049b5e25fSSatish Balay     }
59149b5e25fSSatish Balay   }
59249b5e25fSSatish Balay   PetscFunctionReturn(0);
59349b5e25fSSatish Balay }
59449b5e25fSSatish Balay 
59549b5e25fSSatish Balay 
5964a2ae208SSatish Balay #undef __FUNCT__
5974a2ae208SSatish Balay #define __FUNCT__ "MatSetValuesBlocked_SeqSBAIJ"
59887828ca2SBarry Smith int MatSetValuesBlocked_SeqSBAIJ(Mat A,int m,int *im,int n,int *in,PetscScalar *v,InsertMode is)
59949b5e25fSSatish Balay {
600*0880e062SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
601*0880e062SHong Zhang   int         *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,sorted=a->sorted;
602*0880e062SHong Zhang   int         *imax=a->imax,*ai=a->i,*ailen=a->ilen;
603*0880e062SHong Zhang   int         *aj=a->j,nonew=a->nonew,bs2=a->bs2,bs=a->bs,stepval,ierr;
604*0880e062SHong Zhang   PetscTruth  roworiented=a->roworiented;
605*0880e062SHong Zhang   MatScalar   *value = v,*ap,*aa = a->a,*bap;
606*0880e062SHong Zhang 
60749b5e25fSSatish Balay   PetscFunctionBegin;
608*0880e062SHong Zhang   if (roworiented) {
609*0880e062SHong Zhang     stepval = (n-1)*bs;
610*0880e062SHong Zhang   } else {
611*0880e062SHong Zhang     stepval = (m-1)*bs;
612*0880e062SHong Zhang   }
613*0880e062SHong Zhang   for (k=0; k<m; k++) { /* loop over added rows */
614*0880e062SHong Zhang     row  = im[k];
615*0880e062SHong Zhang     if (row < 0) continue;
616*0880e062SHong Zhang #if defined(PETSC_USE_BOPT_g)
617*0880e062SHong Zhang     if (row >= a->mbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
618*0880e062SHong Zhang #endif
619*0880e062SHong Zhang     rp   = aj + ai[row];
620*0880e062SHong Zhang     ap   = aa + bs2*ai[row];
621*0880e062SHong Zhang     rmax = imax[row];
622*0880e062SHong Zhang     nrow = ailen[row];
623*0880e062SHong Zhang     low  = 0;
624*0880e062SHong Zhang     for (l=0; l<n; l++) { /* loop over added columns */
625*0880e062SHong Zhang       if (in[l] < 0) continue;
626*0880e062SHong Zhang #if defined(PETSC_USE_BOPT_g)
627*0880e062SHong Zhang       if (in[l] >= a->nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
628*0880e062SHong Zhang #endif
629*0880e062SHong Zhang       col = in[l];
630*0880e062SHong Zhang       if (col < row) continue; /* ignore lower triangular block */
631*0880e062SHong Zhang       if (roworiented) {
632*0880e062SHong Zhang         value = v + k*(stepval+bs)*bs + l*bs;
633*0880e062SHong Zhang       } else {
634*0880e062SHong Zhang         value = v + l*(stepval+bs)*bs + k*bs;
635*0880e062SHong Zhang       }
636*0880e062SHong Zhang       if (!sorted) low = 0; high = nrow;
637*0880e062SHong Zhang       while (high-low > 7) {
638*0880e062SHong Zhang         t = (low+high)/2;
639*0880e062SHong Zhang         if (rp[t] > col) high = t;
640*0880e062SHong Zhang         else             low  = t;
641*0880e062SHong Zhang       }
642*0880e062SHong Zhang       for (i=low; i<high; i++) {
643*0880e062SHong Zhang         if (rp[i] > col) break;
644*0880e062SHong Zhang         if (rp[i] == col) {
645*0880e062SHong Zhang           bap  = ap +  bs2*i;
646*0880e062SHong Zhang           if (roworiented) {
647*0880e062SHong Zhang             if (is == ADD_VALUES) {
648*0880e062SHong Zhang               for (ii=0; ii<bs; ii++,value+=stepval) {
649*0880e062SHong Zhang                 for (jj=ii; jj<bs2; jj+=bs) {
650*0880e062SHong Zhang                   bap[jj] += *value++;
651*0880e062SHong Zhang                 }
652*0880e062SHong Zhang               }
653*0880e062SHong Zhang             } else {
654*0880e062SHong Zhang               for (ii=0; ii<bs; ii++,value+=stepval) {
655*0880e062SHong Zhang                 for (jj=ii; jj<bs2; jj+=bs) {
656*0880e062SHong Zhang                   bap[jj] = *value++;
657*0880e062SHong Zhang                 }
658*0880e062SHong Zhang               }
659*0880e062SHong Zhang             }
660*0880e062SHong Zhang           } else {
661*0880e062SHong Zhang             if (is == ADD_VALUES) {
662*0880e062SHong Zhang               for (ii=0; ii<bs; ii++,value+=stepval) {
663*0880e062SHong Zhang                 for (jj=0; jj<bs; jj++) {
664*0880e062SHong Zhang                   *bap++ += *value++;
665*0880e062SHong Zhang                 }
666*0880e062SHong Zhang               }
667*0880e062SHong Zhang             } else {
668*0880e062SHong Zhang               for (ii=0; ii<bs; ii++,value+=stepval) {
669*0880e062SHong Zhang                 for (jj=0; jj<bs; jj++) {
670*0880e062SHong Zhang                   *bap++  = *value++;
671*0880e062SHong Zhang                 }
672*0880e062SHong Zhang               }
673*0880e062SHong Zhang             }
674*0880e062SHong Zhang           }
675*0880e062SHong Zhang           goto noinsert2;
676*0880e062SHong Zhang         }
677*0880e062SHong Zhang       }
678*0880e062SHong Zhang       if (nonew == 1) goto noinsert2;
679*0880e062SHong Zhang       else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix");
680*0880e062SHong Zhang       if (nrow >= rmax) {
681*0880e062SHong Zhang         /* there is no extra room in row, therefore enlarge */
682*0880e062SHong Zhang         int       new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j;
683*0880e062SHong Zhang         MatScalar *new_a;
684*0880e062SHong Zhang 
685*0880e062SHong Zhang         if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix");
686*0880e062SHong Zhang 
687*0880e062SHong Zhang         /* malloc new storage space */
688*0880e062SHong Zhang         len     = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int);
689*0880e062SHong Zhang 	ierr    = PetscMalloc(len,&new_a);CHKERRQ(ierr);
690*0880e062SHong Zhang         new_j   = (int*)(new_a + bs2*new_nz);
691*0880e062SHong Zhang         new_i   = new_j + new_nz;
692*0880e062SHong Zhang 
693*0880e062SHong Zhang         /* copy over old data into new slots */
694*0880e062SHong Zhang         for (ii=0; ii<row+1; ii++) {new_i[ii] = ai[ii];}
695*0880e062SHong Zhang         for (ii=row+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;}
696*0880e062SHong Zhang         ierr = PetscMemcpy(new_j,aj,(ai[row]+nrow)*sizeof(int));CHKERRQ(ierr);
697*0880e062SHong Zhang         len  = (new_nz - CHUNKSIZE - ai[row] - nrow);
698*0880e062SHong Zhang         ierr = PetscMemcpy(new_j+ai[row]+nrow+CHUNKSIZE,aj+ai[row]+nrow,len*sizeof(int));CHKERRQ(ierr);
699*0880e062SHong Zhang         ierr = PetscMemcpy(new_a,aa,(ai[row]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr);
700*0880e062SHong Zhang         ierr = PetscMemzero(new_a+bs2*(ai[row]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr);
701*0880e062SHong Zhang         ierr = PetscMemcpy(new_a+bs2*(ai[row]+nrow+CHUNKSIZE),aa+bs2*(ai[row]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr);
702*0880e062SHong Zhang         /* free up old matrix storage */
703*0880e062SHong Zhang         ierr = PetscFree(a->a);CHKERRQ(ierr);
704*0880e062SHong Zhang         if (!a->singlemalloc) {
705*0880e062SHong Zhang           ierr = PetscFree(a->i);CHKERRQ(ierr);
706*0880e062SHong Zhang           ierr = PetscFree(a->j);CHKERRQ(ierr);
707*0880e062SHong Zhang         }
708*0880e062SHong Zhang         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
709*0880e062SHong Zhang         a->singlemalloc = PETSC_TRUE;
710*0880e062SHong Zhang 
711*0880e062SHong Zhang         rp   = aj + ai[row]; ap = aa + bs2*ai[row];
712*0880e062SHong Zhang         rmax = imax[row] = imax[row] + CHUNKSIZE;
713*0880e062SHong Zhang         PetscLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar)));
714*0880e062SHong Zhang         a->s_maxnz += bs2*CHUNKSIZE;
715*0880e062SHong Zhang         a->reallocs++;
716*0880e062SHong Zhang         a->s_nz++;
717*0880e062SHong Zhang       }
718*0880e062SHong Zhang       N = nrow++ - 1;
719*0880e062SHong Zhang       /* shift up all the later entries in this row */
720*0880e062SHong Zhang       for (ii=N; ii>=i; ii--) {
721*0880e062SHong Zhang         rp[ii+1] = rp[ii];
722*0880e062SHong Zhang         ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr);
723*0880e062SHong Zhang       }
724*0880e062SHong Zhang       if (N >= i) {
725*0880e062SHong Zhang         ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr);
726*0880e062SHong Zhang       }
727*0880e062SHong Zhang       rp[i] = col;
728*0880e062SHong Zhang       bap   = ap +  bs2*i;
729*0880e062SHong Zhang       if (roworiented) {
730*0880e062SHong Zhang         for (ii=0; ii<bs; ii++,value+=stepval) {
731*0880e062SHong Zhang           for (jj=ii; jj<bs2; jj+=bs) {
732*0880e062SHong Zhang             bap[jj] = *value++;
733*0880e062SHong Zhang           }
734*0880e062SHong Zhang         }
735*0880e062SHong Zhang       } else {
736*0880e062SHong Zhang         for (ii=0; ii<bs; ii++,value+=stepval) {
737*0880e062SHong Zhang           for (jj=0; jj<bs; jj++) {
738*0880e062SHong Zhang             *bap++  = *value++;
739*0880e062SHong Zhang           }
740*0880e062SHong Zhang         }
741*0880e062SHong Zhang       }
742*0880e062SHong Zhang       noinsert2:;
743*0880e062SHong Zhang       low = i;
744*0880e062SHong Zhang     }
745*0880e062SHong Zhang     ailen[row] = nrow;
746*0880e062SHong Zhang   }
747*0880e062SHong Zhang   PetscFunctionReturn(0);
74849b5e25fSSatish Balay }
74949b5e25fSSatish Balay 
7504a2ae208SSatish Balay #undef __FUNCT__
7514a2ae208SSatish Balay #define __FUNCT__ "MatAssemblyEnd_SeqSBAIJ"
75249b5e25fSSatish Balay int MatAssemblyEnd_SeqSBAIJ(Mat A,MatAssemblyType mode)
75349b5e25fSSatish Balay {
75449b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
75549b5e25fSSatish Balay   int        fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax;
756c464158bSHong Zhang   int        m = A->m,*ip,N,*ailen = a->ilen;
75749b5e25fSSatish Balay   int        mbs = a->mbs,bs2 = a->bs2,rmax = 0,ierr;
75849b5e25fSSatish Balay   MatScalar  *aa = a->a,*ap;
759cf4441caSHong Zhang #if defined(PETSC_HAVE_SPOOLES)
760cf4441caSHong Zhang   PetscTruth   flag;
761cf4441caSHong Zhang #endif
76249b5e25fSSatish Balay 
76349b5e25fSSatish Balay   PetscFunctionBegin;
76449b5e25fSSatish Balay   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0);
76549b5e25fSSatish Balay 
76649b5e25fSSatish Balay   if (m) rmax = ailen[0];
76749b5e25fSSatish Balay   for (i=1; i<mbs; i++) {
76849b5e25fSSatish Balay     /* move each row back by the amount of empty slots (fshift) before it*/
76949b5e25fSSatish Balay     fshift += imax[i-1] - ailen[i-1];
77049b5e25fSSatish Balay     rmax   = PetscMax(rmax,ailen[i]);
77149b5e25fSSatish Balay     if (fshift) {
77249b5e25fSSatish Balay       ip = aj + ai[i]; ap = aa + bs2*ai[i];
77349b5e25fSSatish Balay       N = ailen[i];
77449b5e25fSSatish Balay       for (j=0; j<N; j++) {
77549b5e25fSSatish Balay         ip[j-fshift] = ip[j];
77649b5e25fSSatish Balay         ierr = PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
77749b5e25fSSatish Balay       }
77849b5e25fSSatish Balay     }
77949b5e25fSSatish Balay     ai[i] = ai[i-1] + ailen[i-1];
78049b5e25fSSatish Balay   }
78149b5e25fSSatish Balay   if (mbs) {
78249b5e25fSSatish Balay     fshift += imax[mbs-1] - ailen[mbs-1];
78349b5e25fSSatish Balay     ai[mbs] = ai[mbs-1] + ailen[mbs-1];
78449b5e25fSSatish Balay   }
78549b5e25fSSatish Balay   /* reset ilen and imax for each row */
78649b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
78749b5e25fSSatish Balay     ailen[i] = imax[i] = ai[i+1] - ai[i];
78849b5e25fSSatish Balay   }
78949b5e25fSSatish Balay   a->s_nz = ai[mbs];
79049b5e25fSSatish Balay 
791b424e231SHong Zhang   /* diagonals may have moved, reset it */
792b424e231SHong Zhang   if (a->diag) {
793b424e231SHong Zhang     ierr = PetscMemcpy(a->diag,ai,(mbs+1)*sizeof(int));CHKERRQ(ierr);
79449b5e25fSSatish Balay   }
795b0a32e0cSBarry Smith   PetscLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Matrix size: %d X %d, block size %d; storage space: %d unneeded, %d used\n",
796c464158bSHong Zhang            m,A->m,a->bs,fshift*bs2,a->s_nz*bs2);
797b0a32e0cSBarry Smith   PetscLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Number of mallocs during MatSetValues is %d\n",
79849b5e25fSSatish Balay            a->reallocs);
799b0a32e0cSBarry Smith   PetscLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Most nonzeros blocks in any row is %d\n",rmax);
80049b5e25fSSatish Balay   a->reallocs          = 0;
80149b5e25fSSatish Balay   A->info.nz_unneeded  = (PetscReal)fshift*bs2;
80249b5e25fSSatish Balay 
803cf4441caSHong Zhang #if defined(PETSC_HAVE_SPOOLES)
8042c535e4dSHong Zhang   ierr = PetscOptionsHasName(A->prefix,"-mat_sbaij_spooles",&flag);CHKERRQ(ierr);
805cf4441caSHong Zhang   if (flag) { ierr = MatUseSpooles_SeqSBAIJ(A);CHKERRQ(ierr); }
806cf4441caSHong Zhang #endif
807cf4441caSHong Zhang 
80849b5e25fSSatish Balay   PetscFunctionReturn(0);
80949b5e25fSSatish Balay }
81049b5e25fSSatish Balay 
81149b5e25fSSatish Balay /*
81249b5e25fSSatish Balay    This function returns an array of flags which indicate the locations of contiguous
81349b5e25fSSatish Balay    blocks that should be zeroed. for eg: if bs = 3  and is = [0,1,2,3,5,6,7,8,9]
81449b5e25fSSatish Balay    then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)]
81549b5e25fSSatish Balay    Assume: sizes should be long enough to hold all the values.
81649b5e25fSSatish Balay */
8174a2ae208SSatish Balay #undef __FUNCT__
8184a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqSBAIJ_Check_Blocks"
81949b5e25fSSatish Balay static int MatZeroRows_SeqSBAIJ_Check_Blocks(int idx[],int n,int bs,int sizes[], int *bs_max)
82049b5e25fSSatish Balay {
82149b5e25fSSatish Balay   int        i,j,k,row;
82249b5e25fSSatish Balay   PetscTruth flg;
82349b5e25fSSatish Balay 
82449b5e25fSSatish Balay   PetscFunctionBegin;
82549b5e25fSSatish Balay   for (i=0,j=0; i<n; j++) {
82649b5e25fSSatish Balay     row = idx[i];
82749b5e25fSSatish Balay     if (row%bs!=0) { /* Not the begining of a block */
82849b5e25fSSatish Balay       sizes[j] = 1;
82949b5e25fSSatish Balay       i++;
83049b5e25fSSatish Balay     } else if (i+bs > n) { /* Beginning of a block, but complete block doesn't exist (at idx end) */
83149b5e25fSSatish Balay       sizes[j] = 1;         /* Also makes sure atleast 'bs' values exist for next else */
83249b5e25fSSatish Balay       i++;
83349b5e25fSSatish Balay     } else { /* Begining of the block, so check if the complete block exists */
83449b5e25fSSatish Balay       flg = PETSC_TRUE;
83549b5e25fSSatish Balay       for (k=1; k<bs; k++) {
83649b5e25fSSatish Balay         if (row+k != idx[i+k]) { /* break in the block */
83749b5e25fSSatish Balay           flg = PETSC_FALSE;
83849b5e25fSSatish Balay           break;
83949b5e25fSSatish Balay         }
84049b5e25fSSatish Balay       }
84149b5e25fSSatish Balay       if (flg == PETSC_TRUE) { /* No break in the bs */
84249b5e25fSSatish Balay         sizes[j] = bs;
84349b5e25fSSatish Balay         i+= bs;
84449b5e25fSSatish Balay       } else {
84549b5e25fSSatish Balay         sizes[j] = 1;
84649b5e25fSSatish Balay         i++;
84749b5e25fSSatish Balay       }
84849b5e25fSSatish Balay     }
84949b5e25fSSatish Balay   }
85049b5e25fSSatish Balay   *bs_max = j;
85149b5e25fSSatish Balay   PetscFunctionReturn(0);
85249b5e25fSSatish Balay }
85349b5e25fSSatish Balay 
8544a2ae208SSatish Balay #undef __FUNCT__
8554a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqSBAIJ"
85687828ca2SBarry Smith int MatZeroRows_SeqSBAIJ(Mat A,IS is,PetscScalar *diag)
85749b5e25fSSatish Balay {
85849b5e25fSSatish Balay   Mat_SeqSBAIJ  *sbaij=(Mat_SeqSBAIJ*)A->data;
85949b5e25fSSatish Balay   int           ierr,i,j,k,count,is_n,*is_idx,*rows;
86049b5e25fSSatish Balay   int           bs=sbaij->bs,bs2=sbaij->bs2,*sizes,row,bs_max;
86187828ca2SBarry Smith   PetscScalar   zero = 0.0;
86249b5e25fSSatish Balay   MatScalar     *aa;
86349b5e25fSSatish Balay 
86449b5e25fSSatish Balay   PetscFunctionBegin;
86549b5e25fSSatish Balay   /* Make a copy of the IS and  sort it */
86649b5e25fSSatish Balay   ierr = ISGetSize(is,&is_n);CHKERRQ(ierr);
86749b5e25fSSatish Balay   ierr = ISGetIndices(is,&is_idx);CHKERRQ(ierr);
86849b5e25fSSatish Balay 
86949b5e25fSSatish Balay   /* allocate memory for rows,sizes */
870b0a32e0cSBarry Smith   ierr = PetscMalloc((3*is_n+1)*sizeof(int),&rows);CHKERRQ(ierr);
87149b5e25fSSatish Balay   sizes = rows + is_n;
87249b5e25fSSatish Balay 
87349b5e25fSSatish Balay   /* initialize copy IS values to rows, and sort them */
87449b5e25fSSatish Balay   for (i=0; i<is_n; i++) { rows[i] = is_idx[i]; }
87549b5e25fSSatish Balay   ierr = PetscSortInt(is_n,rows);CHKERRQ(ierr);
87649b5e25fSSatish Balay   if (sbaij->keepzeroedrows) { /* do not change nonzero structure */
87749b5e25fSSatish Balay     for (i=0; i<is_n; i++) { sizes[i] = 1; } /* sizes: size of blocks, = 1 or bs */
87849b5e25fSSatish Balay     bs_max = is_n;            /* bs_max: num. of contiguous block row in the row */
87949b5e25fSSatish Balay   } else {
88049b5e25fSSatish Balay     ierr = MatZeroRows_SeqSBAIJ_Check_Blocks(rows,is_n,bs,sizes,&bs_max);CHKERRQ(ierr);
88149b5e25fSSatish Balay   }
88249b5e25fSSatish Balay   ierr = ISRestoreIndices(is,&is_idx);CHKERRQ(ierr);
88349b5e25fSSatish Balay 
88449b5e25fSSatish Balay   for (i=0,j=0; i<bs_max; j+=sizes[i],i++) {
88549b5e25fSSatish Balay     row   = rows[j];                  /* row to be zeroed */
886c464158bSHong Zhang     if (row < 0 || row > A->m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"row %d out of range",row);
88749b5e25fSSatish Balay     count = (sbaij->i[row/bs +1] - sbaij->i[row/bs])*bs; /* num. of elements in the row */
88849b5e25fSSatish Balay     aa    = sbaij->a + sbaij->i[row/bs]*bs2 + (row%bs);
88949b5e25fSSatish Balay     if (sizes[i] == bs && !sbaij->keepzeroedrows) {
89049b5e25fSSatish Balay       if (diag) {
89149b5e25fSSatish Balay         if (sbaij->ilen[row/bs] > 0) {
89249b5e25fSSatish Balay           sbaij->ilen[row/bs] = 1;
89349b5e25fSSatish Balay           sbaij->j[sbaij->i[row/bs]] = row/bs;
89449b5e25fSSatish Balay           ierr = PetscMemzero(aa,count*bs*sizeof(MatScalar));CHKERRQ(ierr);
89549b5e25fSSatish Balay         }
89649b5e25fSSatish Balay         /* Now insert all the diagoanl values for this bs */
89749b5e25fSSatish Balay         for (k=0; k<bs; k++) {
89849b5e25fSSatish Balay           ierr = (*A->ops->setvalues)(A,1,rows+j+k,1,rows+j+k,diag,INSERT_VALUES);CHKERRQ(ierr);
89949b5e25fSSatish Balay         }
90049b5e25fSSatish Balay       } else { /* (!diag) */
90149b5e25fSSatish Balay         sbaij->ilen[row/bs] = 0;
90249b5e25fSSatish Balay       } /* end (!diag) */
90349b5e25fSSatish Balay     } else { /* (sizes[i] != bs), broken block */
90449b5e25fSSatish Balay #if defined (PETSC_USE_DEBUG)
90529bbc08cSBarry Smith       if (sizes[i] != 1) SETERRQ(1,"Internal Error. Value should be 1");
90649b5e25fSSatish Balay #endif
90749b5e25fSSatish Balay       for (k=0; k<count; k++) {
90849b5e25fSSatish Balay         aa[0] = zero;
90949b5e25fSSatish Balay         aa+=bs;
91049b5e25fSSatish Balay       }
91149b5e25fSSatish Balay       if (diag) {
91249b5e25fSSatish Balay         ierr = (*A->ops->setvalues)(A,1,rows+j,1,rows+j,diag,INSERT_VALUES);CHKERRQ(ierr);
91349b5e25fSSatish Balay       }
91449b5e25fSSatish Balay     }
91549b5e25fSSatish Balay   }
91649b5e25fSSatish Balay 
91749b5e25fSSatish Balay   ierr = PetscFree(rows);CHKERRQ(ierr);
91849b5e25fSSatish Balay   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
91949b5e25fSSatish Balay   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
92049b5e25fSSatish Balay   PetscFunctionReturn(0);
92149b5e25fSSatish Balay }
92249b5e25fSSatish Balay 
92349b5e25fSSatish Balay /* Only add/insert a(i,j) with i<=j (blocks).
92449b5e25fSSatish Balay    Any a(i,j) with i>j input by user is ingored.
92549b5e25fSSatish Balay */
92649b5e25fSSatish Balay 
9274a2ae208SSatish Balay #undef __FUNCT__
9284a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqSBAIJ"
92987828ca2SBarry Smith int MatSetValues_SeqSBAIJ(Mat A,int m,int *im,int n,int *in,PetscScalar *v,InsertMode is)
93049b5e25fSSatish Balay {
93149b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
93249b5e25fSSatish Balay   int         *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted=a->sorted;
93349b5e25fSSatish Balay   int         *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented;
93449b5e25fSSatish Balay   int         *aj=a->j,nonew=a->nonew,bs=a->bs,brow,bcol;
93549b5e25fSSatish Balay   int         ridx,cidx,bs2=a->bs2,ierr;
93649b5e25fSSatish Balay   MatScalar   *ap,value,*aa=a->a,*bap;
93749b5e25fSSatish Balay 
93849b5e25fSSatish Balay   PetscFunctionBegin;
93949b5e25fSSatish Balay 
94049b5e25fSSatish Balay   for (k=0; k<m; k++) { /* loop over added rows */
94149b5e25fSSatish Balay     row  = im[k];       /* row number */
94249b5e25fSSatish Balay     brow = row/bs;      /* block row number */
94349b5e25fSSatish Balay     if (row < 0) continue;
94449b5e25fSSatish Balay #if defined(PETSC_USE_BOPT_g)
945c464158bSHong Zhang     if (row >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",row,A->m);
94649b5e25fSSatish Balay #endif
94749b5e25fSSatish Balay     rp   = aj + ai[brow]; /*ptr to beginning of column value of the row block*/
94849b5e25fSSatish Balay     ap   = aa + bs2*ai[brow]; /*ptr to beginning of element value of the row block*/
94949b5e25fSSatish Balay     rmax = imax[brow];  /* maximum space allocated for this row */
95049b5e25fSSatish Balay     nrow = ailen[brow]; /* actual length of this row */
95149b5e25fSSatish Balay     low  = 0;
95249b5e25fSSatish Balay 
95349b5e25fSSatish Balay     for (l=0; l<n; l++) { /* loop over added columns */
95449b5e25fSSatish Balay       if (in[l] < 0) continue;
95549b5e25fSSatish Balay #if defined(PETSC_USE_BOPT_g)
956c464158bSHong Zhang       if (in[l] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",in[l],A->m);
95749b5e25fSSatish Balay #endif
95849b5e25fSSatish Balay       col = in[l];
95949b5e25fSSatish Balay       bcol = col/bs;              /* block col number */
96049b5e25fSSatish Balay 
96149b5e25fSSatish Balay       if (brow <= bcol){
96249b5e25fSSatish Balay         ridx = row % bs; cidx = col % bs; /*row and col index inside the block */
9638549e402SHong Zhang         if ((brow==bcol && ridx<=cidx) || (brow<bcol)){
96449b5e25fSSatish Balay           /* element value a(k,l) */
96549b5e25fSSatish Balay           if (roworiented) {
96649b5e25fSSatish Balay             value = v[l + k*n];
96749b5e25fSSatish Balay           } else {
96849b5e25fSSatish Balay             value = v[k + l*m];
96949b5e25fSSatish Balay           }
97049b5e25fSSatish Balay 
97149b5e25fSSatish Balay           /* move pointer bap to a(k,l) quickly and add/insert value */
97249b5e25fSSatish Balay           if (!sorted) low = 0; high = nrow;
97349b5e25fSSatish Balay           while (high-low > 7) {
97449b5e25fSSatish Balay             t = (low+high)/2;
97549b5e25fSSatish Balay             if (rp[t] > bcol) high = t;
97649b5e25fSSatish Balay             else              low  = t;
97749b5e25fSSatish Balay           }
97849b5e25fSSatish Balay           for (i=low; i<high; i++) {
97949b5e25fSSatish Balay             /* printf("The loop of i=low.., rp[%d]=%d\n",i,rp[i]); */
98049b5e25fSSatish Balay             if (rp[i] > bcol) break;
98149b5e25fSSatish Balay             if (rp[i] == bcol) {
98249b5e25fSSatish Balay               bap  = ap +  bs2*i + bs*cidx + ridx;
98349b5e25fSSatish Balay               if (is == ADD_VALUES) *bap += value;
98449b5e25fSSatish Balay               else                  *bap  = value;
9858549e402SHong Zhang               /* for diag block, add/insert its symmetric element a(cidx,ridx) */
9868549e402SHong Zhang               if (brow == bcol && ridx < cidx){
9878549e402SHong Zhang                 bap  = ap +  bs2*i + bs*ridx + cidx;
9888549e402SHong Zhang                 if (is == ADD_VALUES) *bap += value;
9898549e402SHong Zhang                 else                  *bap  = value;
9908549e402SHong Zhang               }
99149b5e25fSSatish Balay               goto noinsert1;
99249b5e25fSSatish Balay             }
99349b5e25fSSatish Balay           }
99449b5e25fSSatish Balay 
99549b5e25fSSatish Balay       if (nonew == 1) goto noinsert1;
99629bbc08cSBarry Smith       else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix");
99749b5e25fSSatish Balay       if (nrow >= rmax) {
99849b5e25fSSatish Balay         /* there is no extra room in row, therefore enlarge */
99949b5e25fSSatish Balay         int       new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j;
100049b5e25fSSatish Balay         MatScalar *new_a;
100149b5e25fSSatish Balay 
100229bbc08cSBarry Smith         if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix");
100349b5e25fSSatish Balay 
100449b5e25fSSatish Balay         /* Malloc new storage space */
100549b5e25fSSatish Balay         len   = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int);
100682502324SSatish Balay         ierr  = PetscMalloc(len,&new_a);CHKERRQ(ierr);
100749b5e25fSSatish Balay         new_j = (int*)(new_a + bs2*new_nz);
100849b5e25fSSatish Balay         new_i = new_j + new_nz;
100949b5e25fSSatish Balay 
101049b5e25fSSatish Balay         /* copy over old data into new slots */
101149b5e25fSSatish Balay         for (ii=0; ii<brow+1; ii++) {new_i[ii] = ai[ii];}
101249b5e25fSSatish Balay         for (ii=brow+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;}
101349b5e25fSSatish Balay         ierr = PetscMemcpy(new_j,aj,(ai[brow]+nrow)*sizeof(int));CHKERRQ(ierr);
101449b5e25fSSatish Balay         len  = (new_nz - CHUNKSIZE - ai[brow] - nrow);
101549b5e25fSSatish Balay         ierr = PetscMemcpy(new_j+ai[brow]+nrow+CHUNKSIZE,aj+ai[brow]+nrow,len*sizeof(int));CHKERRQ(ierr);
101649b5e25fSSatish Balay         ierr = PetscMemcpy(new_a,aa,(ai[brow]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr);
101749b5e25fSSatish Balay         ierr = PetscMemzero(new_a+bs2*(ai[brow]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr);
101849b5e25fSSatish Balay         ierr = PetscMemcpy(new_a+bs2*(ai[brow]+nrow+CHUNKSIZE),aa+bs2*(ai[brow]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr);
101949b5e25fSSatish Balay         /* free up old matrix storage */
102049b5e25fSSatish Balay         ierr = PetscFree(a->a);CHKERRQ(ierr);
102149b5e25fSSatish Balay         if (!a->singlemalloc) {
102249b5e25fSSatish Balay           ierr = PetscFree(a->i);CHKERRQ(ierr);
102349b5e25fSSatish Balay           ierr = PetscFree(a->j);CHKERRQ(ierr);
102449b5e25fSSatish Balay         }
102549b5e25fSSatish Balay         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
102649b5e25fSSatish Balay         a->singlemalloc = PETSC_TRUE;
102749b5e25fSSatish Balay 
102849b5e25fSSatish Balay         rp   = aj + ai[brow]; ap = aa + bs2*ai[brow];
102949b5e25fSSatish Balay         rmax = imax[brow] = imax[brow] + CHUNKSIZE;
1030b0a32e0cSBarry Smith         PetscLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar)));
103149b5e25fSSatish Balay         a->s_maxnz += bs2*CHUNKSIZE;
103249b5e25fSSatish Balay         a->reallocs++;
103349b5e25fSSatish Balay         a->s_nz++;
103449b5e25fSSatish Balay       }
103549b5e25fSSatish Balay 
103649b5e25fSSatish Balay       N = nrow++ - 1;
103749b5e25fSSatish Balay       /* shift up all the later entries in this row */
103849b5e25fSSatish Balay       for (ii=N; ii>=i; ii--) {
103949b5e25fSSatish Balay         rp[ii+1] = rp[ii];
104049b5e25fSSatish Balay         ierr     = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr);
104149b5e25fSSatish Balay       }
104249b5e25fSSatish Balay       if (N>=i) {
104349b5e25fSSatish Balay         ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr);
104449b5e25fSSatish Balay       }
104549b5e25fSSatish Balay       rp[i]                      = bcol;
104649b5e25fSSatish Balay       ap[bs2*i + bs*cidx + ridx] = value;
104749b5e25fSSatish Balay       noinsert1:;
104849b5e25fSSatish Balay       low = i;
104949b5e25fSSatish Balay       /* } */
10508549e402SHong Zhang         }
105149b5e25fSSatish Balay       } /* end of if .. if..  */
105249b5e25fSSatish Balay     }                     /* end of loop over added columns */
105349b5e25fSSatish Balay     ailen[brow] = nrow;
105449b5e25fSSatish Balay   }                       /* end of loop over added rows */
105549b5e25fSSatish Balay 
105649b5e25fSSatish Balay   PetscFunctionReturn(0);
105749b5e25fSSatish Balay }
105849b5e25fSSatish Balay 
1059915354c9SHong Zhang extern int MatCholeskyFactorSymbolic_SeqSBAIJ(Mat,IS,PetscReal,Mat*);
1060915354c9SHong Zhang extern int MatCholeskyFactor_SeqSBAIJ(Mat,IS,PetscReal);
106149b5e25fSSatish Balay extern int MatIncreaseOverlap_SeqSBAIJ(Mat,int,IS*,int);
106249b5e25fSSatish Balay extern int MatGetSubMatrix_SeqSBAIJ(Mat,IS,IS,int,MatReuse,Mat*);
106349b5e25fSSatish Balay extern int MatGetSubMatrices_SeqSBAIJ(Mat,int,IS*,IS*,MatReuse,Mat**);
106449b5e25fSSatish Balay extern int MatMultTranspose_SeqSBAIJ(Mat,Vec,Vec);
106549b5e25fSSatish Balay extern int MatMultTransposeAdd_SeqSBAIJ(Mat,Vec,Vec,Vec);
106687828ca2SBarry Smith extern int MatScale_SeqSBAIJ(PetscScalar*,Mat);
106749b5e25fSSatish Balay extern int MatNorm_SeqSBAIJ(Mat,NormType,PetscReal *);
106849b5e25fSSatish Balay extern int MatEqual_SeqSBAIJ(Mat,Mat,PetscTruth*);
106949b5e25fSSatish Balay extern int MatGetDiagonal_SeqSBAIJ(Mat,Vec);
107049b5e25fSSatish Balay extern int MatDiagonalScale_SeqSBAIJ(Mat,Vec,Vec);
107149b5e25fSSatish Balay extern int MatGetInfo_SeqSBAIJ(Mat,MatInfoType,MatInfo *);
107249b5e25fSSatish Balay extern int MatZeroEntries_SeqSBAIJ(Mat);
107313a02c98SHong Zhang extern int MatGetRowMax_SeqSBAIJ(Mat,Vec);
107449b5e25fSSatish Balay 
107549b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_N(Mat,Vec,Vec);
107649b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_1(Mat,Vec,Vec);
107749b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_2(Mat,Vec,Vec);
107849b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_3(Mat,Vec,Vec);
107949b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_4(Mat,Vec,Vec);
108049b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_5(Mat,Vec,Vec);
108149b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_6(Mat,Vec,Vec);
108249b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_7(Mat,Vec,Vec);
108349b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_7(Mat,Vec,Vec);
108449b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_6(Mat,Vec,Vec);
108549b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_5(Mat,Vec,Vec);
108649b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_4(Mat,Vec,Vec);
108749b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_3(Mat,Vec,Vec);
108849b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_2(Mat,Vec,Vec);
108949b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_1(Mat,Vec,Vec);
109049b5e25fSSatish Balay 
109107f98182SSatish Balay extern int MatSolve_SeqSBAIJ_1_NaturalOrdering(Mat,Vec,Vec);
109207f98182SSatish Balay extern int MatSolve_SeqSBAIJ_2_NaturalOrdering(Mat,Vec,Vec);
109307f98182SSatish Balay extern int MatSolve_SeqSBAIJ_3_NaturalOrdering(Mat,Vec,Vec);
109407f98182SSatish Balay extern int MatSolve_SeqSBAIJ_4_NaturalOrdering(Mat,Vec,Vec);
109507f98182SSatish Balay extern int MatSolve_SeqSBAIJ_5_NaturalOrdering(Mat,Vec,Vec);
109607f98182SSatish Balay extern int MatSolve_SeqSBAIJ_6_NaturalOrdering(Mat,Vec,Vec);
109707f98182SSatish Balay extern int MatSolve_SeqSBAIJ_7_NaturalOrdering(Mat,Vec,Vec);
109807f98182SSatish Balay extern int MatSolve_SeqSBAIJ_N_NaturalOrdering(Mat,Vec,Vec);
109907f98182SSatish Balay 
11005f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat,Mat*);
11015f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_1(Mat,Mat*);
11025f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_2(Mat,Mat*);
11035f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_3(Mat,Mat*);
11045f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_4(Mat,Mat*);
11055f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_5(Mat,Mat*);
11065f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_6(Mat,Mat*);
110757d5e339SHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_7(Mat,Mat*);
110849b5e25fSSatish Balay 
110949b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_1(Mat,Vec,Vec);
111049b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_2(Mat,Vec,Vec);
111149b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_3(Mat,Vec,Vec);
111249b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_4(Mat,Vec,Vec);
111349b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_5(Mat,Vec,Vec);
111449b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_6(Mat,Vec,Vec);
111549b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_7(Mat,Vec,Vec);
111649b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_N(Mat,Vec,Vec);
111749b5e25fSSatish Balay 
111849b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_1(Mat,Vec,Vec,Vec);
111949b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_2(Mat,Vec,Vec,Vec);
112049b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_3(Mat,Vec,Vec,Vec);
112149b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_4(Mat,Vec,Vec,Vec);
112249b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_5(Mat,Vec,Vec,Vec);
112349b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_6(Mat,Vec,Vec,Vec);
112449b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_7(Mat,Vec,Vec,Vec);
112549b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_N(Mat,Vec,Vec,Vec);
112649b5e25fSSatish Balay 
11274d101231SSatish Balay #ifdef HAVE_MatICCFactor
11281a3463dfSHong Zhang /* modefied from MatILUFactor_SeqSBAIJ, needs further work!  */
11294a2ae208SSatish Balay #undef __FUNCT__
11304d101231SSatish Balay #define __FUNCT__ "MatICCFactor_SeqSBAIJ"
11314d101231SSatish Balay int MatICCFactor_SeqSBAIJ(Mat inA,IS row,PetscReal fill,int level)
113249b5e25fSSatish Balay {
11334ccecd49SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)inA->data;
113449b5e25fSSatish Balay   Mat         outA;
113549b5e25fSSatish Balay   int         ierr;
113649b5e25fSSatish Balay   PetscTruth  row_identity,col_identity;
113749b5e25fSSatish Balay 
113849b5e25fSSatish Balay   PetscFunctionBegin;
11391a3463dfSHong Zhang   /*
114029bbc08cSBarry Smith   if (level != 0) SETERRQ(PETSC_ERR_SUP,"Only levels = 0 supported for in-place ILU");
114149b5e25fSSatish Balay   ierr = ISIdentity(row,&row_identity);CHKERRQ(ierr);
114249b5e25fSSatish Balay   ierr = ISIdentity(col,&col_identity);CHKERRQ(ierr);
114349b5e25fSSatish Balay   if (!row_identity || !col_identity) {
114429bbc08cSBarry Smith     SETERRQ(1,"Row and column permutations must be identity for in-place ICC");
114549b5e25fSSatish Balay   }
11461a3463dfSHong Zhang   */
114749b5e25fSSatish Balay 
114849b5e25fSSatish Balay   outA          = inA;
11491260e1f8SHong Zhang   inA->factor   = FACTOR_CHOLESKY;
115049b5e25fSSatish Balay 
115149b5e25fSSatish Balay   if (!a->diag) {
11521a3463dfSHong Zhang     ierr = MatMarkDiagonal_SeqSBAIJ(inA);CHKERRQ(ierr);
115349b5e25fSSatish Balay   }
115449b5e25fSSatish Balay   /*
115549b5e25fSSatish Balay       Blocksize 2, 3, 4, 5, 6 and 7 have a special faster factorization/solver
115649b5e25fSSatish Balay       for ILU(0) factorization with natural ordering
115749b5e25fSSatish Balay   */
115849b5e25fSSatish Balay   switch (a->bs) {
115949b5e25fSSatish Balay   case 1:
11600fe9e5beSBarry Smith     inA->ops->solve            = MatSolve_SeqSBAIJ_1_NaturalOrdering;
116107f98182SSatish Balay     inA->ops->solvetranspose   = MatSolve_SeqSBAIJ_1_NaturalOrdering;
11624d101231SSatish Balay     PetscLoginfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering solvetrans BS=1\n");
116349b5e25fSSatish Balay   case 2:
11641a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering;
11651a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_2_NaturalOrdering;
11661a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_2_NaturalOrdering;
11674d101231SSatish Balay     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=2\n");
116849b5e25fSSatish Balay     break;
116949b5e25fSSatish Balay   case 3:
11701a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering;
11711a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_3_NaturalOrdering;
11721a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_3_NaturalOrdering;
11734d101231SSatish Balay     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=3\n");
117449b5e25fSSatish Balay     break;
117549b5e25fSSatish Balay   case 4:
11761a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering;
11771a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_4_NaturalOrdering;
11781a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_4_NaturalOrdering;
11794d101231SSatish Balay     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=4\n");
118049b5e25fSSatish Balay     break;
118149b5e25fSSatish Balay   case 5:
11821a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering;
11831a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_5_NaturalOrdering;
11841a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_5_NaturalOrdering;
11854d101231SSatish Balay     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=5\n");
118649b5e25fSSatish Balay     break;
118749b5e25fSSatish Balay   case 6:
11881a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering;
11891a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_6_NaturalOrdering;
11901a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_6_NaturalOrdering;
11914d101231SSatish Balay     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=6\n");
119249b5e25fSSatish Balay     break;
119349b5e25fSSatish Balay   case 7:
11941a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering;
11951a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_7_NaturalOrdering;
11961a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_7_NaturalOrdering;
11974d101231SSatish Balay     PetscLogInfo(inA,"MatICCFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=7\n");
119849b5e25fSSatish Balay     break;
119949b5e25fSSatish Balay   default:
120049b5e25fSSatish Balay     a->row        = row;
12011a3463dfSHong Zhang     a->icol       = col;
120249b5e25fSSatish Balay     ierr          = PetscObjectReference((PetscObject)row);CHKERRQ(ierr);
120349b5e25fSSatish Balay     ierr          = PetscObjectReference((PetscObject)col);CHKERRQ(ierr);
120449b5e25fSSatish Balay 
120549b5e25fSSatish Balay     /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */
120649b5e25fSSatish Balay     ierr = ISInvertPermutation(col,PETSC_DECIDE, &(a->icol));CHKERRQ(ierr);
1207b0a32e0cSBarry Smith     PetscLogObjectParent(inA,a->icol);
120849b5e25fSSatish Balay 
120949b5e25fSSatish Balay     if (!a->solve_work) {
121087828ca2SBarry Smith       ierr = PetscMalloc((A->m+a->bs)*sizeof(PetscScalar),&a->solve_work);CHKERRQ(ierr);
121187828ca2SBarry Smith       PetscLogObjectMemory(inA,(A->m+a->bs)*sizeof(PetscScalar));
121249b5e25fSSatish Balay     }
121349b5e25fSSatish Balay   }
121449b5e25fSSatish Balay 
12151a3463dfSHong Zhang   ierr = MatCholeskyFactorNumeric(inA,&outA);CHKERRQ(ierr);
121649b5e25fSSatish Balay 
121749b5e25fSSatish Balay   PetscFunctionReturn(0);
121849b5e25fSSatish Balay }
1219950f1e5bSHong Zhang #endif
1220950f1e5bSHong Zhang 
12214a2ae208SSatish Balay #undef __FUNCT__
12224a2ae208SSatish Balay #define __FUNCT__ "MatPrintHelp_SeqSBAIJ"
122349b5e25fSSatish Balay int MatPrintHelp_SeqSBAIJ(Mat A)
122449b5e25fSSatish Balay {
122549b5e25fSSatish Balay   static PetscTruth called = PETSC_FALSE;
122649b5e25fSSatish Balay   MPI_Comm          comm = A->comm;
122749b5e25fSSatish Balay   int               ierr;
122849b5e25fSSatish Balay 
122949b5e25fSSatish Balay   PetscFunctionBegin;
123049b5e25fSSatish Balay   if (called) {PetscFunctionReturn(0);} else called = PETSC_TRUE;
123149b5e25fSSatish Balay   ierr = (*PetscHelpPrintf)(comm," Options for MATSEQSBAIJ and MATMPISBAIJ matrix formats (the defaults):\n");CHKERRQ(ierr);
123249b5e25fSSatish Balay   ierr = (*PetscHelpPrintf)(comm,"  -mat_block_size <block_size>\n");CHKERRQ(ierr);
123349b5e25fSSatish Balay   PetscFunctionReturn(0);
123449b5e25fSSatish Balay }
123549b5e25fSSatish Balay 
123649b5e25fSSatish Balay EXTERN_C_BEGIN
12374a2ae208SSatish Balay #undef __FUNCT__
12384a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices_SeqSBAIJ"
123949b5e25fSSatish Balay int MatSeqSBAIJSetColumnIndices_SeqSBAIJ(Mat mat,int *indices)
124049b5e25fSSatish Balay {
1241045c9aa0SHong Zhang   Mat_SeqSBAIJ *baij = (Mat_SeqSBAIJ *)mat->data;
124249b5e25fSSatish Balay   int         i,nz,n;
124349b5e25fSSatish Balay 
124449b5e25fSSatish Balay   PetscFunctionBegin;
1245045c9aa0SHong Zhang   nz = baij->s_maxnz;
1246c464158bSHong Zhang   n  = mat->n;
124749b5e25fSSatish Balay   for (i=0; i<nz; i++) {
124849b5e25fSSatish Balay     baij->j[i] = indices[i];
124949b5e25fSSatish Balay   }
1250045c9aa0SHong Zhang   baij->s_nz = nz;
125149b5e25fSSatish Balay   for (i=0; i<n; i++) {
125249b5e25fSSatish Balay     baij->ilen[i] = baij->imax[i];
125349b5e25fSSatish Balay   }
125449b5e25fSSatish Balay 
125549b5e25fSSatish Balay   PetscFunctionReturn(0);
125649b5e25fSSatish Balay }
125749b5e25fSSatish Balay EXTERN_C_END
125849b5e25fSSatish Balay 
12594a2ae208SSatish Balay #undef __FUNCT__
12604a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetColumnIndices"
126149b5e25fSSatish Balay /*@
126219585528SSatish Balay     MatSeqSBAIJSetColumnIndices - Set the column indices for all the rows
126349b5e25fSSatish Balay        in the matrix.
126449b5e25fSSatish Balay 
126549b5e25fSSatish Balay   Input Parameters:
126619585528SSatish Balay +  mat     - the SeqSBAIJ matrix
126749b5e25fSSatish Balay -  indices - the column indices
126849b5e25fSSatish Balay 
126949b5e25fSSatish Balay   Level: advanced
127049b5e25fSSatish Balay 
127149b5e25fSSatish Balay   Notes:
127249b5e25fSSatish Balay     This can be called if you have precomputed the nonzero structure of the
127349b5e25fSSatish Balay   matrix and want to provide it to the matrix object to improve the performance
127449b5e25fSSatish Balay   of the MatSetValues() operation.
127549b5e25fSSatish Balay 
127649b5e25fSSatish Balay     You MUST have set the correct numbers of nonzeros per row in the call to
127719585528SSatish Balay   MatCreateSeqSBAIJ().
127849b5e25fSSatish Balay 
1279ab9f2c04SSatish Balay     MUST be called before any calls to MatSetValues()
128049b5e25fSSatish Balay 
1281ab9f2c04SSatish Balay .seealso: MatCreateSeqSBAIJ
128249b5e25fSSatish Balay @*/
128349b5e25fSSatish Balay int MatSeqSBAIJSetColumnIndices(Mat mat,int *indices)
128449b5e25fSSatish Balay {
128549b5e25fSSatish Balay   int ierr,(*f)(Mat,int *);
128649b5e25fSSatish Balay 
128749b5e25fSSatish Balay   PetscFunctionBegin;
128849b5e25fSSatish Balay   PetscValidHeaderSpecific(mat,MAT_COOKIE);
1289c134de8dSSatish Balay   ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqSBAIJSetColumnIndices_C",(void (**)(void))&f);CHKERRQ(ierr);
129049b5e25fSSatish Balay   if (f) {
129149b5e25fSSatish Balay     ierr = (*f)(mat,indices);CHKERRQ(ierr);
129249b5e25fSSatish Balay   } else {
129329bbc08cSBarry Smith     SETERRQ(1,"Wrong type of matrix to set column indices");
129449b5e25fSSatish Balay   }
129549b5e25fSSatish Balay   PetscFunctionReturn(0);
129649b5e25fSSatish Balay }
129749b5e25fSSatish Balay 
12984a2ae208SSatish Balay #undef __FUNCT__
12994a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqSBAIJ"
1300273d9f13SBarry Smith int MatSetUpPreallocation_SeqSBAIJ(Mat A)
1301273d9f13SBarry Smith {
1302273d9f13SBarry Smith   int        ierr;
1303273d9f13SBarry Smith 
1304273d9f13SBarry Smith   PetscFunctionBegin;
1305273d9f13SBarry Smith   ierr =  MatSeqSBAIJSetPreallocation(A,1,PETSC_DEFAULT,0);CHKERRQ(ierr);
1306273d9f13SBarry Smith   PetscFunctionReturn(0);
1307273d9f13SBarry Smith }
1308273d9f13SBarry Smith 
130942ee4b1aSHong Zhang #include "petscblaslapack.h"
131042ee4b1aSHong Zhang #undef __FUNCT__
131142ee4b1aSHong Zhang #define __FUNCT__ "MatAXPY_SeqSBAIJ"
131242ee4b1aSHong Zhang int MatAXPY_SeqSBAIJ(PetscScalar *a,Mat X,Mat Y,MatStructure str)
131342ee4b1aSHong Zhang {
131442ee4b1aSHong Zhang   int          ierr,one=1;
131542ee4b1aSHong Zhang   Mat_SeqSBAIJ *x  = (Mat_SeqSBAIJ *)X->data,*y = (Mat_SeqSBAIJ *)Y->data;
131642ee4b1aSHong Zhang 
131742ee4b1aSHong Zhang   PetscFunctionBegin;
131842ee4b1aSHong Zhang   if (str == SAME_NONZERO_PATTERN) {
131942ee4b1aSHong Zhang     BLaxpy_(&x->s_nz,a,x->a,&one,y->a,&one);
132042ee4b1aSHong Zhang   } else {
132142ee4b1aSHong Zhang     ierr = MatAXPY_Basic(a,X,Y,str);CHKERRQ(ierr);
132242ee4b1aSHong Zhang   }
132342ee4b1aSHong Zhang   PetscFunctionReturn(0);
132442ee4b1aSHong Zhang }
132542ee4b1aSHong Zhang 
132649b5e25fSSatish Balay /* -------------------------------------------------------------------*/
132749b5e25fSSatish Balay static struct _MatOps MatOps_Values = {MatSetValues_SeqSBAIJ,
132849b5e25fSSatish Balay        MatGetRow_SeqSBAIJ,
132949b5e25fSSatish Balay        MatRestoreRow_SeqSBAIJ,
133049b5e25fSSatish Balay        MatMult_SeqSBAIJ_N,
133149b5e25fSSatish Balay        MatMultAdd_SeqSBAIJ_N,
133249b5e25fSSatish Balay        MatMultTranspose_SeqSBAIJ,
133349b5e25fSSatish Balay        MatMultTransposeAdd_SeqSBAIJ,
133449b5e25fSSatish Balay        MatSolve_SeqSBAIJ_N,
133549b5e25fSSatish Balay        0,
133649b5e25fSSatish Balay        0,
133749b5e25fSSatish Balay        0,
133849b5e25fSSatish Balay        0,
13395f4ef2deSHong Zhang        MatCholeskyFactor_SeqSBAIJ,
1340d06b337dSHong Zhang        MatRelax_SeqSBAIJ,
134149b5e25fSSatish Balay        MatTranspose_SeqSBAIJ,
134249b5e25fSSatish Balay        MatGetInfo_SeqSBAIJ,
134349b5e25fSSatish Balay        MatEqual_SeqSBAIJ,
134449b5e25fSSatish Balay        MatGetDiagonal_SeqSBAIJ,
134549b5e25fSSatish Balay        MatDiagonalScale_SeqSBAIJ,
134649b5e25fSSatish Balay        MatNorm_SeqSBAIJ,
134749b5e25fSSatish Balay        0,
134849b5e25fSSatish Balay        MatAssemblyEnd_SeqSBAIJ,
134949b5e25fSSatish Balay        0,
135049b5e25fSSatish Balay        MatSetOption_SeqSBAIJ,
135149b5e25fSSatish Balay        MatZeroEntries_SeqSBAIJ,
135249b5e25fSSatish Balay        MatZeroRows_SeqSBAIJ,
135349b5e25fSSatish Balay        0,
135449b5e25fSSatish Balay        0,
13555f4ef2deSHong Zhang        MatCholeskyFactorSymbolic_SeqSBAIJ,
13565f4ef2deSHong Zhang        MatCholeskyFactorNumeric_SeqSBAIJ_N,
1357273d9f13SBarry Smith        MatSetUpPreallocation_SeqSBAIJ,
1358c464158bSHong Zhang        0,
13594d101231SSatish Balay        MatICCFactorSymbolic_SeqSBAIJ,
136049b5e25fSSatish Balay        0,
136149b5e25fSSatish Balay        0,
136249b5e25fSSatish Balay        MatDuplicate_SeqSBAIJ,
136349b5e25fSSatish Balay        0,
136449b5e25fSSatish Balay        0,
136549b5e25fSSatish Balay        0,
1366950f1e5bSHong Zhang        0,
136742ee4b1aSHong Zhang        MatAXPY_SeqSBAIJ,
136849b5e25fSSatish Balay        MatGetSubMatrices_SeqSBAIJ,
136949b5e25fSSatish Balay        MatIncreaseOverlap_SeqSBAIJ,
137049b5e25fSSatish Balay        MatGetValues_SeqSBAIJ,
137149b5e25fSSatish Balay        0,
137249b5e25fSSatish Balay        MatPrintHelp_SeqSBAIJ,
137349b5e25fSSatish Balay        MatScale_SeqSBAIJ,
137449b5e25fSSatish Balay        0,
137549b5e25fSSatish Balay        0,
137649b5e25fSSatish Balay        0,
137749b5e25fSSatish Balay        MatGetBlockSize_SeqSBAIJ,
137849b5e25fSSatish Balay        MatGetRowIJ_SeqSBAIJ,
137949b5e25fSSatish Balay        MatRestoreRowIJ_SeqSBAIJ,
138049b5e25fSSatish Balay        0,
138149b5e25fSSatish Balay        0,
138249b5e25fSSatish Balay        0,
138349b5e25fSSatish Balay        0,
138449b5e25fSSatish Balay        0,
138549b5e25fSSatish Balay        0,
138649b5e25fSSatish Balay        MatSetValuesBlocked_SeqSBAIJ,
138749b5e25fSSatish Balay        MatGetSubMatrix_SeqSBAIJ,
138849b5e25fSSatish Balay        0,
138949b5e25fSSatish Balay        0,
13908a124369SBarry Smith        MatGetPetscMaps_Petsc,
1391d959ec07SHong Zhang        0,
1392d959ec07SHong Zhang        0,
1393d959ec07SHong Zhang        0,
1394d959ec07SHong Zhang        0,
1395d959ec07SHong Zhang        0,
1396d959ec07SHong Zhang        0,
139724d5174aSHong Zhang        MatGetRowMax_SeqSBAIJ};
139849b5e25fSSatish Balay 
139949b5e25fSSatish Balay EXTERN_C_BEGIN
14004a2ae208SSatish Balay #undef __FUNCT__
14014a2ae208SSatish Balay #define __FUNCT__ "MatStoreValues_SeqSBAIJ"
140249b5e25fSSatish Balay int MatStoreValues_SeqSBAIJ(Mat mat)
140349b5e25fSSatish Balay {
14044afc71dfSHong Zhang   Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data;
1405c464158bSHong Zhang   int         nz = aij->i[mat->m]*aij->bs*aij->bs2;
140649b5e25fSSatish Balay   int         ierr;
140749b5e25fSSatish Balay 
140849b5e25fSSatish Balay   PetscFunctionBegin;
140949b5e25fSSatish Balay   if (aij->nonew != 1) {
141029bbc08cSBarry Smith     SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
141149b5e25fSSatish Balay   }
141249b5e25fSSatish Balay 
141349b5e25fSSatish Balay   /* allocate space for values if not already there */
141449b5e25fSSatish Balay   if (!aij->saved_values) {
141587828ca2SBarry Smith     ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&aij->saved_values);CHKERRQ(ierr);
141649b5e25fSSatish Balay   }
141749b5e25fSSatish Balay 
141849b5e25fSSatish Balay   /* copy values over */
141987828ca2SBarry Smith   ierr = PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(PetscScalar));CHKERRQ(ierr);
142049b5e25fSSatish Balay   PetscFunctionReturn(0);
142149b5e25fSSatish Balay }
142249b5e25fSSatish Balay EXTERN_C_END
142349b5e25fSSatish Balay 
142449b5e25fSSatish Balay EXTERN_C_BEGIN
14254a2ae208SSatish Balay #undef __FUNCT__
14264a2ae208SSatish Balay #define __FUNCT__ "MatRetrieveValues_SeqSBAIJ"
142749b5e25fSSatish Balay int MatRetrieveValues_SeqSBAIJ(Mat mat)
142849b5e25fSSatish Balay {
14294afc71dfSHong Zhang   Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data;
1430c464158bSHong Zhang   int         nz = aij->i[mat->m]*aij->bs*aij->bs2,ierr;
143149b5e25fSSatish Balay 
143249b5e25fSSatish Balay   PetscFunctionBegin;
143349b5e25fSSatish Balay   if (aij->nonew != 1) {
143429bbc08cSBarry Smith     SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
143549b5e25fSSatish Balay   }
143649b5e25fSSatish Balay   if (!aij->saved_values) {
143729bbc08cSBarry Smith     SETERRQ(1,"Must call MatStoreValues(A);first");
143849b5e25fSSatish Balay   }
143949b5e25fSSatish Balay 
144049b5e25fSSatish Balay   /* copy values over */
144187828ca2SBarry Smith   ierr = PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(PetscScalar));CHKERRQ(ierr);
144249b5e25fSSatish Balay   PetscFunctionReturn(0);
144349b5e25fSSatish Balay }
144449b5e25fSSatish Balay EXTERN_C_END
144549b5e25fSSatish Balay 
14468549e402SHong Zhang EXTERN_C_BEGIN
14474a2ae208SSatish Balay #undef __FUNCT__
14484a2ae208SSatish Balay #define __FUNCT__ "MatCreate_SeqSBAIJ"
1449c464158bSHong Zhang int MatCreate_SeqSBAIJ(Mat B)
1450c464158bSHong Zhang {
1451c464158bSHong Zhang   Mat_SeqSBAIJ *b;
14520c955e93SHong Zhang   int          ierr,size;
1453c464158bSHong Zhang 
1454c464158bSHong Zhang   PetscFunctionBegin;
1455c464158bSHong Zhang   ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr);
1456c464158bSHong Zhang   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1");
1457273d9f13SBarry Smith   B->m = B->M = PetscMax(B->m,B->M);
1458273d9f13SBarry Smith   B->n = B->N = PetscMax(B->n,B->N);
1459c464158bSHong Zhang 
1460b0a32e0cSBarry Smith   ierr    = PetscNew(Mat_SeqSBAIJ,&b);CHKERRQ(ierr);
1461b0a32e0cSBarry Smith   B->data = (void*)b;
1462c464158bSHong Zhang   ierr    = PetscMemzero(b,sizeof(Mat_SeqSBAIJ));CHKERRQ(ierr);
1463c464158bSHong Zhang   ierr    = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
1464c464158bSHong Zhang   B->ops->destroy     = MatDestroy_SeqSBAIJ;
1465c464158bSHong Zhang   B->ops->view        = MatView_SeqSBAIJ;
1466c464158bSHong Zhang   B->factor           = 0;
1467c464158bSHong Zhang   B->lupivotthreshold = 1.0;
1468c464158bSHong Zhang   B->mapping          = 0;
1469c464158bSHong Zhang   b->row              = 0;
1470c464158bSHong Zhang   b->icol             = 0;
1471c464158bSHong Zhang   b->reallocs         = 0;
1472c464158bSHong Zhang   b->saved_values     = 0;
1473c464158bSHong Zhang 
14748a124369SBarry Smith   ierr = PetscMapCreateMPI(B->comm,B->m,B->M,&B->rmap);CHKERRQ(ierr);
14758a124369SBarry Smith   ierr = PetscMapCreateMPI(B->comm,B->n,B->N,&B->cmap);CHKERRQ(ierr);
1476c464158bSHong Zhang 
1477c464158bSHong Zhang   b->sorted           = PETSC_FALSE;
1478c464158bSHong Zhang   b->roworiented      = PETSC_TRUE;
1479c464158bSHong Zhang   b->nonew            = 0;
1480c464158bSHong Zhang   b->diag             = 0;
1481c464158bSHong Zhang   b->solve_work       = 0;
1482c464158bSHong Zhang   b->mult_work        = 0;
14832a1b7f2aSHong Zhang   B->spptr            = 0;
1484c464158bSHong Zhang   b->keepzeroedrows   = PETSC_FALSE;
1485c464158bSHong Zhang 
1486c464158bSHong Zhang   b->inew             = 0;
1487c464158bSHong Zhang   b->jnew             = 0;
1488c464158bSHong Zhang   b->anew             = 0;
1489c464158bSHong Zhang   b->a2anew           = 0;
1490c464158bSHong Zhang   b->permute          = PETSC_FALSE;
1491c464158bSHong Zhang 
1492c464158bSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C",
1493c464158bSHong Zhang                                      "MatStoreValues_SeqSBAIJ",
1494c464158bSHong Zhang                                      (void*)MatStoreValues_SeqSBAIJ);CHKERRQ(ierr);
1495c464158bSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C",
1496c464158bSHong Zhang                                      "MatRetrieveValues_SeqSBAIJ",
1497c464158bSHong Zhang                                      (void*)MatRetrieveValues_SeqSBAIJ);CHKERRQ(ierr);
1498c464158bSHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqSBAIJSetColumnIndices_C",
1499c464158bSHong Zhang                                      "MatSeqSBAIJSetColumnIndices_SeqSBAIJ",
1500c464158bSHong Zhang                                      (void*)MatSeqSBAIJSetColumnIndices_SeqSBAIJ);CHKERRQ(ierr);
1501c464158bSHong Zhang   PetscFunctionReturn(0);
1502c464158bSHong Zhang }
15038549e402SHong Zhang EXTERN_C_END
1504c464158bSHong Zhang 
15054a2ae208SSatish Balay #undef __FUNCT__
15064a2ae208SSatish Balay #define __FUNCT__ "MatSeqSBAIJSetPreallocation"
150749b5e25fSSatish Balay /*@C
1508273d9f13SBarry Smith    MatSeqSBAIJSetPreallocation - Creates a sparse symmetric matrix in block AIJ (block
150949b5e25fSSatish Balay    compressed row) format.  For good matrix assembly performance the
151049b5e25fSSatish Balay    user should preallocate the matrix storage by setting the parameter nz
151149b5e25fSSatish Balay    (or the array nnz).  By setting these parameters accurately, performance
151249b5e25fSSatish Balay    during matrix assembly can be increased by more than a factor of 50.
151349b5e25fSSatish Balay 
1514c464158bSHong Zhang    Collective on Mat
151549b5e25fSSatish Balay 
151649b5e25fSSatish Balay    Input Parameters:
1517c464158bSHong Zhang +  A - the symmetric matrix
151849b5e25fSSatish Balay .  bs - size of block
151949b5e25fSSatish Balay .  nz - number of block nonzeros per block row (same for all rows)
1520744e8345SSatish Balay -  nnz - array containing the number of block nonzeros in the upper triangular plus
1521744e8345SSatish Balay          diagonal portion of each block (possibly different for each block row) or PETSC_NULL
152249b5e25fSSatish Balay 
152349b5e25fSSatish Balay    Options Database Keys:
152449b5e25fSSatish Balay .   -mat_no_unroll - uses code that does not unroll the loops in the
152549b5e25fSSatish Balay                      block calculations (much slower)
152649b5e25fSSatish Balay .    -mat_block_size - size of the blocks to use
152749b5e25fSSatish Balay 
152849b5e25fSSatish Balay    Level: intermediate
152949b5e25fSSatish Balay 
153049b5e25fSSatish Balay    Notes:
153149b5e25fSSatish Balay    Specify the preallocated storage with either nz or nnz (not both).
153249b5e25fSSatish Balay    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
153349b5e25fSSatish Balay    allocation.  For additional details, see the users manual chapter on
153449b5e25fSSatish Balay    matrices.
153549b5e25fSSatish Balay 
1536a209d233SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ()
153749b5e25fSSatish Balay @*/
1538273d9f13SBarry Smith int MatSeqSBAIJSetPreallocation(Mat B,int bs,int nz,int *nnz)
153949b5e25fSSatish Balay {
1540c464158bSHong Zhang   Mat_SeqSBAIJ *b = (Mat_SeqSBAIJ*)B->data;
15410c955e93SHong Zhang   int          i,len,ierr,mbs,bs2;
154249b5e25fSSatish Balay   PetscTruth   flg;
15434afc71dfSHong Zhang   int          s_nz;
154449b5e25fSSatish Balay 
154549b5e25fSSatish Balay   PetscFunctionBegin;
1546273d9f13SBarry Smith   B->preallocated = PETSC_TRUE;
1547e82a3eeeSBarry Smith   ierr = PetscOptionsGetInt(B->prefix,"-mat_block_size",&bs,PETSC_NULL);CHKERRQ(ierr);
1548c464158bSHong Zhang   mbs  = B->m/bs;
154949b5e25fSSatish Balay   bs2  = bs*bs;
155049b5e25fSSatish Balay 
1551c464158bSHong Zhang   if (mbs*bs != B->m) {
155229bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_SIZ,"Number rows, cols must be divisible by blocksize");
155349b5e25fSSatish Balay   }
155449b5e25fSSatish Balay 
1555435da068SBarry Smith   if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 3;
1556435da068SBarry Smith   if (nz < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %d",nz);
155749b5e25fSSatish Balay   if (nnz) {
155849b5e25fSSatish Balay     for (i=0; i<mbs; i++) {
155929bbc08cSBarry Smith       if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %d value %d",i,nnz[i]);
156080fe4e49SBarry 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);
156149b5e25fSSatish Balay     }
156249b5e25fSSatish Balay   }
156349b5e25fSSatish Balay 
1564e82a3eeeSBarry Smith   ierr    = PetscOptionsHasName(B->prefix,"-mat_no_unroll",&flg);CHKERRQ(ierr);
156549b5e25fSSatish Balay   if (!flg) {
156649b5e25fSSatish Balay     switch (bs) {
156749b5e25fSSatish Balay     case 1:
15684ccecd49SHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1;
156949b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_1;
157049b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_1;
157149b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_1;
157249b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_1;
157349b5e25fSSatish Balay       break;
157449b5e25fSSatish Balay     case 2:
15754ccecd49SHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2;
157649b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_2;
157749b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_2;
157849b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_2;
157949b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_2;
158049b5e25fSSatish Balay       break;
158149b5e25fSSatish Balay     case 3:
15825f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3;
158349b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_3;
158449b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_3;
158549b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_3;
158649b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_3;
158749b5e25fSSatish Balay       break;
158849b5e25fSSatish Balay     case 4:
15895f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4;
159049b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_4;
159149b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_4;
159249b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_4;
159349b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_4;
159449b5e25fSSatish Balay       break;
159549b5e25fSSatish Balay     case 5:
15965f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5;
159749b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_5;
159849b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_5;
159949b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_5;
160049b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_5;
160149b5e25fSSatish Balay       break;
160249b5e25fSSatish Balay     case 6:
16035f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6;
160449b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_6;
160549b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_6;
160649b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_6;
160749b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_6;
160849b5e25fSSatish Balay       break;
160949b5e25fSSatish Balay     case 7:
1610de53e5efSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7;
161149b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_7;
161249b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_7;
1613de53e5efSHong Zhang       B->ops->mult            = MatMult_SeqSBAIJ_7;
161449b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_7;
161549b5e25fSSatish Balay       break;
161649b5e25fSSatish Balay     }
161749b5e25fSSatish Balay   }
161849b5e25fSSatish Balay 
161949b5e25fSSatish Balay   b->mbs = mbs;
16204afc71dfSHong Zhang   b->nbs = mbs;
1621b0a32e0cSBarry Smith   ierr   = PetscMalloc((mbs+1)*sizeof(int),&b->imax);CHKERRQ(ierr);
162249b5e25fSSatish Balay   if (!nnz) {
1623435da068SBarry Smith     if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5;
162449b5e25fSSatish Balay     else if (nz <= 0)        nz = 1;
162549b5e25fSSatish Balay     for (i=0; i<mbs; i++) {
16268cef66ccSHong Zhang       b->imax[i] = nz;
162749b5e25fSSatish Balay     }
1628153ea458SHong Zhang     nz = nz*mbs; /* total nz */
162949b5e25fSSatish Balay   } else {
163049b5e25fSSatish Balay     nz = 0;
16318cef66ccSHong Zhang     for (i=0; i<mbs; i++) {b->imax[i] = nnz[i]; nz += nnz[i];}
163249b5e25fSSatish Balay   }
16338cef66ccSHong Zhang   /* s_nz=(nz+mbs)/2; */ /* total diagonal and superdiagonal nonzero blocks */
16348cef66ccSHong Zhang   s_nz = nz;
163549b5e25fSSatish Balay 
163649b5e25fSSatish Balay   /* allocate the matrix space */
1637c464158bSHong Zhang   len  = s_nz*sizeof(int) + s_nz*bs2*sizeof(MatScalar) + (B->m+1)*sizeof(int);
163882502324SSatish Balay   ierr = PetscMalloc(len,&b->a);CHKERRQ(ierr);
163949b5e25fSSatish Balay   ierr = PetscMemzero(b->a,s_nz*bs2*sizeof(MatScalar));CHKERRQ(ierr);
164049b5e25fSSatish Balay   b->j = (int*)(b->a + s_nz*bs2);
164149b5e25fSSatish Balay   ierr = PetscMemzero(b->j,s_nz*sizeof(int));CHKERRQ(ierr);
164249b5e25fSSatish Balay   b->i = b->j + s_nz;
164349b5e25fSSatish Balay   b->singlemalloc = PETSC_TRUE;
164449b5e25fSSatish Balay 
164549b5e25fSSatish Balay   /* pointer to beginning of each row */
164649b5e25fSSatish Balay   b->i[0] = 0;
164749b5e25fSSatish Balay   for (i=1; i<mbs+1; i++) {
164849b5e25fSSatish Balay     b->i[i] = b->i[i-1] + b->imax[i-1];
164949b5e25fSSatish Balay   }
165049b5e25fSSatish Balay 
165149b5e25fSSatish Balay   /* b->ilen will count nonzeros in each block row so far. */
16525033bc48SSatish Balay   ierr = PetscMalloc((mbs+1)*sizeof(int),&b->ilen);CHKERRQ(ierr);
1653b0a32e0cSBarry Smith   PetscLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqSBAIJ));
165449b5e25fSSatish Balay   for (i=0; i<mbs; i++) { b->ilen[i] = 0;}
165549b5e25fSSatish Balay 
165649b5e25fSSatish Balay   b->bs               = bs;
165749b5e25fSSatish Balay   b->bs2              = bs2;
165849b5e25fSSatish Balay   b->s_nz             = 0;
165949b5e25fSSatish Balay   b->s_maxnz          = s_nz*bs2;
1660153ea458SHong Zhang 
166116cdd363SHong Zhang   b->inew             = 0;
166216cdd363SHong Zhang   b->jnew             = 0;
166316cdd363SHong Zhang   b->anew             = 0;
166416cdd363SHong Zhang   b->a2anew           = 0;
16651a3463dfSHong Zhang   b->permute          = PETSC_FALSE;
1666c464158bSHong Zhang   PetscFunctionReturn(0);
1667c464158bSHong Zhang }
1668153ea458SHong Zhang 
166949b5e25fSSatish Balay 
16704a2ae208SSatish Balay #undef __FUNCT__
16714a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqSBAIJ"
1672c464158bSHong Zhang /*@C
1673c464158bSHong Zhang    MatCreateSeqSBAIJ - Creates a sparse symmetric matrix in block AIJ (block
1674c464158bSHong Zhang    compressed row) format.  For good matrix assembly performance the
1675c464158bSHong Zhang    user should preallocate the matrix storage by setting the parameter nz
1676c464158bSHong Zhang    (or the array nnz).  By setting these parameters accurately, performance
1677c464158bSHong Zhang    during matrix assembly can be increased by more than a factor of 50.
167849b5e25fSSatish Balay 
1679c464158bSHong Zhang    Collective on MPI_Comm
1680c464158bSHong Zhang 
1681c464158bSHong Zhang    Input Parameters:
1682c464158bSHong Zhang +  comm - MPI communicator, set to PETSC_COMM_SELF
1683c464158bSHong Zhang .  bs - size of block
1684c464158bSHong Zhang .  m - number of rows, or number of columns
1685c464158bSHong Zhang .  nz - number of block nonzeros per block row (same for all rows)
1686744e8345SSatish Balay -  nnz - array containing the number of block nonzeros in the upper triangular plus
1687744e8345SSatish Balay          diagonal portion of each block (possibly different for each block row) or PETSC_NULL
1688c464158bSHong Zhang 
1689c464158bSHong Zhang    Output Parameter:
1690c464158bSHong Zhang .  A - the symmetric matrix
1691c464158bSHong Zhang 
1692c464158bSHong Zhang    Options Database Keys:
1693c464158bSHong Zhang .   -mat_no_unroll - uses code that does not unroll the loops in the
1694c464158bSHong Zhang                      block calculations (much slower)
1695c464158bSHong Zhang .    -mat_block_size - size of the blocks to use
1696c464158bSHong Zhang 
1697c464158bSHong Zhang    Level: intermediate
1698c464158bSHong Zhang 
1699c464158bSHong Zhang    Notes:
1700c464158bSHong Zhang 
1701c464158bSHong Zhang    Specify the preallocated storage with either nz or nnz (not both).
1702c464158bSHong Zhang    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
1703c464158bSHong Zhang    allocation.  For additional details, see the users manual chapter on
1704c464158bSHong Zhang    matrices.
1705c464158bSHong Zhang 
1706c464158bSHong Zhang .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ()
1707c464158bSHong Zhang @*/
1708c464158bSHong Zhang int MatCreateSeqSBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,int *nnz,Mat *A)
1709c464158bSHong Zhang {
1710c464158bSHong Zhang   int ierr;
1711c464158bSHong Zhang 
1712c464158bSHong Zhang   PetscFunctionBegin;
1713c464158bSHong Zhang   ierr = MatCreate(comm,m,n,m,n,A);CHKERRQ(ierr);
1714c464158bSHong Zhang   ierr = MatSetType(*A,MATSEQSBAIJ);CHKERRQ(ierr);
1715273d9f13SBarry Smith   ierr = MatSeqSBAIJSetPreallocation(*A,bs,nz,nnz);CHKERRQ(ierr);
171649b5e25fSSatish Balay   PetscFunctionReturn(0);
171749b5e25fSSatish Balay }
171849b5e25fSSatish Balay 
17194a2ae208SSatish Balay #undef __FUNCT__
17204a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqSBAIJ"
172149b5e25fSSatish Balay int MatDuplicate_SeqSBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B)
172249b5e25fSSatish Balay {
172349b5e25fSSatish Balay   Mat         C;
172449b5e25fSSatish Balay   Mat_SeqSBAIJ *c,*a = (Mat_SeqSBAIJ*)A->data;
172549b5e25fSSatish Balay   int         i,len,mbs = a->mbs,nz = a->s_nz,bs2 =a->bs2,ierr;
172649b5e25fSSatish Balay 
172749b5e25fSSatish Balay   PetscFunctionBegin;
172829bbc08cSBarry Smith   if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,"Corrupt matrix");
172949b5e25fSSatish Balay 
173049b5e25fSSatish Balay   *B = 0;
1731692f9cbeSHong Zhang   ierr = MatCreate(A->comm,A->m,A->n,A->m,A->n,&C);CHKERRQ(ierr);
1732692f9cbeSHong Zhang   ierr = MatSetType(C,MATSEQSBAIJ);CHKERRQ(ierr);
1733692f9cbeSHong Zhang   c    = (Mat_SeqSBAIJ*)C->data;
1734692f9cbeSHong Zhang 
173549b5e25fSSatish Balay   ierr              = PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
1736273d9f13SBarry Smith   C->preallocated   = PETSC_TRUE;
173749b5e25fSSatish Balay   C->factor         = A->factor;
173849b5e25fSSatish Balay   c->row            = 0;
173949b5e25fSSatish Balay   c->icol           = 0;
174049b5e25fSSatish Balay   c->saved_values   = 0;
174149b5e25fSSatish Balay   c->keepzeroedrows = a->keepzeroedrows;
174249b5e25fSSatish Balay   C->assembled      = PETSC_TRUE;
174349b5e25fSSatish Balay 
174449b5e25fSSatish Balay   c->bs         = a->bs;
174549b5e25fSSatish Balay   c->bs2        = a->bs2;
174649b5e25fSSatish Balay   c->mbs        = a->mbs;
174749b5e25fSSatish Balay   c->nbs        = a->nbs;
174849b5e25fSSatish Balay 
1749b0a32e0cSBarry Smith   ierr = PetscMalloc((mbs+1)*sizeof(int),&c->imax);CHKERRQ(ierr);
1750b0a32e0cSBarry Smith   ierr = PetscMalloc((mbs+1)*sizeof(int),&c->ilen);CHKERRQ(ierr);
175149b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
175249b5e25fSSatish Balay     c->imax[i] = a->imax[i];
175349b5e25fSSatish Balay     c->ilen[i] = a->ilen[i];
175449b5e25fSSatish Balay   }
175549b5e25fSSatish Balay 
175649b5e25fSSatish Balay   /* allocate the matrix space */
175749b5e25fSSatish Balay   c->singlemalloc = PETSC_TRUE;
175849b5e25fSSatish Balay   len  = (mbs+1)*sizeof(int) + nz*(bs2*sizeof(MatScalar) + sizeof(int));
175982502324SSatish Balay   ierr = PetscMalloc(len,&c->a);CHKERRQ(ierr);
176049b5e25fSSatish Balay   c->j = (int*)(c->a + nz*bs2);
176149b5e25fSSatish Balay   c->i = c->j + nz;
176249b5e25fSSatish Balay   ierr = PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int));CHKERRQ(ierr);
176349b5e25fSSatish Balay   if (mbs > 0) {
176449b5e25fSSatish Balay     ierr = PetscMemcpy(c->j,a->j,nz*sizeof(int));CHKERRQ(ierr);
176549b5e25fSSatish Balay     if (cpvalues == MAT_COPY_VALUES) {
176649b5e25fSSatish Balay       ierr = PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr);
176749b5e25fSSatish Balay     } else {
176849b5e25fSSatish Balay       ierr = PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr);
176949b5e25fSSatish Balay     }
177049b5e25fSSatish Balay   }
177149b5e25fSSatish Balay 
1772b0a32e0cSBarry Smith   PetscLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqSBAIJ));
177349b5e25fSSatish Balay   c->sorted      = a->sorted;
177449b5e25fSSatish Balay   c->roworiented = a->roworiented;
177549b5e25fSSatish Balay   c->nonew       = a->nonew;
177649b5e25fSSatish Balay 
177749b5e25fSSatish Balay   if (a->diag) {
1778b0a32e0cSBarry Smith     ierr = PetscMalloc((mbs+1)*sizeof(int),&c->diag);CHKERRQ(ierr);
1779b0a32e0cSBarry Smith     PetscLogObjectMemory(C,(mbs+1)*sizeof(int));
178049b5e25fSSatish Balay     for (i=0; i<mbs; i++) {
178149b5e25fSSatish Balay       c->diag[i] = a->diag[i];
178249b5e25fSSatish Balay     }
178349b5e25fSSatish Balay   } else c->diag        = 0;
178449b5e25fSSatish Balay   c->s_nz               = a->s_nz;
178549b5e25fSSatish Balay   c->s_maxnz            = a->s_maxnz;
178649b5e25fSSatish Balay   c->solve_work         = 0;
17872a1b7f2aSHong Zhang   C->spptr              = 0;      /* Dangerous -I'm throwing away a->spptr */
178849b5e25fSSatish Balay   c->mult_work          = 0;
178949b5e25fSSatish Balay   *B = C;
1790b0a32e0cSBarry Smith   ierr = PetscFListDuplicate(A->qlist,&C->qlist);CHKERRQ(ierr);
179149b5e25fSSatish Balay   PetscFunctionReturn(0);
179249b5e25fSSatish Balay }
179349b5e25fSSatish Balay 
17948549e402SHong Zhang EXTERN_C_BEGIN
17954a2ae208SSatish Balay #undef __FUNCT__
17964a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqSBAIJ"
1797b0a32e0cSBarry Smith int MatLoad_SeqSBAIJ(PetscViewer viewer,MatType type,Mat *A)
179849b5e25fSSatish Balay {
179949b5e25fSSatish Balay   Mat_SeqSBAIJ *a;
180049b5e25fSSatish Balay   Mat          B;
180149b5e25fSSatish Balay   int          i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1;
18023f7326a9SHong Zhang   int          *mask,mbs,*jj,j,rowcount,nzcount,k,*s_browlengths,maskcount;
180349b5e25fSSatish Balay   int          kmax,jcount,block,idx,point,nzcountb,extra_rows;
180449b5e25fSSatish Balay   int          *masked,nmask,tmp,bs2,ishift;
180587828ca2SBarry Smith   PetscScalar  *aa;
180649b5e25fSSatish Balay   MPI_Comm     comm = ((PetscObject)viewer)->comm;
180749b5e25fSSatish Balay 
180849b5e25fSSatish Balay   PetscFunctionBegin;
1809b0a32e0cSBarry Smith   ierr = PetscOptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,PETSC_NULL);CHKERRQ(ierr);
181049b5e25fSSatish Balay   bs2  = bs*bs;
181149b5e25fSSatish Balay 
181249b5e25fSSatish Balay   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
181329bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor");
1814b0a32e0cSBarry Smith   ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
181549b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr);
1816552e946dSBarry Smith   if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not Mat object");
181749b5e25fSSatish Balay   M = header[1]; N = header[2]; nz = header[3];
181849b5e25fSSatish Balay 
181949b5e25fSSatish Balay   if (header[3] < 0) {
182029bbc08cSBarry Smith     SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as SeqSBAIJ");
182149b5e25fSSatish Balay   }
182249b5e25fSSatish Balay 
182329bbc08cSBarry Smith   if (M != N) SETERRQ(PETSC_ERR_SUP,"Can only do square matrices");
182449b5e25fSSatish Balay 
182549b5e25fSSatish Balay   /*
182649b5e25fSSatish Balay      This code adds extra rows to make sure the number of rows is
182749b5e25fSSatish Balay     divisible by the blocksize
182849b5e25fSSatish Balay   */
182949b5e25fSSatish Balay   mbs        = M/bs;
183049b5e25fSSatish Balay   extra_rows = bs - M + bs*(mbs);
183149b5e25fSSatish Balay   if (extra_rows == bs) extra_rows = 0;
183249b5e25fSSatish Balay   else                  mbs++;
183349b5e25fSSatish Balay   if (extra_rows) {
1834b0a32e0cSBarry Smith     PetscLogInfo(0,"MatLoad_SeqSBAIJ:Padding loaded matrix to match blocksize\n");
183549b5e25fSSatish Balay   }
183649b5e25fSSatish Balay 
183749b5e25fSSatish Balay   /* read in row lengths */
1838b0a32e0cSBarry Smith   ierr = PetscMalloc((M+extra_rows)*sizeof(int),&rowlengths);CHKERRQ(ierr);
183949b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
184049b5e25fSSatish Balay   for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1;
184149b5e25fSSatish Balay 
184249b5e25fSSatish Balay   /* read in column indices */
1843b0a32e0cSBarry Smith   ierr = PetscMalloc((nz+extra_rows)*sizeof(int),&jj);CHKERRQ(ierr);
184449b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT);CHKERRQ(ierr);
184549b5e25fSSatish Balay   for (i=0; i<extra_rows; i++) jj[nz+i] = M+i;
184649b5e25fSSatish Balay 
184749b5e25fSSatish Balay   /* loop over row lengths determining block row lengths */
184882502324SSatish Balay   ierr     = PetscMalloc(mbs*sizeof(int),&s_browlengths);CHKERRQ(ierr);
1849a91a614fSBarry Smith   ierr     = PetscMemzero(s_browlengths,mbs*sizeof(int));CHKERRQ(ierr);
185082502324SSatish Balay   ierr     = PetscMalloc(2*mbs*sizeof(int),&mask);CHKERRQ(ierr);
185149b5e25fSSatish Balay   ierr     = PetscMemzero(mask,mbs*sizeof(int));CHKERRQ(ierr);
185249b5e25fSSatish Balay   masked   = mask + mbs;
185349b5e25fSSatish Balay   rowcount = 0; nzcount = 0;
185449b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
185549b5e25fSSatish Balay     nmask = 0;
185649b5e25fSSatish Balay     for (j=0; j<bs; j++) {
185749b5e25fSSatish Balay       kmax = rowlengths[rowcount];
185849b5e25fSSatish Balay       for (k=0; k<kmax; k++) {
18592d703238SHong Zhang         tmp = jj[nzcount++]/bs;   /* block col. index */
186003630b6eSHong Zhang         if (!mask[tmp] && tmp >= i) {masked[nmask++] = tmp; mask[tmp] = 1;}
186149b5e25fSSatish Balay       }
186249b5e25fSSatish Balay       rowcount++;
186349b5e25fSSatish Balay     }
1864574b2666SHong Zhang     s_browlengths[i] += nmask;
1865574b2666SHong Zhang 
186649b5e25fSSatish Balay     /* zero out the mask elements we set */
186749b5e25fSSatish Balay     for (j=0; j<nmask; j++) mask[masked[j]] = 0;
186849b5e25fSSatish Balay   }
186949b5e25fSSatish Balay 
187049b5e25fSSatish Balay   /* create our matrix */
18717fe2be48SHong Zhang   ierr = MatCreateSeqSBAIJ(comm,bs,M+extra_rows,N+extra_rows,0,s_browlengths,A);CHKERRQ(ierr);
187249b5e25fSSatish Balay   B = *A;
187349b5e25fSSatish Balay   a = (Mat_SeqSBAIJ*)B->data;
187449b5e25fSSatish Balay 
187549b5e25fSSatish Balay   /* set matrix "i" values */
187649b5e25fSSatish Balay   a->i[0] = 0;
187749b5e25fSSatish Balay   for (i=1; i<= mbs; i++) {
1878574b2666SHong Zhang     a->i[i]      = a->i[i-1] + s_browlengths[i-1];
1879574b2666SHong Zhang     a->ilen[i-1] = s_browlengths[i-1];
188049b5e25fSSatish Balay   }
18817fe2be48SHong Zhang   a->s_nz = a->i[mbs];
188249b5e25fSSatish Balay 
188349b5e25fSSatish Balay   /* read in nonzero values */
188487828ca2SBarry Smith   ierr = PetscMalloc((nz+extra_rows)*sizeof(PetscScalar),&aa);CHKERRQ(ierr);
188549b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);CHKERRQ(ierr);
188649b5e25fSSatish Balay   for (i=0; i<extra_rows; i++) aa[nz+i] = 1.0;
188749b5e25fSSatish Balay 
188849b5e25fSSatish Balay   /* set "a" and "j" values into matrix */
188949b5e25fSSatish Balay   nzcount = 0; jcount = 0;
189049b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
189149b5e25fSSatish Balay     nzcountb = nzcount;
189249b5e25fSSatish Balay     nmask    = 0;
189349b5e25fSSatish Balay     for (j=0; j<bs; j++) {
189449b5e25fSSatish Balay       kmax = rowlengths[i*bs+j];
189549b5e25fSSatish Balay       for (k=0; k<kmax; k++) {
18962d703238SHong Zhang         tmp = jj[nzcount++]/bs; /* block col. index */
189703630b6eSHong Zhang         if (!mask[tmp] && tmp >= i) { masked[nmask++] = tmp; mask[tmp] = 1;}
18982d703238SHong Zhang       }
18992d703238SHong Zhang     }
19002d703238SHong Zhang     /* sort the masked values */
19012d703238SHong Zhang     ierr = PetscSortInt(nmask,masked);CHKERRQ(ierr);
19022d703238SHong Zhang 
19032d703238SHong Zhang     /* set "j" values into matrix */
19042d703238SHong Zhang     maskcount = 1;
19052d703238SHong Zhang     for (j=0; j<nmask; j++) {
190649b5e25fSSatish Balay       a->j[jcount++]  = masked[j];
190749b5e25fSSatish Balay       mask[masked[j]] = maskcount++;
190849b5e25fSSatish Balay     }
1909574b2666SHong Zhang 
191049b5e25fSSatish Balay     /* set "a" values into matrix */
191149b5e25fSSatish Balay     ishift = bs2*a->i[i];
191249b5e25fSSatish Balay     for (j=0; j<bs; j++) {
191349b5e25fSSatish Balay       kmax = rowlengths[i*bs+j];
191449b5e25fSSatish Balay       for (k=0; k<kmax; k++) {
1915574b2666SHong Zhang         tmp       = jj[nzcountb]/bs ; /* block col. index */
1916574b2666SHong Zhang         if (tmp >= i){
191749b5e25fSSatish Balay           block     = mask[tmp] - 1;
191849b5e25fSSatish Balay           point     = jj[nzcountb] - bs*tmp;
191949b5e25fSSatish Balay           idx       = ishift + bs2*block + j + bs*point;
1920574b2666SHong Zhang           a->a[idx] = aa[nzcountb];
1921574b2666SHong Zhang         }
1922574b2666SHong Zhang         nzcountb++;
192349b5e25fSSatish Balay       }
192449b5e25fSSatish Balay     }
192549b5e25fSSatish Balay     /* zero out the mask elements we set */
192649b5e25fSSatish Balay     for (j=0; j<nmask; j++) mask[masked[j]] = 0;
192749b5e25fSSatish Balay   }
192829bbc08cSBarry Smith   if (jcount != a->s_nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Bad binary matrix");
192949b5e25fSSatish Balay 
193049b5e25fSSatish Balay   ierr = PetscFree(rowlengths);CHKERRQ(ierr);
1931574b2666SHong Zhang   ierr = PetscFree(s_browlengths);CHKERRQ(ierr);
193249b5e25fSSatish Balay   ierr = PetscFree(aa);CHKERRQ(ierr);
193349b5e25fSSatish Balay   ierr = PetscFree(jj);CHKERRQ(ierr);
193449b5e25fSSatish Balay   ierr = PetscFree(mask);CHKERRQ(ierr);
193549b5e25fSSatish Balay 
193649b5e25fSSatish Balay   B->assembled = PETSC_TRUE;
193749b5e25fSSatish Balay   ierr = MatView_Private(B);CHKERRQ(ierr);
193849b5e25fSSatish Balay   PetscFunctionReturn(0);
193949b5e25fSSatish Balay }
19408549e402SHong Zhang EXTERN_C_END
1941574b2666SHong Zhang 
1942d06b337dSHong Zhang #undef __FUNCT__
1943d06b337dSHong Zhang #define __FUNCT__ "MatRelax_SeqSBAIJ"
1944c14dc6b6SHong Zhang int MatRelax_SeqSBAIJ(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,int its,int lits,Vec xx)
1945d06b337dSHong Zhang {
1946d06b337dSHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
1947d06b337dSHong Zhang   MatScalar    *aa=a->a,*v,*v1;
1948d06b337dSHong Zhang   PetscScalar  *x,*b,*t,sum,d;
1949d06b337dSHong Zhang   int          m=a->mbs,bs=a->bs,*ai=a->i,*aj=a->j,ierr;
1950d06b337dSHong Zhang   int          nz,nz1,*vj,*vj1,i;
1951d06b337dSHong Zhang 
1952d06b337dSHong Zhang   PetscFunctionBegin;
1953b965ef7fSBarry Smith   its = its*lits;
195491723122SBarry Smith   if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %d and local its %d both positive",its,lits);
1955d06b337dSHong Zhang 
1956d06b337dSHong Zhang   if (bs > 1)
1957d06b337dSHong Zhang     SETERRQ(PETSC_ERR_SUP,"SSOR for block size > 1 is not yet implemented");
1958d06b337dSHong Zhang 
1959d06b337dSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1960d06b337dSHong Zhang   if (xx != bb) {
1961d06b337dSHong Zhang     ierr = VecGetArray(bb,&b);CHKERRQ(ierr);
1962d06b337dSHong Zhang   } else {
1963d06b337dSHong Zhang     b = x;
1964d06b337dSHong Zhang   }
1965d06b337dSHong Zhang 
1966d06b337dSHong Zhang   ierr = PetscMalloc(m*sizeof(PetscScalar),&t);CHKERRQ(ierr);
1967d06b337dSHong Zhang 
1968d06b337dSHong Zhang   if (flag & SOR_ZERO_INITIAL_GUESS) {
19692798e883SHong Zhang     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
1970d06b337dSHong Zhang       for (i=0; i<m; i++)
1971d06b337dSHong Zhang         t[i] = b[i];
1972d06b337dSHong Zhang 
1973d06b337dSHong Zhang       for (i=0; i<m; i++){
197444706875SHong Zhang         d  = *(aa + ai[i]);  /* diag[i] */
1975d06b337dSHong Zhang         v  = aa + ai[i] + 1;
1976d06b337dSHong Zhang         vj = aj + ai[i] + 1;
1977d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1;
1978d06b337dSHong Zhang         x[i] = omega*t[i]/d;
1979d06b337dSHong Zhang         while (nz--) t[*vj++] -= x[i]*(*v++); /* update rhs */
198044706875SHong Zhang         PetscLogFlops(2*nz-1);
1981d06b337dSHong Zhang       }
1982d06b337dSHong Zhang     }
1983d06b337dSHong Zhang 
19842798e883SHong Zhang     if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
1985d06b337dSHong Zhang       for (i=0; i<m; i++)
1986d06b337dSHong Zhang         t[i] = b[i];
1987d06b337dSHong Zhang 
1988d06b337dSHong Zhang       for (i=0; i<m-1; i++){  /* update rhs */
1989d06b337dSHong Zhang         v  = aa + ai[i] + 1;
1990d06b337dSHong Zhang         vj = aj + ai[i] + 1;
1991d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1;
1992d06b337dSHong Zhang         while (nz--) t[*vj++] -= x[i]*(*v++);
199344706875SHong Zhang         PetscLogFlops(2*nz-1);
1994d06b337dSHong Zhang       }
1995d06b337dSHong Zhang       for (i=m-1; i>=0; i--){
1996d06b337dSHong Zhang         d  = *(aa + ai[i]);
1997d06b337dSHong Zhang         v  = aa + ai[i] + 1;
1998d06b337dSHong Zhang         vj = aj + ai[i] + 1;
1999d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1;
2000d06b337dSHong Zhang         sum = t[i];
2001d06b337dSHong Zhang         while (nz--) sum -= x[*vj++]*(*v++);
200244706875SHong Zhang         PetscLogFlops(2*nz-1);
2003d06b337dSHong Zhang         x[i] =   (1-omega)*x[i] + omega*sum/d;
2004d06b337dSHong Zhang       }
2005d06b337dSHong Zhang     }
2006d06b337dSHong Zhang     its--;
2007d06b337dSHong Zhang   }
2008d06b337dSHong Zhang 
2009d06b337dSHong Zhang   while (its--) {
201044706875SHong Zhang     /*
201144706875SHong Zhang        forward sweep:
201244706875SHong Zhang        for i=0,...,m-1:
201344706875SHong Zhang          sum[i] = (b[i] - U(i,:)x )/d[i];
201444706875SHong Zhang          x[i]   = (1-omega)x[i] + omega*sum[i];
201544706875SHong Zhang          b      = b - x[i]*U^T(i,:);
2016d06b337dSHong Zhang 
201744706875SHong Zhang     */
20182798e883SHong Zhang     if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){
2019d06b337dSHong Zhang       for (i=0; i<m; i++)
2020d06b337dSHong Zhang         t[i] = b[i];
2021d06b337dSHong Zhang 
2022d06b337dSHong Zhang       for (i=0; i<m; i++){
202344706875SHong Zhang         d  = *(aa + ai[i]);  /* diag[i] */
2024d06b337dSHong Zhang         v  = aa + ai[i] + 1; v1=v;
2025d06b337dSHong Zhang         vj = aj + ai[i] + 1; vj1=vj;
2026d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1; nz1=nz;
2027d06b337dSHong Zhang         sum = t[i];
2028d06b337dSHong Zhang         while (nz1--) sum -= (*v1++)*x[*vj1++];
2029d06b337dSHong Zhang         x[i] = (1-omega)*x[i] + omega*sum/d;
2030d06b337dSHong Zhang         while (nz--) t[*vj++] -= x[i]*(*v++);
203144706875SHong Zhang         PetscLogFlops(4*nz-2);
2032d06b337dSHong Zhang       }
2033d06b337dSHong Zhang     }
2034d06b337dSHong Zhang 
20352798e883SHong Zhang   if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){
203644706875SHong Zhang       /*
203744706875SHong Zhang        backward sweep:
203844706875SHong Zhang        b = b - x[i]*U^T(i,:), i=0,...,n-2
203944706875SHong Zhang        for i=m-1,...,0:
204044706875SHong Zhang          sum[i] = (b[i] - U(i,:)x )/d[i];
204144706875SHong Zhang          x[i]   = (1-omega)x[i] + omega*sum[i];
204244706875SHong Zhang       */
2043d06b337dSHong Zhang       for (i=0; i<m; i++)
2044d06b337dSHong Zhang         t[i] = b[i];
2045d06b337dSHong Zhang 
2046d06b337dSHong Zhang       for (i=0; i<m-1; i++){  /* update rhs */
2047d06b337dSHong Zhang         v  = aa + ai[i] + 1;
2048d06b337dSHong Zhang         vj = aj + ai[i] + 1;
2049d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1;
2050d06b337dSHong Zhang         while (nz--) t[*vj++] -= x[i]*(*v++);
205144706875SHong Zhang         PetscLogFlops(2*nz-1);
2052d06b337dSHong Zhang       }
2053d06b337dSHong Zhang       for (i=m-1; i>=0; i--){
2054d06b337dSHong Zhang         d  = *(aa + ai[i]);
2055d06b337dSHong Zhang         v  = aa + ai[i] + 1;
2056d06b337dSHong Zhang         vj = aj + ai[i] + 1;
2057d06b337dSHong Zhang         nz = ai[i+1] - ai[i] - 1;
2058d06b337dSHong Zhang         sum = t[i];
2059d06b337dSHong Zhang         while (nz--) sum -= x[*vj++]*(*v++);
206044706875SHong Zhang         PetscLogFlops(2*nz-1);
2061d06b337dSHong Zhang         x[i] =   (1-omega)*x[i] + omega*sum/d;
2062d06b337dSHong Zhang       }
2063d06b337dSHong Zhang     }
2064d06b337dSHong Zhang   }
2065d06b337dSHong Zhang 
2066d06b337dSHong Zhang   ierr = PetscFree(t); CHKERRQ(ierr);
2067d06b337dSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2068d06b337dSHong Zhang   if (bb != xx) {
2069d06b337dSHong Zhang     ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr);
2070d06b337dSHong Zhang   }
20712798e883SHong Zhang 
2072d06b337dSHong Zhang   PetscFunctionReturn(0);
2073d06b337dSHong Zhang }
2074d06b337dSHong Zhang 
2075d06b337dSHong Zhang 
2076d06b337dSHong Zhang 
2077d06b337dSHong Zhang 
207849b5e25fSSatish Balay 
207949b5e25fSSatish Balay 
2080