xref: /petsc/src/mat/impls/sbaij/seq/sbaij.c (revision ab9f2c04529085dd459171d57e735bfa81c87bde)
1*ab9f2c04SSatish Balay /*$Id: sbaij.c,v 1.12 2000/07/28 21:14:06 hzhang Exp balay $*/
249b5e25fSSatish Balay 
349b5e25fSSatish Balay /*
449b5e25fSSatish Balay     Defines the basic matrix operations for the BAIJ (compressed row)
549b5e25fSSatish Balay   matrix storage format.
649b5e25fSSatish Balay */
7*ab9f2c04SSatish 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 {
2249b5e25fSSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
2349b5e25fSSatish Balay   int         *diag,*jj = a->j,i,ierr;
2449b5e25fSSatish Balay 
2549b5e25fSSatish Balay   PetscFunctionBegin;
2649b5e25fSSatish Balay   ierr = MatMarkDiagonal_SeqBAIJ(A);CHKERRQ(ierr);
2749b5e25fSSatish Balay   diag = a->diag;
2849b5e25fSSatish Balay   for (i=0; i<a->mbs; i++) {
2949b5e25fSSatish Balay     if (jj[diag[i]] != i) {
3049b5e25fSSatish Balay       SETERRQ1(1,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 {
4049b5e25fSSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)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   int         ierr,n = a->mbs,i;
7049b5e25fSSatish Balay 
7149b5e25fSSatish Balay   PetscFunctionBegin;
7249b5e25fSSatish Balay   *nn = n;
7349b5e25fSSatish Balay   if (!ia) PetscFunctionReturn(0);
7449b5e25fSSatish Balay   if (symmetric) {
7549b5e25fSSatish Balay     ierr = MatToSymmetricIJ_SeqAIJ(n,a->i,a->j,0,oshift,ia,ja);CHKERRQ(ierr);
7649b5e25fSSatish Balay   } else if (oshift == 1) {
7749b5e25fSSatish Balay     /* temporarily add 1 to i and j indices */
7849b5e25fSSatish Balay     int nz = a->i[n] + 1;
7949b5e25fSSatish Balay     for (i=0; i<nz; i++) a->j[i]++;
8049b5e25fSSatish Balay     for (i=0; i<n+1; i++) a->i[i]++;
8149b5e25fSSatish Balay     *ia = a->i; *ja = a->j;
8249b5e25fSSatish Balay   } else {
8349b5e25fSSatish Balay     *ia = a->i; *ja = a->j;
8449b5e25fSSatish Balay   }
8549b5e25fSSatish Balay 
8649b5e25fSSatish Balay   PetscFunctionReturn(0);
8749b5e25fSSatish Balay }
8849b5e25fSSatish Balay 
8949b5e25fSSatish Balay #undef __FUNC__
9049b5e25fSSatish Balay #define __FUNC__ "MatRestoreRowIJ_SeqSBAIJ"
9149b5e25fSSatish Balay static int MatRestoreRowIJ_SeqSBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja,
9249b5e25fSSatish Balay                                 PetscTruth *done)
9349b5e25fSSatish Balay {
9449b5e25fSSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
9549b5e25fSSatish Balay   int         i,n = a->mbs,ierr;
9649b5e25fSSatish Balay 
9749b5e25fSSatish Balay   PetscFunctionBegin;
9849b5e25fSSatish Balay   if (!ia) PetscFunctionReturn(0);
9949b5e25fSSatish Balay   if (symmetric) {
10049b5e25fSSatish Balay     ierr = PetscFree(*ia);CHKERRQ(ierr);
10149b5e25fSSatish Balay     ierr = PetscFree(*ja);CHKERRQ(ierr);
10249b5e25fSSatish Balay   } else if (oshift == 1) {
10349b5e25fSSatish Balay     int nz = a->i[n];
10449b5e25fSSatish Balay     for (i=0; i<nz; i++) a->j[i]--;
10549b5e25fSSatish Balay     for (i=0; i<n+1; i++) a->i[i]--;
10649b5e25fSSatish Balay   }
10749b5e25fSSatish Balay   PetscFunctionReturn(0);
10849b5e25fSSatish Balay }
10949b5e25fSSatish Balay 
11049b5e25fSSatish Balay #undef __FUNC__
11149b5e25fSSatish Balay #define __FUNC__ "MatGetBlockSize_SeqSBAIJ"
11249b5e25fSSatish Balay int MatGetBlockSize_SeqSBAIJ(Mat mat,int *bs)
11349b5e25fSSatish Balay {
11449b5e25fSSatish Balay   Mat_SeqSBAIJ *sbaij = (Mat_SeqSBAIJ*)mat->data;
11549b5e25fSSatish Balay 
11649b5e25fSSatish Balay   PetscFunctionBegin;
11749b5e25fSSatish Balay   *bs = sbaij->bs;
11849b5e25fSSatish Balay   PetscFunctionReturn(0);
11949b5e25fSSatish Balay }
12049b5e25fSSatish Balay 
12149b5e25fSSatish Balay #undef __FUNC__
12249b5e25fSSatish Balay #define __FUNC__ "MatDestroy_SeqSBAIJ"
12349b5e25fSSatish Balay int MatDestroy_SeqSBAIJ(Mat A)
12449b5e25fSSatish Balay {
12549b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
12649b5e25fSSatish Balay   int         ierr;
12749b5e25fSSatish Balay 
12849b5e25fSSatish Balay   PetscFunctionBegin;
12949b5e25fSSatish Balay   if (A->mapping) {
13049b5e25fSSatish Balay     ierr = ISLocalToGlobalMappingDestroy(A->mapping);CHKERRQ(ierr);
13149b5e25fSSatish Balay   }
13249b5e25fSSatish Balay   if (A->bmapping) {
13349b5e25fSSatish Balay     ierr = ISLocalToGlobalMappingDestroy(A->bmapping);CHKERRQ(ierr);
13449b5e25fSSatish Balay   }
13549b5e25fSSatish Balay   if (A->rmap) {
13649b5e25fSSatish Balay     ierr = MapDestroy(A->rmap);CHKERRQ(ierr);
13749b5e25fSSatish Balay   }
13849b5e25fSSatish Balay   if (A->cmap) {
13949b5e25fSSatish Balay     ierr = MapDestroy(A->cmap);CHKERRQ(ierr);
14049b5e25fSSatish Balay   }
14149b5e25fSSatish Balay #if defined(PETSC_USE_LOG)
14249b5e25fSSatish Balay   PLogObjectState((PetscObject)A,"Rows=%d, s_NZ=%d",a->m,a->s_nz);
14349b5e25fSSatish Balay #endif
14449b5e25fSSatish Balay   ierr = PetscFree(a->a);CHKERRQ(ierr);
14549b5e25fSSatish Balay   if (!a->singlemalloc) {
14649b5e25fSSatish Balay     ierr = PetscFree(a->i);CHKERRQ(ierr);
14763c8bf9fSHong Zhang     ierr = PetscFree(a->j);CHKERRQ(ierr);
14849b5e25fSSatish Balay   }
14949b5e25fSSatish Balay   if (a->row) {
15049b5e25fSSatish Balay     ierr = ISDestroy(a->row);CHKERRQ(ierr);
15149b5e25fSSatish Balay   }
15249b5e25fSSatish Balay   if (a->col) {
15349b5e25fSSatish Balay     ierr = ISDestroy(a->col);CHKERRQ(ierr);
15449b5e25fSSatish Balay   }
15549b5e25fSSatish Balay   if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);}
15649b5e25fSSatish Balay   if (a->ilen) {ierr = PetscFree(a->ilen);CHKERRQ(ierr);}
15749b5e25fSSatish Balay   if (a->imax) {ierr = PetscFree(a->imax);CHKERRQ(ierr);}
15849b5e25fSSatish Balay   if (a->solve_work) {ierr = PetscFree(a->solve_work);CHKERRQ(ierr);}
15949b5e25fSSatish Balay   if (a->mult_work) {ierr = PetscFree(a->mult_work);CHKERRQ(ierr);}
16049b5e25fSSatish Balay   if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);}
16149b5e25fSSatish Balay   if (a->saved_values) {ierr = PetscFree(a->saved_values);CHKERRQ(ierr);}
16249b5e25fSSatish Balay   ierr = PetscFree(a);CHKERRQ(ierr);
16349b5e25fSSatish Balay   PLogObjectDestroy(A);
16449b5e25fSSatish Balay   PetscHeaderDestroy(A);
16549b5e25fSSatish Balay   PetscFunctionReturn(0);
16649b5e25fSSatish Balay }
16749b5e25fSSatish Balay 
16849b5e25fSSatish Balay #undef __FUNC__
16949b5e25fSSatish Balay #define __FUNC__ "MatSetOption_SeqSBAIJ"
17049b5e25fSSatish Balay int MatSetOption_SeqSBAIJ(Mat A,MatOption op)
17149b5e25fSSatish Balay {
17249b5e25fSSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
17349b5e25fSSatish Balay 
17449b5e25fSSatish Balay   PetscFunctionBegin;
17549b5e25fSSatish Balay   if      (op == MAT_ROW_ORIENTED)                 a->roworiented    = PETSC_TRUE;
17649b5e25fSSatish Balay   else if (op == MAT_COLUMN_ORIENTED)              a->roworiented    = PETSC_FALSE;
17749b5e25fSSatish Balay   else if (op == MAT_COLUMNS_SORTED)               a->sorted         = PETSC_TRUE;
17849b5e25fSSatish Balay   else if (op == MAT_COLUMNS_UNSORTED)             a->sorted         = PETSC_FALSE;
17949b5e25fSSatish Balay   else if (op == MAT_KEEP_ZEROED_ROWS)             a->keepzeroedrows = PETSC_TRUE;
18049b5e25fSSatish Balay   else if (op == MAT_NO_NEW_NONZERO_LOCATIONS)     a->nonew          = 1;
18149b5e25fSSatish Balay   else if (op == MAT_NEW_NONZERO_LOCATION_ERR)     a->nonew          = -1;
18249b5e25fSSatish Balay   else if (op == MAT_NEW_NONZERO_ALLOCATION_ERR)   a->nonew          = -2;
18349b5e25fSSatish Balay   else if (op == MAT_YES_NEW_NONZERO_LOCATIONS)    a->nonew          = 0;
18449b5e25fSSatish Balay   else if (op == MAT_ROWS_SORTED ||
18549b5e25fSSatish Balay            op == MAT_ROWS_UNSORTED ||
18649b5e25fSSatish Balay            op == MAT_SYMMETRIC ||
18749b5e25fSSatish Balay            op == MAT_STRUCTURALLY_SYMMETRIC ||
18849b5e25fSSatish Balay            op == MAT_YES_NEW_DIAGONALS ||
18949b5e25fSSatish Balay            op == MAT_IGNORE_OFF_PROC_ENTRIES ||
19049b5e25fSSatish Balay            op == MAT_USE_HASH_TABLE) {
19149b5e25fSSatish Balay     PLogInfo(A,"MatSetOption_SeqBAIJ:Option ignored\n");
19249b5e25fSSatish Balay   } else if (op == MAT_NO_NEW_DIAGONALS) {
19349b5e25fSSatish Balay     SETERRQ(PETSC_ERR_SUP,0,"MAT_NO_NEW_DIAGONALS");
19449b5e25fSSatish Balay   } else {
19549b5e25fSSatish Balay     SETERRQ(PETSC_ERR_SUP,0,"unknown option");
19649b5e25fSSatish Balay   }
19749b5e25fSSatish Balay   PetscFunctionReturn(0);
19849b5e25fSSatish Balay }
19949b5e25fSSatish Balay 
20049b5e25fSSatish Balay 
20149b5e25fSSatish Balay #undef __FUNC__
20249b5e25fSSatish Balay #define __FUNC__ "MatGetSize_SeqSBAIJ"
20349b5e25fSSatish Balay int MatGetSize_SeqSBAIJ(Mat A,int *m,int *n)
20449b5e25fSSatish Balay {
20549b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
20649b5e25fSSatish Balay 
20749b5e25fSSatish Balay   PetscFunctionBegin;
20849b5e25fSSatish Balay   if (m) *m = a->m;
20949b5e25fSSatish Balay   if (n) *n = a->m;
21049b5e25fSSatish Balay   PetscFunctionReturn(0);
21149b5e25fSSatish Balay }
21249b5e25fSSatish Balay 
21349b5e25fSSatish Balay #undef __FUNC__
21449b5e25fSSatish Balay #define __FUNC__ "MatGetOwnershipRange_SeqSBAIJ"
21549b5e25fSSatish Balay int MatGetOwnershipRange_SeqSBAIJ(Mat A,int *m,int *n)
21649b5e25fSSatish Balay {
21749b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
21849b5e25fSSatish Balay 
21949b5e25fSSatish Balay   PetscFunctionBegin;
22049b5e25fSSatish Balay   *m = 0; *n = a->m;
22149b5e25fSSatish Balay   PetscFunctionReturn(0);
22249b5e25fSSatish Balay }
22349b5e25fSSatish Balay 
22449b5e25fSSatish Balay #undef __FUNC__
22549b5e25fSSatish Balay #define __FUNC__ "MatGetRow_SeqSBAIJ"
22649b5e25fSSatish Balay int MatGetRow_SeqSBAIJ(Mat A,int row,int *ncols,int **cols,Scalar **v)
22749b5e25fSSatish Balay {
22849b5e25fSSatish Balay   Mat_SeqSBAIJ  *a = (Mat_SeqSBAIJ*)A->data;
22949b5e25fSSatish Balay   int          itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*cols_i,bs2;
23049b5e25fSSatish Balay   MatScalar    *aa,*aa_i;
23149b5e25fSSatish Balay   Scalar       *v_i;
23249b5e25fSSatish Balay 
23349b5e25fSSatish Balay   PetscFunctionBegin;
23449b5e25fSSatish Balay   bs  = a->bs;
23549b5e25fSSatish Balay   ai  = a->i;
23649b5e25fSSatish Balay   aj  = a->j;
23749b5e25fSSatish Balay   aa  = a->a;
23849b5e25fSSatish Balay   bs2 = a->bs2;
23949b5e25fSSatish Balay 
24049b5e25fSSatish Balay   if (row < 0 || row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row out of range");
24149b5e25fSSatish Balay 
24249b5e25fSSatish Balay   bn  = row/bs;   /* Block number */
24349b5e25fSSatish Balay   bp  = row % bs; /* Block position */
24449b5e25fSSatish Balay   M   = ai[bn+1] - ai[bn];
24549b5e25fSSatish Balay   *ncols = bs*M;
24649b5e25fSSatish Balay 
24749b5e25fSSatish Balay   if (v) {
24849b5e25fSSatish Balay     *v = 0;
24949b5e25fSSatish Balay     if (*ncols) {
25049b5e25fSSatish Balay       *v = (Scalar*)PetscMalloc((*ncols+row)*sizeof(Scalar));CHKPTRQ(*v);
25149b5e25fSSatish Balay       for (i=0; i<M; i++) { /* for each block in the block row */
25249b5e25fSSatish Balay         v_i  = *v + i*bs;
25349b5e25fSSatish Balay         aa_i = aa + bs2*(ai[bn] + i);
25449b5e25fSSatish Balay         for (j=bp,k=0; j<bs2; j+=bs,k++) {v_i[k] = aa_i[j];}
25549b5e25fSSatish Balay       }
25649b5e25fSSatish Balay     }
25749b5e25fSSatish Balay   }
25849b5e25fSSatish Balay 
25949b5e25fSSatish Balay   if (cols) {
26049b5e25fSSatish Balay     *cols = 0;
26149b5e25fSSatish Balay     if (*ncols) {
26249b5e25fSSatish Balay       *cols = (int*)PetscMalloc((*ncols+row)*sizeof(int));CHKPTRQ(*cols);
26349b5e25fSSatish Balay       for (i=0; i<M; i++) { /* for each block in the block row */
26449b5e25fSSatish Balay         cols_i = *cols + i*bs;
26549b5e25fSSatish Balay         itmp  = bs*aj[ai[bn] + i];
26649b5e25fSSatish Balay         for (j=0; j<bs; j++) {cols_i[j] = itmp++;}
26749b5e25fSSatish Balay       }
26849b5e25fSSatish Balay     }
26949b5e25fSSatish Balay   }
27049b5e25fSSatish Balay 
27149b5e25fSSatish Balay   /*search column A(0:row-1,row) (=A(row,0:row-1)). Could be expensive! */
2725ddb2528SHong Zhang   /* this segment is currently removed, so only entries in the upper triangle are obtained */
2735ddb2528SHong Zhang #ifdef column_search
27449b5e25fSSatish Balay   v_i    = *v    + M*bs;
27549b5e25fSSatish Balay   cols_i = *cols + M*bs;
27649b5e25fSSatish Balay   for (i=0; i<bn; i++){ /* for each block row */
27749b5e25fSSatish Balay     M = ai[i+1] - ai[i];
27849b5e25fSSatish Balay     for (j=0; j<M; j++){
27949b5e25fSSatish Balay       itmp = aj[ai[i] + j];    /* block column value */
28049b5e25fSSatish Balay       if (itmp == bn){
28149b5e25fSSatish Balay         aa_i   = aa    + bs2*(ai[i] + j) + bs*bp;
28249b5e25fSSatish Balay         for (k=0; k<bs; k++) {
28349b5e25fSSatish Balay           *cols_i++ = i*bs+k;
28449b5e25fSSatish Balay           *v_i++    = aa_i[k];
28549b5e25fSSatish Balay         }
28649b5e25fSSatish Balay         *ncols += bs;
28749b5e25fSSatish Balay         break;
28849b5e25fSSatish Balay       }
28949b5e25fSSatish Balay     }
29049b5e25fSSatish Balay   }
2915ddb2528SHong Zhang #endif
29249b5e25fSSatish Balay 
29349b5e25fSSatish Balay   PetscFunctionReturn(0);
29449b5e25fSSatish Balay }
29549b5e25fSSatish Balay 
29649b5e25fSSatish Balay #undef __FUNC__
29749b5e25fSSatish Balay #define __FUNC__ "MatRestoreRow_SeqSBAIJ"
29849b5e25fSSatish Balay int MatRestoreRow_SeqSBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
29949b5e25fSSatish Balay {
30049b5e25fSSatish Balay   int ierr;
30149b5e25fSSatish Balay 
30249b5e25fSSatish Balay   PetscFunctionBegin;
30349b5e25fSSatish Balay   if (idx) {if (*idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}}
30449b5e25fSSatish Balay   if (v)   {if (*v)   {ierr = PetscFree(*v);CHKERRQ(ierr);}}
30549b5e25fSSatish Balay   PetscFunctionReturn(0);
30649b5e25fSSatish Balay }
30749b5e25fSSatish Balay 
30849b5e25fSSatish Balay #undef __FUNC__
30949b5e25fSSatish Balay #define __FUNC__ "MatTranspose_SeqSBAIJ"
31049b5e25fSSatish Balay int MatTranspose_SeqSBAIJ(Mat A,Mat *B)
31149b5e25fSSatish Balay {
31249b5e25fSSatish Balay   Mat_SeqBAIJ *a=(Mat_SeqBAIJ *)A->data;
31349b5e25fSSatish Balay   Mat         C;
31449b5e25fSSatish Balay   int         i,j,k,ierr,*aj=a->j,*ai=a->i,bs=a->bs,mbs=a->mbs,nbs=a->nbs,len,*col;
31549b5e25fSSatish Balay   int         *rows,*cols,bs2=a->bs2,refct;
31649b5e25fSSatish Balay   MatScalar   *array = a->a;
31749b5e25fSSatish Balay 
31849b5e25fSSatish Balay   PetscFunctionBegin;
31949b5e25fSSatish Balay   if (!B && mbs!=nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Square matrix only for in-place");
32049b5e25fSSatish Balay   col  = (int*)PetscMalloc((1+nbs)*sizeof(int));CHKPTRQ(col);
32149b5e25fSSatish Balay   ierr = PetscMemzero(col,(1+nbs)*sizeof(int));CHKERRQ(ierr);
32249b5e25fSSatish Balay 
32349b5e25fSSatish Balay   for (i=0; i<ai[mbs]; i++) col[aj[i]] += 1;
32449b5e25fSSatish Balay   ierr = MatCreateSeqBAIJ(A->comm,bs,a->n,a->m,PETSC_NULL,col,&C);CHKERRQ(ierr);
32549b5e25fSSatish Balay   ierr = PetscFree(col);CHKERRQ(ierr);
32649b5e25fSSatish Balay   rows = (int*)PetscMalloc(2*bs*sizeof(int));CHKPTRQ(rows);
32749b5e25fSSatish Balay   cols = rows + bs;
32849b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
32949b5e25fSSatish Balay     cols[0] = i*bs;
33049b5e25fSSatish Balay     for (k=1; k<bs; k++) cols[k] = cols[k-1] + 1;
33149b5e25fSSatish Balay     len = ai[i+1] - ai[i];
33249b5e25fSSatish Balay     for (j=0; j<len; j++) {
33349b5e25fSSatish Balay       rows[0] = (*aj++)*bs;
33449b5e25fSSatish Balay       for (k=1; k<bs; k++) rows[k] = rows[k-1] + 1;
33549b5e25fSSatish Balay       ierr = MatSetValues(C,bs,rows,bs,cols,array,INSERT_VALUES);CHKERRQ(ierr);
33649b5e25fSSatish Balay       array += bs2;
33749b5e25fSSatish Balay     }
33849b5e25fSSatish Balay   }
33949b5e25fSSatish Balay   ierr = PetscFree(rows);CHKERRQ(ierr);
34049b5e25fSSatish Balay 
34149b5e25fSSatish Balay   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
34249b5e25fSSatish Balay   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
34349b5e25fSSatish Balay 
34449b5e25fSSatish Balay   if (B) {
34549b5e25fSSatish Balay     *B = C;
34649b5e25fSSatish Balay   } else {
34749b5e25fSSatish Balay     PetscOps *Abops;
34849b5e25fSSatish Balay     MatOps   Aops;
34949b5e25fSSatish Balay 
35049b5e25fSSatish Balay     /* This isn't really an in-place transpose */
35149b5e25fSSatish Balay     ierr = PetscFree(a->a);CHKERRQ(ierr);
35249b5e25fSSatish Balay     if (!a->singlemalloc) {
35349b5e25fSSatish Balay       ierr = PetscFree(a->i);CHKERRQ(ierr);
35449b5e25fSSatish Balay       ierr = PetscFree(a->j);CHKERRQ(ierr);
35549b5e25fSSatish Balay     }
35649b5e25fSSatish Balay     if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);}
35749b5e25fSSatish Balay     if (a->ilen) {ierr = PetscFree(a->ilen);CHKERRQ(ierr);}
35849b5e25fSSatish Balay     if (a->imax) {ierr = PetscFree(a->imax);CHKERRQ(ierr);}
35949b5e25fSSatish Balay     if (a->solve_work) {ierr = PetscFree(a->solve_work);CHKERRQ(ierr);}
36049b5e25fSSatish Balay     ierr = PetscFree(a);CHKERRQ(ierr);
36149b5e25fSSatish Balay 
36249b5e25fSSatish Balay 
36349b5e25fSSatish Balay     ierr = MapDestroy(A->rmap);CHKERRQ(ierr);
36449b5e25fSSatish Balay     ierr = MapDestroy(A->cmap);CHKERRQ(ierr);
36549b5e25fSSatish Balay 
36649b5e25fSSatish Balay     /*
36749b5e25fSSatish Balay        This is horrible, horrible code. We need to keep the
36849b5e25fSSatish Balay       A pointers for the bops and ops but copy everything
36949b5e25fSSatish Balay       else from C.
37049b5e25fSSatish Balay     */
37149b5e25fSSatish Balay     Abops    = A->bops;
37249b5e25fSSatish Balay     Aops     = A->ops;
37349b5e25fSSatish Balay     refct    = A->refct;
37449b5e25fSSatish Balay     ierr     = PetscMemcpy(A,C,sizeof(struct _p_Mat));CHKERRQ(ierr);
37549b5e25fSSatish Balay     A->bops  = Abops;
37649b5e25fSSatish Balay     A->ops   = Aops;
37749b5e25fSSatish Balay     A->qlist = 0;
37849b5e25fSSatish Balay     A->refct = refct;
37949b5e25fSSatish Balay 
38049b5e25fSSatish Balay     PetscHeaderDestroy(C);
38149b5e25fSSatish Balay   }
38249b5e25fSSatish Balay   PetscFunctionReturn(0);
38349b5e25fSSatish Balay }
38449b5e25fSSatish Balay 
38549b5e25fSSatish Balay #undef __FUNC__
38649b5e25fSSatish Balay #define __FUNC__ "MatView_SeqSBAIJ_Binary"
38749b5e25fSSatish Balay static int MatView_SeqSBAIJ_Binary(Mat A,Viewer viewer)
38849b5e25fSSatish Balay {
38949b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
39049b5e25fSSatish Balay   int         i,fd,*col_lens,ierr,bs = a->bs,count,*jj,j,k,l,bs2=a->bs2;
39149b5e25fSSatish Balay   Scalar      *aa;
39249b5e25fSSatish Balay   FILE        *file;
39349b5e25fSSatish Balay 
39449b5e25fSSatish Balay   PetscFunctionBegin;
39549b5e25fSSatish Balay   ierr = ViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
39649b5e25fSSatish Balay   col_lens = (int*)PetscMalloc((4+a->m)*sizeof(int));CHKPTRQ(col_lens);
39749b5e25fSSatish Balay   col_lens[0] = MAT_COOKIE;
39849b5e25fSSatish Balay 
39949b5e25fSSatish Balay   col_lens[1] = a->m;
40049b5e25fSSatish Balay   col_lens[2] = a->m;
40149b5e25fSSatish Balay   col_lens[3] = a->s_nz*bs2;
40249b5e25fSSatish Balay 
40349b5e25fSSatish Balay   /* store lengths of each row and write (including header) to file */
40449b5e25fSSatish Balay   count = 0;
40549b5e25fSSatish Balay   for (i=0; i<a->mbs; i++) {
40649b5e25fSSatish Balay     for (j=0; j<bs; j++) {
40749b5e25fSSatish Balay       col_lens[4+count++] = bs*(a->i[i+1] - a->i[i]);
40849b5e25fSSatish Balay     }
40949b5e25fSSatish Balay   }
41049b5e25fSSatish Balay   ierr = PetscBinaryWrite(fd,col_lens,4+a->m,PETSC_INT,1);CHKERRQ(ierr);
41149b5e25fSSatish Balay   ierr = PetscFree(col_lens);CHKERRQ(ierr);
41249b5e25fSSatish Balay 
41349b5e25fSSatish Balay   /* store column indices (zero start index) */
41449b5e25fSSatish Balay   jj = (int*)PetscMalloc((a->s_nz+1)*bs2*sizeof(int));CHKPTRQ(jj);
41549b5e25fSSatish Balay   count = 0;
41649b5e25fSSatish Balay   for (i=0; i<a->mbs; i++) {
41749b5e25fSSatish Balay     for (j=0; j<bs; j++) {
41849b5e25fSSatish Balay       for (k=a->i[i]; k<a->i[i+1]; k++) {
41949b5e25fSSatish Balay         for (l=0; l<bs; l++) {
42049b5e25fSSatish Balay           jj[count++] = bs*a->j[k] + l;
42149b5e25fSSatish Balay         }
42249b5e25fSSatish Balay       }
42349b5e25fSSatish Balay     }
42449b5e25fSSatish Balay   }
42549b5e25fSSatish Balay   ierr = PetscBinaryWrite(fd,jj,bs2*a->s_nz,PETSC_INT,0);CHKERRQ(ierr);
42649b5e25fSSatish Balay   ierr = PetscFree(jj);CHKERRQ(ierr);
42749b5e25fSSatish Balay 
42849b5e25fSSatish Balay   /* store nonzero values */
42949b5e25fSSatish Balay   aa = (Scalar*)PetscMalloc((a->s_nz+1)*bs2*sizeof(Scalar));CHKPTRQ(aa);
43049b5e25fSSatish Balay   count = 0;
43149b5e25fSSatish Balay   for (i=0; i<a->mbs; i++) {
43249b5e25fSSatish Balay     for (j=0; j<bs; j++) {
43349b5e25fSSatish Balay       for (k=a->i[i]; k<a->i[i+1]; k++) {
43449b5e25fSSatish Balay         for (l=0; l<bs; l++) {
43549b5e25fSSatish Balay           aa[count++] = a->a[bs2*k + l*bs + j];
43649b5e25fSSatish Balay         }
43749b5e25fSSatish Balay       }
43849b5e25fSSatish Balay     }
43949b5e25fSSatish Balay   }
44049b5e25fSSatish Balay   ierr = PetscBinaryWrite(fd,aa,bs2*a->s_nz,PETSC_SCALAR,0);CHKERRQ(ierr);
44149b5e25fSSatish Balay   ierr = PetscFree(aa);CHKERRQ(ierr);
44249b5e25fSSatish Balay 
44349b5e25fSSatish Balay   ierr = ViewerBinaryGetInfoPointer(viewer,&file);CHKERRQ(ierr);
44449b5e25fSSatish Balay   if (file) {
44549b5e25fSSatish Balay     fprintf(file,"-matload_block_size %d\n",a->bs);
44649b5e25fSSatish Balay   }
44749b5e25fSSatish Balay   PetscFunctionReturn(0);
44849b5e25fSSatish Balay }
44949b5e25fSSatish Balay 
45049b5e25fSSatish Balay #undef __FUNC__
45149b5e25fSSatish Balay #define __FUNC__ "MatView_SeqSBAIJ_ASCII"
45249b5e25fSSatish Balay static int MatView_SeqSBAIJ_ASCII(Mat A,Viewer viewer)
45349b5e25fSSatish Balay {
45449b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
45549b5e25fSSatish Balay   int         ierr,i,j,format,bs = a->bs,k,l,bs2=a->bs2;
45649b5e25fSSatish Balay   char        *outputname;
45749b5e25fSSatish Balay 
45849b5e25fSSatish Balay   PetscFunctionBegin;
45949b5e25fSSatish Balay   ierr = ViewerGetOutputname(viewer,&outputname);CHKERRQ(ierr);
46049b5e25fSSatish Balay   ierr = ViewerGetFormat(viewer,&format);CHKERRQ(ierr);
46149b5e25fSSatish Balay   if (format == VIEWER_FORMAT_ASCII_INFO || format == VIEWER_FORMAT_ASCII_INFO_LONG) {
46249b5e25fSSatish Balay     ierr = ViewerASCIIPrintf(viewer,"  block size is %d\n",bs);CHKERRQ(ierr);
46349b5e25fSSatish Balay   } else if (format == VIEWER_FORMAT_ASCII_MATLAB) {
46449b5e25fSSatish Balay     SETERRQ(PETSC_ERR_SUP,0,"Matlab format not supported");
46549b5e25fSSatish Balay   } else if (format == VIEWER_FORMAT_ASCII_COMMON) {
46649b5e25fSSatish Balay     ierr = ViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
46749b5e25fSSatish Balay     for (i=0; i<a->mbs; i++) {
46849b5e25fSSatish Balay       for (j=0; j<bs; j++) {
46949b5e25fSSatish Balay         ierr = ViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr);
47049b5e25fSSatish Balay         for (k=a->i[i]; k<a->i[i+1]; k++) {
47149b5e25fSSatish Balay           for (l=0; l<bs; l++) {
47249b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX)
47349b5e25fSSatish Balay             if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
47449b5e25fSSatish Balay               ierr = ViewerASCIIPrintf(viewer," %d %g + %g i",bs*a->j[k]+l,
47549b5e25fSSatish Balay                       PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
47649b5e25fSSatish Balay             } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
47749b5e25fSSatish Balay               ierr = ViewerASCIIPrintf(viewer," %d %g - %g i",bs*a->j[k]+l,
47849b5e25fSSatish Balay                       PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
47949b5e25fSSatish Balay             } else if (PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
48049b5e25fSSatish Balay               ierr = ViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
48149b5e25fSSatish Balay             }
48249b5e25fSSatish Balay #else
48349b5e25fSSatish Balay             if (a->a[bs2*k + l*bs + j] != 0.0) {
48449b5e25fSSatish Balay               ierr = ViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr);
48549b5e25fSSatish Balay             }
48649b5e25fSSatish Balay #endif
48749b5e25fSSatish Balay           }
48849b5e25fSSatish Balay         }
48949b5e25fSSatish Balay         ierr = ViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
49049b5e25fSSatish Balay       }
49149b5e25fSSatish Balay     }
49249b5e25fSSatish Balay     ierr = ViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
49349b5e25fSSatish Balay   } else {
49449b5e25fSSatish Balay     ierr = ViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
49549b5e25fSSatish Balay     for (i=0; i<a->mbs; i++) {
49649b5e25fSSatish Balay       for (j=0; j<bs; j++) {
49749b5e25fSSatish Balay         ierr = ViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr);
49849b5e25fSSatish Balay         for (k=a->i[i]; k<a->i[i+1]; k++) {
49949b5e25fSSatish Balay           for (l=0; l<bs; l++) {
50049b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX)
50149b5e25fSSatish Balay             if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0) {
50249b5e25fSSatish Balay               ierr = ViewerASCIIPrintf(viewer," %d %g + %g i",bs*a->j[k]+l,
50349b5e25fSSatish Balay                 PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
50449b5e25fSSatish Balay             } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0) {
50549b5e25fSSatish Balay               ierr = ViewerASCIIPrintf(viewer," %d %g - %g i",bs*a->j[k]+l,
50649b5e25fSSatish Balay                 PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
50749b5e25fSSatish Balay             } else {
50849b5e25fSSatish Balay               ierr = ViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
50949b5e25fSSatish Balay             }
51049b5e25fSSatish Balay #else
51149b5e25fSSatish Balay             ierr = ViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr);
51249b5e25fSSatish Balay #endif
51349b5e25fSSatish Balay           }
51449b5e25fSSatish Balay         }
51549b5e25fSSatish Balay         ierr = ViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
51649b5e25fSSatish Balay       }
51749b5e25fSSatish Balay     }
51849b5e25fSSatish Balay     ierr = ViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
51949b5e25fSSatish Balay   }
52049b5e25fSSatish Balay   ierr = ViewerFlush(viewer);CHKERRQ(ierr);
52149b5e25fSSatish Balay   PetscFunctionReturn(0);
52249b5e25fSSatish Balay }
52349b5e25fSSatish Balay 
52449b5e25fSSatish Balay #undef __FUNC__
52549b5e25fSSatish Balay #define __FUNC__ "MatView_SeqSBAIJ_Draw_Zoom"
52649b5e25fSSatish Balay static int MatView_SeqSBAIJ_Draw_Zoom(Draw draw,void *Aa)
52749b5e25fSSatish Balay {
52849b5e25fSSatish Balay   Mat          A = (Mat) Aa;
52949b5e25fSSatish Balay   Mat_SeqSBAIJ  *a=(Mat_SeqSBAIJ*)A->data;
53049b5e25fSSatish Balay   int          row,ierr,i,j,k,l,mbs=a->mbs,color,bs=a->bs,bs2=a->bs2,rank;
53149b5e25fSSatish Balay   PetscReal    xl,yl,xr,yr,x_l,x_r,y_l,y_r;
53249b5e25fSSatish Balay   MatScalar    *aa;
53349b5e25fSSatish Balay   MPI_Comm     comm;
53449b5e25fSSatish Balay   Viewer       viewer;
53549b5e25fSSatish Balay 
53649b5e25fSSatish Balay   PetscFunctionBegin;
53749b5e25fSSatish Balay   /*
53849b5e25fSSatish Balay       This is nasty. If this is called from an originally parallel matrix
53949b5e25fSSatish Balay    then all processes call this,but only the first has the matrix so the
54049b5e25fSSatish Balay    rest should return immediately.
54149b5e25fSSatish Balay   */
54249b5e25fSSatish Balay   ierr = PetscObjectGetComm((PetscObject)draw,&comm);CHKERRQ(ierr);
54349b5e25fSSatish Balay   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
54449b5e25fSSatish Balay   if (rank) PetscFunctionReturn(0);
54549b5e25fSSatish Balay 
54649b5e25fSSatish Balay   ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr);
54749b5e25fSSatish Balay 
54849b5e25fSSatish Balay   ierr = DrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr);
54949b5e25fSSatish Balay   DrawString(draw, .3*(xl+xr), .3*(yl+yr), DRAW_BLACK, "symmetric");
55049b5e25fSSatish Balay 
55149b5e25fSSatish Balay   /* loop over matrix elements drawing boxes */
55249b5e25fSSatish Balay   color = DRAW_BLUE;
55349b5e25fSSatish Balay   for (i=0,row=0; i<mbs; i++,row+=bs) {
55449b5e25fSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
55549b5e25fSSatish Balay       y_l = a->m - row - 1.0; y_r = y_l + 1.0;
55649b5e25fSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
55749b5e25fSSatish Balay       aa = a->a + j*bs2;
55849b5e25fSSatish Balay       for (k=0; k<bs; k++) {
55949b5e25fSSatish Balay         for (l=0; l<bs; l++) {
56049b5e25fSSatish Balay           if (PetscRealPart(*aa++) >=  0.) continue;
56149b5e25fSSatish Balay           ierr = DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
56249b5e25fSSatish Balay         }
56349b5e25fSSatish Balay       }
56449b5e25fSSatish Balay     }
56549b5e25fSSatish Balay   }
56649b5e25fSSatish Balay   color = DRAW_CYAN;
56749b5e25fSSatish Balay   for (i=0,row=0; i<mbs; i++,row+=bs) {
56849b5e25fSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
56949b5e25fSSatish Balay       y_l = a->m - row - 1.0; y_r = y_l + 1.0;
57049b5e25fSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
57149b5e25fSSatish Balay       aa = a->a + j*bs2;
57249b5e25fSSatish Balay       for (k=0; k<bs; k++) {
57349b5e25fSSatish Balay         for (l=0; l<bs; l++) {
57449b5e25fSSatish Balay           if (PetscRealPart(*aa++) != 0.) continue;
57549b5e25fSSatish Balay           ierr = DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
57649b5e25fSSatish Balay         }
57749b5e25fSSatish Balay       }
57849b5e25fSSatish Balay     }
57949b5e25fSSatish Balay   }
58049b5e25fSSatish Balay 
58149b5e25fSSatish Balay   color = DRAW_RED;
58249b5e25fSSatish Balay   for (i=0,row=0; i<mbs; i++,row+=bs) {
58349b5e25fSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
58449b5e25fSSatish Balay       y_l = a->m - row - 1.0; y_r = y_l + 1.0;
58549b5e25fSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
58649b5e25fSSatish Balay       aa = a->a + j*bs2;
58749b5e25fSSatish Balay       for (k=0; k<bs; k++) {
58849b5e25fSSatish Balay         for (l=0; l<bs; l++) {
58949b5e25fSSatish Balay           if (PetscRealPart(*aa++) <= 0.) continue;
59049b5e25fSSatish Balay           ierr = DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
59149b5e25fSSatish Balay         }
59249b5e25fSSatish Balay       }
59349b5e25fSSatish Balay     }
59449b5e25fSSatish Balay   }
59549b5e25fSSatish Balay   PetscFunctionReturn(0);
59649b5e25fSSatish Balay }
59749b5e25fSSatish Balay 
59849b5e25fSSatish Balay #undef __FUNC__
59949b5e25fSSatish Balay #define __FUNC__ "MatView_SeqSBAIJ_Draw"
60049b5e25fSSatish Balay static int MatView_SeqSBAIJ_Draw(Mat A,Viewer viewer)
60149b5e25fSSatish Balay {
60249b5e25fSSatish Balay   Mat_SeqSBAIJ  *a=(Mat_SeqSBAIJ*)A->data;
60349b5e25fSSatish Balay   int          ierr;
60449b5e25fSSatish Balay   PetscReal    xl,yl,xr,yr,w,h;
60549b5e25fSSatish Balay   Draw         draw;
60649b5e25fSSatish Balay   PetscTruth   isnull;
60749b5e25fSSatish Balay 
60849b5e25fSSatish Balay   PetscFunctionBegin;
60949b5e25fSSatish Balay 
61049b5e25fSSatish Balay   ierr = ViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
61149b5e25fSSatish Balay   ierr = DrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
61249b5e25fSSatish Balay 
61349b5e25fSSatish Balay   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
61449b5e25fSSatish Balay   xr  = a->m; yr = a->m; h = yr/10.0; w = xr/10.0;
61549b5e25fSSatish Balay   xr += w;    yr += h;  xl = -w;     yl = -h;
61649b5e25fSSatish Balay   ierr = DrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
61749b5e25fSSatish Balay   ierr = DrawZoom(draw,MatView_SeqSBAIJ_Draw_Zoom,A);CHKERRQ(ierr);
61849b5e25fSSatish Balay   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr);
61949b5e25fSSatish Balay   PetscFunctionReturn(0);
62049b5e25fSSatish Balay }
62149b5e25fSSatish Balay 
62249b5e25fSSatish Balay #undef __FUNC__
62349b5e25fSSatish Balay #define __FUNC__ "MatView_SeqSBAIJ"
62449b5e25fSSatish Balay int MatView_SeqSBAIJ(Mat A,Viewer viewer)
62549b5e25fSSatish Balay {
62649b5e25fSSatish Balay   int        ierr;
62749b5e25fSSatish Balay   PetscTruth issocket,isascii,isbinary,isdraw;
62849b5e25fSSatish Balay 
62949b5e25fSSatish Balay   PetscFunctionBegin;
63049b5e25fSSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,SOCKET_VIEWER,&issocket);CHKERRQ(ierr);
63149b5e25fSSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,ASCII_VIEWER,&isascii);CHKERRQ(ierr);
63249b5e25fSSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,BINARY_VIEWER,&isbinary);CHKERRQ(ierr);
63349b5e25fSSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,DRAW_VIEWER,&isdraw);CHKERRQ(ierr);
63449b5e25fSSatish Balay   if (issocket) {
63549b5e25fSSatish Balay     SETERRQ(PETSC_ERR_SUP,0,"Socket viewer not supported");
63649b5e25fSSatish Balay   } else if (isascii){
63749b5e25fSSatish Balay     ierr = MatView_SeqSBAIJ_ASCII(A,viewer);CHKERRQ(ierr);
63849b5e25fSSatish Balay   } else if (isbinary) {
63949b5e25fSSatish Balay     ierr = MatView_SeqSBAIJ_Binary(A,viewer);CHKERRQ(ierr);
64049b5e25fSSatish Balay   } else if (isdraw) {
64149b5e25fSSatish Balay     ierr = MatView_SeqSBAIJ_Draw(A,viewer);CHKERRQ(ierr);
64249b5e25fSSatish Balay   } else {
64349b5e25fSSatish Balay     SETERRQ1(1,1,"Viewer type %s not supported by SeqBAIJ matrices",((PetscObject)viewer)->type_name);
64449b5e25fSSatish Balay   }
64549b5e25fSSatish Balay   PetscFunctionReturn(0);
64649b5e25fSSatish Balay }
64749b5e25fSSatish Balay 
64849b5e25fSSatish Balay 
64949b5e25fSSatish Balay #undef __FUNC__
65049b5e25fSSatish Balay #define __FUNC__ "MatGetValues_SeqSBAIJ"
65149b5e25fSSatish Balay int MatGetValues_SeqSBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v)
65249b5e25fSSatish Balay {
65349b5e25fSSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
65449b5e25fSSatish Balay   int        *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j;
65549b5e25fSSatish Balay   int        *ai = a->i,*ailen = a->ilen;
65649b5e25fSSatish Balay   int        brow,bcol,ridx,cidx,bs=a->bs,bs2=a->bs2;
65749b5e25fSSatish Balay   MatScalar  *ap,*aa = a->a,zero = 0.0;
65849b5e25fSSatish Balay 
65949b5e25fSSatish Balay   PetscFunctionBegin;
66049b5e25fSSatish Balay   for (k=0; k<m; k++) { /* loop over rows */
66149b5e25fSSatish Balay     row  = im[k]; brow = row/bs;
66249b5e25fSSatish Balay     if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row");
66349b5e25fSSatish Balay     if (row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
66449b5e25fSSatish Balay     rp   = aj + ai[brow] ; ap = aa + bs2*ai[brow] ;
66549b5e25fSSatish Balay     nrow = ailen[brow];
66649b5e25fSSatish Balay     for (l=0; l<n; l++) { /* loop over columns */
66749b5e25fSSatish Balay       if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column");
66849b5e25fSSatish Balay       if (in[l] >= a->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");
66949b5e25fSSatish Balay       col  = in[l] ;
67049b5e25fSSatish Balay       bcol = col/bs;
67149b5e25fSSatish Balay       cidx = col%bs;
67249b5e25fSSatish Balay       ridx = row%bs;
67349b5e25fSSatish Balay       high = nrow;
67449b5e25fSSatish Balay       low  = 0; /* assume unsorted */
67549b5e25fSSatish Balay       while (high-low > 5) {
67649b5e25fSSatish Balay         t = (low+high)/2;
67749b5e25fSSatish Balay         if (rp[t] > bcol) high = t;
67849b5e25fSSatish Balay         else             low  = t;
67949b5e25fSSatish Balay       }
68049b5e25fSSatish Balay       for (i=low; i<high; i++) {
68149b5e25fSSatish Balay         if (rp[i] > bcol) break;
68249b5e25fSSatish Balay         if (rp[i] == bcol) {
68349b5e25fSSatish Balay           *v++ = ap[bs2*i+bs*cidx+ridx];
68449b5e25fSSatish Balay           goto finished;
68549b5e25fSSatish Balay         }
68649b5e25fSSatish Balay       }
68749b5e25fSSatish Balay       *v++ = zero;
68849b5e25fSSatish Balay       finished:;
68949b5e25fSSatish Balay     }
69049b5e25fSSatish Balay   }
69149b5e25fSSatish Balay   PetscFunctionReturn(0);
69249b5e25fSSatish Balay }
69349b5e25fSSatish Balay 
69449b5e25fSSatish Balay 
69549b5e25fSSatish Balay #undef __FUNC__
69649b5e25fSSatish Balay #define __FUNC__ "MatSetValuesBlocked_SeqSBAIJ"
69749b5e25fSSatish Balay int MatSetValuesBlocked_SeqSBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is)
69849b5e25fSSatish Balay {
69949b5e25fSSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
70049b5e25fSSatish Balay   int         *rp,k,low,high,t,ii,jj,row,nrow,i,col,l,rmax,N,sorted=a->sorted;
70149b5e25fSSatish Balay   int         *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented;
70249b5e25fSSatish Balay   int         *aj=a->j,nonew=a->nonew,bs2=a->bs2,bs=a->bs,stepval,ierr;
70349b5e25fSSatish Balay   Scalar      *value = v;
70449b5e25fSSatish Balay   MatScalar   *ap,*aa=a->a,*bap;
70549b5e25fSSatish Balay 
70649b5e25fSSatish Balay   PetscFunctionBegin;
70749b5e25fSSatish Balay   if (roworiented) {
70849b5e25fSSatish Balay     stepval = (n-1)*bs;
70949b5e25fSSatish Balay   } else {
71049b5e25fSSatish Balay     stepval = (m-1)*bs;
71149b5e25fSSatish Balay   }
71249b5e25fSSatish Balay   for (k=0; k<m; k++) { /* loop over added rows */
71349b5e25fSSatish Balay     row  = im[k];
71449b5e25fSSatish Balay     if (row < 0) continue;
71549b5e25fSSatish Balay #if defined(PETSC_USE_BOPT_g)
71649b5e25fSSatish Balay     if (row >= a->mbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
71749b5e25fSSatish Balay #endif
71849b5e25fSSatish Balay     rp   = aj + ai[row];
71949b5e25fSSatish Balay     ap   = aa + bs2*ai[row];
72049b5e25fSSatish Balay     rmax = imax[row];
72149b5e25fSSatish Balay     nrow = ailen[row];
72249b5e25fSSatish Balay     low  = 0;
72349b5e25fSSatish Balay     for (l=0; l<n; l++) { /* loop over added columns */
72449b5e25fSSatish Balay       if (in[l] < 0) continue;
72549b5e25fSSatish Balay #if defined(PETSC_USE_BOPT_g)
72649b5e25fSSatish Balay       if (in[l] >= a->nbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");
72749b5e25fSSatish Balay #endif
72849b5e25fSSatish Balay       col = in[l];
72949b5e25fSSatish Balay       if (roworiented) {
73049b5e25fSSatish Balay         value = v + k*(stepval+bs)*bs + l*bs;
73149b5e25fSSatish Balay       } else {
73249b5e25fSSatish Balay         value = v + l*(stepval+bs)*bs + k*bs;
73349b5e25fSSatish Balay       }
73449b5e25fSSatish Balay       if (!sorted) low = 0; high = nrow;
73549b5e25fSSatish Balay       while (high-low > 7) {
73649b5e25fSSatish Balay         t = (low+high)/2;
73749b5e25fSSatish Balay         if (rp[t] > col) high = t;
73849b5e25fSSatish Balay         else             low  = t;
73949b5e25fSSatish Balay       }
74049b5e25fSSatish Balay       for (i=low; i<high; i++) {
74149b5e25fSSatish Balay         if (rp[i] > col) break;
74249b5e25fSSatish Balay         if (rp[i] == col) {
74349b5e25fSSatish Balay           bap  = ap +  bs2*i;
74449b5e25fSSatish Balay           if (roworiented) {
74549b5e25fSSatish Balay             if (is == ADD_VALUES) {
74649b5e25fSSatish Balay               for (ii=0; ii<bs; ii++,value+=stepval) {
74749b5e25fSSatish Balay                 for (jj=ii; jj<bs2; jj+=bs) {
74849b5e25fSSatish Balay                   bap[jj] += *value++;
74949b5e25fSSatish Balay                 }
75049b5e25fSSatish Balay               }
75149b5e25fSSatish Balay             } else {
75249b5e25fSSatish Balay               for (ii=0; ii<bs; ii++,value+=stepval) {
75349b5e25fSSatish Balay                 for (jj=ii; jj<bs2; jj+=bs) {
75449b5e25fSSatish Balay                   bap[jj] = *value++;
75549b5e25fSSatish Balay                 }
75649b5e25fSSatish Balay               }
75749b5e25fSSatish Balay             }
75849b5e25fSSatish Balay           } else {
75949b5e25fSSatish Balay             if (is == ADD_VALUES) {
76049b5e25fSSatish Balay               for (ii=0; ii<bs; ii++,value+=stepval) {
76149b5e25fSSatish Balay                 for (jj=0; jj<bs; jj++) {
76249b5e25fSSatish Balay                   *bap++ += *value++;
76349b5e25fSSatish Balay                 }
76449b5e25fSSatish Balay               }
76549b5e25fSSatish Balay             } else {
76649b5e25fSSatish Balay               for (ii=0; ii<bs; ii++,value+=stepval) {
76749b5e25fSSatish Balay                 for (jj=0; jj<bs; jj++) {
76849b5e25fSSatish Balay                   *bap++  = *value++;
76949b5e25fSSatish Balay                 }
77049b5e25fSSatish Balay               }
77149b5e25fSSatish Balay             }
77249b5e25fSSatish Balay           }
77349b5e25fSSatish Balay           goto noinsert2;
77449b5e25fSSatish Balay         }
77549b5e25fSSatish Balay       }
77649b5e25fSSatish Balay       if (nonew == 1) goto noinsert2;
77749b5e25fSSatish Balay       else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix");
77849b5e25fSSatish Balay       if (nrow >= rmax) {
77949b5e25fSSatish Balay         /* there is no extra room in row, therefore enlarge */
78049b5e25fSSatish Balay         int       new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j;
78149b5e25fSSatish Balay         MatScalar *new_a;
78249b5e25fSSatish Balay 
78349b5e25fSSatish Balay         if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix");
78449b5e25fSSatish Balay 
78549b5e25fSSatish Balay         /* malloc new storage space */
78649b5e25fSSatish Balay         len     = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int);
78749b5e25fSSatish Balay         new_a   = (MatScalar*)PetscMalloc(len);CHKPTRQ(new_a);
78849b5e25fSSatish Balay         new_j   = (int*)(new_a + bs2*new_nz);
78949b5e25fSSatish Balay         new_i   = new_j + new_nz;
79049b5e25fSSatish Balay 
79149b5e25fSSatish Balay         /* copy over old data into new slots */
79249b5e25fSSatish Balay         for (ii=0; ii<row+1; ii++) {new_i[ii] = ai[ii];}
79349b5e25fSSatish Balay         for (ii=row+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;}
79449b5e25fSSatish Balay         ierr = PetscMemcpy(new_j,aj,(ai[row]+nrow)*sizeof(int));CHKERRQ(ierr);
79549b5e25fSSatish Balay         len  = (new_nz - CHUNKSIZE - ai[row] - nrow);
79649b5e25fSSatish Balay         ierr = PetscMemcpy(new_j+ai[row]+nrow+CHUNKSIZE,aj+ai[row]+nrow,len*sizeof(int));CHKERRQ(ierr);
79749b5e25fSSatish Balay         ierr = PetscMemcpy(new_a,aa,(ai[row]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr);
79849b5e25fSSatish Balay         ierr = PetscMemzero(new_a+bs2*(ai[row]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr);
79949b5e25fSSatish Balay         ierr = PetscMemcpy(new_a+bs2*(ai[row]+nrow+CHUNKSIZE),aa+bs2*(ai[row]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr);
80049b5e25fSSatish Balay         /* free up old matrix storage */
80149b5e25fSSatish Balay         ierr = PetscFree(a->a);CHKERRQ(ierr);
80249b5e25fSSatish Balay         if (!a->singlemalloc) {
80349b5e25fSSatish Balay           ierr = PetscFree(a->i);CHKERRQ(ierr);
80449b5e25fSSatish Balay           ierr = PetscFree(a->j);CHKERRQ(ierr);
80549b5e25fSSatish Balay         }
80649b5e25fSSatish Balay         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
80749b5e25fSSatish Balay         a->singlemalloc = PETSC_TRUE;
80849b5e25fSSatish Balay 
80949b5e25fSSatish Balay         rp   = aj + ai[row]; ap = aa + bs2*ai[row];
81049b5e25fSSatish Balay         rmax = imax[row] = imax[row] + CHUNKSIZE;
81149b5e25fSSatish Balay         PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar)));
81249b5e25fSSatish Balay         a->maxnz += bs2*CHUNKSIZE;
81349b5e25fSSatish Balay         a->reallocs++;
81449b5e25fSSatish Balay         a->nz++;
81549b5e25fSSatish Balay       }
81649b5e25fSSatish Balay       N = nrow++ - 1;
81749b5e25fSSatish Balay       /* shift up all the later entries in this row */
81849b5e25fSSatish Balay       for (ii=N; ii>=i; ii--) {
81949b5e25fSSatish Balay         rp[ii+1] = rp[ii];
82049b5e25fSSatish Balay         ierr = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr);
82149b5e25fSSatish Balay       }
82249b5e25fSSatish Balay       if (N >= i) {
82349b5e25fSSatish Balay         ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr);
82449b5e25fSSatish Balay       }
82549b5e25fSSatish Balay       rp[i] = col;
82649b5e25fSSatish Balay       bap   = ap +  bs2*i;
82749b5e25fSSatish Balay       if (roworiented) {
82849b5e25fSSatish Balay         for (ii=0; ii<bs; ii++,value+=stepval) {
82949b5e25fSSatish Balay           for (jj=ii; jj<bs2; jj+=bs) {
83049b5e25fSSatish Balay             bap[jj] = *value++;
83149b5e25fSSatish Balay           }
83249b5e25fSSatish Balay         }
83349b5e25fSSatish Balay       } else {
83449b5e25fSSatish Balay         for (ii=0; ii<bs; ii++,value+=stepval) {
83549b5e25fSSatish Balay           for (jj=0; jj<bs; jj++) {
83649b5e25fSSatish Balay             *bap++  = *value++;
83749b5e25fSSatish Balay           }
83849b5e25fSSatish Balay         }
83949b5e25fSSatish Balay       }
84049b5e25fSSatish Balay       noinsert2:;
84149b5e25fSSatish Balay       low = i;
84249b5e25fSSatish Balay     }
84349b5e25fSSatish Balay     ailen[row] = nrow;
84449b5e25fSSatish Balay   }
84549b5e25fSSatish Balay   PetscFunctionReturn(0);
84649b5e25fSSatish Balay }
84749b5e25fSSatish Balay 
84849b5e25fSSatish Balay #undef __FUNC__
84949b5e25fSSatish Balay #define __FUNC__ "MatAssemblyEnd_SeqSBAIJ"
85049b5e25fSSatish Balay int MatAssemblyEnd_SeqSBAIJ(Mat A,MatAssemblyType mode)
85149b5e25fSSatish Balay {
85249b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
85349b5e25fSSatish Balay   int        fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax;
85449b5e25fSSatish Balay   int        m = a->m,*ip,N,*ailen = a->ilen;
85549b5e25fSSatish Balay   int        mbs = a->mbs,bs2 = a->bs2,rmax = 0,ierr;
85649b5e25fSSatish Balay   MatScalar  *aa = a->a,*ap;
85749b5e25fSSatish Balay 
85849b5e25fSSatish Balay   PetscFunctionBegin;
85949b5e25fSSatish Balay   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0);
86049b5e25fSSatish Balay 
86149b5e25fSSatish Balay   if (m) rmax = ailen[0];
86249b5e25fSSatish Balay   for (i=1; i<mbs; i++) {
86349b5e25fSSatish Balay     /* move each row back by the amount of empty slots (fshift) before it*/
86449b5e25fSSatish Balay     fshift += imax[i-1] - ailen[i-1];
86549b5e25fSSatish Balay     rmax   = PetscMax(rmax,ailen[i]);
86649b5e25fSSatish Balay     if (fshift) {
86749b5e25fSSatish Balay       ip = aj + ai[i]; ap = aa + bs2*ai[i];
86849b5e25fSSatish Balay       N = ailen[i];
86949b5e25fSSatish Balay       for (j=0; j<N; j++) {
87049b5e25fSSatish Balay         ip[j-fshift] = ip[j];
87149b5e25fSSatish Balay         ierr = PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
87249b5e25fSSatish Balay       }
87349b5e25fSSatish Balay     }
87449b5e25fSSatish Balay     ai[i] = ai[i-1] + ailen[i-1];
87549b5e25fSSatish Balay   }
87649b5e25fSSatish Balay   if (mbs) {
87749b5e25fSSatish Balay     fshift += imax[mbs-1] - ailen[mbs-1];
87849b5e25fSSatish Balay     ai[mbs] = ai[mbs-1] + ailen[mbs-1];
87949b5e25fSSatish Balay   }
88049b5e25fSSatish Balay   /* reset ilen and imax for each row */
88149b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
88249b5e25fSSatish Balay     ailen[i] = imax[i] = ai[i+1] - ai[i];
88349b5e25fSSatish Balay   }
88449b5e25fSSatish Balay   a->s_nz = ai[mbs];
88549b5e25fSSatish Balay 
88649b5e25fSSatish Balay   /* diagonals may have moved, so kill the diagonal pointers */
88749b5e25fSSatish Balay   if (fshift && a->diag) {
88849b5e25fSSatish Balay     ierr = PetscFree(a->diag);CHKERRQ(ierr);
88949b5e25fSSatish Balay     PLogObjectMemory(A,-(m+1)*sizeof(int));
89049b5e25fSSatish Balay     a->diag = 0;
89149b5e25fSSatish Balay   }
89249b5e25fSSatish Balay   PLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Matrix size: %d X %d, block size %d; storage space: %d unneeded, %d used\n",
89349b5e25fSSatish Balay            m,a->m,a->bs,fshift*bs2,a->s_nz*bs2);
89449b5e25fSSatish Balay   PLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Number of mallocs during MatSetValues is %d\n",
89549b5e25fSSatish Balay            a->reallocs);
89649b5e25fSSatish Balay   PLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Most nonzeros blocks in any row is %d\n",rmax);
89749b5e25fSSatish Balay   a->reallocs          = 0;
89849b5e25fSSatish Balay   A->info.nz_unneeded  = (PetscReal)fshift*bs2;
89949b5e25fSSatish Balay 
90049b5e25fSSatish Balay   PetscFunctionReturn(0);
90149b5e25fSSatish Balay }
90249b5e25fSSatish Balay 
90349b5e25fSSatish Balay /*
90449b5e25fSSatish Balay    This function returns an array of flags which indicate the locations of contiguous
90549b5e25fSSatish Balay    blocks that should be zeroed. for eg: if bs = 3  and is = [0,1,2,3,5,6,7,8,9]
90649b5e25fSSatish Balay    then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)]
90749b5e25fSSatish Balay    Assume: sizes should be long enough to hold all the values.
90849b5e25fSSatish Balay */
90949b5e25fSSatish Balay #undef __FUNC__
91049b5e25fSSatish Balay #define __FUNC__ "MatZeroRows_SeqSBAIJ_Check_Blocks"
91149b5e25fSSatish Balay static int MatZeroRows_SeqSBAIJ_Check_Blocks(int idx[],int n,int bs,int sizes[], int *bs_max)
91249b5e25fSSatish Balay {
91349b5e25fSSatish Balay   int        i,j,k,row;
91449b5e25fSSatish Balay   PetscTruth flg;
91549b5e25fSSatish Balay 
91649b5e25fSSatish Balay   PetscFunctionBegin;
91749b5e25fSSatish Balay   for (i=0,j=0; i<n; j++) {
91849b5e25fSSatish Balay     row = idx[i];
91949b5e25fSSatish Balay     if (row%bs!=0) { /* Not the begining of a block */
92049b5e25fSSatish Balay       sizes[j] = 1;
92149b5e25fSSatish Balay       i++;
92249b5e25fSSatish Balay     } else if (i+bs > n) { /* Beginning of a block, but complete block doesn't exist (at idx end) */
92349b5e25fSSatish Balay       sizes[j] = 1;         /* Also makes sure atleast 'bs' values exist for next else */
92449b5e25fSSatish Balay       i++;
92549b5e25fSSatish Balay     } else { /* Begining of the block, so check if the complete block exists */
92649b5e25fSSatish Balay       flg = PETSC_TRUE;
92749b5e25fSSatish Balay       for (k=1; k<bs; k++) {
92849b5e25fSSatish Balay         if (row+k != idx[i+k]) { /* break in the block */
92949b5e25fSSatish Balay           flg = PETSC_FALSE;
93049b5e25fSSatish Balay           break;
93149b5e25fSSatish Balay         }
93249b5e25fSSatish Balay       }
93349b5e25fSSatish Balay       if (flg == PETSC_TRUE) { /* No break in the bs */
93449b5e25fSSatish Balay         sizes[j] = bs;
93549b5e25fSSatish Balay         i+= bs;
93649b5e25fSSatish Balay       } else {
93749b5e25fSSatish Balay         sizes[j] = 1;
93849b5e25fSSatish Balay         i++;
93949b5e25fSSatish Balay       }
94049b5e25fSSatish Balay     }
94149b5e25fSSatish Balay   }
94249b5e25fSSatish Balay   *bs_max = j;
94349b5e25fSSatish Balay   PetscFunctionReturn(0);
94449b5e25fSSatish Balay }
94549b5e25fSSatish Balay 
94649b5e25fSSatish Balay #undef __FUNC__
94749b5e25fSSatish Balay #define __FUNC__ "MatZeroRows_SeqSBAIJ"
94849b5e25fSSatish Balay int MatZeroRows_SeqSBAIJ(Mat A,IS is,Scalar *diag)
94949b5e25fSSatish Balay {
95049b5e25fSSatish Balay   Mat_SeqSBAIJ *sbaij=(Mat_SeqSBAIJ*)A->data;
95149b5e25fSSatish Balay   int         ierr,i,j,k,count,is_n,*is_idx,*rows;
95249b5e25fSSatish Balay   int         bs=sbaij->bs,bs2=sbaij->bs2,*sizes,row,bs_max;
95349b5e25fSSatish Balay   Scalar      zero = 0.0;
95449b5e25fSSatish Balay   MatScalar   *aa;
95549b5e25fSSatish Balay 
95649b5e25fSSatish Balay   PetscFunctionBegin;
95749b5e25fSSatish Balay   /* Make a copy of the IS and  sort it */
95849b5e25fSSatish Balay   ierr = ISGetSize(is,&is_n);CHKERRQ(ierr);
95949b5e25fSSatish Balay   ierr = ISGetIndices(is,&is_idx);CHKERRQ(ierr);
96049b5e25fSSatish Balay 
96149b5e25fSSatish Balay   /* allocate memory for rows,sizes */
96249b5e25fSSatish Balay   rows  = (int*)PetscMalloc((3*is_n+1)*sizeof(int));CHKPTRQ(rows);
96349b5e25fSSatish Balay   sizes = rows + is_n;
96449b5e25fSSatish Balay 
96549b5e25fSSatish Balay   /* initialize copy IS values to rows, and sort them */
96649b5e25fSSatish Balay   for (i=0; i<is_n; i++) { rows[i] = is_idx[i]; }
96749b5e25fSSatish Balay   ierr = PetscSortInt(is_n,rows);CHKERRQ(ierr);
96849b5e25fSSatish Balay   if (sbaij->keepzeroedrows) { /* do not change nonzero structure */
96949b5e25fSSatish Balay     for (i=0; i<is_n; i++) { sizes[i] = 1; } /* sizes: size of blocks, = 1 or bs */
97049b5e25fSSatish Balay     bs_max = is_n;            /* bs_max: num. of contiguous block row in the row */
97149b5e25fSSatish Balay   } else {
97249b5e25fSSatish Balay     ierr = MatZeroRows_SeqSBAIJ_Check_Blocks(rows,is_n,bs,sizes,&bs_max);CHKERRQ(ierr);
97349b5e25fSSatish Balay   }
97449b5e25fSSatish Balay   ierr = ISRestoreIndices(is,&is_idx);CHKERRQ(ierr);
97549b5e25fSSatish Balay 
97649b5e25fSSatish Balay   for (i=0,j=0; i<bs_max; j+=sizes[i],i++) {
97749b5e25fSSatish Balay     row   = rows[j];                  /* row to be zeroed */
97849b5e25fSSatish Balay     if (row < 0 || row > sbaij->m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,0,"row %d out of range",row);
97949b5e25fSSatish Balay     count = (sbaij->i[row/bs +1] - sbaij->i[row/bs])*bs; /* num. of elements in the row */
98049b5e25fSSatish Balay     aa    = sbaij->a + sbaij->i[row/bs]*bs2 + (row%bs);
98149b5e25fSSatish Balay     if (sizes[i] == bs && !sbaij->keepzeroedrows) {
98249b5e25fSSatish Balay       if (diag) {
98349b5e25fSSatish Balay         if (sbaij->ilen[row/bs] > 0) {
98449b5e25fSSatish Balay           sbaij->ilen[row/bs] = 1;
98549b5e25fSSatish Balay           sbaij->j[sbaij->i[row/bs]] = row/bs;
98649b5e25fSSatish Balay           ierr = PetscMemzero(aa,count*bs*sizeof(MatScalar));CHKERRQ(ierr);
98749b5e25fSSatish Balay         }
98849b5e25fSSatish Balay         /* Now insert all the diagoanl values for this bs */
98949b5e25fSSatish Balay         for (k=0; k<bs; k++) {
99049b5e25fSSatish Balay           ierr = (*A->ops->setvalues)(A,1,rows+j+k,1,rows+j+k,diag,INSERT_VALUES);CHKERRQ(ierr);
99149b5e25fSSatish Balay         }
99249b5e25fSSatish Balay       } else { /* (!diag) */
99349b5e25fSSatish Balay         sbaij->ilen[row/bs] = 0;
99449b5e25fSSatish Balay       } /* end (!diag) */
99549b5e25fSSatish Balay     } else { /* (sizes[i] != bs), broken block */
99649b5e25fSSatish Balay #if defined (PETSC_USE_DEBUG)
99749b5e25fSSatish Balay       if (sizes[i] != 1) SETERRQ(1,0,"Internal Error. Value should be 1");
99849b5e25fSSatish Balay #endif
99949b5e25fSSatish Balay       for (k=0; k<count; k++) {
100049b5e25fSSatish Balay         aa[0] = zero;
100149b5e25fSSatish Balay         aa+=bs;
100249b5e25fSSatish Balay       }
100349b5e25fSSatish Balay       if (diag) {
100449b5e25fSSatish Balay         ierr = (*A->ops->setvalues)(A,1,rows+j,1,rows+j,diag,INSERT_VALUES);CHKERRQ(ierr);
100549b5e25fSSatish Balay       }
100649b5e25fSSatish Balay     }
100749b5e25fSSatish Balay   }
100849b5e25fSSatish Balay 
100949b5e25fSSatish Balay   ierr = PetscFree(rows);CHKERRQ(ierr);
101049b5e25fSSatish Balay   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
101149b5e25fSSatish Balay   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
101249b5e25fSSatish Balay   PetscFunctionReturn(0);
101349b5e25fSSatish Balay }
101449b5e25fSSatish Balay 
101549b5e25fSSatish Balay /* Only add/insert a(i,j) with i<=j (blocks).
101649b5e25fSSatish Balay    Any a(i,j) with i>j input by user is ingored.
101749b5e25fSSatish Balay */
101849b5e25fSSatish Balay 
101949b5e25fSSatish Balay #undef __FUNC__
102049b5e25fSSatish Balay #define __FUNC__ "MatSetValues_SeqSBAIJ"
102149b5e25fSSatish Balay int MatSetValues_SeqSBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is)
102249b5e25fSSatish Balay {
102349b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
102449b5e25fSSatish Balay   int         *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted=a->sorted;
102549b5e25fSSatish Balay   int         *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented;
102649b5e25fSSatish Balay   int         *aj=a->j,nonew=a->nonew,bs=a->bs,brow,bcol;
102749b5e25fSSatish Balay   int         ridx,cidx,bs2=a->bs2,ierr;
102849b5e25fSSatish Balay   MatScalar   *ap,value,*aa=a->a,*bap;
102949b5e25fSSatish Balay 
103049b5e25fSSatish Balay   PetscFunctionBegin;
103149b5e25fSSatish Balay 
103249b5e25fSSatish Balay   for (k=0; k<m; k++) { /* loop over added rows */
103349b5e25fSSatish Balay     row  = im[k];       /* row number */
103449b5e25fSSatish Balay     brow = row/bs;      /* block row number */
103549b5e25fSSatish Balay     if (row < 0) continue;
103649b5e25fSSatish Balay #if defined(PETSC_USE_BOPT_g)
103749b5e25fSSatish Balay     if (row >= a->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large: row %d max %d",row,a->m);
103849b5e25fSSatish Balay #endif
103949b5e25fSSatish Balay     rp   = aj + ai[brow]; /*ptr to beginning of column value of the row block*/
104049b5e25fSSatish Balay     ap   = aa + bs2*ai[brow]; /*ptr to beginning of element value of the row block*/
104149b5e25fSSatish Balay     rmax = imax[brow];  /* maximum space allocated for this row */
104249b5e25fSSatish Balay     nrow = ailen[brow]; /* actual length of this row */
104349b5e25fSSatish Balay     low  = 0;
104449b5e25fSSatish Balay 
104549b5e25fSSatish Balay     for (l=0; l<n; l++) { /* loop over added columns */
104649b5e25fSSatish Balay       if (in[l] < 0) continue;
104749b5e25fSSatish Balay #if defined(PETSC_USE_BOPT_g)
104849b5e25fSSatish Balay       if (in[l] >= a->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large: col %d max %d",in[l],a->m);
104949b5e25fSSatish Balay #endif
105049b5e25fSSatish Balay       col = in[l];
105149b5e25fSSatish Balay       bcol = col/bs;              /* block col number */
105249b5e25fSSatish Balay 
105349b5e25fSSatish Balay       if (brow <= bcol){
105449b5e25fSSatish Balay         ridx = row % bs; cidx = col % bs; /*row and col index inside the block */
105549b5e25fSSatish Balay         /* if ((brow==bcol && ridx<=cidx) || (brow<bcol)){ */
105649b5e25fSSatish Balay           /* element value a(k,l) */
105749b5e25fSSatish Balay           if (roworiented) {
105849b5e25fSSatish Balay             value = v[l + k*n];
105949b5e25fSSatish Balay           } else {
106049b5e25fSSatish Balay             value = v[k + l*m];
106149b5e25fSSatish Balay           }
106249b5e25fSSatish Balay 
106349b5e25fSSatish Balay           /* move pointer bap to a(k,l) quickly and add/insert value */
106449b5e25fSSatish Balay           if (!sorted) low = 0; high = nrow;
106549b5e25fSSatish Balay           while (high-low > 7) {
106649b5e25fSSatish Balay             t = (low+high)/2;
106749b5e25fSSatish Balay             if (rp[t] > bcol) high = t;
106849b5e25fSSatish Balay             else              low  = t;
106949b5e25fSSatish Balay           }
107049b5e25fSSatish Balay           for (i=low; i<high; i++) {
107149b5e25fSSatish Balay             /* printf("The loop of i=low.., rp[%d]=%d\n",i,rp[i]); */
107249b5e25fSSatish Balay             if (rp[i] > bcol) break;
107349b5e25fSSatish Balay             if (rp[i] == bcol) {
107449b5e25fSSatish Balay               bap  = ap +  bs2*i + bs*cidx + ridx;
107549b5e25fSSatish Balay               if (is == ADD_VALUES) *bap += value;
107649b5e25fSSatish Balay               else                  *bap  = value;
107749b5e25fSSatish Balay               goto noinsert1;
107849b5e25fSSatish Balay             }
107949b5e25fSSatish Balay           }
108049b5e25fSSatish Balay 
108149b5e25fSSatish Balay       if (nonew == 1) goto noinsert1;
108249b5e25fSSatish Balay       else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix");
108349b5e25fSSatish Balay       if (nrow >= rmax) {
108449b5e25fSSatish Balay         /* there is no extra room in row, therefore enlarge */
108549b5e25fSSatish Balay         int       new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j;
108649b5e25fSSatish Balay         MatScalar *new_a;
108749b5e25fSSatish Balay 
108849b5e25fSSatish Balay         if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix");
108949b5e25fSSatish Balay 
109049b5e25fSSatish Balay         /* Malloc new storage space */
109149b5e25fSSatish Balay         len     = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int);
109249b5e25fSSatish Balay         new_a   = (MatScalar*)PetscMalloc(len);CHKPTRQ(new_a);
109349b5e25fSSatish Balay         new_j   = (int*)(new_a + bs2*new_nz);
109449b5e25fSSatish Balay         new_i   = new_j + new_nz;
109549b5e25fSSatish Balay 
109649b5e25fSSatish Balay         /* copy over old data into new slots */
109749b5e25fSSatish Balay         for (ii=0; ii<brow+1; ii++) {new_i[ii] = ai[ii];}
109849b5e25fSSatish Balay         for (ii=brow+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;}
109949b5e25fSSatish Balay         ierr = PetscMemcpy(new_j,aj,(ai[brow]+nrow)*sizeof(int));CHKERRQ(ierr);
110049b5e25fSSatish Balay         len  = (new_nz - CHUNKSIZE - ai[brow] - nrow);
110149b5e25fSSatish Balay         ierr = PetscMemcpy(new_j+ai[brow]+nrow+CHUNKSIZE,aj+ai[brow]+nrow,len*sizeof(int));CHKERRQ(ierr);
110249b5e25fSSatish Balay         ierr = PetscMemcpy(new_a,aa,(ai[brow]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr);
110349b5e25fSSatish Balay         ierr = PetscMemzero(new_a+bs2*(ai[brow]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr);
110449b5e25fSSatish Balay         ierr = PetscMemcpy(new_a+bs2*(ai[brow]+nrow+CHUNKSIZE),aa+bs2*(ai[brow]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr);
110549b5e25fSSatish Balay         /* free up old matrix storage */
110649b5e25fSSatish Balay         ierr = PetscFree(a->a);CHKERRQ(ierr);
110749b5e25fSSatish Balay         if (!a->singlemalloc) {
110849b5e25fSSatish Balay           ierr = PetscFree(a->i);CHKERRQ(ierr);
110949b5e25fSSatish Balay           ierr = PetscFree(a->j);CHKERRQ(ierr);
111049b5e25fSSatish Balay         }
111149b5e25fSSatish Balay         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
111249b5e25fSSatish Balay         a->singlemalloc = PETSC_TRUE;
111349b5e25fSSatish Balay 
111449b5e25fSSatish Balay         rp   = aj + ai[brow]; ap = aa + bs2*ai[brow];
111549b5e25fSSatish Balay         rmax = imax[brow] = imax[brow] + CHUNKSIZE;
111649b5e25fSSatish Balay         PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar)));
111749b5e25fSSatish Balay         a->s_maxnz += bs2*CHUNKSIZE;
111849b5e25fSSatish Balay         a->reallocs++;
111949b5e25fSSatish Balay         a->s_nz++;
112049b5e25fSSatish Balay       }
112149b5e25fSSatish Balay 
112249b5e25fSSatish Balay       N = nrow++ - 1;
112349b5e25fSSatish Balay       /* shift up all the later entries in this row */
112449b5e25fSSatish Balay       for (ii=N; ii>=i; ii--) {
112549b5e25fSSatish Balay         rp[ii+1] = rp[ii];
112649b5e25fSSatish Balay         ierr     = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr);
112749b5e25fSSatish Balay       }
112849b5e25fSSatish Balay       if (N>=i) {
112949b5e25fSSatish Balay         ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr);
113049b5e25fSSatish Balay       }
113149b5e25fSSatish Balay       rp[i]                      = bcol;
113249b5e25fSSatish Balay       ap[bs2*i + bs*cidx + ridx] = value;
113349b5e25fSSatish Balay       noinsert1:;
113449b5e25fSSatish Balay       low = i;
113549b5e25fSSatish Balay       /* } */
113649b5e25fSSatish Balay     } /* end of if .. if..  */
113749b5e25fSSatish Balay     }                     /* end of loop over added columns */
113849b5e25fSSatish Balay     ailen[brow] = nrow;
113949b5e25fSSatish Balay   }                       /* end of loop over added rows */
114049b5e25fSSatish Balay 
114149b5e25fSSatish Balay   PetscFunctionReturn(0);
114249b5e25fSSatish Balay }
114349b5e25fSSatish Balay 
1144915354c9SHong Zhang extern int MatCholeskyFactorSymbolic_SeqSBAIJ(Mat,IS,PetscReal,Mat*);
1145915354c9SHong Zhang extern int MatCholeskyFactor_SeqSBAIJ(Mat,IS,PetscReal);
114649b5e25fSSatish Balay extern int MatIncreaseOverlap_SeqSBAIJ(Mat,int,IS*,int);
114749b5e25fSSatish Balay extern int MatGetSubMatrix_SeqSBAIJ(Mat,IS,IS,int,MatReuse,Mat*);
114849b5e25fSSatish Balay extern int MatGetSubMatrices_SeqSBAIJ(Mat,int,IS*,IS*,MatReuse,Mat**);
114949b5e25fSSatish Balay extern int MatMultTranspose_SeqSBAIJ(Mat,Vec,Vec);
115049b5e25fSSatish Balay extern int MatMultTransposeAdd_SeqSBAIJ(Mat,Vec,Vec,Vec);
115149b5e25fSSatish Balay extern int MatScale_SeqSBAIJ(Scalar*,Mat);
115249b5e25fSSatish Balay extern int MatNorm_SeqSBAIJ(Mat,NormType,PetscReal *);
115349b5e25fSSatish Balay extern int MatEqual_SeqSBAIJ(Mat,Mat,PetscTruth*);
115449b5e25fSSatish Balay extern int MatGetDiagonal_SeqSBAIJ(Mat,Vec);
115549b5e25fSSatish Balay extern int MatDiagonalScale_SeqSBAIJ(Mat,Vec,Vec);
115649b5e25fSSatish Balay extern int MatGetInfo_SeqSBAIJ(Mat,MatInfoType,MatInfo *);
115749b5e25fSSatish Balay extern int MatZeroEntries_SeqSBAIJ(Mat);
115849b5e25fSSatish Balay 
115949b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_N(Mat,Vec,Vec);
116049b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_1(Mat,Vec,Vec);
116149b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_2(Mat,Vec,Vec);
116249b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_3(Mat,Vec,Vec);
116349b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_4(Mat,Vec,Vec);
116449b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_5(Mat,Vec,Vec);
116549b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_6(Mat,Vec,Vec);
116649b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_7(Mat,Vec,Vec);
116749b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_7(Mat,Vec,Vec);
116849b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_6(Mat,Vec,Vec);
116949b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_5(Mat,Vec,Vec);
117049b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_4(Mat,Vec,Vec);
117149b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_3(Mat,Vec,Vec);
117249b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_2(Mat,Vec,Vec);
117349b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_1(Mat,Vec,Vec);
117449b5e25fSSatish Balay 
11755f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat,Mat*);
11765f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_1(Mat,Mat*);
11775f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_2(Mat,Mat*);
11785f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_3(Mat,Mat*);
11795f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_4(Mat,Mat*);
11805f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_5(Mat,Mat*);
11815f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_6(Mat,Mat*);
118249b5e25fSSatish Balay 
118349b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_1(Mat,Vec,Vec);
118449b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_2(Mat,Vec,Vec);
118549b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_3(Mat,Vec,Vec);
118649b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_4(Mat,Vec,Vec);
118749b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_5(Mat,Vec,Vec);
118849b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_6(Mat,Vec,Vec);
118949b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_7(Mat,Vec,Vec);
119049b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_N(Mat,Vec,Vec);
119149b5e25fSSatish Balay 
119249b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_1(Mat,Vec,Vec,Vec);
119349b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_2(Mat,Vec,Vec,Vec);
119449b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_3(Mat,Vec,Vec,Vec);
119549b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_4(Mat,Vec,Vec,Vec);
119649b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_5(Mat,Vec,Vec,Vec);
119749b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_6(Mat,Vec,Vec,Vec);
119849b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_7(Mat,Vec,Vec,Vec);
119949b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_N(Mat,Vec,Vec,Vec);
120049b5e25fSSatish Balay 
1201950f1e5bSHong Zhang #ifdef MatIncompleteCholeskyFactor
1202950f1e5bSHong Zhang /* This function is modified from MatILUFactor_SeqSBAIJ.
1203950f1e5bSHong Zhang    Needs further work! Don't forget to add the function to the matrix table. */
120449b5e25fSSatish Balay #undef __FUNC__
12054ccecd49SHong Zhang #define __FUNC__ "MatIncompleteCholeskyFactor_SeqSBAIJ"
12064ccecd49SHong Zhang int MatIncompleteCholeskyFactor_SeqSBAIJ(Mat inA,IS row,IS col,MatILUInfo *info)
120749b5e25fSSatish Balay {
12084ccecd49SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)inA->data;
120949b5e25fSSatish Balay   Mat         outA;
121049b5e25fSSatish Balay   int         ierr;
121149b5e25fSSatish Balay   PetscTruth  row_identity,col_identity;
121249b5e25fSSatish Balay 
121349b5e25fSSatish Balay   PetscFunctionBegin;
121449b5e25fSSatish Balay   if (info && info->levels != 0) SETERRQ(PETSC_ERR_SUP,0,"Only levels = 0 supported for in-place ILU");
121549b5e25fSSatish Balay   ierr = ISIdentity(row,&row_identity);CHKERRQ(ierr);
121649b5e25fSSatish Balay   ierr = ISIdentity(col,&col_identity);CHKERRQ(ierr);
121749b5e25fSSatish Balay   if (!row_identity || !col_identity) {
121849b5e25fSSatish Balay     SETERRQ(1,1,"Row and column permutations must be identity for in-place ILU");
121949b5e25fSSatish Balay   }
122049b5e25fSSatish Balay 
122149b5e25fSSatish Balay   outA          = inA;
122249b5e25fSSatish Balay   inA->factor   = FACTOR_LU;
122349b5e25fSSatish Balay 
122449b5e25fSSatish Balay   if (!a->diag) {
122549b5e25fSSatish Balay     ierr = MatMarkDiagonal_SeqBAIJ(inA);CHKERRQ(ierr);
122649b5e25fSSatish Balay   }
122749b5e25fSSatish Balay   /*
122849b5e25fSSatish Balay       Blocksize 2, 3, 4, 5, 6 and 7 have a special faster factorization/solver
122949b5e25fSSatish Balay       for ILU(0) factorization with natural ordering
123049b5e25fSSatish Balay   */
123149b5e25fSSatish Balay   switch (a->bs) {
123249b5e25fSSatish Balay   case 1:
123349b5e25fSSatish Balay     inA->ops->solvetranspose   = MatSolveTranspose_SeqBAIJ_2_NaturalOrdering;
123449b5e25fSSatish Balay     PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering solvetrans BS=1\n");
123549b5e25fSSatish Balay   case 2:
123649b5e25fSSatish Balay     inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering;
123749b5e25fSSatish Balay     inA->ops->solve           = MatSolve_SeqBAIJ_2_NaturalOrdering;
123849b5e25fSSatish Balay     inA->ops->solvetranspose  = MatSolveTranspose_SeqBAIJ_2_NaturalOrdering;
123949b5e25fSSatish Balay     PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=2\n");
124049b5e25fSSatish Balay     break;
124149b5e25fSSatish Balay   case 3:
124249b5e25fSSatish Balay     inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering;
124349b5e25fSSatish Balay     inA->ops->solve           = MatSolve_SeqBAIJ_3_NaturalOrdering;
124449b5e25fSSatish Balay     inA->ops->solvetranspose  = MatSolveTranspose_SeqBAIJ_3_NaturalOrdering;
124549b5e25fSSatish Balay     PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=3\n");
124649b5e25fSSatish Balay     break;
124749b5e25fSSatish Balay   case 4:
124849b5e25fSSatish Balay     inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering;
124949b5e25fSSatish Balay     inA->ops->solve           = MatSolve_SeqBAIJ_4_NaturalOrdering;
125049b5e25fSSatish Balay     inA->ops->solvetranspose  = MatSolveTranspose_SeqBAIJ_4_NaturalOrdering;
125149b5e25fSSatish Balay     PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=4\n");
125249b5e25fSSatish Balay     break;
125349b5e25fSSatish Balay   case 5:
125449b5e25fSSatish Balay     inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5_NaturalOrdering;
125549b5e25fSSatish Balay     inA->ops->solve           = MatSolve_SeqBAIJ_5_NaturalOrdering;
125649b5e25fSSatish Balay     inA->ops->solvetranspose  = MatSolveTranspose_SeqBAIJ_5_NaturalOrdering;
125749b5e25fSSatish Balay     PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=5\n");
125849b5e25fSSatish Balay     break;
125949b5e25fSSatish Balay   case 6:
126049b5e25fSSatish Balay     inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6_NaturalOrdering;
126149b5e25fSSatish Balay     inA->ops->solve           = MatSolve_SeqBAIJ_6_NaturalOrdering;
126249b5e25fSSatish Balay     inA->ops->solvetranspose  = MatSolveTranspose_SeqBAIJ_6_NaturalOrdering;
126349b5e25fSSatish Balay     PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=6\n");
126449b5e25fSSatish Balay     break;
126549b5e25fSSatish Balay   case 7:
126649b5e25fSSatish Balay     inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_7_NaturalOrdering;
126749b5e25fSSatish Balay     inA->ops->solvetranspose  = MatSolveTranspose_SeqBAIJ_7_NaturalOrdering;
126849b5e25fSSatish Balay     inA->ops->solve           = MatSolve_SeqBAIJ_7_NaturalOrdering;
126949b5e25fSSatish Balay     PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=7\n");
127049b5e25fSSatish Balay     break;
127149b5e25fSSatish Balay   default:
127249b5e25fSSatish Balay     a->row        = row;
127349b5e25fSSatish Balay     a->col        = col;
127449b5e25fSSatish Balay     ierr          = PetscObjectReference((PetscObject)row);CHKERRQ(ierr);
127549b5e25fSSatish Balay     ierr          = PetscObjectReference((PetscObject)col);CHKERRQ(ierr);
127649b5e25fSSatish Balay 
127749b5e25fSSatish Balay     /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */
127849b5e25fSSatish Balay     ierr = ISInvertPermutation(col,PETSC_DECIDE, &(a->icol));CHKERRQ(ierr);
127949b5e25fSSatish Balay     PLogObjectParent(inA,a->icol);
128049b5e25fSSatish Balay 
128149b5e25fSSatish Balay     if (!a->solve_work) {
128249b5e25fSSatish Balay       a->solve_work = (Scalar*)PetscMalloc((a->m+a->bs)*sizeof(Scalar));CHKPTRQ(a->solve_work);
128349b5e25fSSatish Balay       PLogObjectMemory(inA,(a->m+a->bs)*sizeof(Scalar));
128449b5e25fSSatish Balay     }
128549b5e25fSSatish Balay   }
128649b5e25fSSatish Balay 
128749b5e25fSSatish Balay   ierr = MatLUFactorNumeric(inA,&outA);CHKERRQ(ierr);
128849b5e25fSSatish Balay 
128949b5e25fSSatish Balay   PetscFunctionReturn(0);
129049b5e25fSSatish Balay }
1291950f1e5bSHong Zhang #endif
1292950f1e5bSHong Zhang 
129349b5e25fSSatish Balay #undef __FUNC__
129449b5e25fSSatish Balay #define __FUNC__ "MatPrintHelp_SeqSBAIJ"
129549b5e25fSSatish Balay int MatPrintHelp_SeqSBAIJ(Mat A)
129649b5e25fSSatish Balay {
129749b5e25fSSatish Balay   static PetscTruth called = PETSC_FALSE;
129849b5e25fSSatish Balay   MPI_Comm          comm = A->comm;
129949b5e25fSSatish Balay   int               ierr;
130049b5e25fSSatish Balay 
130149b5e25fSSatish Balay   PetscFunctionBegin;
130249b5e25fSSatish Balay   if (called) {PetscFunctionReturn(0);} else called = PETSC_TRUE;
130349b5e25fSSatish Balay   ierr = (*PetscHelpPrintf)(comm," Options for MATSEQSBAIJ and MATMPISBAIJ matrix formats (the defaults):\n");CHKERRQ(ierr);
130449b5e25fSSatish Balay   ierr = (*PetscHelpPrintf)(comm,"  -mat_block_size <block_size>\n");CHKERRQ(ierr);
130549b5e25fSSatish Balay   PetscFunctionReturn(0);
130649b5e25fSSatish Balay }
130749b5e25fSSatish Balay 
130849b5e25fSSatish Balay EXTERN_C_BEGIN
130949b5e25fSSatish Balay #undef __FUNC__
131049b5e25fSSatish Balay #define __FUNC__ "MatSeqSBAIJSetColumnIndices_SeqSBAIJ"
131149b5e25fSSatish Balay int MatSeqSBAIJSetColumnIndices_SeqSBAIJ(Mat mat,int *indices)
131249b5e25fSSatish Balay {
131349b5e25fSSatish Balay   Mat_SeqBAIJ *baij = (Mat_SeqBAIJ *)mat->data;
131449b5e25fSSatish Balay   int         i,nz,n;
131549b5e25fSSatish Balay 
131649b5e25fSSatish Balay   PetscFunctionBegin;
131749b5e25fSSatish Balay   nz = baij->maxnz;
131849b5e25fSSatish Balay   n  = baij->n;
131949b5e25fSSatish Balay   for (i=0; i<nz; i++) {
132049b5e25fSSatish Balay     baij->j[i] = indices[i];
132149b5e25fSSatish Balay   }
132249b5e25fSSatish Balay   baij->nz = nz;
132349b5e25fSSatish Balay   for (i=0; i<n; i++) {
132449b5e25fSSatish Balay     baij->ilen[i] = baij->imax[i];
132549b5e25fSSatish Balay   }
132649b5e25fSSatish Balay 
132749b5e25fSSatish Balay   PetscFunctionReturn(0);
132849b5e25fSSatish Balay }
132949b5e25fSSatish Balay EXTERN_C_END
133049b5e25fSSatish Balay 
133149b5e25fSSatish Balay #undef __FUNC__
133249b5e25fSSatish Balay #define __FUNC__ "MatSeqSBAIJSetColumnIndices"
133349b5e25fSSatish Balay /*@
133449b5e25fSSatish Balay     MatSeqBAIJSetColumnIndices - Set the column indices for all the rows
133549b5e25fSSatish Balay        in the matrix.
133649b5e25fSSatish Balay 
133749b5e25fSSatish Balay   Input Parameters:
133849b5e25fSSatish Balay +  mat - the SeqBAIJ matrix
133949b5e25fSSatish Balay -  indices - the column indices
134049b5e25fSSatish Balay 
134149b5e25fSSatish Balay   Level: advanced
134249b5e25fSSatish Balay 
134349b5e25fSSatish Balay   Notes:
134449b5e25fSSatish Balay     This can be called if you have precomputed the nonzero structure of the
134549b5e25fSSatish Balay   matrix and want to provide it to the matrix object to improve the performance
134649b5e25fSSatish Balay   of the MatSetValues() operation.
134749b5e25fSSatish Balay 
134849b5e25fSSatish Balay     You MUST have set the correct numbers of nonzeros per row in the call to
134949b5e25fSSatish Balay   MatCreateSeqBAIJ().
135049b5e25fSSatish Balay 
1351*ab9f2c04SSatish Balay     MUST be called before any calls to MatSetValues()
135249b5e25fSSatish Balay 
1353*ab9f2c04SSatish Balay .seealso: MatCreateSeqSBAIJ
135449b5e25fSSatish Balay @*/
135549b5e25fSSatish Balay int MatSeqSBAIJSetColumnIndices(Mat mat,int *indices)
135649b5e25fSSatish Balay {
135749b5e25fSSatish Balay   int ierr,(*f)(Mat,int *);
135849b5e25fSSatish Balay 
135949b5e25fSSatish Balay   PetscFunctionBegin;
136049b5e25fSSatish Balay   PetscValidHeaderSpecific(mat,MAT_COOKIE);
13614afc71dfSHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqSBAIJSetColumnIndices_C",(void **)&f);CHKERRQ(ierr);
136249b5e25fSSatish Balay   if (f) {
136349b5e25fSSatish Balay     ierr = (*f)(mat,indices);CHKERRQ(ierr);
136449b5e25fSSatish Balay   } else {
136549b5e25fSSatish Balay     SETERRQ(1,1,"Wrong type of matrix to set column indices");
136649b5e25fSSatish Balay   }
136749b5e25fSSatish Balay   PetscFunctionReturn(0);
136849b5e25fSSatish Balay }
136949b5e25fSSatish Balay 
137049b5e25fSSatish Balay /* -------------------------------------------------------------------*/
137149b5e25fSSatish Balay static struct _MatOps MatOps_Values = {MatSetValues_SeqSBAIJ,
137249b5e25fSSatish Balay        MatGetRow_SeqSBAIJ,
137349b5e25fSSatish Balay        MatRestoreRow_SeqSBAIJ,
137449b5e25fSSatish Balay        MatMult_SeqSBAIJ_N,
137549b5e25fSSatish Balay        MatMultAdd_SeqSBAIJ_N,
137649b5e25fSSatish Balay        MatMultTranspose_SeqSBAIJ,
137749b5e25fSSatish Balay        MatMultTransposeAdd_SeqSBAIJ,
137849b5e25fSSatish Balay        MatSolve_SeqSBAIJ_N,
137949b5e25fSSatish Balay        0,
138049b5e25fSSatish Balay        0,
138149b5e25fSSatish Balay        0,
138249b5e25fSSatish Balay        0,
13835f4ef2deSHong Zhang        MatCholeskyFactor_SeqSBAIJ,
138449b5e25fSSatish Balay        0,
138549b5e25fSSatish Balay        MatTranspose_SeqSBAIJ,
138649b5e25fSSatish Balay        MatGetInfo_SeqSBAIJ,
138749b5e25fSSatish Balay        MatEqual_SeqSBAIJ,
138849b5e25fSSatish Balay        MatGetDiagonal_SeqSBAIJ,
138949b5e25fSSatish Balay        MatDiagonalScale_SeqSBAIJ,
139049b5e25fSSatish Balay        MatNorm_SeqSBAIJ,
139149b5e25fSSatish Balay        0,
139249b5e25fSSatish Balay        MatAssemblyEnd_SeqSBAIJ,
139349b5e25fSSatish Balay        0,
139449b5e25fSSatish Balay        MatSetOption_SeqSBAIJ,
139549b5e25fSSatish Balay        MatZeroEntries_SeqSBAIJ,
139649b5e25fSSatish Balay        MatZeroRows_SeqSBAIJ,
139749b5e25fSSatish Balay        0,
139849b5e25fSSatish Balay        0,
13995f4ef2deSHong Zhang        MatCholeskyFactorSymbolic_SeqSBAIJ,
14005f4ef2deSHong Zhang        MatCholeskyFactorNumeric_SeqSBAIJ_N,
140149b5e25fSSatish Balay        MatGetSize_SeqSBAIJ,
14024ccecd49SHong Zhang        MatGetSize_SeqSBAIJ,
140349b5e25fSSatish Balay        MatGetOwnershipRange_SeqSBAIJ,
140449b5e25fSSatish Balay        0,
14055f4ef2deSHong Zhang        MatIncompleteCholeskyFactorSymbolic_SeqSBAIJ,
140649b5e25fSSatish Balay        0,
140749b5e25fSSatish Balay        0,
140849b5e25fSSatish Balay        MatDuplicate_SeqSBAIJ,
140949b5e25fSSatish Balay        0,
141049b5e25fSSatish Balay        0,
141149b5e25fSSatish Balay        0,
1412950f1e5bSHong Zhang        0,
141349b5e25fSSatish Balay        0,
141449b5e25fSSatish Balay        MatGetSubMatrices_SeqSBAIJ,
141549b5e25fSSatish Balay        MatIncreaseOverlap_SeqSBAIJ,
141649b5e25fSSatish Balay        MatGetValues_SeqSBAIJ,
141749b5e25fSSatish Balay        0,
141849b5e25fSSatish Balay        MatPrintHelp_SeqSBAIJ,
141949b5e25fSSatish Balay        MatScale_SeqSBAIJ,
142049b5e25fSSatish Balay        0,
142149b5e25fSSatish Balay        0,
142249b5e25fSSatish Balay        0,
142349b5e25fSSatish Balay        MatGetBlockSize_SeqSBAIJ,
142449b5e25fSSatish Balay        MatGetRowIJ_SeqSBAIJ,
142549b5e25fSSatish Balay        MatRestoreRowIJ_SeqSBAIJ,
142649b5e25fSSatish Balay        0,
142749b5e25fSSatish Balay        0,
142849b5e25fSSatish Balay        0,
142949b5e25fSSatish Balay        0,
143049b5e25fSSatish Balay        0,
143149b5e25fSSatish Balay        0,
143249b5e25fSSatish Balay        MatSetValuesBlocked_SeqSBAIJ,
143349b5e25fSSatish Balay        MatGetSubMatrix_SeqSBAIJ,
143449b5e25fSSatish Balay        0,
143549b5e25fSSatish Balay        0,
143649b5e25fSSatish Balay        MatGetMaps_Petsc};
143749b5e25fSSatish Balay 
143849b5e25fSSatish Balay EXTERN_C_BEGIN
143949b5e25fSSatish Balay #undef __FUNC__
144049b5e25fSSatish Balay #define __FUNC__ "MatStoreValues_SeqSBAIJ"
144149b5e25fSSatish Balay int MatStoreValues_SeqSBAIJ(Mat mat)
144249b5e25fSSatish Balay {
14434afc71dfSHong Zhang   Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data;
144449b5e25fSSatish Balay   int         nz = aij->i[aij->m]*aij->bs*aij->bs2;
144549b5e25fSSatish Balay   int         ierr;
144649b5e25fSSatish Balay 
144749b5e25fSSatish Balay   PetscFunctionBegin;
144849b5e25fSSatish Balay   if (aij->nonew != 1) {
144949b5e25fSSatish Balay     SETERRQ(1,1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
145049b5e25fSSatish Balay   }
145149b5e25fSSatish Balay 
145249b5e25fSSatish Balay   /* allocate space for values if not already there */
145349b5e25fSSatish Balay   if (!aij->saved_values) {
145449b5e25fSSatish Balay     aij->saved_values = (Scalar*)PetscMalloc(nz*sizeof(Scalar));CHKPTRQ(aij->saved_values);
145549b5e25fSSatish Balay   }
145649b5e25fSSatish Balay 
145749b5e25fSSatish Balay   /* copy values over */
145849b5e25fSSatish Balay   ierr = PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(Scalar));CHKERRQ(ierr);
145949b5e25fSSatish Balay   PetscFunctionReturn(0);
146049b5e25fSSatish Balay }
146149b5e25fSSatish Balay EXTERN_C_END
146249b5e25fSSatish Balay 
146349b5e25fSSatish Balay EXTERN_C_BEGIN
146449b5e25fSSatish Balay #undef __FUNC__
146549b5e25fSSatish Balay #define __FUNC__ "MatRetrieveValues_SeqSBAIJ"
146649b5e25fSSatish Balay int MatRetrieveValues_SeqSBAIJ(Mat mat)
146749b5e25fSSatish Balay {
14684afc71dfSHong Zhang   Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data;
146949b5e25fSSatish Balay   int         nz = aij->i[aij->m]*aij->bs*aij->bs2,ierr;
147049b5e25fSSatish Balay 
147149b5e25fSSatish Balay   PetscFunctionBegin;
147249b5e25fSSatish Balay   if (aij->nonew != 1) {
147349b5e25fSSatish Balay     SETERRQ(1,1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
147449b5e25fSSatish Balay   }
147549b5e25fSSatish Balay   if (!aij->saved_values) {
147649b5e25fSSatish Balay     SETERRQ(1,1,"Must call MatStoreValues(A);first");
147749b5e25fSSatish Balay   }
147849b5e25fSSatish Balay 
147949b5e25fSSatish Balay   /* copy values over */
148049b5e25fSSatish Balay   ierr = PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(Scalar));CHKERRQ(ierr);
148149b5e25fSSatish Balay   PetscFunctionReturn(0);
148249b5e25fSSatish Balay }
148349b5e25fSSatish Balay EXTERN_C_END
148449b5e25fSSatish Balay 
148549b5e25fSSatish Balay #undef __FUNC__
148649b5e25fSSatish Balay #define __FUNC__ "MatCreateSeqSBAIJ"
148749b5e25fSSatish Balay /*@C
148849b5e25fSSatish Balay    MatCreateSeqSBAIJ - Creates a sparse symmetric matrix in block AIJ (block
148949b5e25fSSatish Balay    compressed row) format.  For good matrix assembly performance the
149049b5e25fSSatish Balay    user should preallocate the matrix storage by setting the parameter nz
149149b5e25fSSatish Balay    (or the array nnz).  By setting these parameters accurately, performance
149249b5e25fSSatish Balay    during matrix assembly can be increased by more than a factor of 50.
149349b5e25fSSatish Balay 
149449b5e25fSSatish Balay    Collective on MPI_Comm
149549b5e25fSSatish Balay 
149649b5e25fSSatish Balay    Input Parameters:
149749b5e25fSSatish Balay +  comm - MPI communicator, set to PETSC_COMM_SELF
149849b5e25fSSatish Balay .  bs - size of block
149949b5e25fSSatish Balay .  m - number of rows, or number of columns
150049b5e25fSSatish Balay .  nz - number of block nonzeros per block row (same for all rows)
150149b5e25fSSatish Balay -  nnz - array containing the number of block nonzeros in the various block rows
150249b5e25fSSatish Balay          (possibly different for each block row) or PETSC_NULL
150349b5e25fSSatish Balay 
150449b5e25fSSatish Balay    Output Parameter:
150549b5e25fSSatish Balay .  A - the symmetric matrix
150649b5e25fSSatish Balay 
150749b5e25fSSatish Balay    Options Database Keys:
150849b5e25fSSatish Balay .   -mat_no_unroll - uses code that does not unroll the loops in the
150949b5e25fSSatish Balay                      block calculations (much slower)
151049b5e25fSSatish Balay .    -mat_block_size - size of the blocks to use
151149b5e25fSSatish Balay 
151249b5e25fSSatish Balay    Level: intermediate
151349b5e25fSSatish Balay 
151449b5e25fSSatish Balay    Notes:
151549b5e25fSSatish Balay    The block AIJ format is fully compatible with standard Fortran 77
151649b5e25fSSatish Balay    storage.  That is, the stored row and column indices can begin at
151749b5e25fSSatish Balay    either one (as in Fortran) or zero.  See the users' manual for details.
151849b5e25fSSatish Balay 
151949b5e25fSSatish Balay    Specify the preallocated storage with either nz or nnz (not both).
152049b5e25fSSatish Balay    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
152149b5e25fSSatish Balay    allocation.  For additional details, see the users manual chapter on
152249b5e25fSSatish Balay    matrices.
152349b5e25fSSatish Balay 
152449b5e25fSSatish Balay .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPIBAIJ()
152549b5e25fSSatish Balay @*/
152649b5e25fSSatish Balay int MatCreateSeqSBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,int *nnz,Mat *A)
152749b5e25fSSatish Balay {
152849b5e25fSSatish Balay   Mat         B;
152949b5e25fSSatish Balay   Mat_SeqSBAIJ *b;
15304afc71dfSHong Zhang   int         i,len,ierr,mbs,bs2,size;
153149b5e25fSSatish Balay   PetscTruth  flg;
15324afc71dfSHong Zhang   int         s_nz;
153349b5e25fSSatish Balay 
153449b5e25fSSatish Balay   PetscFunctionBegin;
153549b5e25fSSatish Balay   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
153649b5e25fSSatish Balay   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"Comm must be of size 1");
153749b5e25fSSatish Balay 
153849b5e25fSSatish Balay   ierr = OptionsGetInt(PETSC_NULL,"-mat_block_size",&bs,PETSC_NULL);CHKERRQ(ierr);
153949b5e25fSSatish Balay   mbs  = m/bs;
154049b5e25fSSatish Balay   bs2  = bs*bs;
154149b5e25fSSatish Balay 
154249b5e25fSSatish Balay   if (mbs*bs!=m) {
154349b5e25fSSatish Balay     SETERRQ(PETSC_ERR_ARG_SIZ,0,"Number rows, cols must be divisible by blocksize");
154449b5e25fSSatish Balay   }
154549b5e25fSSatish Balay 
154649b5e25fSSatish Balay   if (nz < -2) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,0,"nz cannot be less than -2: value %d",nz);
154749b5e25fSSatish Balay   if (nnz) {
154849b5e25fSSatish Balay     for (i=0; i<mbs; i++) {
154949b5e25fSSatish Balay       if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,0,"nnz cannot be less than 0: local row %d value %d",i,nnz[i]);
155049b5e25fSSatish Balay     }
155149b5e25fSSatish Balay   }
155249b5e25fSSatish Balay 
155349b5e25fSSatish Balay   *A = 0;
155449b5e25fSSatish Balay   PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQSBAIJ,"Mat",comm,MatDestroy,MatView);
155549b5e25fSSatish Balay   PLogObjectCreate(B);
155649b5e25fSSatish Balay   B->data = (void*)(b = PetscNew(Mat_SeqSBAIJ));CHKPTRQ(b);
155749b5e25fSSatish Balay   ierr    = PetscMemzero(b,sizeof(Mat_SeqSBAIJ));CHKERRQ(ierr);
155849b5e25fSSatish Balay   ierr    = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
155949b5e25fSSatish Balay   ierr    = OptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg);CHKERRQ(ierr);
156049b5e25fSSatish Balay   if (!flg) {
156149b5e25fSSatish Balay     switch (bs) {
156249b5e25fSSatish Balay     case 1:
15634ccecd49SHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1;
156449b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_1;
156549b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_1;
156649b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_1;
156749b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_1;
156849b5e25fSSatish Balay       break;
156949b5e25fSSatish Balay     case 2:
15704ccecd49SHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2;
157149b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_2;
157249b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_2;
157349b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_2;
157449b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_2;
157549b5e25fSSatish Balay       break;
157649b5e25fSSatish Balay     case 3:
15775f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3;
157849b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_3;
157949b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_3;
158049b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_3;
158149b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_3;
158249b5e25fSSatish Balay       break;
158349b5e25fSSatish Balay     case 4:
15845f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4;
158549b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_4;
158649b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_4;
158749b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_4;
158849b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_4;
158949b5e25fSSatish Balay       break;
159049b5e25fSSatish Balay     case 5:
15915f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5;
159249b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_5;
159349b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_5;
159449b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_5;
159549b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_5;
159649b5e25fSSatish Balay       break;
159749b5e25fSSatish Balay     case 6:
15985f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6;
159949b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_6;
160049b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_6;
160149b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_6;
160249b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_6;
160349b5e25fSSatish Balay       break;
160449b5e25fSSatish Balay     case 7:
160549b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_7;
160649b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_7;
160749b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_7;
160849b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_7;
160949b5e25fSSatish Balay       break;
161049b5e25fSSatish Balay     }
161149b5e25fSSatish Balay   }
161249b5e25fSSatish Balay   B->ops->destroy     = MatDestroy_SeqSBAIJ;
161349b5e25fSSatish Balay   B->ops->view        = MatView_SeqSBAIJ;
161449b5e25fSSatish Balay   B->factor           = 0;
161549b5e25fSSatish Balay   B->lupivotthreshold = 1.0;
161649b5e25fSSatish Balay   B->mapping          = 0;
161749b5e25fSSatish Balay   b->row              = 0;
161849b5e25fSSatish Balay   b->col              = 0;
161949b5e25fSSatish Balay   b->icol             = 0;
162049b5e25fSSatish Balay   b->reallocs         = 0;
162149b5e25fSSatish Balay   b->saved_values     = 0;
162249b5e25fSSatish Balay 
16235f4ef2deSHong Zhang   b->m       = m;
16245f4ef2deSHong Zhang   B->m = m; B->M = m;
16255f4ef2deSHong Zhang   b->n       = m;
162649b5e25fSSatish Balay   B->n = m; B->N = m;
162749b5e25fSSatish Balay 
162849b5e25fSSatish Balay   ierr = MapCreateMPI(comm,m,m,&B->rmap);CHKERRQ(ierr);
162949b5e25fSSatish Balay   ierr = MapCreateMPI(comm,m,m,&B->cmap);CHKERRQ(ierr);
163049b5e25fSSatish Balay 
163149b5e25fSSatish Balay   b->mbs     = mbs;
16324afc71dfSHong Zhang   b->nbs     = mbs;
163349b5e25fSSatish Balay   b->imax    = (int*)PetscMalloc((mbs+1)*sizeof(int));CHKPTRQ(b->imax);
163449b5e25fSSatish Balay   if (!nnz) {
163549b5e25fSSatish Balay     if (nz == PETSC_DEFAULT) nz = 5;
163649b5e25fSSatish Balay     else if (nz <= 0)        nz = 1;
163749b5e25fSSatish Balay     for (i=0; i<mbs; i++) {
163849b5e25fSSatish Balay       b->imax[i] = (nz+1)/2;
163949b5e25fSSatish Balay     }
164049b5e25fSSatish Balay     nz = nz*mbs;
164149b5e25fSSatish Balay   } else {
164249b5e25fSSatish Balay     nz = 0;
164349b5e25fSSatish Balay     for (i=0; i<mbs; i++) {b->imax[i] = (nnz[i]+1)/2; nz += nnz[i];}
164449b5e25fSSatish Balay   }
164549b5e25fSSatish Balay   s_nz=(nz+mbs)/2;  /* total diagonal and superdiagonal nonzero blocks */
164649b5e25fSSatish Balay 
164749b5e25fSSatish Balay   /* allocate the matrix space */
164849b5e25fSSatish Balay   len   = s_nz*sizeof(int) + s_nz*bs2*sizeof(MatScalar) + (b->m+1)*sizeof(int);
164949b5e25fSSatish Balay   b->a  = (MatScalar*)PetscMalloc(len);CHKPTRQ(b->a);
165049b5e25fSSatish Balay   ierr = PetscMemzero(b->a,s_nz*bs2*sizeof(MatScalar));CHKERRQ(ierr);
165149b5e25fSSatish Balay   b->j  = (int*)(b->a + s_nz*bs2);
165249b5e25fSSatish Balay   ierr = PetscMemzero(b->j,s_nz*sizeof(int));CHKERRQ(ierr);
165349b5e25fSSatish Balay   b->i  = b->j + s_nz;
165449b5e25fSSatish Balay   b->singlemalloc = PETSC_TRUE;
165549b5e25fSSatish Balay 
165649b5e25fSSatish Balay   /* pointer to beginning of each row */
165749b5e25fSSatish Balay   b->i[0] = 0;
165849b5e25fSSatish Balay   for (i=1; i<mbs+1; i++) {
165949b5e25fSSatish Balay     b->i[i] = b->i[i-1] + b->imax[i-1];
166049b5e25fSSatish Balay   }
166149b5e25fSSatish Balay 
166249b5e25fSSatish Balay 
166349b5e25fSSatish Balay   /* b->ilen will count nonzeros in each block row so far. */
166449b5e25fSSatish Balay   b->ilen = (int*)PetscMalloc((mbs+1)*sizeof(int));
166549b5e25fSSatish Balay   PLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqSBAIJ));
166649b5e25fSSatish Balay   for (i=0; i<mbs; i++) { b->ilen[i] = 0;}
166749b5e25fSSatish Balay 
166849b5e25fSSatish Balay   b->bs               = bs;
166949b5e25fSSatish Balay   b->bs2              = bs2;
167049b5e25fSSatish Balay   b->mbs              = mbs;
167149b5e25fSSatish Balay   b->s_nz               = 0;
167249b5e25fSSatish Balay   b->s_maxnz            = s_nz*bs2;
167349b5e25fSSatish Balay   b->sorted           = PETSC_FALSE;
167449b5e25fSSatish Balay   b->roworiented      = PETSC_TRUE;
167549b5e25fSSatish Balay   b->nonew            = 0;
167649b5e25fSSatish Balay   b->diag             = 0;
167749b5e25fSSatish Balay   b->solve_work       = 0;
167849b5e25fSSatish Balay   b->mult_work        = 0;
167949b5e25fSSatish Balay   b->spptr            = 0;
168049b5e25fSSatish Balay   B->info.nz_unneeded = (PetscReal)b->s_maxnz;
168149b5e25fSSatish Balay   b->keepzeroedrows   = PETSC_FALSE;
168249b5e25fSSatish Balay 
168349b5e25fSSatish Balay   *A = B;
168449b5e25fSSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-help",&flg);CHKERRQ(ierr);
168549b5e25fSSatish Balay   if (flg) {ierr = MatPrintHelp(B);CHKERRQ(ierr); }
168649b5e25fSSatish Balay 
168749b5e25fSSatish Balay   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C",
1688aec82049SSatish Balay                                      "MatStoreValues_SeqSBAIJ",
1689aec82049SSatish Balay                                      (void*)MatStoreValues_SeqSBAIJ);CHKERRQ(ierr);
169049b5e25fSSatish Balay   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C",
1691aec82049SSatish Balay                                      "MatRetrieveValues_SeqSBAIJ",
1692aec82049SSatish Balay                                      (void*)MatRetrieveValues_SeqSBAIJ);CHKERRQ(ierr);
1693aec82049SSatish Balay   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqSBAIJSetColumnIndices_C",
1694aec82049SSatish Balay                                      "MatSeqSBAIJSetColumnIndices_SeqSBAIJ",
1695aec82049SSatish Balay                                      (void*)MatSeqSBAIJSetColumnIndices_SeqSBAIJ);CHKERRQ(ierr);
169649b5e25fSSatish Balay    /* printf("In MatCreateSeqSBAIJ, \n");
169749b5e25fSSatish Balay    for (i=0; i<mbs; i++){
169849b5e25fSSatish Balay      printf("imax[%d]= %d, ilen[%d]= %d\n", i,b->imax[i], i,b->ilen[i]);
169949b5e25fSSatish Balay    }       */
170049b5e25fSSatish Balay 
170149b5e25fSSatish Balay   PetscFunctionReturn(0);
170249b5e25fSSatish Balay }
170349b5e25fSSatish Balay 
170449b5e25fSSatish Balay #undef __FUNC__
170549b5e25fSSatish Balay #define __FUNC__ "MatDuplicate_SeqSBAIJ"
170649b5e25fSSatish Balay int MatDuplicate_SeqSBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B)
170749b5e25fSSatish Balay {
170849b5e25fSSatish Balay   Mat         C;
170949b5e25fSSatish Balay   Mat_SeqSBAIJ *c,*a = (Mat_SeqSBAIJ*)A->data;
171049b5e25fSSatish Balay   int         i,len,mbs = a->mbs,nz = a->s_nz,bs2 =a->bs2,ierr;
171149b5e25fSSatish Balay 
171249b5e25fSSatish Balay   PetscFunctionBegin;
171349b5e25fSSatish Balay   if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,0,"Corrupt matrix");
171449b5e25fSSatish Balay 
171549b5e25fSSatish Balay   *B = 0;
171649b5e25fSSatish Balay   PetscHeaderCreate(C,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQSBAIJ,"Mat",A->comm,MatDestroy,MatView);
171749b5e25fSSatish Balay   PLogObjectCreate(C);
171849b5e25fSSatish Balay   C->data           = (void*)(c = PetscNew(Mat_SeqSBAIJ));CHKPTRQ(c);
171949b5e25fSSatish Balay   ierr              = PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
172049b5e25fSSatish Balay   C->ops->destroy   = MatDestroy_SeqSBAIJ;
172149b5e25fSSatish Balay   C->ops->view      = MatView_SeqSBAIJ;
172249b5e25fSSatish Balay   C->factor         = A->factor;
172349b5e25fSSatish Balay   c->row            = 0;
172449b5e25fSSatish Balay   c->col            = 0;
172549b5e25fSSatish Balay   c->icol           = 0;
172649b5e25fSSatish Balay   c->saved_values   = 0;
172749b5e25fSSatish Balay   c->keepzeroedrows = a->keepzeroedrows;
172849b5e25fSSatish Balay   C->assembled      = PETSC_TRUE;
172949b5e25fSSatish Balay 
173049b5e25fSSatish Balay   c->m = C->m   = a->m;
173149b5e25fSSatish Balay   c->n = C->n   = a->n;
173249b5e25fSSatish Balay   C->M          = a->m;
173349b5e25fSSatish Balay   C->N          = a->n;
173449b5e25fSSatish Balay 
173549b5e25fSSatish Balay   c->bs         = a->bs;
173649b5e25fSSatish Balay   c->bs2        = a->bs2;
173749b5e25fSSatish Balay   c->mbs        = a->mbs;
173849b5e25fSSatish Balay   c->nbs        = a->nbs;
173949b5e25fSSatish Balay 
174049b5e25fSSatish Balay   c->imax       = (int*)PetscMalloc((mbs+1)*sizeof(int));CHKPTRQ(c->imax);
174149b5e25fSSatish Balay   c->ilen       = (int*)PetscMalloc((mbs+1)*sizeof(int));CHKPTRQ(c->ilen);
174249b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
174349b5e25fSSatish Balay     c->imax[i] = a->imax[i];
174449b5e25fSSatish Balay     c->ilen[i] = a->ilen[i];
174549b5e25fSSatish Balay   }
174649b5e25fSSatish Balay 
174749b5e25fSSatish Balay   /* allocate the matrix space */
174849b5e25fSSatish Balay   c->singlemalloc = PETSC_TRUE;
174949b5e25fSSatish Balay   len  = (mbs+1)*sizeof(int) + nz*(bs2*sizeof(MatScalar) + sizeof(int));
175049b5e25fSSatish Balay   c->a = (MatScalar*)PetscMalloc(len);CHKPTRQ(c->a);
175149b5e25fSSatish Balay   c->j = (int*)(c->a + nz*bs2);
175249b5e25fSSatish Balay   c->i = c->j + nz;
175349b5e25fSSatish Balay   ierr = PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int));CHKERRQ(ierr);
175449b5e25fSSatish Balay   if (mbs > 0) {
175549b5e25fSSatish Balay     ierr = PetscMemcpy(c->j,a->j,nz*sizeof(int));CHKERRQ(ierr);
175649b5e25fSSatish Balay     if (cpvalues == MAT_COPY_VALUES) {
175749b5e25fSSatish Balay       ierr = PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr);
175849b5e25fSSatish Balay     } else {
175949b5e25fSSatish Balay       ierr = PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr);
176049b5e25fSSatish Balay     }
176149b5e25fSSatish Balay   }
176249b5e25fSSatish Balay 
176349b5e25fSSatish Balay   PLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqSBAIJ));
176449b5e25fSSatish Balay   c->sorted      = a->sorted;
176549b5e25fSSatish Balay   c->roworiented = a->roworiented;
176649b5e25fSSatish Balay   c->nonew       = a->nonew;
176749b5e25fSSatish Balay 
176849b5e25fSSatish Balay   if (a->diag) {
176949b5e25fSSatish Balay     c->diag = (int*)PetscMalloc((mbs+1)*sizeof(int));CHKPTRQ(c->diag);
177049b5e25fSSatish Balay     PLogObjectMemory(C,(mbs+1)*sizeof(int));
177149b5e25fSSatish Balay     for (i=0; i<mbs; i++) {
177249b5e25fSSatish Balay       c->diag[i] = a->diag[i];
177349b5e25fSSatish Balay     }
177449b5e25fSSatish Balay   } else c->diag        = 0;
177549b5e25fSSatish Balay   c->s_nz               = a->s_nz;
177649b5e25fSSatish Balay   c->s_maxnz            = a->s_maxnz;
177749b5e25fSSatish Balay   c->solve_work         = 0;
177849b5e25fSSatish Balay   c->spptr              = 0;      /* Dangerous -I'm throwing away a->spptr */
177949b5e25fSSatish Balay   c->mult_work          = 0;
178049b5e25fSSatish Balay   *B = C;
178149b5e25fSSatish Balay   ierr = FListDuplicate(A->qlist,&C->qlist);CHKERRQ(ierr);
178249b5e25fSSatish Balay   PetscFunctionReturn(0);
178349b5e25fSSatish Balay }
178449b5e25fSSatish Balay 
178549b5e25fSSatish Balay #undef __FUNC__
178649b5e25fSSatish Balay #define __FUNC__ "MatLoad_SeqSBAIJ"
178749b5e25fSSatish Balay int MatLoad_SeqSBAIJ(Viewer viewer,MatType type,Mat *A)
178849b5e25fSSatish Balay {
178949b5e25fSSatish Balay   Mat_SeqSBAIJ  *a;
179049b5e25fSSatish Balay   Mat          B;
179149b5e25fSSatish Balay   int          i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1;
179249b5e25fSSatish Balay   int          *mask,mbs,*jj,j,rowcount,nzcount,k,*browlengths,maskcount;
179349b5e25fSSatish Balay   int          kmax,jcount,block,idx,point,nzcountb,extra_rows;
179449b5e25fSSatish Balay   int          *masked,nmask,tmp,bs2,ishift;
179549b5e25fSSatish Balay   Scalar       *aa;
179649b5e25fSSatish Balay   MPI_Comm     comm = ((PetscObject)viewer)->comm;
179749b5e25fSSatish Balay 
179849b5e25fSSatish Balay   PetscFunctionBegin;
179949b5e25fSSatish Balay   ierr = OptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,PETSC_NULL);CHKERRQ(ierr);
180049b5e25fSSatish Balay   bs2  = bs*bs;
180149b5e25fSSatish Balay 
180249b5e25fSSatish Balay   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
180349b5e25fSSatish Balay   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"view must have one processor");
180449b5e25fSSatish Balay   ierr = ViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
180549b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr);
180649b5e25fSSatish Balay   if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"not Mat object");
180749b5e25fSSatish Balay   M = header[1]; N = header[2]; nz = header[3];
180849b5e25fSSatish Balay 
180949b5e25fSSatish Balay   if (header[3] < 0) {
181049b5e25fSSatish Balay     SETERRQ(PETSC_ERR_FILE_UNEXPECTED,1,"Matrix stored in special format, cannot load as SeqSBAIJ");
181149b5e25fSSatish Balay   }
181249b5e25fSSatish Balay 
181349b5e25fSSatish Balay   if (M != N) SETERRQ(PETSC_ERR_SUP,0,"Can only do square matrices");
181449b5e25fSSatish Balay 
181549b5e25fSSatish Balay   /*
181649b5e25fSSatish Balay      This code adds extra rows to make sure the number of rows is
181749b5e25fSSatish Balay     divisible by the blocksize
181849b5e25fSSatish Balay   */
181949b5e25fSSatish Balay   mbs        = M/bs;
182049b5e25fSSatish Balay   extra_rows = bs - M + bs*(mbs);
182149b5e25fSSatish Balay   if (extra_rows == bs) extra_rows = 0;
182249b5e25fSSatish Balay   else                  mbs++;
182349b5e25fSSatish Balay   if (extra_rows) {
182449b5e25fSSatish Balay     PLogInfo(0,"MatLoad_SeqSBAIJ:Padding loaded matrix to match blocksize\n");
182549b5e25fSSatish Balay   }
182649b5e25fSSatish Balay 
182749b5e25fSSatish Balay   /* read in row lengths */
182849b5e25fSSatish Balay   rowlengths = (int*)PetscMalloc((M+extra_rows)*sizeof(int));CHKPTRQ(rowlengths);
182949b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
183049b5e25fSSatish Balay   for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1;
183149b5e25fSSatish Balay 
183249b5e25fSSatish Balay   /* read in column indices */
183349b5e25fSSatish Balay   jj = (int*)PetscMalloc((nz+extra_rows)*sizeof(int));CHKPTRQ(jj);
183449b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT);CHKERRQ(ierr);
183549b5e25fSSatish Balay   for (i=0; i<extra_rows; i++) jj[nz+i] = M+i;
183649b5e25fSSatish Balay 
183749b5e25fSSatish Balay   /* loop over row lengths determining block row lengths */
183849b5e25fSSatish Balay   browlengths = (int*)PetscMalloc(mbs*sizeof(int));CHKPTRQ(browlengths);
183949b5e25fSSatish Balay   ierr        = PetscMemzero(browlengths,mbs*sizeof(int));CHKERRQ(ierr);
184049b5e25fSSatish Balay   mask        = (int*)PetscMalloc(2*mbs*sizeof(int));CHKPTRQ(mask);
184149b5e25fSSatish Balay   ierr        = PetscMemzero(mask,mbs*sizeof(int));CHKERRQ(ierr);
184249b5e25fSSatish Balay   masked      = mask + mbs;
184349b5e25fSSatish Balay   rowcount    = 0; nzcount = 0;
184449b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
184549b5e25fSSatish Balay     nmask = 0;
184649b5e25fSSatish Balay     for (j=0; j<bs; j++) {
184749b5e25fSSatish Balay       kmax = rowlengths[rowcount];
184849b5e25fSSatish Balay       for (k=0; k<kmax; k++) {
184949b5e25fSSatish Balay         tmp = jj[nzcount++]/bs;
185049b5e25fSSatish Balay         if (!mask[tmp]) {masked[nmask++] = tmp; mask[tmp] = 1;}
185149b5e25fSSatish Balay       }
185249b5e25fSSatish Balay       rowcount++;
185349b5e25fSSatish Balay     }
185449b5e25fSSatish Balay     browlengths[i] += nmask;
185549b5e25fSSatish Balay     /* zero out the mask elements we set */
185649b5e25fSSatish Balay     for (j=0; j<nmask; j++) mask[masked[j]] = 0;
185749b5e25fSSatish Balay   }
185849b5e25fSSatish Balay 
185949b5e25fSSatish Balay   /* create our matrix */
186049b5e25fSSatish Balay   ierr = MatCreateSeqSBAIJ(comm,bs,M+extra_rows,N+extra_rows,0,browlengths,A);CHKERRQ(ierr);
186149b5e25fSSatish Balay   B = *A;
186249b5e25fSSatish Balay   a = (Mat_SeqSBAIJ*)B->data;
186349b5e25fSSatish Balay 
186449b5e25fSSatish Balay   /* set matrix "i" values */
186549b5e25fSSatish Balay   a->i[0] = 0;
186649b5e25fSSatish Balay   for (i=1; i<= mbs; i++) {
186749b5e25fSSatish Balay     a->i[i]      = a->i[i-1] + browlengths[i-1];
186849b5e25fSSatish Balay     a->ilen[i-1] = browlengths[i-1];
186949b5e25fSSatish Balay   }
187049b5e25fSSatish Balay   a->s_nz         = 0;
187149b5e25fSSatish Balay   for (i=0; i<mbs; i++) a->s_nz += browlengths[i];
187249b5e25fSSatish Balay 
187349b5e25fSSatish Balay   /* read in nonzero values */
187449b5e25fSSatish Balay   aa = (Scalar*)PetscMalloc((nz+extra_rows)*sizeof(Scalar));CHKPTRQ(aa);
187549b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);CHKERRQ(ierr);
187649b5e25fSSatish Balay   for (i=0; i<extra_rows; i++) aa[nz+i] = 1.0;
187749b5e25fSSatish Balay 
187849b5e25fSSatish Balay   /* set "a" and "j" values into matrix */
187949b5e25fSSatish Balay   nzcount = 0; jcount = 0;
188049b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
188149b5e25fSSatish Balay     nzcountb = nzcount;
188249b5e25fSSatish Balay     nmask    = 0;
188349b5e25fSSatish Balay     for (j=0; j<bs; j++) {
188449b5e25fSSatish Balay       kmax = rowlengths[i*bs+j];
188549b5e25fSSatish Balay       for (k=0; k<kmax; k++) {
188649b5e25fSSatish Balay         tmp = jj[nzcount++]/bs;
188749b5e25fSSatish Balay 	if (!mask[tmp]) { masked[nmask++] = tmp; mask[tmp] = 1;}
188849b5e25fSSatish Balay       }
188949b5e25fSSatish Balay       rowcount++;
189049b5e25fSSatish Balay     }
189149b5e25fSSatish Balay     /* sort the masked values */
189249b5e25fSSatish Balay     ierr = PetscSortInt(nmask,masked);CHKERRQ(ierr);
189349b5e25fSSatish Balay 
189449b5e25fSSatish Balay     /* set "j" values into matrix */
189549b5e25fSSatish Balay     maskcount = 1;
189649b5e25fSSatish Balay     for (j=0; j<nmask; j++) {
189749b5e25fSSatish Balay       a->j[jcount++]  = masked[j];
189849b5e25fSSatish Balay       mask[masked[j]] = maskcount++;
189949b5e25fSSatish Balay     }
190049b5e25fSSatish Balay     /* set "a" values into matrix */
190149b5e25fSSatish Balay     ishift = bs2*a->i[i];
190249b5e25fSSatish Balay     for (j=0; j<bs; j++) {
190349b5e25fSSatish Balay       kmax = rowlengths[i*bs+j];
190449b5e25fSSatish Balay       for (k=0; k<kmax; k++) {
190549b5e25fSSatish Balay         tmp       = jj[nzcountb]/bs ;
190649b5e25fSSatish Balay         block     = mask[tmp] - 1;
190749b5e25fSSatish Balay         point     = jj[nzcountb] - bs*tmp;
190849b5e25fSSatish Balay         idx       = ishift + bs2*block + j + bs*point;
190949b5e25fSSatish Balay         a->a[idx] = aa[nzcountb++];
191049b5e25fSSatish Balay       }
191149b5e25fSSatish Balay     }
191249b5e25fSSatish Balay     /* zero out the mask elements we set */
191349b5e25fSSatish Balay     for (j=0; j<nmask; j++) mask[masked[j]] = 0;
191449b5e25fSSatish Balay   }
191549b5e25fSSatish Balay   if (jcount != a->s_nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"Bad binary matrix");
191649b5e25fSSatish Balay 
191749b5e25fSSatish Balay   ierr = PetscFree(rowlengths);CHKERRQ(ierr);
191849b5e25fSSatish Balay   ierr = PetscFree(browlengths);CHKERRQ(ierr);
191949b5e25fSSatish Balay   ierr = PetscFree(aa);CHKERRQ(ierr);
192049b5e25fSSatish Balay   ierr = PetscFree(jj);CHKERRQ(ierr);
192149b5e25fSSatish Balay   ierr = PetscFree(mask);CHKERRQ(ierr);
192249b5e25fSSatish Balay 
192349b5e25fSSatish Balay   B->assembled = PETSC_TRUE;
192449b5e25fSSatish Balay 
192549b5e25fSSatish Balay   ierr = MatView_Private(B);CHKERRQ(ierr);
192649b5e25fSSatish Balay   PetscFunctionReturn(0);
192749b5e25fSSatish Balay }
192849b5e25fSSatish Balay 
192949b5e25fSSatish Balay 
193049b5e25fSSatish Balay 
1931