xref: /petsc/src/mat/impls/sbaij/seq/sbaij.c (revision 1ce8b1feb7feb99415e7b3f208ae5246e3432011)
1*1ce8b1feSHong Zhang /*$Id: sbaij.c,v 1.25 2000/09/14 20:42:31 hzhang Exp hzhang $*/
249b5e25fSSatish Balay 
349b5e25fSSatish Balay /*
449b5e25fSSatish Balay     Defines the basic matrix operations for the BAIJ (compressed row)
549b5e25fSSatish Balay   matrix storage format.
649b5e25fSSatish Balay */
7ab9f2c04SSatish Balay #include "petscsys.h"                            /*I "petscmat.h" I*/
849b5e25fSSatish Balay #include "src/mat/impls/baij/seq/baij.h"
949b5e25fSSatish Balay #include "src/vec/vecimpl.h"
1049b5e25fSSatish Balay #include "src/inline/spops.h"
11aec82049SSatish Balay #include "src/mat/impls/sbaij/seq/sbaij.h"
1249b5e25fSSatish Balay 
1349b5e25fSSatish Balay #define CHUNKSIZE  10
1449b5e25fSSatish Balay 
1549b5e25fSSatish Balay /*
1649b5e25fSSatish Balay      Checks for missing diagonals
1749b5e25fSSatish Balay */
1849b5e25fSSatish Balay #undef __FUNC__
1949b5e25fSSatish Balay #define __FUNC__ "MatMissingDiagonal_SeqSBAIJ"
2049b5e25fSSatish Balay int MatMissingDiagonal_SeqSBAIJ(Mat A)
2149b5e25fSSatish Balay {
22045c9aa0SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
2349b5e25fSSatish Balay   int         *diag,*jj = a->j,i,ierr;
2449b5e25fSSatish Balay 
2549b5e25fSSatish Balay   PetscFunctionBegin;
26045c9aa0SHong Zhang   ierr = MatMarkDiagonal_SeqSBAIJ(A);CHKERRQ(ierr);
2749b5e25fSSatish Balay   diag = a->diag;
2849b5e25fSSatish Balay   for (i=0; i<a->mbs; i++) {
2949b5e25fSSatish Balay     if (jj[diag[i]] != i) {
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 {
40045c9aa0SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
4149b5e25fSSatish Balay   int          i,j,*diag,m = a->mbs;
4249b5e25fSSatish Balay 
4349b5e25fSSatish Balay   PetscFunctionBegin;
4449b5e25fSSatish Balay   if (a->diag) PetscFunctionReturn(0);
4549b5e25fSSatish Balay 
4649b5e25fSSatish Balay   diag = (int*)PetscMalloc((m+1)*sizeof(int));CHKPTRQ(diag);
4749b5e25fSSatish Balay   PLogObjectMemory(A,(m+1)*sizeof(int));
4849b5e25fSSatish Balay   for (i=0; i<m; i++) {
4949b5e25fSSatish Balay     diag[i] = a->i[i+1];
5049b5e25fSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
5149b5e25fSSatish Balay       if (a->j[j] == i) {
5249b5e25fSSatish Balay         diag[i] = j;
5349b5e25fSSatish Balay         break;
5449b5e25fSSatish Balay       }
5549b5e25fSSatish Balay     }
5649b5e25fSSatish Balay   }
5749b5e25fSSatish Balay   a->diag = diag;
5849b5e25fSSatish Balay   PetscFunctionReturn(0);
5949b5e25fSSatish Balay }
6049b5e25fSSatish Balay 
6149b5e25fSSatish Balay 
6249b5e25fSSatish Balay extern int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int,int**,int**);
6349b5e25fSSatish Balay 
6449b5e25fSSatish Balay #undef __FUNC__
6549b5e25fSSatish Balay #define __FUNC__ "MatGetRowIJ_SeqSBAIJ"
6649b5e25fSSatish Balay static int MatGetRowIJ_SeqSBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja,PetscTruth *done)
6749b5e25fSSatish Balay {
6849b5e25fSSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
6949b5e25fSSatish Balay 
7049b5e25fSSatish Balay   PetscFunctionBegin;
71045c9aa0SHong Zhang   if (ia) {
72045c9aa0SHong Zhang     SETERRQ(1,1,"Function not yet written for SBAIJ format, only supports natural ordering");
7349b5e25fSSatish Balay   }
74045c9aa0SHong Zhang   *nn = a->mbs;
7549b5e25fSSatish Balay   PetscFunctionReturn(0);
7649b5e25fSSatish Balay }
7749b5e25fSSatish Balay 
7849b5e25fSSatish Balay #undef __FUNC__
7949b5e25fSSatish Balay #define __FUNC__ "MatRestoreRowIJ_SeqSBAIJ"
80045c9aa0SHong Zhang static int MatRestoreRowIJ_SeqSBAIJ(Mat A,int oshift,PetscTruth symmetric,int *nn,int **ia,int **ja,PetscTruth *done)
8149b5e25fSSatish Balay {
8249b5e25fSSatish Balay   PetscFunctionBegin;
8349b5e25fSSatish Balay   if (!ia) PetscFunctionReturn(0);
84045c9aa0SHong Zhang   SETERRQ(1,1,"Function not yet written for SBAIJ format, only supports natural ordering");
8516cdd363SHong Zhang   /* PetscFunctionReturn(0); */
8649b5e25fSSatish Balay }
8749b5e25fSSatish Balay 
8849b5e25fSSatish Balay #undef __FUNC__
8949b5e25fSSatish Balay #define __FUNC__ "MatGetBlockSize_SeqSBAIJ"
9049b5e25fSSatish Balay int MatGetBlockSize_SeqSBAIJ(Mat mat,int *bs)
9149b5e25fSSatish Balay {
9249b5e25fSSatish Balay   Mat_SeqSBAIJ *sbaij = (Mat_SeqSBAIJ*)mat->data;
9349b5e25fSSatish Balay 
9449b5e25fSSatish Balay   PetscFunctionBegin;
9549b5e25fSSatish Balay   *bs = sbaij->bs;
9649b5e25fSSatish Balay   PetscFunctionReturn(0);
9749b5e25fSSatish Balay }
9849b5e25fSSatish Balay 
9949b5e25fSSatish Balay #undef __FUNC__
10049b5e25fSSatish Balay #define __FUNC__ "MatDestroy_SeqSBAIJ"
10149b5e25fSSatish Balay int MatDestroy_SeqSBAIJ(Mat A)
10249b5e25fSSatish Balay {
10349b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
10449b5e25fSSatish Balay   int         ierr;
10549b5e25fSSatish Balay 
10649b5e25fSSatish Balay   PetscFunctionBegin;
10749b5e25fSSatish Balay   if (A->mapping) {
10849b5e25fSSatish Balay     ierr = ISLocalToGlobalMappingDestroy(A->mapping);CHKERRQ(ierr);
10949b5e25fSSatish Balay   }
11049b5e25fSSatish Balay   if (A->bmapping) {
11149b5e25fSSatish Balay     ierr = ISLocalToGlobalMappingDestroy(A->bmapping);CHKERRQ(ierr);
11249b5e25fSSatish Balay   }
11349b5e25fSSatish Balay   if (A->rmap) {
11449b5e25fSSatish Balay     ierr = MapDestroy(A->rmap);CHKERRQ(ierr);
11549b5e25fSSatish Balay   }
11649b5e25fSSatish Balay   if (A->cmap) {
11749b5e25fSSatish Balay     ierr = MapDestroy(A->cmap);CHKERRQ(ierr);
11849b5e25fSSatish Balay   }
11949b5e25fSSatish Balay #if defined(PETSC_USE_LOG)
12049b5e25fSSatish Balay   PLogObjectState((PetscObject)A,"Rows=%d, s_NZ=%d",a->m,a->s_nz);
12149b5e25fSSatish Balay #endif
12249b5e25fSSatish Balay   ierr = PetscFree(a->a);CHKERRQ(ierr);
12349b5e25fSSatish Balay   if (!a->singlemalloc) {
12449b5e25fSSatish Balay     ierr = PetscFree(a->i);CHKERRQ(ierr);
12563c8bf9fSHong Zhang     ierr = PetscFree(a->j);CHKERRQ(ierr);
12649b5e25fSSatish Balay   }
12749b5e25fSSatish Balay   if (a->row) {
12849b5e25fSSatish Balay     ierr = ISDestroy(a->row);CHKERRQ(ierr);
12949b5e25fSSatish Balay   }
13049b5e25fSSatish Balay   if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);}
13149b5e25fSSatish Balay   if (a->ilen) {ierr = PetscFree(a->ilen);CHKERRQ(ierr);}
13249b5e25fSSatish Balay   if (a->imax) {ierr = PetscFree(a->imax);CHKERRQ(ierr);}
13349b5e25fSSatish Balay   if (a->solve_work) {ierr = PetscFree(a->solve_work);CHKERRQ(ierr);}
13449b5e25fSSatish Balay   if (a->mult_work) {ierr = PetscFree(a->mult_work);CHKERRQ(ierr);}
13549b5e25fSSatish Balay   if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);}
13649b5e25fSSatish Balay   if (a->saved_values) {ierr = PetscFree(a->saved_values);CHKERRQ(ierr);}
13749b5e25fSSatish Balay   ierr = PetscFree(a);CHKERRQ(ierr);
13849b5e25fSSatish Balay   PLogObjectDestroy(A);
13949b5e25fSSatish Balay   PetscHeaderDestroy(A);
14049b5e25fSSatish Balay   PetscFunctionReturn(0);
14149b5e25fSSatish Balay }
14249b5e25fSSatish Balay 
14349b5e25fSSatish Balay #undef __FUNC__
14449b5e25fSSatish Balay #define __FUNC__ "MatSetOption_SeqSBAIJ"
14549b5e25fSSatish Balay int MatSetOption_SeqSBAIJ(Mat A,MatOption op)
14649b5e25fSSatish Balay {
147045c9aa0SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
14849b5e25fSSatish Balay 
14949b5e25fSSatish Balay   PetscFunctionBegin;
15049b5e25fSSatish Balay   if      (op == MAT_ROW_ORIENTED)                 a->roworiented    = PETSC_TRUE;
15149b5e25fSSatish Balay   else if (op == MAT_COLUMN_ORIENTED)              a->roworiented    = PETSC_FALSE;
15249b5e25fSSatish Balay   else if (op == MAT_COLUMNS_SORTED)               a->sorted         = PETSC_TRUE;
15349b5e25fSSatish Balay   else if (op == MAT_COLUMNS_UNSORTED)             a->sorted         = PETSC_FALSE;
15449b5e25fSSatish Balay   else if (op == MAT_KEEP_ZEROED_ROWS)             a->keepzeroedrows = PETSC_TRUE;
15549b5e25fSSatish Balay   else if (op == MAT_NO_NEW_NONZERO_LOCATIONS)     a->nonew          = 1;
15649b5e25fSSatish Balay   else if (op == MAT_NEW_NONZERO_LOCATION_ERR)     a->nonew          = -1;
15749b5e25fSSatish Balay   else if (op == MAT_NEW_NONZERO_ALLOCATION_ERR)   a->nonew          = -2;
15849b5e25fSSatish Balay   else if (op == MAT_YES_NEW_NONZERO_LOCATIONS)    a->nonew          = 0;
15949b5e25fSSatish Balay   else if (op == MAT_ROWS_SORTED ||
16049b5e25fSSatish Balay            op == MAT_ROWS_UNSORTED ||
16149b5e25fSSatish Balay            op == MAT_SYMMETRIC ||
16249b5e25fSSatish Balay            op == MAT_STRUCTURALLY_SYMMETRIC ||
16349b5e25fSSatish Balay            op == MAT_YES_NEW_DIAGONALS ||
16449b5e25fSSatish Balay            op == MAT_IGNORE_OFF_PROC_ENTRIES ||
16549b5e25fSSatish Balay            op == MAT_USE_HASH_TABLE) {
166045c9aa0SHong Zhang     PLogInfo(A,"MatSetOption_SeqSBAIJ:Option ignored\n");
16749b5e25fSSatish Balay   } else if (op == MAT_NO_NEW_DIAGONALS) {
16849b5e25fSSatish Balay     SETERRQ(PETSC_ERR_SUP,0,"MAT_NO_NEW_DIAGONALS");
16949b5e25fSSatish Balay   } else {
17049b5e25fSSatish Balay     SETERRQ(PETSC_ERR_SUP,0,"unknown option");
17149b5e25fSSatish Balay   }
17249b5e25fSSatish Balay   PetscFunctionReturn(0);
17349b5e25fSSatish Balay }
17449b5e25fSSatish Balay 
17549b5e25fSSatish Balay 
17649b5e25fSSatish Balay #undef __FUNC__
17749b5e25fSSatish Balay #define __FUNC__ "MatGetSize_SeqSBAIJ"
17849b5e25fSSatish Balay int MatGetSize_SeqSBAIJ(Mat A,int *m,int *n)
17949b5e25fSSatish Balay {
18049b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
18149b5e25fSSatish Balay 
18249b5e25fSSatish Balay   PetscFunctionBegin;
18349b5e25fSSatish Balay   if (m) *m = a->m;
184045c9aa0SHong Zhang   if (n) *n = a->n;
18549b5e25fSSatish Balay   PetscFunctionReturn(0);
18649b5e25fSSatish Balay }
18749b5e25fSSatish Balay 
18849b5e25fSSatish Balay #undef __FUNC__
18949b5e25fSSatish Balay #define __FUNC__ "MatGetOwnershipRange_SeqSBAIJ"
19049b5e25fSSatish Balay int MatGetOwnershipRange_SeqSBAIJ(Mat A,int *m,int *n)
19149b5e25fSSatish Balay {
19249b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
19349b5e25fSSatish Balay 
19449b5e25fSSatish Balay   PetscFunctionBegin;
19549b5e25fSSatish Balay   *m = 0; *n = a->m;
19649b5e25fSSatish Balay   PetscFunctionReturn(0);
19749b5e25fSSatish Balay }
19849b5e25fSSatish Balay 
19949b5e25fSSatish Balay #undef __FUNC__
20049b5e25fSSatish Balay #define __FUNC__ "MatGetRow_SeqSBAIJ"
20149b5e25fSSatish Balay int MatGetRow_SeqSBAIJ(Mat A,int row,int *ncols,int **cols,Scalar **v)
20249b5e25fSSatish Balay {
20349b5e25fSSatish Balay   Mat_SeqSBAIJ  *a = (Mat_SeqSBAIJ*)A->data;
20449b5e25fSSatish Balay   int          itmp,i,j,k,M,*ai,*aj,bs,bn,bp,*cols_i,bs2;
20549b5e25fSSatish Balay   MatScalar    *aa,*aa_i;
20649b5e25fSSatish Balay   Scalar       *v_i;
20749b5e25fSSatish Balay 
20849b5e25fSSatish Balay   PetscFunctionBegin;
20949b5e25fSSatish Balay   bs  = a->bs;
21049b5e25fSSatish Balay   ai  = a->i;
21149b5e25fSSatish Balay   aj  = a->j;
21249b5e25fSSatish Balay   aa  = a->a;
21349b5e25fSSatish Balay   bs2 = a->bs2;
21449b5e25fSSatish Balay 
21549b5e25fSSatish Balay   if (row < 0 || row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row out of range");
21649b5e25fSSatish Balay 
21749b5e25fSSatish Balay   bn  = row/bs;   /* Block number */
21849b5e25fSSatish Balay   bp  = row % bs; /* Block position */
21949b5e25fSSatish Balay   M   = ai[bn+1] - ai[bn];
22049b5e25fSSatish Balay   *ncols = bs*M;
22149b5e25fSSatish Balay 
22249b5e25fSSatish Balay   if (v) {
22349b5e25fSSatish Balay     *v = 0;
22449b5e25fSSatish Balay     if (*ncols) {
22549b5e25fSSatish Balay       *v = (Scalar*)PetscMalloc((*ncols+row)*sizeof(Scalar));CHKPTRQ(*v);
22649b5e25fSSatish Balay       for (i=0; i<M; i++) { /* for each block in the block row */
22749b5e25fSSatish Balay         v_i  = *v + i*bs;
22849b5e25fSSatish Balay         aa_i = aa + bs2*(ai[bn] + i);
22949b5e25fSSatish Balay         for (j=bp,k=0; j<bs2; j+=bs,k++) {v_i[k] = aa_i[j];}
23049b5e25fSSatish Balay       }
23149b5e25fSSatish Balay     }
23249b5e25fSSatish Balay   }
23349b5e25fSSatish Balay 
23449b5e25fSSatish Balay   if (cols) {
23549b5e25fSSatish Balay     *cols = 0;
23649b5e25fSSatish Balay     if (*ncols) {
23749b5e25fSSatish Balay       *cols = (int*)PetscMalloc((*ncols+row)*sizeof(int));CHKPTRQ(*cols);
23849b5e25fSSatish Balay       for (i=0; i<M; i++) { /* for each block in the block row */
23949b5e25fSSatish Balay         cols_i = *cols + i*bs;
24049b5e25fSSatish Balay         itmp  = bs*aj[ai[bn] + i];
24149b5e25fSSatish Balay         for (j=0; j<bs; j++) {cols_i[j] = itmp++;}
24249b5e25fSSatish Balay       }
24349b5e25fSSatish Balay     }
24449b5e25fSSatish Balay   }
24549b5e25fSSatish Balay 
24649b5e25fSSatish Balay   /*search column A(0:row-1,row) (=A(row,0:row-1)). Could be expensive! */
2475ddb2528SHong Zhang   /* this segment is currently removed, so only entries in the upper triangle are obtained */
2485ddb2528SHong Zhang #ifdef column_search
24949b5e25fSSatish Balay   v_i    = *v    + M*bs;
25049b5e25fSSatish Balay   cols_i = *cols + M*bs;
25149b5e25fSSatish Balay   for (i=0; i<bn; i++){ /* for each block row */
25249b5e25fSSatish Balay     M = ai[i+1] - ai[i];
25349b5e25fSSatish Balay     for (j=0; j<M; j++){
25449b5e25fSSatish Balay       itmp = aj[ai[i] + j];    /* block column value */
25549b5e25fSSatish Balay       if (itmp == bn){
25649b5e25fSSatish Balay         aa_i   = aa    + bs2*(ai[i] + j) + bs*bp;
25749b5e25fSSatish Balay         for (k=0; k<bs; k++) {
25849b5e25fSSatish Balay           *cols_i++ = i*bs+k;
25949b5e25fSSatish Balay           *v_i++    = aa_i[k];
26049b5e25fSSatish Balay         }
26149b5e25fSSatish Balay         *ncols += bs;
26249b5e25fSSatish Balay         break;
26349b5e25fSSatish Balay       }
26449b5e25fSSatish Balay     }
26549b5e25fSSatish Balay   }
2665ddb2528SHong Zhang #endif
26749b5e25fSSatish Balay 
26849b5e25fSSatish Balay   PetscFunctionReturn(0);
26949b5e25fSSatish Balay }
27049b5e25fSSatish Balay 
27149b5e25fSSatish Balay #undef __FUNC__
27249b5e25fSSatish Balay #define __FUNC__ "MatRestoreRow_SeqSBAIJ"
27349b5e25fSSatish Balay int MatRestoreRow_SeqSBAIJ(Mat A,int row,int *nz,int **idx,Scalar **v)
27449b5e25fSSatish Balay {
27549b5e25fSSatish Balay   int ierr;
27649b5e25fSSatish Balay 
27749b5e25fSSatish Balay   PetscFunctionBegin;
27849b5e25fSSatish Balay   if (idx) {if (*idx) {ierr = PetscFree(*idx);CHKERRQ(ierr);}}
27949b5e25fSSatish Balay   if (v)   {if (*v)   {ierr = PetscFree(*v);CHKERRQ(ierr);}}
28049b5e25fSSatish Balay   PetscFunctionReturn(0);
28149b5e25fSSatish Balay }
28249b5e25fSSatish Balay 
28349b5e25fSSatish Balay #undef __FUNC__
28449b5e25fSSatish Balay #define __FUNC__ "MatTranspose_SeqSBAIJ"
28549b5e25fSSatish Balay int MatTranspose_SeqSBAIJ(Mat A,Mat *B)
28649b5e25fSSatish Balay {
28749b5e25fSSatish Balay   PetscFunctionBegin;
288045c9aa0SHong Zhang   SETERRQ(1,1,"Matrix is symmetric. MatTranspose() should not be called");
28916cdd363SHong Zhang   /* PetscFunctionReturn(0); */
29049b5e25fSSatish Balay }
29149b5e25fSSatish Balay 
29249b5e25fSSatish Balay #undef __FUNC__
29349b5e25fSSatish Balay #define __FUNC__ "MatView_SeqSBAIJ_Binary"
29449b5e25fSSatish Balay static int MatView_SeqSBAIJ_Binary(Mat A,Viewer viewer)
29549b5e25fSSatish Balay {
29649b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
29749b5e25fSSatish Balay   int         i,fd,*col_lens,ierr,bs = a->bs,count,*jj,j,k,l,bs2=a->bs2;
29849b5e25fSSatish Balay   Scalar      *aa;
29949b5e25fSSatish Balay   FILE        *file;
30049b5e25fSSatish Balay 
30149b5e25fSSatish Balay   PetscFunctionBegin;
30249b5e25fSSatish Balay   ierr = ViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
30349b5e25fSSatish Balay   col_lens = (int*)PetscMalloc((4+a->m)*sizeof(int));CHKPTRQ(col_lens);
30449b5e25fSSatish Balay   col_lens[0] = MAT_COOKIE;
30549b5e25fSSatish Balay 
30649b5e25fSSatish Balay   col_lens[1] = a->m;
30749b5e25fSSatish Balay   col_lens[2] = a->m;
30849b5e25fSSatish Balay   col_lens[3] = a->s_nz*bs2;
30949b5e25fSSatish Balay 
31049b5e25fSSatish Balay   /* store lengths of each row and write (including header) to file */
31149b5e25fSSatish Balay   count = 0;
31249b5e25fSSatish Balay   for (i=0; i<a->mbs; i++) {
31349b5e25fSSatish Balay     for (j=0; j<bs; j++) {
31449b5e25fSSatish Balay       col_lens[4+count++] = bs*(a->i[i+1] - a->i[i]);
31549b5e25fSSatish Balay     }
31649b5e25fSSatish Balay   }
31749b5e25fSSatish Balay   ierr = PetscBinaryWrite(fd,col_lens,4+a->m,PETSC_INT,1);CHKERRQ(ierr);
31849b5e25fSSatish Balay   ierr = PetscFree(col_lens);CHKERRQ(ierr);
31949b5e25fSSatish Balay 
32049b5e25fSSatish Balay   /* store column indices (zero start index) */
32149b5e25fSSatish Balay   jj = (int*)PetscMalloc((a->s_nz+1)*bs2*sizeof(int));CHKPTRQ(jj);
32249b5e25fSSatish Balay   count = 0;
32349b5e25fSSatish Balay   for (i=0; i<a->mbs; i++) {
32449b5e25fSSatish Balay     for (j=0; j<bs; j++) {
32549b5e25fSSatish Balay       for (k=a->i[i]; k<a->i[i+1]; k++) {
32649b5e25fSSatish Balay         for (l=0; l<bs; l++) {
32749b5e25fSSatish Balay           jj[count++] = bs*a->j[k] + l;
32849b5e25fSSatish Balay         }
32949b5e25fSSatish Balay       }
33049b5e25fSSatish Balay     }
33149b5e25fSSatish Balay   }
33249b5e25fSSatish Balay   ierr = PetscBinaryWrite(fd,jj,bs2*a->s_nz,PETSC_INT,0);CHKERRQ(ierr);
33349b5e25fSSatish Balay   ierr = PetscFree(jj);CHKERRQ(ierr);
33449b5e25fSSatish Balay 
33549b5e25fSSatish Balay   /* store nonzero values */
33649b5e25fSSatish Balay   aa = (Scalar*)PetscMalloc((a->s_nz+1)*bs2*sizeof(Scalar));CHKPTRQ(aa);
33749b5e25fSSatish Balay   count = 0;
33849b5e25fSSatish Balay   for (i=0; i<a->mbs; i++) {
33949b5e25fSSatish Balay     for (j=0; j<bs; j++) {
34049b5e25fSSatish Balay       for (k=a->i[i]; k<a->i[i+1]; k++) {
34149b5e25fSSatish Balay         for (l=0; l<bs; l++) {
34249b5e25fSSatish Balay           aa[count++] = a->a[bs2*k + l*bs + j];
34349b5e25fSSatish Balay         }
34449b5e25fSSatish Balay       }
34549b5e25fSSatish Balay     }
34649b5e25fSSatish Balay   }
34749b5e25fSSatish Balay   ierr = PetscBinaryWrite(fd,aa,bs2*a->s_nz,PETSC_SCALAR,0);CHKERRQ(ierr);
34849b5e25fSSatish Balay   ierr = PetscFree(aa);CHKERRQ(ierr);
34949b5e25fSSatish Balay 
35049b5e25fSSatish Balay   ierr = ViewerBinaryGetInfoPointer(viewer,&file);CHKERRQ(ierr);
35149b5e25fSSatish Balay   if (file) {
35249b5e25fSSatish Balay     fprintf(file,"-matload_block_size %d\n",a->bs);
35349b5e25fSSatish Balay   }
35449b5e25fSSatish Balay   PetscFunctionReturn(0);
35549b5e25fSSatish Balay }
35649b5e25fSSatish Balay 
35749b5e25fSSatish Balay #undef __FUNC__
35849b5e25fSSatish Balay #define __FUNC__ "MatView_SeqSBAIJ_ASCII"
35949b5e25fSSatish Balay static int MatView_SeqSBAIJ_ASCII(Mat A,Viewer viewer)
36049b5e25fSSatish Balay {
36149b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
36249b5e25fSSatish Balay   int         ierr,i,j,format,bs = a->bs,k,l,bs2=a->bs2;
36349b5e25fSSatish Balay   char        *outputname;
36449b5e25fSSatish Balay 
36549b5e25fSSatish Balay   PetscFunctionBegin;
36649b5e25fSSatish Balay   ierr = ViewerGetOutputname(viewer,&outputname);CHKERRQ(ierr);
36749b5e25fSSatish Balay   ierr = ViewerGetFormat(viewer,&format);CHKERRQ(ierr);
36849b5e25fSSatish Balay   if (format == VIEWER_FORMAT_ASCII_INFO || format == VIEWER_FORMAT_ASCII_INFO_LONG) {
36949b5e25fSSatish Balay     ierr = ViewerASCIIPrintf(viewer,"  block size is %d\n",bs);CHKERRQ(ierr);
37049b5e25fSSatish Balay   } else if (format == VIEWER_FORMAT_ASCII_MATLAB) {
37149b5e25fSSatish Balay     SETERRQ(PETSC_ERR_SUP,0,"Matlab format not supported");
37249b5e25fSSatish Balay   } else if (format == VIEWER_FORMAT_ASCII_COMMON) {
37349b5e25fSSatish Balay     ierr = ViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
37449b5e25fSSatish Balay     for (i=0; i<a->mbs; i++) {
37549b5e25fSSatish Balay       for (j=0; j<bs; j++) {
37649b5e25fSSatish Balay         ierr = ViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr);
37749b5e25fSSatish Balay         for (k=a->i[i]; k<a->i[i+1]; k++) {
37849b5e25fSSatish Balay           for (l=0; l<bs; l++) {
37949b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX)
38049b5e25fSSatish Balay             if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
38149b5e25fSSatish Balay               ierr = ViewerASCIIPrintf(viewer," %d %g + %g i",bs*a->j[k]+l,
38249b5e25fSSatish Balay                       PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
38349b5e25fSSatish Balay             } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0 && PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
38449b5e25fSSatish Balay               ierr = ViewerASCIIPrintf(viewer," %d %g - %g i",bs*a->j[k]+l,
38549b5e25fSSatish Balay                       PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
38649b5e25fSSatish Balay             } else if (PetscRealPart(a->a[bs2*k + l*bs + j]) != 0.0) {
38749b5e25fSSatish Balay               ierr = ViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
38849b5e25fSSatish Balay             }
38949b5e25fSSatish Balay #else
39049b5e25fSSatish Balay             if (a->a[bs2*k + l*bs + j] != 0.0) {
39149b5e25fSSatish Balay               ierr = ViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr);
39249b5e25fSSatish Balay             }
39349b5e25fSSatish Balay #endif
39449b5e25fSSatish Balay           }
39549b5e25fSSatish Balay         }
39649b5e25fSSatish Balay         ierr = ViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
39749b5e25fSSatish Balay       }
39849b5e25fSSatish Balay     }
39949b5e25fSSatish Balay     ierr = ViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
40049b5e25fSSatish Balay   } else {
40149b5e25fSSatish Balay     ierr = ViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr);
40249b5e25fSSatish Balay     for (i=0; i<a->mbs; i++) {
40349b5e25fSSatish Balay       for (j=0; j<bs; j++) {
40449b5e25fSSatish Balay         ierr = ViewerASCIIPrintf(viewer,"row %d:",i*bs+j);CHKERRQ(ierr);
40549b5e25fSSatish Balay         for (k=a->i[i]; k<a->i[i+1]; k++) {
40649b5e25fSSatish Balay           for (l=0; l<bs; l++) {
40749b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX)
40849b5e25fSSatish Balay             if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) > 0.0) {
40949b5e25fSSatish Balay               ierr = ViewerASCIIPrintf(viewer," %d %g + %g i",bs*a->j[k]+l,
41049b5e25fSSatish Balay                 PetscRealPart(a->a[bs2*k + l*bs + j]),PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
41149b5e25fSSatish Balay             } else if (PetscImaginaryPart(a->a[bs2*k + l*bs + j]) < 0.0) {
41249b5e25fSSatish Balay               ierr = ViewerASCIIPrintf(viewer," %d %g - %g i",bs*a->j[k]+l,
41349b5e25fSSatish Balay                 PetscRealPart(a->a[bs2*k + l*bs + j]),-PetscImaginaryPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
41449b5e25fSSatish Balay             } else {
41549b5e25fSSatish Balay               ierr = ViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,PetscRealPart(a->a[bs2*k + l*bs + j]));CHKERRQ(ierr);
41649b5e25fSSatish Balay             }
41749b5e25fSSatish Balay #else
41849b5e25fSSatish Balay             ierr = ViewerASCIIPrintf(viewer," %d %g ",bs*a->j[k]+l,a->a[bs2*k + l*bs + j]);CHKERRQ(ierr);
41949b5e25fSSatish Balay #endif
42049b5e25fSSatish Balay           }
42149b5e25fSSatish Balay         }
42249b5e25fSSatish Balay         ierr = ViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
42349b5e25fSSatish Balay       }
42449b5e25fSSatish Balay     }
42549b5e25fSSatish Balay     ierr = ViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr);
42649b5e25fSSatish Balay   }
42749b5e25fSSatish Balay   ierr = ViewerFlush(viewer);CHKERRQ(ierr);
42849b5e25fSSatish Balay   PetscFunctionReturn(0);
42949b5e25fSSatish Balay }
43049b5e25fSSatish Balay 
43149b5e25fSSatish Balay #undef __FUNC__
43249b5e25fSSatish Balay #define __FUNC__ "MatView_SeqSBAIJ_Draw_Zoom"
43349b5e25fSSatish Balay static int MatView_SeqSBAIJ_Draw_Zoom(Draw draw,void *Aa)
43449b5e25fSSatish Balay {
43549b5e25fSSatish Balay   Mat          A = (Mat) Aa;
43649b5e25fSSatish Balay   Mat_SeqSBAIJ  *a=(Mat_SeqSBAIJ*)A->data;
43749b5e25fSSatish Balay   int          row,ierr,i,j,k,l,mbs=a->mbs,color,bs=a->bs,bs2=a->bs2,rank;
43849b5e25fSSatish Balay   PetscReal    xl,yl,xr,yr,x_l,x_r,y_l,y_r;
43949b5e25fSSatish Balay   MatScalar    *aa;
44049b5e25fSSatish Balay   MPI_Comm     comm;
44149b5e25fSSatish Balay   Viewer       viewer;
44249b5e25fSSatish Balay 
44349b5e25fSSatish Balay   PetscFunctionBegin;
44449b5e25fSSatish Balay   /*
44549b5e25fSSatish Balay       This is nasty. If this is called from an originally parallel matrix
44649b5e25fSSatish Balay    then all processes call this,but only the first has the matrix so the
44749b5e25fSSatish Balay    rest should return immediately.
44849b5e25fSSatish Balay   */
44949b5e25fSSatish Balay   ierr = PetscObjectGetComm((PetscObject)draw,&comm);CHKERRQ(ierr);
45049b5e25fSSatish Balay   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
45149b5e25fSSatish Balay   if (rank) PetscFunctionReturn(0);
45249b5e25fSSatish Balay 
45349b5e25fSSatish Balay   ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr);
45449b5e25fSSatish Balay 
45549b5e25fSSatish Balay   ierr = DrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr);
45649b5e25fSSatish Balay   DrawString(draw, .3*(xl+xr), .3*(yl+yr), DRAW_BLACK, "symmetric");
45749b5e25fSSatish Balay 
45849b5e25fSSatish Balay   /* loop over matrix elements drawing boxes */
45949b5e25fSSatish Balay   color = DRAW_BLUE;
46049b5e25fSSatish Balay   for (i=0,row=0; i<mbs; i++,row+=bs) {
46149b5e25fSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
46249b5e25fSSatish Balay       y_l = a->m - row - 1.0; y_r = y_l + 1.0;
46349b5e25fSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
46449b5e25fSSatish Balay       aa = a->a + j*bs2;
46549b5e25fSSatish Balay       for (k=0; k<bs; k++) {
46649b5e25fSSatish Balay         for (l=0; l<bs; l++) {
46749b5e25fSSatish Balay           if (PetscRealPart(*aa++) >=  0.) continue;
46849b5e25fSSatish Balay           ierr = DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
46949b5e25fSSatish Balay         }
47049b5e25fSSatish Balay       }
47149b5e25fSSatish Balay     }
47249b5e25fSSatish Balay   }
47349b5e25fSSatish Balay   color = DRAW_CYAN;
47449b5e25fSSatish Balay   for (i=0,row=0; i<mbs; i++,row+=bs) {
47549b5e25fSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
47649b5e25fSSatish Balay       y_l = a->m - row - 1.0; y_r = y_l + 1.0;
47749b5e25fSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
47849b5e25fSSatish Balay       aa = a->a + j*bs2;
47949b5e25fSSatish Balay       for (k=0; k<bs; k++) {
48049b5e25fSSatish Balay         for (l=0; l<bs; l++) {
48149b5e25fSSatish Balay           if (PetscRealPart(*aa++) != 0.) continue;
48249b5e25fSSatish Balay           ierr = DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
48349b5e25fSSatish Balay         }
48449b5e25fSSatish Balay       }
48549b5e25fSSatish Balay     }
48649b5e25fSSatish Balay   }
48749b5e25fSSatish Balay 
48849b5e25fSSatish Balay   color = DRAW_RED;
48949b5e25fSSatish Balay   for (i=0,row=0; i<mbs; i++,row+=bs) {
49049b5e25fSSatish Balay     for (j=a->i[i]; j<a->i[i+1]; j++) {
49149b5e25fSSatish Balay       y_l = a->m - row - 1.0; y_r = y_l + 1.0;
49249b5e25fSSatish Balay       x_l = a->j[j]*bs; x_r = x_l + 1.0;
49349b5e25fSSatish Balay       aa = a->a + j*bs2;
49449b5e25fSSatish Balay       for (k=0; k<bs; k++) {
49549b5e25fSSatish Balay         for (l=0; l<bs; l++) {
49649b5e25fSSatish Balay           if (PetscRealPart(*aa++) <= 0.) continue;
49749b5e25fSSatish Balay           ierr = DrawRectangle(draw,x_l+k,y_l-l,x_r+k,y_r-l,color,color,color,color);CHKERRQ(ierr);
49849b5e25fSSatish Balay         }
49949b5e25fSSatish Balay       }
50049b5e25fSSatish Balay     }
50149b5e25fSSatish Balay   }
50249b5e25fSSatish Balay   PetscFunctionReturn(0);
50349b5e25fSSatish Balay }
50449b5e25fSSatish Balay 
50549b5e25fSSatish Balay #undef __FUNC__
50649b5e25fSSatish Balay #define __FUNC__ "MatView_SeqSBAIJ_Draw"
50749b5e25fSSatish Balay static int MatView_SeqSBAIJ_Draw(Mat A,Viewer viewer)
50849b5e25fSSatish Balay {
50949b5e25fSSatish Balay   Mat_SeqSBAIJ  *a=(Mat_SeqSBAIJ*)A->data;
51049b5e25fSSatish Balay   int          ierr;
51149b5e25fSSatish Balay   PetscReal    xl,yl,xr,yr,w,h;
51249b5e25fSSatish Balay   Draw         draw;
51349b5e25fSSatish Balay   PetscTruth   isnull;
51449b5e25fSSatish Balay 
51549b5e25fSSatish Balay   PetscFunctionBegin;
51649b5e25fSSatish Balay 
51749b5e25fSSatish Balay   ierr = ViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
51849b5e25fSSatish Balay   ierr = DrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
51949b5e25fSSatish Balay 
52049b5e25fSSatish Balay   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
52149b5e25fSSatish Balay   xr  = a->m; yr = a->m; h = yr/10.0; w = xr/10.0;
52249b5e25fSSatish Balay   xr += w;    yr += h;  xl = -w;     yl = -h;
52349b5e25fSSatish Balay   ierr = DrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
52449b5e25fSSatish Balay   ierr = DrawZoom(draw,MatView_SeqSBAIJ_Draw_Zoom,A);CHKERRQ(ierr);
52549b5e25fSSatish Balay   ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr);
52649b5e25fSSatish Balay   PetscFunctionReturn(0);
52749b5e25fSSatish Balay }
52849b5e25fSSatish Balay 
52949b5e25fSSatish Balay #undef __FUNC__
53049b5e25fSSatish Balay #define __FUNC__ "MatView_SeqSBAIJ"
53149b5e25fSSatish Balay int MatView_SeqSBAIJ(Mat A,Viewer viewer)
53249b5e25fSSatish Balay {
53349b5e25fSSatish Balay   int        ierr;
53449b5e25fSSatish Balay   PetscTruth issocket,isascii,isbinary,isdraw;
53549b5e25fSSatish Balay 
53649b5e25fSSatish Balay   PetscFunctionBegin;
53749b5e25fSSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,SOCKET_VIEWER,&issocket);CHKERRQ(ierr);
53849b5e25fSSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,ASCII_VIEWER,&isascii);CHKERRQ(ierr);
53949b5e25fSSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,BINARY_VIEWER,&isbinary);CHKERRQ(ierr);
54049b5e25fSSatish Balay   ierr = PetscTypeCompare((PetscObject)viewer,DRAW_VIEWER,&isdraw);CHKERRQ(ierr);
54149b5e25fSSatish Balay   if (issocket) {
54249b5e25fSSatish Balay     SETERRQ(PETSC_ERR_SUP,0,"Socket viewer not supported");
54349b5e25fSSatish Balay   } else if (isascii){
54449b5e25fSSatish Balay     ierr = MatView_SeqSBAIJ_ASCII(A,viewer);CHKERRQ(ierr);
54549b5e25fSSatish Balay   } else if (isbinary) {
54649b5e25fSSatish Balay     ierr = MatView_SeqSBAIJ_Binary(A,viewer);CHKERRQ(ierr);
54749b5e25fSSatish Balay   } else if (isdraw) {
54849b5e25fSSatish Balay     ierr = MatView_SeqSBAIJ_Draw(A,viewer);CHKERRQ(ierr);
54949b5e25fSSatish Balay   } else {
55049b5e25fSSatish Balay     SETERRQ1(1,1,"Viewer type %s not supported by SeqBAIJ matrices",((PetscObject)viewer)->type_name);
55149b5e25fSSatish Balay   }
55249b5e25fSSatish Balay   PetscFunctionReturn(0);
55349b5e25fSSatish Balay }
55449b5e25fSSatish Balay 
55549b5e25fSSatish Balay 
55649b5e25fSSatish Balay #undef __FUNC__
55749b5e25fSSatish Balay #define __FUNC__ "MatGetValues_SeqSBAIJ"
55849b5e25fSSatish Balay int MatGetValues_SeqSBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v)
55949b5e25fSSatish Balay {
560045c9aa0SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
56149b5e25fSSatish Balay   int        *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j;
56249b5e25fSSatish Balay   int        *ai = a->i,*ailen = a->ilen;
56349b5e25fSSatish Balay   int        brow,bcol,ridx,cidx,bs=a->bs,bs2=a->bs2;
56449b5e25fSSatish Balay   MatScalar  *ap,*aa = a->a,zero = 0.0;
56549b5e25fSSatish Balay 
56649b5e25fSSatish Balay   PetscFunctionBegin;
56749b5e25fSSatish Balay   for (k=0; k<m; k++) { /* loop over rows */
56849b5e25fSSatish Balay     row  = im[k]; brow = row/bs;
56949b5e25fSSatish Balay     if (row < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row");
57049b5e25fSSatish Balay     if (row >= a->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large");
57149b5e25fSSatish Balay     rp   = aj + ai[brow] ; ap = aa + bs2*ai[brow] ;
57249b5e25fSSatish Balay     nrow = ailen[brow];
57349b5e25fSSatish Balay     for (l=0; l<n; l++) { /* loop over columns */
57449b5e25fSSatish Balay       if (in[l] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column");
57549b5e25fSSatish Balay       if (in[l] >= a->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large");
57649b5e25fSSatish Balay       col  = in[l] ;
57749b5e25fSSatish Balay       bcol = col/bs;
57849b5e25fSSatish Balay       cidx = col%bs;
57949b5e25fSSatish Balay       ridx = row%bs;
58049b5e25fSSatish Balay       high = nrow;
58149b5e25fSSatish Balay       low  = 0; /* assume unsorted */
58249b5e25fSSatish Balay       while (high-low > 5) {
58349b5e25fSSatish Balay         t = (low+high)/2;
58449b5e25fSSatish Balay         if (rp[t] > bcol) high = t;
58549b5e25fSSatish Balay         else             low  = t;
58649b5e25fSSatish Balay       }
58749b5e25fSSatish Balay       for (i=low; i<high; i++) {
58849b5e25fSSatish Balay         if (rp[i] > bcol) break;
58949b5e25fSSatish Balay         if (rp[i] == bcol) {
59049b5e25fSSatish Balay           *v++ = ap[bs2*i+bs*cidx+ridx];
59149b5e25fSSatish Balay           goto finished;
59249b5e25fSSatish Balay         }
59349b5e25fSSatish Balay       }
59449b5e25fSSatish Balay       *v++ = zero;
59549b5e25fSSatish Balay       finished:;
59649b5e25fSSatish Balay     }
59749b5e25fSSatish Balay   }
59849b5e25fSSatish Balay   PetscFunctionReturn(0);
59949b5e25fSSatish Balay }
60049b5e25fSSatish Balay 
60149b5e25fSSatish Balay 
60249b5e25fSSatish Balay #undef __FUNC__
60349b5e25fSSatish Balay #define __FUNC__ "MatSetValuesBlocked_SeqSBAIJ"
60449b5e25fSSatish Balay int MatSetValuesBlocked_SeqSBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is)
60549b5e25fSSatish Balay {
60649b5e25fSSatish Balay   PetscFunctionBegin;
607045c9aa0SHong Zhang   SETERRQ(1,1,"Function not yet written for SBAIJ format");
60816cdd363SHong Zhang   /* PetscFunctionReturn(0); */
60949b5e25fSSatish Balay }
61049b5e25fSSatish Balay 
61149b5e25fSSatish Balay #undef __FUNC__
61249b5e25fSSatish Balay #define __FUNC__ "MatAssemblyEnd_SeqSBAIJ"
61349b5e25fSSatish Balay int MatAssemblyEnd_SeqSBAIJ(Mat A,MatAssemblyType mode)
61449b5e25fSSatish Balay {
61549b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
61649b5e25fSSatish Balay   int        fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax;
61749b5e25fSSatish Balay   int        m = a->m,*ip,N,*ailen = a->ilen;
61849b5e25fSSatish Balay   int        mbs = a->mbs,bs2 = a->bs2,rmax = 0,ierr;
61949b5e25fSSatish Balay   MatScalar  *aa = a->a,*ap;
62049b5e25fSSatish Balay 
62149b5e25fSSatish Balay   PetscFunctionBegin;
62249b5e25fSSatish Balay   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0);
62349b5e25fSSatish Balay 
62449b5e25fSSatish Balay   if (m) rmax = ailen[0];
62549b5e25fSSatish Balay   for (i=1; i<mbs; i++) {
62649b5e25fSSatish Balay     /* move each row back by the amount of empty slots (fshift) before it*/
62749b5e25fSSatish Balay     fshift += imax[i-1] - ailen[i-1];
62849b5e25fSSatish Balay     rmax   = PetscMax(rmax,ailen[i]);
62949b5e25fSSatish Balay     if (fshift) {
63049b5e25fSSatish Balay       ip = aj + ai[i]; ap = aa + bs2*ai[i];
63149b5e25fSSatish Balay       N = ailen[i];
63249b5e25fSSatish Balay       for (j=0; j<N; j++) {
63349b5e25fSSatish Balay         ip[j-fshift] = ip[j];
63449b5e25fSSatish Balay         ierr = PetscMemcpy(ap+(j-fshift)*bs2,ap+j*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
63549b5e25fSSatish Balay       }
63649b5e25fSSatish Balay     }
63749b5e25fSSatish Balay     ai[i] = ai[i-1] + ailen[i-1];
63849b5e25fSSatish Balay   }
63949b5e25fSSatish Balay   if (mbs) {
64049b5e25fSSatish Balay     fshift += imax[mbs-1] - ailen[mbs-1];
64149b5e25fSSatish Balay     ai[mbs] = ai[mbs-1] + ailen[mbs-1];
64249b5e25fSSatish Balay   }
64349b5e25fSSatish Balay   /* reset ilen and imax for each row */
64449b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
64549b5e25fSSatish Balay     ailen[i] = imax[i] = ai[i+1] - ai[i];
64649b5e25fSSatish Balay   }
64749b5e25fSSatish Balay   a->s_nz = ai[mbs];
64849b5e25fSSatish Balay 
64949b5e25fSSatish Balay   /* diagonals may have moved, so kill the diagonal pointers */
65049b5e25fSSatish Balay   if (fshift && a->diag) {
65149b5e25fSSatish Balay     ierr = PetscFree(a->diag);CHKERRQ(ierr);
65249b5e25fSSatish Balay     PLogObjectMemory(A,-(m+1)*sizeof(int));
65349b5e25fSSatish Balay     a->diag = 0;
65449b5e25fSSatish Balay   }
65549b5e25fSSatish Balay   PLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Matrix size: %d X %d, block size %d; storage space: %d unneeded, %d used\n",
65649b5e25fSSatish Balay            m,a->m,a->bs,fshift*bs2,a->s_nz*bs2);
65749b5e25fSSatish Balay   PLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Number of mallocs during MatSetValues is %d\n",
65849b5e25fSSatish Balay            a->reallocs);
65949b5e25fSSatish Balay   PLogInfo(A,"MatAssemblyEnd_SeqSBAIJ:Most nonzeros blocks in any row is %d\n",rmax);
66049b5e25fSSatish Balay   a->reallocs          = 0;
66149b5e25fSSatish Balay   A->info.nz_unneeded  = (PetscReal)fshift*bs2;
66249b5e25fSSatish Balay 
66349b5e25fSSatish Balay   PetscFunctionReturn(0);
66449b5e25fSSatish Balay }
66549b5e25fSSatish Balay 
66649b5e25fSSatish Balay /*
66749b5e25fSSatish Balay    This function returns an array of flags which indicate the locations of contiguous
66849b5e25fSSatish Balay    blocks that should be zeroed. for eg: if bs = 3  and is = [0,1,2,3,5,6,7,8,9]
66949b5e25fSSatish Balay    then the resulting sizes = [3,1,1,3,1] correspondig to sets [(0,1,2),(3),(5),(6,7,8),(9)]
67049b5e25fSSatish Balay    Assume: sizes should be long enough to hold all the values.
67149b5e25fSSatish Balay */
67249b5e25fSSatish Balay #undef __FUNC__
67349b5e25fSSatish Balay #define __FUNC__ "MatZeroRows_SeqSBAIJ_Check_Blocks"
67449b5e25fSSatish Balay static int MatZeroRows_SeqSBAIJ_Check_Blocks(int idx[],int n,int bs,int sizes[], int *bs_max)
67549b5e25fSSatish Balay {
67649b5e25fSSatish Balay   int        i,j,k,row;
67749b5e25fSSatish Balay   PetscTruth flg;
67849b5e25fSSatish Balay 
67949b5e25fSSatish Balay   PetscFunctionBegin;
68049b5e25fSSatish Balay   for (i=0,j=0; i<n; j++) {
68149b5e25fSSatish Balay     row = idx[i];
68249b5e25fSSatish Balay     if (row%bs!=0) { /* Not the begining of a block */
68349b5e25fSSatish Balay       sizes[j] = 1;
68449b5e25fSSatish Balay       i++;
68549b5e25fSSatish Balay     } else if (i+bs > n) { /* Beginning of a block, but complete block doesn't exist (at idx end) */
68649b5e25fSSatish Balay       sizes[j] = 1;         /* Also makes sure atleast 'bs' values exist for next else */
68749b5e25fSSatish Balay       i++;
68849b5e25fSSatish Balay     } else { /* Begining of the block, so check if the complete block exists */
68949b5e25fSSatish Balay       flg = PETSC_TRUE;
69049b5e25fSSatish Balay       for (k=1; k<bs; k++) {
69149b5e25fSSatish Balay         if (row+k != idx[i+k]) { /* break in the block */
69249b5e25fSSatish Balay           flg = PETSC_FALSE;
69349b5e25fSSatish Balay           break;
69449b5e25fSSatish Balay         }
69549b5e25fSSatish Balay       }
69649b5e25fSSatish Balay       if (flg == PETSC_TRUE) { /* No break in the bs */
69749b5e25fSSatish Balay         sizes[j] = bs;
69849b5e25fSSatish Balay         i+= bs;
69949b5e25fSSatish Balay       } else {
70049b5e25fSSatish Balay         sizes[j] = 1;
70149b5e25fSSatish Balay         i++;
70249b5e25fSSatish Balay       }
70349b5e25fSSatish Balay     }
70449b5e25fSSatish Balay   }
70549b5e25fSSatish Balay   *bs_max = j;
70649b5e25fSSatish Balay   PetscFunctionReturn(0);
70749b5e25fSSatish Balay }
70849b5e25fSSatish Balay 
70949b5e25fSSatish Balay #undef __FUNC__
71049b5e25fSSatish Balay #define __FUNC__ "MatZeroRows_SeqSBAIJ"
71149b5e25fSSatish Balay int MatZeroRows_SeqSBAIJ(Mat A,IS is,Scalar *diag)
71249b5e25fSSatish Balay {
71349b5e25fSSatish Balay   Mat_SeqSBAIJ *sbaij=(Mat_SeqSBAIJ*)A->data;
71449b5e25fSSatish Balay   int         ierr,i,j,k,count,is_n,*is_idx,*rows;
71549b5e25fSSatish Balay   int         bs=sbaij->bs,bs2=sbaij->bs2,*sizes,row,bs_max;
71649b5e25fSSatish Balay   Scalar      zero = 0.0;
71749b5e25fSSatish Balay   MatScalar   *aa;
71849b5e25fSSatish Balay 
71949b5e25fSSatish Balay   PetscFunctionBegin;
72049b5e25fSSatish Balay   /* Make a copy of the IS and  sort it */
72149b5e25fSSatish Balay   ierr = ISGetSize(is,&is_n);CHKERRQ(ierr);
72249b5e25fSSatish Balay   ierr = ISGetIndices(is,&is_idx);CHKERRQ(ierr);
72349b5e25fSSatish Balay 
72449b5e25fSSatish Balay   /* allocate memory for rows,sizes */
72549b5e25fSSatish Balay   rows  = (int*)PetscMalloc((3*is_n+1)*sizeof(int));CHKPTRQ(rows);
72649b5e25fSSatish Balay   sizes = rows + is_n;
72749b5e25fSSatish Balay 
72849b5e25fSSatish Balay   /* initialize copy IS values to rows, and sort them */
72949b5e25fSSatish Balay   for (i=0; i<is_n; i++) { rows[i] = is_idx[i]; }
73049b5e25fSSatish Balay   ierr = PetscSortInt(is_n,rows);CHKERRQ(ierr);
73149b5e25fSSatish Balay   if (sbaij->keepzeroedrows) { /* do not change nonzero structure */
73249b5e25fSSatish Balay     for (i=0; i<is_n; i++) { sizes[i] = 1; } /* sizes: size of blocks, = 1 or bs */
73349b5e25fSSatish Balay     bs_max = is_n;            /* bs_max: num. of contiguous block row in the row */
73449b5e25fSSatish Balay   } else {
73549b5e25fSSatish Balay     ierr = MatZeroRows_SeqSBAIJ_Check_Blocks(rows,is_n,bs,sizes,&bs_max);CHKERRQ(ierr);
73649b5e25fSSatish Balay   }
73749b5e25fSSatish Balay   ierr = ISRestoreIndices(is,&is_idx);CHKERRQ(ierr);
73849b5e25fSSatish Balay 
73949b5e25fSSatish Balay   for (i=0,j=0; i<bs_max; j+=sizes[i],i++) {
74049b5e25fSSatish Balay     row   = rows[j];                  /* row to be zeroed */
74149b5e25fSSatish Balay     if (row < 0 || row > sbaij->m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,0,"row %d out of range",row);
74249b5e25fSSatish Balay     count = (sbaij->i[row/bs +1] - sbaij->i[row/bs])*bs; /* num. of elements in the row */
74349b5e25fSSatish Balay     aa    = sbaij->a + sbaij->i[row/bs]*bs2 + (row%bs);
74449b5e25fSSatish Balay     if (sizes[i] == bs && !sbaij->keepzeroedrows) {
74549b5e25fSSatish Balay       if (diag) {
74649b5e25fSSatish Balay         if (sbaij->ilen[row/bs] > 0) {
74749b5e25fSSatish Balay           sbaij->ilen[row/bs] = 1;
74849b5e25fSSatish Balay           sbaij->j[sbaij->i[row/bs]] = row/bs;
74949b5e25fSSatish Balay           ierr = PetscMemzero(aa,count*bs*sizeof(MatScalar));CHKERRQ(ierr);
75049b5e25fSSatish Balay         }
75149b5e25fSSatish Balay         /* Now insert all the diagoanl values for this bs */
75249b5e25fSSatish Balay         for (k=0; k<bs; k++) {
75349b5e25fSSatish Balay           ierr = (*A->ops->setvalues)(A,1,rows+j+k,1,rows+j+k,diag,INSERT_VALUES);CHKERRQ(ierr);
75449b5e25fSSatish Balay         }
75549b5e25fSSatish Balay       } else { /* (!diag) */
75649b5e25fSSatish Balay         sbaij->ilen[row/bs] = 0;
75749b5e25fSSatish Balay       } /* end (!diag) */
75849b5e25fSSatish Balay     } else { /* (sizes[i] != bs), broken block */
75949b5e25fSSatish Balay #if defined (PETSC_USE_DEBUG)
76049b5e25fSSatish Balay       if (sizes[i] != 1) SETERRQ(1,0,"Internal Error. Value should be 1");
76149b5e25fSSatish Balay #endif
76249b5e25fSSatish Balay       for (k=0; k<count; k++) {
76349b5e25fSSatish Balay         aa[0] = zero;
76449b5e25fSSatish Balay         aa+=bs;
76549b5e25fSSatish Balay       }
76649b5e25fSSatish Balay       if (diag) {
76749b5e25fSSatish Balay         ierr = (*A->ops->setvalues)(A,1,rows+j,1,rows+j,diag,INSERT_VALUES);CHKERRQ(ierr);
76849b5e25fSSatish Balay       }
76949b5e25fSSatish Balay     }
77049b5e25fSSatish Balay   }
77149b5e25fSSatish Balay 
77249b5e25fSSatish Balay   ierr = PetscFree(rows);CHKERRQ(ierr);
77349b5e25fSSatish Balay   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
77449b5e25fSSatish Balay   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
77549b5e25fSSatish Balay   PetscFunctionReturn(0);
77649b5e25fSSatish Balay }
77749b5e25fSSatish Balay 
77849b5e25fSSatish Balay /* Only add/insert a(i,j) with i<=j (blocks).
77949b5e25fSSatish Balay    Any a(i,j) with i>j input by user is ingored.
78049b5e25fSSatish Balay */
78149b5e25fSSatish Balay 
78249b5e25fSSatish Balay #undef __FUNC__
78349b5e25fSSatish Balay #define __FUNC__ "MatSetValues_SeqSBAIJ"
78449b5e25fSSatish Balay int MatSetValues_SeqSBAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is)
78549b5e25fSSatish Balay {
78649b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
78749b5e25fSSatish Balay   int         *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N,sorted=a->sorted;
78849b5e25fSSatish Balay   int         *imax=a->imax,*ai=a->i,*ailen=a->ilen,roworiented=a->roworiented;
78949b5e25fSSatish Balay   int         *aj=a->j,nonew=a->nonew,bs=a->bs,brow,bcol;
79049b5e25fSSatish Balay   int         ridx,cidx,bs2=a->bs2,ierr;
79149b5e25fSSatish Balay   MatScalar   *ap,value,*aa=a->a,*bap;
79249b5e25fSSatish Balay 
79349b5e25fSSatish Balay   PetscFunctionBegin;
79449b5e25fSSatish Balay 
79549b5e25fSSatish Balay   for (k=0; k<m; k++) { /* loop over added rows */
79649b5e25fSSatish Balay     row  = im[k];       /* row number */
79749b5e25fSSatish Balay     brow = row/bs;      /* block row number */
79849b5e25fSSatish Balay     if (row < 0) continue;
79949b5e25fSSatish Balay #if defined(PETSC_USE_BOPT_g)
80049b5e25fSSatish Balay     if (row >= a->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large: row %d max %d",row,a->m);
80149b5e25fSSatish Balay #endif
80249b5e25fSSatish Balay     rp   = aj + ai[brow]; /*ptr to beginning of column value of the row block*/
80349b5e25fSSatish Balay     ap   = aa + bs2*ai[brow]; /*ptr to beginning of element value of the row block*/
80449b5e25fSSatish Balay     rmax = imax[brow];  /* maximum space allocated for this row */
80549b5e25fSSatish Balay     nrow = ailen[brow]; /* actual length of this row */
80649b5e25fSSatish Balay     low  = 0;
80749b5e25fSSatish Balay 
80849b5e25fSSatish Balay     for (l=0; l<n; l++) { /* loop over added columns */
80949b5e25fSSatish Balay       if (in[l] < 0) continue;
81049b5e25fSSatish Balay #if defined(PETSC_USE_BOPT_g)
81149b5e25fSSatish Balay       if (in[l] >= a->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large: col %d max %d",in[l],a->m);
81249b5e25fSSatish Balay #endif
81349b5e25fSSatish Balay       col = in[l];
81449b5e25fSSatish Balay       bcol = col/bs;              /* block col number */
81549b5e25fSSatish Balay 
81649b5e25fSSatish Balay       if (brow <= bcol){
81749b5e25fSSatish Balay         ridx = row % bs; cidx = col % bs; /*row and col index inside the block */
81849b5e25fSSatish Balay         /* if ((brow==bcol && ridx<=cidx) || (brow<bcol)){ */
81949b5e25fSSatish Balay           /* element value a(k,l) */
82049b5e25fSSatish Balay           if (roworiented) {
82149b5e25fSSatish Balay             value = v[l + k*n];
82249b5e25fSSatish Balay           } else {
82349b5e25fSSatish Balay             value = v[k + l*m];
82449b5e25fSSatish Balay           }
82549b5e25fSSatish Balay 
82649b5e25fSSatish Balay           /* move pointer bap to a(k,l) quickly and add/insert value */
82749b5e25fSSatish Balay           if (!sorted) low = 0; high = nrow;
82849b5e25fSSatish Balay           while (high-low > 7) {
82949b5e25fSSatish Balay             t = (low+high)/2;
83049b5e25fSSatish Balay             if (rp[t] > bcol) high = t;
83149b5e25fSSatish Balay             else              low  = t;
83249b5e25fSSatish Balay           }
83349b5e25fSSatish Balay           for (i=low; i<high; i++) {
83449b5e25fSSatish Balay             /* printf("The loop of i=low.., rp[%d]=%d\n",i,rp[i]); */
83549b5e25fSSatish Balay             if (rp[i] > bcol) break;
83649b5e25fSSatish Balay             if (rp[i] == bcol) {
83749b5e25fSSatish Balay               bap  = ap +  bs2*i + bs*cidx + ridx;
83849b5e25fSSatish Balay               if (is == ADD_VALUES) *bap += value;
83949b5e25fSSatish Balay               else                  *bap  = value;
84049b5e25fSSatish Balay               goto noinsert1;
84149b5e25fSSatish Balay             }
84249b5e25fSSatish Balay           }
84349b5e25fSSatish Balay 
84449b5e25fSSatish Balay       if (nonew == 1) goto noinsert1;
84549b5e25fSSatish Balay       else if (nonew == -1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix");
84649b5e25fSSatish Balay       if (nrow >= rmax) {
84749b5e25fSSatish Balay         /* there is no extra room in row, therefore enlarge */
84849b5e25fSSatish Balay         int       new_nz = ai[a->mbs] + CHUNKSIZE,len,*new_i,*new_j;
84949b5e25fSSatish Balay         MatScalar *new_a;
85049b5e25fSSatish Balay 
85149b5e25fSSatish Balay         if (nonew == -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Inserting a new nonzero in the matrix");
85249b5e25fSSatish Balay 
85349b5e25fSSatish Balay         /* Malloc new storage space */
85449b5e25fSSatish Balay         len     = new_nz*(sizeof(int)+bs2*sizeof(MatScalar))+(a->mbs+1)*sizeof(int);
85549b5e25fSSatish Balay         new_a   = (MatScalar*)PetscMalloc(len);CHKPTRQ(new_a);
85649b5e25fSSatish Balay         new_j   = (int*)(new_a + bs2*new_nz);
85749b5e25fSSatish Balay         new_i   = new_j + new_nz;
85849b5e25fSSatish Balay 
85949b5e25fSSatish Balay         /* copy over old data into new slots */
86049b5e25fSSatish Balay         for (ii=0; ii<brow+1; ii++) {new_i[ii] = ai[ii];}
86149b5e25fSSatish Balay         for (ii=brow+1; ii<a->mbs+1; ii++) {new_i[ii] = ai[ii]+CHUNKSIZE;}
86249b5e25fSSatish Balay         ierr = PetscMemcpy(new_j,aj,(ai[brow]+nrow)*sizeof(int));CHKERRQ(ierr);
86349b5e25fSSatish Balay         len  = (new_nz - CHUNKSIZE - ai[brow] - nrow);
86449b5e25fSSatish Balay         ierr = PetscMemcpy(new_j+ai[brow]+nrow+CHUNKSIZE,aj+ai[brow]+nrow,len*sizeof(int));CHKERRQ(ierr);
86549b5e25fSSatish Balay         ierr = PetscMemcpy(new_a,aa,(ai[brow]+nrow)*bs2*sizeof(MatScalar));CHKERRQ(ierr);
86649b5e25fSSatish Balay         ierr = PetscMemzero(new_a+bs2*(ai[brow]+nrow),bs2*CHUNKSIZE*sizeof(MatScalar));CHKERRQ(ierr);
86749b5e25fSSatish Balay         ierr = PetscMemcpy(new_a+bs2*(ai[brow]+nrow+CHUNKSIZE),aa+bs2*(ai[brow]+nrow),bs2*len*sizeof(MatScalar));CHKERRQ(ierr);
86849b5e25fSSatish Balay         /* free up old matrix storage */
86949b5e25fSSatish Balay         ierr = PetscFree(a->a);CHKERRQ(ierr);
87049b5e25fSSatish Balay         if (!a->singlemalloc) {
87149b5e25fSSatish Balay           ierr = PetscFree(a->i);CHKERRQ(ierr);
87249b5e25fSSatish Balay           ierr = PetscFree(a->j);CHKERRQ(ierr);
87349b5e25fSSatish Balay         }
87449b5e25fSSatish Balay         aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j;
87549b5e25fSSatish Balay         a->singlemalloc = PETSC_TRUE;
87649b5e25fSSatish Balay 
87749b5e25fSSatish Balay         rp   = aj + ai[brow]; ap = aa + bs2*ai[brow];
87849b5e25fSSatish Balay         rmax = imax[brow] = imax[brow] + CHUNKSIZE;
87949b5e25fSSatish Balay         PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + bs2*sizeof(MatScalar)));
88049b5e25fSSatish Balay         a->s_maxnz += bs2*CHUNKSIZE;
88149b5e25fSSatish Balay         a->reallocs++;
88249b5e25fSSatish Balay         a->s_nz++;
88349b5e25fSSatish Balay       }
88449b5e25fSSatish Balay 
88549b5e25fSSatish Balay       N = nrow++ - 1;
88649b5e25fSSatish Balay       /* shift up all the later entries in this row */
88749b5e25fSSatish Balay       for (ii=N; ii>=i; ii--) {
88849b5e25fSSatish Balay         rp[ii+1] = rp[ii];
88949b5e25fSSatish Balay         ierr     = PetscMemcpy(ap+bs2*(ii+1),ap+bs2*(ii),bs2*sizeof(MatScalar));CHKERRQ(ierr);
89049b5e25fSSatish Balay       }
89149b5e25fSSatish Balay       if (N>=i) {
89249b5e25fSSatish Balay         ierr = PetscMemzero(ap+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr);
89349b5e25fSSatish Balay       }
89449b5e25fSSatish Balay       rp[i]                      = bcol;
89549b5e25fSSatish Balay       ap[bs2*i + bs*cidx + ridx] = value;
89649b5e25fSSatish Balay       noinsert1:;
89749b5e25fSSatish Balay       low = i;
89849b5e25fSSatish Balay       /* } */
89949b5e25fSSatish Balay     } /* end of if .. if..  */
90049b5e25fSSatish Balay     }                     /* end of loop over added columns */
90149b5e25fSSatish Balay     ailen[brow] = nrow;
90249b5e25fSSatish Balay   }                       /* end of loop over added rows */
90349b5e25fSSatish Balay 
90449b5e25fSSatish Balay   PetscFunctionReturn(0);
90549b5e25fSSatish Balay }
90649b5e25fSSatish Balay 
907915354c9SHong Zhang extern int MatCholeskyFactorSymbolic_SeqSBAIJ(Mat,IS,PetscReal,Mat*);
908915354c9SHong Zhang extern int MatCholeskyFactor_SeqSBAIJ(Mat,IS,PetscReal);
90949b5e25fSSatish Balay extern int MatIncreaseOverlap_SeqSBAIJ(Mat,int,IS*,int);
91049b5e25fSSatish Balay extern int MatGetSubMatrix_SeqSBAIJ(Mat,IS,IS,int,MatReuse,Mat*);
91149b5e25fSSatish Balay extern int MatGetSubMatrices_SeqSBAIJ(Mat,int,IS*,IS*,MatReuse,Mat**);
91249b5e25fSSatish Balay extern int MatMultTranspose_SeqSBAIJ(Mat,Vec,Vec);
91349b5e25fSSatish Balay extern int MatMultTransposeAdd_SeqSBAIJ(Mat,Vec,Vec,Vec);
91449b5e25fSSatish Balay extern int MatScale_SeqSBAIJ(Scalar*,Mat);
91549b5e25fSSatish Balay extern int MatNorm_SeqSBAIJ(Mat,NormType,PetscReal *);
91649b5e25fSSatish Balay extern int MatEqual_SeqSBAIJ(Mat,Mat,PetscTruth*);
91749b5e25fSSatish Balay extern int MatGetDiagonal_SeqSBAIJ(Mat,Vec);
91849b5e25fSSatish Balay extern int MatDiagonalScale_SeqSBAIJ(Mat,Vec,Vec);
91949b5e25fSSatish Balay extern int MatGetInfo_SeqSBAIJ(Mat,MatInfoType,MatInfo *);
92049b5e25fSSatish Balay extern int MatZeroEntries_SeqSBAIJ(Mat);
92149b5e25fSSatish Balay 
92249b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_N(Mat,Vec,Vec);
92349b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_1(Mat,Vec,Vec);
92449b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_2(Mat,Vec,Vec);
92549b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_3(Mat,Vec,Vec);
92649b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_4(Mat,Vec,Vec);
92749b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_5(Mat,Vec,Vec);
92849b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_6(Mat,Vec,Vec);
92949b5e25fSSatish Balay extern int MatSolve_SeqSBAIJ_7(Mat,Vec,Vec);
93049b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_7(Mat,Vec,Vec);
93149b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_6(Mat,Vec,Vec);
93249b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_5(Mat,Vec,Vec);
93349b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_4(Mat,Vec,Vec);
93449b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_3(Mat,Vec,Vec);
93549b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_2(Mat,Vec,Vec);
93649b5e25fSSatish Balay extern int MatSolveTranspose_SeqSBAIJ_1(Mat,Vec,Vec);
93749b5e25fSSatish Balay 
9385f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat,Mat*);
9395f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_1(Mat,Mat*);
9405f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_2(Mat,Mat*);
9415f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_3(Mat,Mat*);
9425f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_4(Mat,Mat*);
9435f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_5(Mat,Mat*);
9445f19b6beSHong Zhang extern int MatCholeskyFactorNumeric_SeqSBAIJ_6(Mat,Mat*);
94549b5e25fSSatish Balay 
94649b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_1(Mat,Vec,Vec);
94749b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_2(Mat,Vec,Vec);
94849b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_3(Mat,Vec,Vec);
94949b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_4(Mat,Vec,Vec);
95049b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_5(Mat,Vec,Vec);
95149b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_6(Mat,Vec,Vec);
95249b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_7(Mat,Vec,Vec);
95349b5e25fSSatish Balay extern int MatMult_SeqSBAIJ_N(Mat,Vec,Vec);
95449b5e25fSSatish Balay 
95549b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_1(Mat,Vec,Vec,Vec);
95649b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_2(Mat,Vec,Vec,Vec);
95749b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_3(Mat,Vec,Vec,Vec);
95849b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_4(Mat,Vec,Vec,Vec);
95949b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_5(Mat,Vec,Vec,Vec);
96049b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_6(Mat,Vec,Vec,Vec);
96149b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_7(Mat,Vec,Vec,Vec);
96249b5e25fSSatish Balay extern int MatMultAdd_SeqSBAIJ_N(Mat,Vec,Vec,Vec);
96349b5e25fSSatish Balay 
964950f1e5bSHong Zhang #ifdef MatIncompleteCholeskyFactor
965950f1e5bSHong Zhang /* This function is modified from MatILUFactor_SeqSBAIJ.
966950f1e5bSHong Zhang    Needs further work! Don't forget to add the function to the matrix table. */
96749b5e25fSSatish Balay #undef __FUNC__
9684ccecd49SHong Zhang #define __FUNC__ "MatIncompleteCholeskyFactor_SeqSBAIJ"
9694ccecd49SHong Zhang int MatIncompleteCholeskyFactor_SeqSBAIJ(Mat inA,IS row,IS col,MatILUInfo *info)
97049b5e25fSSatish Balay {
9714ccecd49SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)inA->data;
97249b5e25fSSatish Balay   Mat         outA;
97349b5e25fSSatish Balay   int         ierr;
97449b5e25fSSatish Balay   PetscTruth  row_identity,col_identity;
97549b5e25fSSatish Balay 
97649b5e25fSSatish Balay   PetscFunctionBegin;
97749b5e25fSSatish Balay   if (info && info->levels != 0) SETERRQ(PETSC_ERR_SUP,0,"Only levels = 0 supported for in-place ILU");
97849b5e25fSSatish Balay   ierr = ISIdentity(row,&row_identity);CHKERRQ(ierr);
97949b5e25fSSatish Balay   ierr = ISIdentity(col,&col_identity);CHKERRQ(ierr);
98049b5e25fSSatish Balay   if (!row_identity || !col_identity) {
98149b5e25fSSatish Balay     SETERRQ(1,1,"Row and column permutations must be identity for in-place ILU");
98249b5e25fSSatish Balay   }
98349b5e25fSSatish Balay 
98449b5e25fSSatish Balay   outA          = inA;
98549b5e25fSSatish Balay   inA->factor   = FACTOR_LU;
98649b5e25fSSatish Balay 
98749b5e25fSSatish Balay   if (!a->diag) {
98849b5e25fSSatish Balay     ierr = MatMarkDiagonal_SeqBAIJ(inA);CHKERRQ(ierr);
98949b5e25fSSatish Balay   }
99049b5e25fSSatish Balay   /*
99149b5e25fSSatish Balay       Blocksize 2, 3, 4, 5, 6 and 7 have a special faster factorization/solver
99249b5e25fSSatish Balay       for ILU(0) factorization with natural ordering
99349b5e25fSSatish Balay   */
99449b5e25fSSatish Balay   switch (a->bs) {
99549b5e25fSSatish Balay   case 1:
99649b5e25fSSatish Balay     inA->ops->solvetranspose   = MatSolveTranspose_SeqBAIJ_2_NaturalOrdering;
99749b5e25fSSatish Balay     PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering solvetrans BS=1\n");
99849b5e25fSSatish Balay   case 2:
99949b5e25fSSatish Balay     inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering;
100049b5e25fSSatish Balay     inA->ops->solve           = MatSolve_SeqBAIJ_2_NaturalOrdering;
100149b5e25fSSatish Balay     inA->ops->solvetranspose  = MatSolveTranspose_SeqBAIJ_2_NaturalOrdering;
100249b5e25fSSatish Balay     PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=2\n");
100349b5e25fSSatish Balay     break;
100449b5e25fSSatish Balay   case 3:
100549b5e25fSSatish Balay     inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering;
100649b5e25fSSatish Balay     inA->ops->solve           = MatSolve_SeqBAIJ_3_NaturalOrdering;
100749b5e25fSSatish Balay     inA->ops->solvetranspose  = MatSolveTranspose_SeqBAIJ_3_NaturalOrdering;
100849b5e25fSSatish Balay     PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=3\n");
100949b5e25fSSatish Balay     break;
101049b5e25fSSatish Balay   case 4:
101149b5e25fSSatish Balay     inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_4_NaturalOrdering;
101249b5e25fSSatish Balay     inA->ops->solve           = MatSolve_SeqBAIJ_4_NaturalOrdering;
101349b5e25fSSatish Balay     inA->ops->solvetranspose  = MatSolveTranspose_SeqBAIJ_4_NaturalOrdering;
101449b5e25fSSatish Balay     PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=4\n");
101549b5e25fSSatish Balay     break;
101649b5e25fSSatish Balay   case 5:
101749b5e25fSSatish Balay     inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_5_NaturalOrdering;
101849b5e25fSSatish Balay     inA->ops->solve           = MatSolve_SeqBAIJ_5_NaturalOrdering;
101949b5e25fSSatish Balay     inA->ops->solvetranspose  = MatSolveTranspose_SeqBAIJ_5_NaturalOrdering;
102049b5e25fSSatish Balay     PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=5\n");
102149b5e25fSSatish Balay     break;
102249b5e25fSSatish Balay   case 6:
102349b5e25fSSatish Balay     inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_6_NaturalOrdering;
102449b5e25fSSatish Balay     inA->ops->solve           = MatSolve_SeqBAIJ_6_NaturalOrdering;
102549b5e25fSSatish Balay     inA->ops->solvetranspose  = MatSolveTranspose_SeqBAIJ_6_NaturalOrdering;
102649b5e25fSSatish Balay     PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=6\n");
102749b5e25fSSatish Balay     break;
102849b5e25fSSatish Balay   case 7:
102949b5e25fSSatish Balay     inA->ops->lufactornumeric = MatLUFactorNumeric_SeqBAIJ_7_NaturalOrdering;
103049b5e25fSSatish Balay     inA->ops->solvetranspose  = MatSolveTranspose_SeqBAIJ_7_NaturalOrdering;
103149b5e25fSSatish Balay     inA->ops->solve           = MatSolve_SeqBAIJ_7_NaturalOrdering;
103249b5e25fSSatish Balay     PLogInfo(inA,"MatILUFactor_SeqBAIJ:Using special in-place natural ordering factor and solve BS=7\n");
103349b5e25fSSatish Balay     break;
103449b5e25fSSatish Balay   default:
103549b5e25fSSatish Balay     a->row        = row;
103649b5e25fSSatish Balay     a->col        = col;
103749b5e25fSSatish Balay     ierr          = PetscObjectReference((PetscObject)row);CHKERRQ(ierr);
103849b5e25fSSatish Balay     ierr          = PetscObjectReference((PetscObject)col);CHKERRQ(ierr);
103949b5e25fSSatish Balay 
104049b5e25fSSatish Balay     /* Create the invert permutation so that it can be used in MatLUFactorNumeric() */
104149b5e25fSSatish Balay     ierr = ISInvertPermutation(col,PETSC_DECIDE, &(a->icol));CHKERRQ(ierr);
104249b5e25fSSatish Balay     PLogObjectParent(inA,a->icol);
104349b5e25fSSatish Balay 
104449b5e25fSSatish Balay     if (!a->solve_work) {
104549b5e25fSSatish Balay       a->solve_work = (Scalar*)PetscMalloc((a->m+a->bs)*sizeof(Scalar));CHKPTRQ(a->solve_work);
104649b5e25fSSatish Balay       PLogObjectMemory(inA,(a->m+a->bs)*sizeof(Scalar));
104749b5e25fSSatish Balay     }
104849b5e25fSSatish Balay   }
104949b5e25fSSatish Balay 
105049b5e25fSSatish Balay   ierr = MatLUFactorNumeric(inA,&outA);CHKERRQ(ierr);
105149b5e25fSSatish Balay 
105249b5e25fSSatish Balay   PetscFunctionReturn(0);
105349b5e25fSSatish Balay }
1054950f1e5bSHong Zhang #endif
1055950f1e5bSHong Zhang 
105649b5e25fSSatish Balay #undef __FUNC__
105749b5e25fSSatish Balay #define __FUNC__ "MatPrintHelp_SeqSBAIJ"
105849b5e25fSSatish Balay int MatPrintHelp_SeqSBAIJ(Mat A)
105949b5e25fSSatish Balay {
106049b5e25fSSatish Balay   static PetscTruth called = PETSC_FALSE;
106149b5e25fSSatish Balay   MPI_Comm          comm = A->comm;
106249b5e25fSSatish Balay   int               ierr;
106349b5e25fSSatish Balay 
106449b5e25fSSatish Balay   PetscFunctionBegin;
106549b5e25fSSatish Balay   if (called) {PetscFunctionReturn(0);} else called = PETSC_TRUE;
106649b5e25fSSatish Balay   ierr = (*PetscHelpPrintf)(comm," Options for MATSEQSBAIJ and MATMPISBAIJ matrix formats (the defaults):\n");CHKERRQ(ierr);
106749b5e25fSSatish Balay   ierr = (*PetscHelpPrintf)(comm,"  -mat_block_size <block_size>\n");CHKERRQ(ierr);
106849b5e25fSSatish Balay   PetscFunctionReturn(0);
106949b5e25fSSatish Balay }
107049b5e25fSSatish Balay 
107149b5e25fSSatish Balay EXTERN_C_BEGIN
107249b5e25fSSatish Balay #undef __FUNC__
107349b5e25fSSatish Balay #define __FUNC__ "MatSeqSBAIJSetColumnIndices_SeqSBAIJ"
107449b5e25fSSatish Balay int MatSeqSBAIJSetColumnIndices_SeqSBAIJ(Mat mat,int *indices)
107549b5e25fSSatish Balay {
1076045c9aa0SHong Zhang   Mat_SeqSBAIJ *baij = (Mat_SeqSBAIJ *)mat->data;
107749b5e25fSSatish Balay   int         i,nz,n;
107849b5e25fSSatish Balay 
107949b5e25fSSatish Balay   PetscFunctionBegin;
1080045c9aa0SHong Zhang   nz = baij->s_maxnz;
108149b5e25fSSatish Balay   n  = baij->n;
108249b5e25fSSatish Balay   for (i=0; i<nz; i++) {
108349b5e25fSSatish Balay     baij->j[i] = indices[i];
108449b5e25fSSatish Balay   }
1085045c9aa0SHong Zhang   baij->s_nz = nz;
108649b5e25fSSatish Balay   for (i=0; i<n; i++) {
108749b5e25fSSatish Balay     baij->ilen[i] = baij->imax[i];
108849b5e25fSSatish Balay   }
108949b5e25fSSatish Balay 
109049b5e25fSSatish Balay   PetscFunctionReturn(0);
109149b5e25fSSatish Balay }
109249b5e25fSSatish Balay EXTERN_C_END
109349b5e25fSSatish Balay 
109449b5e25fSSatish Balay #undef __FUNC__
109549b5e25fSSatish Balay #define __FUNC__ "MatSeqSBAIJSetColumnIndices"
109649b5e25fSSatish Balay /*@
109719585528SSatish Balay     MatSeqSBAIJSetColumnIndices - Set the column indices for all the rows
109849b5e25fSSatish Balay        in the matrix.
109949b5e25fSSatish Balay 
110049b5e25fSSatish Balay   Input Parameters:
110119585528SSatish Balay +  mat     - the SeqSBAIJ matrix
110249b5e25fSSatish Balay -  indices - the column indices
110349b5e25fSSatish Balay 
110449b5e25fSSatish Balay   Level: advanced
110549b5e25fSSatish Balay 
110649b5e25fSSatish Balay   Notes:
110749b5e25fSSatish Balay     This can be called if you have precomputed the nonzero structure of the
110849b5e25fSSatish Balay   matrix and want to provide it to the matrix object to improve the performance
110949b5e25fSSatish Balay   of the MatSetValues() operation.
111049b5e25fSSatish Balay 
111149b5e25fSSatish Balay     You MUST have set the correct numbers of nonzeros per row in the call to
111219585528SSatish Balay   MatCreateSeqSBAIJ().
111349b5e25fSSatish Balay 
1114ab9f2c04SSatish Balay     MUST be called before any calls to MatSetValues()
111549b5e25fSSatish Balay 
1116ab9f2c04SSatish Balay .seealso: MatCreateSeqSBAIJ
111749b5e25fSSatish Balay @*/
111849b5e25fSSatish Balay int MatSeqSBAIJSetColumnIndices(Mat mat,int *indices)
111949b5e25fSSatish Balay {
112049b5e25fSSatish Balay   int ierr,(*f)(Mat,int *);
112149b5e25fSSatish Balay 
112249b5e25fSSatish Balay   PetscFunctionBegin;
112349b5e25fSSatish Balay   PetscValidHeaderSpecific(mat,MAT_COOKIE);
11244afc71dfSHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqSBAIJSetColumnIndices_C",(void **)&f);CHKERRQ(ierr);
112549b5e25fSSatish Balay   if (f) {
112649b5e25fSSatish Balay     ierr = (*f)(mat,indices);CHKERRQ(ierr);
112749b5e25fSSatish Balay   } else {
112849b5e25fSSatish Balay     SETERRQ(1,1,"Wrong type of matrix to set column indices");
112949b5e25fSSatish Balay   }
113049b5e25fSSatish Balay   PetscFunctionReturn(0);
113149b5e25fSSatish Balay }
113249b5e25fSSatish Balay 
113349b5e25fSSatish Balay /* -------------------------------------------------------------------*/
113449b5e25fSSatish Balay static struct _MatOps MatOps_Values = {MatSetValues_SeqSBAIJ,
113549b5e25fSSatish Balay        MatGetRow_SeqSBAIJ,
113649b5e25fSSatish Balay        MatRestoreRow_SeqSBAIJ,
113749b5e25fSSatish Balay        MatMult_SeqSBAIJ_N,
113849b5e25fSSatish Balay        MatMultAdd_SeqSBAIJ_N,
113949b5e25fSSatish Balay        MatMultTranspose_SeqSBAIJ,
114049b5e25fSSatish Balay        MatMultTransposeAdd_SeqSBAIJ,
114149b5e25fSSatish Balay        MatSolve_SeqSBAIJ_N,
114249b5e25fSSatish Balay        0,
114349b5e25fSSatish Balay        0,
114449b5e25fSSatish Balay        0,
114549b5e25fSSatish Balay        0,
11465f4ef2deSHong Zhang        MatCholeskyFactor_SeqSBAIJ,
114749b5e25fSSatish Balay        0,
114849b5e25fSSatish Balay        MatTranspose_SeqSBAIJ,
114949b5e25fSSatish Balay        MatGetInfo_SeqSBAIJ,
115049b5e25fSSatish Balay        MatEqual_SeqSBAIJ,
115149b5e25fSSatish Balay        MatGetDiagonal_SeqSBAIJ,
115249b5e25fSSatish Balay        MatDiagonalScale_SeqSBAIJ,
115349b5e25fSSatish Balay        MatNorm_SeqSBAIJ,
115449b5e25fSSatish Balay        0,
115549b5e25fSSatish Balay        MatAssemblyEnd_SeqSBAIJ,
115649b5e25fSSatish Balay        0,
115749b5e25fSSatish Balay        MatSetOption_SeqSBAIJ,
115849b5e25fSSatish Balay        MatZeroEntries_SeqSBAIJ,
115949b5e25fSSatish Balay        MatZeroRows_SeqSBAIJ,
116049b5e25fSSatish Balay        0,
116149b5e25fSSatish Balay        0,
11625f4ef2deSHong Zhang        MatCholeskyFactorSymbolic_SeqSBAIJ,
11635f4ef2deSHong Zhang        MatCholeskyFactorNumeric_SeqSBAIJ_N,
116449b5e25fSSatish Balay        MatGetSize_SeqSBAIJ,
11654ccecd49SHong Zhang        MatGetSize_SeqSBAIJ,
116649b5e25fSSatish Balay        MatGetOwnershipRange_SeqSBAIJ,
116749b5e25fSSatish Balay        0,
11685f4ef2deSHong Zhang        MatIncompleteCholeskyFactorSymbolic_SeqSBAIJ,
116949b5e25fSSatish Balay        0,
117049b5e25fSSatish Balay        0,
117149b5e25fSSatish Balay        MatDuplicate_SeqSBAIJ,
117249b5e25fSSatish Balay        0,
117349b5e25fSSatish Balay        0,
117449b5e25fSSatish Balay        0,
1175950f1e5bSHong Zhang        0,
117649b5e25fSSatish Balay        0,
117749b5e25fSSatish Balay        MatGetSubMatrices_SeqSBAIJ,
117849b5e25fSSatish Balay        MatIncreaseOverlap_SeqSBAIJ,
117949b5e25fSSatish Balay        MatGetValues_SeqSBAIJ,
118049b5e25fSSatish Balay        0,
118149b5e25fSSatish Balay        MatPrintHelp_SeqSBAIJ,
118249b5e25fSSatish Balay        MatScale_SeqSBAIJ,
118349b5e25fSSatish Balay        0,
118449b5e25fSSatish Balay        0,
118549b5e25fSSatish Balay        0,
118649b5e25fSSatish Balay        MatGetBlockSize_SeqSBAIJ,
118749b5e25fSSatish Balay        MatGetRowIJ_SeqSBAIJ,
118849b5e25fSSatish Balay        MatRestoreRowIJ_SeqSBAIJ,
118949b5e25fSSatish Balay        0,
119049b5e25fSSatish Balay        0,
119149b5e25fSSatish Balay        0,
119249b5e25fSSatish Balay        0,
119349b5e25fSSatish Balay        0,
119449b5e25fSSatish Balay        0,
119549b5e25fSSatish Balay        MatSetValuesBlocked_SeqSBAIJ,
119649b5e25fSSatish Balay        MatGetSubMatrix_SeqSBAIJ,
119749b5e25fSSatish Balay        0,
119849b5e25fSSatish Balay        0,
119949b5e25fSSatish Balay        MatGetMaps_Petsc};
120049b5e25fSSatish Balay 
120149b5e25fSSatish Balay EXTERN_C_BEGIN
120249b5e25fSSatish Balay #undef __FUNC__
120349b5e25fSSatish Balay #define __FUNC__ "MatStoreValues_SeqSBAIJ"
120449b5e25fSSatish Balay int MatStoreValues_SeqSBAIJ(Mat mat)
120549b5e25fSSatish Balay {
12064afc71dfSHong Zhang   Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data;
120749b5e25fSSatish Balay   int         nz = aij->i[aij->m]*aij->bs*aij->bs2;
120849b5e25fSSatish Balay   int         ierr;
120949b5e25fSSatish Balay 
121049b5e25fSSatish Balay   PetscFunctionBegin;
121149b5e25fSSatish Balay   if (aij->nonew != 1) {
121249b5e25fSSatish Balay     SETERRQ(1,1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
121349b5e25fSSatish Balay   }
121449b5e25fSSatish Balay 
121549b5e25fSSatish Balay   /* allocate space for values if not already there */
121649b5e25fSSatish Balay   if (!aij->saved_values) {
121749b5e25fSSatish Balay     aij->saved_values = (Scalar*)PetscMalloc(nz*sizeof(Scalar));CHKPTRQ(aij->saved_values);
121849b5e25fSSatish Balay   }
121949b5e25fSSatish Balay 
122049b5e25fSSatish Balay   /* copy values over */
122149b5e25fSSatish Balay   ierr = PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(Scalar));CHKERRQ(ierr);
122249b5e25fSSatish Balay   PetscFunctionReturn(0);
122349b5e25fSSatish Balay }
122449b5e25fSSatish Balay EXTERN_C_END
122549b5e25fSSatish Balay 
122649b5e25fSSatish Balay EXTERN_C_BEGIN
122749b5e25fSSatish Balay #undef __FUNC__
122849b5e25fSSatish Balay #define __FUNC__ "MatRetrieveValues_SeqSBAIJ"
122949b5e25fSSatish Balay int MatRetrieveValues_SeqSBAIJ(Mat mat)
123049b5e25fSSatish Balay {
12314afc71dfSHong Zhang   Mat_SeqSBAIJ *aij = (Mat_SeqSBAIJ *)mat->data;
123249b5e25fSSatish Balay   int         nz = aij->i[aij->m]*aij->bs*aij->bs2,ierr;
123349b5e25fSSatish Balay 
123449b5e25fSSatish Balay   PetscFunctionBegin;
123549b5e25fSSatish Balay   if (aij->nonew != 1) {
123649b5e25fSSatish Balay     SETERRQ(1,1,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first");
123749b5e25fSSatish Balay   }
123849b5e25fSSatish Balay   if (!aij->saved_values) {
123949b5e25fSSatish Balay     SETERRQ(1,1,"Must call MatStoreValues(A);first");
124049b5e25fSSatish Balay   }
124149b5e25fSSatish Balay 
124249b5e25fSSatish Balay   /* copy values over */
124349b5e25fSSatish Balay   ierr = PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(Scalar));CHKERRQ(ierr);
124449b5e25fSSatish Balay   PetscFunctionReturn(0);
124549b5e25fSSatish Balay }
124649b5e25fSSatish Balay EXTERN_C_END
124749b5e25fSSatish Balay 
124849b5e25fSSatish Balay #undef __FUNC__
124949b5e25fSSatish Balay #define __FUNC__ "MatCreateSeqSBAIJ"
125049b5e25fSSatish Balay /*@C
125149b5e25fSSatish Balay    MatCreateSeqSBAIJ - Creates a sparse symmetric matrix in block AIJ (block
125249b5e25fSSatish Balay    compressed row) format.  For good matrix assembly performance the
125349b5e25fSSatish Balay    user should preallocate the matrix storage by setting the parameter nz
125449b5e25fSSatish Balay    (or the array nnz).  By setting these parameters accurately, performance
125549b5e25fSSatish Balay    during matrix assembly can be increased by more than a factor of 50.
125649b5e25fSSatish Balay 
125749b5e25fSSatish Balay    Collective on MPI_Comm
125849b5e25fSSatish Balay 
125949b5e25fSSatish Balay    Input Parameters:
126049b5e25fSSatish Balay +  comm - MPI communicator, set to PETSC_COMM_SELF
126149b5e25fSSatish Balay .  bs - size of block
126249b5e25fSSatish Balay .  m - number of rows, or number of columns
126349b5e25fSSatish Balay .  nz - number of block nonzeros per block row (same for all rows)
126449b5e25fSSatish Balay -  nnz - array containing the number of block nonzeros in the various block rows
126549b5e25fSSatish Balay          (possibly different for each block row) or PETSC_NULL
126649b5e25fSSatish Balay 
126749b5e25fSSatish Balay    Output Parameter:
126849b5e25fSSatish Balay .  A - the symmetric matrix
126949b5e25fSSatish Balay 
127049b5e25fSSatish Balay    Options Database Keys:
127149b5e25fSSatish Balay .   -mat_no_unroll - uses code that does not unroll the loops in the
127249b5e25fSSatish Balay                      block calculations (much slower)
127349b5e25fSSatish Balay .    -mat_block_size - size of the blocks to use
127449b5e25fSSatish Balay 
127549b5e25fSSatish Balay    Level: intermediate
127649b5e25fSSatish Balay 
127749b5e25fSSatish Balay    Notes:
127849b5e25fSSatish Balay    The block AIJ format is fully compatible with standard Fortran 77
127949b5e25fSSatish Balay    storage.  That is, the stored row and column indices can begin at
128049b5e25fSSatish Balay    either one (as in Fortran) or zero.  See the users' manual for details.
128149b5e25fSSatish Balay 
128249b5e25fSSatish Balay    Specify the preallocated storage with either nz or nnz (not both).
128349b5e25fSSatish Balay    Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory
128449b5e25fSSatish Balay    allocation.  For additional details, see the users manual chapter on
128549b5e25fSSatish Balay    matrices.
128649b5e25fSSatish Balay 
128749b5e25fSSatish Balay .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatCreateMPIBAIJ()
128849b5e25fSSatish Balay @*/
128949b5e25fSSatish Balay int MatCreateSeqSBAIJ(MPI_Comm comm,int bs,int m,int n,int nz,int *nnz,Mat *A)
129049b5e25fSSatish Balay {
129149b5e25fSSatish Balay   Mat         B;
129249b5e25fSSatish Balay   Mat_SeqSBAIJ *b;
12934afc71dfSHong Zhang   int         i,len,ierr,mbs,bs2,size;
129449b5e25fSSatish Balay   PetscTruth  flg;
12954afc71dfSHong Zhang   int         s_nz;
129649b5e25fSSatish Balay 
129749b5e25fSSatish Balay   PetscFunctionBegin;
129849b5e25fSSatish Balay   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
129949b5e25fSSatish Balay   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"Comm must be of size 1");
130049b5e25fSSatish Balay 
130149b5e25fSSatish Balay   ierr = OptionsGetInt(PETSC_NULL,"-mat_block_size",&bs,PETSC_NULL);CHKERRQ(ierr);
130249b5e25fSSatish Balay   mbs  = m/bs;
130349b5e25fSSatish Balay   bs2  = bs*bs;
130449b5e25fSSatish Balay 
130549b5e25fSSatish Balay   if (mbs*bs!=m) {
130649b5e25fSSatish Balay     SETERRQ(PETSC_ERR_ARG_SIZ,0,"Number rows, cols must be divisible by blocksize");
130749b5e25fSSatish Balay   }
130849b5e25fSSatish Balay 
130949b5e25fSSatish Balay   if (nz < -2) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,0,"nz cannot be less than -2: value %d",nz);
131049b5e25fSSatish Balay   if (nnz) {
131149b5e25fSSatish Balay     for (i=0; i<mbs; i++) {
131249b5e25fSSatish 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]);
131349b5e25fSSatish Balay     }
131449b5e25fSSatish Balay   }
131549b5e25fSSatish Balay 
131649b5e25fSSatish Balay   *A = 0;
131749b5e25fSSatish Balay   PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQSBAIJ,"Mat",comm,MatDestroy,MatView);
131849b5e25fSSatish Balay   PLogObjectCreate(B);
131949b5e25fSSatish Balay   B->data = (void*)(b = PetscNew(Mat_SeqSBAIJ));CHKPTRQ(b);
132049b5e25fSSatish Balay   ierr    = PetscMemzero(b,sizeof(Mat_SeqSBAIJ));CHKERRQ(ierr);
132149b5e25fSSatish Balay   ierr    = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr);
132249b5e25fSSatish Balay   ierr    = OptionsHasName(PETSC_NULL,"-mat_no_unroll",&flg);CHKERRQ(ierr);
132349b5e25fSSatish Balay   if (!flg) {
132449b5e25fSSatish Balay     switch (bs) {
132549b5e25fSSatish Balay     case 1:
13264ccecd49SHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1;
132749b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_1;
132849b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_1;
132949b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_1;
133049b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_1;
133149b5e25fSSatish Balay       break;
133249b5e25fSSatish Balay     case 2:
13334ccecd49SHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2;
133449b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_2;
133549b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_2;
133649b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_2;
133749b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_2;
133849b5e25fSSatish Balay       break;
133949b5e25fSSatish Balay     case 3:
13405f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3;
134149b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_3;
134249b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_3;
134349b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_3;
134449b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_3;
134549b5e25fSSatish Balay       break;
134649b5e25fSSatish Balay     case 4:
13475f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4;
134849b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_4;
134949b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_4;
135049b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_4;
135149b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_4;
135249b5e25fSSatish Balay       break;
135349b5e25fSSatish Balay     case 5:
13545f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5;
135549b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_5;
135649b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_5;
135749b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_5;
135849b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_5;
135949b5e25fSSatish Balay       break;
136049b5e25fSSatish Balay     case 6:
13615f4ef2deSHong Zhang       B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6;
136249b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_6;
136349b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_6;
136449b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_6;
136549b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_6;
136649b5e25fSSatish Balay       break;
136749b5e25fSSatish Balay     case 7:
136849b5e25fSSatish Balay       B->ops->mult            = MatMult_SeqSBAIJ_7;
136949b5e25fSSatish Balay       B->ops->solve           = MatSolve_SeqSBAIJ_7;
137049b5e25fSSatish Balay       B->ops->solvetranspose  = MatSolveTranspose_SeqSBAIJ_7;
137149b5e25fSSatish Balay       B->ops->multadd         = MatMultAdd_SeqSBAIJ_7;
137249b5e25fSSatish Balay       break;
137349b5e25fSSatish Balay     }
137449b5e25fSSatish Balay   }
137549b5e25fSSatish Balay   B->ops->destroy     = MatDestroy_SeqSBAIJ;
137649b5e25fSSatish Balay   B->ops->view        = MatView_SeqSBAIJ;
137749b5e25fSSatish Balay   B->factor           = 0;
137849b5e25fSSatish Balay   B->lupivotthreshold = 1.0;
137949b5e25fSSatish Balay   B->mapping          = 0;
138049b5e25fSSatish Balay   b->row              = 0;
138149b5e25fSSatish Balay   b->icol             = 0;
138249b5e25fSSatish Balay   b->reallocs         = 0;
138349b5e25fSSatish Balay   b->saved_values     = 0;
138449b5e25fSSatish Balay 
13855f4ef2deSHong Zhang   b->m       = m;
13865f4ef2deSHong Zhang   B->m = m; B->M = m;
13875f4ef2deSHong Zhang   b->n       = m;
138849b5e25fSSatish Balay   B->n = m; B->N = m;
138949b5e25fSSatish Balay 
139049b5e25fSSatish Balay   ierr = MapCreateMPI(comm,m,m,&B->rmap);CHKERRQ(ierr);
139149b5e25fSSatish Balay   ierr = MapCreateMPI(comm,m,m,&B->cmap);CHKERRQ(ierr);
139249b5e25fSSatish Balay 
139349b5e25fSSatish Balay   b->mbs     = mbs;
13944afc71dfSHong Zhang   b->nbs     = mbs;
139549b5e25fSSatish Balay   b->imax    = (int*)PetscMalloc((mbs+1)*sizeof(int));CHKPTRQ(b->imax);
139649b5e25fSSatish Balay   if (!nnz) {
139749b5e25fSSatish Balay     if (nz == PETSC_DEFAULT) nz = 5;
139849b5e25fSSatish Balay     else if (nz <= 0)        nz = 1;
139949b5e25fSSatish Balay     for (i=0; i<mbs; i++) {
14008cef66ccSHong Zhang       b->imax[i] = nz;
140149b5e25fSSatish Balay     }
1402153ea458SHong Zhang     nz = nz*mbs; /* total nz */
140349b5e25fSSatish Balay   } else {
140449b5e25fSSatish Balay     nz = 0;
14058cef66ccSHong Zhang     for (i=0; i<mbs; i++) {b->imax[i] = nnz[i]; nz += nnz[i];}
140649b5e25fSSatish Balay   }
14078cef66ccSHong Zhang   /* s_nz=(nz+mbs)/2; */ /* total diagonal and superdiagonal nonzero blocks */
14088cef66ccSHong Zhang   s_nz = nz;
140949b5e25fSSatish Balay 
141049b5e25fSSatish Balay   /* allocate the matrix space */
141149b5e25fSSatish Balay   len   = s_nz*sizeof(int) + s_nz*bs2*sizeof(MatScalar) + (b->m+1)*sizeof(int);
141249b5e25fSSatish Balay   b->a  = (MatScalar*)PetscMalloc(len);CHKPTRQ(b->a);
141349b5e25fSSatish Balay   ierr = PetscMemzero(b->a,s_nz*bs2*sizeof(MatScalar));CHKERRQ(ierr);
141449b5e25fSSatish Balay   b->j  = (int*)(b->a + s_nz*bs2);
141549b5e25fSSatish Balay   ierr = PetscMemzero(b->j,s_nz*sizeof(int));CHKERRQ(ierr);
141649b5e25fSSatish Balay   b->i  = b->j + s_nz;
141749b5e25fSSatish Balay   b->singlemalloc = PETSC_TRUE;
141849b5e25fSSatish Balay 
141949b5e25fSSatish Balay   /* pointer to beginning of each row */
142049b5e25fSSatish Balay   b->i[0] = 0;
142149b5e25fSSatish Balay   for (i=1; i<mbs+1; i++) {
142249b5e25fSSatish Balay     b->i[i] = b->i[i-1] + b->imax[i-1];
142349b5e25fSSatish Balay   }
142449b5e25fSSatish Balay 
142549b5e25fSSatish Balay   /* b->ilen will count nonzeros in each block row so far. */
142649b5e25fSSatish Balay   b->ilen = (int*)PetscMalloc((mbs+1)*sizeof(int));
142749b5e25fSSatish Balay   PLogObjectMemory(B,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqSBAIJ));
142849b5e25fSSatish Balay   for (i=0; i<mbs; i++) { b->ilen[i] = 0;}
142949b5e25fSSatish Balay 
143049b5e25fSSatish Balay   b->bs               = bs;
143149b5e25fSSatish Balay   b->bs2              = bs2;
1432*1ce8b1feSHong Zhang   /* b->mbs              = mbs; redundant */
143349b5e25fSSatish Balay   b->s_nz               = 0;
143449b5e25fSSatish Balay   b->s_maxnz            = s_nz*bs2;
143549b5e25fSSatish Balay   b->sorted           = PETSC_FALSE;
143649b5e25fSSatish Balay   b->roworiented      = PETSC_TRUE;
143749b5e25fSSatish Balay   b->nonew            = 0;
143849b5e25fSSatish Balay   b->diag             = 0;
143949b5e25fSSatish Balay   b->solve_work       = 0;
144049b5e25fSSatish Balay   b->mult_work        = 0;
144149b5e25fSSatish Balay   b->spptr            = 0;
144249b5e25fSSatish Balay   B->info.nz_unneeded = (PetscReal)b->s_maxnz;
144349b5e25fSSatish Balay   b->keepzeroedrows   = PETSC_FALSE;
1444153ea458SHong Zhang 
144516cdd363SHong Zhang   b->inew             = 0;
144616cdd363SHong Zhang   b->jnew             = 0;
144716cdd363SHong Zhang   b->anew             = 0;
144816cdd363SHong Zhang   b->a2anew           = 0;
1449153ea458SHong Zhang 
145049b5e25fSSatish Balay   *A = B;
145149b5e25fSSatish Balay   ierr = OptionsHasName(PETSC_NULL,"-help",&flg);CHKERRQ(ierr);
145249b5e25fSSatish Balay   if (flg) {ierr = MatPrintHelp(B);CHKERRQ(ierr); }
145349b5e25fSSatish Balay 
145449b5e25fSSatish Balay   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C",
1455aec82049SSatish Balay                                      "MatStoreValues_SeqSBAIJ",
1456aec82049SSatish Balay                                      (void*)MatStoreValues_SeqSBAIJ);CHKERRQ(ierr);
145749b5e25fSSatish Balay   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C",
1458aec82049SSatish Balay                                      "MatRetrieveValues_SeqSBAIJ",
1459aec82049SSatish Balay                                      (void*)MatRetrieveValues_SeqSBAIJ);CHKERRQ(ierr);
1460aec82049SSatish Balay   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqSBAIJSetColumnIndices_C",
1461aec82049SSatish Balay                                      "MatSeqSBAIJSetColumnIndices_SeqSBAIJ",
1462*1ce8b1feSHong Zhang                                      (void*)MatSeqSBAIJSetColumnIndices_SeqSBAIJ);CHKERRQ(ierr);
1463*1ce8b1feSHong Zhang   /* printf("In MatCreateSeqSBAIJ, s_nz = %d, bs2=%d, m=%d, mbs=%d \n", s_nz,bs2,m,mbs); */
1464153ea458SHong Zhang    /*
146549b5e25fSSatish Balay    for (i=0; i<mbs; i++){
146649b5e25fSSatish Balay      printf("imax[%d]= %d, ilen[%d]= %d\n", i,b->imax[i], i,b->ilen[i]);
146749b5e25fSSatish Balay    }       */
146849b5e25fSSatish Balay 
146949b5e25fSSatish Balay   PetscFunctionReturn(0);
147049b5e25fSSatish Balay }
147149b5e25fSSatish Balay 
147249b5e25fSSatish Balay #undef __FUNC__
147349b5e25fSSatish Balay #define __FUNC__ "MatDuplicate_SeqSBAIJ"
147449b5e25fSSatish Balay int MatDuplicate_SeqSBAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B)
147549b5e25fSSatish Balay {
147649b5e25fSSatish Balay   Mat         C;
147749b5e25fSSatish Balay   Mat_SeqSBAIJ *c,*a = (Mat_SeqSBAIJ*)A->data;
147849b5e25fSSatish Balay   int         i,len,mbs = a->mbs,nz = a->s_nz,bs2 =a->bs2,ierr;
147949b5e25fSSatish Balay 
148049b5e25fSSatish Balay   PetscFunctionBegin;
148149b5e25fSSatish Balay   if (a->i[mbs] != nz) SETERRQ(PETSC_ERR_PLIB,0,"Corrupt matrix");
148249b5e25fSSatish Balay 
148349b5e25fSSatish Balay   *B = 0;
148449b5e25fSSatish Balay   PetscHeaderCreate(C,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQSBAIJ,"Mat",A->comm,MatDestroy,MatView);
148549b5e25fSSatish Balay   PLogObjectCreate(C);
148649b5e25fSSatish Balay   C->data           = (void*)(c = PetscNew(Mat_SeqSBAIJ));CHKPTRQ(c);
148749b5e25fSSatish Balay   ierr              = PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
148849b5e25fSSatish Balay   C->ops->destroy   = MatDestroy_SeqSBAIJ;
148949b5e25fSSatish Balay   C->ops->view      = MatView_SeqSBAIJ;
149049b5e25fSSatish Balay   C->factor         = A->factor;
149149b5e25fSSatish Balay   c->row            = 0;
149249b5e25fSSatish Balay   c->icol           = 0;
149349b5e25fSSatish Balay   c->saved_values   = 0;
149449b5e25fSSatish Balay   c->keepzeroedrows = a->keepzeroedrows;
149549b5e25fSSatish Balay   C->assembled      = PETSC_TRUE;
149649b5e25fSSatish Balay 
149749b5e25fSSatish Balay   c->m = C->m   = a->m;
149849b5e25fSSatish Balay   c->n = C->n   = a->n;
149949b5e25fSSatish Balay   C->M          = a->m;
150049b5e25fSSatish Balay   C->N          = a->n;
150149b5e25fSSatish Balay 
150249b5e25fSSatish Balay   c->bs         = a->bs;
150349b5e25fSSatish Balay   c->bs2        = a->bs2;
150449b5e25fSSatish Balay   c->mbs        = a->mbs;
150549b5e25fSSatish Balay   c->nbs        = a->nbs;
150649b5e25fSSatish Balay 
150749b5e25fSSatish Balay   c->imax       = (int*)PetscMalloc((mbs+1)*sizeof(int));CHKPTRQ(c->imax);
150849b5e25fSSatish Balay   c->ilen       = (int*)PetscMalloc((mbs+1)*sizeof(int));CHKPTRQ(c->ilen);
150949b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
151049b5e25fSSatish Balay     c->imax[i] = a->imax[i];
151149b5e25fSSatish Balay     c->ilen[i] = a->ilen[i];
151249b5e25fSSatish Balay   }
151349b5e25fSSatish Balay 
151449b5e25fSSatish Balay   /* allocate the matrix space */
151549b5e25fSSatish Balay   c->singlemalloc = PETSC_TRUE;
151649b5e25fSSatish Balay   len  = (mbs+1)*sizeof(int) + nz*(bs2*sizeof(MatScalar) + sizeof(int));
151749b5e25fSSatish Balay   c->a = (MatScalar*)PetscMalloc(len);CHKPTRQ(c->a);
151849b5e25fSSatish Balay   c->j = (int*)(c->a + nz*bs2);
151949b5e25fSSatish Balay   c->i = c->j + nz;
152049b5e25fSSatish Balay   ierr = PetscMemcpy(c->i,a->i,(mbs+1)*sizeof(int));CHKERRQ(ierr);
152149b5e25fSSatish Balay   if (mbs > 0) {
152249b5e25fSSatish Balay     ierr = PetscMemcpy(c->j,a->j,nz*sizeof(int));CHKERRQ(ierr);
152349b5e25fSSatish Balay     if (cpvalues == MAT_COPY_VALUES) {
152449b5e25fSSatish Balay       ierr = PetscMemcpy(c->a,a->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr);
152549b5e25fSSatish Balay     } else {
152649b5e25fSSatish Balay       ierr = PetscMemzero(c->a,bs2*nz*sizeof(MatScalar));CHKERRQ(ierr);
152749b5e25fSSatish Balay     }
152849b5e25fSSatish Balay   }
152949b5e25fSSatish Balay 
153049b5e25fSSatish Balay   PLogObjectMemory(C,len+2*(mbs+1)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_SeqSBAIJ));
153149b5e25fSSatish Balay   c->sorted      = a->sorted;
153249b5e25fSSatish Balay   c->roworiented = a->roworiented;
153349b5e25fSSatish Balay   c->nonew       = a->nonew;
153449b5e25fSSatish Balay 
153549b5e25fSSatish Balay   if (a->diag) {
153649b5e25fSSatish Balay     c->diag = (int*)PetscMalloc((mbs+1)*sizeof(int));CHKPTRQ(c->diag);
153749b5e25fSSatish Balay     PLogObjectMemory(C,(mbs+1)*sizeof(int));
153849b5e25fSSatish Balay     for (i=0; i<mbs; i++) {
153949b5e25fSSatish Balay       c->diag[i] = a->diag[i];
154049b5e25fSSatish Balay     }
154149b5e25fSSatish Balay   } else c->diag        = 0;
154249b5e25fSSatish Balay   c->s_nz               = a->s_nz;
154349b5e25fSSatish Balay   c->s_maxnz            = a->s_maxnz;
154449b5e25fSSatish Balay   c->solve_work         = 0;
154549b5e25fSSatish Balay   c->spptr              = 0;      /* Dangerous -I'm throwing away a->spptr */
154649b5e25fSSatish Balay   c->mult_work          = 0;
154749b5e25fSSatish Balay   *B = C;
154849b5e25fSSatish Balay   ierr = FListDuplicate(A->qlist,&C->qlist);CHKERRQ(ierr);
154949b5e25fSSatish Balay   PetscFunctionReturn(0);
155049b5e25fSSatish Balay }
155149b5e25fSSatish Balay 
155249b5e25fSSatish Balay #undef __FUNC__
155349b5e25fSSatish Balay #define __FUNC__ "MatLoad_SeqSBAIJ"
155449b5e25fSSatish Balay int MatLoad_SeqSBAIJ(Viewer viewer,MatType type,Mat *A)
155549b5e25fSSatish Balay {
155649b5e25fSSatish Balay   Mat_SeqSBAIJ  *a;
155749b5e25fSSatish Balay   Mat          B;
155849b5e25fSSatish Balay   int          i,nz,ierr,fd,header[4],size,*rowlengths=0,M,N,bs=1;
1559574b2666SHong Zhang   int          *mask,mbs,*jj,j,rowcount,nzcount,k,*browlengths,*s_browlengths,maskcount;
156049b5e25fSSatish Balay   int          kmax,jcount,block,idx,point,nzcountb,extra_rows;
156149b5e25fSSatish Balay   int          *masked,nmask,tmp,bs2,ishift;
156249b5e25fSSatish Balay   Scalar       *aa;
156349b5e25fSSatish Balay   MPI_Comm     comm = ((PetscObject)viewer)->comm;
156449b5e25fSSatish Balay 
156549b5e25fSSatish Balay   PetscFunctionBegin;
156649b5e25fSSatish Balay   ierr = OptionsGetInt(PETSC_NULL,"-matload_block_size",&bs,PETSC_NULL);CHKERRQ(ierr);
156749b5e25fSSatish Balay   bs2  = bs*bs;
156849b5e25fSSatish Balay 
156949b5e25fSSatish Balay   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
157049b5e25fSSatish Balay   if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"view must have one processor");
157149b5e25fSSatish Balay   ierr = ViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr);
157249b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr);
157349b5e25fSSatish Balay   if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"not Mat object");
157449b5e25fSSatish Balay   M = header[1]; N = header[2]; nz = header[3];
157549b5e25fSSatish Balay 
157649b5e25fSSatish Balay   if (header[3] < 0) {
157749b5e25fSSatish Balay     SETERRQ(PETSC_ERR_FILE_UNEXPECTED,1,"Matrix stored in special format, cannot load as SeqSBAIJ");
157849b5e25fSSatish Balay   }
157949b5e25fSSatish Balay 
158049b5e25fSSatish Balay   if (M != N) SETERRQ(PETSC_ERR_SUP,0,"Can only do square matrices");
158149b5e25fSSatish Balay 
158249b5e25fSSatish Balay   /*
158349b5e25fSSatish Balay      This code adds extra rows to make sure the number of rows is
158449b5e25fSSatish Balay     divisible by the blocksize
158549b5e25fSSatish Balay   */
158649b5e25fSSatish Balay   mbs        = M/bs;
158749b5e25fSSatish Balay   extra_rows = bs - M + bs*(mbs);
158849b5e25fSSatish Balay   if (extra_rows == bs) extra_rows = 0;
158949b5e25fSSatish Balay   else                  mbs++;
159049b5e25fSSatish Balay   if (extra_rows) {
159149b5e25fSSatish Balay     PLogInfo(0,"MatLoad_SeqSBAIJ:Padding loaded matrix to match blocksize\n");
159249b5e25fSSatish Balay   }
159349b5e25fSSatish Balay 
159449b5e25fSSatish Balay   /* read in row lengths */
159549b5e25fSSatish Balay   rowlengths = (int*)PetscMalloc((M+extra_rows)*sizeof(int));CHKPTRQ(rowlengths);
159649b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr);
159749b5e25fSSatish Balay   for (i=0; i<extra_rows; i++) rowlengths[M+i] = 1;
159849b5e25fSSatish Balay 
159949b5e25fSSatish Balay   /* read in column indices */
160049b5e25fSSatish Balay   jj = (int*)PetscMalloc((nz+extra_rows)*sizeof(int));CHKPTRQ(jj);
160149b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,jj,nz,PETSC_INT);CHKERRQ(ierr);
160249b5e25fSSatish Balay   for (i=0; i<extra_rows; i++) jj[nz+i] = M+i;
160349b5e25fSSatish Balay 
160449b5e25fSSatish Balay   /* loop over row lengths determining block row lengths */
160549b5e25fSSatish Balay   browlengths = (int*)PetscMalloc(mbs*sizeof(int));CHKPTRQ(browlengths);
1606574b2666SHong Zhang   s_browlengths = (int*)PetscMalloc(mbs*sizeof(int));CHKPTRQ(s_browlengths);
1607574b2666SHong Zhang   ierr        = PetscMemzero(s_browlengths,mbs*sizeof(int));CHKERRQ(ierr);
160849b5e25fSSatish Balay   mask        = (int*)PetscMalloc(2*mbs*sizeof(int));CHKPTRQ(mask);
160949b5e25fSSatish Balay   ierr        = PetscMemzero(mask,mbs*sizeof(int));CHKERRQ(ierr);
161049b5e25fSSatish Balay   masked      = mask + mbs;
161149b5e25fSSatish Balay   rowcount    = 0; nzcount = 0;
161249b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
161349b5e25fSSatish Balay     nmask = 0;
161449b5e25fSSatish Balay     for (j=0; j<bs; j++) {
161549b5e25fSSatish Balay       kmax = rowlengths[rowcount];
161649b5e25fSSatish Balay       for (k=0; k<kmax; k++) {
16172d703238SHong Zhang         tmp = jj[nzcount++]/bs;   /* block col. index */
161803630b6eSHong Zhang         if (!mask[tmp] && tmp >= i) {masked[nmask++] = tmp; mask[tmp] = 1;}
161949b5e25fSSatish Balay       }
162049b5e25fSSatish Balay       rowcount++;
162149b5e25fSSatish Balay     }
1622574b2666SHong Zhang     s_browlengths[i] += nmask;
1623574b2666SHong Zhang     browlengths[i]    = 2*s_browlengths[i];
1624574b2666SHong Zhang 
162549b5e25fSSatish Balay     /* zero out the mask elements we set */
162649b5e25fSSatish Balay     for (j=0; j<nmask; j++) mask[masked[j]] = 0;
162749b5e25fSSatish Balay   }
162849b5e25fSSatish Balay 
162949b5e25fSSatish Balay   /* create our matrix */
163049b5e25fSSatish Balay   ierr = MatCreateSeqSBAIJ(comm,bs,M+extra_rows,N+extra_rows,0,browlengths,A);CHKERRQ(ierr);
163149b5e25fSSatish Balay   B = *A;
163249b5e25fSSatish Balay   a = (Mat_SeqSBAIJ*)B->data;
163349b5e25fSSatish Balay 
163449b5e25fSSatish Balay   /* set matrix "i" values */
163549b5e25fSSatish Balay   a->i[0] = 0;
163649b5e25fSSatish Balay   for (i=1; i<= mbs; i++) {
1637574b2666SHong Zhang     a->i[i]      = a->i[i-1] + s_browlengths[i-1];
1638574b2666SHong Zhang     a->ilen[i-1] = s_browlengths[i-1];
163949b5e25fSSatish Balay   }
164049b5e25fSSatish Balay   a->s_nz         = 0;
1641574b2666SHong Zhang   for (i=0; i<mbs; i++) a->s_nz += s_browlengths[i];
164249b5e25fSSatish Balay 
164349b5e25fSSatish Balay   /* read in nonzero values */
164449b5e25fSSatish Balay   aa = (Scalar*)PetscMalloc((nz+extra_rows)*sizeof(Scalar));CHKPTRQ(aa);
164549b5e25fSSatish Balay   ierr = PetscBinaryRead(fd,aa,nz,PETSC_SCALAR);CHKERRQ(ierr);
164649b5e25fSSatish Balay   for (i=0; i<extra_rows; i++) aa[nz+i] = 1.0;
164749b5e25fSSatish Balay 
164849b5e25fSSatish Balay   /* set "a" and "j" values into matrix */
164949b5e25fSSatish Balay   nzcount = 0; jcount = 0;
165049b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
165149b5e25fSSatish Balay     nzcountb = nzcount;
165249b5e25fSSatish Balay     nmask    = 0;
165349b5e25fSSatish Balay     for (j=0; j<bs; j++) {
165449b5e25fSSatish Balay       kmax = rowlengths[i*bs+j];
165549b5e25fSSatish Balay       for (k=0; k<kmax; k++) {
16562d703238SHong Zhang         tmp = jj[nzcount++]/bs; /* block col. index */
165703630b6eSHong Zhang         if (!mask[tmp] && tmp >= i) { masked[nmask++] = tmp; mask[tmp] = 1;}
16582d703238SHong Zhang       }
16592d703238SHong Zhang       rowcount++;
16602d703238SHong Zhang     }
16612d703238SHong Zhang     /* sort the masked values */
16622d703238SHong Zhang     ierr = PetscSortInt(nmask,masked);CHKERRQ(ierr);
16632d703238SHong Zhang 
16642d703238SHong Zhang     /* set "j" values into matrix */
16652d703238SHong Zhang     maskcount = 1;
16662d703238SHong Zhang     for (j=0; j<nmask; j++) {
166749b5e25fSSatish Balay       a->j[jcount++]  = masked[j];
166849b5e25fSSatish Balay       mask[masked[j]] = maskcount++;
166949b5e25fSSatish Balay     }
1670574b2666SHong Zhang 
167149b5e25fSSatish Balay     /* set "a" values into matrix */
167249b5e25fSSatish Balay     ishift = bs2*a->i[i];
167349b5e25fSSatish Balay     for (j=0; j<bs; j++) {
167449b5e25fSSatish Balay       kmax = rowlengths[i*bs+j];
167549b5e25fSSatish Balay       for (k=0; k<kmax; k++) {
1676574b2666SHong Zhang         tmp       = jj[nzcountb]/bs ; /* block col. index */
1677574b2666SHong Zhang         if (tmp >= i){
167849b5e25fSSatish Balay           block     = mask[tmp] - 1;
167949b5e25fSSatish Balay           point     = jj[nzcountb] - bs*tmp;
168049b5e25fSSatish Balay           idx       = ishift + bs2*block + j + bs*point;
1681574b2666SHong Zhang           a->a[idx] = aa[nzcountb];
1682574b2666SHong Zhang         }
1683574b2666SHong Zhang         nzcountb++;
168449b5e25fSSatish Balay       }
168549b5e25fSSatish Balay     }
168649b5e25fSSatish Balay     /* zero out the mask elements we set */
168749b5e25fSSatish Balay     for (j=0; j<nmask; j++) mask[masked[j]] = 0;
168849b5e25fSSatish Balay   }
1689574b2666SHong Zhang   if (jcount != a->s_nz) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"Bad binary matrix");
169049b5e25fSSatish Balay 
169149b5e25fSSatish Balay   ierr = PetscFree(rowlengths);CHKERRQ(ierr);
169249b5e25fSSatish Balay   ierr = PetscFree(browlengths);CHKERRQ(ierr);
1693574b2666SHong Zhang   ierr = PetscFree(s_browlengths);CHKERRQ(ierr);
169449b5e25fSSatish Balay   ierr = PetscFree(aa);CHKERRQ(ierr);
169549b5e25fSSatish Balay   ierr = PetscFree(jj);CHKERRQ(ierr);
169649b5e25fSSatish Balay   ierr = PetscFree(mask);CHKERRQ(ierr);
169749b5e25fSSatish Balay 
169849b5e25fSSatish Balay   B->assembled = PETSC_TRUE;
169949b5e25fSSatish Balay   ierr = MatView_Private(B);CHKERRQ(ierr);
170049b5e25fSSatish Balay   PetscFunctionReturn(0);
170149b5e25fSSatish Balay }
1702574b2666SHong Zhang 
1703574b2666SHong Zhang 
170449b5e25fSSatish Balay 
170549b5e25fSSatish Balay 
1706