xref: /petsc/src/mat/impls/sbaij/seq/sbaij2.c (revision dbe03f88ef651362992d50b70cdfc3144f10fda2)
173f4d377SMatthew Knepley /*$Id: sbaij2.c,v 1.32 2001/08/07 03:03:01 balay Exp $*/
249b5e25fSSatish Balay 
349b5e25fSSatish Balay #include "src/mat/impls/baij/seq/baij.h"
449b5e25fSSatish Balay #include "src/inline/spops.h"
549b5e25fSSatish Balay #include "src/inline/ilu.h"
649b5e25fSSatish Balay #include "petscbt.h"
73a7fca6bSBarry Smith #include "src/mat/impls/sbaij/seq/sbaij.h"
849b5e25fSSatish Balay 
94a2ae208SSatish Balay #undef __FUNCT__
104a2ae208SSatish Balay #define __FUNCT__ "MatIncreaseOverlap_SeqSBAIJ"
11268466fbSBarry Smith int MatIncreaseOverlap_SeqSBAIJ(Mat A,int is_max,IS is[],int ov)
1249b5e25fSSatish Balay {
135eee224dSHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
14*dbe03f88SHong Zhang   int          brow,i,j,k,l,mbs,n,*idx,ierr,*nidx,isz,bcol,bcol_max,
15d94109b8SHong Zhang                start,end,*ai,*aj,bs,*nidx2;
16d94109b8SHong Zhang   PetscBT      table;
17d94109b8SHong Zhang   PetscBT      table0;
18d94109b8SHong Zhang 
19d94109b8SHong Zhang   PetscFunctionBegin;
20d94109b8SHong Zhang   mbs = a->mbs;
21d94109b8SHong Zhang   ai  = a->i;
22d94109b8SHong Zhang   aj  = a->j;
23d94109b8SHong Zhang   bs  = a->bs;
24d94109b8SHong Zhang 
25d94109b8SHong Zhang   if (ov < 0)  SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap specified");
26d94109b8SHong Zhang 
27d94109b8SHong Zhang   ierr = PetscBTCreate(mbs,table);CHKERRQ(ierr);
28d94109b8SHong Zhang   ierr = PetscMalloc((mbs+1)*sizeof(int),&nidx);CHKERRQ(ierr);
29d94109b8SHong Zhang   ierr = PetscMalloc((A->m+1)*sizeof(int),&nidx2);CHKERRQ(ierr);
30d94109b8SHong Zhang   ierr = PetscBTCreate(mbs,table0);CHKERRQ(ierr);
31d94109b8SHong Zhang 
32d94109b8SHong Zhang   for (i=0; i<is_max; i++) { /* for each is */
33d94109b8SHong Zhang     isz  = 0;
34d94109b8SHong Zhang     ierr = PetscBTMemzero(mbs,table);CHKERRQ(ierr);
35d94109b8SHong Zhang 
36d94109b8SHong Zhang     /* Extract the indices, assume there can be duplicate entries */
37d94109b8SHong Zhang     ierr = ISGetIndices(is[i],&idx);CHKERRQ(ierr);
38d94109b8SHong Zhang     ierr = ISGetLocalSize(is[i],&n);CHKERRQ(ierr);
39d94109b8SHong Zhang 
40d94109b8SHong Zhang     /* Enter these into the temp arrays i.e mark table[brow], enter brow into new index */
41*dbe03f88SHong Zhang     bcol_max = 0;
42d94109b8SHong Zhang     for (j=0; j<n ; ++j){
43d94109b8SHong Zhang       brow = idx[j]/bs; /* convert the indices into block indices */
44d94109b8SHong Zhang       if (brow >= mbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"index greater than mat-dim");
45*dbe03f88SHong Zhang       if(!PetscBTLookupSet(table,brow)) {
46*dbe03f88SHong Zhang         nidx[isz++] = brow;
47*dbe03f88SHong Zhang         if (bcol_max < brow) bcol_max = brow;
48*dbe03f88SHong Zhang       }
49d94109b8SHong Zhang     }
50d94109b8SHong Zhang     ierr = ISRestoreIndices(is[i],&idx);CHKERRQ(ierr);
51d94109b8SHong Zhang     ierr = ISDestroy(is[i]);CHKERRQ(ierr);
52d94109b8SHong Zhang 
53d94109b8SHong Zhang     k = 0;
54d94109b8SHong Zhang     for (j=0; j<ov; j++){ /* for each overlap */
55d94109b8SHong Zhang       /* set table0 for lookup - only mark entries that are added onto nidx in (j-1)-th overlap */
56d94109b8SHong Zhang       ierr = PetscBTMemzero(mbs,table0);CHKERRQ(ierr);
57d94109b8SHong Zhang       for (l=k; l<isz; l++) PetscBTSet(table0,nidx[l]);
58d94109b8SHong Zhang 
59d94109b8SHong Zhang       n = isz;  /* length of the updated is[i] */
60d94109b8SHong Zhang       for (brow=0; brow<mbs; brow++){
61d94109b8SHong Zhang         start = ai[brow]; end   = ai[brow+1];
62d94109b8SHong Zhang         if (PetscBTLookup(table0,brow)){ /* brow is on nidx - row search: collect all bcol in this brow */
63d94109b8SHong Zhang           for (l = start; l<end ; l++){
64d94109b8SHong Zhang             bcol = aj[l];
65d94109b8SHong Zhang             if (!PetscBTLookupSet(table,bcol)) {nidx[isz++] = bcol;}
66d94109b8SHong Zhang           }
67d94109b8SHong Zhang           k++;
68d94109b8SHong Zhang           if (k >= n) break; /* for (brow=0; brow<mbs; brow++) */
69d94109b8SHong Zhang         } else { /* brow is not on nidx - col serach: add brow onto nidx if there is a bcol in nidx */
70d94109b8SHong Zhang           for (l = start; l<end ; l++){
71d94109b8SHong Zhang             bcol = aj[l];
72*dbe03f88SHong Zhang             if (bcol > bcol_max) break;
73d94109b8SHong Zhang             if (PetscBTLookup(table0,bcol)){
74d94109b8SHong Zhang               if (!PetscBTLookupSet(table,brow)) {nidx[isz++] = brow;}
75d94109b8SHong Zhang               break; /* for l = start; l<end ; l++) */
76d94109b8SHong Zhang             }
77d94109b8SHong Zhang           }
78d94109b8SHong Zhang         }
79d94109b8SHong Zhang       }
80d94109b8SHong Zhang     } /* for each overlap */
81d94109b8SHong Zhang 
82d94109b8SHong Zhang     /* expand the Index Set */
83d94109b8SHong Zhang     for (j=0; j<isz; j++) {
84d94109b8SHong Zhang       for (k=0; k<bs; k++)
85d94109b8SHong Zhang         nidx2[j*bs+k] = nidx[j]*bs+k;
86d94109b8SHong Zhang     }
87d94109b8SHong Zhang     ierr = ISCreateGeneral(PETSC_COMM_SELF,isz*bs,nidx2,is+i);CHKERRQ(ierr);
88d94109b8SHong Zhang   }
89d94109b8SHong Zhang   ierr = PetscBTDestroy(table);CHKERRQ(ierr);
90d94109b8SHong Zhang   ierr = PetscFree(nidx);CHKERRQ(ierr);
91d94109b8SHong Zhang   ierr = PetscFree(nidx2);CHKERRQ(ierr);
92d94109b8SHong Zhang   ierr = PetscBTDestroy(table0);CHKERRQ(ierr);
935eee224dSHong Zhang   PetscFunctionReturn(0);
9449b5e25fSSatish Balay }
9549b5e25fSSatish Balay 
964a2ae208SSatish Balay #undef __FUNCT__
974a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqSBAIJ_Private"
9849b5e25fSSatish Balay int MatGetSubMatrix_SeqSBAIJ_Private(Mat A,IS isrow,IS iscol,int cs,MatReuse scall,Mat *B)
9949b5e25fSSatish Balay {
10049b5e25fSSatish Balay   Mat_SeqSBAIJ  *a = (Mat_SeqSBAIJ*)A->data,*c;
10149b5e25fSSatish Balay   int          *smap,i,k,kstart,kend,ierr,oldcols = a->mbs,*lens;
10249b5e25fSSatish Balay   int          row,mat_i,*mat_j,tcol,*mat_ilen;
10349b5e25fSSatish Balay   int          *irow,nrows,*ssmap,bs=a->bs,bs2=a->bs2;
10449b5e25fSSatish Balay   int          *aj = a->j,*ai = a->i;
10549b5e25fSSatish Balay   MatScalar    *mat_a;
10649b5e25fSSatish Balay   Mat          C;
10749b5e25fSSatish Balay   PetscTruth   flag;
10849b5e25fSSatish Balay 
10949b5e25fSSatish Balay   PetscFunctionBegin;
11049b5e25fSSatish Balay 
111ac355199SBarry Smith   if (isrow != iscol) SETERRQ(1,"MatGetSubmatrices_SeqSBAIJ: For symm. format, iscol must equal isro");
11249b5e25fSSatish Balay   ierr = ISSorted(iscol,(PetscTruth*)&i);CHKERRQ(ierr);
113347d480fSBarry Smith   if (!i) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"IS is not sorted");
11449b5e25fSSatish Balay 
11549b5e25fSSatish Balay   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
11649b5e25fSSatish Balay   ierr = ISGetSize(isrow,&nrows);CHKERRQ(ierr);
11749b5e25fSSatish Balay 
118b0a32e0cSBarry Smith   ierr  = PetscMalloc((1+oldcols)*sizeof(int),&smap);CHKERRQ(ierr);
11949b5e25fSSatish Balay   ssmap = smap;
120b0a32e0cSBarry Smith   ierr  = PetscMalloc((1+nrows)*sizeof(int),&lens);CHKERRQ(ierr);
12149b5e25fSSatish Balay   ierr  = PetscMemzero(smap,oldcols*sizeof(int));CHKERRQ(ierr);
12249b5e25fSSatish Balay   for (i=0; i<nrows; i++) smap[irow[i]] = i+1; /* nrows = ncols */
12349b5e25fSSatish Balay   /* determine lens of each row */
12449b5e25fSSatish Balay   for (i=0; i<nrows; i++) {
12549b5e25fSSatish Balay     kstart  = ai[irow[i]];
12649b5e25fSSatish Balay     kend    = kstart + a->ilen[irow[i]];
12749b5e25fSSatish Balay     lens[i] = 0;
12849b5e25fSSatish Balay       for (k=kstart; k<kend; k++) {
12949b5e25fSSatish Balay         if (ssmap[aj[k]]) {
13049b5e25fSSatish Balay           lens[i]++;
13149b5e25fSSatish Balay         }
13249b5e25fSSatish Balay       }
13349b5e25fSSatish Balay     }
13449b5e25fSSatish Balay   /* Create and fill new matrix */
13549b5e25fSSatish Balay   if (scall == MAT_REUSE_MATRIX) {
13649b5e25fSSatish Balay     c = (Mat_SeqSBAIJ *)((*B)->data);
13749b5e25fSSatish Balay 
138347d480fSBarry Smith     if (c->mbs!=nrows || c->bs!=bs) SETERRQ(PETSC_ERR_ARG_SIZ,"Submatrix wrong size");
13949b5e25fSSatish Balay     ierr = PetscMemcmp(c->ilen,lens,c->mbs *sizeof(int),&flag);CHKERRQ(ierr);
14049b5e25fSSatish Balay     if (flag == PETSC_FALSE) {
141347d480fSBarry Smith       SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong no of nonzeros");
14249b5e25fSSatish Balay     }
14349b5e25fSSatish Balay     ierr = PetscMemzero(c->ilen,c->mbs*sizeof(int));CHKERRQ(ierr);
14449b5e25fSSatish Balay     C = *B;
14549b5e25fSSatish Balay   } else {
146e2d9671bSKris Buschelman     ierr = MatCreate(A->comm,nrows*bs,nrows*bs,PETSC_DETERMINE,PETSC_DETERMINE,&C);CHKERRQ(ierr);
147e2d9671bSKris Buschelman     ierr = MatSetType(C,A->type_name);CHKERRQ(ierr);
148e2d9671bSKris Buschelman     ierr = MatSeqSBAIJSetPreallocation(C,bs,0,lens);CHKERRQ(ierr);
14949b5e25fSSatish Balay   }
15049b5e25fSSatish Balay   c = (Mat_SeqSBAIJ *)(C->data);
15149b5e25fSSatish Balay   for (i=0; i<nrows; i++) {
15249b5e25fSSatish Balay     row    = irow[i];
15349b5e25fSSatish Balay     kstart = ai[row];
15449b5e25fSSatish Balay     kend   = kstart + a->ilen[row];
15549b5e25fSSatish Balay     mat_i  = c->i[i];
15649b5e25fSSatish Balay     mat_j  = c->j + mat_i;
15749b5e25fSSatish Balay     mat_a  = c->a + mat_i*bs2;
15849b5e25fSSatish Balay     mat_ilen = c->ilen + i;
15949b5e25fSSatish Balay     for (k=kstart; k<kend; k++) {
16049b5e25fSSatish Balay       if ((tcol=ssmap[a->j[k]])) {
16149b5e25fSSatish Balay         *mat_j++ = tcol - 1;
16249b5e25fSSatish Balay         ierr     = PetscMemcpy(mat_a,a->a+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
16349b5e25fSSatish Balay         mat_a   += bs2;
16449b5e25fSSatish Balay         (*mat_ilen)++;
16549b5e25fSSatish Balay       }
16649b5e25fSSatish Balay     }
16749b5e25fSSatish Balay   }
16849b5e25fSSatish Balay 
16949b5e25fSSatish Balay   /* Free work space */
17049b5e25fSSatish Balay   ierr = PetscFree(smap);CHKERRQ(ierr);
17149b5e25fSSatish Balay   ierr = PetscFree(lens);CHKERRQ(ierr);
17249b5e25fSSatish Balay   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
17349b5e25fSSatish Balay   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
17449b5e25fSSatish Balay 
17549b5e25fSSatish Balay   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
17649b5e25fSSatish Balay   *B = C;
17749b5e25fSSatish Balay   PetscFunctionReturn(0);
17849b5e25fSSatish Balay }
17949b5e25fSSatish Balay 
1804a2ae208SSatish Balay #undef __FUNCT__
1814a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqSBAIJ"
18249b5e25fSSatish Balay int MatGetSubMatrix_SeqSBAIJ(Mat A,IS isrow,IS iscol,int cs,MatReuse scall,Mat *B)
18349b5e25fSSatish Balay {
18449b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
18549b5e25fSSatish Balay   IS          is1;
18649b5e25fSSatish Balay   int         *vary,*iary,*irow,nrows,i,ierr,bs=a->bs,count;
18749b5e25fSSatish Balay 
18849b5e25fSSatish Balay   PetscFunctionBegin;
189ac355199SBarry Smith   if (isrow != iscol) SETERRQ(1,"MatGetSubmatrices_SeqSBAIJ: For symm. format, iscol must equal isro");
19049b5e25fSSatish Balay 
19149b5e25fSSatish Balay   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
19249b5e25fSSatish Balay   ierr = ISGetSize(isrow,&nrows);CHKERRQ(ierr);
19349b5e25fSSatish Balay 
19449b5e25fSSatish Balay   /* Verify if the indices corespond to each element in a block
19549b5e25fSSatish Balay    and form the IS with compressed IS */
19682502324SSatish Balay   ierr = PetscMalloc(2*(a->mbs+1)*sizeof(int),&vary);CHKERRQ(ierr);
19749b5e25fSSatish Balay   iary = vary + a->mbs;
19849b5e25fSSatish Balay   ierr = PetscMemzero(vary,(a->mbs)*sizeof(int));CHKERRQ(ierr);
19949b5e25fSSatish Balay   for (i=0; i<nrows; i++) vary[irow[i]/bs]++;
20049b5e25fSSatish Balay 
20149b5e25fSSatish Balay   count = 0;
20249b5e25fSSatish Balay   for (i=0; i<a->mbs; i++) {
203ac355199SBarry Smith     if (vary[i]!=0 && vary[i]!=bs) SETERRQ(1,"Index set does not match blocks");
20449b5e25fSSatish Balay     if (vary[i]==bs) iary[count++] = i;
20549b5e25fSSatish Balay   }
20649b5e25fSSatish Balay   ierr = ISCreateGeneral(PETSC_COMM_SELF,count,iary,&is1);CHKERRQ(ierr);
20749b5e25fSSatish Balay 
20849b5e25fSSatish Balay   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
20949b5e25fSSatish Balay   ierr = PetscFree(vary);CHKERRQ(ierr);
21049b5e25fSSatish Balay 
21149b5e25fSSatish Balay   ierr = MatGetSubMatrix_SeqSBAIJ_Private(A,is1,is1,cs,scall,B);CHKERRQ(ierr);
21249b5e25fSSatish Balay   ISDestroy(is1);
21349b5e25fSSatish Balay   PetscFunctionReturn(0);
21449b5e25fSSatish Balay }
21549b5e25fSSatish Balay 
2164a2ae208SSatish Balay #undef __FUNCT__
2174a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_SeqSBAIJ"
218268466fbSBarry Smith int MatGetSubMatrices_SeqSBAIJ(Mat A,int n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
21949b5e25fSSatish Balay {
22049b5e25fSSatish Balay   int ierr,i;
22149b5e25fSSatish Balay 
22249b5e25fSSatish Balay   PetscFunctionBegin;
22349b5e25fSSatish Balay   if (scall == MAT_INITIAL_MATRIX) {
22482502324SSatish Balay     ierr = PetscMalloc((n+1)*sizeof(Mat),B);CHKERRQ(ierr);
22549b5e25fSSatish Balay   }
22649b5e25fSSatish Balay 
22749b5e25fSSatish Balay   for (i=0; i<n; i++) {
22849b5e25fSSatish Balay     ierr = MatGetSubMatrix_SeqSBAIJ(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr);
22949b5e25fSSatish Balay   }
23049b5e25fSSatish Balay   PetscFunctionReturn(0);
23149b5e25fSSatish Balay }
23249b5e25fSSatish Balay 
23349b5e25fSSatish Balay /* -------------------------------------------------------*/
23449b5e25fSSatish Balay /* Should check that shapes of vectors and matrices match */
23549b5e25fSSatish Balay /* -------------------------------------------------------*/
236d9eff348SSatish Balay #include "petscblaslapack.h"
23749b5e25fSSatish Balay 
2384a2ae208SSatish Balay #undef __FUNCT__
2394a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqSBAIJ_1"
24049b5e25fSSatish Balay int MatMult_SeqSBAIJ_1(Mat A,Vec xx,Vec zz)
24149b5e25fSSatish Balay {
24249b5e25fSSatish Balay   Mat_SeqSBAIJ    *a = (Mat_SeqSBAIJ*)A->data;
24387828ca2SBarry Smith   PetscScalar     *x,*z,*xb,x1,zero=0.0;
24449b5e25fSSatish Balay   MatScalar       *v;
245831a3094SHong Zhang   int             mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,ierr,*ib,cval,j,jmin;
24649b5e25fSSatish Balay 
24749b5e25fSSatish Balay   PetscFunctionBegin;
24849b5e25fSSatish Balay   ierr = VecSet(&zero,zz);CHKERRQ(ierr);
2491ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2501ebc52fbSHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
25149b5e25fSSatish Balay 
25249b5e25fSSatish Balay   v  = a->a;
25349b5e25fSSatish Balay   xb = x;
25449b5e25fSSatish Balay 
25549b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
25649b5e25fSSatish Balay     n  = ai[1] - ai[0];  /* length of i_th row of A */
25749b5e25fSSatish Balay     x1 = xb[0];
25849b5e25fSSatish Balay     ib = aj + *ai;
259831a3094SHong Zhang     jmin = 0;
260831a3094SHong Zhang     if (*ib == i) {      /* (diag of A)*x */
261831a3094SHong Zhang       z[i] += *v++ * x[*ib++];
262831a3094SHong Zhang       jmin++;
263831a3094SHong Zhang     }
264831a3094SHong Zhang     for (j=jmin; j<n; j++) {
26549b5e25fSSatish Balay       cval    = *ib;
26649b5e25fSSatish Balay       z[cval] += *v * x1;      /* (strict lower triangular part of A)*x  */
26749b5e25fSSatish Balay       z[i] += *v++ * x[*ib++]; /* (strict upper triangular part of A)*x  */
26849b5e25fSSatish Balay     }
26949b5e25fSSatish Balay     xb++; ai++;
27049b5e25fSSatish Balay   }
27149b5e25fSSatish Balay 
2721ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2731ebc52fbSHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
2746c6c5352SBarry Smith   PetscLogFlops(2*(a->nz*2 - A->m) - A->m);  /* nz = (nz+m)/2 */
27549b5e25fSSatish Balay   PetscFunctionReturn(0);
27649b5e25fSSatish Balay }
27749b5e25fSSatish Balay 
2784a2ae208SSatish Balay #undef __FUNCT__
2794a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqSBAIJ_2"
28049b5e25fSSatish Balay int MatMult_SeqSBAIJ_2(Mat A,Vec xx,Vec zz)
28149b5e25fSSatish Balay {
28249b5e25fSSatish Balay   Mat_SeqSBAIJ     *a = (Mat_SeqSBAIJ*)A->data;
28387828ca2SBarry Smith   PetscScalar      *x,*z,*xb,x1,x2,zero=0.0;
28449b5e25fSSatish Balay   MatScalar        *v;
285831a3094SHong Zhang   int              mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,ierr,*ib,cval,j,jmin;
28649b5e25fSSatish Balay 
28749b5e25fSSatish Balay 
28849b5e25fSSatish Balay   PetscFunctionBegin;
28949b5e25fSSatish Balay   ierr = VecSet(&zero,zz);CHKERRQ(ierr);
2901ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2911ebc52fbSHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
29249b5e25fSSatish Balay 
29349b5e25fSSatish Balay   v     = a->a;
29449b5e25fSSatish Balay   xb = x;
29549b5e25fSSatish Balay 
29649b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
29749b5e25fSSatish Balay     n  = ai[1] - ai[0]; /* length of i_th block row of A */
29849b5e25fSSatish Balay     x1 = xb[0]; x2 = xb[1];
29949b5e25fSSatish Balay     ib = aj + *ai;
300831a3094SHong Zhang     jmin = 0;
3017fbae186SHong Zhang     if (*ib == i){     /* (diag of A)*x */
30249b5e25fSSatish Balay       z[2*i]   += v[0]*x1 + v[2]*x2;
30349b5e25fSSatish Balay       z[2*i+1] += v[2]*x1 + v[3]*x2;
304831a3094SHong Zhang       v += 4; jmin++;
3057fbae186SHong Zhang     }
306831a3094SHong Zhang     for (j=jmin; j<n; j++) {
30749b5e25fSSatish Balay       /* (strict lower triangular part of A)*x  */
30849b5e25fSSatish Balay       cval       = ib[j]*2;
30949b5e25fSSatish Balay       z[cval]     += v[0]*x1 + v[1]*x2;
31049b5e25fSSatish Balay       z[cval+1]   += v[2]*x1 + v[3]*x2;
31149b5e25fSSatish Balay       /* (strict upper triangular part of A)*x  */
31249b5e25fSSatish Balay       z[2*i]   += v[0]*x[cval] + v[2]*x[cval+1];
31349b5e25fSSatish Balay       z[2*i+1] += v[1]*x[cval] + v[3]*x[cval+1];
31449b5e25fSSatish Balay       v  += 4;
31549b5e25fSSatish Balay     }
31649b5e25fSSatish Balay     xb +=2; ai++;
31749b5e25fSSatish Balay   }
31849b5e25fSSatish Balay 
3191ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
3201ebc52fbSHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
3216c6c5352SBarry Smith   PetscLogFlops(8*(a->nz*2 - A->m) - A->m);
32249b5e25fSSatish Balay   PetscFunctionReturn(0);
32349b5e25fSSatish Balay }
32449b5e25fSSatish Balay 
3254a2ae208SSatish Balay #undef __FUNCT__
3264a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqSBAIJ_3"
32749b5e25fSSatish Balay int MatMult_SeqSBAIJ_3(Mat A,Vec xx,Vec zz)
32849b5e25fSSatish Balay {
32949b5e25fSSatish Balay   Mat_SeqSBAIJ  *a = (Mat_SeqSBAIJ*)A->data;
33087828ca2SBarry Smith   PetscScalar   *x,*z,*xb,x1,x2,x3,zero=0.0;
33149b5e25fSSatish Balay   MatScalar     *v;
332831a3094SHong Zhang   int           mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,ierr,*ib,cval,j,jmin;
33349b5e25fSSatish Balay 
33449b5e25fSSatish Balay 
33549b5e25fSSatish Balay   PetscFunctionBegin;
33649b5e25fSSatish Balay   ierr = VecSet(&zero,zz);CHKERRQ(ierr);
3371ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
3381ebc52fbSHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
33949b5e25fSSatish Balay 
34049b5e25fSSatish Balay   v     = a->a;
34149b5e25fSSatish Balay   xb = x;
34249b5e25fSSatish Balay 
34349b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
34449b5e25fSSatish Balay     n  = ai[1] - ai[0]; /* length of i_th block row of A */
34549b5e25fSSatish Balay     x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
34649b5e25fSSatish Balay     ib = aj + *ai;
347831a3094SHong Zhang     jmin = 0;
3487fbae186SHong Zhang     if (*ib == i){     /* (diag of A)*x */
34949b5e25fSSatish Balay       z[3*i]   += v[0]*x1 + v[3]*x2 + v[6]*x3;
35049b5e25fSSatish Balay       z[3*i+1] += v[3]*x1 + v[4]*x2 + v[7]*x3;
35149b5e25fSSatish Balay       z[3*i+2] += v[6]*x1 + v[7]*x2 + v[8]*x3;
352831a3094SHong Zhang       v += 9; jmin++;
3537fbae186SHong Zhang     }
354831a3094SHong Zhang     for (j=jmin; j<n; j++) {
35549b5e25fSSatish Balay       /* (strict lower triangular part of A)*x  */
35649b5e25fSSatish Balay       cval       = ib[j]*3;
35749b5e25fSSatish Balay       z[cval]     += v[0]*x1 + v[1]*x2 + v[2]*x3;
35849b5e25fSSatish Balay       z[cval+1]   += v[3]*x1 + v[4]*x2 + v[5]*x3;
35949b5e25fSSatish Balay       z[cval+2]   += v[6]*x1 + v[7]*x2 + v[8]*x3;
36049b5e25fSSatish Balay       /* (strict upper triangular part of A)*x  */
36149b5e25fSSatish Balay       z[3*i]   += v[0]*x[cval] + v[3]*x[cval+1]+ v[6]*x[cval+2];
36249b5e25fSSatish Balay       z[3*i+1] += v[1]*x[cval] + v[4]*x[cval+1]+ v[7]*x[cval+2];
36349b5e25fSSatish Balay       z[3*i+2] += v[2]*x[cval] + v[5]*x[cval+1]+ v[8]*x[cval+2];
36449b5e25fSSatish Balay       v  += 9;
36549b5e25fSSatish Balay     }
36649b5e25fSSatish Balay     xb +=3; ai++;
36749b5e25fSSatish Balay   }
36849b5e25fSSatish Balay 
3691ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
3701ebc52fbSHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
3716c6c5352SBarry Smith   PetscLogFlops(18*(a->nz*2 - A->m) - A->m);
37249b5e25fSSatish Balay   PetscFunctionReturn(0);
37349b5e25fSSatish Balay }
37449b5e25fSSatish Balay 
3754a2ae208SSatish Balay #undef __FUNCT__
3764a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqSBAIJ_4"
37749b5e25fSSatish Balay int MatMult_SeqSBAIJ_4(Mat A,Vec xx,Vec zz)
37849b5e25fSSatish Balay {
37949b5e25fSSatish Balay   Mat_SeqSBAIJ     *a = (Mat_SeqSBAIJ*)A->data;
38087828ca2SBarry Smith   PetscScalar      *x,*z,*xb,x1,x2,x3,x4,zero=0.0;
38149b5e25fSSatish Balay   MatScalar        *v;
382831a3094SHong Zhang   int              mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,ierr,*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);
4236c6c5352SBarry Smith   PetscLogFlops(32*(a->nz*2 - A->m) - A->m);
42449b5e25fSSatish Balay   PetscFunctionReturn(0);
42549b5e25fSSatish Balay }
42649b5e25fSSatish Balay 
4274a2ae208SSatish Balay #undef __FUNCT__
4284a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqSBAIJ_5"
42949b5e25fSSatish Balay int 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;
434831a3094SHong Zhang   int             mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,ierr,*ib,cval,j,jmin;
43549b5e25fSSatish Balay 
43649b5e25fSSatish Balay   PetscFunctionBegin;
43749b5e25fSSatish Balay   ierr = VecSet(&zero,zz);CHKERRQ(ierr);
4381ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
4391ebc52fbSHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
44049b5e25fSSatish Balay 
44149b5e25fSSatish Balay   v     = a->a;
44249b5e25fSSatish Balay   xb = x;
44349b5e25fSSatish Balay 
44449b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
44549b5e25fSSatish Balay     n  = ai[1] - ai[0]; /* length of i_th block row of A */
44649b5e25fSSatish Balay     x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5=xb[4];
44749b5e25fSSatish Balay     ib = aj + *ai;
448831a3094SHong Zhang     jmin = 0;
4497fbae186SHong Zhang     if (*ib == i){      /* (diag of A)*x */
45049b5e25fSSatish Balay       z[5*i]   += v[0]*x1  + v[5]*x2 + v[10]*x3 + v[15]*x4+ v[20]*x5;
45149b5e25fSSatish Balay       z[5*i+1] += v[5]*x1  + v[6]*x2 + v[11]*x3 + v[16]*x4+ v[21]*x5;
45249b5e25fSSatish Balay       z[5*i+2] += v[10]*x1 +v[11]*x2 + v[12]*x3 + v[17]*x4+ v[22]*x5;
45349b5e25fSSatish Balay       z[5*i+3] += v[15]*x1 +v[16]*x2 + v[17]*x3 + v[18]*x4+ v[23]*x5;
45449b5e25fSSatish Balay       z[5*i+4] += v[20]*x1 +v[21]*x2 + v[22]*x3 + v[23]*x4+ v[24]*x5;
455831a3094SHong Zhang       v += 25; jmin++;
4567fbae186SHong Zhang     }
457831a3094SHong Zhang     for (j=jmin; j<n; j++) {
45849b5e25fSSatish Balay       /* (strict lower triangular part of A)*x  */
45949b5e25fSSatish Balay       cval       = ib[j]*5;
46049b5e25fSSatish Balay       z[cval]     += v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4 + v[4]*x5;
46149b5e25fSSatish Balay       z[cval+1]   += v[5]*x1 + v[6]*x2 + v[7]*x3 + v[8]*x4 + v[9]*x5;
46249b5e25fSSatish Balay       z[cval+2]   += v[10]*x1 + v[11]*x2 + v[12]*x3 + v[13]*x4+ v[14]*x5;
46349b5e25fSSatish Balay       z[cval+3]   += v[15]*x1 + v[16]*x2 + v[17]*x3 + v[18]*x4+ v[19]*x5;
46449b5e25fSSatish Balay       z[cval+4]   += v[20]*x1 + v[21]*x2 + v[22]*x3 + v[23]*x4+ v[24]*x5;
46549b5e25fSSatish Balay       /* (strict upper triangular part of A)*x  */
46649b5e25fSSatish 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];
46749b5e25fSSatish 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];
46849b5e25fSSatish 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];
46949b5e25fSSatish 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];
47049b5e25fSSatish 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];
47149b5e25fSSatish Balay       v  += 25;
47249b5e25fSSatish Balay     }
47349b5e25fSSatish Balay     xb +=5; ai++;
47449b5e25fSSatish Balay   }
47549b5e25fSSatish Balay 
4761ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
4771ebc52fbSHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
4786c6c5352SBarry Smith   PetscLogFlops(50*(a->nz*2 - A->m) - A->m);
47949b5e25fSSatish Balay   PetscFunctionReturn(0);
48049b5e25fSSatish Balay }
48149b5e25fSSatish Balay 
48249b5e25fSSatish Balay 
4834a2ae208SSatish Balay #undef __FUNCT__
4844a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqSBAIJ_6"
48549b5e25fSSatish Balay int MatMult_SeqSBAIJ_6(Mat A,Vec xx,Vec zz)
48649b5e25fSSatish Balay {
48749b5e25fSSatish Balay   Mat_SeqSBAIJ    *a = (Mat_SeqSBAIJ*)A->data;
48887828ca2SBarry Smith   PetscScalar     *x,*z,*xb,x1,x2,x3,x4,x5,x6,zero=0.0;
48949b5e25fSSatish Balay   MatScalar       *v;
490831a3094SHong Zhang   int             mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,ierr,*ib,cval,j,jmin;
49149b5e25fSSatish Balay 
49249b5e25fSSatish Balay   PetscFunctionBegin;
49349b5e25fSSatish Balay   ierr = VecSet(&zero,zz);CHKERRQ(ierr);
4941ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
4951ebc52fbSHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
49649b5e25fSSatish Balay 
49749b5e25fSSatish Balay   v     = a->a;
49849b5e25fSSatish Balay   xb = x;
49949b5e25fSSatish Balay 
50049b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
50149b5e25fSSatish Balay     n  = ai[1] - ai[0]; /* length of i_th block row of A */
50249b5e25fSSatish Balay     x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5=xb[4]; x6=xb[5];
50349b5e25fSSatish Balay     ib = aj + *ai;
504831a3094SHong Zhang     jmin = 0;
5057fbae186SHong Zhang     if (*ib == i){      /* (diag of A)*x */
50649b5e25fSSatish Balay       z[6*i]   += v[0]*x1  + v[6]*x2 + v[12]*x3 + v[18]*x4+ v[24]*x5 + v[30]*x6;
50749b5e25fSSatish Balay       z[6*i+1] += v[6]*x1  + v[7]*x2 + v[13]*x3 + v[19]*x4+ v[25]*x5 + v[31]*x6;
50849b5e25fSSatish Balay       z[6*i+2] += v[12]*x1 +v[13]*x2 + v[14]*x3 + v[20]*x4+ v[26]*x5 + v[32]*x6;
50949b5e25fSSatish Balay       z[6*i+3] += v[18]*x1 +v[19]*x2 + v[20]*x3 + v[21]*x4+ v[27]*x5 + v[33]*x6;
51049b5e25fSSatish Balay       z[6*i+4] += v[24]*x1 +v[25]*x2 + v[26]*x3 + v[27]*x4+ v[28]*x5 + v[34]*x6;
51149b5e25fSSatish Balay       z[6*i+5] += v[30]*x1 +v[31]*x2 + v[32]*x3 + v[33]*x4+ v[34]*x5 + v[35]*x6;
512831a3094SHong Zhang       v += 36; jmin++;
5137fbae186SHong Zhang     }
514831a3094SHong Zhang     for (j=jmin; j<n; j++) {
51549b5e25fSSatish Balay       /* (strict lower triangular part of A)*x  */
51649b5e25fSSatish Balay       cval       = ib[j]*6;
51749b5e25fSSatish Balay       z[cval]   += v[0]*x1  + v[1]*x2 + v[2]*x3 + v[3]*x4+ v[4]*x5 + v[5]*x6;
51849b5e25fSSatish Balay       z[cval+1] += v[6]*x1  + v[7]*x2 + v[8]*x3 + v[9]*x4+ v[10]*x5 + v[11]*x6;
51949b5e25fSSatish Balay       z[cval+2] += v[12]*x1  + v[13]*x2 + v[14]*x3 + v[15]*x4+ v[16]*x5 + v[17]*x6;
52049b5e25fSSatish Balay       z[cval+3] += v[18]*x1  + v[19]*x2 + v[20]*x3 + v[21]*x4+ v[22]*x5 + v[23]*x6;
52149b5e25fSSatish Balay       z[cval+4] += v[24]*x1  + v[25]*x2 + v[26]*x3 + v[27]*x4+ v[28]*x5 + v[29]*x6;
52249b5e25fSSatish Balay       z[cval+5] += v[30]*x1  + v[31]*x2 + v[32]*x3 + v[33]*x4+ v[34]*x5 + v[35]*x6;
52349b5e25fSSatish Balay       /* (strict upper triangular part of A)*x  */
52449b5e25fSSatish 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];
52549b5e25fSSatish 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];
52649b5e25fSSatish 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];
52749b5e25fSSatish 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];
52849b5e25fSSatish 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];
52949b5e25fSSatish 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];
53049b5e25fSSatish Balay       v  += 36;
53149b5e25fSSatish Balay     }
53249b5e25fSSatish Balay     xb +=6; ai++;
53349b5e25fSSatish Balay   }
53449b5e25fSSatish Balay 
5351ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
5361ebc52fbSHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
5376c6c5352SBarry Smith   PetscLogFlops(72*(a->nz*2 - A->m) - A->m);
53849b5e25fSSatish Balay   PetscFunctionReturn(0);
53949b5e25fSSatish Balay }
5404a2ae208SSatish Balay #undef __FUNCT__
5414a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqSBAIJ_7"
54249b5e25fSSatish Balay int MatMult_SeqSBAIJ_7(Mat A,Vec xx,Vec zz)
54349b5e25fSSatish Balay {
54449b5e25fSSatish Balay   Mat_SeqSBAIJ    *a = (Mat_SeqSBAIJ*)A->data;
54587828ca2SBarry Smith   PetscScalar     *x,*z,*xb,x1,x2,x3,x4,x5,x6,x7,zero=0.0;
54649b5e25fSSatish Balay   MatScalar       *v;
547831a3094SHong Zhang   int             mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,ierr,*ib,cval,j,jmin;
54849b5e25fSSatish Balay 
54949b5e25fSSatish Balay   PetscFunctionBegin;
55049b5e25fSSatish Balay   ierr = VecSet(&zero,zz);CHKERRQ(ierr);
5511ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
5521ebc52fbSHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
55349b5e25fSSatish Balay 
55449b5e25fSSatish Balay   v     = a->a;
55549b5e25fSSatish Balay   xb = x;
55649b5e25fSSatish Balay 
55749b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
55849b5e25fSSatish Balay     n  = ai[1] - ai[0]; /* length of i_th block row of A */
55949b5e25fSSatish Balay     x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5=xb[4]; x6=xb[5]; x7=xb[6];
56049b5e25fSSatish Balay     ib = aj + *ai;
561831a3094SHong Zhang     jmin = 0;
5627fbae186SHong Zhang     if (*ib == i){      /* (diag of A)*x */
56349b5e25fSSatish 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;
56449b5e25fSSatish 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;
56549b5e25fSSatish 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;
56649b5e25fSSatish 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;
56749b5e25fSSatish 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;
56849b5e25fSSatish 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;
56949b5e25fSSatish 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;
570831a3094SHong Zhang       v += 49; jmin++;
5717fbae186SHong Zhang     }
572831a3094SHong Zhang     for (j=jmin; j<n; j++) {
57349b5e25fSSatish Balay       /* (strict lower triangular part of A)*x  */
57449b5e25fSSatish Balay       cval       = ib[j]*7;
57549b5e25fSSatish Balay       z[cval]   += v[0]*x1  + v[1]*x2 + v[2]*x3 + v[3]*x4+ v[4]*x5 + v[5]*x6+ v[6]*x7;
57649b5e25fSSatish 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;
57749b5e25fSSatish 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;
57849b5e25fSSatish 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;
57949b5e25fSSatish 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;
58049b5e25fSSatish 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;
58149b5e25fSSatish 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;
58249b5e25fSSatish Balay       /* (strict upper triangular part of A)*x  */
58349b5e25fSSatish 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];
58449b5e25fSSatish 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];
58549b5e25fSSatish 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];
58649b5e25fSSatish 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];
58749b5e25fSSatish 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];
58849b5e25fSSatish 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];
58949b5e25fSSatish 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];
59049b5e25fSSatish Balay       v  += 49;
59149b5e25fSSatish Balay     }
59249b5e25fSSatish Balay     xb +=7; ai++;
59349b5e25fSSatish Balay   }
5941ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
5951ebc52fbSHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
5966c6c5352SBarry Smith   PetscLogFlops(98*(a->nz*2 - A->m) - A->m);
59749b5e25fSSatish Balay   PetscFunctionReturn(0);
59849b5e25fSSatish Balay }
59949b5e25fSSatish Balay 
60049b5e25fSSatish Balay /*
60149b5e25fSSatish Balay     This will not work with MatScalar == float because it calls the BLAS
60249b5e25fSSatish Balay */
6034a2ae208SSatish Balay #undef __FUNCT__
6044a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqSBAIJ_N"
60549b5e25fSSatish Balay int MatMult_SeqSBAIJ_N(Mat A,Vec xx,Vec zz)
60649b5e25fSSatish Balay {
60749b5e25fSSatish Balay   Mat_SeqSBAIJ    *a = (Mat_SeqSBAIJ*)A->data;
60887828ca2SBarry Smith   PetscScalar     *x,*x_ptr,*z,*z_ptr,*xb,*zb,*work,*workt,zero=0.0;
6090b60a74dSHong Zhang   MatScalar       *v;
610066653e3SSatish Balay   int             ierr,mbs=a->mbs,i,*idx,*aj,*ii,bs=a->bs,j,n,bs2=a->bs2;
6110b60a74dSHong Zhang   int             ncols,k;
61249b5e25fSSatish Balay 
61349b5e25fSSatish Balay   PetscFunctionBegin;
61449b5e25fSSatish Balay   ierr = VecSet(&zero,zz);CHKERRQ(ierr);
6151ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr); x_ptr=x;
6161ebc52fbSHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr); z_ptr=z;
61749b5e25fSSatish Balay 
61849b5e25fSSatish Balay   aj   = a->j;
61949b5e25fSSatish Balay   v    = a->a;
62049b5e25fSSatish Balay   ii   = a->i;
62149b5e25fSSatish Balay 
62249b5e25fSSatish Balay   if (!a->mult_work) {
62387828ca2SBarry Smith     ierr = PetscMalloc((A->m+1)*sizeof(PetscScalar),&a->mult_work);CHKERRQ(ierr);
62449b5e25fSSatish Balay   }
62549b5e25fSSatish Balay   work = a->mult_work;
62649b5e25fSSatish Balay 
62749b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
62849b5e25fSSatish Balay     n     = ii[1] - ii[0]; ncols = n*bs;
62949b5e25fSSatish Balay     workt = work; idx=aj+ii[0];
63049b5e25fSSatish Balay 
63149b5e25fSSatish Balay     /* upper triangular part */
63249b5e25fSSatish Balay     for (j=0; j<n; j++) {
63349b5e25fSSatish Balay       xb = x_ptr + bs*(*idx++);
63449b5e25fSSatish Balay       for (k=0; k<bs; k++) workt[k] = xb[k];
63549b5e25fSSatish Balay       workt += bs;
63649b5e25fSSatish Balay     }
63749b5e25fSSatish Balay     /* z(i*bs:(i+1)*bs-1) += A(i,:)*x */
63849b5e25fSSatish Balay     Kernel_w_gets_w_plus_Ar_times_v(bs,ncols,work,v,z);
63949b5e25fSSatish Balay 
64049b5e25fSSatish Balay     /* strict lower triangular part */
641831a3094SHong Zhang     idx = aj+ii[0];
642831a3094SHong Zhang     if (*idx == i){
64396b9376eSHong Zhang       ncols -= bs; v += bs2; idx++; n--;
644831a3094SHong Zhang     }
64596b9376eSHong Zhang 
64649b5e25fSSatish Balay     if (ncols > 0){
64749b5e25fSSatish Balay       workt = work;
64887828ca2SBarry Smith       ierr  = PetscMemzero(workt,ncols*sizeof(PetscScalar));CHKERRQ(ierr);
649831a3094SHong Zhang       Kernel_w_gets_w_plus_trans_Ar_times_v(bs,ncols,x,v,workt);
650831a3094SHong Zhang       for (j=0; j<n; j++) {
651831a3094SHong Zhang         zb = z_ptr + bs*(*idx++);
65249b5e25fSSatish Balay         for (k=0; k<bs; k++) zb[k] += workt[k] ;
65349b5e25fSSatish Balay         workt += bs;
65449b5e25fSSatish Balay       }
65549b5e25fSSatish Balay     }
65649b5e25fSSatish Balay     x += bs; v += n*bs2; z += bs; ii++;
65749b5e25fSSatish Balay   }
65849b5e25fSSatish Balay 
6591ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
6601ebc52fbSHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
6616c6c5352SBarry Smith   PetscLogFlops(2*(a->nz*2 - A->m)*bs2 - A->m);
66249b5e25fSSatish Balay   PetscFunctionReturn(0);
66349b5e25fSSatish Balay }
66449b5e25fSSatish Balay 
6654a2ae208SSatish Balay #undef __FUNCT__
6664a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqSBAIJ_1"
66749b5e25fSSatish Balay int MatMultAdd_SeqSBAIJ_1(Mat A,Vec xx,Vec yy,Vec zz)
66849b5e25fSSatish Balay {
66949b5e25fSSatish Balay   Mat_SeqSBAIJ    *a = (Mat_SeqSBAIJ*)A->data;
67087828ca2SBarry Smith   PetscScalar     *x,*y,*z,*xb,x1;
67149b5e25fSSatish Balay   MatScalar       *v;
672831a3094SHong Zhang   int             mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,ierr,*ib,cval,j,jmin;
67349b5e25fSSatish Balay 
67449b5e25fSSatish Balay   PetscFunctionBegin;
6751ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
67649b5e25fSSatish Balay   if (yy != xx) {
6771ebc52fbSHong Zhang     ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
67849b5e25fSSatish Balay   } else {
67949b5e25fSSatish Balay     y = x;
68049b5e25fSSatish Balay   }
68149b5e25fSSatish Balay   if (zz != yy) {
682a9d4b620SHong Zhang     /* ierr = VecCopy(yy,zz);CHKERRQ(ierr); */
6831ebc52fbSHong Zhang     ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
684a9d4b620SHong Zhang     ierr = PetscMemcpy(z,y,yy->n*sizeof(PetscScalar));CHKERRQ(ierr);
68549b5e25fSSatish Balay   } else {
68649b5e25fSSatish Balay     z = y;
68749b5e25fSSatish Balay   }
68849b5e25fSSatish Balay 
68949b5e25fSSatish Balay   v  = a->a;
69049b5e25fSSatish Balay   xb = x;
69149b5e25fSSatish Balay 
69249b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
69349b5e25fSSatish Balay     n  = ai[1] - ai[0];  /* length of i_th row of A */
69449b5e25fSSatish Balay     x1 = xb[0];
69549b5e25fSSatish Balay     ib = aj + *ai;
696831a3094SHong Zhang     jmin = 0;
697831a3094SHong Zhang     if (*ib == i) {            /* (diag of A)*x */
698831a3094SHong Zhang       z[i] += *v++ * x[*ib++]; jmin++;
699831a3094SHong Zhang     }
700831a3094SHong Zhang     for (j=jmin; j<n; j++) {
70149b5e25fSSatish Balay       cval    = *ib;
70249b5e25fSSatish Balay       z[cval] += *v * x1;      /* (strict lower triangular part of A)*x  */
70349b5e25fSSatish Balay       z[i] += *v++ * x[*ib++]; /* (strict upper triangular part of A)*x  */
70449b5e25fSSatish Balay     }
70549b5e25fSSatish Balay     xb++; ai++;
70649b5e25fSSatish Balay   }
70749b5e25fSSatish Balay 
7081ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
7091ebc52fbSHong Zhang   if (yy != xx) ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
7101ebc52fbSHong Zhang   if (zz != yy) ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
71149b5e25fSSatish Balay 
7126c6c5352SBarry Smith   PetscLogFlops(2*(a->nz*2 - A->m));
71349b5e25fSSatish Balay   PetscFunctionReturn(0);
71449b5e25fSSatish Balay }
71549b5e25fSSatish Balay 
7164a2ae208SSatish Balay #undef __FUNCT__
7174a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqSBAIJ_2"
71849b5e25fSSatish Balay int MatMultAdd_SeqSBAIJ_2(Mat A,Vec xx,Vec yy,Vec zz)
71949b5e25fSSatish Balay {
72049b5e25fSSatish Balay   Mat_SeqSBAIJ    *a = (Mat_SeqSBAIJ*)A->data;
72187828ca2SBarry Smith   PetscScalar     *x,*y,*z,*xb,x1,x2;
72249b5e25fSSatish Balay   MatScalar       *v;
723831a3094SHong Zhang   int             mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,ierr,*ib,cval,j,jmin;
72449b5e25fSSatish Balay 
72549b5e25fSSatish Balay   PetscFunctionBegin;
7261ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
72749b5e25fSSatish Balay   if (yy != xx) {
7281ebc52fbSHong Zhang     ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
72949b5e25fSSatish Balay   } else {
73049b5e25fSSatish Balay     y = x;
73149b5e25fSSatish Balay   }
73249b5e25fSSatish Balay   if (zz != yy) {
733a9d4b620SHong Zhang     /* ierr = VecCopy(yy,zz);CHKERRQ(ierr); */
7341ebc52fbSHong Zhang     ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
735a9d4b620SHong Zhang     ierr = PetscMemcpy(z,y,yy->n*sizeof(PetscScalar));CHKERRQ(ierr);
73649b5e25fSSatish Balay   } else {
73749b5e25fSSatish Balay     z = y;
73849b5e25fSSatish Balay   }
73949b5e25fSSatish Balay 
74049b5e25fSSatish Balay   v     = a->a;
74149b5e25fSSatish Balay   xb = x;
74249b5e25fSSatish Balay 
74349b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
74449b5e25fSSatish Balay     n  = ai[1] - ai[0]; /* length of i_th block row of A */
74549b5e25fSSatish Balay     x1 = xb[0]; x2 = xb[1];
74649b5e25fSSatish Balay     ib = aj + *ai;
747831a3094SHong Zhang     jmin = 0;
7487fbae186SHong Zhang     if (*ib == i){      /* (diag of A)*x */
74949b5e25fSSatish Balay       z[2*i]   += v[0]*x1 + v[2]*x2;
75049b5e25fSSatish Balay       z[2*i+1] += v[2]*x1 + v[3]*x2;
751831a3094SHong Zhang       v += 4; jmin++;
7527fbae186SHong Zhang     }
753831a3094SHong Zhang     for (j=jmin; j<n; j++) {
75449b5e25fSSatish Balay       /* (strict lower triangular part of A)*x  */
75549b5e25fSSatish Balay       cval       = ib[j]*2;
75649b5e25fSSatish Balay       z[cval]     += v[0]*x1 + v[1]*x2;
75749b5e25fSSatish Balay       z[cval+1]   += v[2]*x1 + v[3]*x2;
75849b5e25fSSatish Balay       /* (strict upper triangular part of A)*x  */
75949b5e25fSSatish Balay       z[2*i]   += v[0]*x[cval] + v[2]*x[cval+1];
76049b5e25fSSatish Balay       z[2*i+1] += v[1]*x[cval] + v[3]*x[cval+1];
76149b5e25fSSatish Balay       v  += 4;
76249b5e25fSSatish Balay     }
76349b5e25fSSatish Balay     xb +=2; ai++;
76449b5e25fSSatish Balay   }
76549b5e25fSSatish Balay 
7661ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
7671ebc52fbSHong Zhang   if (yy != xx) ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
7681ebc52fbSHong Zhang   if (zz != yy) ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
76949b5e25fSSatish Balay 
7706c6c5352SBarry Smith   PetscLogFlops(4*(a->nz*2 - A->m));
77149b5e25fSSatish Balay   PetscFunctionReturn(0);
77249b5e25fSSatish Balay }
77349b5e25fSSatish Balay 
7744a2ae208SSatish Balay #undef __FUNCT__
7754a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqSBAIJ_3"
77649b5e25fSSatish Balay int MatMultAdd_SeqSBAIJ_3(Mat A,Vec xx,Vec yy,Vec zz)
77749b5e25fSSatish Balay {
77849b5e25fSSatish Balay   Mat_SeqSBAIJ    *a = (Mat_SeqSBAIJ*)A->data;
77987828ca2SBarry Smith   PetscScalar     *x,*y,*z,*xb,x1,x2,x3;
78049b5e25fSSatish Balay   MatScalar       *v;
781831a3094SHong Zhang   int             mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,ierr,*ib,cval,j,jmin;
78249b5e25fSSatish Balay 
78349b5e25fSSatish Balay   PetscFunctionBegin;
7841ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
78549b5e25fSSatish Balay   if (yy != xx) {
7861ebc52fbSHong Zhang     ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
78749b5e25fSSatish Balay   } else {
78849b5e25fSSatish Balay     y = x;
78949b5e25fSSatish Balay   }
79049b5e25fSSatish Balay   if (zz != yy) {
791a9d4b620SHong Zhang     /* ierr = VecCopy(yy,zz);CHKERRQ(ierr); */
7921ebc52fbSHong Zhang     ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
793a9d4b620SHong Zhang     ierr = PetscMemcpy(z,y,yy->n*sizeof(PetscScalar));CHKERRQ(ierr);
79449b5e25fSSatish Balay   } else {
79549b5e25fSSatish Balay     z = y;
79649b5e25fSSatish Balay   }
79749b5e25fSSatish Balay 
79849b5e25fSSatish Balay   v     = a->a;
79949b5e25fSSatish Balay   xb = x;
80049b5e25fSSatish Balay 
80149b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
80249b5e25fSSatish Balay     n  = ai[1] - ai[0]; /* length of i_th block row of A */
80349b5e25fSSatish Balay     x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
80449b5e25fSSatish Balay     ib = aj + *ai;
805831a3094SHong Zhang     jmin = 0;
8067fbae186SHong Zhang     if (*ib == i){     /* (diag of A)*x */
80749b5e25fSSatish Balay      z[3*i]   += v[0]*x1 + v[3]*x2 + v[6]*x3;
80849b5e25fSSatish Balay      z[3*i+1] += v[3]*x1 + v[4]*x2 + v[7]*x3;
80949b5e25fSSatish Balay      z[3*i+2] += v[6]*x1 + v[7]*x2 + v[8]*x3;
810831a3094SHong Zhang      v += 9; jmin++;
8117fbae186SHong Zhang     }
812831a3094SHong Zhang     for (j=jmin; j<n; j++) {
81349b5e25fSSatish Balay       /* (strict lower triangular part of A)*x  */
81449b5e25fSSatish Balay       cval       = ib[j]*3;
81549b5e25fSSatish Balay       z[cval]     += v[0]*x1 + v[1]*x2 + v[2]*x3;
81649b5e25fSSatish Balay       z[cval+1]   += v[3]*x1 + v[4]*x2 + v[5]*x3;
81749b5e25fSSatish Balay       z[cval+2]   += v[6]*x1 + v[7]*x2 + v[8]*x3;
81849b5e25fSSatish Balay       /* (strict upper triangular part of A)*x  */
81949b5e25fSSatish Balay       z[3*i]   += v[0]*x[cval] + v[3]*x[cval+1]+ v[6]*x[cval+2];
82049b5e25fSSatish Balay       z[3*i+1] += v[1]*x[cval] + v[4]*x[cval+1]+ v[7]*x[cval+2];
82149b5e25fSSatish Balay       z[3*i+2] += v[2]*x[cval] + v[5]*x[cval+1]+ v[8]*x[cval+2];
82249b5e25fSSatish Balay       v  += 9;
82349b5e25fSSatish Balay     }
82449b5e25fSSatish Balay     xb +=3; ai++;
82549b5e25fSSatish Balay   }
82649b5e25fSSatish Balay 
8271ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
8281ebc52fbSHong Zhang   if (yy != xx) ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
8291ebc52fbSHong Zhang   if (zz != yy) ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
83049b5e25fSSatish Balay 
8316c6c5352SBarry Smith   PetscLogFlops(18*(a->nz*2 - A->m));
83249b5e25fSSatish Balay   PetscFunctionReturn(0);
83349b5e25fSSatish Balay }
83449b5e25fSSatish Balay 
8354a2ae208SSatish Balay #undef __FUNCT__
8364a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqSBAIJ_4"
83749b5e25fSSatish Balay int MatMultAdd_SeqSBAIJ_4(Mat A,Vec xx,Vec yy,Vec zz)
83849b5e25fSSatish Balay {
83949b5e25fSSatish Balay   Mat_SeqSBAIJ    *a = (Mat_SeqSBAIJ*)A->data;
84087828ca2SBarry Smith   PetscScalar     *x,*y,*z,*xb,x1,x2,x3,x4;
84149b5e25fSSatish Balay   MatScalar       *v;
842831a3094SHong Zhang   int             mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,ierr,*ib,cval,j,jmin;
84349b5e25fSSatish Balay 
84449b5e25fSSatish Balay   PetscFunctionBegin;
8451ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
84649b5e25fSSatish Balay   if (yy != xx) {
8471ebc52fbSHong Zhang     ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
84849b5e25fSSatish Balay   } else {
84949b5e25fSSatish Balay     y = x;
85049b5e25fSSatish Balay   }
85149b5e25fSSatish Balay   if (zz != yy) {
852a9d4b620SHong Zhang     /* ierr = VecCopy(yy,zz);CHKERRQ(ierr); */
8531ebc52fbSHong Zhang     ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
854a9d4b620SHong Zhang     ierr = PetscMemcpy(z,y,yy->n*sizeof(PetscScalar));CHKERRQ(ierr);
85549b5e25fSSatish Balay   } else {
85649b5e25fSSatish Balay     z = y;
85749b5e25fSSatish Balay   }
85849b5e25fSSatish Balay 
85949b5e25fSSatish Balay   v     = a->a;
86049b5e25fSSatish Balay   xb = x;
86149b5e25fSSatish Balay 
86249b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
86349b5e25fSSatish Balay     n  = ai[1] - ai[0]; /* length of i_th block row of A */
86449b5e25fSSatish Balay     x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
86549b5e25fSSatish Balay     ib = aj + *ai;
866831a3094SHong Zhang     jmin = 0;
8677fbae186SHong Zhang     if (*ib == i){      /* (diag of A)*x */
86849b5e25fSSatish Balay       z[4*i]   += v[0]*x1 + v[4]*x2 +  v[8]*x3 + v[12]*x4;
86949b5e25fSSatish Balay       z[4*i+1] += v[4]*x1 + v[5]*x2 +  v[9]*x3 + v[13]*x4;
87049b5e25fSSatish Balay       z[4*i+2] += v[8]*x1 + v[9]*x2 + v[10]*x3 + v[14]*x4;
87149b5e25fSSatish Balay       z[4*i+3] += v[12]*x1+ v[13]*x2+ v[14]*x3 + v[15]*x4;
872831a3094SHong Zhang       v += 16; jmin++;
8737fbae186SHong Zhang     }
874831a3094SHong Zhang     for (j=jmin; j<n; j++) {
87549b5e25fSSatish Balay       /* (strict lower triangular part of A)*x  */
87649b5e25fSSatish Balay       cval       = ib[j]*4;
87749b5e25fSSatish Balay       z[cval]     += v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4;
87849b5e25fSSatish Balay       z[cval+1]   += v[4]*x1 + v[5]*x2 + v[6]*x3 + v[7]*x4;
87949b5e25fSSatish Balay       z[cval+2]   += v[8]*x1 + v[9]*x2 + v[10]*x3 + v[11]*x4;
88049b5e25fSSatish Balay       z[cval+3]   += v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4;
88149b5e25fSSatish Balay       /* (strict upper triangular part of A)*x  */
88249b5e25fSSatish Balay       z[4*i]   += v[0]*x[cval] + v[4]*x[cval+1]+ v[8]*x[cval+2] + v[12]*x[cval+3];
88349b5e25fSSatish Balay       z[4*i+1] += v[1]*x[cval] + v[5]*x[cval+1]+ v[9]*x[cval+2] + v[13]*x[cval+3];
88449b5e25fSSatish Balay       z[4*i+2] += v[2]*x[cval] + v[6]*x[cval+1]+ v[10]*x[cval+2]+ v[14]*x[cval+3];
88549b5e25fSSatish Balay       z[4*i+3] += v[3]*x[cval] + v[7]*x[cval+1]+ v[11]*x[cval+2]+ v[15]*x[cval+3];
88649b5e25fSSatish Balay       v  += 16;
88749b5e25fSSatish Balay     }
88849b5e25fSSatish Balay     xb +=4; ai++;
88949b5e25fSSatish Balay   }
89049b5e25fSSatish Balay 
8911ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
8921ebc52fbSHong Zhang   if (yy != xx) ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
8931ebc52fbSHong Zhang   if (zz != yy) ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
89449b5e25fSSatish Balay 
8956c6c5352SBarry Smith   PetscLogFlops(32*(a->nz*2 - A->m));
89649b5e25fSSatish Balay   PetscFunctionReturn(0);
89749b5e25fSSatish Balay }
89849b5e25fSSatish Balay 
8994a2ae208SSatish Balay #undef __FUNCT__
9004a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqSBAIJ_5"
90149b5e25fSSatish Balay int MatMultAdd_SeqSBAIJ_5(Mat A,Vec xx,Vec yy,Vec zz)
90249b5e25fSSatish Balay {
90349b5e25fSSatish Balay   Mat_SeqSBAIJ    *a = (Mat_SeqSBAIJ*)A->data;
90487828ca2SBarry Smith   PetscScalar     *x,*y,*z,*xb,x1,x2,x3,x4,x5;
90549b5e25fSSatish Balay   MatScalar       *v;
906831a3094SHong Zhang   int             mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,ierr,*ib,cval,j,jmin;
90749b5e25fSSatish Balay 
90849b5e25fSSatish Balay   PetscFunctionBegin;
9091ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
91049b5e25fSSatish Balay   if (yy != xx) {
9111ebc52fbSHong Zhang     ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
91249b5e25fSSatish Balay   } else {
91349b5e25fSSatish Balay     y = x;
91449b5e25fSSatish Balay   }
91549b5e25fSSatish Balay   if (zz != yy) {
916a9d4b620SHong Zhang     /* ierr = VecCopy(yy,zz);CHKERRQ(ierr); */
9171ebc52fbSHong Zhang     ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
918a9d4b620SHong Zhang     ierr = PetscMemcpy(z,y,yy->n*sizeof(PetscScalar));CHKERRQ(ierr);
91949b5e25fSSatish Balay   } else {
92049b5e25fSSatish Balay     z = y;
92149b5e25fSSatish Balay   }
92249b5e25fSSatish Balay 
92349b5e25fSSatish Balay   v     = a->a;
92449b5e25fSSatish Balay   xb = x;
92549b5e25fSSatish Balay 
92649b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
92749b5e25fSSatish Balay     n  = ai[1] - ai[0]; /* length of i_th block row of A */
92849b5e25fSSatish Balay     x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5=xb[4];
92949b5e25fSSatish Balay     ib = aj + *ai;
930831a3094SHong Zhang     jmin = 0;
9317fbae186SHong Zhang     if (*ib == i){      /* (diag of A)*x */
93249b5e25fSSatish Balay       z[5*i]   += v[0]*x1  + v[5]*x2 + v[10]*x3 + v[15]*x4+ v[20]*x5;
93349b5e25fSSatish Balay       z[5*i+1] += v[5]*x1  + v[6]*x2 + v[11]*x3 + v[16]*x4+ v[21]*x5;
93449b5e25fSSatish Balay       z[5*i+2] += v[10]*x1 +v[11]*x2 + v[12]*x3 + v[17]*x4+ v[22]*x5;
93549b5e25fSSatish Balay       z[5*i+3] += v[15]*x1 +v[16]*x2 + v[17]*x3 + v[18]*x4+ v[23]*x5;
93649b5e25fSSatish Balay       z[5*i+4] += v[20]*x1 +v[21]*x2 + v[22]*x3 + v[23]*x4+ v[24]*x5;
937831a3094SHong Zhang       v += 25; jmin++;
9387fbae186SHong Zhang     }
939831a3094SHong Zhang     for (j=jmin; j<n; j++) {
94049b5e25fSSatish Balay       /* (strict lower triangular part of A)*x  */
94149b5e25fSSatish Balay       cval       = ib[j]*5;
94249b5e25fSSatish Balay       z[cval]     += v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4 + v[4]*x5;
94349b5e25fSSatish Balay       z[cval+1]   += v[5]*x1 + v[6]*x2 + v[7]*x3 + v[8]*x4 + v[9]*x5;
94449b5e25fSSatish Balay       z[cval+2]   += v[10]*x1 + v[11]*x2 + v[12]*x3 + v[13]*x4+ v[14]*x5;
94549b5e25fSSatish Balay       z[cval+3]   += v[15]*x1 + v[16]*x2 + v[17]*x3 + v[18]*x4+ v[19]*x5;
94649b5e25fSSatish Balay       z[cval+4]   += v[20]*x1 + v[21]*x2 + v[22]*x3 + v[23]*x4+ v[24]*x5;
94749b5e25fSSatish Balay       /* (strict upper triangular part of A)*x  */
94849b5e25fSSatish 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];
94949b5e25fSSatish 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];
95049b5e25fSSatish 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];
95149b5e25fSSatish 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];
95249b5e25fSSatish 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];
95349b5e25fSSatish Balay       v  += 25;
95449b5e25fSSatish Balay     }
95549b5e25fSSatish Balay     xb +=5; ai++;
95649b5e25fSSatish Balay   }
95749b5e25fSSatish Balay 
9581ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
9591ebc52fbSHong Zhang   if (yy != xx) ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
9601ebc52fbSHong Zhang   if (zz != yy) ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
96149b5e25fSSatish Balay 
9626c6c5352SBarry Smith   PetscLogFlops(50*(a->nz*2 - A->m));
96349b5e25fSSatish Balay   PetscFunctionReturn(0);
96449b5e25fSSatish Balay }
9654a2ae208SSatish Balay #undef __FUNCT__
9664a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqSBAIJ_6"
96749b5e25fSSatish Balay int MatMultAdd_SeqSBAIJ_6(Mat A,Vec xx,Vec yy,Vec zz)
96849b5e25fSSatish Balay {
96949b5e25fSSatish Balay   Mat_SeqSBAIJ    *a = (Mat_SeqSBAIJ*)A->data;
97087828ca2SBarry Smith   PetscScalar     *x,*y,*z,*xb,x1,x2,x3,x4,x5,x6;
97149b5e25fSSatish Balay   MatScalar       *v;
972831a3094SHong Zhang   int             mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,ierr,*ib,cval,j,jmin;
97349b5e25fSSatish Balay 
97449b5e25fSSatish Balay   PetscFunctionBegin;
9751ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
97649b5e25fSSatish Balay   if (yy != xx) {
9771ebc52fbSHong Zhang     ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
97849b5e25fSSatish Balay   } else {
97949b5e25fSSatish Balay     y = x;
98049b5e25fSSatish Balay   }
98149b5e25fSSatish Balay   if (zz != yy) {
982a9d4b620SHong Zhang     /* ierr = VecCopy(yy,zz);CHKERRQ(ierr); */
9831ebc52fbSHong Zhang     ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
984a9d4b620SHong Zhang     ierr = PetscMemcpy(z,y,yy->n*sizeof(PetscScalar));CHKERRQ(ierr);
98549b5e25fSSatish Balay   } else {
98649b5e25fSSatish Balay     z = y;
98749b5e25fSSatish Balay   }
98849b5e25fSSatish Balay 
98949b5e25fSSatish Balay   v     = a->a;
99049b5e25fSSatish Balay   xb = x;
99149b5e25fSSatish Balay 
99249b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
99349b5e25fSSatish Balay     n  = ai[1] - ai[0]; /* length of i_th block row of A */
99449b5e25fSSatish Balay     x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5=xb[4]; x6=xb[5];
99549b5e25fSSatish Balay     ib = aj + *ai;
996831a3094SHong Zhang     jmin = 0;
9977fbae186SHong Zhang     if (*ib == i){     /* (diag of A)*x */
99849b5e25fSSatish Balay       z[6*i]   += v[0]*x1  + v[6]*x2 + v[12]*x3 + v[18]*x4+ v[24]*x5 + v[30]*x6;
99949b5e25fSSatish Balay       z[6*i+1] += v[6]*x1  + v[7]*x2 + v[13]*x3 + v[19]*x4+ v[25]*x5 + v[31]*x6;
100049b5e25fSSatish Balay       z[6*i+2] += v[12]*x1 +v[13]*x2 + v[14]*x3 + v[20]*x4+ v[26]*x5 + v[32]*x6;
100149b5e25fSSatish Balay       z[6*i+3] += v[18]*x1 +v[19]*x2 + v[20]*x3 + v[21]*x4+ v[27]*x5 + v[33]*x6;
100249b5e25fSSatish Balay       z[6*i+4] += v[24]*x1 +v[25]*x2 + v[26]*x3 + v[27]*x4+ v[28]*x5 + v[34]*x6;
100349b5e25fSSatish Balay       z[6*i+5] += v[30]*x1 +v[31]*x2 + v[32]*x3 + v[33]*x4+ v[34]*x5 + v[35]*x6;
1004831a3094SHong Zhang       v += 36; jmin++;
10057fbae186SHong Zhang     }
1006831a3094SHong Zhang     for (j=jmin; j<n; j++) {
100749b5e25fSSatish Balay       /* (strict lower triangular part of A)*x  */
100849b5e25fSSatish Balay       cval       = ib[j]*6;
100949b5e25fSSatish Balay       z[cval]   += v[0]*x1  + v[1]*x2 + v[2]*x3 + v[3]*x4+ v[4]*x5 + v[5]*x6;
101049b5e25fSSatish Balay       z[cval+1] += v[6]*x1  + v[7]*x2 + v[8]*x3 + v[9]*x4+ v[10]*x5 + v[11]*x6;
101149b5e25fSSatish Balay       z[cval+2] += v[12]*x1  + v[13]*x2 + v[14]*x3 + v[15]*x4+ v[16]*x5 + v[17]*x6;
101249b5e25fSSatish Balay       z[cval+3] += v[18]*x1  + v[19]*x2 + v[20]*x3 + v[21]*x4+ v[22]*x5 + v[23]*x6;
101349b5e25fSSatish Balay       z[cval+4] += v[24]*x1  + v[25]*x2 + v[26]*x3 + v[27]*x4+ v[28]*x5 + v[29]*x6;
101449b5e25fSSatish Balay       z[cval+5] += v[30]*x1  + v[31]*x2 + v[32]*x3 + v[33]*x4+ v[34]*x5 + v[35]*x6;
101549b5e25fSSatish Balay       /* (strict upper triangular part of A)*x  */
101649b5e25fSSatish 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];
101749b5e25fSSatish 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];
101849b5e25fSSatish 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];
101949b5e25fSSatish 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];
102049b5e25fSSatish 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];
102149b5e25fSSatish 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];
102249b5e25fSSatish Balay       v  += 36;
102349b5e25fSSatish Balay     }
102449b5e25fSSatish Balay     xb +=6; ai++;
102549b5e25fSSatish Balay   }
102649b5e25fSSatish Balay 
10271ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
10281ebc52fbSHong Zhang   if (yy != xx) ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
10291ebc52fbSHong Zhang   if (zz != yy) ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
103049b5e25fSSatish Balay 
10316c6c5352SBarry Smith   PetscLogFlops(72*(a->nz*2 - A->m));
103249b5e25fSSatish Balay   PetscFunctionReturn(0);
103349b5e25fSSatish Balay }
103449b5e25fSSatish Balay 
10354a2ae208SSatish Balay #undef __FUNCT__
10364a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqSBAIJ_7"
103749b5e25fSSatish Balay int MatMultAdd_SeqSBAIJ_7(Mat A,Vec xx,Vec yy,Vec zz)
103849b5e25fSSatish Balay {
103949b5e25fSSatish Balay   Mat_SeqSBAIJ    *a = (Mat_SeqSBAIJ*)A->data;
104087828ca2SBarry Smith   PetscScalar     *x,*y,*z,*xb,x1,x2,x3,x4,x5,x6,x7;
104149b5e25fSSatish Balay   MatScalar       *v;
1042831a3094SHong Zhang   int             mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,ierr,*ib,cval,j,jmin;
104349b5e25fSSatish Balay 
104449b5e25fSSatish Balay   PetscFunctionBegin;
10451ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
104649b5e25fSSatish Balay   if (yy != xx) {
10471ebc52fbSHong Zhang     ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
104849b5e25fSSatish Balay   } else {
104949b5e25fSSatish Balay     y = x;
105049b5e25fSSatish Balay   }
105149b5e25fSSatish Balay   if (zz != yy) {
1052a9d4b620SHong Zhang     /* ierr = VecCopy(yy,zz);CHKERRQ(ierr); */
10531ebc52fbSHong Zhang     ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
1054a9d4b620SHong Zhang     ierr = PetscMemcpy(z,y,yy->n*sizeof(PetscScalar));CHKERRQ(ierr);
105549b5e25fSSatish Balay   } else {
105649b5e25fSSatish Balay     z = y;
105749b5e25fSSatish Balay   }
105849b5e25fSSatish Balay 
105949b5e25fSSatish Balay   v     = a->a;
106049b5e25fSSatish Balay   xb = x;
106149b5e25fSSatish Balay 
106249b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
106349b5e25fSSatish Balay     n  = ai[1] - ai[0]; /* length of i_th block row of A */
106449b5e25fSSatish Balay     x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5=xb[4]; x6=xb[5]; x7=xb[6];
106549b5e25fSSatish Balay     ib = aj + *ai;
1066831a3094SHong Zhang     jmin = 0;
10677fbae186SHong Zhang     if (*ib == i){     /* (diag of A)*x */
106849b5e25fSSatish 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;
106949b5e25fSSatish 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;
107049b5e25fSSatish 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;
107149b5e25fSSatish 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;
107249b5e25fSSatish 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;
107349b5e25fSSatish 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;
107449b5e25fSSatish 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;
1075831a3094SHong Zhang       v += 49; jmin++;
10767fbae186SHong Zhang     }
1077831a3094SHong Zhang     for (j=jmin; j<n; j++) {
107849b5e25fSSatish Balay       /* (strict lower triangular part of A)*x  */
107949b5e25fSSatish Balay       cval       = ib[j]*7;
108049b5e25fSSatish Balay       z[cval]   += v[0]*x1  + v[1]*x2 + v[2]*x3 + v[3]*x4+ v[4]*x5 + v[5]*x6+ v[6]*x7;
108149b5e25fSSatish 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;
108249b5e25fSSatish 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;
108349b5e25fSSatish 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;
108449b5e25fSSatish 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;
108549b5e25fSSatish 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;
108649b5e25fSSatish 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;
108749b5e25fSSatish Balay       /* (strict upper triangular part of A)*x  */
108849b5e25fSSatish 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];
108949b5e25fSSatish 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];
109049b5e25fSSatish 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];
109149b5e25fSSatish 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];
109249b5e25fSSatish 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];
109349b5e25fSSatish 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];
109449b5e25fSSatish 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];
109549b5e25fSSatish Balay       v  += 49;
109649b5e25fSSatish Balay     }
109749b5e25fSSatish Balay     xb +=7; ai++;
109849b5e25fSSatish Balay   }
109949b5e25fSSatish Balay 
11001ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
11011ebc52fbSHong Zhang   if (yy != xx) ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
11021ebc52fbSHong Zhang   if (zz != yy) ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
110349b5e25fSSatish Balay 
11046c6c5352SBarry Smith   PetscLogFlops(98*(a->nz*2 - A->m));
110549b5e25fSSatish Balay   PetscFunctionReturn(0);
110649b5e25fSSatish Balay }
110749b5e25fSSatish Balay 
11084a2ae208SSatish Balay #undef __FUNCT__
11094a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqSBAIJ_N"
111049b5e25fSSatish Balay int MatMultAdd_SeqSBAIJ_N(Mat A,Vec xx,Vec yy,Vec zz)
111149b5e25fSSatish Balay {
111249b5e25fSSatish Balay   Mat_SeqSBAIJ    *a = (Mat_SeqSBAIJ*)A->data;
111387828ca2SBarry Smith   PetscScalar     *x,*x_ptr,*y,*z,*z_ptr=0,*xb,*zb,*work,*workt;
1114066653e3SSatish Balay   MatScalar       *v;
1115066653e3SSatish Balay   int             ierr,mbs=a->mbs,i,*idx,*aj,*ii,bs=a->bs,j,n,bs2=a->bs2;
1116066653e3SSatish Balay   int             ncols,k;
111749b5e25fSSatish Balay 
111849b5e25fSSatish Balay   PetscFunctionBegin;
11191ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr); x_ptr=x;
112049b5e25fSSatish Balay   if (yy != xx) {
11211ebc52fbSHong Zhang     ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
112249b5e25fSSatish Balay   } else {
112349b5e25fSSatish Balay     y = x;
112449b5e25fSSatish Balay   }
112549b5e25fSSatish Balay   if (zz != yy) {
1126a9d4b620SHong Zhang     /* ierr = VecCopy(yy,zz);CHKERRQ(ierr); */
11271ebc52fbSHong Zhang     ierr = VecGetArray(zz,&z);CHKERRQ(ierr); z_ptr=z;
1128a9d4b620SHong Zhang     ierr = PetscMemcpy(z,y,yy->n*sizeof(PetscScalar));CHKERRQ(ierr);
112949b5e25fSSatish Balay   } else {
113049b5e25fSSatish Balay     z = y;
113149b5e25fSSatish Balay   }
113249b5e25fSSatish Balay 
113349b5e25fSSatish Balay   aj   = a->j;
113449b5e25fSSatish Balay   v    = a->a;
113549b5e25fSSatish Balay   ii   = a->i;
113649b5e25fSSatish Balay 
113749b5e25fSSatish Balay   if (!a->mult_work) {
113887828ca2SBarry Smith     ierr = PetscMalloc((A->m+1)*sizeof(PetscScalar),&a->mult_work);CHKERRQ(ierr);
113949b5e25fSSatish Balay   }
114049b5e25fSSatish Balay   work = a->mult_work;
114149b5e25fSSatish Balay 
114249b5e25fSSatish Balay 
114349b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
114449b5e25fSSatish Balay     n     = ii[1] - ii[0]; ncols = n*bs;
114549b5e25fSSatish Balay     workt = work; idx=aj+ii[0];
114649b5e25fSSatish Balay 
114749b5e25fSSatish Balay     /* upper triangular part */
114849b5e25fSSatish Balay     for (j=0; j<n; j++) {
114949b5e25fSSatish Balay       xb = x_ptr + bs*(*idx++);
115049b5e25fSSatish Balay       for (k=0; k<bs; k++) workt[k] = xb[k];
115149b5e25fSSatish Balay       workt += bs;
115249b5e25fSSatish Balay     }
115349b5e25fSSatish Balay     /* z(i*bs:(i+1)*bs-1) += A(i,:)*x */
115449b5e25fSSatish Balay     Kernel_w_gets_w_plus_Ar_times_v(bs,ncols,work,v,z);
115549b5e25fSSatish Balay 
115649b5e25fSSatish Balay     /* strict lower triangular part */
1157831a3094SHong Zhang     idx = aj+ii[0];
1158831a3094SHong Zhang     if (*idx == i){
115996b9376eSHong Zhang       ncols -= bs; v += bs2; idx++; n--;
1160831a3094SHong Zhang     }
116149b5e25fSSatish Balay     if (ncols > 0){
116249b5e25fSSatish Balay       workt = work;
116387828ca2SBarry Smith       ierr  = PetscMemzero(workt,ncols*sizeof(PetscScalar));CHKERRQ(ierr);
1164831a3094SHong Zhang       Kernel_w_gets_w_plus_trans_Ar_times_v(bs,ncols,x,v,workt);
1165831a3094SHong Zhang       for (j=0; j<n; j++) {
1166831a3094SHong Zhang         zb = z_ptr + bs*(*idx++);
1167831a3094SHong Zhang         /* idx++; */
116849b5e25fSSatish Balay         for (k=0; k<bs; k++) zb[k] += workt[k] ;
116949b5e25fSSatish Balay         workt += bs;
117049b5e25fSSatish Balay       }
117149b5e25fSSatish Balay     }
117249b5e25fSSatish Balay 
117349b5e25fSSatish Balay     x += bs; v += n*bs2; z += bs; ii++;
117449b5e25fSSatish Balay   }
117549b5e25fSSatish Balay 
11761ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
11771ebc52fbSHong Zhang   if (yy != xx) ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
11781ebc52fbSHong Zhang   if (zz != yy) ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
117949b5e25fSSatish Balay 
11806c6c5352SBarry Smith   PetscLogFlops(2*(a->nz*2 - A->m));
118149b5e25fSSatish Balay   PetscFunctionReturn(0);
118249b5e25fSSatish Balay }
118349b5e25fSSatish Balay 
11844a2ae208SSatish Balay #undef __FUNCT__
11854a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_SeqSBAIJ"
118649b5e25fSSatish Balay int MatMultTranspose_SeqSBAIJ(Mat A,Vec xx,Vec zz)
118749b5e25fSSatish Balay {
1188eeeff2ecSHong Zhang   int ierr;
1189eeeff2ecSHong Zhang 
119049b5e25fSSatish Balay   PetscFunctionBegin;
1191eeeff2ecSHong Zhang   ierr = MatMult(A,xx,zz);CHKERRQ(ierr);
1192eeeff2ecSHong Zhang   PetscFunctionReturn(0);
119349b5e25fSSatish Balay }
119449b5e25fSSatish Balay 
11954a2ae208SSatish Balay #undef __FUNCT__
11964a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_SeqSBAIJ"
119749b5e25fSSatish Balay int MatMultTransposeAdd_SeqSBAIJ(Mat A,Vec xx,Vec yy,Vec zz)
119849b5e25fSSatish Balay 
119949b5e25fSSatish Balay {
1200eeeff2ecSHong Zhang   int ierr;
1201eeeff2ecSHong Zhang 
120249b5e25fSSatish Balay   PetscFunctionBegin;
1203eeeff2ecSHong Zhang   ierr = MatMultAdd(A,xx,yy,zz);CHKERRQ(ierr);
1204eeeff2ecSHong Zhang   PetscFunctionReturn(0);
120549b5e25fSSatish Balay }
120649b5e25fSSatish Balay 
12074a2ae208SSatish Balay #undef __FUNCT__
12084a2ae208SSatish Balay #define __FUNCT__ "MatScale_SeqSBAIJ"
1209268466fbSBarry Smith int MatScale_SeqSBAIJ(const PetscScalar *alpha,Mat inA)
121049b5e25fSSatish Balay {
121149b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)inA->data;
12126c6c5352SBarry Smith   int         one = 1,totalnz = a->bs2*a->nz;
121349b5e25fSSatish Balay 
121449b5e25fSSatish Balay   PetscFunctionBegin;
1215268466fbSBarry Smith   BLscal_(&totalnz,(PetscScalar*)alpha,a->a,&one);
1216b0a32e0cSBarry Smith   PetscLogFlops(totalnz);
121749b5e25fSSatish Balay   PetscFunctionReturn(0);
121849b5e25fSSatish Balay }
121949b5e25fSSatish Balay 
12204a2ae208SSatish Balay #undef __FUNCT__
12214a2ae208SSatish Balay #define __FUNCT__ "MatNorm_SeqSBAIJ"
122249b5e25fSSatish Balay int MatNorm_SeqSBAIJ(Mat A,NormType type,PetscReal *norm)
122349b5e25fSSatish Balay {
122449b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
122549b5e25fSSatish Balay   MatScalar   *v = a->a;
122649b5e25fSSatish Balay   PetscReal   sum_diag = 0.0, sum_off = 0.0, *sum;
1227831a3094SHong Zhang   int         i,j,k,bs = a->bs,bs2=a->bs2,k1,mbs=a->mbs,*aj=a->j;
1228831a3094SHong Zhang   int         *jl,*il,jmin,jmax,ierr,nexti,ik,*col;
122949b5e25fSSatish Balay 
123049b5e25fSSatish Balay   PetscFunctionBegin;
123149b5e25fSSatish Balay   if (type == NORM_FROBENIUS) {
123249b5e25fSSatish Balay     for (k=0; k<mbs; k++){
123349b5e25fSSatish Balay       jmin = a->i[k]; jmax = a->i[k+1];
1234831a3094SHong Zhang       col  = aj + jmin;
1235831a3094SHong Zhang       if (*col == k){         /* diagonal block */
123649b5e25fSSatish Balay         for (i=0; i<bs2; i++){
123749b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX)
123849b5e25fSSatish Balay           sum_diag += PetscRealPart(PetscConj(*v)*(*v)); v++;
123949b5e25fSSatish Balay #else
124049b5e25fSSatish Balay           sum_diag += (*v)*(*v); v++;
124149b5e25fSSatish Balay #endif
124249b5e25fSSatish Balay         }
1243831a3094SHong Zhang         jmin++;
1244831a3094SHong Zhang       }
1245831a3094SHong Zhang       for (j=jmin; j<jmax; j++){  /* off-diagonal blocks */
124649b5e25fSSatish Balay         for (i=0; i<bs2; i++){
124749b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX)
124849b5e25fSSatish Balay           sum_off += PetscRealPart(PetscConj(*v)*(*v)); v++;
124949b5e25fSSatish Balay #else
125049b5e25fSSatish Balay           sum_off += (*v)*(*v); v++;
125149b5e25fSSatish Balay #endif
125249b5e25fSSatish Balay         }
125349b5e25fSSatish Balay       }
125449b5e25fSSatish Balay     }
125549b5e25fSSatish Balay     *norm = sqrt(sum_diag + 2*sum_off);
125649b5e25fSSatish Balay 
125749b5e25fSSatish Balay   }  else if (type == NORM_INFINITY) { /* maximum row sum */
125882502324SSatish Balay     ierr = PetscMalloc(mbs*sizeof(int),&il);CHKERRQ(ierr);
125982502324SSatish Balay     ierr = PetscMalloc(mbs*sizeof(int),&jl);CHKERRQ(ierr);
126082502324SSatish Balay     ierr = PetscMalloc(bs*sizeof(PetscReal),&sum);CHKERRQ(ierr);
126149b5e25fSSatish Balay     for (i=0; i<mbs; i++) {
126249b5e25fSSatish Balay       jl[i] = mbs; il[0] = 0;
126349b5e25fSSatish Balay     }
126449b5e25fSSatish Balay 
126549b5e25fSSatish Balay     *norm = 0.0;
126649b5e25fSSatish Balay     for (k=0; k<mbs; k++) { /* k_th block row */
126749b5e25fSSatish Balay       for (j=0; j<bs; j++) sum[j]=0.0;
126849b5e25fSSatish Balay 
126949b5e25fSSatish Balay       /*-- col sum --*/
127049b5e25fSSatish Balay       i = jl[k]; /* first |A(i,k)| to be added */
1271831a3094SHong Zhang       /* jl[k]=i: first nozero element in row i for submatrix A(1:k,k:n) (active window)
1272831a3094SHong Zhang                   at step k */
127349b5e25fSSatish Balay       while (i<mbs){
127449b5e25fSSatish Balay         nexti = jl[i];  /* next block row to be added */
127549b5e25fSSatish Balay         ik    = il[i];  /* block index of A(i,k) in the array a */
127649b5e25fSSatish Balay         for (j=0; j<bs; j++){
127749b5e25fSSatish Balay           v = a->a + ik*bs2 + j*bs;
127849b5e25fSSatish Balay           for (k1=0; k1<bs; k1++) {
127949b5e25fSSatish Balay             sum[j] += PetscAbsScalar(*v); v++;
128049b5e25fSSatish Balay           }
128149b5e25fSSatish Balay         }
128249b5e25fSSatish Balay         /* update il, jl */
1283831a3094SHong Zhang         jmin = ik + 1; /* block index of array a: points to the next nonzero of A in row i */
1284831a3094SHong Zhang         jmax = a->i[i+1];
128549b5e25fSSatish Balay         if (jmin < jmax){
128649b5e25fSSatish Balay           il[i] = jmin;
128749b5e25fSSatish Balay           j   = a->j[jmin];
128849b5e25fSSatish Balay           jl[i] = jl[j]; jl[j]=i;
128949b5e25fSSatish Balay         }
129049b5e25fSSatish Balay         i = nexti;
129149b5e25fSSatish Balay       }
129249b5e25fSSatish Balay 
129349b5e25fSSatish Balay       /*-- row sum --*/
129449b5e25fSSatish Balay       jmin = a->i[k]; jmax = a->i[k+1];
129549b5e25fSSatish Balay       for (i=jmin; i<jmax; i++) {
129649b5e25fSSatish Balay         for (j=0; j<bs; j++){
129749b5e25fSSatish Balay           v = a->a + i*bs2 + j;
129849b5e25fSSatish Balay           for (k1=0; k1<bs; k1++){
129949b5e25fSSatish Balay             sum[j] += PetscAbsScalar(*v);
130049b5e25fSSatish Balay             v   += bs;
130149b5e25fSSatish Balay           }
130249b5e25fSSatish Balay         }
130349b5e25fSSatish Balay       }
130449b5e25fSSatish Balay       /* add k_th block row to il, jl */
1305831a3094SHong Zhang       col = aj+jmin;
1306831a3094SHong Zhang       if (*col == k) jmin++;
130749b5e25fSSatish Balay       if (jmin < jmax){
130849b5e25fSSatish Balay         il[k] = jmin;
130949b5e25fSSatish Balay         j   = a->j[jmin];
131049b5e25fSSatish Balay         jl[k] = jl[j]; jl[j] = k;
131149b5e25fSSatish Balay       }
131249b5e25fSSatish Balay       for (j=0; j<bs; j++){
131349b5e25fSSatish Balay         if (sum[j] > *norm) *norm = sum[j];
131449b5e25fSSatish Balay       }
131549b5e25fSSatish Balay     }
131649b5e25fSSatish Balay     ierr = PetscFree(il);CHKERRQ(ierr);
131749b5e25fSSatish Balay     ierr = PetscFree(jl);CHKERRQ(ierr);
131849b5e25fSSatish Balay     ierr = PetscFree(sum);CHKERRQ(ierr);
131949b5e25fSSatish Balay   } else {
1320347d480fSBarry Smith     SETERRQ(PETSC_ERR_SUP,"No support for this norm yet");
132149b5e25fSSatish Balay   }
132249b5e25fSSatish Balay   PetscFunctionReturn(0);
132349b5e25fSSatish Balay }
132449b5e25fSSatish Balay 
13254a2ae208SSatish Balay #undef __FUNCT__
13264a2ae208SSatish Balay #define __FUNCT__ "MatEqual_SeqSBAIJ"
132749b5e25fSSatish Balay int MatEqual_SeqSBAIJ(Mat A,Mat B,PetscTruth* flg)
132849b5e25fSSatish Balay {
132949b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ *)A->data,*b = (Mat_SeqSBAIJ *)B->data;
133049b5e25fSSatish Balay   int          ierr;
133149b5e25fSSatish Balay 
133249b5e25fSSatish Balay   PetscFunctionBegin;
133349b5e25fSSatish Balay 
133449b5e25fSSatish Balay   /* If the  matrix/block dimensions are not equal, or no of nonzeros or shift */
13356c6c5352SBarry Smith   if ((A->m != B->m) || (A->n != B->n) || (a->bs != b->bs)|| (a->nz != b->nz)) {
1336ef511fbeSHong Zhang     *flg = PETSC_FALSE;
1337ef511fbeSHong Zhang     PetscFunctionReturn(0);
133849b5e25fSSatish Balay   }
133949b5e25fSSatish Balay 
134049b5e25fSSatish Balay   /* if the a->i are the same */
134149b5e25fSSatish Balay   ierr = PetscMemcmp(a->i,b->i,(a->mbs+1)*sizeof(int),flg);CHKERRQ(ierr);
134249b5e25fSSatish Balay   if (*flg == PETSC_FALSE) {
134349b5e25fSSatish Balay     PetscFunctionReturn(0);
134449b5e25fSSatish Balay   }
134549b5e25fSSatish Balay 
134649b5e25fSSatish Balay   /* if a->j are the same */
13476c6c5352SBarry Smith   ierr = PetscMemcmp(a->j,b->j,(a->nz)*sizeof(int),flg);CHKERRQ(ierr);
134849b5e25fSSatish Balay   if (*flg == PETSC_FALSE) {
134949b5e25fSSatish Balay     PetscFunctionReturn(0);
135049b5e25fSSatish Balay   }
135149b5e25fSSatish Balay   /* if a->a are the same */
13526c6c5352SBarry Smith   ierr = PetscMemcmp(a->a,b->a,(a->nz)*(a->bs)*(a->bs)*sizeof(PetscScalar),flg);CHKERRQ(ierr);
135349b5e25fSSatish Balay 
1354935af2e7SHong Zhang   PetscFunctionReturn(0);
135549b5e25fSSatish Balay }
135649b5e25fSSatish Balay 
13574a2ae208SSatish Balay #undef __FUNCT__
13584a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_SeqSBAIJ"
135949b5e25fSSatish Balay int MatGetDiagonal_SeqSBAIJ(Mat A,Vec v)
136049b5e25fSSatish Balay {
136149b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
136249b5e25fSSatish Balay   int          ierr,i,j,k,n,row,bs,*ai,*aj,ambs,bs2;
136387828ca2SBarry Smith   PetscScalar  *x,zero = 0.0;
136449b5e25fSSatish Balay   MatScalar    *aa,*aa_j;
136549b5e25fSSatish Balay 
136649b5e25fSSatish Balay   PetscFunctionBegin;
136749b5e25fSSatish Balay   bs   = a->bs;
136882799104SHong Zhang   if (A->factor && bs>1) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix with bs>1");
136982799104SHong Zhang 
137049b5e25fSSatish Balay   aa   = a->a;
137149b5e25fSSatish Balay   ai   = a->i;
137249b5e25fSSatish Balay   aj   = a->j;
137349b5e25fSSatish Balay   ambs = a->mbs;
137449b5e25fSSatish Balay   bs2  = a->bs2;
137549b5e25fSSatish Balay 
137649b5e25fSSatish Balay   ierr = VecSet(&zero,v);CHKERRQ(ierr);
13771ebc52fbSHong Zhang   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
137849b5e25fSSatish Balay   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
1379ef511fbeSHong Zhang   if (n != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
138049b5e25fSSatish Balay   for (i=0; i<ambs; i++) {
138149b5e25fSSatish Balay     j=ai[i];
138249b5e25fSSatish Balay     if (aj[j] == i) {             /* if this is a diagonal element */
138349b5e25fSSatish Balay       row  = i*bs;
138449b5e25fSSatish Balay       aa_j = aa + j*bs2;
138582799104SHong Zhang       if (A->factor && bs==1){
138682799104SHong Zhang         for (k=0; k<bs2; k+=(bs+1),row++) x[row] = 1.0/aa_j[k];
138782799104SHong Zhang       } else {
138849b5e25fSSatish Balay         for (k=0; k<bs2; k+=(bs+1),row++) x[row] = aa_j[k];
138949b5e25fSSatish Balay       }
139049b5e25fSSatish Balay     }
139182799104SHong Zhang   }
139282799104SHong Zhang 
13931ebc52fbSHong Zhang   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
139449b5e25fSSatish Balay   PetscFunctionReturn(0);
139549b5e25fSSatish Balay }
139649b5e25fSSatish Balay 
13974a2ae208SSatish Balay #undef __FUNCT__
13984a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_SeqSBAIJ"
139949b5e25fSSatish Balay int MatDiagonalScale_SeqSBAIJ(Mat A,Vec ll,Vec rr)
140049b5e25fSSatish Balay {
140149b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
140287828ca2SBarry Smith   PetscScalar  *l,*r,x,*li,*ri;
140349b5e25fSSatish Balay   MatScalar    *aa,*v;
1404066653e3SSatish Balay   int          ierr,i,j,k,lm,rn,M,m,*ai,*aj,mbs,tmp,bs,bs2;
140549b5e25fSSatish Balay 
140649b5e25fSSatish Balay   PetscFunctionBegin;
140749b5e25fSSatish Balay   ai  = a->i;
140849b5e25fSSatish Balay   aj  = a->j;
140949b5e25fSSatish Balay   aa  = a->a;
1410ef511fbeSHong Zhang   m   = A->m;
141149b5e25fSSatish Balay   bs  = a->bs;
141249b5e25fSSatish Balay   mbs = a->mbs;
141349b5e25fSSatish Balay   bs2 = a->bs2;
141449b5e25fSSatish Balay 
141549b5e25fSSatish Balay   if (ll != rr) {
1416347d480fSBarry Smith     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"For symmetric format, left and right scaling vectors must be same\n");
141749b5e25fSSatish Balay   }
141849b5e25fSSatish Balay   if (ll) {
14191ebc52fbSHong Zhang     ierr = VecGetArray(ll,&l);CHKERRQ(ierr);
142049b5e25fSSatish Balay     ierr = VecGetLocalSize(ll,&lm);CHKERRQ(ierr);
1421347d480fSBarry Smith     if (lm != m) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vector wrong length");
142249b5e25fSSatish Balay     for (i=0; i<mbs; i++) { /* for each block row */
142349b5e25fSSatish Balay       M  = ai[i+1] - ai[i];
142449b5e25fSSatish Balay       li = l + i*bs;
142549b5e25fSSatish Balay       v  = aa + bs2*ai[i];
142649b5e25fSSatish Balay       for (j=0; j<M; j++) { /* for each block */
142749b5e25fSSatish Balay         for (k=0; k<bs2; k++) {
142849b5e25fSSatish Balay           (*v++) *= li[k%bs];
142949b5e25fSSatish Balay         }
143049b5e25fSSatish Balay #ifdef CONT
143149b5e25fSSatish Balay         /* will be used to replace the above loop */
143249b5e25fSSatish Balay         ri = l + bs*aj[ai[i]+j];
143349b5e25fSSatish Balay         for (k=0; k<bs; k++) { /* column value */
143449b5e25fSSatish Balay           x = ri[k];
143549b5e25fSSatish Balay           for (tmp=0; tmp<bs; tmp++) (*v++) *= li[tmp]*x;
143649b5e25fSSatish Balay         }
143749b5e25fSSatish Balay #endif
143849b5e25fSSatish Balay 
143949b5e25fSSatish Balay       }
144049b5e25fSSatish Balay     }
14411ebc52fbSHong Zhang     ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr);
14426c6c5352SBarry Smith     PetscLogFlops(2*a->nz);
144349b5e25fSSatish Balay   }
144449b5e25fSSatish Balay   /* will be deleted */
144549b5e25fSSatish Balay   if (rr) {
14461ebc52fbSHong Zhang     ierr = VecGetArray(rr,&r);CHKERRQ(ierr);
144749b5e25fSSatish Balay     ierr = VecGetLocalSize(rr,&rn);CHKERRQ(ierr);
1448347d480fSBarry Smith     if (rn != m) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vector wrong length");
144949b5e25fSSatish Balay     for (i=0; i<mbs; i++) { /* for each block row */
145049b5e25fSSatish Balay       M  = ai[i+1] - ai[i];
145149b5e25fSSatish Balay       v  = aa + bs2*ai[i];
145249b5e25fSSatish Balay       for (j=0; j<M; j++) { /* for each block */
145349b5e25fSSatish Balay         ri = r + bs*aj[ai[i]+j];
145449b5e25fSSatish Balay         for (k=0; k<bs; k++) {
145549b5e25fSSatish Balay           x = ri[k];
145649b5e25fSSatish Balay           for (tmp=0; tmp<bs; tmp++) (*v++) *= x;
145749b5e25fSSatish Balay         }
145849b5e25fSSatish Balay       }
145949b5e25fSSatish Balay     }
14601ebc52fbSHong Zhang     ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr);
14616c6c5352SBarry Smith     PetscLogFlops(a->nz);
146249b5e25fSSatish Balay   }
146349b5e25fSSatish Balay   PetscFunctionReturn(0);
146449b5e25fSSatish Balay }
146549b5e25fSSatish Balay 
14664a2ae208SSatish Balay #undef __FUNCT__
14674a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_SeqSBAIJ"
146849b5e25fSSatish Balay int MatGetInfo_SeqSBAIJ(Mat A,MatInfoType flag,MatInfo *info)
146949b5e25fSSatish Balay {
147049b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
147149b5e25fSSatish Balay 
147249b5e25fSSatish Balay   PetscFunctionBegin;
1473ef511fbeSHong Zhang   info->rows_global    = (double)A->m;
1474ef511fbeSHong Zhang   info->columns_global = (double)A->m;
1475ef511fbeSHong Zhang   info->rows_local     = (double)A->m;
1476ef511fbeSHong Zhang   info->columns_local  = (double)A->m;
147749b5e25fSSatish Balay   info->block_size     = a->bs2;
14786c6c5352SBarry Smith   info->nz_allocated   = a->maxnz; /*num. of nonzeros in upper triangular part */
14796c6c5352SBarry Smith   info->nz_used        = a->bs2*a->nz; /*num. of nonzeros in upper triangular part */
148049b5e25fSSatish Balay   info->nz_unneeded    = (double)(info->nz_allocated - info->nz_used);
148149b5e25fSSatish Balay   info->assemblies   = A->num_ass;
148249b5e25fSSatish Balay   info->mallocs      = a->reallocs;
148349b5e25fSSatish Balay   info->memory       = A->mem;
148449b5e25fSSatish Balay   if (A->factor) {
148549b5e25fSSatish Balay     info->fill_ratio_given  = A->info.fill_ratio_given;
148649b5e25fSSatish Balay     info->fill_ratio_needed = A->info.fill_ratio_needed;
148749b5e25fSSatish Balay     info->factor_mallocs    = A->info.factor_mallocs;
148849b5e25fSSatish Balay   } else {
148949b5e25fSSatish Balay     info->fill_ratio_given  = 0;
149049b5e25fSSatish Balay     info->fill_ratio_needed = 0;
149149b5e25fSSatish Balay     info->factor_mallocs    = 0;
149249b5e25fSSatish Balay   }
149349b5e25fSSatish Balay   PetscFunctionReturn(0);
149449b5e25fSSatish Balay }
149549b5e25fSSatish Balay 
149649b5e25fSSatish Balay 
14974a2ae208SSatish Balay #undef __FUNCT__
14984a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_SeqSBAIJ"
149949b5e25fSSatish Balay int MatZeroEntries_SeqSBAIJ(Mat A)
150049b5e25fSSatish Balay {
150149b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
150249b5e25fSSatish Balay   int         ierr;
150349b5e25fSSatish Balay 
150449b5e25fSSatish Balay   PetscFunctionBegin;
150549b5e25fSSatish Balay   ierr = PetscMemzero(a->a,a->bs2*a->i[a->mbs]*sizeof(MatScalar));CHKERRQ(ierr);
150649b5e25fSSatish Balay   PetscFunctionReturn(0);
150749b5e25fSSatish Balay }
1508dc354874SHong Zhang 
15094a2ae208SSatish Balay #undef __FUNCT__
15104a2ae208SSatish Balay #define __FUNCT__ "MatGetRowMax_SeqSBAIJ"
1511dc354874SHong Zhang int MatGetRowMax_SeqSBAIJ(Mat A,Vec v)
1512dc354874SHong Zhang {
1513dc354874SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
1514d0f6400bSHong Zhang   int          ierr,i,j,n,row,col,bs,*ai,*aj,mbs;
1515c3fca9a7SHong Zhang   PetscReal    atmp;
1516273d9f13SBarry Smith   MatScalar    *aa;
151787828ca2SBarry Smith   PetscScalar  zero = 0.0,*x;
1518d0f6400bSHong Zhang   int          ncols,brow,bcol,krow,kcol;
1519dc354874SHong Zhang 
1520dc354874SHong Zhang   PetscFunctionBegin;
1521dc354874SHong Zhang   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1522dc354874SHong Zhang   bs   = a->bs;
1523dc354874SHong Zhang   aa   = a->a;
1524dc354874SHong Zhang   ai   = a->i;
1525dc354874SHong Zhang   aj   = a->j;
152644117c81SHong Zhang   mbs = a->mbs;
1527dc354874SHong Zhang 
1528dc354874SHong Zhang   ierr = VecSet(&zero,v);CHKERRQ(ierr);
15291ebc52fbSHong Zhang   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1530dc354874SHong Zhang   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
1531ef511fbeSHong Zhang   if (n != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
153244117c81SHong Zhang   for (i=0; i<mbs; i++) {
1533d0f6400bSHong Zhang     ncols = ai[1] - ai[0]; ai++;
1534d0f6400bSHong Zhang     brow  = bs*i;
153544117c81SHong Zhang     for (j=0; j<ncols; j++){
1536d0f6400bSHong Zhang       bcol = bs*(*aj);
153744117c81SHong Zhang       for (kcol=0; kcol<bs; kcol++){
1538d0f6400bSHong Zhang         col = bcol + kcol;      /* col index */
153944117c81SHong Zhang         for (krow=0; krow<bs; krow++){
1540d0f6400bSHong Zhang           atmp = PetscAbsScalar(*aa); aa++;
1541d0f6400bSHong Zhang           row = brow + krow;    /* row index */
154244117c81SHong Zhang           /* printf("val[%d,%d]: %g\n",row,col,atmp); */
1543c3fca9a7SHong Zhang           if (PetscRealPart(x[row]) < atmp) x[row] = atmp;
1544c3fca9a7SHong Zhang           if (*aj > i && PetscRealPart(x[col]) < atmp) x[col] = atmp;
154544117c81SHong Zhang         }
154644117c81SHong Zhang       }
1547d0f6400bSHong Zhang       aj++;
1548dc354874SHong Zhang     }
1549dc354874SHong Zhang   }
15501ebc52fbSHong Zhang   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1551dc354874SHong Zhang   PetscFunctionReturn(0);
1552dc354874SHong Zhang }
1553