xref: /petsc/src/mat/impls/sbaij/seq/sbaij2.c (revision 0b8dc8d2aad2e8b6c929528b6bafed2bccee6aac)
149b5e25fSSatish Balay 
249b5e25fSSatish Balay #include "src/mat/impls/baij/seq/baij.h"
349b5e25fSSatish Balay #include "src/inline/spops.h"
449b5e25fSSatish Balay #include "src/inline/ilu.h"
549b5e25fSSatish Balay #include "petscbt.h"
63a7fca6bSBarry Smith #include "src/mat/impls/sbaij/seq/sbaij.h"
749b5e25fSSatish Balay 
84a2ae208SSatish Balay #undef __FUNCT__
94a2ae208SSatish Balay #define __FUNCT__ "MatIncreaseOverlap_SeqSBAIJ"
1013f74950SBarry Smith PetscErrorCode MatIncreaseOverlap_SeqSBAIJ(Mat A,PetscInt is_max,IS is[],PetscInt ov)
1149b5e25fSSatish Balay {
125eee224dSHong Zhang   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
136849ba73SBarry Smith   PetscErrorCode ierr;
14521d7252SBarry Smith   PetscInt       brow,i,j,k,l,mbs,n,*idx,*nidx,isz,bcol,bcol_max,start,end,*ai,*aj,bs,*nidx2;
15521d7252SBarry Smith   PetscBT        table,table0;
16d94109b8SHong Zhang 
17d94109b8SHong Zhang   PetscFunctionBegin;
18b3bf805bSHong Zhang   if (ov < 0)  SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap specified");
19d94109b8SHong Zhang   mbs = a->mbs;
20d94109b8SHong Zhang   ai  = a->i;
21d94109b8SHong Zhang   aj  = a->j;
22521d7252SBarry Smith   bs  = A->bs;
23d94109b8SHong Zhang   ierr = PetscBTCreate(mbs,table);CHKERRQ(ierr);
2413f74950SBarry Smith   ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&nidx);CHKERRQ(ierr);
2513f74950SBarry Smith   ierr = PetscMalloc((A->m+1)*sizeof(PetscInt),&nidx2);CHKERRQ(ierr);
26d94109b8SHong Zhang   ierr = PetscBTCreate(mbs,table0);CHKERRQ(ierr);
27d94109b8SHong Zhang 
28d94109b8SHong Zhang   for (i=0; i<is_max; i++) { /* for each is */
29d94109b8SHong Zhang     isz  = 0;
30d94109b8SHong Zhang     ierr = PetscBTMemzero(mbs,table);CHKERRQ(ierr);
31d94109b8SHong Zhang 
32d94109b8SHong Zhang     /* Extract the indices, assume there can be duplicate entries */
33d94109b8SHong Zhang     ierr = ISGetIndices(is[i],&idx);CHKERRQ(ierr);
34d94109b8SHong Zhang     ierr = ISGetLocalSize(is[i],&n);CHKERRQ(ierr);
35d94109b8SHong Zhang 
36d94109b8SHong Zhang     /* Enter these into the temp arrays i.e mark table[brow], enter brow into new index */
37dbe03f88SHong Zhang     bcol_max = 0;
38d94109b8SHong Zhang     for (j=0; j<n ; ++j){
39d94109b8SHong Zhang       brow = idx[j]/bs; /* convert the indices into block indices */
40d94109b8SHong Zhang       if (brow >= mbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"index greater than mat-dim");
41dbe03f88SHong Zhang       if(!PetscBTLookupSet(table,brow)) {
42dbe03f88SHong Zhang         nidx[isz++] = brow;
43dbe03f88SHong Zhang         if (bcol_max < brow) bcol_max = brow;
44dbe03f88SHong Zhang       }
45d94109b8SHong Zhang     }
46d94109b8SHong Zhang     ierr = ISRestoreIndices(is[i],&idx);CHKERRQ(ierr);
47d94109b8SHong Zhang     ierr = ISDestroy(is[i]);CHKERRQ(ierr);
48d94109b8SHong Zhang 
49d94109b8SHong Zhang     k = 0;
50d94109b8SHong Zhang     for (j=0; j<ov; j++){ /* for each overlap */
51d94109b8SHong Zhang       /* set table0 for lookup - only mark entries that are added onto nidx in (j-1)-th overlap */
52d94109b8SHong Zhang       ierr = PetscBTMemzero(mbs,table0);CHKERRQ(ierr);
53efee365bSSatish Balay       for (l=k; l<isz; l++) { ierr = PetscBTSet(table0,nidx[l]);CHKERRQ(ierr); }
54d94109b8SHong Zhang 
55d94109b8SHong Zhang       n = isz;  /* length of the updated is[i] */
56d94109b8SHong Zhang       for (brow=0; brow<mbs; brow++){
57d94109b8SHong Zhang         start = ai[brow]; end   = ai[brow+1];
58d94109b8SHong Zhang         if (PetscBTLookup(table0,brow)){ /* brow is on nidx - row search: collect all bcol in this brow */
59d94109b8SHong Zhang           for (l = start; l<end ; l++){
60d94109b8SHong Zhang             bcol = aj[l];
61d94109b8SHong Zhang             if (!PetscBTLookupSet(table,bcol)) {nidx[isz++] = bcol;}
62d94109b8SHong Zhang           }
63d94109b8SHong Zhang           k++;
64d94109b8SHong Zhang           if (k >= n) break; /* for (brow=0; brow<mbs; brow++) */
65d94109b8SHong Zhang         } else { /* brow is not on nidx - col serach: add brow onto nidx if there is a bcol in nidx */
66d94109b8SHong Zhang           for (l = start; l<end ; l++){
67d94109b8SHong Zhang             bcol = aj[l];
68dbe03f88SHong Zhang             if (bcol > bcol_max) break;
69d94109b8SHong Zhang             if (PetscBTLookup(table0,bcol)){
70d94109b8SHong Zhang               if (!PetscBTLookupSet(table,brow)) {nidx[isz++] = brow;}
71d94109b8SHong Zhang               break; /* for l = start; l<end ; l++) */
72d94109b8SHong Zhang             }
73d94109b8SHong Zhang           }
74d94109b8SHong Zhang         }
75d94109b8SHong Zhang       }
76d94109b8SHong Zhang     } /* for each overlap */
77d94109b8SHong Zhang 
78d94109b8SHong Zhang     /* expand the Index Set */
79d94109b8SHong Zhang     for (j=0; j<isz; j++) {
80d94109b8SHong Zhang       for (k=0; k<bs; k++)
81d94109b8SHong Zhang         nidx2[j*bs+k] = nidx[j]*bs+k;
82d94109b8SHong Zhang     }
83d94109b8SHong Zhang     ierr = ISCreateGeneral(PETSC_COMM_SELF,isz*bs,nidx2,is+i);CHKERRQ(ierr);
84d94109b8SHong Zhang   }
85d94109b8SHong Zhang   ierr = PetscBTDestroy(table);CHKERRQ(ierr);
86d94109b8SHong Zhang   ierr = PetscFree(nidx);CHKERRQ(ierr);
87d94109b8SHong Zhang   ierr = PetscFree(nidx2);CHKERRQ(ierr);
88d94109b8SHong Zhang   ierr = PetscBTDestroy(table0);CHKERRQ(ierr);
895eee224dSHong Zhang   PetscFunctionReturn(0);
9049b5e25fSSatish Balay }
9149b5e25fSSatish Balay 
924a2ae208SSatish Balay #undef __FUNCT__
934a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqSBAIJ_Private"
9413f74950SBarry Smith PetscErrorCode MatGetSubMatrix_SeqSBAIJ_Private(Mat A,IS isrow,IS iscol,PetscInt cs,MatReuse scall,Mat *B)
9549b5e25fSSatish Balay {
9649b5e25fSSatish Balay   Mat_SeqSBAIJ  *a = (Mat_SeqSBAIJ*)A->data,*c;
976849ba73SBarry Smith   PetscErrorCode ierr;
9813f74950SBarry Smith   PetscInt       *smap,i,k,kstart,kend,oldcols = a->mbs,*lens;
9913f74950SBarry Smith   PetscInt       row,mat_i,*mat_j,tcol,*mat_ilen;
10013f74950SBarry Smith   PetscInt       *irow,nrows,*ssmap,bs=A->bs,bs2=a->bs2;
10113f74950SBarry Smith   PetscInt       *aj = a->j,*ai = a->i;
10249b5e25fSSatish Balay   MatScalar      *mat_a;
10349b5e25fSSatish Balay   Mat            C;
10449b5e25fSSatish Balay   PetscTruth     flag;
10549b5e25fSSatish Balay 
10649b5e25fSSatish Balay   PetscFunctionBegin;
107e005ede5SBarry Smith   if (isrow != iscol) SETERRQ(PETSC_ERR_ARG_INCOMP,"For symmetric format, iscol must equal isro");
10849b5e25fSSatish Balay   ierr = ISSorted(iscol,(PetscTruth*)&i);CHKERRQ(ierr);
109347d480fSBarry Smith   if (!i) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"IS is not sorted");
11049b5e25fSSatish Balay 
11149b5e25fSSatish Balay   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
11249b5e25fSSatish Balay   ierr = ISGetSize(isrow,&nrows);CHKERRQ(ierr);
11349b5e25fSSatish Balay 
11413f74950SBarry Smith   ierr  = PetscMalloc((1+oldcols)*sizeof(PetscInt),&smap);CHKERRQ(ierr);
11549b5e25fSSatish Balay   ssmap = smap;
11613f74950SBarry Smith   ierr  = PetscMalloc((1+nrows)*sizeof(PetscInt),&lens);CHKERRQ(ierr);
11713f74950SBarry Smith   ierr  = PetscMemzero(smap,oldcols*sizeof(PetscInt));CHKERRQ(ierr);
11849b5e25fSSatish Balay   for (i=0; i<nrows; i++) smap[irow[i]] = i+1; /* nrows = ncols */
11949b5e25fSSatish Balay   /* determine lens of each row */
12049b5e25fSSatish Balay   for (i=0; i<nrows; i++) {
12149b5e25fSSatish Balay     kstart  = ai[irow[i]];
12249b5e25fSSatish Balay     kend    = kstart + a->ilen[irow[i]];
12349b5e25fSSatish Balay     lens[i] = 0;
12449b5e25fSSatish Balay       for (k=kstart; k<kend; k++) {
12549b5e25fSSatish Balay         if (ssmap[aj[k]]) {
12649b5e25fSSatish Balay           lens[i]++;
12749b5e25fSSatish Balay         }
12849b5e25fSSatish Balay       }
12949b5e25fSSatish Balay     }
13049b5e25fSSatish Balay   /* Create and fill new matrix */
13149b5e25fSSatish Balay   if (scall == MAT_REUSE_MATRIX) {
13249b5e25fSSatish Balay     c = (Mat_SeqSBAIJ *)((*B)->data);
13349b5e25fSSatish Balay 
134521d7252SBarry Smith     if (c->mbs!=nrows || (*B)->bs!=bs) SETERRQ(PETSC_ERR_ARG_SIZ,"Submatrix wrong size");
13513f74950SBarry Smith     ierr = PetscMemcmp(c->ilen,lens,c->mbs *sizeof(PetscInt),&flag);CHKERRQ(ierr);
136abc0a331SBarry Smith     if (!flag) {
137347d480fSBarry Smith       SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong no of nonzeros");
13849b5e25fSSatish Balay     }
13913f74950SBarry Smith     ierr = PetscMemzero(c->ilen,c->mbs*sizeof(PetscInt));CHKERRQ(ierr);
14049b5e25fSSatish Balay     C = *B;
14149b5e25fSSatish Balay   } else {
142e2d9671bSKris Buschelman     ierr = MatCreate(A->comm,nrows*bs,nrows*bs,PETSC_DETERMINE,PETSC_DETERMINE,&C);CHKERRQ(ierr);
143e2d9671bSKris Buschelman     ierr = MatSetType(C,A->type_name);CHKERRQ(ierr);
144e2d9671bSKris Buschelman     ierr = MatSeqSBAIJSetPreallocation(C,bs,0,lens);CHKERRQ(ierr);
14549b5e25fSSatish Balay   }
14649b5e25fSSatish Balay   c = (Mat_SeqSBAIJ *)(C->data);
14749b5e25fSSatish Balay   for (i=0; i<nrows; i++) {
14849b5e25fSSatish Balay     row    = irow[i];
14949b5e25fSSatish Balay     kstart = ai[row];
15049b5e25fSSatish Balay     kend   = kstart + a->ilen[row];
15149b5e25fSSatish Balay     mat_i  = c->i[i];
15249b5e25fSSatish Balay     mat_j  = c->j + mat_i;
15349b5e25fSSatish Balay     mat_a  = c->a + mat_i*bs2;
15449b5e25fSSatish Balay     mat_ilen = c->ilen + i;
15549b5e25fSSatish Balay     for (k=kstart; k<kend; k++) {
15649b5e25fSSatish Balay       if ((tcol=ssmap[a->j[k]])) {
15749b5e25fSSatish Balay         *mat_j++ = tcol - 1;
15849b5e25fSSatish Balay         ierr     = PetscMemcpy(mat_a,a->a+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
15949b5e25fSSatish Balay         mat_a   += bs2;
16049b5e25fSSatish Balay         (*mat_ilen)++;
16149b5e25fSSatish Balay       }
16249b5e25fSSatish Balay     }
16349b5e25fSSatish Balay   }
16449b5e25fSSatish Balay 
16549b5e25fSSatish Balay   /* Free work space */
16649b5e25fSSatish Balay   ierr = PetscFree(smap);CHKERRQ(ierr);
16749b5e25fSSatish Balay   ierr = PetscFree(lens);CHKERRQ(ierr);
16849b5e25fSSatish Balay   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
16949b5e25fSSatish Balay   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
17049b5e25fSSatish Balay 
17149b5e25fSSatish Balay   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
17249b5e25fSSatish Balay   *B = C;
17349b5e25fSSatish Balay   PetscFunctionReturn(0);
17449b5e25fSSatish Balay }
17549b5e25fSSatish Balay 
1764a2ae208SSatish Balay #undef __FUNCT__
1774a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqSBAIJ"
17813f74950SBarry Smith PetscErrorCode MatGetSubMatrix_SeqSBAIJ(Mat A,IS isrow,IS iscol,PetscInt cs,MatReuse scall,Mat *B)
17949b5e25fSSatish Balay {
18049b5e25fSSatish Balay   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
18149b5e25fSSatish Balay   IS             is1;
1826849ba73SBarry Smith   PetscErrorCode ierr;
18313f74950SBarry Smith   PetscInt       *vary,*iary,*irow,nrows,i,bs=A->bs,count;
18449b5e25fSSatish Balay 
18549b5e25fSSatish Balay   PetscFunctionBegin;
186e005ede5SBarry Smith   if (isrow != iscol) SETERRQ(PETSC_ERR_ARG_INCOMP,"For symmetric format, iscol must equal isro");
18749b5e25fSSatish Balay 
18849b5e25fSSatish Balay   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
18949b5e25fSSatish Balay   ierr = ISGetSize(isrow,&nrows);CHKERRQ(ierr);
19049b5e25fSSatish Balay 
19149b5e25fSSatish Balay   /* Verify if the indices corespond to each element in a block
19249b5e25fSSatish Balay    and form the IS with compressed IS */
19313f74950SBarry Smith   ierr = PetscMalloc(2*(a->mbs+1)*sizeof(PetscInt),&vary);CHKERRQ(ierr);
19449b5e25fSSatish Balay   iary = vary + a->mbs;
19513f74950SBarry Smith   ierr = PetscMemzero(vary,(a->mbs)*sizeof(PetscInt));CHKERRQ(ierr);
19649b5e25fSSatish Balay   for (i=0; i<nrows; i++) vary[irow[i]/bs]++;
19749b5e25fSSatish Balay 
19849b5e25fSSatish Balay   count = 0;
19949b5e25fSSatish Balay   for (i=0; i<a->mbs; i++) {
200e005ede5SBarry Smith     if (vary[i]!=0 && vary[i]!=bs) SETERRQ(PETSC_ERR_ARG_INCOMP,"Index set does not match blocks");
20149b5e25fSSatish Balay     if (vary[i]==bs) iary[count++] = i;
20249b5e25fSSatish Balay   }
20349b5e25fSSatish Balay   ierr = ISCreateGeneral(PETSC_COMM_SELF,count,iary,&is1);CHKERRQ(ierr);
20449b5e25fSSatish Balay 
20549b5e25fSSatish Balay   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
20649b5e25fSSatish Balay   ierr = PetscFree(vary);CHKERRQ(ierr);
20749b5e25fSSatish Balay 
20849b5e25fSSatish Balay   ierr = MatGetSubMatrix_SeqSBAIJ_Private(A,is1,is1,cs,scall,B);CHKERRQ(ierr);
20949b5e25fSSatish Balay   ISDestroy(is1);
21049b5e25fSSatish Balay   PetscFunctionReturn(0);
21149b5e25fSSatish Balay }
21249b5e25fSSatish Balay 
2134a2ae208SSatish Balay #undef __FUNCT__
2144a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_SeqSBAIJ"
21513f74950SBarry Smith PetscErrorCode MatGetSubMatrices_SeqSBAIJ(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
21649b5e25fSSatish Balay {
2176849ba73SBarry Smith   PetscErrorCode ierr;
21813f74950SBarry Smith   PetscInt       i;
21949b5e25fSSatish Balay 
22049b5e25fSSatish Balay   PetscFunctionBegin;
22149b5e25fSSatish Balay   if (scall == MAT_INITIAL_MATRIX) {
22282502324SSatish Balay     ierr = PetscMalloc((n+1)*sizeof(Mat),B);CHKERRQ(ierr);
22349b5e25fSSatish Balay   }
22449b5e25fSSatish Balay 
22549b5e25fSSatish Balay   for (i=0; i<n; i++) {
22649b5e25fSSatish Balay     ierr = MatGetSubMatrix_SeqSBAIJ(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr);
22749b5e25fSSatish Balay   }
22849b5e25fSSatish Balay   PetscFunctionReturn(0);
22949b5e25fSSatish Balay }
23049b5e25fSSatish Balay 
23149b5e25fSSatish Balay /* -------------------------------------------------------*/
23249b5e25fSSatish Balay /* Should check that shapes of vectors and matrices match */
23349b5e25fSSatish Balay /* -------------------------------------------------------*/
234d9eff348SSatish Balay #include "petscblaslapack.h"
23549b5e25fSSatish Balay 
2364a2ae208SSatish Balay #undef __FUNCT__
2374a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqSBAIJ_1"
238dfbe8321SBarry Smith PetscErrorCode MatMult_SeqSBAIJ_1(Mat A,Vec xx,Vec zz)
23949b5e25fSSatish Balay {
24049b5e25fSSatish Balay   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
24187828ca2SBarry Smith   PetscScalar    *x,*z,*xb,x1,zero=0.0;
24249b5e25fSSatish Balay   MatScalar      *v;
2436849ba73SBarry Smith   PetscErrorCode ierr;
24413f74950SBarry Smith   PetscInt       mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,*ib,cval,j,jmin;
24549b5e25fSSatish Balay 
24649b5e25fSSatish Balay   PetscFunctionBegin;
24749b5e25fSSatish Balay   ierr = VecSet(&zero,zz);CHKERRQ(ierr);
2481ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2491ebc52fbSHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
25049b5e25fSSatish Balay 
25149b5e25fSSatish Balay   v  = a->a;
25249b5e25fSSatish Balay   xb = x;
25349b5e25fSSatish Balay 
25449b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
25549b5e25fSSatish Balay     n  = ai[1] - ai[0];  /* length of i_th row of A */
25649b5e25fSSatish Balay     x1 = xb[0];
25749b5e25fSSatish Balay     ib = aj + *ai;
258831a3094SHong Zhang     jmin = 0;
259831a3094SHong Zhang     if (*ib == i) {      /* (diag of A)*x */
260831a3094SHong Zhang       z[i] += *v++ * x[*ib++];
261831a3094SHong Zhang       jmin++;
262831a3094SHong Zhang     }
263831a3094SHong Zhang     for (j=jmin; j<n; j++) {
26449b5e25fSSatish Balay       cval    = *ib;
26549b5e25fSSatish Balay       z[cval] += *v * x1;      /* (strict lower triangular part of A)*x  */
26649b5e25fSSatish Balay       z[i] += *v++ * x[*ib++]; /* (strict upper triangular part of A)*x  */
26749b5e25fSSatish Balay     }
26849b5e25fSSatish Balay     xb++; ai++;
26949b5e25fSSatish Balay   }
27049b5e25fSSatish Balay 
2711ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2721ebc52fbSHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
273efee365bSSatish Balay   ierr = PetscLogFlops(2*(a->nz*2 - A->m) - A->m);CHKERRQ(ierr);  /* nz = (nz+m)/2 */
27449b5e25fSSatish Balay   PetscFunctionReturn(0);
27549b5e25fSSatish Balay }
27649b5e25fSSatish Balay 
2774a2ae208SSatish Balay #undef __FUNCT__
2784a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqSBAIJ_2"
279dfbe8321SBarry Smith PetscErrorCode MatMult_SeqSBAIJ_2(Mat A,Vec xx,Vec zz)
28049b5e25fSSatish Balay {
28149b5e25fSSatish Balay   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
28287828ca2SBarry Smith   PetscScalar    *x,*z,*xb,x1,x2,zero=0.0;
28349b5e25fSSatish Balay   MatScalar      *v;
2846849ba73SBarry Smith   PetscErrorCode ierr;
28513f74950SBarry Smith   PetscInt       mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,*ib,cval,j,jmin;
28649b5e25fSSatish Balay 
28749b5e25fSSatish Balay   PetscFunctionBegin;
28849b5e25fSSatish Balay   ierr = VecSet(&zero,zz);CHKERRQ(ierr);
2891ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2901ebc52fbSHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
29149b5e25fSSatish Balay 
29249b5e25fSSatish Balay   v     = a->a;
29349b5e25fSSatish Balay   xb = x;
29449b5e25fSSatish Balay 
29549b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
29649b5e25fSSatish Balay     n  = ai[1] - ai[0]; /* length of i_th block row of A */
29749b5e25fSSatish Balay     x1 = xb[0]; x2 = xb[1];
29849b5e25fSSatish Balay     ib = aj + *ai;
299831a3094SHong Zhang     jmin = 0;
3007fbae186SHong Zhang     if (*ib == i){     /* (diag of A)*x */
30149b5e25fSSatish Balay       z[2*i]   += v[0]*x1 + v[2]*x2;
30249b5e25fSSatish Balay       z[2*i+1] += v[2]*x1 + v[3]*x2;
303831a3094SHong Zhang       v += 4; jmin++;
3047fbae186SHong Zhang     }
305831a3094SHong Zhang     for (j=jmin; j<n; j++) {
30649b5e25fSSatish Balay       /* (strict lower triangular part of A)*x  */
30749b5e25fSSatish Balay       cval       = ib[j]*2;
30849b5e25fSSatish Balay       z[cval]     += v[0]*x1 + v[1]*x2;
30949b5e25fSSatish Balay       z[cval+1]   += v[2]*x1 + v[3]*x2;
31049b5e25fSSatish Balay       /* (strict upper triangular part of A)*x  */
31149b5e25fSSatish Balay       z[2*i]   += v[0]*x[cval] + v[2]*x[cval+1];
31249b5e25fSSatish Balay       z[2*i+1] += v[1]*x[cval] + v[3]*x[cval+1];
31349b5e25fSSatish Balay       v  += 4;
31449b5e25fSSatish Balay     }
31549b5e25fSSatish Balay     xb +=2; ai++;
31649b5e25fSSatish Balay   }
31749b5e25fSSatish Balay 
3181ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
3191ebc52fbSHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
320efee365bSSatish Balay   ierr = PetscLogFlops(8*(a->nz*2 - A->m) - A->m);CHKERRQ(ierr);
32149b5e25fSSatish Balay   PetscFunctionReturn(0);
32249b5e25fSSatish Balay }
32349b5e25fSSatish Balay 
3244a2ae208SSatish Balay #undef __FUNCT__
3254a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqSBAIJ_3"
326dfbe8321SBarry Smith PetscErrorCode MatMult_SeqSBAIJ_3(Mat A,Vec xx,Vec zz)
32749b5e25fSSatish Balay {
32849b5e25fSSatish Balay   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
32987828ca2SBarry Smith   PetscScalar    *x,*z,*xb,x1,x2,x3,zero=0.0;
33049b5e25fSSatish Balay   MatScalar      *v;
3316849ba73SBarry Smith   PetscErrorCode ierr;
33213f74950SBarry Smith   PetscInt       mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,*ib,cval,j,jmin;
33349b5e25fSSatish Balay 
33449b5e25fSSatish Balay   PetscFunctionBegin;
33549b5e25fSSatish Balay   ierr = VecSet(&zero,zz);CHKERRQ(ierr);
3361ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
3371ebc52fbSHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
33849b5e25fSSatish Balay 
33949b5e25fSSatish Balay   v    = a->a;
34049b5e25fSSatish Balay   xb   = x;
34149b5e25fSSatish Balay 
34249b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
34349b5e25fSSatish Balay     n  = ai[1] - ai[0]; /* length of i_th block row of A */
34449b5e25fSSatish Balay     x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
34549b5e25fSSatish Balay     ib = aj + *ai;
346831a3094SHong Zhang     jmin = 0;
3477fbae186SHong Zhang     if (*ib == i){     /* (diag of A)*x */
34849b5e25fSSatish Balay       z[3*i]   += v[0]*x1 + v[3]*x2 + v[6]*x3;
34949b5e25fSSatish Balay       z[3*i+1] += v[3]*x1 + v[4]*x2 + v[7]*x3;
35049b5e25fSSatish Balay       z[3*i+2] += v[6]*x1 + v[7]*x2 + v[8]*x3;
351831a3094SHong Zhang       v += 9; jmin++;
3527fbae186SHong Zhang     }
353831a3094SHong Zhang     for (j=jmin; j<n; j++) {
35449b5e25fSSatish Balay       /* (strict lower triangular part of A)*x  */
35549b5e25fSSatish Balay       cval       = ib[j]*3;
35649b5e25fSSatish Balay       z[cval]     += v[0]*x1 + v[1]*x2 + v[2]*x3;
35749b5e25fSSatish Balay       z[cval+1]   += v[3]*x1 + v[4]*x2 + v[5]*x3;
35849b5e25fSSatish Balay       z[cval+2]   += v[6]*x1 + v[7]*x2 + v[8]*x3;
35949b5e25fSSatish Balay       /* (strict upper triangular part of A)*x  */
36049b5e25fSSatish Balay       z[3*i]   += v[0]*x[cval] + v[3]*x[cval+1]+ v[6]*x[cval+2];
36149b5e25fSSatish Balay       z[3*i+1] += v[1]*x[cval] + v[4]*x[cval+1]+ v[7]*x[cval+2];
36249b5e25fSSatish Balay       z[3*i+2] += v[2]*x[cval] + v[5]*x[cval+1]+ v[8]*x[cval+2];
36349b5e25fSSatish Balay       v  += 9;
36449b5e25fSSatish Balay     }
36549b5e25fSSatish Balay     xb +=3; ai++;
36649b5e25fSSatish Balay   }
36749b5e25fSSatish Balay 
3681ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
3691ebc52fbSHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
370efee365bSSatish Balay   ierr = PetscLogFlops(18*(a->nz*2 - A->m) - A->m);CHKERRQ(ierr);
37149b5e25fSSatish Balay   PetscFunctionReturn(0);
37249b5e25fSSatish Balay }
37349b5e25fSSatish Balay 
3744a2ae208SSatish Balay #undef __FUNCT__
3754a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqSBAIJ_4"
376dfbe8321SBarry Smith PetscErrorCode MatMult_SeqSBAIJ_4(Mat A,Vec xx,Vec zz)
37749b5e25fSSatish Balay {
37849b5e25fSSatish Balay   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
37987828ca2SBarry Smith   PetscScalar    *x,*z,*xb,x1,x2,x3,x4,zero=0.0;
38049b5e25fSSatish Balay   MatScalar      *v;
3816849ba73SBarry Smith   PetscErrorCode ierr;
38213f74950SBarry Smith   PetscInt       mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,*ib,cval,j,jmin;
38349b5e25fSSatish Balay 
38449b5e25fSSatish Balay   PetscFunctionBegin;
38549b5e25fSSatish Balay   ierr = VecSet(&zero,zz);CHKERRQ(ierr);
3861ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
3871ebc52fbSHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
38849b5e25fSSatish Balay 
38949b5e25fSSatish Balay   v     = a->a;
39049b5e25fSSatish Balay   xb = x;
39149b5e25fSSatish Balay 
39249b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
39349b5e25fSSatish Balay     n  = ai[1] - ai[0]; /* length of i_th block row of A */
39449b5e25fSSatish Balay     x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
39549b5e25fSSatish Balay     ib = aj + *ai;
396831a3094SHong Zhang     jmin = 0;
3977fbae186SHong Zhang     if (*ib == i){     /* (diag of A)*x */
39849b5e25fSSatish Balay       z[4*i]   += v[0]*x1 + v[4]*x2 +  v[8]*x3 + v[12]*x4;
39949b5e25fSSatish Balay       z[4*i+1] += v[4]*x1 + v[5]*x2 +  v[9]*x3 + v[13]*x4;
40049b5e25fSSatish Balay       z[4*i+2] += v[8]*x1 + v[9]*x2 + v[10]*x3 + v[14]*x4;
40149b5e25fSSatish Balay       z[4*i+3] += v[12]*x1+ v[13]*x2+ v[14]*x3 + v[15]*x4;
402831a3094SHong Zhang       v += 16; jmin++;
4037fbae186SHong Zhang     }
404831a3094SHong Zhang     for (j=jmin; j<n; j++) {
40549b5e25fSSatish Balay       /* (strict lower triangular part of A)*x  */
40649b5e25fSSatish Balay       cval       = ib[j]*4;
40749b5e25fSSatish Balay       z[cval]     += v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4;
40849b5e25fSSatish Balay       z[cval+1]   += v[4]*x1 + v[5]*x2 + v[6]*x3 + v[7]*x4;
40949b5e25fSSatish Balay       z[cval+2]   += v[8]*x1 + v[9]*x2 + v[10]*x3 + v[11]*x4;
41049b5e25fSSatish Balay       z[cval+3]   += v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4;
41149b5e25fSSatish Balay       /* (strict upper triangular part of A)*x  */
41249b5e25fSSatish Balay       z[4*i]   += v[0]*x[cval] + v[4]*x[cval+1]+ v[8]*x[cval+2] + v[12]*x[cval+3];
41349b5e25fSSatish Balay       z[4*i+1] += v[1]*x[cval] + v[5]*x[cval+1]+ v[9]*x[cval+2] + v[13]*x[cval+3];
41449b5e25fSSatish Balay       z[4*i+2] += v[2]*x[cval] + v[6]*x[cval+1]+ v[10]*x[cval+2]+ v[14]*x[cval+3];
41549b5e25fSSatish Balay       z[4*i+3] += v[3]*x[cval] + v[7]*x[cval+1]+ v[11]*x[cval+2]+ v[15]*x[cval+3];
41649b5e25fSSatish Balay       v  += 16;
41749b5e25fSSatish Balay     }
41849b5e25fSSatish Balay     xb +=4; ai++;
41949b5e25fSSatish Balay   }
42049b5e25fSSatish Balay 
4211ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
4221ebc52fbSHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
423efee365bSSatish Balay   ierr = PetscLogFlops(32*(a->nz*2 - A->m) - A->m);CHKERRQ(ierr);
42449b5e25fSSatish Balay   PetscFunctionReturn(0);
42549b5e25fSSatish Balay }
42649b5e25fSSatish Balay 
4274a2ae208SSatish Balay #undef __FUNCT__
4284a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqSBAIJ_5"
429dfbe8321SBarry Smith PetscErrorCode MatMult_SeqSBAIJ_5(Mat A,Vec xx,Vec zz)
43049b5e25fSSatish Balay {
43149b5e25fSSatish Balay   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
43287828ca2SBarry Smith   PetscScalar    *x,*z,*xb,x1,x2,x3,x4,x5,zero=0.0;
43349b5e25fSSatish Balay   MatScalar      *v;
4346849ba73SBarry Smith   PetscErrorCode ierr;
43513f74950SBarry Smith   PetscInt       mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,*ib,cval,j,jmin;
43649b5e25fSSatish Balay 
43749b5e25fSSatish Balay   PetscFunctionBegin;
43849b5e25fSSatish Balay   ierr = VecSet(&zero,zz);CHKERRQ(ierr);
4391ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
4401ebc52fbSHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
44149b5e25fSSatish Balay 
44249b5e25fSSatish Balay   v     = a->a;
44349b5e25fSSatish Balay   xb = x;
44449b5e25fSSatish Balay 
44549b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
44649b5e25fSSatish Balay     n  = ai[1] - ai[0]; /* length of i_th block row of A */
44749b5e25fSSatish Balay     x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5=xb[4];
44849b5e25fSSatish Balay     ib = aj + *ai;
449831a3094SHong Zhang     jmin = 0;
4507fbae186SHong Zhang     if (*ib == i){      /* (diag of A)*x */
45149b5e25fSSatish Balay       z[5*i]   += v[0]*x1  + v[5]*x2 + v[10]*x3 + v[15]*x4+ v[20]*x5;
45249b5e25fSSatish Balay       z[5*i+1] += v[5]*x1  + v[6]*x2 + v[11]*x3 + v[16]*x4+ v[21]*x5;
45349b5e25fSSatish Balay       z[5*i+2] += v[10]*x1 +v[11]*x2 + v[12]*x3 + v[17]*x4+ v[22]*x5;
45449b5e25fSSatish Balay       z[5*i+3] += v[15]*x1 +v[16]*x2 + v[17]*x3 + v[18]*x4+ v[23]*x5;
45549b5e25fSSatish Balay       z[5*i+4] += v[20]*x1 +v[21]*x2 + v[22]*x3 + v[23]*x4+ v[24]*x5;
456831a3094SHong Zhang       v += 25; jmin++;
4577fbae186SHong Zhang     }
458831a3094SHong Zhang     for (j=jmin; j<n; j++) {
45949b5e25fSSatish Balay       /* (strict lower triangular part of A)*x  */
46049b5e25fSSatish Balay       cval       = ib[j]*5;
46149b5e25fSSatish Balay       z[cval]     += v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4 + v[4]*x5;
46249b5e25fSSatish Balay       z[cval+1]   += v[5]*x1 + v[6]*x2 + v[7]*x3 + v[8]*x4 + v[9]*x5;
46349b5e25fSSatish Balay       z[cval+2]   += v[10]*x1 + v[11]*x2 + v[12]*x3 + v[13]*x4+ v[14]*x5;
46449b5e25fSSatish Balay       z[cval+3]   += v[15]*x1 + v[16]*x2 + v[17]*x3 + v[18]*x4+ v[19]*x5;
46549b5e25fSSatish Balay       z[cval+4]   += v[20]*x1 + v[21]*x2 + v[22]*x3 + v[23]*x4+ v[24]*x5;
46649b5e25fSSatish Balay       /* (strict upper triangular part of A)*x  */
46749b5e25fSSatish Balay       z[5*i]   +=v[0]*x[cval]+v[5]*x[cval+1]+v[10]*x[cval+2]+v[15]*x[cval+3]+v[20]*x[cval+4];
46849b5e25fSSatish Balay       z[5*i+1] +=v[1]*x[cval]+v[6]*x[cval+1]+v[11]*x[cval+2]+v[16]*x[cval+3]+v[21]*x[cval+4];
46949b5e25fSSatish Balay       z[5*i+2] +=v[2]*x[cval]+v[7]*x[cval+1]+v[12]*x[cval+2]+v[17]*x[cval+3]+v[22]*x[cval+4];
47049b5e25fSSatish Balay       z[5*i+3] +=v[3]*x[cval]+v[8]*x[cval+1]+v[13]*x[cval+2]+v[18]*x[cval+3]+v[23]*x[cval+4];
47149b5e25fSSatish Balay       z[5*i+4] +=v[4]*x[cval]+v[9]*x[cval+1]+v[14]*x[cval+2]+v[19]*x[cval+3]+v[24]*x[cval+4];
47249b5e25fSSatish Balay       v  += 25;
47349b5e25fSSatish Balay     }
47449b5e25fSSatish Balay     xb +=5; ai++;
47549b5e25fSSatish Balay   }
47649b5e25fSSatish Balay 
4771ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
4781ebc52fbSHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
479efee365bSSatish Balay   ierr = PetscLogFlops(50*(a->nz*2 - A->m) - A->m);CHKERRQ(ierr);
48049b5e25fSSatish Balay   PetscFunctionReturn(0);
48149b5e25fSSatish Balay }
48249b5e25fSSatish Balay 
48349b5e25fSSatish Balay 
4844a2ae208SSatish Balay #undef __FUNCT__
4854a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqSBAIJ_6"
486dfbe8321SBarry Smith PetscErrorCode MatMult_SeqSBAIJ_6(Mat A,Vec xx,Vec zz)
48749b5e25fSSatish Balay {
48849b5e25fSSatish Balay   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
48987828ca2SBarry Smith   PetscScalar    *x,*z,*xb,x1,x2,x3,x4,x5,x6,zero=0.0;
49049b5e25fSSatish Balay   MatScalar      *v;
4916849ba73SBarry Smith   PetscErrorCode ierr;
49213f74950SBarry Smith   PetscInt       mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,*ib,cval,j,jmin;
49349b5e25fSSatish Balay 
49449b5e25fSSatish Balay   PetscFunctionBegin;
49549b5e25fSSatish Balay   ierr = VecSet(&zero,zz);CHKERRQ(ierr);
4961ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
4971ebc52fbSHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
49849b5e25fSSatish Balay 
49949b5e25fSSatish Balay   v     = a->a;
50049b5e25fSSatish Balay   xb = x;
50149b5e25fSSatish Balay 
50249b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
50349b5e25fSSatish Balay     n  = ai[1] - ai[0]; /* length of i_th block row of A */
50449b5e25fSSatish Balay     x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5=xb[4]; x6=xb[5];
50549b5e25fSSatish Balay     ib = aj + *ai;
506831a3094SHong Zhang     jmin = 0;
5077fbae186SHong Zhang     if (*ib == i){      /* (diag of A)*x */
50849b5e25fSSatish Balay       z[6*i]   += v[0]*x1  + v[6]*x2 + v[12]*x3 + v[18]*x4+ v[24]*x5 + v[30]*x6;
50949b5e25fSSatish Balay       z[6*i+1] += v[6]*x1  + v[7]*x2 + v[13]*x3 + v[19]*x4+ v[25]*x5 + v[31]*x6;
51049b5e25fSSatish Balay       z[6*i+2] += v[12]*x1 +v[13]*x2 + v[14]*x3 + v[20]*x4+ v[26]*x5 + v[32]*x6;
51149b5e25fSSatish Balay       z[6*i+3] += v[18]*x1 +v[19]*x2 + v[20]*x3 + v[21]*x4+ v[27]*x5 + v[33]*x6;
51249b5e25fSSatish Balay       z[6*i+4] += v[24]*x1 +v[25]*x2 + v[26]*x3 + v[27]*x4+ v[28]*x5 + v[34]*x6;
51349b5e25fSSatish Balay       z[6*i+5] += v[30]*x1 +v[31]*x2 + v[32]*x3 + v[33]*x4+ v[34]*x5 + v[35]*x6;
514831a3094SHong Zhang       v += 36; jmin++;
5157fbae186SHong Zhang     }
516831a3094SHong Zhang     for (j=jmin; j<n; j++) {
51749b5e25fSSatish Balay       /* (strict lower triangular part of A)*x  */
51849b5e25fSSatish Balay       cval       = ib[j]*6;
51949b5e25fSSatish Balay       z[cval]   += v[0]*x1  + v[1]*x2 + v[2]*x3 + v[3]*x4+ v[4]*x5 + v[5]*x6;
52049b5e25fSSatish Balay       z[cval+1] += v[6]*x1  + v[7]*x2 + v[8]*x3 + v[9]*x4+ v[10]*x5 + v[11]*x6;
52149b5e25fSSatish Balay       z[cval+2] += v[12]*x1  + v[13]*x2 + v[14]*x3 + v[15]*x4+ v[16]*x5 + v[17]*x6;
52249b5e25fSSatish Balay       z[cval+3] += v[18]*x1  + v[19]*x2 + v[20]*x3 + v[21]*x4+ v[22]*x5 + v[23]*x6;
52349b5e25fSSatish Balay       z[cval+4] += v[24]*x1  + v[25]*x2 + v[26]*x3 + v[27]*x4+ v[28]*x5 + v[29]*x6;
52449b5e25fSSatish Balay       z[cval+5] += v[30]*x1  + v[31]*x2 + v[32]*x3 + v[33]*x4+ v[34]*x5 + v[35]*x6;
52549b5e25fSSatish Balay       /* (strict upper triangular part of A)*x  */
52649b5e25fSSatish Balay       z[6*i]   +=v[0]*x[cval]+v[6]*x[cval+1]+v[12]*x[cval+2]+v[18]*x[cval+3]+v[24]*x[cval+4]+v[30]*x[cval+5];
52749b5e25fSSatish Balay       z[6*i+1] +=v[1]*x[cval]+v[7]*x[cval+1]+v[13]*x[cval+2]+v[19]*x[cval+3]+v[25]*x[cval+4]+v[31]*x[cval+5];
52849b5e25fSSatish Balay       z[6*i+2] +=v[2]*x[cval]+v[8]*x[cval+1]+v[14]*x[cval+2]+v[20]*x[cval+3]+v[26]*x[cval+4]+v[32]*x[cval+5];
52949b5e25fSSatish Balay       z[6*i+3] +=v[3]*x[cval]+v[9]*x[cval+1]+v[15]*x[cval+2]+v[21]*x[cval+3]+v[27]*x[cval+4]+v[33]*x[cval+5];
53049b5e25fSSatish Balay       z[6*i+4] +=v[4]*x[cval]+v[10]*x[cval+1]+v[16]*x[cval+2]+v[22]*x[cval+3]+v[28]*x[cval+4]+v[34]*x[cval+5];
53149b5e25fSSatish Balay       z[6*i+5] +=v[5]*x[cval]+v[11]*x[cval+1]+v[17]*x[cval+2]+v[23]*x[cval+3]+v[29]*x[cval+4]+v[35]*x[cval+5];
53249b5e25fSSatish Balay       v  += 36;
53349b5e25fSSatish Balay     }
53449b5e25fSSatish Balay     xb +=6; ai++;
53549b5e25fSSatish Balay   }
53649b5e25fSSatish Balay 
5371ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
5381ebc52fbSHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
539efee365bSSatish Balay   ierr = PetscLogFlops(72*(a->nz*2 - A->m) - A->m);CHKERRQ(ierr);
54049b5e25fSSatish Balay   PetscFunctionReturn(0);
54149b5e25fSSatish Balay }
5424a2ae208SSatish Balay #undef __FUNCT__
5434a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqSBAIJ_7"
544dfbe8321SBarry Smith PetscErrorCode MatMult_SeqSBAIJ_7(Mat A,Vec xx,Vec zz)
54549b5e25fSSatish Balay {
54649b5e25fSSatish Balay   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
54787828ca2SBarry Smith   PetscScalar    *x,*z,*xb,x1,x2,x3,x4,x5,x6,x7,zero=0.0;
54849b5e25fSSatish Balay   MatScalar      *v;
5496849ba73SBarry Smith   PetscErrorCode ierr;
55013f74950SBarry Smith   PetscInt       mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,*ib,cval,j,jmin;
55149b5e25fSSatish Balay 
55249b5e25fSSatish Balay   PetscFunctionBegin;
55349b5e25fSSatish Balay   ierr = VecSet(&zero,zz);CHKERRQ(ierr);
5541ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
5551ebc52fbSHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
55649b5e25fSSatish Balay 
55749b5e25fSSatish Balay   v     = a->a;
55849b5e25fSSatish Balay   xb = x;
55949b5e25fSSatish Balay 
56049b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
56149b5e25fSSatish Balay     n  = ai[1] - ai[0]; /* length of i_th block row of A */
56249b5e25fSSatish Balay     x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5=xb[4]; x6=xb[5]; x7=xb[6];
56349b5e25fSSatish Balay     ib = aj + *ai;
564831a3094SHong Zhang     jmin = 0;
5657fbae186SHong Zhang     if (*ib == i){      /* (diag of A)*x */
56649b5e25fSSatish Balay       z[7*i]   += v[0]*x1 + v[7]*x2 + v[14]*x3 + v[21]*x4+ v[28]*x5 + v[35]*x6+ v[42]*x7;
56749b5e25fSSatish Balay       z[7*i+1] += v[7]*x1 + v[8]*x2 + v[15]*x3 + v[22]*x4+ v[29]*x5 + v[36]*x6+ v[43]*x7;
56849b5e25fSSatish Balay       z[7*i+2] += v[14]*x1+ v[15]*x2 +v[16]*x3 + v[23]*x4+ v[30]*x5 + v[37]*x6+ v[44]*x7;
56949b5e25fSSatish Balay       z[7*i+3] += v[21]*x1+ v[22]*x2 +v[23]*x3 + v[24]*x4+ v[31]*x5 + v[38]*x6+ v[45]*x7;
57049b5e25fSSatish Balay       z[7*i+4] += v[28]*x1+ v[29]*x2 +v[30]*x3 + v[31]*x4+ v[32]*x5 + v[39]*x6+ v[46]*x7;
57149b5e25fSSatish Balay       z[7*i+5] += v[35]*x1+ v[36]*x2 +v[37]*x3 + v[38]*x4+ v[39]*x5 + v[40]*x6+ v[47]*x7;
57249b5e25fSSatish Balay       z[7*i+6] += v[42]*x1+ v[43]*x2 +v[44]*x3 + v[45]*x4+ v[46]*x5 + v[47]*x6+ v[48]*x7;
573831a3094SHong Zhang       v += 49; jmin++;
5747fbae186SHong Zhang     }
575831a3094SHong Zhang     for (j=jmin; j<n; j++) {
57649b5e25fSSatish Balay       /* (strict lower triangular part of A)*x  */
57749b5e25fSSatish Balay       cval       = ib[j]*7;
57849b5e25fSSatish Balay       z[cval]   += v[0]*x1  + v[1]*x2 + v[2]*x3 + v[3]*x4+ v[4]*x5 + v[5]*x6+ v[6]*x7;
57949b5e25fSSatish Balay       z[cval+1] += v[7]*x1  + v[8]*x2 + v[9]*x3 + v[10]*x4+ v[11]*x5 + v[12]*x6+ v[13]*x7;
58049b5e25fSSatish Balay       z[cval+2] += v[14]*x1  + v[15]*x2 + v[16]*x3 + v[17]*x4+ v[18]*x5 + v[19]*x6+ v[20]*x7;
58149b5e25fSSatish Balay       z[cval+3] += v[21]*x1  + v[22]*x2 + v[23]*x3 + v[24]*x4+ v[25]*x5 + v[26]*x6+ v[27]*x7;
58249b5e25fSSatish Balay       z[cval+4] += v[28]*x1  + v[29]*x2 + v[30]*x3 + v[31]*x4+ v[32]*x5 + v[33]*x6+ v[34]*x7;
58349b5e25fSSatish Balay       z[cval+5] += v[35]*x1  + v[36]*x2 + v[37]*x3 + v[38]*x4+ v[39]*x5 + v[40]*x6+ v[41]*x7;
58449b5e25fSSatish Balay       z[cval+6] += v[42]*x1  + v[43]*x2 + v[44]*x3 + v[45]*x4+ v[46]*x5 + v[47]*x6+ v[48]*x7;
58549b5e25fSSatish Balay       /* (strict upper triangular part of A)*x  */
58649b5e25fSSatish Balay       z[7*i]  +=v[0]*x[cval]+v[7]*x[cval+1]+v[14]*x[cval+2]+v[21]*x[cval+3]+v[28]*x[cval+4]+v[35]*x[cval+5]+v[42]*x[cval+6];
58749b5e25fSSatish Balay       z[7*i+1]+=v[1]*x[cval]+v[8]*x[cval+1]+v[15]*x[cval+2]+v[22]*x[cval+3]+v[29]*x[cval+4]+v[36]*x[cval+5]+v[43]*x[cval+6];
58849b5e25fSSatish Balay       z[7*i+2]+=v[2]*x[cval]+v[9]*x[cval+1]+v[16]*x[cval+2]+v[23]*x[cval+3]+v[30]*x[cval+4]+v[37]*x[cval+5]+v[44]*x[cval+6];
58949b5e25fSSatish Balay       z[7*i+3]+=v[3]*x[cval]+v[10]*x[cval+1]+v[17]*x[cval+2]+v[24]*x[cval+3]+v[31]*x[cval+4]+v[38]*x[cval+5]+v[45]*x[cval+6];
59049b5e25fSSatish Balay       z[7*i+4]+=v[4]*x[cval]+v[11]*x[cval+1]+v[18]*x[cval+2]+v[25]*x[cval+3]+v[32]*x[cval+4]+v[39]*x[cval+5]+v[46]*x[cval+6];
59149b5e25fSSatish Balay       z[7*i+5]+=v[5]*x[cval]+v[12]*x[cval+1]+v[19]*x[cval+2]+v[26]*x[cval+3]+v[33]*x[cval+4]+v[40]*x[cval+5]+v[47]*x[cval+6];
59249b5e25fSSatish Balay       z[7*i+6]+=v[6]*x[cval]+v[13]*x[cval+1]+v[20]*x[cval+2]+v[27]*x[cval+3]+v[34]*x[cval+4]+v[41]*x[cval+5]+v[48]*x[cval+6];
59349b5e25fSSatish Balay       v  += 49;
59449b5e25fSSatish Balay     }
59549b5e25fSSatish Balay     xb +=7; ai++;
59649b5e25fSSatish Balay   }
5971ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
5981ebc52fbSHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
599efee365bSSatish Balay   ierr = PetscLogFlops(98*(a->nz*2 - A->m) - A->m);CHKERRQ(ierr);
60049b5e25fSSatish Balay   PetscFunctionReturn(0);
60149b5e25fSSatish Balay }
60249b5e25fSSatish Balay 
60349b5e25fSSatish Balay /*
60449b5e25fSSatish Balay     This will not work with MatScalar == float because it calls the BLAS
60549b5e25fSSatish Balay */
6064a2ae208SSatish Balay #undef __FUNCT__
6074a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqSBAIJ_N"
608dfbe8321SBarry Smith PetscErrorCode MatMult_SeqSBAIJ_N(Mat A,Vec xx,Vec zz)
60949b5e25fSSatish Balay {
61049b5e25fSSatish Balay   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
61187828ca2SBarry Smith   PetscScalar    *x,*x_ptr,*z,*z_ptr,*xb,*zb,*work,*workt,zero=0.0;
6120b60a74dSHong Zhang   MatScalar      *v;
613dfbe8321SBarry Smith   PetscErrorCode ierr;
614521d7252SBarry Smith   PetscInt       mbs=a->mbs,i,*idx,*aj,*ii,bs=A->bs,j,n,bs2=a->bs2,ncols,k;
61549b5e25fSSatish Balay 
61649b5e25fSSatish Balay   PetscFunctionBegin;
61749b5e25fSSatish Balay   ierr = VecSet(&zero,zz);CHKERRQ(ierr);
6181ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr); x_ptr=x;
6191ebc52fbSHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr); z_ptr=z;
62049b5e25fSSatish Balay 
62149b5e25fSSatish Balay   aj   = a->j;
62249b5e25fSSatish Balay   v    = a->a;
62349b5e25fSSatish Balay   ii   = a->i;
62449b5e25fSSatish Balay 
62549b5e25fSSatish Balay   if (!a->mult_work) {
62687828ca2SBarry Smith     ierr = PetscMalloc((A->m+1)*sizeof(PetscScalar),&a->mult_work);CHKERRQ(ierr);
62749b5e25fSSatish Balay   }
62849b5e25fSSatish Balay   work = a->mult_work;
62949b5e25fSSatish Balay 
63049b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
63149b5e25fSSatish Balay     n     = ii[1] - ii[0]; ncols = n*bs;
63249b5e25fSSatish Balay     workt = work; idx=aj+ii[0];
63349b5e25fSSatish Balay 
63449b5e25fSSatish Balay     /* upper triangular part */
63549b5e25fSSatish Balay     for (j=0; j<n; j++) {
63649b5e25fSSatish Balay       xb = x_ptr + bs*(*idx++);
63749b5e25fSSatish Balay       for (k=0; k<bs; k++) workt[k] = xb[k];
63849b5e25fSSatish Balay       workt += bs;
63949b5e25fSSatish Balay     }
64049b5e25fSSatish Balay     /* z(i*bs:(i+1)*bs-1) += A(i,:)*x */
64149b5e25fSSatish Balay     Kernel_w_gets_w_plus_Ar_times_v(bs,ncols,work,v,z);
64249b5e25fSSatish Balay 
64349b5e25fSSatish Balay     /* strict lower triangular part */
644831a3094SHong Zhang     idx = aj+ii[0];
645831a3094SHong Zhang     if (*idx == i){
64696b9376eSHong Zhang       ncols -= bs; v += bs2; idx++; n--;
647831a3094SHong Zhang     }
64896b9376eSHong Zhang 
64949b5e25fSSatish Balay     if (ncols > 0){
65049b5e25fSSatish Balay       workt = work;
65187828ca2SBarry Smith       ierr  = PetscMemzero(workt,ncols*sizeof(PetscScalar));CHKERRQ(ierr);
652831a3094SHong Zhang       Kernel_w_gets_w_plus_trans_Ar_times_v(bs,ncols,x,v,workt);
653831a3094SHong Zhang       for (j=0; j<n; j++) {
654831a3094SHong Zhang         zb = z_ptr + bs*(*idx++);
65549b5e25fSSatish Balay         for (k=0; k<bs; k++) zb[k] += workt[k] ;
65649b5e25fSSatish Balay         workt += bs;
65749b5e25fSSatish Balay       }
65849b5e25fSSatish Balay     }
65949b5e25fSSatish Balay     x += bs; v += n*bs2; z += bs; ii++;
66049b5e25fSSatish Balay   }
66149b5e25fSSatish Balay 
6621ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
6631ebc52fbSHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
664efee365bSSatish Balay   ierr = PetscLogFlops(2*(a->nz*2 - A->m)*bs2 - A->m);CHKERRQ(ierr);
66549b5e25fSSatish Balay   PetscFunctionReturn(0);
66649b5e25fSSatish Balay }
66749b5e25fSSatish Balay 
6684a2ae208SSatish Balay #undef __FUNCT__
6694a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqSBAIJ_1"
670dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqSBAIJ_1(Mat A,Vec xx,Vec yy,Vec zz)
67149b5e25fSSatish Balay {
67249b5e25fSSatish Balay   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
67387828ca2SBarry Smith   PetscScalar    *x,*y,*z,*xb,x1;
67449b5e25fSSatish Balay   MatScalar      *v;
6756849ba73SBarry Smith   PetscErrorCode ierr;
67613f74950SBarry Smith   PetscInt       mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,*ib,cval,j,jmin;
67749b5e25fSSatish Balay 
67849b5e25fSSatish Balay   PetscFunctionBegin;
6791ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
68049b5e25fSSatish Balay   if (yy != xx) {
6811ebc52fbSHong Zhang     ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
68249b5e25fSSatish Balay   } else {
68349b5e25fSSatish Balay     y = x;
68449b5e25fSSatish Balay   }
68549b5e25fSSatish Balay   if (zz != yy) {
686a9d4b620SHong Zhang     /* ierr = VecCopy(yy,zz);CHKERRQ(ierr); */
6871ebc52fbSHong Zhang     ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
688a9d4b620SHong Zhang     ierr = PetscMemcpy(z,y,yy->n*sizeof(PetscScalar));CHKERRQ(ierr);
68949b5e25fSSatish Balay   } else {
69049b5e25fSSatish Balay     z = y;
69149b5e25fSSatish Balay   }
69249b5e25fSSatish Balay 
69349b5e25fSSatish Balay   v  = a->a;
69449b5e25fSSatish Balay   xb = x;
69549b5e25fSSatish Balay 
69649b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
69749b5e25fSSatish Balay     n  = ai[1] - ai[0];  /* length of i_th row of A */
69849b5e25fSSatish Balay     x1 = xb[0];
69949b5e25fSSatish Balay     ib = aj + *ai;
700831a3094SHong Zhang     jmin = 0;
701831a3094SHong Zhang     if (*ib == i) {            /* (diag of A)*x */
702831a3094SHong Zhang       z[i] += *v++ * x[*ib++]; jmin++;
703831a3094SHong Zhang     }
704831a3094SHong Zhang     for (j=jmin; j<n; j++) {
70549b5e25fSSatish Balay       cval    = *ib;
70649b5e25fSSatish Balay       z[cval] += *v * x1;      /* (strict lower triangular part of A)*x  */
70749b5e25fSSatish Balay       z[i] += *v++ * x[*ib++]; /* (strict upper triangular part of A)*x  */
70849b5e25fSSatish Balay     }
70949b5e25fSSatish Balay     xb++; ai++;
71049b5e25fSSatish Balay   }
71149b5e25fSSatish Balay 
7121ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
7131ebc52fbSHong Zhang   if (yy != xx) ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
7141ebc52fbSHong Zhang   if (zz != yy) ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
71549b5e25fSSatish Balay 
716efee365bSSatish Balay   ierr = PetscLogFlops(2*(a->nz*2 - A->m));CHKERRQ(ierr);
71749b5e25fSSatish Balay   PetscFunctionReturn(0);
71849b5e25fSSatish Balay }
71949b5e25fSSatish Balay 
7204a2ae208SSatish Balay #undef __FUNCT__
7214a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqSBAIJ_2"
722dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqSBAIJ_2(Mat A,Vec xx,Vec yy,Vec zz)
72349b5e25fSSatish Balay {
72449b5e25fSSatish Balay   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
72587828ca2SBarry Smith   PetscScalar    *x,*y,*z,*xb,x1,x2;
72649b5e25fSSatish Balay   MatScalar      *v;
7276849ba73SBarry Smith   PetscErrorCode ierr;
72813f74950SBarry Smith   PetscInt       mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,*ib,cval,j,jmin;
72949b5e25fSSatish Balay 
73049b5e25fSSatish Balay   PetscFunctionBegin;
7311ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
73249b5e25fSSatish Balay   if (yy != xx) {
7331ebc52fbSHong Zhang     ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
73449b5e25fSSatish Balay   } else {
73549b5e25fSSatish Balay     y = x;
73649b5e25fSSatish Balay   }
73749b5e25fSSatish Balay   if (zz != yy) {
738a9d4b620SHong Zhang     /* ierr = VecCopy(yy,zz);CHKERRQ(ierr); */
7391ebc52fbSHong Zhang     ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
740a9d4b620SHong Zhang     ierr = PetscMemcpy(z,y,yy->n*sizeof(PetscScalar));CHKERRQ(ierr);
74149b5e25fSSatish Balay   } else {
74249b5e25fSSatish Balay     z = y;
74349b5e25fSSatish Balay   }
74449b5e25fSSatish Balay 
74549b5e25fSSatish Balay   v     = a->a;
74649b5e25fSSatish Balay   xb = x;
74749b5e25fSSatish Balay 
74849b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
74949b5e25fSSatish Balay     n  = ai[1] - ai[0]; /* length of i_th block row of A */
75049b5e25fSSatish Balay     x1 = xb[0]; x2 = xb[1];
75149b5e25fSSatish Balay     ib = aj + *ai;
752831a3094SHong Zhang     jmin = 0;
7537fbae186SHong Zhang     if (*ib == i){      /* (diag of A)*x */
75449b5e25fSSatish Balay       z[2*i]   += v[0]*x1 + v[2]*x2;
75549b5e25fSSatish Balay       z[2*i+1] += v[2]*x1 + v[3]*x2;
756831a3094SHong Zhang       v += 4; jmin++;
7577fbae186SHong Zhang     }
758831a3094SHong Zhang     for (j=jmin; j<n; j++) {
75949b5e25fSSatish Balay       /* (strict lower triangular part of A)*x  */
76049b5e25fSSatish Balay       cval       = ib[j]*2;
76149b5e25fSSatish Balay       z[cval]     += v[0]*x1 + v[1]*x2;
76249b5e25fSSatish Balay       z[cval+1]   += v[2]*x1 + v[3]*x2;
76349b5e25fSSatish Balay       /* (strict upper triangular part of A)*x  */
76449b5e25fSSatish Balay       z[2*i]   += v[0]*x[cval] + v[2]*x[cval+1];
76549b5e25fSSatish Balay       z[2*i+1] += v[1]*x[cval] + v[3]*x[cval+1];
76649b5e25fSSatish Balay       v  += 4;
76749b5e25fSSatish Balay     }
76849b5e25fSSatish Balay     xb +=2; ai++;
76949b5e25fSSatish Balay   }
77049b5e25fSSatish Balay 
7711ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
7721ebc52fbSHong Zhang   if (yy != xx) ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
7731ebc52fbSHong Zhang   if (zz != yy) ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
77449b5e25fSSatish Balay 
775efee365bSSatish Balay   ierr = PetscLogFlops(4*(a->nz*2 - A->m));CHKERRQ(ierr);
77649b5e25fSSatish Balay   PetscFunctionReturn(0);
77749b5e25fSSatish Balay }
77849b5e25fSSatish Balay 
7794a2ae208SSatish Balay #undef __FUNCT__
7804a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqSBAIJ_3"
781dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqSBAIJ_3(Mat A,Vec xx,Vec yy,Vec zz)
78249b5e25fSSatish Balay {
78349b5e25fSSatish Balay   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
78487828ca2SBarry Smith   PetscScalar    *x,*y,*z,*xb,x1,x2,x3;
78549b5e25fSSatish Balay   MatScalar      *v;
7866849ba73SBarry Smith   PetscErrorCode ierr;
78713f74950SBarry Smith   PetscInt       mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,*ib,cval,j,jmin;
78849b5e25fSSatish Balay 
78949b5e25fSSatish Balay   PetscFunctionBegin;
7901ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
79149b5e25fSSatish Balay   if (yy != xx) {
7921ebc52fbSHong Zhang     ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
79349b5e25fSSatish Balay   } else {
79449b5e25fSSatish Balay     y = x;
79549b5e25fSSatish Balay   }
79649b5e25fSSatish Balay   if (zz != yy) {
797a9d4b620SHong Zhang     /* ierr = VecCopy(yy,zz);CHKERRQ(ierr); */
7981ebc52fbSHong Zhang     ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
799a9d4b620SHong Zhang     ierr = PetscMemcpy(z,y,yy->n*sizeof(PetscScalar));CHKERRQ(ierr);
80049b5e25fSSatish Balay   } else {
80149b5e25fSSatish Balay     z = y;
80249b5e25fSSatish Balay   }
80349b5e25fSSatish Balay 
80449b5e25fSSatish Balay   v     = a->a;
80549b5e25fSSatish Balay   xb = x;
80649b5e25fSSatish Balay 
80749b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
80849b5e25fSSatish Balay     n  = ai[1] - ai[0]; /* length of i_th block row of A */
80949b5e25fSSatish Balay     x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
81049b5e25fSSatish Balay     ib = aj + *ai;
811831a3094SHong Zhang     jmin = 0;
8127fbae186SHong Zhang     if (*ib == i){     /* (diag of A)*x */
81349b5e25fSSatish Balay      z[3*i]   += v[0]*x1 + v[3]*x2 + v[6]*x3;
81449b5e25fSSatish Balay      z[3*i+1] += v[3]*x1 + v[4]*x2 + v[7]*x3;
81549b5e25fSSatish Balay      z[3*i+2] += v[6]*x1 + v[7]*x2 + v[8]*x3;
816831a3094SHong Zhang      v += 9; jmin++;
8177fbae186SHong Zhang     }
818831a3094SHong Zhang     for (j=jmin; j<n; j++) {
81949b5e25fSSatish Balay       /* (strict lower triangular part of A)*x  */
82049b5e25fSSatish Balay       cval       = ib[j]*3;
82149b5e25fSSatish Balay       z[cval]     += v[0]*x1 + v[1]*x2 + v[2]*x3;
82249b5e25fSSatish Balay       z[cval+1]   += v[3]*x1 + v[4]*x2 + v[5]*x3;
82349b5e25fSSatish Balay       z[cval+2]   += v[6]*x1 + v[7]*x2 + v[8]*x3;
82449b5e25fSSatish Balay       /* (strict upper triangular part of A)*x  */
82549b5e25fSSatish Balay       z[3*i]   += v[0]*x[cval] + v[3]*x[cval+1]+ v[6]*x[cval+2];
82649b5e25fSSatish Balay       z[3*i+1] += v[1]*x[cval] + v[4]*x[cval+1]+ v[7]*x[cval+2];
82749b5e25fSSatish Balay       z[3*i+2] += v[2]*x[cval] + v[5]*x[cval+1]+ v[8]*x[cval+2];
82849b5e25fSSatish Balay       v  += 9;
82949b5e25fSSatish Balay     }
83049b5e25fSSatish Balay     xb +=3; ai++;
83149b5e25fSSatish Balay   }
83249b5e25fSSatish Balay 
8331ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
8341ebc52fbSHong Zhang   if (yy != xx) ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
8351ebc52fbSHong Zhang   if (zz != yy) ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
83649b5e25fSSatish Balay 
837efee365bSSatish Balay   ierr = PetscLogFlops(18*(a->nz*2 - A->m));CHKERRQ(ierr);
83849b5e25fSSatish Balay   PetscFunctionReturn(0);
83949b5e25fSSatish Balay }
84049b5e25fSSatish Balay 
8414a2ae208SSatish Balay #undef __FUNCT__
8424a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqSBAIJ_4"
843dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqSBAIJ_4(Mat A,Vec xx,Vec yy,Vec zz)
84449b5e25fSSatish Balay {
84549b5e25fSSatish Balay   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
84687828ca2SBarry Smith   PetscScalar    *x,*y,*z,*xb,x1,x2,x3,x4;
84749b5e25fSSatish Balay   MatScalar      *v;
8486849ba73SBarry Smith   PetscErrorCode ierr;
84913f74950SBarry Smith   PetscInt       mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,*ib,cval,j,jmin;
85049b5e25fSSatish Balay 
85149b5e25fSSatish Balay   PetscFunctionBegin;
8521ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
85349b5e25fSSatish Balay   if (yy != xx) {
8541ebc52fbSHong Zhang     ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
85549b5e25fSSatish Balay   } else {
85649b5e25fSSatish Balay     y = x;
85749b5e25fSSatish Balay   }
85849b5e25fSSatish Balay   if (zz != yy) {
859a9d4b620SHong Zhang     /* ierr = VecCopy(yy,zz);CHKERRQ(ierr); */
8601ebc52fbSHong Zhang     ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
861a9d4b620SHong Zhang     ierr = PetscMemcpy(z,y,yy->n*sizeof(PetscScalar));CHKERRQ(ierr);
86249b5e25fSSatish Balay   } else {
86349b5e25fSSatish Balay     z = y;
86449b5e25fSSatish Balay   }
86549b5e25fSSatish Balay 
86649b5e25fSSatish Balay   v     = a->a;
86749b5e25fSSatish Balay   xb = x;
86849b5e25fSSatish Balay 
86949b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
87049b5e25fSSatish Balay     n  = ai[1] - ai[0]; /* length of i_th block row of A */
87149b5e25fSSatish Balay     x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
87249b5e25fSSatish Balay     ib = aj + *ai;
873831a3094SHong Zhang     jmin = 0;
8747fbae186SHong Zhang     if (*ib == i){      /* (diag of A)*x */
87549b5e25fSSatish Balay       z[4*i]   += v[0]*x1 + v[4]*x2 +  v[8]*x3 + v[12]*x4;
87649b5e25fSSatish Balay       z[4*i+1] += v[4]*x1 + v[5]*x2 +  v[9]*x3 + v[13]*x4;
87749b5e25fSSatish Balay       z[4*i+2] += v[8]*x1 + v[9]*x2 + v[10]*x3 + v[14]*x4;
87849b5e25fSSatish Balay       z[4*i+3] += v[12]*x1+ v[13]*x2+ v[14]*x3 + v[15]*x4;
879831a3094SHong Zhang       v += 16; jmin++;
8807fbae186SHong Zhang     }
881831a3094SHong Zhang     for (j=jmin; j<n; j++) {
88249b5e25fSSatish Balay       /* (strict lower triangular part of A)*x  */
88349b5e25fSSatish Balay       cval       = ib[j]*4;
88449b5e25fSSatish Balay       z[cval]     += v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4;
88549b5e25fSSatish Balay       z[cval+1]   += v[4]*x1 + v[5]*x2 + v[6]*x3 + v[7]*x4;
88649b5e25fSSatish Balay       z[cval+2]   += v[8]*x1 + v[9]*x2 + v[10]*x3 + v[11]*x4;
88749b5e25fSSatish Balay       z[cval+3]   += v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4;
88849b5e25fSSatish Balay       /* (strict upper triangular part of A)*x  */
88949b5e25fSSatish Balay       z[4*i]   += v[0]*x[cval] + v[4]*x[cval+1]+ v[8]*x[cval+2] + v[12]*x[cval+3];
89049b5e25fSSatish Balay       z[4*i+1] += v[1]*x[cval] + v[5]*x[cval+1]+ v[9]*x[cval+2] + v[13]*x[cval+3];
89149b5e25fSSatish Balay       z[4*i+2] += v[2]*x[cval] + v[6]*x[cval+1]+ v[10]*x[cval+2]+ v[14]*x[cval+3];
89249b5e25fSSatish Balay       z[4*i+3] += v[3]*x[cval] + v[7]*x[cval+1]+ v[11]*x[cval+2]+ v[15]*x[cval+3];
89349b5e25fSSatish Balay       v  += 16;
89449b5e25fSSatish Balay     }
89549b5e25fSSatish Balay     xb +=4; ai++;
89649b5e25fSSatish Balay   }
89749b5e25fSSatish Balay 
8981ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
8991ebc52fbSHong Zhang   if (yy != xx) ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
9001ebc52fbSHong Zhang   if (zz != yy) ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
90149b5e25fSSatish Balay 
902efee365bSSatish Balay   ierr = PetscLogFlops(32*(a->nz*2 - A->m));CHKERRQ(ierr);
90349b5e25fSSatish Balay   PetscFunctionReturn(0);
90449b5e25fSSatish Balay }
90549b5e25fSSatish Balay 
9064a2ae208SSatish Balay #undef __FUNCT__
9074a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqSBAIJ_5"
908dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqSBAIJ_5(Mat A,Vec xx,Vec yy,Vec zz)
90949b5e25fSSatish Balay {
91049b5e25fSSatish Balay   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
91187828ca2SBarry Smith   PetscScalar    *x,*y,*z,*xb,x1,x2,x3,x4,x5;
91249b5e25fSSatish Balay   MatScalar      *v;
9136849ba73SBarry Smith   PetscErrorCode ierr;
91413f74950SBarry Smith   PetscInt       mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,*ib,cval,j,jmin;
91549b5e25fSSatish Balay 
91649b5e25fSSatish Balay   PetscFunctionBegin;
9171ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
91849b5e25fSSatish Balay   if (yy != xx) {
9191ebc52fbSHong Zhang     ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
92049b5e25fSSatish Balay   } else {
92149b5e25fSSatish Balay     y = x;
92249b5e25fSSatish Balay   }
92349b5e25fSSatish Balay   if (zz != yy) {
924a9d4b620SHong Zhang     /* ierr = VecCopy(yy,zz);CHKERRQ(ierr); */
9251ebc52fbSHong Zhang     ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
926a9d4b620SHong Zhang     ierr = PetscMemcpy(z,y,yy->n*sizeof(PetscScalar));CHKERRQ(ierr);
92749b5e25fSSatish Balay   } else {
92849b5e25fSSatish Balay     z = y;
92949b5e25fSSatish Balay   }
93049b5e25fSSatish Balay 
93149b5e25fSSatish Balay   v     = a->a;
93249b5e25fSSatish Balay   xb = x;
93349b5e25fSSatish Balay 
93449b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
93549b5e25fSSatish Balay     n  = ai[1] - ai[0]; /* length of i_th block row of A */
93649b5e25fSSatish Balay     x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5=xb[4];
93749b5e25fSSatish Balay     ib = aj + *ai;
938831a3094SHong Zhang     jmin = 0;
9397fbae186SHong Zhang     if (*ib == i){      /* (diag of A)*x */
94049b5e25fSSatish Balay       z[5*i]   += v[0]*x1  + v[5]*x2 + v[10]*x3 + v[15]*x4+ v[20]*x5;
94149b5e25fSSatish Balay       z[5*i+1] += v[5]*x1  + v[6]*x2 + v[11]*x3 + v[16]*x4+ v[21]*x5;
94249b5e25fSSatish Balay       z[5*i+2] += v[10]*x1 +v[11]*x2 + v[12]*x3 + v[17]*x4+ v[22]*x5;
94349b5e25fSSatish Balay       z[5*i+3] += v[15]*x1 +v[16]*x2 + v[17]*x3 + v[18]*x4+ v[23]*x5;
94449b5e25fSSatish Balay       z[5*i+4] += v[20]*x1 +v[21]*x2 + v[22]*x3 + v[23]*x4+ v[24]*x5;
945831a3094SHong Zhang       v += 25; jmin++;
9467fbae186SHong Zhang     }
947831a3094SHong Zhang     for (j=jmin; j<n; j++) {
94849b5e25fSSatish Balay       /* (strict lower triangular part of A)*x  */
94949b5e25fSSatish Balay       cval       = ib[j]*5;
95049b5e25fSSatish Balay       z[cval]     += v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4 + v[4]*x5;
95149b5e25fSSatish Balay       z[cval+1]   += v[5]*x1 + v[6]*x2 + v[7]*x3 + v[8]*x4 + v[9]*x5;
95249b5e25fSSatish Balay       z[cval+2]   += v[10]*x1 + v[11]*x2 + v[12]*x3 + v[13]*x4+ v[14]*x5;
95349b5e25fSSatish Balay       z[cval+3]   += v[15]*x1 + v[16]*x2 + v[17]*x3 + v[18]*x4+ v[19]*x5;
95449b5e25fSSatish Balay       z[cval+4]   += v[20]*x1 + v[21]*x2 + v[22]*x3 + v[23]*x4+ v[24]*x5;
95549b5e25fSSatish Balay       /* (strict upper triangular part of A)*x  */
95649b5e25fSSatish Balay       z[5*i]   +=v[0]*x[cval]+v[5]*x[cval+1]+v[10]*x[cval+2]+v[15]*x[cval+3]+v[20]*x[cval+4];
95749b5e25fSSatish Balay       z[5*i+1] +=v[1]*x[cval]+v[6]*x[cval+1]+v[11]*x[cval+2]+v[16]*x[cval+3]+v[21]*x[cval+4];
95849b5e25fSSatish Balay       z[5*i+2] +=v[2]*x[cval]+v[7]*x[cval+1]+v[12]*x[cval+2]+v[17]*x[cval+3]+v[22]*x[cval+4];
95949b5e25fSSatish Balay       z[5*i+3] +=v[3]*x[cval]+v[8]*x[cval+1]+v[13]*x[cval+2]+v[18]*x[cval+3]+v[23]*x[cval+4];
96049b5e25fSSatish Balay       z[5*i+4] +=v[4]*x[cval]+v[9]*x[cval+1]+v[14]*x[cval+2]+v[19]*x[cval+3]+v[24]*x[cval+4];
96149b5e25fSSatish Balay       v  += 25;
96249b5e25fSSatish Balay     }
96349b5e25fSSatish Balay     xb +=5; ai++;
96449b5e25fSSatish Balay   }
96549b5e25fSSatish Balay 
9661ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
9671ebc52fbSHong Zhang   if (yy != xx) ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
9681ebc52fbSHong Zhang   if (zz != yy) ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
96949b5e25fSSatish Balay 
970efee365bSSatish Balay   ierr = PetscLogFlops(50*(a->nz*2 - A->m));CHKERRQ(ierr);
97149b5e25fSSatish Balay   PetscFunctionReturn(0);
97249b5e25fSSatish Balay }
9734a2ae208SSatish Balay #undef __FUNCT__
9744a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqSBAIJ_6"
975dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqSBAIJ_6(Mat A,Vec xx,Vec yy,Vec zz)
97649b5e25fSSatish Balay {
97749b5e25fSSatish Balay   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
97887828ca2SBarry Smith   PetscScalar    *x,*y,*z,*xb,x1,x2,x3,x4,x5,x6;
97949b5e25fSSatish Balay   MatScalar      *v;
9806849ba73SBarry Smith   PetscErrorCode ierr;
98113f74950SBarry Smith   PetscInt       mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,*ib,cval,j,jmin;
98249b5e25fSSatish Balay 
98349b5e25fSSatish Balay   PetscFunctionBegin;
9841ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
98549b5e25fSSatish Balay   if (yy != xx) {
9861ebc52fbSHong Zhang     ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
98749b5e25fSSatish Balay   } else {
98849b5e25fSSatish Balay     y = x;
98949b5e25fSSatish Balay   }
99049b5e25fSSatish Balay   if (zz != yy) {
9911ebc52fbSHong Zhang     ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
992a9d4b620SHong Zhang     ierr = PetscMemcpy(z,y,yy->n*sizeof(PetscScalar));CHKERRQ(ierr);
99349b5e25fSSatish Balay   } else {
99449b5e25fSSatish Balay     z = y;
99549b5e25fSSatish Balay   }
99649b5e25fSSatish Balay 
99749b5e25fSSatish Balay   v     = a->a;
99849b5e25fSSatish Balay   xb = x;
99949b5e25fSSatish Balay 
100049b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
100149b5e25fSSatish Balay     n  = ai[1] - ai[0]; /* length of i_th block row of A */
100249b5e25fSSatish Balay     x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5=xb[4]; x6=xb[5];
100349b5e25fSSatish Balay     ib = aj + *ai;
1004831a3094SHong Zhang     jmin = 0;
10057fbae186SHong Zhang     if (*ib == i){     /* (diag of A)*x */
100649b5e25fSSatish Balay       z[6*i]   += v[0]*x1  + v[6]*x2 + v[12]*x3 + v[18]*x4+ v[24]*x5 + v[30]*x6;
100749b5e25fSSatish Balay       z[6*i+1] += v[6]*x1  + v[7]*x2 + v[13]*x3 + v[19]*x4+ v[25]*x5 + v[31]*x6;
100849b5e25fSSatish Balay       z[6*i+2] += v[12]*x1 +v[13]*x2 + v[14]*x3 + v[20]*x4+ v[26]*x5 + v[32]*x6;
100949b5e25fSSatish Balay       z[6*i+3] += v[18]*x1 +v[19]*x2 + v[20]*x3 + v[21]*x4+ v[27]*x5 + v[33]*x6;
101049b5e25fSSatish Balay       z[6*i+4] += v[24]*x1 +v[25]*x2 + v[26]*x3 + v[27]*x4+ v[28]*x5 + v[34]*x6;
101149b5e25fSSatish Balay       z[6*i+5] += v[30]*x1 +v[31]*x2 + v[32]*x3 + v[33]*x4+ v[34]*x5 + v[35]*x6;
1012831a3094SHong Zhang       v += 36; jmin++;
10137fbae186SHong Zhang     }
1014831a3094SHong Zhang     for (j=jmin; j<n; j++) {
101549b5e25fSSatish Balay       /* (strict lower triangular part of A)*x  */
101649b5e25fSSatish Balay       cval       = ib[j]*6;
101749b5e25fSSatish Balay       z[cval]   += v[0]*x1  + v[1]*x2 + v[2]*x3 + v[3]*x4+ v[4]*x5 + v[5]*x6;
101849b5e25fSSatish Balay       z[cval+1] += v[6]*x1  + v[7]*x2 + v[8]*x3 + v[9]*x4+ v[10]*x5 + v[11]*x6;
101949b5e25fSSatish Balay       z[cval+2] += v[12]*x1  + v[13]*x2 + v[14]*x3 + v[15]*x4+ v[16]*x5 + v[17]*x6;
102049b5e25fSSatish Balay       z[cval+3] += v[18]*x1  + v[19]*x2 + v[20]*x3 + v[21]*x4+ v[22]*x5 + v[23]*x6;
102149b5e25fSSatish Balay       z[cval+4] += v[24]*x1  + v[25]*x2 + v[26]*x3 + v[27]*x4+ v[28]*x5 + v[29]*x6;
102249b5e25fSSatish Balay       z[cval+5] += v[30]*x1  + v[31]*x2 + v[32]*x3 + v[33]*x4+ v[34]*x5 + v[35]*x6;
102349b5e25fSSatish Balay       /* (strict upper triangular part of A)*x  */
102449b5e25fSSatish Balay       z[6*i]   +=v[0]*x[cval]+v[6]*x[cval+1]+v[12]*x[cval+2]+v[18]*x[cval+3]+v[24]*x[cval+4]+v[30]*x[cval+5];
102549b5e25fSSatish Balay       z[6*i+1] +=v[1]*x[cval]+v[7]*x[cval+1]+v[13]*x[cval+2]+v[19]*x[cval+3]+v[25]*x[cval+4]+v[31]*x[cval+5];
102649b5e25fSSatish Balay       z[6*i+2] +=v[2]*x[cval]+v[8]*x[cval+1]+v[14]*x[cval+2]+v[20]*x[cval+3]+v[26]*x[cval+4]+v[32]*x[cval+5];
102749b5e25fSSatish Balay       z[6*i+3] +=v[3]*x[cval]+v[9]*x[cval+1]+v[15]*x[cval+2]+v[21]*x[cval+3]+v[27]*x[cval+4]+v[33]*x[cval+5];
102849b5e25fSSatish Balay       z[6*i+4] +=v[4]*x[cval]+v[10]*x[cval+1]+v[16]*x[cval+2]+v[22]*x[cval+3]+v[28]*x[cval+4]+v[34]*x[cval+5];
102949b5e25fSSatish Balay       z[6*i+5] +=v[5]*x[cval]+v[11]*x[cval+1]+v[17]*x[cval+2]+v[23]*x[cval+3]+v[29]*x[cval+4]+v[35]*x[cval+5];
103049b5e25fSSatish Balay       v  += 36;
103149b5e25fSSatish Balay     }
103249b5e25fSSatish Balay     xb +=6; ai++;
103349b5e25fSSatish Balay   }
103449b5e25fSSatish Balay 
10351ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
10361ebc52fbSHong Zhang   if (yy != xx) ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
10371ebc52fbSHong Zhang   if (zz != yy) ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
103849b5e25fSSatish Balay 
1039efee365bSSatish Balay   ierr = PetscLogFlops(72*(a->nz*2 - A->m));CHKERRQ(ierr);
104049b5e25fSSatish Balay   PetscFunctionReturn(0);
104149b5e25fSSatish Balay }
104249b5e25fSSatish Balay 
10434a2ae208SSatish Balay #undef __FUNCT__
10444a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqSBAIJ_7"
1045dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqSBAIJ_7(Mat A,Vec xx,Vec yy,Vec zz)
104649b5e25fSSatish Balay {
104749b5e25fSSatish Balay   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
104887828ca2SBarry Smith   PetscScalar    *x,*y,*z,*xb,x1,x2,x3,x4,x5,x6,x7;
104949b5e25fSSatish Balay   MatScalar      *v;
10506849ba73SBarry Smith   PetscErrorCode ierr;
105113f74950SBarry Smith   PetscInt       mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,*ib,cval,j,jmin;
105249b5e25fSSatish Balay 
105349b5e25fSSatish Balay   PetscFunctionBegin;
10541ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
105549b5e25fSSatish Balay   if (yy != xx) {
10561ebc52fbSHong Zhang     ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
105749b5e25fSSatish Balay   } else {
105849b5e25fSSatish Balay     y = x;
105949b5e25fSSatish Balay   }
106049b5e25fSSatish Balay   if (zz != yy) {
10611ebc52fbSHong Zhang     ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
1062a9d4b620SHong Zhang     ierr = PetscMemcpy(z,y,yy->n*sizeof(PetscScalar));CHKERRQ(ierr);
106349b5e25fSSatish Balay   } else {
106449b5e25fSSatish Balay     z = y;
106549b5e25fSSatish Balay   }
106649b5e25fSSatish Balay 
106749b5e25fSSatish Balay   v     = a->a;
106849b5e25fSSatish Balay   xb = x;
106949b5e25fSSatish Balay 
107049b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
107149b5e25fSSatish Balay     n  = ai[1] - ai[0]; /* length of i_th block row of A */
107249b5e25fSSatish Balay     x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5=xb[4]; x6=xb[5]; x7=xb[6];
107349b5e25fSSatish Balay     ib = aj + *ai;
1074831a3094SHong Zhang     jmin = 0;
10757fbae186SHong Zhang     if (*ib == i){     /* (diag of A)*x */
107649b5e25fSSatish Balay       z[7*i]   += v[0]*x1 + v[7]*x2 + v[14]*x3 + v[21]*x4+ v[28]*x5 + v[35]*x6+ v[42]*x7;
107749b5e25fSSatish Balay       z[7*i+1] += v[7]*x1 + v[8]*x2 + v[15]*x3 + v[22]*x4+ v[29]*x5 + v[36]*x6+ v[43]*x7;
107849b5e25fSSatish Balay       z[7*i+2] += v[14]*x1+ v[15]*x2 +v[16]*x3 + v[23]*x4+ v[30]*x5 + v[37]*x6+ v[44]*x7;
107949b5e25fSSatish Balay       z[7*i+3] += v[21]*x1+ v[22]*x2 +v[23]*x3 + v[24]*x4+ v[31]*x5 + v[38]*x6+ v[45]*x7;
108049b5e25fSSatish Balay       z[7*i+4] += v[28]*x1+ v[29]*x2 +v[30]*x3 + v[31]*x4+ v[32]*x5 + v[39]*x6+ v[46]*x7;
108149b5e25fSSatish Balay       z[7*i+5] += v[35]*x1+ v[36]*x2 +v[37]*x3 + v[38]*x4+ v[39]*x5 + v[40]*x6+ v[47]*x7;
108249b5e25fSSatish Balay       z[7*i+6] += v[42]*x1+ v[43]*x2 +v[44]*x3 + v[45]*x4+ v[46]*x5 + v[47]*x6+ v[48]*x7;
1083831a3094SHong Zhang       v += 49; jmin++;
10847fbae186SHong Zhang     }
1085831a3094SHong Zhang     for (j=jmin; j<n; j++) {
108649b5e25fSSatish Balay       /* (strict lower triangular part of A)*x  */
108749b5e25fSSatish Balay       cval       = ib[j]*7;
108849b5e25fSSatish Balay       z[cval]   += v[0]*x1  + v[1]*x2 + v[2]*x3 + v[3]*x4+ v[4]*x5 + v[5]*x6+ v[6]*x7;
108949b5e25fSSatish Balay       z[cval+1] += v[7]*x1  + v[8]*x2 + v[9]*x3 + v[10]*x4+ v[11]*x5 + v[12]*x6+ v[13]*x7;
109049b5e25fSSatish Balay       z[cval+2] += v[14]*x1  + v[15]*x2 + v[16]*x3 + v[17]*x4+ v[18]*x5 + v[19]*x6+ v[20]*x7;
109149b5e25fSSatish Balay       z[cval+3] += v[21]*x1  + v[22]*x2 + v[23]*x3 + v[24]*x4+ v[25]*x5 + v[26]*x6+ v[27]*x7;
109249b5e25fSSatish Balay       z[cval+4] += v[28]*x1  + v[29]*x2 + v[30]*x3 + v[31]*x4+ v[32]*x5 + v[33]*x6+ v[34]*x7;
109349b5e25fSSatish Balay       z[cval+5] += v[35]*x1  + v[36]*x2 + v[37]*x3 + v[38]*x4+ v[39]*x5 + v[40]*x6+ v[41]*x7;
109449b5e25fSSatish Balay       z[cval+6] += v[42]*x1  + v[43]*x2 + v[44]*x3 + v[45]*x4+ v[46]*x5 + v[47]*x6+ v[48]*x7;
109549b5e25fSSatish Balay       /* (strict upper triangular part of A)*x  */
109649b5e25fSSatish Balay       z[7*i]  +=v[0]*x[cval]+v[7]*x[cval+1]+v[14]*x[cval+2]+v[21]*x[cval+3]+v[28]*x[cval+4]+v[35]*x[cval+5]+v[42]*x[cval+6];
109749b5e25fSSatish Balay       z[7*i+1]+=v[1]*x[cval]+v[8]*x[cval+1]+v[15]*x[cval+2]+v[22]*x[cval+3]+v[29]*x[cval+4]+v[36]*x[cval+5]+v[43]*x[cval+6];
109849b5e25fSSatish Balay       z[7*i+2]+=v[2]*x[cval]+v[9]*x[cval+1]+v[16]*x[cval+2]+v[23]*x[cval+3]+v[30]*x[cval+4]+v[37]*x[cval+5]+v[44]*x[cval+6];
109949b5e25fSSatish Balay       z[7*i+3]+=v[3]*x[cval]+v[10]*x[cval+1]+v[17]*x[cval+2]+v[24]*x[cval+3]+v[31]*x[cval+4]+v[38]*x[cval+5]+v[45]*x[cval+6];
110049b5e25fSSatish Balay       z[7*i+4]+=v[4]*x[cval]+v[11]*x[cval+1]+v[18]*x[cval+2]+v[25]*x[cval+3]+v[32]*x[cval+4]+v[39]*x[cval+5]+v[46]*x[cval+6];
110149b5e25fSSatish Balay       z[7*i+5]+=v[5]*x[cval]+v[12]*x[cval+1]+v[19]*x[cval+2]+v[26]*x[cval+3]+v[33]*x[cval+4]+v[40]*x[cval+5]+v[47]*x[cval+6];
110249b5e25fSSatish Balay       z[7*i+6]+=v[6]*x[cval]+v[13]*x[cval+1]+v[20]*x[cval+2]+v[27]*x[cval+3]+v[34]*x[cval+4]+v[41]*x[cval+5]+v[48]*x[cval+6];
110349b5e25fSSatish Balay       v  += 49;
110449b5e25fSSatish Balay     }
110549b5e25fSSatish Balay     xb +=7; ai++;
110649b5e25fSSatish Balay   }
110749b5e25fSSatish Balay 
11081ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
11091ebc52fbSHong Zhang   if (yy != xx) ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
11101ebc52fbSHong Zhang   if (zz != yy) ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
111149b5e25fSSatish Balay 
1112efee365bSSatish Balay   ierr = PetscLogFlops(98*(a->nz*2 - A->m));CHKERRQ(ierr);
111349b5e25fSSatish Balay   PetscFunctionReturn(0);
111449b5e25fSSatish Balay }
111549b5e25fSSatish Balay 
11164a2ae208SSatish Balay #undef __FUNCT__
11174a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqSBAIJ_N"
1118dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqSBAIJ_N(Mat A,Vec xx,Vec yy,Vec zz)
111949b5e25fSSatish Balay {
112049b5e25fSSatish Balay   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
112187828ca2SBarry Smith   PetscScalar    *x,*x_ptr,*y,*z,*z_ptr=0,*xb,*zb,*work,*workt;
1122066653e3SSatish Balay   MatScalar      *v;
1123dfbe8321SBarry Smith   PetscErrorCode ierr;
1124521d7252SBarry Smith   PetscInt       mbs=a->mbs,i,*idx,*aj,*ii,bs=A->bs,j,n,bs2=a->bs2,ncols,k;
112549b5e25fSSatish Balay 
112649b5e25fSSatish Balay   PetscFunctionBegin;
11271ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr); x_ptr=x;
112849b5e25fSSatish Balay   if (yy != xx) {
11291ebc52fbSHong Zhang     ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
113049b5e25fSSatish Balay   } else {
113149b5e25fSSatish Balay     y = x;
113249b5e25fSSatish Balay   }
113349b5e25fSSatish Balay   if (zz != yy) {
11341ebc52fbSHong Zhang     ierr = VecGetArray(zz,&z);CHKERRQ(ierr); z_ptr=z;
1135a9d4b620SHong Zhang     ierr = PetscMemcpy(z,y,yy->n*sizeof(PetscScalar));CHKERRQ(ierr);
113649b5e25fSSatish Balay   } else {
113749b5e25fSSatish Balay     z = y;
113849b5e25fSSatish Balay   }
113949b5e25fSSatish Balay 
114049b5e25fSSatish Balay   aj   = a->j;
114149b5e25fSSatish Balay   v    = a->a;
114249b5e25fSSatish Balay   ii   = a->i;
114349b5e25fSSatish Balay 
114449b5e25fSSatish Balay   if (!a->mult_work) {
114587828ca2SBarry Smith     ierr = PetscMalloc((A->m+1)*sizeof(PetscScalar),&a->mult_work);CHKERRQ(ierr);
114649b5e25fSSatish Balay   }
114749b5e25fSSatish Balay   work = a->mult_work;
114849b5e25fSSatish Balay 
114949b5e25fSSatish Balay 
115049b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
115149b5e25fSSatish Balay     n     = ii[1] - ii[0]; ncols = n*bs;
115249b5e25fSSatish Balay     workt = work; idx=aj+ii[0];
115349b5e25fSSatish Balay 
115449b5e25fSSatish Balay     /* upper triangular part */
115549b5e25fSSatish Balay     for (j=0; j<n; j++) {
115649b5e25fSSatish Balay       xb = x_ptr + bs*(*idx++);
115749b5e25fSSatish Balay       for (k=0; k<bs; k++) workt[k] = xb[k];
115849b5e25fSSatish Balay       workt += bs;
115949b5e25fSSatish Balay     }
116049b5e25fSSatish Balay     /* z(i*bs:(i+1)*bs-1) += A(i,:)*x */
116149b5e25fSSatish Balay     Kernel_w_gets_w_plus_Ar_times_v(bs,ncols,work,v,z);
116249b5e25fSSatish Balay 
116349b5e25fSSatish Balay     /* strict lower triangular part */
1164831a3094SHong Zhang     idx = aj+ii[0];
1165831a3094SHong Zhang     if (*idx == i){
116696b9376eSHong Zhang       ncols -= bs; v += bs2; idx++; n--;
1167831a3094SHong Zhang     }
116849b5e25fSSatish Balay     if (ncols > 0){
116949b5e25fSSatish Balay       workt = work;
117087828ca2SBarry Smith       ierr  = PetscMemzero(workt,ncols*sizeof(PetscScalar));CHKERRQ(ierr);
1171831a3094SHong Zhang       Kernel_w_gets_w_plus_trans_Ar_times_v(bs,ncols,x,v,workt);
1172831a3094SHong Zhang       for (j=0; j<n; j++) {
1173831a3094SHong Zhang         zb = z_ptr + bs*(*idx++);
117449b5e25fSSatish Balay         for (k=0; k<bs; k++) zb[k] += workt[k] ;
117549b5e25fSSatish Balay         workt += bs;
117649b5e25fSSatish Balay       }
117749b5e25fSSatish Balay     }
117849b5e25fSSatish Balay 
117949b5e25fSSatish Balay     x += bs; v += n*bs2; z += bs; ii++;
118049b5e25fSSatish Balay   }
118149b5e25fSSatish Balay 
11821ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
11831ebc52fbSHong Zhang   if (yy != xx) ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
11841ebc52fbSHong Zhang   if (zz != yy) ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
118549b5e25fSSatish Balay 
1186efee365bSSatish Balay   ierr = PetscLogFlops(2*(a->nz*2 - A->m));CHKERRQ(ierr);
118749b5e25fSSatish Balay   PetscFunctionReturn(0);
118849b5e25fSSatish Balay }
118949b5e25fSSatish Balay 
11904a2ae208SSatish Balay #undef __FUNCT__
11914a2ae208SSatish Balay #define __FUNCT__ "MatScale_SeqSBAIJ"
1192dfbe8321SBarry Smith PetscErrorCode MatScale_SeqSBAIJ(const PetscScalar *alpha,Mat inA)
119349b5e25fSSatish Balay {
119449b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)inA->data;
11954ce68768SBarry Smith   PetscBLASInt one = 1,totalnz = (PetscBLASInt)a->bs2*a->nz;
1196efee365bSSatish Balay   PetscErrorCode ierr;
119749b5e25fSSatish Balay 
119849b5e25fSSatish Balay   PetscFunctionBegin;
119971044d3cSBarry Smith   BLASscal_(&totalnz,(PetscScalar*)alpha,a->a,&one);
1200efee365bSSatish Balay   ierr = PetscLogFlops(totalnz);CHKERRQ(ierr);
120149b5e25fSSatish Balay   PetscFunctionReturn(0);
120249b5e25fSSatish Balay }
120349b5e25fSSatish Balay 
12044a2ae208SSatish Balay #undef __FUNCT__
12054a2ae208SSatish Balay #define __FUNCT__ "MatNorm_SeqSBAIJ"
1206dfbe8321SBarry Smith PetscErrorCode MatNorm_SeqSBAIJ(Mat A,NormType type,PetscReal *norm)
120749b5e25fSSatish Balay {
120849b5e25fSSatish Balay   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
120949b5e25fSSatish Balay   MatScalar      *v = a->a;
121049b5e25fSSatish Balay   PetscReal      sum_diag = 0.0, sum_off = 0.0, *sum;
121113f74950SBarry Smith   PetscInt       i,j,k,bs = A->bs,bs2=a->bs2,k1,mbs=a->mbs,*aj=a->j;
12126849ba73SBarry Smith   PetscErrorCode ierr;
121313f74950SBarry Smith   PetscInt       *jl,*il,jmin,jmax,nexti,ik,*col;
121449b5e25fSSatish Balay 
121549b5e25fSSatish Balay   PetscFunctionBegin;
121649b5e25fSSatish Balay   if (type == NORM_FROBENIUS) {
121749b5e25fSSatish Balay     for (k=0; k<mbs; k++){
121849b5e25fSSatish Balay       jmin = a->i[k]; jmax = a->i[k+1];
1219831a3094SHong Zhang       col  = aj + jmin;
1220831a3094SHong Zhang       if (*col == k){         /* diagonal block */
122149b5e25fSSatish Balay         for (i=0; i<bs2; i++){
122249b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX)
122349b5e25fSSatish Balay           sum_diag += PetscRealPart(PetscConj(*v)*(*v)); v++;
122449b5e25fSSatish Balay #else
122549b5e25fSSatish Balay           sum_diag += (*v)*(*v); v++;
122649b5e25fSSatish Balay #endif
122749b5e25fSSatish Balay         }
1228831a3094SHong Zhang         jmin++;
1229831a3094SHong Zhang       }
1230831a3094SHong Zhang       for (j=jmin; j<jmax; j++){  /* off-diagonal blocks */
123149b5e25fSSatish Balay         for (i=0; i<bs2; i++){
123249b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX)
123349b5e25fSSatish Balay           sum_off += PetscRealPart(PetscConj(*v)*(*v)); v++;
123449b5e25fSSatish Balay #else
123549b5e25fSSatish Balay           sum_off += (*v)*(*v); v++;
123649b5e25fSSatish Balay #endif
123749b5e25fSSatish Balay         }
123849b5e25fSSatish Balay       }
123949b5e25fSSatish Balay     }
124049b5e25fSSatish Balay     *norm = sqrt(sum_diag + 2*sum_off);
1241*0b8dc8d2SHong Zhang   }  else if (type == NORM_INFINITY || type == NORM_1) { /* maximum row/column sum */
1242*0b8dc8d2SHong Zhang     ierr = PetscMalloc((2*mbs+1)*sizeof(PetscInt)+bs*sizeof(PetscReal),&il);CHKERRQ(ierr);
1243*0b8dc8d2SHong Zhang     jl   = il + mbs;
1244*0b8dc8d2SHong Zhang     sum  = (PetscReal*)(jl + mbs);
1245*0b8dc8d2SHong Zhang     for (i=0; i<mbs; i++) jl[i] = mbs;
1246*0b8dc8d2SHong Zhang     il[0] = 0;
124749b5e25fSSatish Balay 
124849b5e25fSSatish Balay     *norm = 0.0;
124949b5e25fSSatish Balay     for (k=0; k<mbs; k++) { /* k_th block row */
125049b5e25fSSatish Balay       for (j=0; j<bs; j++) sum[j]=0.0;
125149b5e25fSSatish Balay       /*-- col sum --*/
125249b5e25fSSatish Balay       i = jl[k]; /* first |A(i,k)| to be added */
1253831a3094SHong Zhang       /* jl[k]=i: first nozero element in row i for submatrix A(1:k,k:n) (active window)
1254831a3094SHong Zhang                   at step k */
125549b5e25fSSatish Balay       while (i<mbs){
125649b5e25fSSatish Balay         nexti = jl[i];  /* next block row to be added */
125749b5e25fSSatish Balay         ik    = il[i];  /* block index of A(i,k) in the array a */
125849b5e25fSSatish Balay         for (j=0; j<bs; j++){
125949b5e25fSSatish Balay           v = a->a + ik*bs2 + j*bs;
126049b5e25fSSatish Balay           for (k1=0; k1<bs; k1++) {
126149b5e25fSSatish Balay             sum[j] += PetscAbsScalar(*v); v++;
126249b5e25fSSatish Balay           }
126349b5e25fSSatish Balay         }
126449b5e25fSSatish Balay         /* update il, jl */
1265831a3094SHong Zhang         jmin = ik + 1; /* block index of array a: points to the next nonzero of A in row i */
1266831a3094SHong Zhang         jmax = a->i[i+1];
126749b5e25fSSatish Balay         if (jmin < jmax){
126849b5e25fSSatish Balay           il[i] = jmin;
126949b5e25fSSatish Balay           j   = a->j[jmin];
127049b5e25fSSatish Balay           jl[i] = jl[j]; jl[j]=i;
127149b5e25fSSatish Balay         }
127249b5e25fSSatish Balay         i = nexti;
127349b5e25fSSatish Balay       }
127449b5e25fSSatish Balay       /*-- row sum --*/
127549b5e25fSSatish Balay       jmin = a->i[k]; jmax = a->i[k+1];
127649b5e25fSSatish Balay       for (i=jmin; i<jmax; i++) {
127749b5e25fSSatish Balay         for (j=0; j<bs; j++){
127849b5e25fSSatish Balay           v = a->a + i*bs2 + j;
127949b5e25fSSatish Balay           for (k1=0; k1<bs; k1++){
1280*0b8dc8d2SHong Zhang             sum[j] += PetscAbsScalar(*v); v += bs;
128149b5e25fSSatish Balay           }
128249b5e25fSSatish Balay         }
128349b5e25fSSatish Balay       }
128449b5e25fSSatish Balay       /* add k_th block row to il, jl */
1285831a3094SHong Zhang       col = aj+jmin;
1286831a3094SHong Zhang       if (*col == k) jmin++;
128749b5e25fSSatish Balay       if (jmin < jmax){
128849b5e25fSSatish Balay         il[k] = jmin;
1289*0b8dc8d2SHong Zhang         j = a->j[jmin]; jl[k] = jl[j]; jl[j] = k;
129049b5e25fSSatish Balay       }
129149b5e25fSSatish Balay       for (j=0; j<bs; j++){
129249b5e25fSSatish Balay         if (sum[j] > *norm) *norm = sum[j];
129349b5e25fSSatish Balay       }
129449b5e25fSSatish Balay     }
129549b5e25fSSatish Balay     ierr = PetscFree(il);CHKERRQ(ierr);
129649b5e25fSSatish Balay   } else {
1297347d480fSBarry Smith     SETERRQ(PETSC_ERR_SUP,"No support for this norm yet");
129849b5e25fSSatish Balay   }
129949b5e25fSSatish Balay   PetscFunctionReturn(0);
130049b5e25fSSatish Balay }
130149b5e25fSSatish Balay 
13024a2ae208SSatish Balay #undef __FUNCT__
13034a2ae208SSatish Balay #define __FUNCT__ "MatEqual_SeqSBAIJ"
1304dfbe8321SBarry Smith PetscErrorCode MatEqual_SeqSBAIJ(Mat A,Mat B,PetscTruth* flg)
130549b5e25fSSatish Balay {
130649b5e25fSSatish Balay   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ *)A->data,*b = (Mat_SeqSBAIJ *)B->data;
1307dfbe8321SBarry Smith   PetscErrorCode ierr;
130849b5e25fSSatish Balay 
130949b5e25fSSatish Balay   PetscFunctionBegin;
131049b5e25fSSatish Balay 
131149b5e25fSSatish Balay   /* If the  matrix/block dimensions are not equal, or no of nonzeros or shift */
1312521d7252SBarry Smith   if ((A->m != B->m) || (A->n != B->n) || (A->bs != B->bs)|| (a->nz != b->nz)) {
1313ef511fbeSHong Zhang     *flg = PETSC_FALSE;
1314ef511fbeSHong Zhang     PetscFunctionReturn(0);
131549b5e25fSSatish Balay   }
131649b5e25fSSatish Balay 
131749b5e25fSSatish Balay   /* if the a->i are the same */
131813f74950SBarry Smith   ierr = PetscMemcmp(a->i,b->i,(a->mbs+1)*sizeof(PetscInt),flg);CHKERRQ(ierr);
1319abc0a331SBarry Smith   if (!*flg) {
132049b5e25fSSatish Balay     PetscFunctionReturn(0);
132149b5e25fSSatish Balay   }
132249b5e25fSSatish Balay 
132349b5e25fSSatish Balay   /* if a->j are the same */
132413f74950SBarry Smith   ierr = PetscMemcmp(a->j,b->j,(a->nz)*sizeof(PetscInt),flg);CHKERRQ(ierr);
1325abc0a331SBarry Smith   if (!*flg) {
132649b5e25fSSatish Balay     PetscFunctionReturn(0);
132749b5e25fSSatish Balay   }
132849b5e25fSSatish Balay   /* if a->a are the same */
1329521d7252SBarry Smith   ierr = PetscMemcmp(a->a,b->a,(a->nz)*(A->bs)*(A->bs)*sizeof(PetscScalar),flg);CHKERRQ(ierr);
1330935af2e7SHong Zhang   PetscFunctionReturn(0);
133149b5e25fSSatish Balay }
133249b5e25fSSatish Balay 
13334a2ae208SSatish Balay #undef __FUNCT__
13344a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_SeqSBAIJ"
1335dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_SeqSBAIJ(Mat A,Vec v)
133649b5e25fSSatish Balay {
133749b5e25fSSatish Balay   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
1338dfbe8321SBarry Smith   PetscErrorCode ierr;
133913f74950SBarry Smith   PetscInt       i,j,k,n,row,bs,*ai,*aj,ambs,bs2;
134087828ca2SBarry Smith   PetscScalar    *x,zero = 0.0;
134149b5e25fSSatish Balay   MatScalar      *aa,*aa_j;
134249b5e25fSSatish Balay 
134349b5e25fSSatish Balay   PetscFunctionBegin;
1344521d7252SBarry Smith   bs   = A->bs;
134582799104SHong Zhang   if (A->factor && bs>1) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix with bs>1");
134682799104SHong Zhang 
134749b5e25fSSatish Balay   aa   = a->a;
134849b5e25fSSatish Balay   ai   = a->i;
134949b5e25fSSatish Balay   aj   = a->j;
135049b5e25fSSatish Balay   ambs = a->mbs;
135149b5e25fSSatish Balay   bs2  = a->bs2;
135249b5e25fSSatish Balay 
135349b5e25fSSatish Balay   ierr = VecSet(&zero,v);CHKERRQ(ierr);
13541ebc52fbSHong Zhang   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
135549b5e25fSSatish Balay   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
1356ef511fbeSHong Zhang   if (n != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
135749b5e25fSSatish Balay   for (i=0; i<ambs; i++) {
135849b5e25fSSatish Balay     j=ai[i];
135949b5e25fSSatish Balay     if (aj[j] == i) {             /* if this is a diagonal element */
136049b5e25fSSatish Balay       row  = i*bs;
136149b5e25fSSatish Balay       aa_j = aa + j*bs2;
136282799104SHong Zhang       if (A->factor && bs==1){
136382799104SHong Zhang         for (k=0; k<bs2; k+=(bs+1),row++) x[row] = 1.0/aa_j[k];
136482799104SHong Zhang       } else {
136549b5e25fSSatish Balay         for (k=0; k<bs2; k+=(bs+1),row++) x[row] = aa_j[k];
136649b5e25fSSatish Balay       }
136749b5e25fSSatish Balay     }
136882799104SHong Zhang   }
136982799104SHong Zhang 
13701ebc52fbSHong Zhang   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
137149b5e25fSSatish Balay   PetscFunctionReturn(0);
137249b5e25fSSatish Balay }
137349b5e25fSSatish Balay 
13744a2ae208SSatish Balay #undef __FUNCT__
13754a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_SeqSBAIJ"
1376dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_SeqSBAIJ(Mat A,Vec ll,Vec rr)
137749b5e25fSSatish Balay {
137849b5e25fSSatish Balay   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
13795e90f9d9SHong Zhang   PetscScalar    *l,x,*li,*ri;
138049b5e25fSSatish Balay   MatScalar      *aa,*v;
1381dfbe8321SBarry Smith   PetscErrorCode ierr;
13825e90f9d9SHong Zhang   PetscInt       i,j,k,lm,M,m,*ai,*aj,mbs,tmp,bs,bs2;
1383b3bf805bSHong Zhang   PetscTruth     flg;
138449b5e25fSSatish Balay 
138549b5e25fSSatish Balay   PetscFunctionBegin;
1386b3bf805bSHong Zhang   if (ll != rr){
1387b3bf805bSHong Zhang     ierr = VecEqual(ll,rr,&flg);CHKERRQ(ierr);
1388b3bf805bSHong Zhang     if (!flg)
1389b3bf805bSHong Zhang       SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"For symmetric format, left and right scaling vectors must be same\n");
1390b3bf805bSHong Zhang   }
1391b3bf805bSHong Zhang   if (!ll) PetscFunctionReturn(0);
139249b5e25fSSatish Balay   ai  = a->i;
139349b5e25fSSatish Balay   aj  = a->j;
139449b5e25fSSatish Balay   aa  = a->a;
1395ef511fbeSHong Zhang   m   = A->m;
1396521d7252SBarry Smith   bs  = A->bs;
139749b5e25fSSatish Balay   mbs = a->mbs;
139849b5e25fSSatish Balay   bs2 = a->bs2;
139949b5e25fSSatish Balay 
14001ebc52fbSHong Zhang   ierr = VecGetArray(ll,&l);CHKERRQ(ierr);
140149b5e25fSSatish Balay   ierr = VecGetLocalSize(ll,&lm);CHKERRQ(ierr);
1402347d480fSBarry Smith   if (lm != m) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vector wrong length");
140349b5e25fSSatish Balay   for (i=0; i<mbs; i++) { /* for each block row */
140449b5e25fSSatish Balay     M  = ai[i+1] - ai[i];
140549b5e25fSSatish Balay     li = l + i*bs;
140649b5e25fSSatish Balay     v  = aa + bs2*ai[i];
140749b5e25fSSatish Balay     for (j=0; j<M; j++) { /* for each block */
140849b5e25fSSatish Balay       ri = l + bs*aj[ai[i]+j];
14095e90f9d9SHong Zhang       for (k=0; k<bs; k++) {
141049b5e25fSSatish Balay         x = ri[k];
141149b5e25fSSatish Balay         for (tmp=0; tmp<bs; tmp++) (*v++) *= li[tmp]*x;
141249b5e25fSSatish Balay       }
141349b5e25fSSatish Balay     }
141449b5e25fSSatish Balay   }
14151ebc52fbSHong Zhang   ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr);
1416efee365bSSatish Balay   ierr = PetscLogFlops(2*a->nz);CHKERRQ(ierr);
141749b5e25fSSatish Balay   PetscFunctionReturn(0);
141849b5e25fSSatish Balay }
141949b5e25fSSatish Balay 
14204a2ae208SSatish Balay #undef __FUNCT__
14214a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_SeqSBAIJ"
1422dfbe8321SBarry Smith PetscErrorCode MatGetInfo_SeqSBAIJ(Mat A,MatInfoType flag,MatInfo *info)
142349b5e25fSSatish Balay {
142449b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
142549b5e25fSSatish Balay 
142649b5e25fSSatish Balay   PetscFunctionBegin;
1427ef511fbeSHong Zhang   info->rows_global    = (double)A->m;
1428ef511fbeSHong Zhang   info->columns_global = (double)A->m;
1429ef511fbeSHong Zhang   info->rows_local     = (double)A->m;
1430ef511fbeSHong Zhang   info->columns_local  = (double)A->m;
143149b5e25fSSatish Balay   info->block_size     = a->bs2;
14326c6c5352SBarry Smith   info->nz_allocated   = a->maxnz; /*num. of nonzeros in upper triangular part */
14336c6c5352SBarry Smith   info->nz_used        = a->bs2*a->nz; /*num. of nonzeros in upper triangular part */
143449b5e25fSSatish Balay   info->nz_unneeded    = (double)(info->nz_allocated - info->nz_used);
143549b5e25fSSatish Balay   info->assemblies   = A->num_ass;
143649b5e25fSSatish Balay   info->mallocs      = a->reallocs;
143749b5e25fSSatish Balay   info->memory       = A->mem;
143849b5e25fSSatish Balay   if (A->factor) {
143949b5e25fSSatish Balay     info->fill_ratio_given  = A->info.fill_ratio_given;
144049b5e25fSSatish Balay     info->fill_ratio_needed = A->info.fill_ratio_needed;
144149b5e25fSSatish Balay     info->factor_mallocs    = A->info.factor_mallocs;
144249b5e25fSSatish Balay   } else {
144349b5e25fSSatish Balay     info->fill_ratio_given  = 0;
144449b5e25fSSatish Balay     info->fill_ratio_needed = 0;
144549b5e25fSSatish Balay     info->factor_mallocs    = 0;
144649b5e25fSSatish Balay   }
144749b5e25fSSatish Balay   PetscFunctionReturn(0);
144849b5e25fSSatish Balay }
144949b5e25fSSatish Balay 
145049b5e25fSSatish Balay 
14514a2ae208SSatish Balay #undef __FUNCT__
14524a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_SeqSBAIJ"
1453dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_SeqSBAIJ(Mat A)
145449b5e25fSSatish Balay {
145549b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
1456dfbe8321SBarry Smith   PetscErrorCode ierr;
145749b5e25fSSatish Balay 
145849b5e25fSSatish Balay   PetscFunctionBegin;
145949b5e25fSSatish Balay   ierr = PetscMemzero(a->a,a->bs2*a->i[a->mbs]*sizeof(MatScalar));CHKERRQ(ierr);
146049b5e25fSSatish Balay   PetscFunctionReturn(0);
146149b5e25fSSatish Balay }
1462dc354874SHong Zhang 
14634a2ae208SSatish Balay #undef __FUNCT__
14644a2ae208SSatish Balay #define __FUNCT__ "MatGetRowMax_SeqSBAIJ"
1465dfbe8321SBarry Smith PetscErrorCode MatGetRowMax_SeqSBAIJ(Mat A,Vec v)
1466dc354874SHong Zhang {
1467dc354874SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
1468dfbe8321SBarry Smith   PetscErrorCode ierr;
146913f74950SBarry Smith   PetscInt i,j,n,row,col,bs,*ai,*aj,mbs;
1470c3fca9a7SHong Zhang   PetscReal    atmp;
1471273d9f13SBarry Smith   MatScalar    *aa;
147287828ca2SBarry Smith   PetscScalar  zero = 0.0,*x;
147313f74950SBarry Smith   PetscInt          ncols,brow,bcol,krow,kcol;
1474dc354874SHong Zhang 
1475dc354874SHong Zhang   PetscFunctionBegin;
1476dc354874SHong Zhang   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1477521d7252SBarry Smith   bs   = A->bs;
1478dc354874SHong Zhang   aa   = a->a;
1479dc354874SHong Zhang   ai   = a->i;
1480dc354874SHong Zhang   aj   = a->j;
148144117c81SHong Zhang   mbs = a->mbs;
1482dc354874SHong Zhang 
1483dc354874SHong Zhang   ierr = VecSet(&zero,v);CHKERRQ(ierr);
14841ebc52fbSHong Zhang   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1485dc354874SHong Zhang   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
1486ef511fbeSHong Zhang   if (n != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
148744117c81SHong Zhang   for (i=0; i<mbs; i++) {
1488d0f6400bSHong Zhang     ncols = ai[1] - ai[0]; ai++;
1489d0f6400bSHong Zhang     brow  = bs*i;
149044117c81SHong Zhang     for (j=0; j<ncols; j++){
1491d0f6400bSHong Zhang       bcol = bs*(*aj);
149244117c81SHong Zhang       for (kcol=0; kcol<bs; kcol++){
1493d0f6400bSHong Zhang         col = bcol + kcol;      /* col index */
149444117c81SHong Zhang         for (krow=0; krow<bs; krow++){
1495d0f6400bSHong Zhang           atmp = PetscAbsScalar(*aa); aa++;
1496d0f6400bSHong Zhang           row = brow + krow;    /* row index */
1497c3fca9a7SHong Zhang           if (PetscRealPart(x[row]) < atmp) x[row] = atmp;
1498c3fca9a7SHong Zhang           if (*aj > i && PetscRealPart(x[col]) < atmp) x[col] = atmp;
149944117c81SHong Zhang         }
150044117c81SHong Zhang       }
1501d0f6400bSHong Zhang       aj++;
1502dc354874SHong Zhang     }
1503dc354874SHong Zhang   }
15041ebc52fbSHong Zhang   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1505dc354874SHong Zhang   PetscFunctionReturn(0);
1506dc354874SHong Zhang }
1507