xref: /petsc/src/mat/impls/sbaij/seq/sbaij.c (revision 29bbc08cd5461c366aba6645924e8ff42acd1de0)
1*29bbc08cSBarry Smith /*$Id: sbaij.c,v 1.29 2000/09/25 19:53:56 curfman Exp bsmith $*/
249b5e25fSSatish Balay 
349b5e25fSSatish Balay /*
449b5e25fSSatish Balay     Defines the basic matrix operations for the BAIJ (compressed row)
549b5e25fSSatish Balay   matrix storage format.
649b5e25fSSatish Balay */
7ab9f2c04SSatish Balay #include "petscsys.h"                            /*I "petscmat.h" I*/
849b5e25fSSatish Balay #include "src/mat/impls/baij/seq/baij.h"
949b5e25fSSatish Balay #include "src/vec/vecimpl.h"
1049b5e25fSSatish Balay #include "src/inline/spops.h"
11aec82049SSatish Balay #include "src/mat/impls/sbaij/seq/sbaij.h"
1249b5e25fSSatish Balay 
1349b5e25fSSatish Balay #define CHUNKSIZE  10
1449b5e25fSSatish Balay 
1549b5e25fSSatish Balay /*
1649b5e25fSSatish Balay      Checks for missing diagonals
1749b5e25fSSatish Balay */
1849b5e25fSSatish Balay #undef __FUNC__
1949b5e25fSSatish Balay #define __FUNC__ "MatMissingDiagonal_SeqSBAIJ"
2049b5e25fSSatish Balay int MatMissingDiagonal_SeqSBAIJ(Mat A)
2149b5e25fSSatish Balay {
22045c9aa0SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
2349b5e25fSSatish Balay   int         *diag,*jj = a->j,i,ierr;
2449b5e25fSSatish Balay 
2549b5e25fSSatish Balay   PetscFunctionBegin;
26045c9aa0SHong Zhang   ierr = MatMarkDiagonal_SeqSBAIJ(A);CHKERRQ(ierr);
2749b5e25fSSatish Balay   diag = a->diag;
2849b5e25fSSatish Balay   for (i=0; i<a->mbs; i++) {
2949b5e25fSSatish Balay     if (jj[diag[i]] != i) {
30*29bbc08cSBarry Smith       SETERRQ1(1,"Matrix is missing diagonal number %d",i);
3149b5e25fSSatish Balay     }
3249b5e25fSSatish Balay   }
3349b5e25fSSatish Balay   PetscFunctionReturn(0);
3449b5e25fSSatish Balay }
3549b5e25fSSatish Balay 
3649b5e25fSSatish Balay #undef __FUNC__
3749b5e25fSSatish Balay #define __FUNC__ "MatMarkDiagonal_SeqSBAIJ"
3849b5e25fSSatish Balay int MatMarkDiagonal_SeqSBAIJ(Mat A)
3949b5e25fSSatish Balay {
40045c9aa0SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
4149b5e25fSSatish Balay   int          i,j,*diag,m = a->mbs;
4249b5e25fSSatish Balay 
4349b5e25fSSatish Balay   PetscFunctionBegin;
4449b5e25fSSatish Balay   if (a->diag) PetscFunctionReturn(0);
4549b5e25fSSatish Balay 
4649b5e25fSSatish Balay   diag = (int*)PetscMalloc((m+1)*sizeof(int));CHKPTRQ(diag);
4749b5e25fSSatish Balay   PLogObjectMemory(A,(m+1)*sizeof(int));
4849b5e25fSSatish Balay   for (i=0; i<m; i++) {
4949b5e25fSSatish Balay     diag[i] = a->i[i+1];
5049b5e25fSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
5149b5e25fSSatish Balay       if (a->j[j] == i) {
5249b5e25fSSatish Balay         diag[i] = j;
5349b5e25fSSatish Balay         break;
5449b5e25fSSatish Balay       }
5549b5e25fSSatish Balay     }
5649b5e25fSSatish Balay   }
5749b5e25fSSatish Balay   a->diag = diag;
5849b5e25fSSatish Balay   PetscFunctionReturn(0);
5949b5e25fSSatish Balay }
6049b5e25fSSatish Balay 
6149b5e25fSSatish Balay 
6249b5e25fSSatish Balay extern int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int,int**,int**);
6349b5e25fSSatish Balay 
6449b5e25fSSatish Balay #undef __FUNC__
6549b5e25fSSatish Balay #define __FUNC__ "MatGetRowIJ_SeqSBAIJ"
6649b5e25fSSatish Balay static int MatGetRowIJ_SeqSBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja,PetscTruth *done)
6749b5e25fSSatish Balay {
6849b5e25fSSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
6949b5e25fSSatish Balay 
7049b5e25fSSatish Balay   PetscFunctionBegin;
71045c9aa0SHong Zhang   if (ia) {
72*29bbc08cSBarry Smith     SETERRQ(1,"Function not yet written for SBAIJ format, only supports natural ordering");
7349b5e25fSSatish Balay   }
74045c9aa0SHong Zhang   *nn = a->mbs;
7549b5e25fSSatish Balay   PetscFunctionReturn(0);
7649b5e25fSSatish Balay }
7749b5e25fSSatish Balay 
7849b5e25fSSatish Balay #undef __FUNC__
7949b5e25fSSatish Balay #define __FUNC__ "MatRestoreRowIJ_SeqSBAIJ"
80045c9aa0SHong Zhang static int MatRestoreRowIJ_SeqSBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja,PetscTruth *done)
8149b5e25fSSatish Balay {
8249b5e25fSSatish Balay   PetscFunctionBegin;
8349b5e25fSSatish Balay   if (!ia) PetscFunctionReturn(0);
84*29bbc08cSBarry Smith   SETERRQ(1,"Function not yet written for SBAIJ format, only supports natural ordering");
8516cdd363SHong Zhang   /* PetscFunctionReturn(0); */
8649b5e25fSSatish Balay }
8749b5e25fSSatish Balay 
8849b5e25fSSatish Balay #undef __FUNC__
8949b5e25fSSatish Balay #define __FUNC__ "MatGetBlockSize_SeqSBAIJ"
9049b5e25fSSatish Balay int MatGetBlockSize_SeqSBAIJ(Mat mat,int *bs)
9149b5e25fSSatish Balay {
9249b5e25fSSatish Balay   Mat_SeqSBAIJ *sbaij = (Mat_SeqSBAIJ*)mat->data;
9349b5e25fSSatish Balay 
9449b5e25fSSatish Balay   PetscFunctionBegin;
9549b5e25fSSatish Balay   *bs = sbaij->bs;
9649b5e25fSSatish Balay   PetscFunctionReturn(0);
9749b5e25fSSatish Balay }
9849b5e25fSSatish Balay 
9949b5e25fSSatish Balay #undef __FUNC__
10049b5e25fSSatish Balay #define __FUNC__ "MatDestroy_SeqSBAIJ"
10149b5e25fSSatish Balay int MatDestroy_SeqSBAIJ(Mat A)
10249b5e25fSSatish Balay {
10349b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
10449b5e25fSSatish Balay   int         ierr;
10549b5e25fSSatish Balay 
10649b5e25fSSatish Balay   PetscFunctionBegin;
10749b5e25fSSatish Balay   if (A->mapping) {
10849b5e25fSSatish Balay     ierr = ISLocalToGlobalMappingDestroy(A->mapping);CHKERRQ(ierr);
10949b5e25fSSatish Balay   }
11049b5e25fSSatish Balay   if (A->bmapping) {
11149b5e25fSSatish Balay     ierr = ISLocalToGlobalMappingDestroy(A->bmapping);CHKERRQ(ierr);
11249b5e25fSSatish Balay   }
11349b5e25fSSatish Balay   if (A->rmap) {
11449b5e25fSSatish Balay     ierr = MapDestroy(A->rmap);CHKERRQ(ierr);
11549b5e25fSSatish Balay   }
11649b5e25fSSatish Balay   if (A->cmap) {
11749b5e25fSSatish Balay     ierr = MapDestroy(A->cmap);CHKERRQ(ierr);
11849b5e25fSSatish Balay   }
11949b5e25fSSatish Balay #if defined(PETSC_USE_LOG)
12049b5e25fSSatish Balay   PLogObjectState((PetscObject)A,"Rows=%d, s_NZ=%d",a->m,a->s_nz);
12149b5e25fSSatish Balay #endif
12249b5e25fSSatish Balay   ierr = PetscFree(a->a);CHKERRQ(ierr);
12349b5e25fSSatish Balay   if (!a->singlemalloc) {
12449b5e25fSSatish Balay     ierr = PetscFree(a->i);CHKERRQ(ierr);
12563c8bf9fSHong Zhang     ierr = PetscFree(a->j);CHKERRQ(ierr);
12649b5e25fSSatish Balay   }
12749b5e25fSSatish Balay   if (a->row) {
12849b5e25fSSatish Balay     ierr = ISDestroy(a->row);CHKERRQ(ierr);
12949b5e25fSSatish Balay   }
13049b5e25fSSatish Balay   if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);}
13149b5e25fSSatish Balay   if (a->ilen) {ierr = PetscFree(a->ilen);CHKERRQ(ierr);}
13249b5e25fSSatish Balay   if (a->imax) {ierr = PetscFree(a->imax);CHKERRQ(ierr);}
13349b5e25fSSatish Balay   if (a->solve_work) {ierr = PetscFree(a->solve_work);CHKERRQ(ierr);}
13449b5e25fSSatish Balay   if (a->mult_work) {ierr = PetscFree(a->mult_work);CHKERRQ(ierr);}
13549b5e25fSSatish Balay   if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);}
13649b5e25fSSatish Balay   if (a->saved_values) {ierr = PetscFree(a->saved_values);CHKERRQ(ierr);}
1371a3463dfSHong Zhang 
1381a3463dfSHong Zhang   if (a->inew){
1391a3463dfSHong Zhang     ierr = PetscFree(a->inew);CHKERRQ(ierr);
1401a3463dfSHong Zhang     a->inew = 0;
1411a3463dfSHong Zhang   }
14249b5e25fSSatish Balay   ierr = PetscFree(a);CHKERRQ(ierr);
14349b5e25fSSatish Balay   PLogObjectDestroy(A);
14449b5e25fSSatish Balay   PetscHeaderDestroy(A);
14549b5e25fSSatish Balay   PetscFunctionReturn(0);
14649b5e25fSSatish Balay }
14749b5e25fSSatish Balay 
14849b5e25fSSatish Balay #undef __FUNC__
14949b5e25fSSatish Balay #define __FUNC__ "MatSetOption_SeqSBAIJ"
15049b5e25fSSatish Balay int MatSetOption_SeqSBAIJ(Mat A,MatOption op)
15149b5e25fSSatish Balay {
152045c9aa0SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
15349b5e25fSSatish Balay 
15449b5e25fSSatish Balay   PetscFunctionBegin;
15549b5e25fSSatish Balay   if      (op == MAT_ROW_ORIENTED)                 a->roworiented    = PETSC_TRUE;
15649b5e25fSSatish Balay   else if (op == MAT_COLUMN_ORIENTED)              a->roworiented    = PETSC_FALSE;
15749b5e25fSSatish Balay   else if (op == MAT_COLUMNS_SORTED)               a->sorted         = PETSC_TRUE;
15849b5e25fSSatish Balay   else if (op == MAT_COLUMNS_UNSORTED)             a->sorted         = PETSC_FALSE;
15949b5e25fSSatish Balay   else if (op == MAT_KEEP_ZEROED_ROWS)             a->keepzeroedrows = PETSC_TRUE;
16049b5e25fSSatish Balay   else if (op == MAT_NO_NEW_NONZERO_LOCATIONS)     a->nonew          = 1;
16149b5e25fSSatish Balay   else if (op == MAT_NEW_NONZERO_LOCATION_ERR)     a->nonew          = -1;
16249b5e25fSSatish Balay   else if (op == MAT_NEW_NONZERO_ALLOCATION_ERR)   a->nonew          = -2;
16349b5e25fSSatish Balay   else if (op == MAT_YES_NEW_NONZERO_LOCATIONS)    a->nonew          = 0;
16449b5e25fSSatish Balay   else if (op == MAT_ROWS_SORTED ||
16549b5e25fSSatish Balay            op == MAT_ROWS_UNSORTED ||
16649b5e25fSSatish Balay            op == MAT_SYMMETRIC ||
16749b5e25fSSatish Balay            op == MAT_STRUCTURALLY_SYMMETRIC ||
16849b5e25fSSatish Balay            op == MAT_YES_NEW_DIAGONALS ||
16949b5e25fSSatish Balay            op == MAT_IGNORE_OFF_PROC_ENTRIES ||
17049b5e25fSSatish Balay            op == MAT_USE_HASH_TABLE) {
171045c9aa0SHong Zhang     PLogInfo(A,"MatSetOption_SeqSBAIJ:Option ignored\n");
17249b5e25fSSatish Balay   } else if (op == MAT_NO_NEW_DIAGONALS) {
173*29bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS");
17449b5e25fSSatish Balay   } else {
175*29bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"unknown option");
17649b5e25fSSatish Balay   }
17749b5e25fSSatish Balay   PetscFunctionReturn(0);
17849b5e25fSSatish Balay }
17949b5e25fSSatish Balay 
18049b5e25fSSatish Balay 
18149b5e25fSSatish Balay #undef __FUNC__
18249b5e25fSSatish Balay #define __FUNC__ "MatGetSize_SeqSBAIJ"
18349b5e25fSSatish Balay int MatGetSize_SeqSBAIJ(Mat A,int *m,int *n)
18449b5e25fSSatish Balay {
18549b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
18649b5e25fSSatish Balay 
18749b5e25fSSatish Balay   PetscFunctionBegin;
18849b5e25fSSatish Balay   if (m) *m = a->m;
189045c9aa0SHong Zhang   if (n) *n = a->n;
19049b5e25fSSatish Balay   PetscFunctionReturn(0);
19149b5e25fSSatish Balay }
19249b5e25fSSatish Balay 
19349b5e25fSSatish Balay #undef __FUNC__
19449b5e25fSSatish Balay #define __FUNC__ "MatGetOwnershipRange_SeqSBAIJ"
19549b5e25fSSatish Balay int MatGetOwnershipRange_SeqSBAIJ(Mat A,int *m,int *n)
19649b5e25fSSatish Balay {
19749b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
19849b5e25fSSatish Balay 
19949b5e25fSSatish Balay   PetscFunctionBegin;
20049b5e25fSSatish Balay   *m = 0; *n = a->m;
20149b5e25fSSatish Balay   PetscFunctionReturn(0);
20249b5e25fSSatish Balay }
20349b5e25fSSatish Balay 
20449b5e25fSSatish Balay #undef __FUNC__
20549b5e25fSSatish Balay #define __FUNC__ "MatGetRow_SeqSBAIJ"
20649b5e25fSSatish Balay int MatGetRow_SeqSBAIJ(Mat A,int row,int *ncols,int **cols,Scalar **v)
20749b5e25fSSatish Balay {
20849b5e25fSSatish Balay   Mat_SeqSBAIJ  *a = (Mat_SeqSBAIJ*)A->data;
20949b5e25fSSatish Balay   int          itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*cols_i,bs2;
21049b5e25fSSatish Balay   MatScalar    *aa,*aa_i;
21149b5e25fSSatish Balay   Scalar       *v_i;
21249b5e25fSSatish Balay 
21349b5e25fSSatish Balay   PetscFunctionBegin;
21449b5e25fSSatish Balay   bs  = a->bs;
21549b5e25fSSatish Balay   ai  = a->i;
21649b5e25fSSatish Balay   aj  = a->j;
21749b5e25fSSatish Balay   aa  = a->a;
21849b5e25fSSatish Balay   bs2 = a->bs2;
21949b5e25fSSatish Balay 
220*29bbc08cSBarry Smith   if (row < 0 || row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row out of range");
22149b5e25fSSatish Balay 
22249b5e25fSSatish Balay   bn  = row/bs;   /* Block number */
22349b5e25fSSatish Balay   bp  = row % bs; /* Block position */
22449b5e25fSSatish Balay   M   = ai[bn+1] - ai[bn];
22549b5e25fSSatish Balay   *ncols = bs*M;
22649b5e25fSSatish Balay 
22749b5e25fSSatish Balay   if (v) {
22849b5e25fSSatish Balay     *v = 0;
22949b5e25fSSatish Balay     if (*ncols) {
23049b5e25fSSatish Balay       *v = (Scalar*)PetscMalloc((*ncols+row)*sizeof(Scalar));CHKPTRQ(*v);
23149b5e25fSSatish Balay       for (i=0; i<M; i++) { /* for each block in the block row */
23249b5e25fSSatish Balay         v_i  = *v + i*bs;
23349b5e25fSSatish Balay         aa_i = aa + bs2*(ai[bn] + i);
23449b5e25fSSatish Balay         for (j=bp,k=0; j<bs2; j+=bs,k++) {v_i[k] = aa_i[j];}
23549b5e25fSSatish Balay       }
23649b5e25fSSatish Balay     }
23749b5e25fSSatish Balay   }
23849b5e25fSSatish Balay 
23949b5e25fSSatish Balay   if (cols) {
24049b5e25fSSatish Balay     *cols = 0;
24149b5e25fSSatish Balay     if (*ncols) {
24249b5e25fSSatish Balay       *cols = (int*)PetscMalloc((*ncols+row)*sizeof(int));CHKPTRQ(*cols);
24349b5e25fSSatish Balay       for (i=0; i<M; i++) { /* for each block in the block row */
24449b5e25fSSatish Balay         cols_i = *cols + i*bs;
24549b5e25fSSatish Balay         itmp  = bs*aj[ai[bn] + i];
24649b5e25fSSatish Balay         for (j=0; j<bs; j++) {cols_i[j] = itmp++;}
24749b5e25fSSatish Balay       }
24849b5e25fSSatish Balay     }
24949b5e25fSSatish Balay   }
25049b5e25fSSatish Balay 
25149b5e25fSSatish Balay   /*search column A(0:row-1,row) (=A(row,0:row-1)). Could be expensive! */
2525ddb2528SHong Zhang   /* this segment is currently removed, so only entries in the upper triangle are obtained */
2535ddb2528SHong Zhang #ifdef column_search
25449b5e25fSSatish Balay   v_i    = *v    + M*bs;
25549b5e25fSSatish Balay   cols_i = *cols + M*bs;
25649b5e25fSSatish Balay   for (i=0; i<bn; i++){ /* for each block row */
25749b5e25fSSatish Balay     M = ai[i+1] - ai[i];
25849b5e25fSSatish Balay     for (j=0; j<M; j++){
25949b5e25fSSatish Balay       itmp = aj[ai[i] + j];    /* block column value */
26049b5e25fSSatish Balay       if (itmp == bn){
26149b5e25fSSatish Balay         aa_i   = aa    + bs2*(ai[i] + j) + bs*bp;
26249b5e25fSSatish Balay         for (k=0; k<bs; k++) {
26349b5e25fSSatish Balay           *cols_i++ = i*bs+k;
26449b5e25fSSatish Balay           *v_i++    = aa_i[k];
26549b5e25fSSatish Balay         }
26649b5e25fSSatish Balay         *ncols += bs;
26749b5e25fSSatish Balay         break;
26849b5e25fSSatish Balay       }
26949b5e25fSSatish Balay     }
27049b5e25fSSatish Balay   }
2715ddb2528SHong Zhang #endif
27249b5e25fSSatish Balay 
27349b5e25fSSatish Balay   PetscFunctionReturn(0);
27449b5e25fSSatish Balay }
27549b5e25fSSatish Balay 
27649b5e25fSSatish Balay #undef __FUNC__
27749b5e25fSSatish Balay #define __FUNC__ "MatRestoreRow_SeqSBAIJ"
27849b5e25fSSatish Balay int MatRestoreRow_SeqSBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
27949b5e25fSSatish Balay {
28049b5e25fSSatish Balay   int ierr;
28149b5e25fSSatish Balay 
28249b5e25fSSatish Balay   PetscFunctionBegin;
28349b5e25fSSatish Balay   if (idx) {if (*idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}}
28449b5e25fSSatish Balay   if (v)   {if (*v)   {ierr = PetscFree(*v);CHKERRQ(ierr);}}
28549b5e25fSSatish Balay   PetscFunctionReturn(0);
28649b5e25fSSatish Balay }
28749b5e25fSSatish Balay 
28849b5e25fSSatish Balay #undef __FUNC__
28949b5e25fSSatish Balay #define __FUNC__ "MatTranspose_SeqSBAIJ"
29049b5e25fSSatish Balay int MatTranspose_SeqSBAIJ(Mat A,Mat *B)
29149b5e25fSSatish Balay {
29249b5e25fSSatish Balay   PetscFunctionBegin;
293*29bbc08cSBarry Smith   SETERRQ(1,"Matrix is symmetric. MatTranspose() should not be called");
29416cdd363SHong Zhang   /* PetscFunctionReturn(0); */
29549b5e25fSSatish Balay }
29649b5e25fSSatish Balay 
29749b5e25fSSatish Balay #undef __FUNC__
29849b5e25fSSatish Balay #define __FUNC__ "MatView_SeqSBAIJ_Binary"
29949b5e25fSSatish Balay static int MatView_SeqSBAIJ_Binary(Mat A,Viewer viewer)
30049b5e25fSSatish Balay {
30149b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
30249b5e25fSSatish Balay   int         i,fd,*col_lens,ierr,bs = a->bs,count,*jj,j,k,l,bs2=a->bs2;
30349b5e25fSSatish Balay   Scalar      *aa;
30449b5e25fSSatish Balay   FILE        *file;
30549b5e25fSSatish Balay 
30649b5e25fSSatish Balay   PetscFunctionBegin;
30749b5e25fSSatish Balay   ierr = ViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
30849b5e25fSSatish Balay   col_lens = (int*)PetscMalloc((4+a->m)*sizeof(int));CHKPTRQ(col_lens);
30949b5e25fSSatish Balay   col_lens[0] = MAT_COOKIE;
31049b5e25fSSatish Balay 
31149b5e25fSSatish Balay   col_lens[1] = a->m;
31249b5e25fSSatish Balay   col_lens[2] = a->m;
31349b5e25fSSatish Balay   col_lens[3] = a->s_nz*bs2;
31449b5e25fSSatish Balay 
31549b5e25fSSatish Balay   /* store lengths of each row and write (including header) to file */
31649b5e25fSSatish Balay   count = 0;
31749b5e25fSSatish Balay   for (i=0; i<a->mbs; i++) {
31849b5e25fSSatish Balay     for (j=0; j<bs; j++) {
31949b5e25fSSatish Balay       col_lens[4+count++] = bs*(a->i[i+1] - a->i[i]);
32049b5e25fSSatish Balay     }
32149b5e25fSSatish Balay   }
32249b5e25fSSatish Balay   ierr = PetscBinaryWrite(fd,col_lens,4+a->m,PETSC_INT,1);CHKERRQ(ierr);
32349b5e25fSSatish Balay   ierr = PetscFree(col_lens);CHKERRQ(ierr);
32449b5e25fSSatish Balay 
32549b5e25fSSatish Balay   /* store column indices (zero start index) */
32649b5e25fSSatish Balay   jj = (int*)PetscMalloc((a->s_nz+1)*bs2*sizeof(int));CHKPTRQ(jj);
32749b5e25fSSatish Balay   count = 0;
32849b5e25fSSatish Balay   for (i=0; i<a->mbs; i++) {
32949b5e25fSSatish Balay     for (j=0; j<bs; j++) {
33049b5e25fSSatish Balay       for (k=a->i[i]; k<a->i[i+1]; k++) {
33149b5e25fSSatish Balay         for (l=0; l<bs; l++) {
33249b5e25fSSatish Balay           jj[count++] = bs*a->j[k] + l;
33349b5e25fSSatish Balay         }
33449b5e25fSSatish Balay       }
33549b5e25fSSatish Balay     }
33649b5e25fSSatish Balay   }
33749b5e25fSSatish Balay   ierr = PetscBinaryWrite(fd,jj,bs2*a->s_nz,PETSC_INT,0);CHKERRQ(ierr);
33849b5e25fSSatish Balay   ierr = PetscFree(jj);CHKERRQ(ierr);
33949b5e25fSSatish Balay 
34049b5e25fSSatish Balay   /* store nonzero values */
34149b5e25fSSatish Balay   aa = (Scalar*)PetscMalloc((a->s_nz+1)*bs2*sizeof(Scalar));CHKPTRQ(aa);
34249b5e25fSSatish Balay   count = 0;
34349b5e25fSSatish Balay   for (i=0; i<a->mbs; i++) {
34449b5e25fSSatish Balay     for (j=0; j<bs; j++) {
34549b5e25fSSatish Balay       for (k=a->i[i]; k<a->i[i+1]; k++) {
34649b5e25fSSatish Balay         for (l=0; l<bs; l++) {
34749b5e25fSSatish Balay           aa[count++] = a->a[bs2*k + l*bs + j];
34849b5e25fSSatish Balay         }
34949b5e25fSSatish Balay       }
35049b5e25fSSatish Balay     }
35149b5e25fSSatish Balay   }
35249b5e25fSSatish Balay   ierr = PetscBinaryWrite(fd,aa,bs2*a->s_nz,PETSC_SCALAR,0);CHKERRQ(ierr);
35349b5e25fSSatish Balay   ierr = PetscFree(aa);CHKERRQ(ierr);
35449b5e25fSSatish Balay 
35549b5e25fSSatish Balay   ierr = ViewerBinaryGetInfoPointer(viewer,&file);CHKERRQ(ierr);
35649b5e25fSSatish Balay   if (file) {
35749b5e25fSSatish Balay     fprintf(file,"-matload_block_size %d\n",a->bs);
35849b5e25fSSatish Balay   }
35949b5e25fSSatish Balay   PetscFunctionReturn(0);
36049b5e25fSSatish Balay }
36149b5e25fSSatish Balay 
36249b5e25fSSatish Balay #undef __FUNC__
36349b5e25fSSatish Balay #define __FUNC__ "MatView_SeqSBAIJ_ASCII"
36449b5e25fSSatish Balay static int MatView_SeqSBAIJ_ASCII(Mat A,Viewer viewer)
36549b5e25fSSatish Balay {
36649b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
36749b5e25fSSatish Balay   int         ierr,i,j,format,bs = a->bs,k,l,bs2=a->bs2;
36849b5e25fSSatish Balay   char        *outputname;
36949b5e25fSSatish Balay 
37049b5e25fSSatish Balay   PetscFunctionBegin;
37149b5e25fSSatish Balay   ierr = ViewerGetOutputname(viewer,&outputname);CHKERRQ(ierr);
37249b5e25fSSatish Balay   ierr = ViewerGetFormat(viewer,&format);CHKERRQ(ierr);
37349b5e25fSSatish Balay   if (format == VIEWER_FORMAT_ASCII_INFO || format == VIEWER_FORMAT_ASCII_INFO_LONG) {
37449b5e25fSSatish Balay     ierr = ViewerASCIIPrintf(viewer,"  block size is %d\n",bs);CHKERRQ(ierr);
37549b5e25fSSatish Balay   } else if (format == VIEWER_FORMAT_ASCII_MATLAB) {
376*29bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Matlab format not supported");
37749b5e25fSSatish Balay   } else if (format == VIEWER_FORMAT_ASCII_COMMON) {
37849b5e25fSSatish Balay     ierr = ViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
37949b5e25fSSatish Balay     for (i=0; i<a->mbs; i++) {
38049b5e25fSSatish Balay       for (j=0; j<bs; j++) {
38149b5e25fSSatish Balay         ierr = ViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr);
38249b5e25fSSatish Balay         for (k=a->i[i]; k<a->i[i+1]; k++) {
38349b5e25fSSatish Balay           for (l=0; l<bs; l++) {
38449b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX)
38549b5e25fSSatish Balay             if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
38649b5e25fSSatish Balay               ierr = ViewerASCIIPrintf(viewer," %d %g + %g i",bs*a->j[k]+l,
38749b5e25fSSatish Balay                       PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
38849b5e25fSSatish Balay             } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
38949b5e25fSSatish Balay               ierr = ViewerASCIIPrintf(viewer," %d %g - %g i",bs*a->j[k]+l,
39049b5e25fSSatish Balay                       PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
39149b5e25fSSatish Balay             } else if (PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
39249b5e25fSSatish Balay               ierr = ViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
39349b5e25fSSatish Balay             }
39449b5e25fSSatish Balay #else
39549b5e25fSSatish Balay             if (a->a[bs2*k + l*bs + j] != 0.0) {
39649b5e25fSSatish Balay               ierr = ViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr);
39749b5e25fSSatish Balay             }
39849b5e25fSSatish Balay #endif
39949b5e25fSSatish Balay           }
40049b5e25fSSatish Balay         }
40149b5e25fSSatish Balay         ierr = ViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
40249b5e25fSSatish Balay       }
40349b5e25fSSatish Balay     }
40449b5e25fSSatish Balay     ierr = ViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
40549b5e25fSSatish Balay   } else {
40649b5e25fSSatish Balay     ierr = ViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
40749b5e25fSSatish Balay     for (i=0; i<a->mbs; i++) {
40849b5e25fSSatish Balay       for (j=0; j<bs; j++) {
40949b5e25fSSatish Balay         ierr = ViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr);
41049b5e25fSSatish Balay         for (k=a->i[i]; k<a->i[i+1]; k++) {
41149b5e25fSSatish Balay           for (l=0; l<bs; l++) {
41249b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX)
41349b5e25fSSatish Balay             if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0) {
41449b5e25fSSatish Balay               ierr = ViewerASCIIPrintf(viewer," %d %g + %g i",bs*a->j[k]+l,
41549b5e25fSSatish Balay                 PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
41649b5e25fSSatish Balay             } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0) {
41749b5e25fSSatish Balay               ierr = ViewerASCIIPrintf(viewer," %d %g - %g i",bs*a->j[k]+l,
41849b5e25fSSatish Balay                 PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
41949b5e25fSSatish Balay             } else {
42049b5e25fSSatish Balay               ierr = ViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
42149b5e25fSSatish Balay             }
42249b5e25fSSatish Balay #else
42349b5e25fSSatish Balay             ierr = ViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr);
42449b5e25fSSatish Balay #endif
42549b5e25fSSatish Balay           }
42649b5e25fSSatish Balay         }
42749b5e25fSSatish Balay         ierr = ViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
42849b5e25fSSatish Balay       }
42949b5e25fSSatish Balay     }
43049b5e25fSSatish Balay     ierr = ViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
43149b5e25fSSatish Balay   }
43249b5e25fSSatish Balay   ierr = ViewerFlush(viewer);CHKERRQ(ierr);
43349b5e25fSSatish Balay   PetscFunctionReturn(0);
43449b5e25fSSatish Balay }
43549b5e25fSSatish Balay 
43649b5e25fSSatish Balay #undef __FUNC__
43749b5e25fSSatish Balay #define __FUNC__ "MatView_SeqSBAIJ_Draw_Zoom"
43849b5e25fSSatish Balay static int MatView_SeqSBAIJ_Draw_Zoom(Draw draw,void *Aa)
43949b5e25fSSatish Balay {
44049b5e25fSSatish Balay   Mat          A = (Mat) Aa;
44149b5e25fSSatish Balay   Mat_SeqSBAIJ  *a=(Mat_SeqSBAIJ*)A->data;
44249b5e25fSSatish Balay   int          row,ierr,i,j,k,l,mbs=a->mbs,color,bs=a->bs,bs2=a->bs2,rank;
44349b5e25fSSatish Balay   PetscReal    xl,yl,xr,yr,x_l,x_r,y_l,y_r;
44449b5e25fSSatish Balay   MatScalar    *aa;
44549b5e25fSSatish Balay   MPI_Comm     comm;
44649b5e25fSSatish Balay   Viewer       viewer;
44749b5e25fSSatish Balay 
44849b5e25fSSatish Balay   PetscFunctionBegin;
44949b5e25fSSatish Balay   /*
45049b5e25fSSatish Balay       This is nasty. If this is called from an originally parallel matrix
45149b5e25fSSatish Balay    then all processes call this,but only the first has the matrix so the
45249b5e25fSSatish Balay    rest should return immediately.
45349b5e25fSSatish Balay   */
45449b5e25fSSatish Balay   ierr = PetscObjectGetComm((PetscObject)draw,&comm);CHKERRQ(ierr);
45549b5e25fSSatish Balay   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
45649b5e25fSSatish Balay   if (rank) PetscFunctionReturn(0);
45749b5e25fSSatish Balay 
45849b5e25fSSatish Balay   ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr);
45949b5e25fSSatish Balay 
46049b5e25fSSatish Balay   ierr = DrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr);
46149b5e25fSSatish Balay   DrawString(draw, .3*(xl+xr), .3*(yl+yr), DRAW_BLACK, "symmetric");
46249b5e25fSSatish Balay 
46349b5e25fSSatish Balay   /* loop over matrix elements drawing boxes */
46449b5e25fSSatish Balay   color = DRAW_BLUE;
46549b5e25fSSatish Balay   for (i=0,row=0; i<mbs; i++,row+=bs) {
46649b5e25fSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
46749b5e25fSSatish Balay       y_l = a->m - row - 1.0; y_r = y_l + 1.0;
46849b5e25fSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
46949b5e25fSSatish Balay       aa = a->a + j*bs2;
47049b5e25fSSatish Balay       for (k=0; k<bs; k++) {
47149b5e25fSSatish Balay         for (l=0; l<bs; l++) {
47249b5e25fSSatish Balay           if (PetscRealPart(*aa++) >=  0.) continue;
47349b5e25fSSatish Balay           ierr = DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
47449b5e25fSSatish Balay         }
47549b5e25fSSatish Balay       }
47649b5e25fSSatish Balay     }
47749b5e25fSSatish Balay   }
47849b5e25fSSatish Balay   color = DRAW_CYAN;
47949b5e25fSSatish Balay   for (i=0,row=0; i<mbs; i++,row+=bs) {
48049b5e25fSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
48149b5e25fSSatish Balay       y_l = a->m - row - 1.0; y_r = y_l + 1.0;
48249b5e25fSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
48349b5e25fSSatish Balay       aa = a->a + j*bs2;
48449b5e25fSSatish Balay       for (k=0; k<bs; k++) {
48549b5e25fSSatish Balay         for (l=0; l<bs; l++) {
48649b5e25fSSatish Balay           if (PetscRealPart(*aa++) != 0.) continue;
48749b5e25fSSatish Balay           ierr = DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
48849b5e25fSSatish Balay         }
48949b5e25fSSatish Balay       }
49049b5e25fSSatish Balay     }
49149b5e25fSSatish Balay   }
49249b5e25fSSatish Balay 
49349b5e25fSSatish Balay   color = DRAW_RED;
49449b5e25fSSatish Balay   for (i=0,row=0; i<mbs; i++,row+=bs) {
49549b5e25fSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
49649b5e25fSSatish Balay       y_l = a->m - row - 1.0; y_r = y_l + 1.0;
49749b5e25fSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
49849b5e25fSSatish Balay       aa = a->a + j*bs2;
49949b5e25fSSatish Balay       for (k=0; k<bs; k++) {
50049b5e25fSSatish Balay         for (l=0; l<bs; l++) {
50149b5e25fSSatish Balay           if (PetscRealPart(*aa++) <= 0.) continue;
50249b5e25fSSatish Balay           ierr = DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
50349b5e25fSSatish Balay         }
50449b5e25fSSatish Balay       }
50549b5e25fSSatish Balay     }
50649b5e25fSSatish Balay   }
50749b5e25fSSatish Balay   PetscFunctionReturn(0);
50849b5e25fSSatish Balay }
50949b5e25fSSatish Balay 
51049b5e25fSSatish Balay #undef __FUNC__
51149b5e25fSSatish Balay #define __FUNC__ "MatView_SeqSBAIJ_Draw"
51249b5e25fSSatish Balay static int MatView_SeqSBAIJ_Draw(Mat A,Viewer viewer)
51349b5e25fSSatish Balay {
51449b5e25fSSatish Balay   Mat_SeqSBAIJ  *a=(Mat_SeqSBAIJ*)A->data;
51549b5e25fSSatish Balay   int          ierr;
51649b5e25fSSatish Balay   PetscReal    xl,yl,xr,yr,w,h;
51749b5e25fSSatish Balay   Draw         draw;
51849b5e25fSSatish Balay   PetscTruth   isnull;
51949b5e25fSSatish Balay 
52049b5e25fSSatish Balay   PetscFunctionBegin;
52149b5e25fSSatish Balay 
52249b5e25fSSatish Balay   ierr = ViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
52349b5e25fSSatish Balay   ierr = DrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
52449b5e25fSSatish Balay 
52549b5e25fSSatish Balay   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
52649b5e25fSSatish Balay   xr  = a->m; yr = a->m; h = yr/10.0; w = xr/10.0;
52749b5e25fSSatish Balay   xr += w;    yr += h;  xl = -w;     yl = -h;
52849b5e25fSSatish Balay   ierr = DrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
52949b5e25fSSatish Balay   ierr = DrawZoom(draw,MatView_SeqSBAIJ_Draw_Zoom,A);CHKERRQ(ierr);
53049b5e25fSSatish Balay   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr);
53149b5e25fSSatish Balay   PetscFunctionReturn(0);
53249b5e25fSSatish Balay }
53349b5e25fSSatish Balay 
53449b5e25fSSatish Balay #undef __FUNC__
53549b5e25fSSatish Balay #define __FUNC__ "MatView_SeqSBAIJ"
53649b5e25fSSatish Balay int MatView_SeqSBAIJ(Mat A,Viewer viewer)
53749b5e25fSSatish Balay {
53849b5e25fSSatish Balay   int        ierr;
53949b5e25fSSatish Balay   PetscTruth issocket,isascii,isbinary,isdraw;
54049b5e25fSSatish Balay 
54149b5e25fSSatish Balay   PetscFunctionBegin;
54249b5e25fSSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,SOCKET_VIEWER,&issocket);CHKERRQ(ierr);
54349b5e25fSSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,ASCII_VIEWER,&isascii);CHKERRQ(ierr);
54449b5e25fSSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,BINARY_VIEWER,&isbinary);CHKERRQ(ierr);
54549b5e25fSSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,DRAW_VIEWER,&isdraw);CHKERRQ(ierr);
54649b5e25fSSatish Balay   if (issocket) {
547*29bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Socket viewer not supported");
54849b5e25fSSatish Balay   } else if (isascii){
54949b5e25fSSatish Balay     ierr = MatView_SeqSBAIJ_ASCII(A,viewer);CHKERRQ(ierr);
55049b5e25fSSatish Balay   } else if (isbinary) {
55149b5e25fSSatish Balay     ierr = MatView_SeqSBAIJ_Binary(A,viewer);CHKERRQ(ierr);
55249b5e25fSSatish Balay   } else if (isdraw) {
55349b5e25fSSatish Balay     ierr = MatView_SeqSBAIJ_Draw(A,viewer);CHKERRQ(ierr);
55449b5e25fSSatish Balay   } else {
555*29bbc08cSBarry Smith     SETERRQ1(1,"Viewer type %s not supported by SeqBAIJ matrices",((PetscObject)viewer)->type_name);
55649b5e25fSSatish Balay   }
55749b5e25fSSatish Balay   PetscFunctionReturn(0);
55849b5e25fSSatish Balay }
55949b5e25fSSatish Balay 
56049b5e25fSSatish Balay 
56149b5e25fSSatish Balay #undef __FUNC__
56249b5e25fSSatish Balay #define __FUNC__ "MatGetValues_SeqSBAIJ"
56349b5e25fSSatish Balay int MatGetValues_SeqSBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v)
56449b5e25fSSatish Balay {
565045c9aa0SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
56649b5e25fSSatish Balay   int        *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j;
56749b5e25fSSatish Balay   int        *ai = a->i,*ailen = a->ilen;
56849b5e25fSSatish Balay   int        brow,bcol,ridx,cidx,bs=a->bs,bs2=a->bs2;
56949b5e25fSSatish Balay   MatScalar  *ap,*aa = a->a,zero = 0.0;
57049b5e25fSSatish Balay 
57149b5e25fSSatish Balay   PetscFunctionBegin;
57249b5e25fSSatish Balay   for (k=0; k<m; k++) { /* loop over rows */
57349b5e25fSSatish Balay     row  = im[k]; brow = row/bs;
574*29bbc08cSBarry Smith     if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative row");
575*29bbc08cSBarry Smith     if (row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row too large");
57649b5e25fSSatish Balay     rp   = aj + ai[brow] ; ap = aa + bs2*ai[brow] ;
57749b5e25fSSatish Balay     nrow = ailen[brow];
57849b5e25fSSatish Balay     for (l=0; l<n; l++) { /* loop over columns */
579*29bbc08cSBarry Smith       if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative column");
580*29bbc08cSBarry Smith       if (in[l] >= a->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Column too large");
58149b5e25fSSatish Balay       col  = in[l] ;
58249b5e25fSSatish Balay       bcol = col/bs;
58349b5e25fSSatish Balay       cidx = col%bs;
58449b5e25fSSatish Balay       ridx = row%bs;
58549b5e25fSSatish Balay       high = nrow;
58649b5e25fSSatish Balay       low  = 0; /* assume unsorted */
58749b5e25fSSatish Balay       while (high-low > 5) {
58849b5e25fSSatish Balay         t = (low+high)/2;
58949b5e25fSSatish Balay         if (rp[t] > bcol) high = t;
59049b5e25fSSatish Balay         else             low  = t;
59149b5e25fSSatish Balay       }
59249b5e25fSSatish Balay       for (i=low; i<high; i++) {
59349b5e25fSSatish Balay         if (rp[i] > bcol) break;
59449b5e25fSSatish Balay         if (rp[i] == bcol) {
59549b5e25fSSatish Balay           *v++ = ap[bs2*i+bs*cidx+ridx];
59649b5e25fSSatish Balay           goto finished;
59749b5e25fSSatish Balay         }
59849b5e25fSSatish Balay       }
59949b5e25fSSatish Balay       *v++ = zero;
60049b5e25fSSatish Balay       finished:;
60149b5e25fSSatish Balay     }
60249b5e25fSSatish Balay   }
60349b5e25fSSatish Balay   PetscFunctionReturn(0);
60449b5e25fSSatish Balay }
60549b5e25fSSatish Balay 
60649b5e25fSSatish Balay 
60749b5e25fSSatish Balay #undef __FUNC__
60849b5e25fSSatish Balay #define __FUNC__ "MatSetValuesBlocked_SeqSBAIJ"
60949b5e25fSSatish Balay int MatSetValuesBlocked_SeqSBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is)
61049b5e25fSSatish Balay {
61149b5e25fSSatish Balay   PetscFunctionBegin;
612*29bbc08cSBarry Smith   SETERRQ(1,"Function not yet written for SBAIJ format");
61316cdd363SHong Zhang   /* PetscFunctionReturn(0); */
61449b5e25fSSatish Balay }
61549b5e25fSSatish Balay 
61649b5e25fSSatish Balay #undef __FUNC__
61749b5e25fSSatish Balay #define __FUNC__ "MatAssemblyEnd_SeqSBAIJ"
61849b5e25fSSatish Balay int MatAssemblyEnd_SeqSBAIJ(Mat A,MatAssemblyType mode)
61949b5e25fSSatish Balay {
62049b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
62149b5e25fSSatish Balay   int        fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax;
62249b5e25fSSatish Balay   int        m = a->m,*ip,N,*ailen = a->ilen;
62349b5e25fSSatish Balay   int        mbs = a->mbs,bs2 = a->bs2,rmax = 0,ierr;
62449b5e25fSSatish Balay   MatScalar  *aa = a->a,*ap;
62549b5e25fSSatish Balay 
62649b5e25fSSatish Balay   PetscFunctionBegin;
62749b5e25fSSatish Balay   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0);
62849b5e25fSSatish Balay 
62949b5e25fSSatish Balay   if (m) rmax = ailen[0];
63049b5e25fSSatish Balay   for (i=1; i<mbs; i++) {
63149b5e25fSSatish Balay     /* move each row back by the amount of empty slots (fshift) before it*/
63249b5e25fSSatish Balay     fshift += imax[i-1] - ailen[i-1];
63349b5e25fSSatish Balay     rmax   = PetscMax(rmax,ailen[i]);
63449b5e25fSSatish Balay     if (fshift) {
63549b5e25fSSatish Balay       ip = aj + ai[i]; ap = aa + bs2*ai[i];
63649b5e25fSSatish Balay       N = ailen[i];
63749b5e25fSSatish Balay       for (j=0; j<N; j++) {
63849b5e25fSSatish Balay         ip[j-fshift] = ip[j];
63949b5e25fSSatish Balay         ierr = PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
64049b5e25fSSatish Balay       }
64149b5e25fSSatish Balay     }
64249b5e25fSSatish Balay     ai[i] = ai[i-1] + ailen[i-1];
64349b5e25fSSatish Balay   }
64449b5e25fSSatish Balay   if (mbs) {
64549b5e25fSSatish Balay     fshift += imax[mbs-1] - ailen[mbs-1];
64649b5e25fSSatish Balay     ai[mbs] = ai[mbs-1] + ailen[mbs-1];
64749b5e25fSSatish Balay   }
64849b5e25fSSatish Balay   /* reset ilen and imax for each row */
64949b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
65049b5e25fSSatish Balay     ailen[i] = imax[i] = ai[i+1] - ai[i];
65149b5e25fSSatish Balay   }
65249b5e25fSSatish Balay   a->s_nz = ai[mbs];
65349b5e25fSSatish Balay 
65449b5e25fSSatish Balay   /* diagonals may have moved, so kill the diagonal pointers */
65549b5e25fSSatish Balay   if (fshift && a->diag) {
65649b5e25fSSatish Balay     ierr = PetscFree(a->diag);CHKERRQ(ierr);
65749b5e25fSSatish Balay     PLogObjectMemory(A,-(m+1)*sizeof(int));
65849b5e25fSSatish Balay     a->diag = 0;
65949b5e25fSSatish Balay   }
66049b5e25fSSatish Balay   PLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Matrix size: %d X %d, block size %d; storage space: %d unneeded, %d used\n",
66149b5e25fSSatish Balay            m,a->m,a->bs,fshift*bs2,a->s_nz*bs2);
66249b5e25fSSatish Balay   PLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Number of mallocs during MatSetValues is %d\n",
66349b5e25fSSatish Balay            a->reallocs);
66449b5e25fSSatish Balay   PLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Most nonzeros blocks in any row is %d\n",rmax);
66549b5e25fSSatish Balay   a->reallocs          = 0;
66649b5e25fSSatish Balay   A->info.nz_unneeded  = (PetscReal)fshift*bs2;
66749b5e25fSSatish Balay 
66849b5e25fSSatish Balay   PetscFunctionReturn(0);
66949b5e25fSSatish Balay }
67049b5e25fSSatish Balay 
67149b5e25fSSatish Balay /*
67249b5e25fSSatish Balay    This function returns an array of flags which indicate the locations of contiguous
67349b5e25fSSatish Balay    blocks that should be zeroed. for eg: if bs = 3  and is = [0,1,2,3,5,6,7,8,9]
67449b5e25fSSatish Balay    then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)]
67549b5e25fSSatish Balay    Assume: sizes should be long enough to hold all the values.
67649b5e25fSSatish Balay */
67749b5e25fSSatish Balay #undef __FUNC__
67849b5e25fSSatish Balay #define __FUNC__ "MatZeroRows_SeqSBAIJ_Check_Blocks"
67949b5e25fSSatish Balay static int MatZeroRows_SeqSBAIJ_Check_Blocks(int idx[],int n,int bs,int sizes[], int *bs_max)
68049b5e25fSSatish Balay {
68149b5e25fSSatish Balay   int        i,j,k,row;
68249b5e25fSSatish Balay   PetscTruth flg;
68349b5e25fSSatish Balay 
68449b5e25fSSatish Balay   PetscFunctionBegin;
68549b5e25fSSatish Balay   for (i=0,j=0; i<n; j++) {
68649b5e25fSSatish Balay     row = idx[i];
68749b5e25fSSatish Balay     if (row%bs!=0) { /* Not the begining of a block */
68849b5e25fSSatish Balay       sizes[j] = 1;
68949b5e25fSSatish Balay       i++;
69049b5e25fSSatish Balay     } else if (i+bs > n) { /* Beginning of a block, but complete block doesn't exist (at idx end) */
69149b5e25fSSatish Balay       sizes[j] = 1;         /* Also makes sure atleast 'bs' values exist for next else */
69249b5e25fSSatish Balay       i++;
69349b5e25fSSatish Balay     } else { /* Begining of the block, so check if the complete block exists */
69449b5e25fSSatish Balay       flg = PETSC_TRUE;
69549b5e25fSSatish Balay       for (k=1; k<bs; k++) {
69649b5e25fSSatish Balay         if (row+k != idx[i+k]) { /* break in the block */
69749b5e25fSSatish Balay           flg = PETSC_FALSE;
69849b5e25fSSatish Balay           break;
69949b5e25fSSatish Balay         }
70049b5e25fSSatish Balay       }
70149b5e25fSSatish Balay       if (flg == PETSC_TRUE) { /* No break in the bs */
70249b5e25fSSatish Balay         sizes[j] = bs;
70349b5e25fSSatish Balay         i+= bs;
70449b5e25fSSatish Balay       } else {
70549b5e25fSSatish Balay         sizes[j] = 1;
70649b5e25fSSatish Balay         i++;
70749b5e25fSSatish Balay       }
70849b5e25fSSatish Balay     }
70949b5e25fSSatish Balay   }
71049b5e25fSSatish Balay   *bs_max = j;
71149b5e25fSSatish Balay   PetscFunctionReturn(0);
71249b5e25fSSatish Balay }
71349b5e25fSSatish Balay 
71449b5e25fSSatish Balay #undef __FUNC__
71549b5e25fSSatish Balay #define __FUNC__ "MatZeroRows_SeqSBAIJ"
71649b5e25fSSatish Balay int MatZeroRows_SeqSBAIJ(Mat A,IS is,Scalar *diag)
71749b5e25fSSatish Balay {
71849b5e25fSSatish Balay   Mat_SeqSBAIJ *sbaij=(Mat_SeqSBAIJ*)A->data;
71949b5e25fSSatish Balay   int         ierr,i,j,k,count,is_n,*is_idx,*rows;
72049b5e25fSSatish Balay   int         bs=sbaij->bs,bs2=sbaij->bs2,*sizes,row,bs_max;
72149b5e25fSSatish Balay   Scalar      zero = 0.0;
72249b5e25fSSatish Balay   MatScalar   *aa;
72349b5e25fSSatish Balay 
72449b5e25fSSatish Balay   PetscFunctionBegin;
72549b5e25fSSatish Balay   /* Make a copy of the IS and  sort it */
72649b5e25fSSatish Balay   ierr = ISGetSize(is,&is_n);CHKERRQ(ierr);
72749b5e25fSSatish Balay   ierr = ISGetIndices(is,&is_idx);CHKERRQ(ierr);
72849b5e25fSSatish Balay 
72949b5e25fSSatish Balay   /* allocate memory for rows,sizes */
73049b5e25fSSatish Balay   rows  = (int*)PetscMalloc((3*is_n+1)*sizeof(int));CHKPTRQ(rows);
73149b5e25fSSatish Balay   sizes = rows + is_n;
73249b5e25fSSatish Balay 
73349b5e25fSSatish Balay   /* initialize copy IS values to rows, and sort them */
73449b5e25fSSatish Balay   for (i=0; i<is_n; i++) { rows[i] = is_idx[i]; }
73549b5e25fSSatish Balay   ierr = PetscSortInt(is_n,rows);CHKERRQ(ierr);
73649b5e25fSSatish Balay   if (sbaij->keepzeroedrows) { /* do not change nonzero structure */
73749b5e25fSSatish Balay     for (i=0; i<is_n; i++) { sizes[i] = 1; } /* sizes: size of blocks, = 1 or bs */
73849b5e25fSSatish Balay     bs_max = is_n;            /* bs_max: num. of contiguous block row in the row */
73949b5e25fSSatish Balay   } else {
74049b5e25fSSatish Balay     ierr = MatZeroRows_SeqSBAIJ_Check_Blocks(rows,is_n,bs,sizes,&bs_max);CHKERRQ(ierr);
74149b5e25fSSatish Balay   }
74249b5e25fSSatish Balay   ierr = ISRestoreIndices(is,&is_idx);CHKERRQ(ierr);
74349b5e25fSSatish Balay 
74449b5e25fSSatish Balay   for (i=0,j=0; i<bs_max; j+=sizes[i],i++) {
74549b5e25fSSatish Balay     row   = rows[j];                  /* row to be zeroed */
746*29bbc08cSBarry Smith     if (row < 0 || row > sbaij->m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"row %d out of range",row);
74749b5e25fSSatish Balay     count = (sbaij->i[row/bs +1] - sbaij->i[row/bs])*bs; /* num. of elements in the row */
74849b5e25fSSatish Balay     aa    = sbaij->a + sbaij->i[row/bs]*bs2 + (row%bs);
74949b5e25fSSatish Balay     if (sizes[i] == bs && !sbaij->keepzeroedrows) {
75049b5e25fSSatish Balay       if (diag) {
75149b5e25fSSatish Balay         if (sbaij->ilen[row/bs] > 0) {
75249b5e25fSSatish Balay           sbaij->ilen[row/bs] = 1;
75349b5e25fSSatish Balay           sbaij->j[sbaij->i[row/bs]] = row/bs;
75449b5e25fSSatish Balay           ierr = PetscMemzero(aa,count*bs*sizeof(MatScalar));CHKERRQ(ierr);
75549b5e25fSSatish Balay         }
75649b5e25fSSatish Balay         /* Now insert all the diagoanl values for this bs */
75749b5e25fSSatish Balay         for (k=0; k<bs; k++) {
75849b5e25fSSatish Balay           ierr = (*A->ops->setvalues)(A,1,rows+j+k,1,rows+j+k,diag,INSERT_VALUES);CHKERRQ(ierr);
75949b5e25fSSatish Balay         }
76049b5e25fSSatish Balay       } else { /* (!diag) */
76149b5e25fSSatish Balay         sbaij->ilen[row/bs] = 0;
76249b5e25fSSatish Balay       } /* end (!diag) */
76349b5e25fSSatish Balay     } else { /* (sizes[i] != bs), broken block */
76449b5e25fSSatish Balay #if defined (PETSC_USE_DEBUG)
765*29bbc08cSBarry Smith       if (sizes[i] != 1) SETERRQ(1,"Internal Error. Value should be 1");
76649b5e25fSSatish Balay #endif
76749b5e25fSSatish Balay       for (k=0; k<count; k++) {
76849b5e25fSSatish Balay         aa[0] = zero;
76949b5e25fSSatish Balay         aa+=bs;
77049b5e25fSSatish Balay       }
77149b5e25fSSatish Balay       if (diag) {
77249b5e25fSSatish Balay         ierr = (*A->ops->setvalues)(A,1,rows+j,1,rows+j,diag,INSERT_VALUES);CHKERRQ(ierr);
77349b5e25fSSatish Balay       }
77449b5e25fSSatish Balay     }
77549b5e25fSSatish Balay   }
77649b5e25fSSatish Balay 
77749b5e25fSSatish Balay   ierr = PetscFree(rows);CHKERRQ(ierr);
77849b5e25fSSatish Balay   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
77949b5e25fSSatish Balay   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
78049b5e25fSSatish Balay   PetscFunctionReturn(0);
78149b5e25fSSatish Balay }
78249b5e25fSSatish Balay 
78349b5e25fSSatish Balay /* Only add/insert a(i,j) with i<=j (blocks).
78449b5e25fSSatish Balay    Any a(i,j) with i>j input by user is ingored.
78549b5e25fSSatish Balay */
78649b5e25fSSatish Balay 
78749b5e25fSSatish Balay #undef __FUNC__
78849b5e25fSSatish Balay #define __FUNC__ "MatSetValues_SeqSBAIJ"
78949b5e25fSSatish Balay int MatSetValues_SeqSBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is)
79049b5e25fSSatish Balay {
79149b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
79249b5e25fSSatish Balay   int         *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted=a->sorted;
79349b5e25fSSatish Balay   int         *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented;
79449b5e25fSSatish Balay   int         *aj=a->j,nonew=a->nonew,bs=a->bs,brow,bcol;
79549b5e25fSSatish Balay   int         ridx,cidx,bs2=a->bs2,ierr;
79649b5e25fSSatish Balay   MatScalar   *ap,value,*aa=a->a,*bap;
79749b5e25fSSatish Balay 
79849b5e25fSSatish Balay   PetscFunctionBegin;
79949b5e25fSSatish Balay 
80049b5e25fSSatish Balay   for (k=0; k<m; k++) { /* loop over added rows */
80149b5e25fSSatish Balay     row  = im[k];       /* row number */
80249b5e25fSSatish Balay     brow = row/bs;      /* block row number */
80349b5e25fSSatish Balay     if (row < 0) continue;
80449b5e25fSSatish Balay #if defined(PETSC_USE_BOPT_g)
805*29bbc08cSBarry Smith     if (row >= a->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",row,a->m);
80649b5e25fSSatish Balay #endif
80749b5e25fSSatish Balay     rp   = aj + ai[brow]; /*ptr to beginning of column value of the row block*/
80849b5e25fSSatish Balay     ap   = aa + bs2*ai[brow]; /*ptr to beginning of element value of the row block*/
80949b5e25fSSatish Balay     rmax = imax[brow];  /* maximum space allocated for this row */
81049b5e25fSSatish Balay     nrow = ailen[brow]; /* actual length of this row */
81149b5e25fSSatish Balay     low  = 0;
81249b5e25fSSatish Balay 
81349b5e25fSSatish Balay     for (l=0; l<n; l++) { /* loop over added columns */
81449b5e25fSSatish Balay       if (in[l] < 0) continue;
81549b5e25fSSatish Balay #if defined(PETSC_USE_BOPT_g)
816*29bbc08cSBarry Smith       if (in[l] >= a->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",in[l],a->m);
81749b5e25fSSatish Balay #endif
81849b5e25fSSatish Balay       col = in[l];
81949b5e25fSSatish Balay       bcol = col/bs;              /* block col number */
82049b5e25fSSatish Balay 
82149b5e25fSSatish Balay       if (brow <= bcol){
82249b5e25fSSatish Balay         ridx = row % bs; cidx = col % bs; /*row and col index inside the block */
82349b5e25fSSatish Balay         /* if ((brow==bcol && ridx<=cidx) || (brow<bcol)){ */
82449b5e25fSSatish Balay           /* element value a(k,l) */
82549b5e25fSSatish Balay           if (roworiented) {
82649b5e25fSSatish Balay             value = v[l + k*n];
82749b5e25fSSatish Balay           } else {
82849b5e25fSSatish Balay             value = v[k + l*m];
82949b5e25fSSatish Balay           }
83049b5e25fSSatish Balay 
83149b5e25fSSatish Balay           /* move pointer bap to a(k,l) quickly and add/insert value */
83249b5e25fSSatish Balay           if (!sorted) low = 0; high = nrow;
83349b5e25fSSatish Balay           while (high-low > 7) {
83449b5e25fSSatish Balay             t = (low+high)/2;
83549b5e25fSSatish Balay             if (rp[t] > bcol) high = t;
83649b5e25fSSatish Balay             else              low  = t;
83749b5e25fSSatish Balay           }
83849b5e25fSSatish Balay           for (i=low; i<high; i++) {
83949b5e25fSSatish Balay             /* printf("The loop of i=low.., rp[%d]=%d\n",i,rp[i]); */
84049b5e25fSSatish Balay             if (rp[i] > bcol) break;
84149b5e25fSSatish Balay             if (rp[i] == bcol) {
84249b5e25fSSatish Balay               bap  = ap +  bs2*i + bs*cidx + ridx;
84349b5e25fSSatish Balay               if (is == ADD_VALUES) *bap += value;
84449b5e25fSSatish Balay               else                  *bap  = value;
84549b5e25fSSatish Balay               goto noinsert1;
84649b5e25fSSatish Balay             }
84749b5e25fSSatish Balay           }
84849b5e25fSSatish Balay 
84949b5e25fSSatish Balay       if (nonew == 1) goto noinsert1;
850*29bbc08cSBarry Smith       else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix");
85149b5e25fSSatish Balay       if (nrow >= rmax) {
85249b5e25fSSatish Balay         /* there is no extra room in row, therefore enlarge */
85349b5e25fSSatish Balay         int       new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j;
85449b5e25fSSatish Balay         MatScalar *new_a;
85549b5e25fSSatish Balay 
856*29bbc08cSBarry Smith         if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix");
85749b5e25fSSatish Balay 
85849b5e25fSSatish Balay         /* Malloc new storage space */
85949b5e25fSSatish Balay         len     = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int);
86049b5e25fSSatish Balay         new_a   = (MatScalar*)PetscMalloc(len);CHKPTRQ(new_a);
86149b5e25fSSatish Balay         new_j   = (int*)(new_a + bs2*new_nz);
86249b5e25fSSatish Balay         new_i   = new_j + new_nz;
86349b5e25fSSatish Balay 
86449b5e25fSSatish Balay         /* copy over old data into new slots */
86549b5e25fSSatish Balay         for (ii=0; ii<brow+1; ii++) {new_i[ii] = ai[ii];}
86649b5e25fSSatish Balay         for (ii=brow+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;}
86749b5e25fSSatish Balay         ierr = PetscMemcpy(new_j,aj,(ai[brow]+nrow)*sizeof(int));CHKERRQ(ierr);
86849b5e25fSSatish Balay         len  = (new_nz - CHUNKSIZE - ai[brow] - nrow);
86949b5e25fSSatish Balay         ierr = PetscMemcpy(new_j+ai[brow]+nrow+CHUNKSIZE,aj+ai[brow]+nrow,len*sizeof(int));CHKERRQ(ierr);
87049b5e25fSSatish Balay         ierr = PetscMemcpy(new_a,aa,(ai[brow]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr);
87149b5e25fSSatish Balay         ierr = PetscMemzero(new_a+bs2*(ai[brow]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr);
87249b5e25fSSatish Balay         ierr = PetscMemcpy(new_a+bs2*(ai[brow]+nrow+CHUNKSIZE),aa+bs2*(ai[brow]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr);
87349b5e25fSSatish Balay         /* free up old matrix storage */
87449b5e25fSSatish Balay         ierr = PetscFree(a->a);CHKERRQ(ierr);
87549b5e25fSSatish Balay         if (!a->singlemalloc) {
87649b5e25fSSatish Balay           ierr = PetscFree(a->i);CHKERRQ(ierr);
87749b5e25fSSatish Balay           ierr = PetscFree(a->j);CHKERRQ(ierr);
87849b5e25fSSatish Balay         }
87949b5e25fSSatish Balay         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
88049b5e25fSSatish Balay         a->singlemalloc = PETSC_TRUE;
88149b5e25fSSatish Balay 
88249b5e25fSSatish Balay         rp   = aj + ai[brow]; ap = aa + bs2*ai[brow];
88349b5e25fSSatish Balay         rmax = imax[brow] = imax[brow] + CHUNKSIZE;
88449b5e25fSSatish Balay         PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar)));
88549b5e25fSSatish Balay         a->s_maxnz += bs2*CHUNKSIZE;
88649b5e25fSSatish Balay         a->reallocs++;
88749b5e25fSSatish Balay         a->s_nz++;
88849b5e25fSSatish Balay       }
88949b5e25fSSatish Balay 
89049b5e25fSSatish Balay       N = nrow++ - 1;
89149b5e25fSSatish Balay       /* shift up all the later entries in this row */
89249b5e25fSSatish Balay       for (ii=N; ii>=i; ii--) {
89349b5e25fSSatish Balay         rp[ii+1] = rp[ii];
89449b5e25fSSatish Balay         ierr     = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr);
89549b5e25fSSatish Balay       }
89649b5e25fSSatish Balay       if (N>=i) {
89749b5e25fSSatish Balay         ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr);
89849b5e25fSSatish Balay       }
89949b5e25fSSatish Balay       rp[i]                      = bcol;
90049b5e25fSSatish Balay       ap[bs2*i + bs*cidx + ridx] = value;
90149b5e25fSSatish Balay       noinsert1:;
90249b5e25fSSatish Balay       low = i;
90349b5e25fSSatish Balay       /* } */
90449b5e25fSSatish Balay     } /* end of if .. if..  */
90549b5e25fSSatish Balay     }                     /* end of loop over added columns */
90649b5e25fSSatish Balay     ailen[brow] = nrow;
90749b5e25fSSatish Balay   }                       /* end of loop over added rows */
90849b5e25fSSatish Balay 
90949b5e25fSSatish Balay   PetscFunctionReturn(0);
91049b5e25fSSatish Balay }
91149b5e25fSSatish Balay 
912915354c9SHong Zhang extern int MatCholeskyFactorSymbolic_SeqSBAIJ(Mat,IS,PetscReal,Mat*);
913915354c9SHong Zhang extern int MatCholeskyFactor_SeqSBAIJ(Mat,IS,PetscReal);
91449b5e25fSSatish Balay extern int MatIncreaseOverlap_SeqSBAIJ(Mat,int,IS*,int);
91549b5e25fSSatish Balay extern int MatGetSubMatrix_SeqSBAIJ(Mat,IS,IS,int,MatReuse,Mat*);
91649b5e25fSSatish Balay extern int MatGetSubMatrices_SeqSBAIJ(Mat,int,IS*,IS*,MatReuse,Mat**);
91749b5e25fSSatish Balay extern int MatMultTranspose_SeqSBAIJ(Mat,Vec,Vec);
91849b5e25fSSatish Balay extern int MatMultTransposeAdd_SeqSBAIJ(Mat,Vec,Vec,Vec);
91949b5e25fSSatish Balay extern int MatScale_SeqSBAIJ(Scalar*,Mat);
92049b5e25fSSatish Balay extern int MatNorm_SeqSBAIJ(Mat,NormType,PetscReal *);
92149b5e25fSSatish Balay extern int MatEqual_SeqSBAIJ(Mat,Mat,PetscTruth*);
92249b5e25fSSatish Balay extern int MatGetDiagonal_SeqSBAIJ(Mat,Vec);
92349b5e25fSSatish Balay extern int MatDiagonalScale_SeqSBAIJ(Mat,Vec,Vec);
92449b5e25fSSatish Balay extern int MatGetInfo_SeqSBAIJ(Mat,MatInfoType,MatInfo *);
92549b5e25fSSatish Balay extern int MatZeroEntries_SeqSBAIJ(Mat);
92649b5e25fSSatish Balay 
92749b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_N(Mat,Vec,Vec);
92849b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_1(Mat,Vec,Vec);
92949b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_2(Mat,Vec,Vec);
93049b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_3(Mat,Vec,Vec);
93149b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_4(Mat,Vec,Vec);
93249b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_5(Mat,Vec,Vec);
93349b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_6(Mat,Vec,Vec);
93449b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_7(Mat,Vec,Vec);
93549b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_7(Mat,Vec,Vec);
93649b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_6(Mat,Vec,Vec);
93749b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_5(Mat,Vec,Vec);
93849b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_4(Mat,Vec,Vec);
93949b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_3(Mat,Vec,Vec);
94049b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_2(Mat,Vec,Vec);
94149b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_1(Mat,Vec,Vec);
94249b5e25fSSatish Balay 
9435f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat,Mat*);
9445f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_1(Mat,Mat*);
9455f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_2(Mat,Mat*);
9465f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_3(Mat,Mat*);
9475f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_4(Mat,Mat*);
9485f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_5(Mat,Mat*);
9495f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_6(Mat,Mat*);
95049b5e25fSSatish Balay 
95149b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_1(Mat,Vec,Vec);
95249b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_2(Mat,Vec,Vec);
95349b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_3(Mat,Vec,Vec);
95449b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_4(Mat,Vec,Vec);
95549b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_5(Mat,Vec,Vec);
95649b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_6(Mat,Vec,Vec);
95749b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_7(Mat,Vec,Vec);
95849b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_N(Mat,Vec,Vec);
95949b5e25fSSatish Balay 
96049b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_1(Mat,Vec,Vec,Vec);
96149b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_2(Mat,Vec,Vec,Vec);
96249b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_3(Mat,Vec,Vec,Vec);
96349b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_4(Mat,Vec,Vec,Vec);
96449b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_5(Mat,Vec,Vec,Vec);
96549b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_6(Mat,Vec,Vec,Vec);
96649b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_7(Mat,Vec,Vec,Vec);
96749b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_N(Mat,Vec,Vec,Vec);
96849b5e25fSSatish Balay 
9691a3463dfSHong Zhang #ifdef HAVE_MatIncompleteCholeskyFactor
9701a3463dfSHong Zhang /* modefied from MatILUFactor_SeqSBAIJ, needs further work!  */
97149b5e25fSSatish Balay #undef __FUNC__
9724ccecd49SHong Zhang #define __FUNC__ "MatIncompleteCholeskyFactor_SeqSBAIJ"
9731a3463dfSHong Zhang int MatIncompleteCholeskyFactor_SeqSBAIJ(Mat inA,IS row,PetscReal fill,int level)
97449b5e25fSSatish Balay {
9754ccecd49SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)inA->data;
97649b5e25fSSatish Balay   Mat         outA;
97749b5e25fSSatish Balay   int         ierr;
97849b5e25fSSatish Balay   PetscTruth  row_identity,col_identity;
97949b5e25fSSatish Balay 
98049b5e25fSSatish Balay   PetscFunctionBegin;
9811a3463dfSHong Zhang   /*
982*29bbc08cSBarry Smith   if (level != 0) SETERRQ(PETSC_ERR_SUP,"Only levels = 0 supported for in-place ILU");
98349b5e25fSSatish Balay   ierr = ISIdentity(row,&row_identity);CHKERRQ(ierr);
98449b5e25fSSatish Balay   ierr = ISIdentity(col,&col_identity);CHKERRQ(ierr);
98549b5e25fSSatish Balay   if (!row_identity || !col_identity) {
986*29bbc08cSBarry Smith     SETERRQ(1,"Row and column permutations must be identity for in-place ICC");
98749b5e25fSSatish Balay   }
9881a3463dfSHong Zhang   */
98949b5e25fSSatish Balay 
99049b5e25fSSatish Balay   outA          = inA;
99149b5e25fSSatish Balay   inA->factor   = FACTOR_LU;
99249b5e25fSSatish Balay 
99349b5e25fSSatish Balay   if (!a->diag) {
9941a3463dfSHong Zhang     ierr = MatMarkDiagonal_SeqSBAIJ(inA);CHKERRQ(ierr);
99549b5e25fSSatish Balay   }
99649b5e25fSSatish Balay   /*
99749b5e25fSSatish Balay       Blocksize 2, 3, 4, 5, 6 and 7 have a special faster factorization/solver
99849b5e25fSSatish Balay       for ILU(0) factorization with natural ordering
99949b5e25fSSatish Balay   */
100049b5e25fSSatish Balay   switch (a->bs) {
100149b5e25fSSatish Balay   case 1:
10021a3463dfSHong Zhang     inA->ops->solvetranspose   = MatSolve_SeqSBAIJ_2_NaturalOrdering; /* 2? */
10031a3463dfSHong Zhang     PLogInfo(inA,"MatIncompleteCholeskyFactor_SeqSBAIJ:Using special in-place natural ordering solvetrans BS=1\n");
100449b5e25fSSatish Balay   case 2:
10051a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering;
10061a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_2_NaturalOrdering;
10071a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_2_NaturalOrdering;
10081a3463dfSHong Zhang     PLogInfo(inA,"MatIncompleteCholeskyFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=2\n");
100949b5e25fSSatish Balay     break;
101049b5e25fSSatish Balay   case 3:
10111a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering;
10121a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_3_NaturalOrdering;
10131a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_3_NaturalOrdering;
10141a3463dfSHong Zhang     PLogInfo(inA,"MatIncompleteCholeskyFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=3\n");
101549b5e25fSSatish Balay     break;
101649b5e25fSSatish Balay   case 4:
10171a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering;
10181a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_4_NaturalOrdering;
10191a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_4_NaturalOrdering;
10201a3463dfSHong Zhang     PLogInfo(inA,"MatIncompleteCholeskyFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=4\n");
102149b5e25fSSatish Balay     break;
102249b5e25fSSatish Balay   case 5:
10231a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering;
10241a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_5_NaturalOrdering;
10251a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_5_NaturalOrdering;
10261a3463dfSHong Zhang     PLogInfo(inA,"MatIncompleteCholeskyFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=5\n");
102749b5e25fSSatish Balay     break;
102849b5e25fSSatish Balay   case 6:
10291a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering;
10301a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_6_NaturalOrdering;
10311a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_6_NaturalOrdering;
10321a3463dfSHong Zhang     PLogInfo(inA,"MatIncompleteCholeskyFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=6\n");
103349b5e25fSSatish Balay     break;
103449b5e25fSSatish Balay   case 7:
10351a3463dfSHong Zhang     inA->ops->lufactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering;
10361a3463dfSHong Zhang     inA->ops->solvetranspose  = MatSolve_SeqSBAIJ_7_NaturalOrdering;
10371a3463dfSHong Zhang     inA->ops->solve           = MatSolve_SeqSBAIJ_7_NaturalOrdering;
10381a3463dfSHong Zhang     PLogInfo(inA,"MatIncompleteCholeskyFactor_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=7\n");
103949b5e25fSSatish Balay     break;
104049b5e25fSSatish Balay   default:
104149b5e25fSSatish Balay     a->row        = row;
10421a3463dfSHong Zhang     a->icol       = col;
104349b5e25fSSatish Balay     ierr          = PetscObjectReference((PetscObject)row);CHKERRQ(ierr);
104449b5e25fSSatish Balay     ierr          = PetscObjectReference((PetscObject)col);CHKERRQ(ierr);
104549b5e25fSSatish Balay 
104649b5e25fSSatish Balay     /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */
104749b5e25fSSatish Balay     ierr = ISInvertPermutation(col,PETSC_DECIDE, &(a->icol));CHKERRQ(ierr);
104849b5e25fSSatish Balay     PLogObjectParent(inA,a->icol);
104949b5e25fSSatish Balay 
105049b5e25fSSatish Balay     if (!a->solve_work) {
105149b5e25fSSatish Balay       a->solve_work = (Scalar*)PetscMalloc((a->m+a->bs)*sizeof(Scalar));CHKPTRQ(a->solve_work);
105249b5e25fSSatish Balay       PLogObjectMemory(inA,(a->m+a->bs)*sizeof(Scalar));
105349b5e25fSSatish Balay     }
105449b5e25fSSatish Balay   }
105549b5e25fSSatish Balay 
10561a3463dfSHong Zhang   ierr = MatCholeskyFactorNumeric(inA,&outA);CHKERRQ(ierr);
105749b5e25fSSatish Balay 
105849b5e25fSSatish Balay   PetscFunctionReturn(0);
105949b5e25fSSatish Balay }
1060950f1e5bSHong Zhang #endif
1061950f1e5bSHong Zhang 
106249b5e25fSSatish Balay #undef __FUNC__
106349b5e25fSSatish Balay #define __FUNC__ "MatPrintHelp_SeqSBAIJ"
106449b5e25fSSatish Balay int MatPrintHelp_SeqSBAIJ(Mat A)
106549b5e25fSSatish Balay {
106649b5e25fSSatish Balay   static PetscTruth called = PETSC_FALSE;
106749b5e25fSSatish Balay   MPI_Comm          comm = A->comm;
106849b5e25fSSatish Balay   int               ierr;
106949b5e25fSSatish Balay 
107049b5e25fSSatish Balay   PetscFunctionBegin;
107149b5e25fSSatish Balay   if (called) {PetscFunctionReturn(0);} else called = PETSC_TRUE;
107249b5e25fSSatish Balay   ierr = (*PetscHelpPrintf)(comm," Options for MATSEQSBAIJ and MATMPISBAIJ matrix formats (the defaults):\n");CHKERRQ(ierr);
107349b5e25fSSatish Balay   ierr = (*PetscHelpPrintf)(comm,"  -mat_block_size <block_size>\n");CHKERRQ(ierr);
107449b5e25fSSatish Balay   PetscFunctionReturn(0);
107549b5e25fSSatish Balay }
107649b5e25fSSatish Balay 
107749b5e25fSSatish Balay EXTERN_C_BEGIN
107849b5e25fSSatish Balay #undef __FUNC__
107949b5e25fSSatish Balay #define __FUNC__ "MatSeqSBAIJSetColumnIndices_SeqSBAIJ"
108049b5e25fSSatish Balay int MatSeqSBAIJSetColumnIndices_SeqSBAIJ(Mat mat,int *indices)
108149b5e25fSSatish Balay {
1082045c9aa0SHong Zhang   Mat_SeqSBAIJ *baij = (Mat_SeqSBAIJ *)mat->data;
108349b5e25fSSatish Balay   int         i,nz,n;
108449b5e25fSSatish Balay 
108549b5e25fSSatish Balay   PetscFunctionBegin;
1086045c9aa0SHong Zhang   nz = baij->s_maxnz;
108749b5e25fSSatish Balay   n  = baij->n;
108849b5e25fSSatish Balay   for (i=0; i<nz; i++) {
108949b5e25fSSatish Balay     baij->j[i] = indices[i];
109049b5e25fSSatish Balay   }
1091045c9aa0SHong Zhang   baij->s_nz = nz;
109249b5e25fSSatish Balay   for (i=0; i<n; i++) {
109349b5e25fSSatish Balay     baij->ilen[i] = baij->imax[i];
109449b5e25fSSatish Balay   }
109549b5e25fSSatish Balay 
109649b5e25fSSatish Balay   PetscFunctionReturn(0);
109749b5e25fSSatish Balay }
109849b5e25fSSatish Balay EXTERN_C_END
109949b5e25fSSatish Balay 
110049b5e25fSSatish Balay #undef __FUNC__
110149b5e25fSSatish Balay #define __FUNC__ "MatSeqSBAIJSetColumnIndices"
110249b5e25fSSatish Balay /*@
110319585528SSatish Balay     MatSeqSBAIJSetColumnIndices - Set the column indices for all the rows
110449b5e25fSSatish Balay        in the matrix.
110549b5e25fSSatish Balay 
110649b5e25fSSatish Balay   Input Parameters:
110719585528SSatish Balay +  mat     - the SeqSBAIJ matrix
110849b5e25fSSatish Balay -  indices - the column indices
110949b5e25fSSatish Balay 
111049b5e25fSSatish Balay   Level: advanced
111149b5e25fSSatish Balay 
111249b5e25fSSatish Balay   Notes:
111349b5e25fSSatish Balay     This can be called if you have precomputed the nonzero structure of the
111449b5e25fSSatish Balay   matrix and want to provide it to the matrix object to improve the performance
111549b5e25fSSatish Balay   of the MatSetValues() operation.
111649b5e25fSSatish Balay 
111749b5e25fSSatish Balay     You MUST have set the correct numbers of nonzeros per row in the call to
111819585528SSatish Balay   MatCreateSeqSBAIJ().
111949b5e25fSSatish Balay 
1120ab9f2c04SSatish Balay     MUST be called before any calls to MatSetValues()
112149b5e25fSSatish Balay 
1122ab9f2c04SSatish Balay .seealso: MatCreateSeqSBAIJ
112349b5e25fSSatish Balay @*/
112449b5e25fSSatish Balay int MatSeqSBAIJSetColumnIndices(Mat mat,int *indices)
112549b5e25fSSatish Balay {
112649b5e25fSSatish Balay   int ierr,(*f)(Mat,int *);
112749b5e25fSSatish Balay 
112849b5e25fSSatish Balay   PetscFunctionBegin;
112949b5e25fSSatish Balay   PetscValidHeaderSpecific(mat,MAT_COOKIE);
11304afc71dfSHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqSBAIJSetColumnIndices_C",(void **)&f);CHKERRQ(ierr);
113149b5e25fSSatish Balay   if (f) {
113249b5e25fSSatish Balay     ierr = (*f)(mat,indices);CHKERRQ(ierr);
113349b5e25fSSatish Balay   } else {
1134*29bbc08cSBarry Smith     SETERRQ(1,"Wrong type of matrix to set column indices");
113549b5e25fSSatish Balay   }
113649b5e25fSSatish Balay   PetscFunctionReturn(0);
113749b5e25fSSatish Balay }
113849b5e25fSSatish Balay 
113949b5e25fSSatish Balay /* -------------------------------------------------------------------*/
114049b5e25fSSatish Balay static struct _MatOps MatOps_Values = {MatSetValues_SeqSBAIJ,
114149b5e25fSSatish Balay        MatGetRow_SeqSBAIJ,
114249b5e25fSSatish Balay        MatRestoreRow_SeqSBAIJ,
114349b5e25fSSatish Balay        MatMult_SeqSBAIJ_N,
114449b5e25fSSatish Balay        MatMultAdd_SeqSBAIJ_N,
114549b5e25fSSatish Balay        MatMultTranspose_SeqSBAIJ,
114649b5e25fSSatish Balay        MatMultTransposeAdd_SeqSBAIJ,
114749b5e25fSSatish Balay        MatSolve_SeqSBAIJ_N,
114849b5e25fSSatish Balay        0,
114949b5e25fSSatish Balay        0,
115049b5e25fSSatish Balay        0,
115149b5e25fSSatish Balay        0,
11525f4ef2deSHong Zhang        MatCholeskyFactor_SeqSBAIJ,
115349b5e25fSSatish Balay        0,
115449b5e25fSSatish Balay        MatTranspose_SeqSBAIJ,
115549b5e25fSSatish Balay        MatGetInfo_SeqSBAIJ,
115649b5e25fSSatish Balay        MatEqual_SeqSBAIJ,
115749b5e25fSSatish Balay        MatGetDiagonal_SeqSBAIJ,
115849b5e25fSSatish Balay        MatDiagonalScale_SeqSBAIJ,
115949b5e25fSSatish Balay        MatNorm_SeqSBAIJ,
116049b5e25fSSatish Balay        0,
116149b5e25fSSatish Balay        MatAssemblyEnd_SeqSBAIJ,
116249b5e25fSSatish Balay        0,
116349b5e25fSSatish Balay        MatSetOption_SeqSBAIJ,
116449b5e25fSSatish Balay        MatZeroEntries_SeqSBAIJ,
116549b5e25fSSatish Balay        MatZeroRows_SeqSBAIJ,
116649b5e25fSSatish Balay        0,
116749b5e25fSSatish Balay        0,
11685f4ef2deSHong Zhang        MatCholeskyFactorSymbolic_SeqSBAIJ,
11695f4ef2deSHong Zhang        MatCholeskyFactorNumeric_SeqSBAIJ_N,
117049b5e25fSSatish Balay        MatGetSize_SeqSBAIJ,
11714ccecd49SHong Zhang        MatGetSize_SeqSBAIJ,
117249b5e25fSSatish Balay        MatGetOwnershipRange_SeqSBAIJ,
117349b5e25fSSatish Balay        0,
11745f4ef2deSHong Zhang        MatIncompleteCholeskyFactorSymbolic_SeqSBAIJ,
117549b5e25fSSatish Balay        0,
117649b5e25fSSatish Balay        0,
117749b5e25fSSatish Balay        MatDuplicate_SeqSBAIJ,
117849b5e25fSSatish Balay        0,
117949b5e25fSSatish Balay        0,
118049b5e25fSSatish Balay        0,
1181950f1e5bSHong Zhang        0,
118249b5e25fSSatish Balay        0,
118349b5e25fSSatish Balay        MatGetSubMatrices_SeqSBAIJ,
118449b5e25fSSatish Balay        MatIncreaseOverlap_SeqSBAIJ,
118549b5e25fSSatish Balay        MatGetValues_SeqSBAIJ,
118649b5e25fSSatish Balay        0,
118749b5e25fSSatish Balay        MatPrintHelp_SeqSBAIJ,
118849b5e25fSSatish Balay        MatScale_SeqSBAIJ,
118949b5e25fSSatish Balay        0,
119049b5e25fSSatish Balay        0,
119149b5e25fSSatish Balay        0,
119249b5e25fSSatish Balay        MatGetBlockSize_SeqSBAIJ,
119349b5e25fSSatish Balay        MatGetRowIJ_SeqSBAIJ,
119449b5e25fSSatish Balay        MatRestoreRowIJ_SeqSBAIJ,
119549b5e25fSSatish Balay        0,
119649b5e25fSSatish Balay        0,
119749b5e25fSSatish Balay        0,
119849b5e25fSSatish Balay        0,
119949b5e25fSSatish Balay        0,
120049b5e25fSSatish Balay        0,
120149b5e25fSSatish Balay        MatSetValuesBlocked_SeqSBAIJ,
120249b5e25fSSatish Balay        MatGetSubMatrix_SeqSBAIJ,
120349b5e25fSSatish Balay        0,
120449b5e25fSSatish Balay        0,
120549b5e25fSSatish Balay        MatGetMaps_Petsc};
120649b5e25fSSatish Balay 
120749b5e25fSSatish Balay EXTERN_C_BEGIN
120849b5e25fSSatish Balay #undef __FUNC__
120949b5e25fSSatish Balay #define __FUNC__ "MatStoreValues_SeqSBAIJ"
121049b5e25fSSatish Balay int MatStoreValues_SeqSBAIJ(Mat mat)
121149b5e25fSSatish Balay {
12124afc71dfSHong Zhang   Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data;
121349b5e25fSSatish Balay   int         nz = aij->i[aij->m]*aij->bs*aij->bs2;
121449b5e25fSSatish Balay   int         ierr;
121549b5e25fSSatish Balay 
121649b5e25fSSatish Balay   PetscFunctionBegin;
121749b5e25fSSatish Balay   if (aij->nonew != 1) {
1218*29bbc08cSBarry Smith     SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
121949b5e25fSSatish Balay   }
122049b5e25fSSatish Balay 
122149b5e25fSSatish Balay   /* allocate space for values if not already there */
122249b5e25fSSatish Balay   if (!aij->saved_values) {
122349b5e25fSSatish Balay     aij->saved_values = (Scalar*)PetscMalloc(nz*sizeof(Scalar));CHKPTRQ(aij->saved_values);
122449b5e25fSSatish Balay   }
122549b5e25fSSatish Balay 
122649b5e25fSSatish Balay   /* copy values over */
122749b5e25fSSatish Balay   ierr = PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(Scalar));CHKERRQ(ierr);
122849b5e25fSSatish Balay   PetscFunctionReturn(0);
122949b5e25fSSatish Balay }
123049b5e25fSSatish Balay EXTERN_C_END
123149b5e25fSSatish Balay 
123249b5e25fSSatish Balay EXTERN_C_BEGIN
123349b5e25fSSatish Balay #undef __FUNC__
123449b5e25fSSatish Balay #define __FUNC__ "MatRetrieveValues_SeqSBAIJ"
123549b5e25fSSatish Balay int MatRetrieveValues_SeqSBAIJ(Mat mat)
123649b5e25fSSatish Balay {
12374afc71dfSHong Zhang   Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data;
123849b5e25fSSatish Balay   int         nz = aij->i[aij->m]*aij->bs*aij->bs2,ierr;
123949b5e25fSSatish Balay 
124049b5e25fSSatish Balay   PetscFunctionBegin;
124149b5e25fSSatish Balay   if (aij->nonew != 1) {
1242*29bbc08cSBarry Smith     SETERRQ(1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
124349b5e25fSSatish Balay   }
124449b5e25fSSatish Balay   if (!aij->saved_values) {
1245*29bbc08cSBarry Smith     SETERRQ(1,"Must call MatStoreValues(A);first");
124649b5e25fSSatish Balay   }
124749b5e25fSSatish Balay 
124849b5e25fSSatish Balay   /* copy values over */
124949b5e25fSSatish Balay   ierr = PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(Scalar));CHKERRQ(ierr);
125049b5e25fSSatish Balay   PetscFunctionReturn(0);
125149b5e25fSSatish Balay }
125249b5e25fSSatish Balay EXTERN_C_END
125349b5e25fSSatish Balay 
125449b5e25fSSatish Balay #undef __FUNC__
125549b5e25fSSatish Balay #define __FUNC__ "MatCreateSeqSBAIJ"
125649b5e25fSSatish Balay /*@C
125749b5e25fSSatish Balay    MatCreateSeqSBAIJ - Creates a sparse symmetric matrix in block AIJ (block
125849b5e25fSSatish Balay    compressed row) format.  For good matrix assembly performance the
125949b5e25fSSatish Balay    user should preallocate the matrix storage by setting the parameter nz
126049b5e25fSSatish Balay    (or the array nnz).  By setting these parameters accurately, performance
126149b5e25fSSatish Balay    during matrix assembly can be increased by more than a factor of 50.
126249b5e25fSSatish Balay 
126349b5e25fSSatish Balay    Collective on MPI_Comm
126449b5e25fSSatish Balay 
126549b5e25fSSatish Balay    Input Parameters:
126649b5e25fSSatish Balay +  comm - MPI communicator, set to PETSC_COMM_SELF
126749b5e25fSSatish Balay .  bs - size of block
126849b5e25fSSatish Balay .  m - number of rows, or number of columns
126949b5e25fSSatish Balay .  nz - number of block nonzeros per block row (same for all rows)
127049b5e25fSSatish Balay -  nnz - array containing the number of block nonzeros in the various block rows
127149b5e25fSSatish Balay          (possibly different for each block row) or PETSC_NULL
127249b5e25fSSatish Balay 
127349b5e25fSSatish Balay    Output Parameter:
127449b5e25fSSatish Balay .  A - the symmetric matrix
127549b5e25fSSatish Balay 
127649b5e25fSSatish Balay    Options Database Keys:
127749b5e25fSSatish Balay .   -mat_no_unroll - uses code that does not unroll the loops in the
127849b5e25fSSatish Balay                      block calculations (much slower)
127949b5e25fSSatish Balay .    -mat_block_size - size of the blocks to use
128049b5e25fSSatish Balay 
128149b5e25fSSatish Balay    Level: intermediate
128249b5e25fSSatish Balay 
128349b5e25fSSatish Balay    Notes:
128449b5e25fSSatish Balay    The block AIJ format is fully compatible with standard Fortran 77
128549b5e25fSSatish Balay    storage.  That is, the stored row and column indices can begin at
128649b5e25fSSatish Balay    either one (as in Fortran) or zero.  See the users' manual for details.
128749b5e25fSSatish Balay 
128849b5e25fSSatish Balay    Specify the preallocated storage with either nz or nnz (not both).
128949b5e25fSSatish Balay    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
129049b5e25fSSatish Balay    allocation.  For additional details, see the users manual chapter on
129149b5e25fSSatish Balay    matrices.
129249b5e25fSSatish Balay 
1293a209d233SLois Curfman McInnes .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPISBAIJ()
129449b5e25fSSatish Balay @*/
129549b5e25fSSatish Balay int MatCreateSeqSBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,int *nnz,Mat *A)
129649b5e25fSSatish Balay {
129749b5e25fSSatish Balay   Mat         B;
129849b5e25fSSatish Balay   Mat_SeqSBAIJ *b;
12994afc71dfSHong Zhang   int         i,len,ierr,mbs,bs2,size;
130049b5e25fSSatish Balay   PetscTruth  flg;
13014afc71dfSHong Zhang   int         s_nz;
130249b5e25fSSatish Balay 
130349b5e25fSSatish Balay   PetscFunctionBegin;
130449b5e25fSSatish Balay   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1305*29bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1");
130649b5e25fSSatish Balay 
130749b5e25fSSatish Balay   ierr = OptionsGetInt(PETSC_NULL,"-mat_block_size",&bs,PETSC_NULL);CHKERRQ(ierr);
130849b5e25fSSatish Balay   mbs  = m/bs;
130949b5e25fSSatish Balay   bs2  = bs*bs;
131049b5e25fSSatish Balay 
131149b5e25fSSatish Balay   if (mbs*bs!=m) {
1312*29bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_SIZ,"Number rows, cols must be divisible by blocksize");
131349b5e25fSSatish Balay   }
131449b5e25fSSatish Balay 
1315*29bbc08cSBarry Smith   if (nz < -2) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than -2: value %d",nz);
131649b5e25fSSatish Balay   if (nnz) {
131749b5e25fSSatish Balay     for (i=0; i<mbs; i++) {
1318*29bbc08cSBarry Smith       if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %d value %d",i,nnz[i]);
131949b5e25fSSatish Balay     }
132049b5e25fSSatish Balay   }
132149b5e25fSSatish Balay 
132249b5e25fSSatish Balay   *A = 0;
132349b5e25fSSatish Balay   PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQSBAIJ,"Mat",comm,MatDestroy,MatView);
132449b5e25fSSatish Balay   PLogObjectCreate(B);
132549b5e25fSSatish Balay   B->data = (void*)(b = PetscNew(Mat_SeqSBAIJ));CHKPTRQ(b);
132649b5e25fSSatish Balay   ierr    = PetscMemzero(b,sizeof(Mat_SeqSBAIJ));CHKERRQ(ierr);
132749b5e25fSSatish Balay   ierr    = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
132849b5e25fSSatish Balay   ierr    = OptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg);CHKERRQ(ierr);
132949b5e25fSSatish Balay   if (!flg) {
133049b5e25fSSatish Balay     switch (bs) {
133149b5e25fSSatish Balay     case 1:
13324ccecd49SHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1;
133349b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_1;
133449b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_1;
133549b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_1;
133649b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_1;
133749b5e25fSSatish Balay       break;
133849b5e25fSSatish Balay     case 2:
13394ccecd49SHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2;
134049b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_2;
134149b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_2;
134249b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_2;
134349b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_2;
134449b5e25fSSatish Balay       break;
134549b5e25fSSatish Balay     case 3:
13465f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3;
134749b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_3;
134849b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_3;
134949b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_3;
135049b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_3;
135149b5e25fSSatish Balay       break;
135249b5e25fSSatish Balay     case 4:
13535f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4;
135449b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_4;
135549b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_4;
135649b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_4;
135749b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_4;
135849b5e25fSSatish Balay       break;
135949b5e25fSSatish Balay     case 5:
13605f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5;
136149b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_5;
136249b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_5;
136349b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_5;
136449b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_5;
136549b5e25fSSatish Balay       break;
136649b5e25fSSatish Balay     case 6:
13675f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6;
136849b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_6;
136949b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_6;
137049b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_6;
137149b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_6;
137249b5e25fSSatish Balay       break;
137349b5e25fSSatish Balay     case 7:
137449b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_7;
137549b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_7;
137649b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_7;
137749b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_7;
137849b5e25fSSatish Balay       break;
137949b5e25fSSatish Balay     }
138049b5e25fSSatish Balay   }
138149b5e25fSSatish Balay   B->ops->destroy     = MatDestroy_SeqSBAIJ;
138249b5e25fSSatish Balay   B->ops->view        = MatView_SeqSBAIJ;
138349b5e25fSSatish Balay   B->factor           = 0;
138449b5e25fSSatish Balay   B->lupivotthreshold = 1.0;
138549b5e25fSSatish Balay   B->mapping          = 0;
138649b5e25fSSatish Balay   b->row              = 0;
138749b5e25fSSatish Balay   b->icol             = 0;
138849b5e25fSSatish Balay   b->reallocs         = 0;
138949b5e25fSSatish Balay   b->saved_values     = 0;
139049b5e25fSSatish Balay 
13915f4ef2deSHong Zhang   b->m       = m;
13925f4ef2deSHong Zhang   B->m = m; B->M = m;
13935f4ef2deSHong Zhang   b->n       = m;
139449b5e25fSSatish Balay   B->n = m; B->N = m;
139549b5e25fSSatish Balay 
139649b5e25fSSatish Balay   ierr = MapCreateMPI(comm,m,m,&B->rmap);CHKERRQ(ierr);
139749b5e25fSSatish Balay   ierr = MapCreateMPI(comm,m,m,&B->cmap);CHKERRQ(ierr);
139849b5e25fSSatish Balay 
139949b5e25fSSatish Balay   b->mbs     = mbs;
14004afc71dfSHong Zhang   b->nbs     = mbs;
140149b5e25fSSatish Balay   b->imax    = (int*)PetscMalloc((mbs+1)*sizeof(int));CHKPTRQ(b->imax);
140249b5e25fSSatish Balay   if (!nnz) {
140349b5e25fSSatish Balay     if (nz == PETSC_DEFAULT) nz = 5;
140449b5e25fSSatish Balay     else if (nz <= 0)        nz = 1;
140549b5e25fSSatish Balay     for (i=0; i<mbs; i++) {
14068cef66ccSHong Zhang       b->imax[i] = nz;
140749b5e25fSSatish Balay     }
1408153ea458SHong Zhang     nz = nz*mbs; /* total nz */
140949b5e25fSSatish Balay   } else {
141049b5e25fSSatish Balay     nz = 0;
14118cef66ccSHong Zhang     for (i=0; i<mbs; i++) {b->imax[i] = nnz[i]; nz += nnz[i];}
141249b5e25fSSatish Balay   }
14138cef66ccSHong Zhang   /* s_nz=(nz+mbs)/2; */ /* total diagonal and superdiagonal nonzero blocks */
14148cef66ccSHong Zhang   s_nz = nz;
141549b5e25fSSatish Balay 
141649b5e25fSSatish Balay   /* allocate the matrix space */
141749b5e25fSSatish Balay   len   = s_nz*sizeof(int) + s_nz*bs2*sizeof(MatScalar) + (b->m+1)*sizeof(int);
141849b5e25fSSatish Balay   b->a  = (MatScalar*)PetscMalloc(len);CHKPTRQ(b->a);
141949b5e25fSSatish Balay   ierr = PetscMemzero(b->a,s_nz*bs2*sizeof(MatScalar));CHKERRQ(ierr);
142049b5e25fSSatish Balay   b->j  = (int*)(b->a + s_nz*bs2);
142149b5e25fSSatish Balay   ierr = PetscMemzero(b->j,s_nz*sizeof(int));CHKERRQ(ierr);
142249b5e25fSSatish Balay   b->i  = b->j + s_nz;
142349b5e25fSSatish Balay   b->singlemalloc = PETSC_TRUE;
142449b5e25fSSatish Balay 
142549b5e25fSSatish Balay   /* pointer to beginning of each row */
142649b5e25fSSatish Balay   b->i[0] = 0;
142749b5e25fSSatish Balay   for (i=1; i<mbs+1; i++) {
142849b5e25fSSatish Balay     b->i[i] = b->i[i-1] + b->imax[i-1];
142949b5e25fSSatish Balay   }
143049b5e25fSSatish Balay 
143149b5e25fSSatish Balay   /* b->ilen will count nonzeros in each block row so far. */
143249b5e25fSSatish Balay   b->ilen = (int*)PetscMalloc((mbs+1)*sizeof(int));
143349b5e25fSSatish Balay   PLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqSBAIJ));
143449b5e25fSSatish Balay   for (i=0; i<mbs; i++) { b->ilen[i] = 0;}
143549b5e25fSSatish Balay 
143649b5e25fSSatish Balay   b->bs               = bs;
143749b5e25fSSatish Balay   b->bs2              = bs2;
14381ce8b1feSHong Zhang   /* b->mbs              = mbs; redundant */
143949b5e25fSSatish Balay   b->s_nz               = 0;
144049b5e25fSSatish Balay   b->s_maxnz            = s_nz*bs2;
144149b5e25fSSatish Balay   b->sorted           = PETSC_FALSE;
144249b5e25fSSatish Balay   b->roworiented      = PETSC_TRUE;
144349b5e25fSSatish Balay   b->nonew            = 0;
144449b5e25fSSatish Balay   b->diag             = 0;
144549b5e25fSSatish Balay   b->solve_work       = 0;
144649b5e25fSSatish Balay   b->mult_work        = 0;
144749b5e25fSSatish Balay   b->spptr            = 0;
144849b5e25fSSatish Balay   B->info.nz_unneeded = (PetscReal)b->s_maxnz;
144949b5e25fSSatish Balay   b->keepzeroedrows   = PETSC_FALSE;
1450153ea458SHong Zhang 
145116cdd363SHong Zhang   b->inew             = 0;
145216cdd363SHong Zhang   b->jnew             = 0;
145316cdd363SHong Zhang   b->anew             = 0;
145416cdd363SHong Zhang   b->a2anew           = 0;
14551a3463dfSHong Zhang   b->permute          = PETSC_FALSE;
1456153ea458SHong Zhang 
145749b5e25fSSatish Balay   *A = B;
145849b5e25fSSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-help",&flg);CHKERRQ(ierr);
145949b5e25fSSatish Balay   if (flg) {ierr = MatPrintHelp(B);CHKERRQ(ierr); }
146049b5e25fSSatish Balay 
146149b5e25fSSatish Balay   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C",
1462aec82049SSatish Balay                                      "MatStoreValues_SeqSBAIJ",
1463aec82049SSatish Balay                                      (void*)MatStoreValues_SeqSBAIJ);CHKERRQ(ierr);
146449b5e25fSSatish Balay   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C",
1465aec82049SSatish Balay                                      "MatRetrieveValues_SeqSBAIJ",
1466aec82049SSatish Balay                                      (void*)MatRetrieveValues_SeqSBAIJ);CHKERRQ(ierr);
1467aec82049SSatish Balay   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqSBAIJSetColumnIndices_C",
1468aec82049SSatish Balay                                      "MatSeqSBAIJSetColumnIndices_SeqSBAIJ",
14691ce8b1feSHong Zhang                                      (void*)MatSeqSBAIJSetColumnIndices_SeqSBAIJ);CHKERRQ(ierr);
14701ce8b1feSHong Zhang   /* printf("In MatCreateSeqSBAIJ, s_nz = %d, bs2=%d, m=%d, mbs=%d \n", s_nz,bs2,m,mbs); */
1471153ea458SHong Zhang    /*
147249b5e25fSSatish Balay    for (i=0; i<mbs; i++){
147349b5e25fSSatish Balay      printf("imax[%d]= %d, ilen[%d]= %d\n", i,b->imax[i], i,b->ilen[i]);
147449b5e25fSSatish Balay    }       */
147549b5e25fSSatish Balay 
147649b5e25fSSatish Balay   PetscFunctionReturn(0);
147749b5e25fSSatish Balay }
147849b5e25fSSatish Balay 
147949b5e25fSSatish Balay #undef __FUNC__
148049b5e25fSSatish Balay #define __FUNC__ "MatDuplicate_SeqSBAIJ"
148149b5e25fSSatish Balay int MatDuplicate_SeqSBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B)
148249b5e25fSSatish Balay {
148349b5e25fSSatish Balay   Mat         C;
148449b5e25fSSatish Balay   Mat_SeqSBAIJ *c,*a = (Mat_SeqSBAIJ*)A->data;
148549b5e25fSSatish Balay   int         i,len,mbs = a->mbs,nz = a->s_nz,bs2 =a->bs2,ierr;
148649b5e25fSSatish Balay 
148749b5e25fSSatish Balay   PetscFunctionBegin;
1488*29bbc08cSBarry Smith   if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,"Corrupt matrix");
148949b5e25fSSatish Balay 
149049b5e25fSSatish Balay   *B = 0;
149149b5e25fSSatish Balay   PetscHeaderCreate(C,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQSBAIJ,"Mat",A->comm,MatDestroy,MatView);
149249b5e25fSSatish Balay   PLogObjectCreate(C);
149349b5e25fSSatish Balay   C->data           = (void*)(c = PetscNew(Mat_SeqSBAIJ));CHKPTRQ(c);
149449b5e25fSSatish Balay   ierr              = PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
149549b5e25fSSatish Balay   C->ops->destroy   = MatDestroy_SeqSBAIJ;
149649b5e25fSSatish Balay   C->ops->view      = MatView_SeqSBAIJ;
149749b5e25fSSatish Balay   C->factor         = A->factor;
149849b5e25fSSatish Balay   c->row            = 0;
149949b5e25fSSatish Balay   c->icol           = 0;
150049b5e25fSSatish Balay   c->saved_values   = 0;
150149b5e25fSSatish Balay   c->keepzeroedrows = a->keepzeroedrows;
150249b5e25fSSatish Balay   C->assembled      = PETSC_TRUE;
150349b5e25fSSatish Balay 
150449b5e25fSSatish Balay   c->m = C->m   = a->m;
150549b5e25fSSatish Balay   c->n = C->n   = a->n;
150649b5e25fSSatish Balay   C->M          = a->m;
150749b5e25fSSatish Balay   C->N          = a->n;
150849b5e25fSSatish Balay 
150949b5e25fSSatish Balay   c->bs         = a->bs;
151049b5e25fSSatish Balay   c->bs2        = a->bs2;
151149b5e25fSSatish Balay   c->mbs        = a->mbs;
151249b5e25fSSatish Balay   c->nbs        = a->nbs;
151349b5e25fSSatish Balay 
151449b5e25fSSatish Balay   c->imax       = (int*)PetscMalloc((mbs+1)*sizeof(int));CHKPTRQ(c->imax);
151549b5e25fSSatish Balay   c->ilen       = (int*)PetscMalloc((mbs+1)*sizeof(int));CHKPTRQ(c->ilen);
151649b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
151749b5e25fSSatish Balay     c->imax[i] = a->imax[i];
151849b5e25fSSatish Balay     c->ilen[i] = a->ilen[i];
151949b5e25fSSatish Balay   }
152049b5e25fSSatish Balay 
152149b5e25fSSatish Balay   /* allocate the matrix space */
152249b5e25fSSatish Balay   c->singlemalloc = PETSC_TRUE;
152349b5e25fSSatish Balay   len  = (mbs+1)*sizeof(int) + nz*(bs2*sizeof(MatScalar) + sizeof(int));
152449b5e25fSSatish Balay   c->a = (MatScalar*)PetscMalloc(len);CHKPTRQ(c->a);
152549b5e25fSSatish Balay   c->j = (int*)(c->a + nz*bs2);
152649b5e25fSSatish Balay   c->i = c->j + nz;
152749b5e25fSSatish Balay   ierr = PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int));CHKERRQ(ierr);
152849b5e25fSSatish Balay   if (mbs > 0) {
152949b5e25fSSatish Balay     ierr = PetscMemcpy(c->j,a->j,nz*sizeof(int));CHKERRQ(ierr);
153049b5e25fSSatish Balay     if (cpvalues == MAT_COPY_VALUES) {
153149b5e25fSSatish Balay       ierr = PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr);
153249b5e25fSSatish Balay     } else {
153349b5e25fSSatish Balay       ierr = PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr);
153449b5e25fSSatish Balay     }
153549b5e25fSSatish Balay   }
153649b5e25fSSatish Balay 
153749b5e25fSSatish Balay   PLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqSBAIJ));
153849b5e25fSSatish Balay   c->sorted      = a->sorted;
153949b5e25fSSatish Balay   c->roworiented = a->roworiented;
154049b5e25fSSatish Balay   c->nonew       = a->nonew;
154149b5e25fSSatish Balay 
154249b5e25fSSatish Balay   if (a->diag) {
154349b5e25fSSatish Balay     c->diag = (int*)PetscMalloc((mbs+1)*sizeof(int));CHKPTRQ(c->diag);
154449b5e25fSSatish Balay     PLogObjectMemory(C,(mbs+1)*sizeof(int));
154549b5e25fSSatish Balay     for (i=0; i<mbs; i++) {
154649b5e25fSSatish Balay       c->diag[i] = a->diag[i];
154749b5e25fSSatish Balay     }
154849b5e25fSSatish Balay   } else c->diag        = 0;
154949b5e25fSSatish Balay   c->s_nz               = a->s_nz;
155049b5e25fSSatish Balay   c->s_maxnz            = a->s_maxnz;
155149b5e25fSSatish Balay   c->solve_work         = 0;
155249b5e25fSSatish Balay   c->spptr              = 0;      /* Dangerous -I'm throwing away a->spptr */
155349b5e25fSSatish Balay   c->mult_work          = 0;
155449b5e25fSSatish Balay   *B = C;
155549b5e25fSSatish Balay   ierr = FListDuplicate(A->qlist,&C->qlist);CHKERRQ(ierr);
155649b5e25fSSatish Balay   PetscFunctionReturn(0);
155749b5e25fSSatish Balay }
155849b5e25fSSatish Balay 
155949b5e25fSSatish Balay #undef __FUNC__
156049b5e25fSSatish Balay #define __FUNC__ "MatLoad_SeqSBAIJ"
156149b5e25fSSatish Balay int MatLoad_SeqSBAIJ(Viewer viewer,MatType type,Mat *A)
156249b5e25fSSatish Balay {
156349b5e25fSSatish Balay   Mat_SeqSBAIJ  *a;
156449b5e25fSSatish Balay   Mat          B;
156549b5e25fSSatish Balay   int          i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1;
1566574b2666SHong Zhang   int          *mask,mbs,*jj,j,rowcount,nzcount,k,*browlengths,*s_browlengths,maskcount;
156749b5e25fSSatish Balay   int          kmax,jcount,block,idx,point,nzcountb,extra_rows;
156849b5e25fSSatish Balay   int          *masked,nmask,tmp,bs2,ishift;
156949b5e25fSSatish Balay   Scalar       *aa;
157049b5e25fSSatish Balay   MPI_Comm     comm = ((PetscObject)viewer)->comm;
157149b5e25fSSatish Balay 
157249b5e25fSSatish Balay   PetscFunctionBegin;
157349b5e25fSSatish Balay   ierr = OptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,PETSC_NULL);CHKERRQ(ierr);
157449b5e25fSSatish Balay   bs2  = bs*bs;
157549b5e25fSSatish Balay 
157649b5e25fSSatish Balay   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1577*29bbc08cSBarry Smith   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor");
157849b5e25fSSatish Balay   ierr = ViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
157949b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr);
1580*29bbc08cSBarry Smith   if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not Mat object");
158149b5e25fSSatish Balay   M = header[1]; N = header[2]; nz = header[3];
158249b5e25fSSatish Balay 
158349b5e25fSSatish Balay   if (header[3] < 0) {
1584*29bbc08cSBarry Smith     SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format, cannot load as SeqSBAIJ");
158549b5e25fSSatish Balay   }
158649b5e25fSSatish Balay 
1587*29bbc08cSBarry Smith   if (M != N) SETERRQ(PETSC_ERR_SUP,"Can only do square matrices");
158849b5e25fSSatish Balay 
158949b5e25fSSatish Balay   /*
159049b5e25fSSatish Balay      This code adds extra rows to make sure the number of rows is
159149b5e25fSSatish Balay     divisible by the blocksize
159249b5e25fSSatish Balay   */
159349b5e25fSSatish Balay   mbs        = M/bs;
159449b5e25fSSatish Balay   extra_rows = bs - M + bs*(mbs);
159549b5e25fSSatish Balay   if (extra_rows == bs) extra_rows = 0;
159649b5e25fSSatish Balay   else                  mbs++;
159749b5e25fSSatish Balay   if (extra_rows) {
159849b5e25fSSatish Balay     PLogInfo(0,"MatLoad_SeqSBAIJ:Padding loaded matrix to match blocksize\n");
159949b5e25fSSatish Balay   }
160049b5e25fSSatish Balay 
160149b5e25fSSatish Balay   /* read in row lengths */
160249b5e25fSSatish Balay   rowlengths = (int*)PetscMalloc((M+extra_rows)*sizeof(int));CHKPTRQ(rowlengths);
160349b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
160449b5e25fSSatish Balay   for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1;
160549b5e25fSSatish Balay 
160649b5e25fSSatish Balay   /* read in column indices */
160749b5e25fSSatish Balay   jj = (int*)PetscMalloc((nz+extra_rows)*sizeof(int));CHKPTRQ(jj);
160849b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT);CHKERRQ(ierr);
160949b5e25fSSatish Balay   for (i=0; i<extra_rows; i++) jj[nz+i] = M+i;
161049b5e25fSSatish Balay 
161149b5e25fSSatish Balay   /* loop over row lengths determining block row lengths */
161249b5e25fSSatish Balay   browlengths = (int*)PetscMalloc(mbs*sizeof(int));CHKPTRQ(browlengths);
1613574b2666SHong Zhang   s_browlengths = (int*)PetscMalloc(mbs*sizeof(int));CHKPTRQ(s_browlengths);
1614574b2666SHong Zhang   ierr        = PetscMemzero(s_browlengths,mbs*sizeof(int));CHKERRQ(ierr);
161549b5e25fSSatish Balay   mask        = (int*)PetscMalloc(2*mbs*sizeof(int));CHKPTRQ(mask);
161649b5e25fSSatish Balay   ierr        = PetscMemzero(mask,mbs*sizeof(int));CHKERRQ(ierr);
161749b5e25fSSatish Balay   masked      = mask + mbs;
161849b5e25fSSatish Balay   rowcount    = 0; nzcount = 0;
161949b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
162049b5e25fSSatish Balay     nmask = 0;
162149b5e25fSSatish Balay     for (j=0; j<bs; j++) {
162249b5e25fSSatish Balay       kmax = rowlengths[rowcount];
162349b5e25fSSatish Balay       for (k=0; k<kmax; k++) {
16242d703238SHong Zhang         tmp = jj[nzcount++]/bs;   /* block col. index */
162503630b6eSHong Zhang         if (!mask[tmp] && tmp >= i) {masked[nmask++] = tmp; mask[tmp] = 1;}
162649b5e25fSSatish Balay       }
162749b5e25fSSatish Balay       rowcount++;
162849b5e25fSSatish Balay     }
1629574b2666SHong Zhang     s_browlengths[i] += nmask;
1630574b2666SHong Zhang     browlengths[i]    = 2*s_browlengths[i];
1631574b2666SHong Zhang 
163249b5e25fSSatish Balay     /* zero out the mask elements we set */
163349b5e25fSSatish Balay     for (j=0; j<nmask; j++) mask[masked[j]] = 0;
163449b5e25fSSatish Balay   }
163549b5e25fSSatish Balay 
163649b5e25fSSatish Balay   /* create our matrix */
163749b5e25fSSatish Balay   ierr = MatCreateSeqSBAIJ(comm,bs,M+extra_rows,N+extra_rows,0,browlengths,A);CHKERRQ(ierr);
163849b5e25fSSatish Balay   B = *A;
163949b5e25fSSatish Balay   a = (Mat_SeqSBAIJ*)B->data;
164049b5e25fSSatish Balay 
164149b5e25fSSatish Balay   /* set matrix "i" values */
164249b5e25fSSatish Balay   a->i[0] = 0;
164349b5e25fSSatish Balay   for (i=1; i<= mbs; i++) {
1644574b2666SHong Zhang     a->i[i]      = a->i[i-1] + s_browlengths[i-1];
1645574b2666SHong Zhang     a->ilen[i-1] = s_browlengths[i-1];
164649b5e25fSSatish Balay   }
164749b5e25fSSatish Balay   a->s_nz         = 0;
1648574b2666SHong Zhang   for (i=0; i<mbs; i++) a->s_nz += s_browlengths[i];
164949b5e25fSSatish Balay 
165049b5e25fSSatish Balay   /* read in nonzero values */
165149b5e25fSSatish Balay   aa = (Scalar*)PetscMalloc((nz+extra_rows)*sizeof(Scalar));CHKPTRQ(aa);
165249b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);CHKERRQ(ierr);
165349b5e25fSSatish Balay   for (i=0; i<extra_rows; i++) aa[nz+i] = 1.0;
165449b5e25fSSatish Balay 
165549b5e25fSSatish Balay   /* set "a" and "j" values into matrix */
165649b5e25fSSatish Balay   nzcount = 0; jcount = 0;
165749b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
165849b5e25fSSatish Balay     nzcountb = nzcount;
165949b5e25fSSatish Balay     nmask    = 0;
166049b5e25fSSatish Balay     for (j=0; j<bs; j++) {
166149b5e25fSSatish Balay       kmax = rowlengths[i*bs+j];
166249b5e25fSSatish Balay       for (k=0; k<kmax; k++) {
16632d703238SHong Zhang         tmp = jj[nzcount++]/bs; /* block col. index */
166403630b6eSHong Zhang         if (!mask[tmp] && tmp >= i) { masked[nmask++] = tmp; mask[tmp] = 1;}
16652d703238SHong Zhang       }
16662d703238SHong Zhang       rowcount++;
16672d703238SHong Zhang     }
16682d703238SHong Zhang     /* sort the masked values */
16692d703238SHong Zhang     ierr = PetscSortInt(nmask,masked);CHKERRQ(ierr);
16702d703238SHong Zhang 
16712d703238SHong Zhang     /* set "j" values into matrix */
16722d703238SHong Zhang     maskcount = 1;
16732d703238SHong Zhang     for (j=0; j<nmask; j++) {
167449b5e25fSSatish Balay       a->j[jcount++]  = masked[j];
167549b5e25fSSatish Balay       mask[masked[j]] = maskcount++;
167649b5e25fSSatish Balay     }
1677574b2666SHong Zhang 
167849b5e25fSSatish Balay     /* set "a" values into matrix */
167949b5e25fSSatish Balay     ishift = bs2*a->i[i];
168049b5e25fSSatish Balay     for (j=0; j<bs; j++) {
168149b5e25fSSatish Balay       kmax = rowlengths[i*bs+j];
168249b5e25fSSatish Balay       for (k=0; k<kmax; k++) {
1683574b2666SHong Zhang         tmp       = jj[nzcountb]/bs ; /* block col. index */
1684574b2666SHong Zhang         if (tmp >= i){
168549b5e25fSSatish Balay           block     = mask[tmp] - 1;
168649b5e25fSSatish Balay           point     = jj[nzcountb] - bs*tmp;
168749b5e25fSSatish Balay           idx       = ishift + bs2*block + j + bs*point;
1688574b2666SHong Zhang           a->a[idx] = aa[nzcountb];
1689574b2666SHong Zhang         }
1690574b2666SHong Zhang         nzcountb++;
169149b5e25fSSatish Balay       }
169249b5e25fSSatish Balay     }
169349b5e25fSSatish Balay     /* zero out the mask elements we set */
169449b5e25fSSatish Balay     for (j=0; j<nmask; j++) mask[masked[j]] = 0;
169549b5e25fSSatish Balay   }
1696*29bbc08cSBarry Smith   if (jcount != a->s_nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Bad binary matrix");
169749b5e25fSSatish Balay 
169849b5e25fSSatish Balay   ierr = PetscFree(rowlengths);CHKERRQ(ierr);
169949b5e25fSSatish Balay   ierr = PetscFree(browlengths);CHKERRQ(ierr);
1700574b2666SHong Zhang   ierr = PetscFree(s_browlengths);CHKERRQ(ierr);
170149b5e25fSSatish Balay   ierr = PetscFree(aa);CHKERRQ(ierr);
170249b5e25fSSatish Balay   ierr = PetscFree(jj);CHKERRQ(ierr);
170349b5e25fSSatish Balay   ierr = PetscFree(mask);CHKERRQ(ierr);
170449b5e25fSSatish Balay 
170549b5e25fSSatish Balay   B->assembled = PETSC_TRUE;
170649b5e25fSSatish Balay   ierr = MatView_Private(B);CHKERRQ(ierr);
170749b5e25fSSatish Balay   PetscFunctionReturn(0);
170849b5e25fSSatish Balay }
1709574b2666SHong Zhang 
1710574b2666SHong Zhang 
171149b5e25fSSatish Balay 
171249b5e25fSSatish Balay 
1713