xref: /petsc/src/mat/impls/sbaij/seq/sbaij2.c (revision 7fbae186bacbe2f615b0bc86cd5234b9970bdcb7)
1*7fbae186SHong Zhang /*$Id: sbaij2.c,v 1.13 2000/09/28 20:51:54 bsmith Exp hzhang $*/
249b5e25fSSatish Balay 
349b5e25fSSatish Balay #include "petscsys.h"
449b5e25fSSatish Balay #include "src/mat/impls/baij/seq/baij.h"
549b5e25fSSatish Balay #include "src/vec/vecimpl.h"
649b5e25fSSatish Balay #include "src/inline/spops.h"
749b5e25fSSatish Balay #include "src/inline/ilu.h"
849b5e25fSSatish Balay #include "petscbt.h"
949b5e25fSSatish Balay #include "sbaij.h"
1049b5e25fSSatish Balay 
1149b5e25fSSatish Balay #undef __FUNC__
1249b5e25fSSatish Balay #define __FUNC__ "MatIncreaseOverlap_SeqSBAIJ"
1349b5e25fSSatish Balay int MatIncreaseOverlap_SeqSBAIJ(Mat A,int is_max,IS *is,int ov)
1449b5e25fSSatish Balay {
1549b5e25fSSatish Balay   PetscFunctionBegin;
160eaf17bfSBarry Smith   SETERRQ(1,"Function not yet written for SBAIJ format");
175535bfb1SHong Zhang   /* PetscFunctionReturn(0); */
1849b5e25fSSatish Balay }
1949b5e25fSSatish Balay 
2049b5e25fSSatish Balay #undef __FUNC__
2149b5e25fSSatish Balay #define __FUNC__ "MatGetSubMatrix_SeqSBAIJ_Private"
2249b5e25fSSatish Balay int MatGetSubMatrix_SeqSBAIJ_Private(Mat A,IS isrow,IS iscol,int cs,MatReuse scall,Mat *B)
2349b5e25fSSatish Balay {
2449b5e25fSSatish Balay   Mat_SeqSBAIJ  *a = (Mat_SeqSBAIJ*)A->data,*c;
2549b5e25fSSatish Balay   int          *smap,i,k,kstart,kend,ierr,oldcols = a->mbs,*lens;
2649b5e25fSSatish Balay   int          row,mat_i,*mat_j,tcol,*mat_ilen;
2749b5e25fSSatish Balay   int          *irow,nrows,*ssmap,bs=a->bs,bs2=a->bs2;
2849b5e25fSSatish Balay   int          *aj = a->j,*ai = a->i;
2949b5e25fSSatish Balay   MatScalar    *mat_a;
3049b5e25fSSatish Balay   Mat          C;
3149b5e25fSSatish Balay   PetscTruth   flag;
3249b5e25fSSatish Balay 
3349b5e25fSSatish Balay   PetscFunctionBegin;
3449b5e25fSSatish Balay 
35347d480fSBarry Smith   if (isrow != iscol) SETERRA(1,"MatGetSubmatrices_SeqSBAIJ: For symm. format, iscol must equal isro");
3649b5e25fSSatish Balay   ierr = ISSorted(iscol,(PetscTruth*)&i);CHKERRQ(ierr);
37347d480fSBarry Smith   if (!i) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"IS is not sorted");
3849b5e25fSSatish Balay 
3949b5e25fSSatish Balay   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
4049b5e25fSSatish Balay   ierr = ISGetSize(isrow,&nrows);CHKERRQ(ierr);
4149b5e25fSSatish Balay 
4249b5e25fSSatish Balay   smap  = (int*)PetscMalloc((1+oldcols)*sizeof(int));CHKPTRQ(smap);
4349b5e25fSSatish Balay   ssmap = smap;
4449b5e25fSSatish Balay   lens  = (int*)PetscMalloc((1+nrows)*sizeof(int));CHKPTRQ(lens);
4549b5e25fSSatish Balay   ierr  = PetscMemzero(smap,oldcols*sizeof(int));CHKERRQ(ierr);
4649b5e25fSSatish Balay   for (i=0; i<nrows; i++) smap[irow[i]] = i+1; /* nrows = ncols */
4749b5e25fSSatish Balay   /* determine lens of each row */
4849b5e25fSSatish Balay   for (i=0; i<nrows; i++) {
4949b5e25fSSatish Balay     kstart  = ai[irow[i]];
5049b5e25fSSatish Balay     kend    = kstart + a->ilen[irow[i]];
5149b5e25fSSatish Balay     lens[i] = 0;
5249b5e25fSSatish Balay       for (k=kstart; k<kend; k++) {
5349b5e25fSSatish Balay         if (ssmap[aj[k]]) {
5449b5e25fSSatish Balay           lens[i]++;
5549b5e25fSSatish Balay         }
5649b5e25fSSatish Balay       }
5749b5e25fSSatish Balay     }
5849b5e25fSSatish Balay   /* Create and fill new matrix */
5949b5e25fSSatish Balay   if (scall == MAT_REUSE_MATRIX) {
6049b5e25fSSatish Balay     c = (Mat_SeqSBAIJ *)((*B)->data);
6149b5e25fSSatish Balay 
62347d480fSBarry Smith     if (c->mbs!=nrows || c->bs!=bs) SETERRQ(PETSC_ERR_ARG_SIZ,"Submatrix wrong size");
6349b5e25fSSatish Balay     ierr = PetscMemcmp(c->ilen,lens,c->mbs *sizeof(int),&flag);CHKERRQ(ierr);
6449b5e25fSSatish Balay     if (flag == PETSC_FALSE) {
65347d480fSBarry Smith       SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong no of nonzeros");
6649b5e25fSSatish Balay     }
6749b5e25fSSatish Balay     ierr = PetscMemzero(c->ilen,c->mbs*sizeof(int));CHKERRQ(ierr);
6849b5e25fSSatish Balay     C = *B;
6949b5e25fSSatish Balay   } else {
7049b5e25fSSatish Balay     ierr = MatCreateSeqSBAIJ(A->comm,bs,nrows*bs,nrows*bs,0,lens,&C);CHKERRQ(ierr);
7149b5e25fSSatish Balay   }
7249b5e25fSSatish Balay   c = (Mat_SeqSBAIJ *)(C->data);
7349b5e25fSSatish Balay   for (i=0; i<nrows; i++) {
7449b5e25fSSatish Balay     row    = irow[i];
7549b5e25fSSatish Balay     kstart = ai[row];
7649b5e25fSSatish Balay     kend   = kstart + a->ilen[row];
7749b5e25fSSatish Balay     mat_i  = c->i[i];
7849b5e25fSSatish Balay     mat_j  = c->j + mat_i;
7949b5e25fSSatish Balay     mat_a  = c->a + mat_i*bs2;
8049b5e25fSSatish Balay     mat_ilen = c->ilen + i;
8149b5e25fSSatish Balay     for (k=kstart; k<kend; k++) {
8249b5e25fSSatish Balay       if ((tcol=ssmap[a->j[k]])) {
8349b5e25fSSatish Balay         *mat_j++ = tcol - 1;
8449b5e25fSSatish Balay         ierr     = PetscMemcpy(mat_a,a->a+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
8549b5e25fSSatish Balay         mat_a   += bs2;
8649b5e25fSSatish Balay         (*mat_ilen)++;
8749b5e25fSSatish Balay       }
8849b5e25fSSatish Balay     }
8949b5e25fSSatish Balay   }
9049b5e25fSSatish Balay 
9149b5e25fSSatish Balay   /* Free work space */
9249b5e25fSSatish Balay   ierr = PetscFree(smap);CHKERRQ(ierr);
9349b5e25fSSatish Balay   ierr = PetscFree(lens);CHKERRQ(ierr);
9449b5e25fSSatish Balay   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
9549b5e25fSSatish Balay   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
9649b5e25fSSatish Balay 
9749b5e25fSSatish Balay   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
9849b5e25fSSatish Balay   *B = C;
9949b5e25fSSatish Balay   PetscFunctionReturn(0);
10049b5e25fSSatish Balay }
10149b5e25fSSatish Balay 
10249b5e25fSSatish Balay #undef __FUNC__
10349b5e25fSSatish Balay #define __FUNC__ "MatGetSubMatrix_SeqSBAIJ"
10449b5e25fSSatish Balay int MatGetSubMatrix_SeqSBAIJ(Mat A,IS isrow,IS iscol,int cs,MatReuse scall,Mat *B)
10549b5e25fSSatish Balay {
10649b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
10749b5e25fSSatish Balay   IS          is1;
10849b5e25fSSatish Balay   int         *vary,*iary,*irow,nrows,i,ierr,bs=a->bs,count;
10949b5e25fSSatish Balay 
11049b5e25fSSatish Balay   PetscFunctionBegin;
111347d480fSBarry Smith   if (isrow != iscol) SETERRA(1,"MatGetSubmatrices_SeqSBAIJ: For symm. format, iscol must equal isro");
11249b5e25fSSatish Balay 
11349b5e25fSSatish Balay   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
11449b5e25fSSatish Balay   ierr = ISGetSize(isrow,&nrows);CHKERRQ(ierr);
11549b5e25fSSatish Balay 
11649b5e25fSSatish Balay   /* Verify if the indices corespond to each element in a block
11749b5e25fSSatish Balay    and form the IS with compressed IS */
11849b5e25fSSatish Balay   vary = (int*)PetscMalloc(2*(a->mbs+1)*sizeof(int));CHKPTRQ(vary);
11949b5e25fSSatish Balay   iary = vary + a->mbs;
12049b5e25fSSatish Balay   ierr = PetscMemzero(vary,(a->mbs)*sizeof(int));CHKERRQ(ierr);
12149b5e25fSSatish Balay   for (i=0; i<nrows; i++) vary[irow[i]/bs]++;
12249b5e25fSSatish Balay 
12349b5e25fSSatish Balay   count = 0;
12449b5e25fSSatish Balay   for (i=0; i<a->mbs; i++) {
125347d480fSBarry Smith     if (vary[i]!=0 && vary[i]!=bs) SETERRA(1,"Index set does not match blocks");
12649b5e25fSSatish Balay     if (vary[i]==bs) iary[count++] = i;
12749b5e25fSSatish Balay   }
12849b5e25fSSatish Balay   ierr = ISCreateGeneral(PETSC_COMM_SELF,count,iary,&is1);CHKERRQ(ierr);
12949b5e25fSSatish Balay 
13049b5e25fSSatish Balay   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
13149b5e25fSSatish Balay   ierr = PetscFree(vary);CHKERRQ(ierr);
13249b5e25fSSatish Balay 
13349b5e25fSSatish Balay   ierr = MatGetSubMatrix_SeqSBAIJ_Private(A,is1,is1,cs,scall,B);CHKERRQ(ierr);
13449b5e25fSSatish Balay   ISDestroy(is1);
13549b5e25fSSatish Balay   PetscFunctionReturn(0);
13649b5e25fSSatish Balay }
13749b5e25fSSatish Balay 
13849b5e25fSSatish Balay #undef __FUNC__
13949b5e25fSSatish Balay #define __FUNC__ "MatGetSubMatrices_SeqSBAIJ"
14049b5e25fSSatish Balay int MatGetSubMatrices_SeqSBAIJ(Mat A,int n,IS *irow,IS *icol,MatReuse scall,Mat **B)
14149b5e25fSSatish Balay {
14249b5e25fSSatish Balay   int ierr,i;
14349b5e25fSSatish Balay 
14449b5e25fSSatish Balay   PetscFunctionBegin;
14549b5e25fSSatish Balay   if (scall == MAT_INITIAL_MATRIX) {
14649b5e25fSSatish Balay     *B = (Mat*)PetscMalloc((n+1)*sizeof(Mat));CHKPTRQ(*B);
14749b5e25fSSatish Balay   }
14849b5e25fSSatish Balay 
14949b5e25fSSatish Balay   for (i=0; i<n; i++) {
15049b5e25fSSatish Balay     ierr = MatGetSubMatrix_SeqSBAIJ(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr);
15149b5e25fSSatish Balay   }
15249b5e25fSSatish Balay   PetscFunctionReturn(0);
15349b5e25fSSatish Balay }
15449b5e25fSSatish Balay 
15549b5e25fSSatish Balay /* -------------------------------------------------------*/
15649b5e25fSSatish Balay /* Should check that shapes of vectors and matrices match */
15749b5e25fSSatish Balay /* -------------------------------------------------------*/
158d9eff348SSatish Balay #include "petscblaslapack.h"
15949b5e25fSSatish Balay 
16049b5e25fSSatish Balay #undef __FUNC__
16149b5e25fSSatish Balay #define __FUNC__ "MatMult_SeqSBAIJ_1"
16249b5e25fSSatish Balay int MatMult_SeqSBAIJ_1(Mat A,Vec xx,Vec zz)
16349b5e25fSSatish Balay {
16449b5e25fSSatish Balay   Mat_SeqSBAIJ    *a = (Mat_SeqSBAIJ*)A->data;
16549b5e25fSSatish Balay   Scalar          *x,*z,*xb,x1,zero=0.0;
16649b5e25fSSatish Balay   MatScalar       *v;
16749b5e25fSSatish Balay   int             mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,ierr,*ib,cval,j;
16849b5e25fSSatish Balay 
16949b5e25fSSatish Balay   PetscFunctionBegin;
17049b5e25fSSatish Balay   ierr = VecSet(&zero,zz);CHKERRQ(ierr);
17149b5e25fSSatish Balay   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
17249b5e25fSSatish Balay   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
17349b5e25fSSatish Balay 
17449b5e25fSSatish Balay   v  = a->a;
17549b5e25fSSatish Balay   xb = x;
17649b5e25fSSatish Balay 
17749b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
17849b5e25fSSatish Balay     n  = ai[1] - ai[0];  /* length of i_th row of A */
17949b5e25fSSatish Balay     x1 = xb[0];
18049b5e25fSSatish Balay     ib = aj + *ai;
181*7fbae186SHong Zhang     if (*ib == i) z[i] += *v++ * x[*ib++];   /* (diag of A)*x */
18249b5e25fSSatish Balay     for (j=1; j<n; j++) {
18349b5e25fSSatish Balay       cval    = *ib;
18449b5e25fSSatish Balay       z[cval] += *v * x1;      /* (strict lower triangular part of A)*x  */
18549b5e25fSSatish Balay       z[i] += *v++ * x[*ib++]; /* (strict upper triangular part of A)*x  */
18649b5e25fSSatish Balay     }
18749b5e25fSSatish Balay     xb++; ai++;
18849b5e25fSSatish Balay   }
18949b5e25fSSatish Balay 
19049b5e25fSSatish Balay   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
19149b5e25fSSatish Balay   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
19249b5e25fSSatish Balay   PLogFlops(2*(a->s_nz*2 - a->m) - a->m);  /* s_nz = (nz+m)/2 */
19349b5e25fSSatish Balay   PetscFunctionReturn(0);
19449b5e25fSSatish Balay }
19549b5e25fSSatish Balay 
19649b5e25fSSatish Balay #undef __FUNC__
19749b5e25fSSatish Balay #define __FUNC__ "MatMult_SeqSBAIJ_2"
19849b5e25fSSatish Balay int MatMult_SeqSBAIJ_2(Mat A,Vec xx,Vec zz)
19949b5e25fSSatish Balay {
20049b5e25fSSatish Balay   Mat_SeqSBAIJ     *a = (Mat_SeqSBAIJ*)A->data;
20149b5e25fSSatish Balay   Scalar          *x,*z,*xb,x1,x2,zero=0.0;
20249b5e25fSSatish Balay   MatScalar       *v;
20349b5e25fSSatish Balay   int             mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,ierr,*ib,cval,j;
20449b5e25fSSatish Balay 
20549b5e25fSSatish Balay 
20649b5e25fSSatish Balay   PetscFunctionBegin;
20749b5e25fSSatish Balay   ierr = VecSet(&zero,zz);CHKERRQ(ierr);
20849b5e25fSSatish Balay   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
20949b5e25fSSatish Balay   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
21049b5e25fSSatish Balay 
21149b5e25fSSatish Balay   v     = a->a;
21249b5e25fSSatish Balay   xb = x;
21349b5e25fSSatish Balay 
21449b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
21549b5e25fSSatish Balay     n  = ai[1] - ai[0]; /* length of i_th block row of A */
21649b5e25fSSatish Balay     x1 = xb[0]; x2 = xb[1];
21749b5e25fSSatish Balay     ib = aj + *ai;
218*7fbae186SHong Zhang     if (*ib == i){     /* (diag of A)*x */
21949b5e25fSSatish Balay       z[2*i]   += v[0]*x1 + v[2]*x2;
22049b5e25fSSatish Balay       z[2*i+1] += v[2]*x1 + v[3]*x2;
22149b5e25fSSatish Balay       v += 4;
222*7fbae186SHong Zhang     }
22349b5e25fSSatish Balay     for (j=1; j<n; j++) {
22449b5e25fSSatish Balay       /* (strict lower triangular part of A)*x  */
22549b5e25fSSatish Balay       cval       = ib[j]*2;
22649b5e25fSSatish Balay       z[cval]     += v[0]*x1 + v[1]*x2;
22749b5e25fSSatish Balay       z[cval+1]   += v[2]*x1 + v[3]*x2;
22849b5e25fSSatish Balay       /* (strict upper triangular part of A)*x  */
22949b5e25fSSatish Balay       z[2*i]   += v[0]*x[cval] + v[2]*x[cval+1];
23049b5e25fSSatish Balay       z[2*i+1] += v[1]*x[cval] + v[3]*x[cval+1];
23149b5e25fSSatish Balay       v  += 4;
23249b5e25fSSatish Balay     }
23349b5e25fSSatish Balay     xb +=2; ai++;
23449b5e25fSSatish Balay   }
23549b5e25fSSatish Balay 
23649b5e25fSSatish Balay   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
23749b5e25fSSatish Balay   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
23849b5e25fSSatish Balay   PLogFlops(8*(a->s_nz*2 - a->m) - a->m);
23949b5e25fSSatish Balay   PetscFunctionReturn(0);
24049b5e25fSSatish Balay }
24149b5e25fSSatish Balay 
24249b5e25fSSatish Balay #undef __FUNC__
24349b5e25fSSatish Balay #define __FUNC__ "MatMult_SeqSBAIJ_3"
24449b5e25fSSatish Balay int MatMult_SeqSBAIJ_3(Mat A,Vec xx,Vec zz)
24549b5e25fSSatish Balay {
24649b5e25fSSatish Balay   Mat_SeqSBAIJ  *a = (Mat_SeqSBAIJ*)A->data;
24749b5e25fSSatish Balay   Scalar          *x,*z,*xb,x1,x2,x3,zero=0.0;
24849b5e25fSSatish Balay   MatScalar       *v;
24949b5e25fSSatish Balay   int             mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,ierr,*ib,cval,j;
25049b5e25fSSatish Balay 
25149b5e25fSSatish Balay 
25249b5e25fSSatish Balay   PetscFunctionBegin;
25349b5e25fSSatish Balay   ierr = VecSet(&zero,zz);CHKERRQ(ierr);
25449b5e25fSSatish Balay   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
25549b5e25fSSatish Balay   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
25649b5e25fSSatish Balay 
25749b5e25fSSatish Balay   v     = a->a;
25849b5e25fSSatish Balay   xb = x;
25949b5e25fSSatish Balay 
26049b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
26149b5e25fSSatish Balay     n  = ai[1] - ai[0]; /* length of i_th block row of A */
26249b5e25fSSatish Balay     x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
26349b5e25fSSatish Balay     ib = aj + *ai;
264*7fbae186SHong Zhang     if (*ib == i){     /* (diag of A)*x */
26549b5e25fSSatish Balay       z[3*i]   += v[0]*x1 + v[3]*x2 + v[6]*x3;
26649b5e25fSSatish Balay       z[3*i+1] += v[3]*x1 + v[4]*x2 + v[7]*x3;
26749b5e25fSSatish Balay       z[3*i+2] += v[6]*x1 + v[7]*x2 + v[8]*x3;
26849b5e25fSSatish Balay       v += 9;
269*7fbae186SHong Zhang     }
27049b5e25fSSatish Balay     for (j=1; j<n; j++) {
27149b5e25fSSatish Balay       /* (strict lower triangular part of A)*x  */
27249b5e25fSSatish Balay       cval       = ib[j]*3;
27349b5e25fSSatish Balay       z[cval]     += v[0]*x1 + v[1]*x2 + v[2]*x3;
27449b5e25fSSatish Balay       z[cval+1]   += v[3]*x1 + v[4]*x2 + v[5]*x3;
27549b5e25fSSatish Balay       z[cval+2]   += v[6]*x1 + v[7]*x2 + v[8]*x3;
27649b5e25fSSatish Balay       /* (strict upper triangular part of A)*x  */
27749b5e25fSSatish Balay       z[3*i]   += v[0]*x[cval] + v[3]*x[cval+1]+ v[6]*x[cval+2];
27849b5e25fSSatish Balay       z[3*i+1] += v[1]*x[cval] + v[4]*x[cval+1]+ v[7]*x[cval+2];
27949b5e25fSSatish Balay       z[3*i+2] += v[2]*x[cval] + v[5]*x[cval+1]+ v[8]*x[cval+2];
28049b5e25fSSatish Balay       v  += 9;
28149b5e25fSSatish Balay     }
28249b5e25fSSatish Balay     xb +=3; ai++;
28349b5e25fSSatish Balay   }
28449b5e25fSSatish Balay 
28549b5e25fSSatish Balay   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
28649b5e25fSSatish Balay   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
28749b5e25fSSatish Balay   PLogFlops(18*(a->s_nz*2 - a->m) - a->m);
28849b5e25fSSatish Balay   PetscFunctionReturn(0);
28949b5e25fSSatish Balay }
29049b5e25fSSatish Balay 
29149b5e25fSSatish Balay #undef __FUNC__
29249b5e25fSSatish Balay #define __FUNC__ "MatMult_SeqSBAIJ_4"
29349b5e25fSSatish Balay int MatMult_SeqSBAIJ_4(Mat A,Vec xx,Vec zz)
29449b5e25fSSatish Balay {
29549b5e25fSSatish Balay   Mat_SeqSBAIJ     *a = (Mat_SeqSBAIJ*)A->data;
29649b5e25fSSatish Balay   Scalar          *x,*z,*xb,x1,x2,x3,x4,zero=0.0;
29749b5e25fSSatish Balay   MatScalar       *v;
29849b5e25fSSatish Balay   int             mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,ierr,*ib,cval,j;
29949b5e25fSSatish Balay 
30049b5e25fSSatish Balay   PetscFunctionBegin;
30149b5e25fSSatish Balay   ierr = VecSet(&zero,zz);CHKERRQ(ierr);
30249b5e25fSSatish Balay   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
30349b5e25fSSatish Balay   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
30449b5e25fSSatish Balay 
30549b5e25fSSatish Balay   v     = a->a;
30649b5e25fSSatish Balay   xb = x;
30749b5e25fSSatish Balay 
30849b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
30949b5e25fSSatish Balay     n  = ai[1] - ai[0]; /* length of i_th block row of A */
31049b5e25fSSatish Balay     x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
31149b5e25fSSatish Balay     ib = aj + *ai;
312*7fbae186SHong Zhang     if (*ib == i){     /* (diag of A)*x */
31349b5e25fSSatish Balay       z[4*i]   += v[0]*x1 + v[4]*x2 +  v[8]*x3 + v[12]*x4;
31449b5e25fSSatish Balay       z[4*i+1] += v[4]*x1 + v[5]*x2 +  v[9]*x3 + v[13]*x4;
31549b5e25fSSatish Balay       z[4*i+2] += v[8]*x1 + v[9]*x2 + v[10]*x3 + v[14]*x4;
31649b5e25fSSatish Balay       z[4*i+3] += v[12]*x1+ v[13]*x2+ v[14]*x3 + v[15]*x4;
31749b5e25fSSatish Balay       v += 16;
318*7fbae186SHong Zhang     }
31949b5e25fSSatish Balay     for (j=1; j<n; j++) {
32049b5e25fSSatish Balay       /* (strict lower triangular part of A)*x  */
32149b5e25fSSatish Balay       cval       = ib[j]*4;
32249b5e25fSSatish Balay       z[cval]     += v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4;
32349b5e25fSSatish Balay       z[cval+1]   += v[4]*x1 + v[5]*x2 + v[6]*x3 + v[7]*x4;
32449b5e25fSSatish Balay       z[cval+2]   += v[8]*x1 + v[9]*x2 + v[10]*x3 + v[11]*x4;
32549b5e25fSSatish Balay       z[cval+3]   += v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4;
32649b5e25fSSatish Balay       /* (strict upper triangular part of A)*x  */
32749b5e25fSSatish Balay       z[4*i]   += v[0]*x[cval] + v[4]*x[cval+1]+ v[8]*x[cval+2] + v[12]*x[cval+3];
32849b5e25fSSatish Balay       z[4*i+1] += v[1]*x[cval] + v[5]*x[cval+1]+ v[9]*x[cval+2] + v[13]*x[cval+3];
32949b5e25fSSatish Balay       z[4*i+2] += v[2]*x[cval] + v[6]*x[cval+1]+ v[10]*x[cval+2]+ v[14]*x[cval+3];
33049b5e25fSSatish Balay       z[4*i+3] += v[3]*x[cval] + v[7]*x[cval+1]+ v[11]*x[cval+2]+ v[15]*x[cval+3];
33149b5e25fSSatish Balay       v  += 16;
33249b5e25fSSatish Balay     }
33349b5e25fSSatish Balay     xb +=4; ai++;
33449b5e25fSSatish Balay   }
33549b5e25fSSatish Balay 
33649b5e25fSSatish Balay   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
33749b5e25fSSatish Balay   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
33849b5e25fSSatish Balay   PLogFlops(32*(a->s_nz*2 - a->m) - a->m);
33949b5e25fSSatish Balay   PetscFunctionReturn(0);
34049b5e25fSSatish Balay }
34149b5e25fSSatish Balay 
34249b5e25fSSatish Balay #undef __FUNC__
34349b5e25fSSatish Balay #define __FUNC__ "MatMult_SeqSBAIJ_5"
34449b5e25fSSatish Balay int MatMult_SeqSBAIJ_5(Mat A,Vec xx,Vec zz)
34549b5e25fSSatish Balay {
34649b5e25fSSatish Balay   Mat_SeqSBAIJ     *a = (Mat_SeqSBAIJ*)A->data;
34749b5e25fSSatish Balay   Scalar          *x,*z,*xb,x1,x2,x3,x4,x5,zero=0.0;
34849b5e25fSSatish Balay   MatScalar       *v;
34949b5e25fSSatish Balay   int             mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,ierr,*ib,cval,j;
35049b5e25fSSatish Balay 
35149b5e25fSSatish Balay   PetscFunctionBegin;
35249b5e25fSSatish Balay   ierr = VecSet(&zero,zz);CHKERRQ(ierr);
35349b5e25fSSatish Balay   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
35449b5e25fSSatish Balay   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
35549b5e25fSSatish Balay 
35649b5e25fSSatish Balay   v     = a->a;
35749b5e25fSSatish Balay   xb = x;
35849b5e25fSSatish Balay 
35949b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
36049b5e25fSSatish Balay     n  = ai[1] - ai[0]; /* length of i_th block row of A */
36149b5e25fSSatish Balay     x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5=xb[4];
36249b5e25fSSatish Balay     ib = aj + *ai;
363*7fbae186SHong Zhang     if (*ib == i){      /* (diag of A)*x */
36449b5e25fSSatish Balay       z[5*i]   += v[0]*x1  + v[5]*x2 + v[10]*x3 + v[15]*x4+ v[20]*x5;
36549b5e25fSSatish Balay       z[5*i+1] += v[5]*x1  + v[6]*x2 + v[11]*x3 + v[16]*x4+ v[21]*x5;
36649b5e25fSSatish Balay       z[5*i+2] += v[10]*x1 +v[11]*x2 + v[12]*x3 + v[17]*x4+ v[22]*x5;
36749b5e25fSSatish Balay       z[5*i+3] += v[15]*x1 +v[16]*x2 + v[17]*x3 + v[18]*x4+ v[23]*x5;
36849b5e25fSSatish Balay       z[5*i+4] += v[20]*x1 +v[21]*x2 + v[22]*x3 + v[23]*x4+ v[24]*x5;
36949b5e25fSSatish Balay       v += 25;
370*7fbae186SHong Zhang     }
37149b5e25fSSatish Balay     for (j=1; j<n; j++) {
37249b5e25fSSatish Balay       /* (strict lower triangular part of A)*x  */
37349b5e25fSSatish Balay       cval       = ib[j]*5;
37449b5e25fSSatish Balay       z[cval]     += v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4 + v[4]*x5;
37549b5e25fSSatish Balay       z[cval+1]   += v[5]*x1 + v[6]*x2 + v[7]*x3 + v[8]*x4 + v[9]*x5;
37649b5e25fSSatish Balay       z[cval+2]   += v[10]*x1 + v[11]*x2 + v[12]*x3 + v[13]*x4+ v[14]*x5;
37749b5e25fSSatish Balay       z[cval+3]   += v[15]*x1 + v[16]*x2 + v[17]*x3 + v[18]*x4+ v[19]*x5;
37849b5e25fSSatish Balay       z[cval+4]   += v[20]*x1 + v[21]*x2 + v[22]*x3 + v[23]*x4+ v[24]*x5;
37949b5e25fSSatish Balay       /* (strict upper triangular part of A)*x  */
38049b5e25fSSatish 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];
38149b5e25fSSatish 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];
38249b5e25fSSatish 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];
38349b5e25fSSatish 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];
38449b5e25fSSatish 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];
38549b5e25fSSatish Balay       v  += 25;
38649b5e25fSSatish Balay     }
38749b5e25fSSatish Balay     xb +=5; ai++;
38849b5e25fSSatish Balay   }
38949b5e25fSSatish Balay 
39049b5e25fSSatish Balay   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
39149b5e25fSSatish Balay   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
39249b5e25fSSatish Balay   PLogFlops(50*(a->s_nz*2 - a->m) - a->m);
39349b5e25fSSatish Balay   PetscFunctionReturn(0);
39449b5e25fSSatish Balay }
39549b5e25fSSatish Balay 
39649b5e25fSSatish Balay 
39749b5e25fSSatish Balay #undef __FUNC__
39849b5e25fSSatish Balay #define __FUNC__ "MatMult_SeqSBAIJ_6"
39949b5e25fSSatish Balay int MatMult_SeqSBAIJ_6(Mat A,Vec xx,Vec zz)
40049b5e25fSSatish Balay {
40149b5e25fSSatish Balay   Mat_SeqSBAIJ     *a = (Mat_SeqSBAIJ*)A->data;
40249b5e25fSSatish Balay   Scalar          *x,*z,*xb,x1,x2,x3,x4,x5,x6,zero=0.0;
40349b5e25fSSatish Balay   MatScalar       *v;
40449b5e25fSSatish Balay   int             mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,ierr,*ib,cval,j;
40549b5e25fSSatish Balay 
40649b5e25fSSatish Balay   PetscFunctionBegin;
40749b5e25fSSatish Balay   ierr = VecSet(&zero,zz);CHKERRQ(ierr);
40849b5e25fSSatish Balay   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
40949b5e25fSSatish Balay   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
41049b5e25fSSatish Balay 
41149b5e25fSSatish Balay   v     = a->a;
41249b5e25fSSatish Balay   xb = x;
41349b5e25fSSatish Balay 
41449b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
41549b5e25fSSatish Balay     n  = ai[1] - ai[0]; /* length of i_th block row of A */
41649b5e25fSSatish Balay     x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5=xb[4]; x6=xb[5];
41749b5e25fSSatish Balay     ib = aj + *ai;
418*7fbae186SHong Zhang     if (*ib == i){      /* (diag of A)*x */
41949b5e25fSSatish Balay       z[6*i]   += v[0]*x1  + v[6]*x2 + v[12]*x3 + v[18]*x4+ v[24]*x5 + v[30]*x6;
42049b5e25fSSatish Balay       z[6*i+1] += v[6]*x1  + v[7]*x2 + v[13]*x3 + v[19]*x4+ v[25]*x5 + v[31]*x6;
42149b5e25fSSatish Balay       z[6*i+2] += v[12]*x1 +v[13]*x2 + v[14]*x3 + v[20]*x4+ v[26]*x5 + v[32]*x6;
42249b5e25fSSatish Balay       z[6*i+3] += v[18]*x1 +v[19]*x2 + v[20]*x3 + v[21]*x4+ v[27]*x5 + v[33]*x6;
42349b5e25fSSatish Balay       z[6*i+4] += v[24]*x1 +v[25]*x2 + v[26]*x3 + v[27]*x4+ v[28]*x5 + v[34]*x6;
42449b5e25fSSatish Balay       z[6*i+5] += v[30]*x1 +v[31]*x2 + v[32]*x3 + v[33]*x4+ v[34]*x5 + v[35]*x6;
42549b5e25fSSatish Balay       v += 36;
426*7fbae186SHong Zhang     }
42749b5e25fSSatish Balay     for (j=1; j<n; j++) {
42849b5e25fSSatish Balay       /* (strict lower triangular part of A)*x  */
42949b5e25fSSatish Balay       cval       = ib[j]*6;
43049b5e25fSSatish Balay       z[cval]   += v[0]*x1  + v[1]*x2 + v[2]*x3 + v[3]*x4+ v[4]*x5 + v[5]*x6;
43149b5e25fSSatish Balay       z[cval+1] += v[6]*x1  + v[7]*x2 + v[8]*x3 + v[9]*x4+ v[10]*x5 + v[11]*x6;
43249b5e25fSSatish Balay       z[cval+2] += v[12]*x1  + v[13]*x2 + v[14]*x3 + v[15]*x4+ v[16]*x5 + v[17]*x6;
43349b5e25fSSatish Balay       z[cval+3] += v[18]*x1  + v[19]*x2 + v[20]*x3 + v[21]*x4+ v[22]*x5 + v[23]*x6;
43449b5e25fSSatish Balay       z[cval+4] += v[24]*x1  + v[25]*x2 + v[26]*x3 + v[27]*x4+ v[28]*x5 + v[29]*x6;
43549b5e25fSSatish Balay       z[cval+5] += v[30]*x1  + v[31]*x2 + v[32]*x3 + v[33]*x4+ v[34]*x5 + v[35]*x6;
43649b5e25fSSatish Balay       /* (strict upper triangular part of A)*x  */
43749b5e25fSSatish 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];
43849b5e25fSSatish 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];
43949b5e25fSSatish 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];
44049b5e25fSSatish 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];
44149b5e25fSSatish 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];
44249b5e25fSSatish 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];
44349b5e25fSSatish Balay       v  += 36;
44449b5e25fSSatish Balay     }
44549b5e25fSSatish Balay     xb +=6; ai++;
44649b5e25fSSatish Balay   }
44749b5e25fSSatish Balay 
44849b5e25fSSatish Balay   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
44949b5e25fSSatish Balay   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
45049b5e25fSSatish Balay   PLogFlops(72*(a->s_nz*2 - a->m) - a->m);
45149b5e25fSSatish Balay   PetscFunctionReturn(0);
45249b5e25fSSatish Balay }
45349b5e25fSSatish Balay #undef __FUNC__
45449b5e25fSSatish Balay #define __FUNC__ "MatMult_SeqSBAIJ_7"
45549b5e25fSSatish Balay int MatMult_SeqSBAIJ_7(Mat A,Vec xx,Vec zz)
45649b5e25fSSatish Balay {
45749b5e25fSSatish Balay   Mat_SeqSBAIJ     *a = (Mat_SeqSBAIJ*)A->data;
45849b5e25fSSatish Balay   Scalar          *x,*z,*xb,x1,x2,x3,x4,x5,x6,x7,zero=0.0;
45949b5e25fSSatish Balay   MatScalar       *v;
46049b5e25fSSatish Balay   int             mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,ierr,*ib,cval,j;
46149b5e25fSSatish Balay 
46249b5e25fSSatish Balay   PetscFunctionBegin;
46349b5e25fSSatish Balay   ierr = VecSet(&zero,zz);CHKERRQ(ierr);
46449b5e25fSSatish Balay   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
46549b5e25fSSatish Balay   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
46649b5e25fSSatish Balay 
46749b5e25fSSatish Balay   v     = a->a;
46849b5e25fSSatish Balay   xb = x;
46949b5e25fSSatish Balay 
47049b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
47149b5e25fSSatish Balay     n  = ai[1] - ai[0]; /* length of i_th block row of A */
47249b5e25fSSatish Balay     x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5=xb[4]; x6=xb[5]; x7=xb[6];
47349b5e25fSSatish Balay     ib = aj + *ai;
474*7fbae186SHong Zhang     if (*ib == i){      /* (diag of A)*x */
47549b5e25fSSatish 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;
47649b5e25fSSatish 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;
47749b5e25fSSatish 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;
47849b5e25fSSatish 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;
47949b5e25fSSatish 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;
48049b5e25fSSatish 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;
48149b5e25fSSatish 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;
48249b5e25fSSatish Balay       v += 49;
483*7fbae186SHong Zhang     }
48449b5e25fSSatish Balay     for (j=1; j<n; j++) {
48549b5e25fSSatish Balay       /* (strict lower triangular part of A)*x  */
48649b5e25fSSatish Balay       cval       = ib[j]*7;
48749b5e25fSSatish Balay       z[cval]   += v[0]*x1  + v[1]*x2 + v[2]*x3 + v[3]*x4+ v[4]*x5 + v[5]*x6+ v[6]*x7;
48849b5e25fSSatish 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;
48949b5e25fSSatish 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;
49049b5e25fSSatish 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;
49149b5e25fSSatish 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;
49249b5e25fSSatish 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;
49349b5e25fSSatish 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;
49449b5e25fSSatish Balay       /* (strict upper triangular part of A)*x  */
49549b5e25fSSatish 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];
49649b5e25fSSatish 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];
49749b5e25fSSatish 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];
49849b5e25fSSatish 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];
49949b5e25fSSatish 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];
50049b5e25fSSatish 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];
50149b5e25fSSatish 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];
50249b5e25fSSatish Balay       v  += 49;
50349b5e25fSSatish Balay     }
50449b5e25fSSatish Balay     xb +=7; ai++;
50549b5e25fSSatish Balay   }
50649b5e25fSSatish Balay   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
50749b5e25fSSatish Balay   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
50849b5e25fSSatish Balay   PLogFlops(98*(a->s_nz*2 - a->m) - a->m);
50949b5e25fSSatish Balay   PetscFunctionReturn(0);
51049b5e25fSSatish Balay }
51149b5e25fSSatish Balay 
51249b5e25fSSatish Balay /*
51349b5e25fSSatish Balay     This will not work with MatScalar == float because it calls the BLAS
51449b5e25fSSatish Balay */
51549b5e25fSSatish Balay #undef __FUNC__
51649b5e25fSSatish Balay #define __FUNC__ "MatMult_SeqSBAIJ_N"
51749b5e25fSSatish Balay int MatMult_SeqSBAIJ_N(Mat A,Vec xx,Vec zz)
51849b5e25fSSatish Balay {
51949b5e25fSSatish Balay   Mat_SeqSBAIJ     *a = (Mat_SeqSBAIJ*)A->data;
520066653e3SSatish Balay   Scalar          *x,*x_ptr,*z,*z_ptr,*xb,*zb,*work,*workt,zero=0.0;
521066653e3SSatish Balay   MatScalar       *v;
522066653e3SSatish Balay   int             ierr,mbs=a->mbs,i,*idx,*aj,*ii,bs=a->bs,j,n,bs2=a->bs2;
523066653e3SSatish Balay   int             ncols,k;
52449b5e25fSSatish Balay 
52549b5e25fSSatish Balay   PetscFunctionBegin;
52649b5e25fSSatish Balay   ierr = VecSet(&zero,zz);CHKERRQ(ierr);
52749b5e25fSSatish Balay   ierr = VecGetArray(xx,&x);CHKERRQ(ierr); x_ptr=x;
52849b5e25fSSatish Balay   ierr = VecGetArray(zz,&z);CHKERRQ(ierr); z_ptr=z;
52949b5e25fSSatish Balay 
53049b5e25fSSatish Balay   aj   = a->j;
53149b5e25fSSatish Balay   v    = a->a;
53249b5e25fSSatish Balay   ii   = a->i;
53349b5e25fSSatish Balay 
53449b5e25fSSatish Balay   if (!a->mult_work) {
53549b5e25fSSatish Balay     k = a->m;
53649b5e25fSSatish Balay     a->mult_work = (Scalar*)PetscMalloc((k+1)*sizeof(Scalar));CHKPTRQ(a->mult_work);
53749b5e25fSSatish Balay   }
53849b5e25fSSatish Balay   work = a->mult_work;
53949b5e25fSSatish Balay 
54049b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
54149b5e25fSSatish Balay     n     = ii[1] - ii[0]; ncols = n*bs;
54249b5e25fSSatish Balay     workt = work; idx=aj+ii[0];
54349b5e25fSSatish Balay 
54449b5e25fSSatish Balay     /* upper triangular part */
54549b5e25fSSatish Balay     for (j=0; j<n; j++) {
546066653e3SSatish Balay       /* col = bs*(*idx); */
54749b5e25fSSatish Balay       xb = x_ptr + bs*(*idx++);
54849b5e25fSSatish Balay       for (k=0; k<bs; k++) workt[k] = xb[k];
54949b5e25fSSatish Balay       workt += bs;
55049b5e25fSSatish Balay     }
55149b5e25fSSatish Balay     /* z(i*bs:(i+1)*bs-1) += A(i,:)*x */
55249b5e25fSSatish Balay     Kernel_w_gets_w_plus_Ar_times_v(bs,ncols,work,v,z);
55349b5e25fSSatish Balay 
55449b5e25fSSatish Balay     /* strict lower triangular part */
55549b5e25fSSatish Balay     ncols -= bs;
55649b5e25fSSatish Balay     if (ncols > 0){
55749b5e25fSSatish Balay       workt = work;
55849b5e25fSSatish Balay       ierr  = PetscMemzero(workt,ncols*sizeof(Scalar));CHKERRQ(ierr);
55949b5e25fSSatish Balay       Kernel_w_gets_w_plus_trans_Ar_times_v(bs,ncols,x,v+bs2,workt);
56049b5e25fSSatish Balay 
56149b5e25fSSatish Balay       idx=aj+ii[0]+1;
56249b5e25fSSatish Balay       for (j=1; j<n; j++) {
56349b5e25fSSatish Balay         zb = z_ptr + bs*(*idx);
56449b5e25fSSatish Balay         idx++;
56549b5e25fSSatish Balay         for (k=0; k<bs; k++) zb[k] += workt[k] ;
56649b5e25fSSatish Balay         workt += bs;
56749b5e25fSSatish Balay       }
56849b5e25fSSatish Balay     }
56949b5e25fSSatish Balay 
57049b5e25fSSatish Balay     x += bs; v += n*bs2; z += bs; ii++;
57149b5e25fSSatish Balay   }
57249b5e25fSSatish Balay 
57349b5e25fSSatish Balay   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
57449b5e25fSSatish Balay   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
57549b5e25fSSatish Balay   PLogFlops(2*(a->s_nz*2 - a->m)*bs2 - a->m);
57649b5e25fSSatish Balay   PetscFunctionReturn(0);
57749b5e25fSSatish Balay }
57849b5e25fSSatish Balay 
57949b5e25fSSatish Balay #undef __FUNC__
58049b5e25fSSatish Balay #define __FUNC__ "MatMultAdd_SeqSBAIJ_1"
58149b5e25fSSatish Balay int MatMultAdd_SeqSBAIJ_1(Mat A,Vec xx,Vec yy,Vec zz)
58249b5e25fSSatish Balay {
58349b5e25fSSatish Balay   Mat_SeqSBAIJ    *a = (Mat_SeqSBAIJ*)A->data;
58449b5e25fSSatish Balay   Scalar          *x,*y,*z,*xb,x1;
58549b5e25fSSatish Balay   MatScalar       *v;
58649b5e25fSSatish Balay   int             mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,ierr,*ib,cval,j;
58749b5e25fSSatish Balay 
58849b5e25fSSatish Balay   PetscFunctionBegin;
58949b5e25fSSatish Balay   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
59049b5e25fSSatish Balay   if (yy != xx) {
59149b5e25fSSatish Balay     ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
59249b5e25fSSatish Balay   } else {
59349b5e25fSSatish Balay     y = x;
59449b5e25fSSatish Balay   }
59549b5e25fSSatish Balay   if (zz != yy) {
59649b5e25fSSatish Balay     ierr = VecCopy(yy,zz);CHKERRQ(ierr);
59749b5e25fSSatish Balay     ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
59849b5e25fSSatish Balay   } else {
59949b5e25fSSatish Balay     z = y;
60049b5e25fSSatish Balay   }
60149b5e25fSSatish Balay 
60249b5e25fSSatish Balay   v  = a->a;
60349b5e25fSSatish Balay   xb = x;
60449b5e25fSSatish Balay 
60549b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
60649b5e25fSSatish Balay     n  = ai[1] - ai[0];  /* length of i_th row of A */
60749b5e25fSSatish Balay     x1 = xb[0];
60849b5e25fSSatish Balay     ib = aj + *ai;
609*7fbae186SHong Zhang     if (*ib == i) z[i] += *v++ * x[*ib++];   /* (diag of A)*x */
61049b5e25fSSatish Balay     for (j=1; j<n; j++) {
61149b5e25fSSatish Balay       cval    = *ib;
61249b5e25fSSatish Balay       z[cval] += *v * x1;      /* (strict lower triangular part of A)*x  */
61349b5e25fSSatish Balay       z[i] += *v++ * x[*ib++]; /* (strict upper triangular part of A)*x  */
61449b5e25fSSatish Balay     }
61549b5e25fSSatish Balay     xb++; ai++;
61649b5e25fSSatish Balay   }
61749b5e25fSSatish Balay 
61849b5e25fSSatish Balay   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
61949b5e25fSSatish Balay   if (yy != xx) ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
62049b5e25fSSatish Balay   if (zz != yy) ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
62149b5e25fSSatish Balay 
62249b5e25fSSatish Balay   PLogFlops(2*(a->s_nz*2 - a->m));
62349b5e25fSSatish Balay   PetscFunctionReturn(0);
62449b5e25fSSatish Balay }
62549b5e25fSSatish Balay 
62649b5e25fSSatish Balay #undef __FUNC__
62749b5e25fSSatish Balay #define __FUNC__ "MatMultAdd_SeqSBAIJ_2"
62849b5e25fSSatish Balay int MatMultAdd_SeqSBAIJ_2(Mat A,Vec xx,Vec yy,Vec zz)
62949b5e25fSSatish Balay {
63049b5e25fSSatish Balay   Mat_SeqSBAIJ     *a = (Mat_SeqSBAIJ*)A->data;
63149b5e25fSSatish Balay   Scalar          *x,*y,*z,*xb,x1,x2;
63249b5e25fSSatish Balay   MatScalar       *v;
63349b5e25fSSatish Balay   int             mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,ierr,*ib,cval,j;
63449b5e25fSSatish Balay 
63549b5e25fSSatish Balay   PetscFunctionBegin;
63649b5e25fSSatish Balay   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
63749b5e25fSSatish Balay   if (yy != xx) {
63849b5e25fSSatish Balay     ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
63949b5e25fSSatish Balay   } else {
64049b5e25fSSatish Balay     y = x;
64149b5e25fSSatish Balay   }
64249b5e25fSSatish Balay   if (zz != yy) {
64349b5e25fSSatish Balay     ierr = VecCopy(yy,zz);CHKERRQ(ierr);
64449b5e25fSSatish Balay     ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
64549b5e25fSSatish Balay   } else {
64649b5e25fSSatish Balay     z = y;
64749b5e25fSSatish Balay   }
64849b5e25fSSatish Balay 
64949b5e25fSSatish Balay   v     = a->a;
65049b5e25fSSatish Balay   xb = x;
65149b5e25fSSatish Balay 
65249b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
65349b5e25fSSatish Balay     n  = ai[1] - ai[0]; /* length of i_th block row of A */
65449b5e25fSSatish Balay     x1 = xb[0]; x2 = xb[1];
65549b5e25fSSatish Balay     ib = aj + *ai;
656*7fbae186SHong Zhang     if (*ib == i){      /* (diag of A)*x */
65749b5e25fSSatish Balay       z[2*i]   += v[0]*x1 + v[2]*x2;
65849b5e25fSSatish Balay       z[2*i+1] += v[2]*x1 + v[3]*x2;
65949b5e25fSSatish Balay       v += 4;
660*7fbae186SHong Zhang     }
66149b5e25fSSatish Balay     for (j=1; j<n; j++) {
66249b5e25fSSatish Balay       /* (strict lower triangular part of A)*x  */
66349b5e25fSSatish Balay       cval       = ib[j]*2;
66449b5e25fSSatish Balay       z[cval]     += v[0]*x1 + v[1]*x2;
66549b5e25fSSatish Balay       z[cval+1]   += v[2]*x1 + v[3]*x2;
66649b5e25fSSatish Balay       /* (strict upper triangular part of A)*x  */
66749b5e25fSSatish Balay       z[2*i]   += v[0]*x[cval] + v[2]*x[cval+1];
66849b5e25fSSatish Balay       z[2*i+1] += v[1]*x[cval] + v[3]*x[cval+1];
66949b5e25fSSatish Balay       v  += 4;
67049b5e25fSSatish Balay     }
67149b5e25fSSatish Balay     xb +=2; ai++;
67249b5e25fSSatish Balay   }
67349b5e25fSSatish Balay 
67449b5e25fSSatish Balay   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
67549b5e25fSSatish Balay   if (yy != xx) ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
67649b5e25fSSatish Balay   if (zz != yy) ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
67749b5e25fSSatish Balay 
67849b5e25fSSatish Balay   PLogFlops(4*(a->s_nz*2 - a->m));
67949b5e25fSSatish Balay   PetscFunctionReturn(0);
68049b5e25fSSatish Balay }
68149b5e25fSSatish Balay 
68249b5e25fSSatish Balay #undef __FUNC__
68349b5e25fSSatish Balay #define __FUNC__ "MatMultAdd_SeqSBAIJ_3"
68449b5e25fSSatish Balay int MatMultAdd_SeqSBAIJ_3(Mat A,Vec xx,Vec yy,Vec zz)
68549b5e25fSSatish Balay {
68649b5e25fSSatish Balay   Mat_SeqSBAIJ  *a = (Mat_SeqSBAIJ*)A->data;
68749b5e25fSSatish Balay   Scalar          *x,*y,*z,*xb,x1,x2,x3;
68849b5e25fSSatish Balay   MatScalar       *v;
68949b5e25fSSatish Balay   int             mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,ierr,*ib,cval,j;
69049b5e25fSSatish Balay 
69149b5e25fSSatish Balay   PetscFunctionBegin;
69249b5e25fSSatish Balay   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
69349b5e25fSSatish Balay   if (yy != xx) {
69449b5e25fSSatish Balay     ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
69549b5e25fSSatish Balay   } else {
69649b5e25fSSatish Balay     y = x;
69749b5e25fSSatish Balay   }
69849b5e25fSSatish Balay   if (zz != yy) {
69949b5e25fSSatish Balay     ierr = VecCopy(yy,zz);CHKERRQ(ierr);
70049b5e25fSSatish Balay     ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
70149b5e25fSSatish Balay   } else {
70249b5e25fSSatish Balay     z = y;
70349b5e25fSSatish Balay   }
70449b5e25fSSatish Balay 
70549b5e25fSSatish Balay   v     = a->a;
70649b5e25fSSatish Balay   xb = x;
70749b5e25fSSatish Balay 
70849b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
70949b5e25fSSatish Balay     n  = ai[1] - ai[0]; /* length of i_th block row of A */
71049b5e25fSSatish Balay     x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
71149b5e25fSSatish Balay     ib = aj + *ai;
712*7fbae186SHong Zhang     if (*ib == i){     /* (diag of A)*x */
71349b5e25fSSatish Balay      z[3*i]   += v[0]*x1 + v[3]*x2 + v[6]*x3;
71449b5e25fSSatish Balay      z[3*i+1] += v[3]*x1 + v[4]*x2 + v[7]*x3;
71549b5e25fSSatish Balay      z[3*i+2] += v[6]*x1 + v[7]*x2 + v[8]*x3;
71649b5e25fSSatish Balay      v += 9;
717*7fbae186SHong Zhang     }
71849b5e25fSSatish Balay     for (j=1; j<n; j++) {
71949b5e25fSSatish Balay       /* (strict lower triangular part of A)*x  */
72049b5e25fSSatish Balay       cval       = ib[j]*3;
72149b5e25fSSatish Balay       z[cval]     += v[0]*x1 + v[1]*x2 + v[2]*x3;
72249b5e25fSSatish Balay       z[cval+1]   += v[3]*x1 + v[4]*x2 + v[5]*x3;
72349b5e25fSSatish Balay       z[cval+2]   += v[6]*x1 + v[7]*x2 + v[8]*x3;
72449b5e25fSSatish Balay       /* (strict upper triangular part of A)*x  */
72549b5e25fSSatish Balay       z[3*i]   += v[0]*x[cval] + v[3]*x[cval+1]+ v[6]*x[cval+2];
72649b5e25fSSatish Balay       z[3*i+1] += v[1]*x[cval] + v[4]*x[cval+1]+ v[7]*x[cval+2];
72749b5e25fSSatish Balay       z[3*i+2] += v[2]*x[cval] + v[5]*x[cval+1]+ v[8]*x[cval+2];
72849b5e25fSSatish Balay       v  += 9;
72949b5e25fSSatish Balay     }
73049b5e25fSSatish Balay     xb +=3; ai++;
73149b5e25fSSatish Balay   }
73249b5e25fSSatish Balay 
73349b5e25fSSatish Balay   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
73449b5e25fSSatish Balay   if (yy != xx) ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
73549b5e25fSSatish Balay   if (zz != yy) ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
73649b5e25fSSatish Balay 
73749b5e25fSSatish Balay   PLogFlops(18*(a->s_nz*2 - a->m));
73849b5e25fSSatish Balay   PetscFunctionReturn(0);
73949b5e25fSSatish Balay }
74049b5e25fSSatish Balay 
74149b5e25fSSatish Balay #undef __FUNC__
74249b5e25fSSatish Balay #define __FUNC__ "MatMultAdd_SeqSBAIJ_4"
74349b5e25fSSatish Balay int MatMultAdd_SeqSBAIJ_4(Mat A,Vec xx,Vec yy,Vec zz)
74449b5e25fSSatish Balay {
74549b5e25fSSatish Balay   Mat_SeqSBAIJ     *a = (Mat_SeqSBAIJ*)A->data;
74649b5e25fSSatish Balay   Scalar          *x,*y,*z,*xb,x1,x2,x3,x4;
74749b5e25fSSatish Balay   MatScalar       *v;
74849b5e25fSSatish Balay   int             mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,ierr,*ib,cval,j;
74949b5e25fSSatish Balay 
75049b5e25fSSatish Balay   PetscFunctionBegin;
75149b5e25fSSatish Balay   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
75249b5e25fSSatish Balay   if (yy != xx) {
75349b5e25fSSatish Balay     ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
75449b5e25fSSatish Balay   } else {
75549b5e25fSSatish Balay     y = x;
75649b5e25fSSatish Balay   }
75749b5e25fSSatish Balay   if (zz != yy) {
75849b5e25fSSatish Balay     ierr = VecCopy(yy,zz);CHKERRQ(ierr);
75949b5e25fSSatish Balay     ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
76049b5e25fSSatish Balay   } else {
76149b5e25fSSatish Balay     z = y;
76249b5e25fSSatish Balay   }
76349b5e25fSSatish Balay 
76449b5e25fSSatish Balay   v     = a->a;
76549b5e25fSSatish Balay   xb = x;
76649b5e25fSSatish Balay 
76749b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
76849b5e25fSSatish Balay     n  = ai[1] - ai[0]; /* length of i_th block row of A */
76949b5e25fSSatish Balay     x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
77049b5e25fSSatish Balay     ib = aj + *ai;
771*7fbae186SHong Zhang     if (*ib == i){      /* (diag of A)*x */
77249b5e25fSSatish Balay       z[4*i]   += v[0]*x1 + v[4]*x2 +  v[8]*x3 + v[12]*x4;
77349b5e25fSSatish Balay       z[4*i+1] += v[4]*x1 + v[5]*x2 +  v[9]*x3 + v[13]*x4;
77449b5e25fSSatish Balay       z[4*i+2] += v[8]*x1 + v[9]*x2 + v[10]*x3 + v[14]*x4;
77549b5e25fSSatish Balay       z[4*i+3] += v[12]*x1+ v[13]*x2+ v[14]*x3 + v[15]*x4;
77649b5e25fSSatish Balay       v += 16;
777*7fbae186SHong Zhang     }
77849b5e25fSSatish Balay     for (j=1; j<n; j++) {
77949b5e25fSSatish Balay       /* (strict lower triangular part of A)*x  */
78049b5e25fSSatish Balay       cval       = ib[j]*4;
78149b5e25fSSatish Balay       z[cval]     += v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4;
78249b5e25fSSatish Balay       z[cval+1]   += v[4]*x1 + v[5]*x2 + v[6]*x3 + v[7]*x4;
78349b5e25fSSatish Balay       z[cval+2]   += v[8]*x1 + v[9]*x2 + v[10]*x3 + v[11]*x4;
78449b5e25fSSatish Balay       z[cval+3]   += v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4;
78549b5e25fSSatish Balay       /* (strict upper triangular part of A)*x  */
78649b5e25fSSatish Balay       z[4*i]   += v[0]*x[cval] + v[4]*x[cval+1]+ v[8]*x[cval+2] + v[12]*x[cval+3];
78749b5e25fSSatish Balay       z[4*i+1] += v[1]*x[cval] + v[5]*x[cval+1]+ v[9]*x[cval+2] + v[13]*x[cval+3];
78849b5e25fSSatish Balay       z[4*i+2] += v[2]*x[cval] + v[6]*x[cval+1]+ v[10]*x[cval+2]+ v[14]*x[cval+3];
78949b5e25fSSatish Balay       z[4*i+3] += v[3]*x[cval] + v[7]*x[cval+1]+ v[11]*x[cval+2]+ v[15]*x[cval+3];
79049b5e25fSSatish Balay       v  += 16;
79149b5e25fSSatish Balay     }
79249b5e25fSSatish Balay     xb +=4; ai++;
79349b5e25fSSatish Balay   }
79449b5e25fSSatish Balay 
79549b5e25fSSatish Balay   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
79649b5e25fSSatish Balay   if (yy != xx) ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
79749b5e25fSSatish Balay   if (zz != yy) ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
79849b5e25fSSatish Balay 
79949b5e25fSSatish Balay   PLogFlops(32*(a->s_nz*2 - a->m));
80049b5e25fSSatish Balay   PetscFunctionReturn(0);
80149b5e25fSSatish Balay }
80249b5e25fSSatish Balay 
80349b5e25fSSatish Balay #undef __FUNC__
80449b5e25fSSatish Balay #define __FUNC__ "MatMultAdd_SeqSBAIJ_5"
80549b5e25fSSatish Balay int MatMultAdd_SeqSBAIJ_5(Mat A,Vec xx,Vec yy,Vec zz)
80649b5e25fSSatish Balay {
80749b5e25fSSatish Balay   Mat_SeqSBAIJ     *a = (Mat_SeqSBAIJ*)A->data;
80849b5e25fSSatish Balay   Scalar          *x,*y,*z,*xb,x1,x2,x3,x4,x5;
80949b5e25fSSatish Balay   MatScalar       *v;
81049b5e25fSSatish Balay   int             mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,ierr,*ib,cval,j;
81149b5e25fSSatish Balay 
81249b5e25fSSatish Balay   PetscFunctionBegin;
81349b5e25fSSatish Balay   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
81449b5e25fSSatish Balay   if (yy != xx) {
81549b5e25fSSatish Balay     ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
81649b5e25fSSatish Balay   } else {
81749b5e25fSSatish Balay     y = x;
81849b5e25fSSatish Balay   }
81949b5e25fSSatish Balay   if (zz != yy) {
82049b5e25fSSatish Balay     ierr = VecCopy(yy,zz);CHKERRQ(ierr);
82149b5e25fSSatish Balay     ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
82249b5e25fSSatish Balay   } else {
82349b5e25fSSatish Balay     z = y;
82449b5e25fSSatish Balay   }
82549b5e25fSSatish Balay 
82649b5e25fSSatish Balay   v     = a->a;
82749b5e25fSSatish Balay   xb = x;
82849b5e25fSSatish Balay 
82949b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
83049b5e25fSSatish Balay     n  = ai[1] - ai[0]; /* length of i_th block row of A */
83149b5e25fSSatish Balay     x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5=xb[4];
83249b5e25fSSatish Balay     ib = aj + *ai;
833*7fbae186SHong Zhang     if (*ib == i){      /* (diag of A)*x */
83449b5e25fSSatish Balay       z[5*i]   += v[0]*x1  + v[5]*x2 + v[10]*x3 + v[15]*x4+ v[20]*x5;
83549b5e25fSSatish Balay       z[5*i+1] += v[5]*x1  + v[6]*x2 + v[11]*x3 + v[16]*x4+ v[21]*x5;
83649b5e25fSSatish Balay       z[5*i+2] += v[10]*x1 +v[11]*x2 + v[12]*x3 + v[17]*x4+ v[22]*x5;
83749b5e25fSSatish Balay       z[5*i+3] += v[15]*x1 +v[16]*x2 + v[17]*x3 + v[18]*x4+ v[23]*x5;
83849b5e25fSSatish Balay       z[5*i+4] += v[20]*x1 +v[21]*x2 + v[22]*x3 + v[23]*x4+ v[24]*x5;
83949b5e25fSSatish Balay       v += 25;
840*7fbae186SHong Zhang     }
84149b5e25fSSatish Balay     for (j=1; j<n; j++) {
84249b5e25fSSatish Balay       /* (strict lower triangular part of A)*x  */
84349b5e25fSSatish Balay       cval       = ib[j]*5;
84449b5e25fSSatish Balay       z[cval]     += v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4 + v[4]*x5;
84549b5e25fSSatish Balay       z[cval+1]   += v[5]*x1 + v[6]*x2 + v[7]*x3 + v[8]*x4 + v[9]*x5;
84649b5e25fSSatish Balay       z[cval+2]   += v[10]*x1 + v[11]*x2 + v[12]*x3 + v[13]*x4+ v[14]*x5;
84749b5e25fSSatish Balay       z[cval+3]   += v[15]*x1 + v[16]*x2 + v[17]*x3 + v[18]*x4+ v[19]*x5;
84849b5e25fSSatish Balay       z[cval+4]   += v[20]*x1 + v[21]*x2 + v[22]*x3 + v[23]*x4+ v[24]*x5;
84949b5e25fSSatish Balay       /* (strict upper triangular part of A)*x  */
85049b5e25fSSatish 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];
85149b5e25fSSatish 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];
85249b5e25fSSatish 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];
85349b5e25fSSatish 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];
85449b5e25fSSatish 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];
85549b5e25fSSatish Balay       v  += 25;
85649b5e25fSSatish Balay     }
85749b5e25fSSatish Balay     xb +=5; ai++;
85849b5e25fSSatish Balay   }
85949b5e25fSSatish Balay 
86049b5e25fSSatish Balay   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
86149b5e25fSSatish Balay   if (yy != xx) ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
86249b5e25fSSatish Balay   if (zz != yy) ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
86349b5e25fSSatish Balay 
86449b5e25fSSatish Balay   PLogFlops(50*(a->s_nz*2 - a->m));
86549b5e25fSSatish Balay   PetscFunctionReturn(0);
86649b5e25fSSatish Balay }
86749b5e25fSSatish Balay #undef __FUNC__
86849b5e25fSSatish Balay #define __FUNC__ "MatMultAdd_SeqSBAIJ_6"
86949b5e25fSSatish Balay int MatMultAdd_SeqSBAIJ_6(Mat A,Vec xx,Vec yy,Vec zz)
87049b5e25fSSatish Balay {
87149b5e25fSSatish Balay   Mat_SeqSBAIJ     *a = (Mat_SeqSBAIJ*)A->data;
87249b5e25fSSatish Balay   Scalar          *x,*y,*z,*xb,x1,x2,x3,x4,x5,x6;
87349b5e25fSSatish Balay   MatScalar       *v;
87449b5e25fSSatish Balay   int             mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,ierr,*ib,cval,j;
87549b5e25fSSatish Balay 
87649b5e25fSSatish Balay   PetscFunctionBegin;
87749b5e25fSSatish Balay   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
87849b5e25fSSatish Balay   if (yy != xx) {
87949b5e25fSSatish Balay     ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
88049b5e25fSSatish Balay   } else {
88149b5e25fSSatish Balay     y = x;
88249b5e25fSSatish Balay   }
88349b5e25fSSatish Balay   if (zz != yy) {
88449b5e25fSSatish Balay     ierr = VecCopy(yy,zz);CHKERRQ(ierr);
88549b5e25fSSatish Balay     ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
88649b5e25fSSatish Balay   } else {
88749b5e25fSSatish Balay     z = y;
88849b5e25fSSatish Balay   }
88949b5e25fSSatish Balay 
89049b5e25fSSatish Balay   v     = a->a;
89149b5e25fSSatish Balay   xb = x;
89249b5e25fSSatish Balay 
89349b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
89449b5e25fSSatish Balay     n  = ai[1] - ai[0]; /* length of i_th block row of A */
89549b5e25fSSatish Balay     x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5=xb[4]; x6=xb[5];
89649b5e25fSSatish Balay     ib = aj + *ai;
897*7fbae186SHong Zhang     if (*ib == i){     /* (diag of A)*x */
89849b5e25fSSatish Balay       z[6*i]   += v[0]*x1  + v[6]*x2 + v[12]*x3 + v[18]*x4+ v[24]*x5 + v[30]*x6;
89949b5e25fSSatish Balay       z[6*i+1] += v[6]*x1  + v[7]*x2 + v[13]*x3 + v[19]*x4+ v[25]*x5 + v[31]*x6;
90049b5e25fSSatish Balay       z[6*i+2] += v[12]*x1 +v[13]*x2 + v[14]*x3 + v[20]*x4+ v[26]*x5 + v[32]*x6;
90149b5e25fSSatish Balay       z[6*i+3] += v[18]*x1 +v[19]*x2 + v[20]*x3 + v[21]*x4+ v[27]*x5 + v[33]*x6;
90249b5e25fSSatish Balay       z[6*i+4] += v[24]*x1 +v[25]*x2 + v[26]*x3 + v[27]*x4+ v[28]*x5 + v[34]*x6;
90349b5e25fSSatish Balay       z[6*i+5] += v[30]*x1 +v[31]*x2 + v[32]*x3 + v[33]*x4+ v[34]*x5 + v[35]*x6;
90449b5e25fSSatish Balay       v += 36;
905*7fbae186SHong Zhang     }
90649b5e25fSSatish Balay     for (j=1; j<n; j++) {
90749b5e25fSSatish Balay       /* (strict lower triangular part of A)*x  */
90849b5e25fSSatish Balay       cval       = ib[j]*6;
90949b5e25fSSatish Balay       z[cval]   += v[0]*x1  + v[1]*x2 + v[2]*x3 + v[3]*x4+ v[4]*x5 + v[5]*x6;
91049b5e25fSSatish Balay       z[cval+1] += v[6]*x1  + v[7]*x2 + v[8]*x3 + v[9]*x4+ v[10]*x5 + v[11]*x6;
91149b5e25fSSatish Balay       z[cval+2] += v[12]*x1  + v[13]*x2 + v[14]*x3 + v[15]*x4+ v[16]*x5 + v[17]*x6;
91249b5e25fSSatish Balay       z[cval+3] += v[18]*x1  + v[19]*x2 + v[20]*x3 + v[21]*x4+ v[22]*x5 + v[23]*x6;
91349b5e25fSSatish Balay       z[cval+4] += v[24]*x1  + v[25]*x2 + v[26]*x3 + v[27]*x4+ v[28]*x5 + v[29]*x6;
91449b5e25fSSatish Balay       z[cval+5] += v[30]*x1  + v[31]*x2 + v[32]*x3 + v[33]*x4+ v[34]*x5 + v[35]*x6;
91549b5e25fSSatish Balay       /* (strict upper triangular part of A)*x  */
91649b5e25fSSatish 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];
91749b5e25fSSatish 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];
91849b5e25fSSatish 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];
91949b5e25fSSatish 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];
92049b5e25fSSatish 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];
92149b5e25fSSatish 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];
92249b5e25fSSatish Balay       v  += 36;
92349b5e25fSSatish Balay     }
92449b5e25fSSatish Balay     xb +=6; ai++;
92549b5e25fSSatish Balay   }
92649b5e25fSSatish Balay 
92749b5e25fSSatish Balay   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
92849b5e25fSSatish Balay   if (yy != xx) ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
92949b5e25fSSatish Balay   if (zz != yy) ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
93049b5e25fSSatish Balay 
93149b5e25fSSatish Balay   PLogFlops(72*(a->s_nz*2 - a->m));
93249b5e25fSSatish Balay   PetscFunctionReturn(0);
93349b5e25fSSatish Balay }
93449b5e25fSSatish Balay 
93549b5e25fSSatish Balay #undef __FUNC__
93649b5e25fSSatish Balay #define __FUNC__ "MatMultAdd_SeqSBAIJ_7"
93749b5e25fSSatish Balay int MatMultAdd_SeqSBAIJ_7(Mat A,Vec xx,Vec yy,Vec zz)
93849b5e25fSSatish Balay {
93949b5e25fSSatish Balay   Mat_SeqSBAIJ     *a = (Mat_SeqSBAIJ*)A->data;
94049b5e25fSSatish Balay   Scalar          *x,*y,*z,*xb,x1,x2,x3,x4,x5,x6,x7;
94149b5e25fSSatish Balay   MatScalar       *v;
94249b5e25fSSatish Balay   int             mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,ierr,*ib,cval,j;
94349b5e25fSSatish Balay 
94449b5e25fSSatish Balay   PetscFunctionBegin;
94549b5e25fSSatish Balay   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
94649b5e25fSSatish Balay   if (yy != xx) {
94749b5e25fSSatish Balay     ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
94849b5e25fSSatish Balay   } else {
94949b5e25fSSatish Balay     y = x;
95049b5e25fSSatish Balay   }
95149b5e25fSSatish Balay   if (zz != yy) {
95249b5e25fSSatish Balay     ierr = VecCopy(yy,zz);CHKERRQ(ierr);
95349b5e25fSSatish Balay     ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
95449b5e25fSSatish Balay   } else {
95549b5e25fSSatish Balay     z = y;
95649b5e25fSSatish Balay   }
95749b5e25fSSatish Balay 
95849b5e25fSSatish Balay   v     = a->a;
95949b5e25fSSatish Balay   xb = x;
96049b5e25fSSatish Balay 
96149b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
96249b5e25fSSatish Balay     n  = ai[1] - ai[0]; /* length of i_th block row of A */
96349b5e25fSSatish Balay     x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5=xb[4]; x6=xb[5]; x7=xb[6];
96449b5e25fSSatish Balay     ib = aj + *ai;
965*7fbae186SHong Zhang     if (*ib == i){     /* (diag of A)*x */
96649b5e25fSSatish 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;
96749b5e25fSSatish 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;
96849b5e25fSSatish 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;
96949b5e25fSSatish 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;
97049b5e25fSSatish 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;
97149b5e25fSSatish 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;
97249b5e25fSSatish 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;
97349b5e25fSSatish Balay       v += 49;
974*7fbae186SHong Zhang     }
97549b5e25fSSatish Balay     for (j=1; j<n; j++) {
97649b5e25fSSatish Balay       /* (strict lower triangular part of A)*x  */
97749b5e25fSSatish Balay       cval       = ib[j]*7;
97849b5e25fSSatish Balay       z[cval]   += v[0]*x1  + v[1]*x2 + v[2]*x3 + v[3]*x4+ v[4]*x5 + v[5]*x6+ v[6]*x7;
97949b5e25fSSatish 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;
98049b5e25fSSatish 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;
98149b5e25fSSatish 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;
98249b5e25fSSatish 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;
98349b5e25fSSatish 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;
98449b5e25fSSatish 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;
98549b5e25fSSatish Balay       /* (strict upper triangular part of A)*x  */
98649b5e25fSSatish 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];
98749b5e25fSSatish 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];
98849b5e25fSSatish 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];
98949b5e25fSSatish 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];
99049b5e25fSSatish 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];
99149b5e25fSSatish 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];
99249b5e25fSSatish 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];
99349b5e25fSSatish Balay       v  += 49;
99449b5e25fSSatish Balay     }
99549b5e25fSSatish Balay     xb +=7; ai++;
99649b5e25fSSatish Balay   }
99749b5e25fSSatish Balay 
99849b5e25fSSatish Balay   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
99949b5e25fSSatish Balay   if (yy != xx) ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
100049b5e25fSSatish Balay   if (zz != yy) ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
100149b5e25fSSatish Balay 
100249b5e25fSSatish Balay   PLogFlops(98*(a->s_nz*2 - a->m));
100349b5e25fSSatish Balay   PetscFunctionReturn(0);
100449b5e25fSSatish Balay }
100549b5e25fSSatish Balay 
100649b5e25fSSatish Balay #undef __FUNC__
100749b5e25fSSatish Balay #define __FUNC__ "MatMultAdd_SeqSBAIJ_N"
100849b5e25fSSatish Balay int MatMultAdd_SeqSBAIJ_N(Mat A,Vec xx,Vec yy,Vec zz)
100949b5e25fSSatish Balay {
101049b5e25fSSatish Balay   Mat_SeqSBAIJ     *a = (Mat_SeqSBAIJ*)A->data;
1011ecfa4421SSatish Balay   Scalar          *x,*x_ptr,*y,*z,*z_ptr=0,*xb,*zb,*work,*workt;
1012066653e3SSatish Balay   MatScalar       *v;
1013066653e3SSatish Balay   int             ierr,mbs=a->mbs,i,*idx,*aj,*ii,bs=a->bs,j,n,bs2=a->bs2;
1014066653e3SSatish Balay   int             ncols,k;
101549b5e25fSSatish Balay 
101649b5e25fSSatish Balay   PetscFunctionBegin;
101749b5e25fSSatish Balay   ierr = VecGetArray(xx,&x);CHKERRQ(ierr); x_ptr=x;
101849b5e25fSSatish Balay   if (yy != xx) {
101949b5e25fSSatish Balay     ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
102049b5e25fSSatish Balay   } else {
102149b5e25fSSatish Balay     y = x;
102249b5e25fSSatish Balay   }
102349b5e25fSSatish Balay   if (zz != yy) {
102449b5e25fSSatish Balay     ierr = VecCopy(yy,zz);CHKERRQ(ierr);
102549b5e25fSSatish Balay     ierr = VecGetArray(zz,&z);CHKERRQ(ierr); z_ptr=z;
102649b5e25fSSatish Balay   } else {
102749b5e25fSSatish Balay     z = y;
102849b5e25fSSatish Balay   }
102949b5e25fSSatish Balay 
103049b5e25fSSatish Balay   aj   = a->j;
103149b5e25fSSatish Balay   v    = a->a;
103249b5e25fSSatish Balay   ii   = a->i;
103349b5e25fSSatish Balay 
103449b5e25fSSatish Balay   if (!a->mult_work) {
103549b5e25fSSatish Balay     k = a->m;
103649b5e25fSSatish Balay     a->mult_work = (Scalar*)PetscMalloc((k+1)*sizeof(Scalar));CHKPTRQ(a->mult_work);
103749b5e25fSSatish Balay   }
103849b5e25fSSatish Balay   work = a->mult_work;
103949b5e25fSSatish Balay 
104049b5e25fSSatish Balay 
104149b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
104249b5e25fSSatish Balay     n     = ii[1] - ii[0]; ncols = n*bs;
104349b5e25fSSatish Balay     workt = work; idx=aj+ii[0];
104449b5e25fSSatish Balay 
104549b5e25fSSatish Balay     /* upper triangular part */
104649b5e25fSSatish Balay     for (j=0; j<n; j++) {
1047066653e3SSatish Balay       /* col = bs*(*idx); */
104849b5e25fSSatish Balay       xb = x_ptr + bs*(*idx++);
104949b5e25fSSatish Balay       for (k=0; k<bs; k++) workt[k] = xb[k];
105049b5e25fSSatish Balay       workt += bs;
105149b5e25fSSatish Balay     }
105249b5e25fSSatish Balay     /* z(i*bs:(i+1)*bs-1) += A(i,:)*x */
105349b5e25fSSatish Balay     Kernel_w_gets_w_plus_Ar_times_v(bs,ncols,work,v,z);
105449b5e25fSSatish Balay 
105549b5e25fSSatish Balay     /* strict lower triangular part */
105649b5e25fSSatish Balay     ncols -= bs;
105749b5e25fSSatish Balay     if (ncols > 0){
105849b5e25fSSatish Balay       workt = work;
105949b5e25fSSatish Balay       ierr  = PetscMemzero(workt,ncols*sizeof(Scalar));CHKERRQ(ierr);
106049b5e25fSSatish Balay       Kernel_w_gets_w_plus_trans_Ar_times_v(bs,ncols,x,v+bs2,workt);
106149b5e25fSSatish Balay 
106249b5e25fSSatish Balay       idx=aj+ii[0]+1;
106349b5e25fSSatish Balay       for (j=1; j<n; j++) {
106449b5e25fSSatish Balay         zb = z_ptr + bs*(*idx);
106549b5e25fSSatish Balay         idx++;
106649b5e25fSSatish Balay         for (k=0; k<bs; k++) zb[k] += workt[k] ;
106749b5e25fSSatish Balay         workt += bs;
106849b5e25fSSatish Balay       }
106949b5e25fSSatish Balay     }
107049b5e25fSSatish Balay 
107149b5e25fSSatish Balay     x += bs; v += n*bs2; z += bs; ii++;
107249b5e25fSSatish Balay   }
107349b5e25fSSatish Balay 
107449b5e25fSSatish Balay   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
107549b5e25fSSatish Balay   if (yy != xx) ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
107649b5e25fSSatish Balay   if (zz != yy) ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
107749b5e25fSSatish Balay 
107849b5e25fSSatish Balay   PLogFlops(2*(a->s_nz*2 - a->m));
107949b5e25fSSatish Balay   PetscFunctionReturn(0);
108049b5e25fSSatish Balay }
108149b5e25fSSatish Balay 
108249b5e25fSSatish Balay #undef __FUNC__
108349b5e25fSSatish Balay #define __FUNC__ "MatMultTranspose_SeqSBAIJ"
108449b5e25fSSatish Balay int MatMultTranspose_SeqSBAIJ(Mat A,Vec xx,Vec zz)
108549b5e25fSSatish Balay {
108649b5e25fSSatish Balay   PetscFunctionBegin;
10870eaf17bfSBarry Smith   SETERRQ(1,"Matrix is symmetric. Call MatMult().");
10885535bfb1SHong Zhang   /* PetscFunctionReturn(0); */
108949b5e25fSSatish Balay }
109049b5e25fSSatish Balay 
109149b5e25fSSatish Balay #undef __FUNC__
109249b5e25fSSatish Balay #define __FUNC__ "MatMultTransposeAdd_SeqSBAIJ"
109349b5e25fSSatish Balay int MatMultTransposeAdd_SeqSBAIJ(Mat A,Vec xx,Vec yy,Vec zz)
109449b5e25fSSatish Balay 
109549b5e25fSSatish Balay {
109649b5e25fSSatish Balay   PetscFunctionBegin;
10970eaf17bfSBarry Smith   SETERRQ(1,"Matrix is symmetric. Call MatMultAdd().");
10985535bfb1SHong Zhang   /* PetscFunctionReturn(0); */
109949b5e25fSSatish Balay }
110049b5e25fSSatish Balay 
110149b5e25fSSatish Balay #undef __FUNC__
110249b5e25fSSatish Balay #define __FUNC__ "MatScale_SeqSBAIJ"
110349b5e25fSSatish Balay int MatScale_SeqSBAIJ(Scalar *alpha,Mat inA)
110449b5e25fSSatish Balay {
110549b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)inA->data;
110649b5e25fSSatish Balay   int         one = 1,totalnz = a->bs2*a->s_nz;
110749b5e25fSSatish Balay 
110849b5e25fSSatish Balay   PetscFunctionBegin;
110949b5e25fSSatish Balay   BLscal_(&totalnz,alpha,a->a,&one);
111049b5e25fSSatish Balay   PLogFlops(totalnz);
111149b5e25fSSatish Balay   PetscFunctionReturn(0);
111249b5e25fSSatish Balay }
111349b5e25fSSatish Balay 
111449b5e25fSSatish Balay #undef __FUNC__
111549b5e25fSSatish Balay #define __FUNC__ "MatNorm_SeqSBAIJ"
111649b5e25fSSatish Balay int MatNorm_SeqSBAIJ(Mat A,NormType type,PetscReal *norm)
111749b5e25fSSatish Balay {
111849b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
111949b5e25fSSatish Balay   MatScalar   *v = a->a;
112049b5e25fSSatish Balay   PetscReal   sum_diag = 0.0, sum_off = 0.0, *sum;
1121066653e3SSatish Balay   int         i,j,k,bs = a->bs,bs2=a->bs2,k1,mbs=a->mbs;
112249b5e25fSSatish Balay   int         *jl,*il,jmin,jmax,ierr,nexti,ik;
112349b5e25fSSatish Balay 
112449b5e25fSSatish Balay   PetscFunctionBegin;
112549b5e25fSSatish Balay   if (type == NORM_FROBENIUS) {
112649b5e25fSSatish Balay     for (k=0; k<mbs; k++){
112749b5e25fSSatish Balay       jmin = a->i[k]; jmax = a->i[k+1];
112849b5e25fSSatish Balay       /* diagonal block */
112949b5e25fSSatish Balay       for (i=0; i<bs2; i++){
113049b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX)
113149b5e25fSSatish Balay         sum_diag += PetscRealPart(PetscConj(*v)*(*v)); v++;
113249b5e25fSSatish Balay #else
113349b5e25fSSatish Balay         sum_diag += (*v)*(*v); v++;
113449b5e25fSSatish Balay #endif
113549b5e25fSSatish Balay       }
113649b5e25fSSatish Balay       for (j=jmin+1; j<jmax; j++){  /* off-diagonal blocks */
113749b5e25fSSatish Balay         for (i=0; i<bs2; i++){
113849b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX)
113949b5e25fSSatish Balay           sum_off += PetscRealPart(PetscConj(*v)*(*v)); v++;
114049b5e25fSSatish Balay #else
114149b5e25fSSatish Balay           sum_off += (*v)*(*v); v++;
114249b5e25fSSatish Balay #endif
114349b5e25fSSatish Balay         }
114449b5e25fSSatish Balay       }
114549b5e25fSSatish Balay     }
114649b5e25fSSatish Balay     *norm = sqrt(sum_diag + 2*sum_off);
114749b5e25fSSatish Balay 
114849b5e25fSSatish Balay   }  else if (type == NORM_INFINITY) { /* maximum row sum */
114949b5e25fSSatish Balay     il = (int*)PetscMalloc(mbs*sizeof(int));CHKPTRQ(il);
115049b5e25fSSatish Balay     jl = (int*)PetscMalloc(mbs*sizeof(int));CHKPTRQ(jl);
1151f3dd7b2dSKris Buschelman     sum = (PetscReal*)PetscMalloc(bs*sizeof(PetscReal));CHKPTRQ(sum);
115249b5e25fSSatish Balay     for (i=0; i<mbs; i++) {
115349b5e25fSSatish Balay       jl[i] = mbs; il[0] = 0;
115449b5e25fSSatish Balay     }
115549b5e25fSSatish Balay 
115649b5e25fSSatish Balay     *norm = 0.0;
115749b5e25fSSatish Balay     for (k=0; k<mbs; k++) { /* k_th block row */
115849b5e25fSSatish Balay       for (j=0; j<bs; j++) sum[j]=0.0;
115949b5e25fSSatish Balay 
116049b5e25fSSatish Balay       /*-- col sum --*/
116149b5e25fSSatish Balay       i = jl[k]; /* first |A(i,k)| to be added */
116249b5e25fSSatish Balay       while (i<mbs){
116349b5e25fSSatish Balay         nexti = jl[i];  /* next block row to be added */
116449b5e25fSSatish Balay         ik    = il[i];  /* block index of A(i,k) in the array a */
116549b5e25fSSatish Balay         for (j=0; j<bs; j++){
116649b5e25fSSatish Balay           v = a->a + ik*bs2 + j*bs;
116749b5e25fSSatish Balay           for (k1=0; k1<bs; k1++) {
116849b5e25fSSatish Balay             sum[j] += PetscAbsScalar(*v); v++;
116949b5e25fSSatish Balay           }
117049b5e25fSSatish Balay         }
117149b5e25fSSatish Balay         /* update il, jl */
117249b5e25fSSatish Balay         jmin = ik + 1; jmax = a->i[i+1];
117349b5e25fSSatish Balay         if (jmin < jmax){
117449b5e25fSSatish Balay           il[i] = jmin;
117549b5e25fSSatish Balay           j   = a->j[jmin];
117649b5e25fSSatish Balay           jl[i] = jl[j]; jl[j]=i;
117749b5e25fSSatish Balay         }
117849b5e25fSSatish Balay         i = nexti;
117949b5e25fSSatish Balay       }
118049b5e25fSSatish Balay 
118149b5e25fSSatish Balay       /*-- row sum --*/
118249b5e25fSSatish Balay       jmin = a->i[k]; jmax = a->i[k+1];
118349b5e25fSSatish Balay       for (i=jmin; i<jmax; i++) {
118449b5e25fSSatish Balay         for (j=0; j<bs; j++){
118549b5e25fSSatish Balay           v = a->a + i*bs2 + j;
118649b5e25fSSatish Balay           for (k1=0; k1<bs; k1++){
118749b5e25fSSatish Balay             sum[j] += PetscAbsScalar(*v);
118849b5e25fSSatish Balay             v   += bs;
118949b5e25fSSatish Balay           }
119049b5e25fSSatish Balay         }
119149b5e25fSSatish Balay       }
119249b5e25fSSatish Balay       /* add k_th block row to il, jl */
119349b5e25fSSatish Balay       jmin++;
119449b5e25fSSatish Balay       if (jmin < jmax){
119549b5e25fSSatish Balay         il[k] = jmin;
119649b5e25fSSatish Balay         j   = a->j[jmin];
119749b5e25fSSatish Balay         jl[k] = jl[j]; jl[j] = k;
119849b5e25fSSatish Balay       }
119949b5e25fSSatish Balay       for (j=0; j<bs; j++){
120049b5e25fSSatish Balay         if (sum[j] > *norm) *norm = sum[j];
120149b5e25fSSatish Balay       }
120249b5e25fSSatish Balay     }
120349b5e25fSSatish Balay     ierr = PetscFree(il);CHKERRQ(ierr);
120449b5e25fSSatish Balay     ierr = PetscFree(jl);CHKERRQ(ierr);
120549b5e25fSSatish Balay     ierr = PetscFree(sum);CHKERRQ(ierr);
120649b5e25fSSatish Balay   } else {
1207347d480fSBarry Smith     SETERRQ(PETSC_ERR_SUP,"No support for this norm yet");
120849b5e25fSSatish Balay   }
120949b5e25fSSatish Balay   PetscFunctionReturn(0);
121049b5e25fSSatish Balay }
121149b5e25fSSatish Balay 
121249b5e25fSSatish Balay #undef __FUNC__
121349b5e25fSSatish Balay #define __FUNC__ "MatEqual_SeqSBAIJ"
121449b5e25fSSatish Balay int MatEqual_SeqSBAIJ(Mat A,Mat B,PetscTruth* flg)
121549b5e25fSSatish Balay {
121649b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ *)A->data,*b = (Mat_SeqSBAIJ *)B->data;
121749b5e25fSSatish Balay   int         ierr;
121849b5e25fSSatish Balay 
121949b5e25fSSatish Balay   PetscFunctionBegin;
1220347d480fSBarry Smith   if (B->type !=MATSEQSBAIJ) SETERRQ(PETSC_ERR_ARG_INCOMP,"Matrices must be same type");
122149b5e25fSSatish Balay 
122249b5e25fSSatish Balay   /* If the  matrix/block dimensions are not equal, or no of nonzeros or shift */
122349b5e25fSSatish Balay   if ((a->m != b->m) || (a->n !=b->n) || (a->bs != b->bs)|| (a->s_nz != b->s_nz)) {
122449b5e25fSSatish Balay     *flg = PETSC_FALSE; PetscFunctionReturn(0);
122549b5e25fSSatish Balay   }
122649b5e25fSSatish Balay 
122749b5e25fSSatish Balay   /* if the a->i are the same */
122849b5e25fSSatish Balay   ierr = PetscMemcmp(a->i,b->i,(a->mbs+1)*sizeof(int),flg);CHKERRQ(ierr);
122949b5e25fSSatish Balay   if (*flg == PETSC_FALSE) {
123049b5e25fSSatish Balay     PetscFunctionReturn(0);
123149b5e25fSSatish Balay   }
123249b5e25fSSatish Balay 
123349b5e25fSSatish Balay   /* if a->j are the same */
123449b5e25fSSatish Balay   ierr = PetscMemcmp(a->j,b->j,(a->s_nz)*sizeof(int),flg);CHKERRQ(ierr);
123549b5e25fSSatish Balay   if (*flg == PETSC_FALSE) {
123649b5e25fSSatish Balay     PetscFunctionReturn(0);
123749b5e25fSSatish Balay   }
123849b5e25fSSatish Balay   /* if a->a are the same */
123949b5e25fSSatish Balay   ierr = PetscMemcmp(a->a,b->a,(a->s_nz)*(a->bs)*(a->bs)*sizeof(Scalar),flg);CHKERRQ(ierr);
124049b5e25fSSatish Balay   PetscFunctionReturn(0);
124149b5e25fSSatish Balay 
124249b5e25fSSatish Balay }
124349b5e25fSSatish Balay 
124449b5e25fSSatish Balay #undef __FUNC__
124549b5e25fSSatish Balay #define __FUNC__ "MatGetDiagonal_SeqSBAIJ"
124649b5e25fSSatish Balay int MatGetDiagonal_SeqSBAIJ(Mat A,Vec v)
124749b5e25fSSatish Balay {
124849b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
124949b5e25fSSatish Balay   int         ierr,i,j,k,n,row,bs,*ai,*aj,ambs,bs2;
125049b5e25fSSatish Balay   Scalar      *x,zero = 0.0;
125149b5e25fSSatish Balay   MatScalar   *aa,*aa_j;
125249b5e25fSSatish Balay 
125349b5e25fSSatish Balay   PetscFunctionBegin;
1254347d480fSBarry Smith   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
125549b5e25fSSatish Balay   bs   = a->bs;
125649b5e25fSSatish Balay   aa   = a->a;
125749b5e25fSSatish Balay   ai   = a->i;
125849b5e25fSSatish Balay   aj   = a->j;
125949b5e25fSSatish Balay   ambs = a->mbs;
126049b5e25fSSatish Balay   bs2  = a->bs2;
126149b5e25fSSatish Balay 
126249b5e25fSSatish Balay   ierr = VecSet(&zero,v);CHKERRQ(ierr);
126349b5e25fSSatish Balay   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
126449b5e25fSSatish Balay   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
1265347d480fSBarry Smith   /* if (n != a->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");*/
126649b5e25fSSatish Balay   for (i=0; i<ambs; i++) {
126749b5e25fSSatish Balay     j=ai[i];
126849b5e25fSSatish Balay     if (aj[j] == i) {             /* if this is a diagonal element */
126949b5e25fSSatish Balay       row  = i*bs;
127049b5e25fSSatish Balay       aa_j = aa + j*bs2;
127149b5e25fSSatish Balay       for (k=0; k<bs2; k+=(bs+1),row++) x[row] = aa_j[k];
127249b5e25fSSatish Balay     }
127349b5e25fSSatish Balay   }
127449b5e25fSSatish Balay   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
127549b5e25fSSatish Balay   PetscFunctionReturn(0);
127649b5e25fSSatish Balay }
127749b5e25fSSatish Balay 
127849b5e25fSSatish Balay #undef __FUNC__
127949b5e25fSSatish Balay #define __FUNC__ "MatDiagonalScale_SeqSBAIJ"
128049b5e25fSSatish Balay int MatDiagonalScale_SeqSBAIJ(Mat A,Vec ll,Vec rr)
128149b5e25fSSatish Balay {
128249b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
128349b5e25fSSatish Balay   Scalar      *l,*r,x,*li,*ri;
128449b5e25fSSatish Balay   MatScalar   *aa,*v;
1285066653e3SSatish Balay   int         ierr,i,j,k,lm,rn,M,m,*ai,*aj,mbs,tmp,bs,bs2;
128649b5e25fSSatish Balay 
128749b5e25fSSatish Balay   PetscFunctionBegin;
128849b5e25fSSatish Balay   ai  = a->i;
128949b5e25fSSatish Balay   aj  = a->j;
129049b5e25fSSatish Balay   aa  = a->a;
129149b5e25fSSatish Balay   m   = a->m;
129249b5e25fSSatish Balay   bs  = a->bs;
129349b5e25fSSatish Balay   mbs = a->mbs;
129449b5e25fSSatish Balay   bs2 = a->bs2;
129549b5e25fSSatish Balay 
129649b5e25fSSatish Balay   if (ll != rr) {
1297347d480fSBarry Smith     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"For symmetric format, left and right scaling vectors must be same\n");
129849b5e25fSSatish Balay   }
129949b5e25fSSatish Balay   if (ll) {
130049b5e25fSSatish Balay     ierr = VecGetArray(ll,&l);CHKERRQ(ierr);
130149b5e25fSSatish Balay     ierr = VecGetLocalSize(ll,&lm);CHKERRQ(ierr);
1302347d480fSBarry Smith     if (lm != m) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vector wrong length");
130349b5e25fSSatish Balay     for (i=0; i<mbs; i++) { /* for each block row */
130449b5e25fSSatish Balay       M  = ai[i+1] - ai[i];
130549b5e25fSSatish Balay       li = l + i*bs;
130649b5e25fSSatish Balay       v  = aa + bs2*ai[i];
130749b5e25fSSatish Balay       for (j=0; j<M; j++) { /* for each block */
130849b5e25fSSatish Balay         for (k=0; k<bs2; k++) {
130949b5e25fSSatish Balay           (*v++) *= li[k%bs];
131049b5e25fSSatish Balay         }
131149b5e25fSSatish Balay #ifdef CONT
131249b5e25fSSatish Balay         /* will be used to replace the above loop */
131349b5e25fSSatish Balay         ri = l + bs*aj[ai[i]+j];
131449b5e25fSSatish Balay         for (k=0; k<bs; k++) { /* column value */
131549b5e25fSSatish Balay           x = ri[k];
131649b5e25fSSatish Balay           for (tmp=0; tmp<bs; tmp++) (*v++) *= li[tmp]*x;
131749b5e25fSSatish Balay         }
131849b5e25fSSatish Balay #endif
131949b5e25fSSatish Balay 
132049b5e25fSSatish Balay       }
132149b5e25fSSatish Balay     }
132249b5e25fSSatish Balay     ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr);
132349b5e25fSSatish Balay     PLogFlops(2*a->s_nz);
132449b5e25fSSatish Balay   }
132549b5e25fSSatish Balay   /* will be deleted */
132649b5e25fSSatish Balay   if (rr) {
132749b5e25fSSatish Balay     ierr = VecGetArray(rr,&r);CHKERRQ(ierr);
132849b5e25fSSatish Balay     ierr = VecGetLocalSize(rr,&rn);CHKERRQ(ierr);
1329347d480fSBarry Smith     if (rn != m) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vector wrong length");
133049b5e25fSSatish Balay     for (i=0; i<mbs; i++) { /* for each block row */
133149b5e25fSSatish Balay       M  = ai[i+1] - ai[i];
133249b5e25fSSatish Balay       v  = aa + bs2*ai[i];
133349b5e25fSSatish Balay       for (j=0; j<M; j++) { /* for each block */
133449b5e25fSSatish Balay         ri = r + bs*aj[ai[i]+j];
133549b5e25fSSatish Balay         for (k=0; k<bs; k++) {
133649b5e25fSSatish Balay           x = ri[k];
133749b5e25fSSatish Balay           for (tmp=0; tmp<bs; tmp++) (*v++) *= x;
133849b5e25fSSatish Balay         }
133949b5e25fSSatish Balay       }
134049b5e25fSSatish Balay     }
134149b5e25fSSatish Balay     ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr);
134249b5e25fSSatish Balay     PLogFlops(a->s_nz);
134349b5e25fSSatish Balay   }
134449b5e25fSSatish Balay   PetscFunctionReturn(0);
134549b5e25fSSatish Balay }
134649b5e25fSSatish Balay 
134749b5e25fSSatish Balay #undef __FUNC__
134849b5e25fSSatish Balay #define __FUNC__ "MatGetInfo_SeqSBAIJ"
134949b5e25fSSatish Balay int MatGetInfo_SeqSBAIJ(Mat A,MatInfoType flag,MatInfo *info)
135049b5e25fSSatish Balay {
135149b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
135249b5e25fSSatish Balay 
135349b5e25fSSatish Balay   PetscFunctionBegin;
135449b5e25fSSatish Balay   info->rows_global    = (double)a->m;
135549b5e25fSSatish Balay   info->columns_global = (double)a->m;
135649b5e25fSSatish Balay   info->rows_local     = (double)a->m;
135749b5e25fSSatish Balay   info->columns_local  = (double)a->m;
135849b5e25fSSatish Balay   info->block_size     = a->bs2;
135949b5e25fSSatish Balay   info->nz_allocated   = a->s_maxnz; /*num. of nonzeros in upper triangular part */
136049b5e25fSSatish Balay   info->nz_used        = a->bs2*a->s_nz; /*num. of nonzeros in upper triangular part */
136149b5e25fSSatish Balay   info->nz_unneeded    = (double)(info->nz_allocated - info->nz_used);
136249b5e25fSSatish Balay   info->assemblies   = A->num_ass;
136349b5e25fSSatish Balay   info->mallocs      = a->reallocs;
136449b5e25fSSatish Balay   info->memory       = A->mem;
136549b5e25fSSatish Balay   if (A->factor) {
136649b5e25fSSatish Balay     info->fill_ratio_given  = A->info.fill_ratio_given;
136749b5e25fSSatish Balay     info->fill_ratio_needed = A->info.fill_ratio_needed;
136849b5e25fSSatish Balay     info->factor_mallocs    = A->info.factor_mallocs;
136949b5e25fSSatish Balay   } else {
137049b5e25fSSatish Balay     info->fill_ratio_given  = 0;
137149b5e25fSSatish Balay     info->fill_ratio_needed = 0;
137249b5e25fSSatish Balay     info->factor_mallocs    = 0;
137349b5e25fSSatish Balay   }
137449b5e25fSSatish Balay   PetscFunctionReturn(0);
137549b5e25fSSatish Balay }
137649b5e25fSSatish Balay 
137749b5e25fSSatish Balay 
137849b5e25fSSatish Balay #undef __FUNC__
137949b5e25fSSatish Balay #define __FUNC__ "MatZeroEntries_SeqSBAIJ"
138049b5e25fSSatish Balay int MatZeroEntries_SeqSBAIJ(Mat A)
138149b5e25fSSatish Balay {
138249b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
138349b5e25fSSatish Balay   int         ierr;
138449b5e25fSSatish Balay 
138549b5e25fSSatish Balay   PetscFunctionBegin;
138649b5e25fSSatish Balay   ierr = PetscMemzero(a->a,a->bs2*a->i[a->mbs]*sizeof(MatScalar));CHKERRQ(ierr);
138749b5e25fSSatish Balay   PetscFunctionReturn(0);
138849b5e25fSSatish Balay }
1389