xref: /petsc/src/mat/impls/sbaij/seq/sbaij2.c (revision f69a0ea3504314d029ee423e0de2950ece2dab86)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
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"
1113f74950SBarry Smith PetscErrorCode MatIncreaseOverlap_SeqSBAIJ(Mat A,PetscInt is_max,IS is[],PetscInt ov)
1249b5e25fSSatish Balay {
135eee224dSHong Zhang   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
146849ba73SBarry Smith   PetscErrorCode ierr;
15521d7252SBarry Smith   PetscInt       brow,i,j,k,l,mbs,n,*idx,*nidx,isz,bcol,bcol_max,start,end,*ai,*aj,bs,*nidx2;
16521d7252SBarry Smith   PetscBT        table,table0;
17d94109b8SHong Zhang 
18d94109b8SHong Zhang   PetscFunctionBegin;
19b3bf805bSHong Zhang   if (ov < 0)  SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap specified");
20d94109b8SHong Zhang   mbs = a->mbs;
21d94109b8SHong Zhang   ai  = a->i;
22d94109b8SHong Zhang   aj  = a->j;
23521d7252SBarry Smith   bs  = A->bs;
24d94109b8SHong Zhang   ierr = PetscBTCreate(mbs,table);CHKERRQ(ierr);
2513f74950SBarry Smith   ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&nidx);CHKERRQ(ierr);
2613f74950SBarry Smith   ierr = PetscMalloc((A->m+1)*sizeof(PetscInt),&nidx2);CHKERRQ(ierr);
27d94109b8SHong Zhang   ierr = PetscBTCreate(mbs,table0);CHKERRQ(ierr);
28d94109b8SHong Zhang 
29d94109b8SHong Zhang   for (i=0; i<is_max; i++) { /* for each is */
30d94109b8SHong Zhang     isz  = 0;
31d94109b8SHong Zhang     ierr = PetscBTMemzero(mbs,table);CHKERRQ(ierr);
32d94109b8SHong Zhang 
33d94109b8SHong Zhang     /* Extract the indices, assume there can be duplicate entries */
34d94109b8SHong Zhang     ierr = ISGetIndices(is[i],&idx);CHKERRQ(ierr);
35d94109b8SHong Zhang     ierr = ISGetLocalSize(is[i],&n);CHKERRQ(ierr);
36d94109b8SHong Zhang 
37d94109b8SHong Zhang     /* Enter these into the temp arrays i.e mark table[brow], enter brow into new index */
38dbe03f88SHong Zhang     bcol_max = 0;
39d94109b8SHong Zhang     for (j=0; j<n ; ++j){
40d94109b8SHong Zhang       brow = idx[j]/bs; /* convert the indices into block indices */
41d94109b8SHong Zhang       if (brow >= mbs) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"index greater than mat-dim");
42dbe03f88SHong Zhang       if(!PetscBTLookupSet(table,brow)) {
43dbe03f88SHong Zhang         nidx[isz++] = brow;
44dbe03f88SHong Zhang         if (bcol_max < brow) bcol_max = brow;
45dbe03f88SHong Zhang       }
46d94109b8SHong Zhang     }
47d94109b8SHong Zhang     ierr = ISRestoreIndices(is[i],&idx);CHKERRQ(ierr);
48d94109b8SHong Zhang     ierr = ISDestroy(is[i]);CHKERRQ(ierr);
49d94109b8SHong Zhang 
50d94109b8SHong Zhang     k = 0;
51d94109b8SHong Zhang     for (j=0; j<ov; j++){ /* for each overlap */
52d94109b8SHong Zhang       /* set table0 for lookup - only mark entries that are added onto nidx in (j-1)-th overlap */
53d94109b8SHong Zhang       ierr = PetscBTMemzero(mbs,table0);CHKERRQ(ierr);
54efee365bSSatish Balay       for (l=k; l<isz; l++) { ierr = PetscBTSet(table0,nidx[l]);CHKERRQ(ierr); }
55d94109b8SHong Zhang 
56d94109b8SHong Zhang       n = isz;  /* length of the updated is[i] */
57d94109b8SHong Zhang       for (brow=0; brow<mbs; brow++){
58d94109b8SHong Zhang         start = ai[brow]; end   = ai[brow+1];
59d94109b8SHong Zhang         if (PetscBTLookup(table0,brow)){ /* brow is on nidx - row search: collect all bcol in this brow */
60d94109b8SHong Zhang           for (l = start; l<end ; l++){
61d94109b8SHong Zhang             bcol = aj[l];
62d94109b8SHong Zhang             if (!PetscBTLookupSet(table,bcol)) {nidx[isz++] = bcol;}
63d94109b8SHong Zhang           }
64d94109b8SHong Zhang           k++;
65d94109b8SHong Zhang           if (k >= n) break; /* for (brow=0; brow<mbs; brow++) */
66d94109b8SHong Zhang         } else { /* brow is not on nidx - col serach: add brow onto nidx if there is a bcol in nidx */
67d94109b8SHong Zhang           for (l = start; l<end ; l++){
68d94109b8SHong Zhang             bcol = aj[l];
69dbe03f88SHong Zhang             if (bcol > bcol_max) break;
70d94109b8SHong Zhang             if (PetscBTLookup(table0,bcol)){
71d94109b8SHong Zhang               if (!PetscBTLookupSet(table,brow)) {nidx[isz++] = brow;}
72d94109b8SHong Zhang               break; /* for l = start; l<end ; l++) */
73d94109b8SHong Zhang             }
74d94109b8SHong Zhang           }
75d94109b8SHong Zhang         }
76d94109b8SHong Zhang       }
77d94109b8SHong Zhang     } /* for each overlap */
78d94109b8SHong Zhang 
79d94109b8SHong Zhang     /* expand the Index Set */
80d94109b8SHong Zhang     for (j=0; j<isz; j++) {
81d94109b8SHong Zhang       for (k=0; k<bs; k++)
82d94109b8SHong Zhang         nidx2[j*bs+k] = nidx[j]*bs+k;
83d94109b8SHong Zhang     }
84d94109b8SHong Zhang     ierr = ISCreateGeneral(PETSC_COMM_SELF,isz*bs,nidx2,is+i);CHKERRQ(ierr);
85d94109b8SHong Zhang   }
86d94109b8SHong Zhang   ierr = PetscBTDestroy(table);CHKERRQ(ierr);
87d94109b8SHong Zhang   ierr = PetscFree(nidx);CHKERRQ(ierr);
88d94109b8SHong Zhang   ierr = PetscFree(nidx2);CHKERRQ(ierr);
89d94109b8SHong Zhang   ierr = PetscBTDestroy(table0);CHKERRQ(ierr);
905eee224dSHong Zhang   PetscFunctionReturn(0);
9149b5e25fSSatish Balay }
9249b5e25fSSatish Balay 
934a2ae208SSatish Balay #undef __FUNCT__
944a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqSBAIJ_Private"
9513f74950SBarry Smith PetscErrorCode MatGetSubMatrix_SeqSBAIJ_Private(Mat A,IS isrow,IS iscol,PetscInt cs,MatReuse scall,Mat *B)
9649b5e25fSSatish Balay {
9749b5e25fSSatish Balay   Mat_SeqSBAIJ  *a = (Mat_SeqSBAIJ*)A->data,*c;
986849ba73SBarry Smith   PetscErrorCode ierr;
9913f74950SBarry Smith   PetscInt       *smap,i,k,kstart,kend,oldcols = a->mbs,*lens;
10013f74950SBarry Smith   PetscInt       row,mat_i,*mat_j,tcol,*mat_ilen;
10113f74950SBarry Smith   PetscInt       *irow,nrows,*ssmap,bs=A->bs,bs2=a->bs2;
10213f74950SBarry Smith   PetscInt       *aj = a->j,*ai = a->i;
10349b5e25fSSatish Balay   MatScalar      *mat_a;
10449b5e25fSSatish Balay   Mat            C;
10549b5e25fSSatish Balay   PetscTruth     flag;
10649b5e25fSSatish Balay 
10749b5e25fSSatish Balay   PetscFunctionBegin;
108e005ede5SBarry Smith   if (isrow != iscol) SETERRQ(PETSC_ERR_ARG_INCOMP,"For symmetric format, iscol must equal isro");
10949b5e25fSSatish Balay   ierr = ISSorted(iscol,(PetscTruth*)&i);CHKERRQ(ierr);
110347d480fSBarry Smith   if (!i) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"IS is not sorted");
11149b5e25fSSatish Balay 
11249b5e25fSSatish Balay   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
11349b5e25fSSatish Balay   ierr = ISGetSize(isrow,&nrows);CHKERRQ(ierr);
11449b5e25fSSatish Balay 
11513f74950SBarry Smith   ierr  = PetscMalloc((1+oldcols)*sizeof(PetscInt),&smap);CHKERRQ(ierr);
11649b5e25fSSatish Balay   ssmap = smap;
11713f74950SBarry Smith   ierr  = PetscMalloc((1+nrows)*sizeof(PetscInt),&lens);CHKERRQ(ierr);
11813f74950SBarry Smith   ierr  = PetscMemzero(smap,oldcols*sizeof(PetscInt));CHKERRQ(ierr);
11949b5e25fSSatish Balay   for (i=0; i<nrows; i++) smap[irow[i]] = i+1; /* nrows = ncols */
12049b5e25fSSatish Balay   /* determine lens of each row */
12149b5e25fSSatish Balay   for (i=0; i<nrows; i++) {
12249b5e25fSSatish Balay     kstart  = ai[irow[i]];
12349b5e25fSSatish Balay     kend    = kstart + a->ilen[irow[i]];
12449b5e25fSSatish Balay     lens[i] = 0;
12549b5e25fSSatish Balay       for (k=kstart; k<kend; k++) {
12649b5e25fSSatish Balay         if (ssmap[aj[k]]) {
12749b5e25fSSatish Balay           lens[i]++;
12849b5e25fSSatish Balay         }
12949b5e25fSSatish Balay       }
13049b5e25fSSatish Balay     }
13149b5e25fSSatish Balay   /* Create and fill new matrix */
13249b5e25fSSatish Balay   if (scall == MAT_REUSE_MATRIX) {
13349b5e25fSSatish Balay     c = (Mat_SeqSBAIJ *)((*B)->data);
13449b5e25fSSatish Balay 
135521d7252SBarry Smith     if (c->mbs!=nrows || (*B)->bs!=bs) SETERRQ(PETSC_ERR_ARG_SIZ,"Submatrix wrong size");
13613f74950SBarry Smith     ierr = PetscMemcmp(c->ilen,lens,c->mbs *sizeof(PetscInt),&flag);CHKERRQ(ierr);
137abc0a331SBarry Smith     if (!flag) {
138347d480fSBarry Smith       SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong no of nonzeros");
13949b5e25fSSatish Balay     }
14013f74950SBarry Smith     ierr = PetscMemzero(c->ilen,c->mbs*sizeof(PetscInt));CHKERRQ(ierr);
14149b5e25fSSatish Balay     C = *B;
14249b5e25fSSatish Balay   } else {
143*f69a0ea3SMatthew Knepley     ierr = MatCreate(A->comm,&C);CHKERRQ(ierr);
144*f69a0ea3SMatthew Knepley     ierr = MatSetSizes(C,nrows*bs,nrows*bs,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
145e2d9671bSKris Buschelman     ierr = MatSetType(C,A->type_name);CHKERRQ(ierr);
146ab93d7beSBarry Smith     ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(C,bs,0,lens);CHKERRQ(ierr);
14749b5e25fSSatish Balay   }
14849b5e25fSSatish Balay   c = (Mat_SeqSBAIJ *)(C->data);
14949b5e25fSSatish Balay   for (i=0; i<nrows; i++) {
15049b5e25fSSatish Balay     row    = irow[i];
15149b5e25fSSatish Balay     kstart = ai[row];
15249b5e25fSSatish Balay     kend   = kstart + a->ilen[row];
15349b5e25fSSatish Balay     mat_i  = c->i[i];
15449b5e25fSSatish Balay     mat_j  = c->j + mat_i;
15549b5e25fSSatish Balay     mat_a  = c->a + mat_i*bs2;
15649b5e25fSSatish Balay     mat_ilen = c->ilen + i;
15749b5e25fSSatish Balay     for (k=kstart; k<kend; k++) {
15849b5e25fSSatish Balay       if ((tcol=ssmap[a->j[k]])) {
15949b5e25fSSatish Balay         *mat_j++ = tcol - 1;
16049b5e25fSSatish Balay         ierr     = PetscMemcpy(mat_a,a->a+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
16149b5e25fSSatish Balay         mat_a   += bs2;
16249b5e25fSSatish Balay         (*mat_ilen)++;
16349b5e25fSSatish Balay       }
16449b5e25fSSatish Balay     }
16549b5e25fSSatish Balay   }
16649b5e25fSSatish Balay 
16749b5e25fSSatish Balay   /* Free work space */
16849b5e25fSSatish Balay   ierr = PetscFree(smap);CHKERRQ(ierr);
16949b5e25fSSatish Balay   ierr = PetscFree(lens);CHKERRQ(ierr);
17049b5e25fSSatish Balay   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
17149b5e25fSSatish Balay   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
17249b5e25fSSatish Balay 
17349b5e25fSSatish Balay   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
17449b5e25fSSatish Balay   *B = C;
17549b5e25fSSatish Balay   PetscFunctionReturn(0);
17649b5e25fSSatish Balay }
17749b5e25fSSatish Balay 
1784a2ae208SSatish Balay #undef __FUNCT__
1794a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqSBAIJ"
18013f74950SBarry Smith PetscErrorCode MatGetSubMatrix_SeqSBAIJ(Mat A,IS isrow,IS iscol,PetscInt cs,MatReuse scall,Mat *B)
18149b5e25fSSatish Balay {
18249b5e25fSSatish Balay   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
18349b5e25fSSatish Balay   IS             is1;
1846849ba73SBarry Smith   PetscErrorCode ierr;
18513f74950SBarry Smith   PetscInt       *vary,*iary,*irow,nrows,i,bs=A->bs,count;
18649b5e25fSSatish Balay 
18749b5e25fSSatish Balay   PetscFunctionBegin;
188e005ede5SBarry Smith   if (isrow != iscol) SETERRQ(PETSC_ERR_ARG_INCOMP,"For symmetric format, iscol must equal isro");
18949b5e25fSSatish Balay 
19049b5e25fSSatish Balay   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
19149b5e25fSSatish Balay   ierr = ISGetSize(isrow,&nrows);CHKERRQ(ierr);
19249b5e25fSSatish Balay 
19349b5e25fSSatish Balay   /* Verify if the indices corespond to each element in a block
19449b5e25fSSatish Balay    and form the IS with compressed IS */
19513f74950SBarry Smith   ierr = PetscMalloc(2*(a->mbs+1)*sizeof(PetscInt),&vary);CHKERRQ(ierr);
19649b5e25fSSatish Balay   iary = vary + a->mbs;
19713f74950SBarry Smith   ierr = PetscMemzero(vary,(a->mbs)*sizeof(PetscInt));CHKERRQ(ierr);
19849b5e25fSSatish Balay   for (i=0; i<nrows; i++) vary[irow[i]/bs]++;
19949b5e25fSSatish Balay 
20049b5e25fSSatish Balay   count = 0;
20149b5e25fSSatish Balay   for (i=0; i<a->mbs; i++) {
202e005ede5SBarry Smith     if (vary[i]!=0 && vary[i]!=bs) SETERRQ(PETSC_ERR_ARG_INCOMP,"Index set does not match blocks");
20349b5e25fSSatish Balay     if (vary[i]==bs) iary[count++] = i;
20449b5e25fSSatish Balay   }
20549b5e25fSSatish Balay   ierr = ISCreateGeneral(PETSC_COMM_SELF,count,iary,&is1);CHKERRQ(ierr);
20649b5e25fSSatish Balay 
20749b5e25fSSatish Balay   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
20849b5e25fSSatish Balay   ierr = PetscFree(vary);CHKERRQ(ierr);
20949b5e25fSSatish Balay 
21049b5e25fSSatish Balay   ierr = MatGetSubMatrix_SeqSBAIJ_Private(A,is1,is1,cs,scall,B);CHKERRQ(ierr);
21149b5e25fSSatish Balay   ISDestroy(is1);
21249b5e25fSSatish Balay   PetscFunctionReturn(0);
21349b5e25fSSatish Balay }
21449b5e25fSSatish Balay 
2154a2ae208SSatish Balay #undef __FUNCT__
2164a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_SeqSBAIJ"
21713f74950SBarry Smith PetscErrorCode MatGetSubMatrices_SeqSBAIJ(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
21849b5e25fSSatish Balay {
2196849ba73SBarry Smith   PetscErrorCode ierr;
22013f74950SBarry Smith   PetscInt       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"
240dfbe8321SBarry Smith PetscErrorCode 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;
2456849ba73SBarry Smith   PetscErrorCode ierr;
24613f74950SBarry Smith   PetscInt       mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,*ib,cval,j,jmin;
24749b5e25fSSatish Balay 
24849b5e25fSSatish Balay   PetscFunctionBegin;
2492dcb1b2aSMatthew Knepley   ierr = VecSet(zz,zero);CHKERRQ(ierr);
2501ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2511ebc52fbSHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
25249b5e25fSSatish Balay 
25349b5e25fSSatish Balay   v  = a->a;
25449b5e25fSSatish Balay   xb = x;
25549b5e25fSSatish Balay 
25649b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
25749b5e25fSSatish Balay     n  = ai[1] - ai[0];  /* length of i_th row of A */
25849b5e25fSSatish Balay     x1 = xb[0];
25949b5e25fSSatish Balay     ib = aj + *ai;
260831a3094SHong Zhang     jmin = 0;
261831a3094SHong Zhang     if (*ib == i) {      /* (diag of A)*x */
262831a3094SHong Zhang       z[i] += *v++ * x[*ib++];
263831a3094SHong Zhang       jmin++;
264831a3094SHong Zhang     }
265831a3094SHong Zhang     for (j=jmin; j<n; j++) {
26649b5e25fSSatish Balay       cval    = *ib;
26749b5e25fSSatish Balay       z[cval] += *v * x1;      /* (strict lower triangular part of A)*x  */
26849b5e25fSSatish Balay       z[i] += *v++ * x[*ib++]; /* (strict upper triangular part of A)*x  */
26949b5e25fSSatish Balay     }
27049b5e25fSSatish Balay     xb++; ai++;
27149b5e25fSSatish Balay   }
27249b5e25fSSatish Balay 
2731ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2741ebc52fbSHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
275efee365bSSatish Balay   ierr = PetscLogFlops(2*(a->nz*2 - A->m) - A->m);CHKERRQ(ierr);  /* nz = (nz+m)/2 */
27649b5e25fSSatish Balay   PetscFunctionReturn(0);
27749b5e25fSSatish Balay }
27849b5e25fSSatish Balay 
2794a2ae208SSatish Balay #undef __FUNCT__
2804a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqSBAIJ_2"
281dfbe8321SBarry Smith PetscErrorCode MatMult_SeqSBAIJ_2(Mat A,Vec xx,Vec zz)
28249b5e25fSSatish Balay {
28349b5e25fSSatish Balay   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
28487828ca2SBarry Smith   PetscScalar    *x,*z,*xb,x1,x2,zero=0.0;
28549b5e25fSSatish Balay   MatScalar      *v;
2866849ba73SBarry Smith   PetscErrorCode ierr;
28713f74950SBarry Smith   PetscInt       mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,*ib,cval,j,jmin;
28849b5e25fSSatish Balay 
28949b5e25fSSatish Balay   PetscFunctionBegin;
2902dcb1b2aSMatthew Knepley   ierr = VecSet(zz,zero);CHKERRQ(ierr);
2911ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2921ebc52fbSHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
29349b5e25fSSatish Balay 
29449b5e25fSSatish Balay   v     = a->a;
29549b5e25fSSatish Balay   xb = x;
29649b5e25fSSatish Balay 
29749b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
29849b5e25fSSatish Balay     n  = ai[1] - ai[0]; /* length of i_th block row of A */
29949b5e25fSSatish Balay     x1 = xb[0]; x2 = xb[1];
30049b5e25fSSatish Balay     ib = aj + *ai;
301831a3094SHong Zhang     jmin = 0;
3027fbae186SHong Zhang     if (*ib == i){     /* (diag of A)*x */
30349b5e25fSSatish Balay       z[2*i]   += v[0]*x1 + v[2]*x2;
30449b5e25fSSatish Balay       z[2*i+1] += v[2]*x1 + v[3]*x2;
305831a3094SHong Zhang       v += 4; jmin++;
3067fbae186SHong Zhang     }
307831a3094SHong Zhang     for (j=jmin; j<n; j++) {
30849b5e25fSSatish Balay       /* (strict lower triangular part of A)*x  */
30949b5e25fSSatish Balay       cval       = ib[j]*2;
31049b5e25fSSatish Balay       z[cval]     += v[0]*x1 + v[1]*x2;
31149b5e25fSSatish Balay       z[cval+1]   += v[2]*x1 + v[3]*x2;
31249b5e25fSSatish Balay       /* (strict upper triangular part of A)*x  */
31349b5e25fSSatish Balay       z[2*i]   += v[0]*x[cval] + v[2]*x[cval+1];
31449b5e25fSSatish Balay       z[2*i+1] += v[1]*x[cval] + v[3]*x[cval+1];
31549b5e25fSSatish Balay       v  += 4;
31649b5e25fSSatish Balay     }
31749b5e25fSSatish Balay     xb +=2; ai++;
31849b5e25fSSatish Balay   }
31949b5e25fSSatish Balay 
3201ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
3211ebc52fbSHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
322efee365bSSatish Balay   ierr = PetscLogFlops(8*(a->nz*2 - A->m) - A->m);CHKERRQ(ierr);
32349b5e25fSSatish Balay   PetscFunctionReturn(0);
32449b5e25fSSatish Balay }
32549b5e25fSSatish Balay 
3264a2ae208SSatish Balay #undef __FUNCT__
3274a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqSBAIJ_3"
328dfbe8321SBarry Smith PetscErrorCode MatMult_SeqSBAIJ_3(Mat A,Vec xx,Vec zz)
32949b5e25fSSatish Balay {
33049b5e25fSSatish Balay   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
33187828ca2SBarry Smith   PetscScalar    *x,*z,*xb,x1,x2,x3,zero=0.0;
33249b5e25fSSatish Balay   MatScalar      *v;
3336849ba73SBarry Smith   PetscErrorCode ierr;
33413f74950SBarry Smith   PetscInt       mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,*ib,cval,j,jmin;
33549b5e25fSSatish Balay 
33649b5e25fSSatish Balay   PetscFunctionBegin;
3372dcb1b2aSMatthew Knepley   ierr = VecSet(zz,zero);CHKERRQ(ierr);
3381ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
3391ebc52fbSHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
34049b5e25fSSatish Balay 
34149b5e25fSSatish Balay   v    = a->a;
34249b5e25fSSatish Balay   xb   = x;
34349b5e25fSSatish Balay 
34449b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
34549b5e25fSSatish Balay     n  = ai[1] - ai[0]; /* length of i_th block row of A */
34649b5e25fSSatish Balay     x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
34749b5e25fSSatish Balay     ib = aj + *ai;
348831a3094SHong Zhang     jmin = 0;
3497fbae186SHong Zhang     if (*ib == i){     /* (diag of A)*x */
35049b5e25fSSatish Balay       z[3*i]   += v[0]*x1 + v[3]*x2 + v[6]*x3;
35149b5e25fSSatish Balay       z[3*i+1] += v[3]*x1 + v[4]*x2 + v[7]*x3;
35249b5e25fSSatish Balay       z[3*i+2] += v[6]*x1 + v[7]*x2 + v[8]*x3;
353831a3094SHong Zhang       v += 9; jmin++;
3547fbae186SHong Zhang     }
355831a3094SHong Zhang     for (j=jmin; j<n; j++) {
35649b5e25fSSatish Balay       /* (strict lower triangular part of A)*x  */
35749b5e25fSSatish Balay       cval       = ib[j]*3;
35849b5e25fSSatish Balay       z[cval]     += v[0]*x1 + v[1]*x2 + v[2]*x3;
35949b5e25fSSatish Balay       z[cval+1]   += v[3]*x1 + v[4]*x2 + v[5]*x3;
36049b5e25fSSatish Balay       z[cval+2]   += v[6]*x1 + v[7]*x2 + v[8]*x3;
36149b5e25fSSatish Balay       /* (strict upper triangular part of A)*x  */
36249b5e25fSSatish Balay       z[3*i]   += v[0]*x[cval] + v[3]*x[cval+1]+ v[6]*x[cval+2];
36349b5e25fSSatish Balay       z[3*i+1] += v[1]*x[cval] + v[4]*x[cval+1]+ v[7]*x[cval+2];
36449b5e25fSSatish Balay       z[3*i+2] += v[2]*x[cval] + v[5]*x[cval+1]+ v[8]*x[cval+2];
36549b5e25fSSatish Balay       v  += 9;
36649b5e25fSSatish Balay     }
36749b5e25fSSatish Balay     xb +=3; ai++;
36849b5e25fSSatish Balay   }
36949b5e25fSSatish Balay 
3701ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
3711ebc52fbSHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
372efee365bSSatish Balay   ierr = PetscLogFlops(18*(a->nz*2 - A->m) - A->m);CHKERRQ(ierr);
37349b5e25fSSatish Balay   PetscFunctionReturn(0);
37449b5e25fSSatish Balay }
37549b5e25fSSatish Balay 
3764a2ae208SSatish Balay #undef __FUNCT__
3774a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqSBAIJ_4"
378dfbe8321SBarry Smith PetscErrorCode MatMult_SeqSBAIJ_4(Mat A,Vec xx,Vec zz)
37949b5e25fSSatish Balay {
38049b5e25fSSatish Balay   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
38187828ca2SBarry Smith   PetscScalar    *x,*z,*xb,x1,x2,x3,x4,zero=0.0;
38249b5e25fSSatish Balay   MatScalar      *v;
3836849ba73SBarry Smith   PetscErrorCode ierr;
38413f74950SBarry Smith   PetscInt       mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,*ib,cval,j,jmin;
38549b5e25fSSatish Balay 
38649b5e25fSSatish Balay   PetscFunctionBegin;
3872dcb1b2aSMatthew Knepley   ierr = VecSet(zz,zero);CHKERRQ(ierr);
3881ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
3891ebc52fbSHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
39049b5e25fSSatish Balay 
39149b5e25fSSatish Balay   v     = a->a;
39249b5e25fSSatish Balay   xb = x;
39349b5e25fSSatish Balay 
39449b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
39549b5e25fSSatish Balay     n  = ai[1] - ai[0]; /* length of i_th block row of A */
39649b5e25fSSatish Balay     x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
39749b5e25fSSatish Balay     ib = aj + *ai;
398831a3094SHong Zhang     jmin = 0;
3997fbae186SHong Zhang     if (*ib == i){     /* (diag of A)*x */
40049b5e25fSSatish Balay       z[4*i]   += v[0]*x1 + v[4]*x2 +  v[8]*x3 + v[12]*x4;
40149b5e25fSSatish Balay       z[4*i+1] += v[4]*x1 + v[5]*x2 +  v[9]*x3 + v[13]*x4;
40249b5e25fSSatish Balay       z[4*i+2] += v[8]*x1 + v[9]*x2 + v[10]*x3 + v[14]*x4;
40349b5e25fSSatish Balay       z[4*i+3] += v[12]*x1+ v[13]*x2+ v[14]*x3 + v[15]*x4;
404831a3094SHong Zhang       v += 16; jmin++;
4057fbae186SHong Zhang     }
406831a3094SHong Zhang     for (j=jmin; j<n; j++) {
40749b5e25fSSatish Balay       /* (strict lower triangular part of A)*x  */
40849b5e25fSSatish Balay       cval       = ib[j]*4;
40949b5e25fSSatish Balay       z[cval]     += v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4;
41049b5e25fSSatish Balay       z[cval+1]   += v[4]*x1 + v[5]*x2 + v[6]*x3 + v[7]*x4;
41149b5e25fSSatish Balay       z[cval+2]   += v[8]*x1 + v[9]*x2 + v[10]*x3 + v[11]*x4;
41249b5e25fSSatish Balay       z[cval+3]   += v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4;
41349b5e25fSSatish Balay       /* (strict upper triangular part of A)*x  */
41449b5e25fSSatish Balay       z[4*i]   += v[0]*x[cval] + v[4]*x[cval+1]+ v[8]*x[cval+2] + v[12]*x[cval+3];
41549b5e25fSSatish Balay       z[4*i+1] += v[1]*x[cval] + v[5]*x[cval+1]+ v[9]*x[cval+2] + v[13]*x[cval+3];
41649b5e25fSSatish Balay       z[4*i+2] += v[2]*x[cval] + v[6]*x[cval+1]+ v[10]*x[cval+2]+ v[14]*x[cval+3];
41749b5e25fSSatish Balay       z[4*i+3] += v[3]*x[cval] + v[7]*x[cval+1]+ v[11]*x[cval+2]+ v[15]*x[cval+3];
41849b5e25fSSatish Balay       v  += 16;
41949b5e25fSSatish Balay     }
42049b5e25fSSatish Balay     xb +=4; ai++;
42149b5e25fSSatish Balay   }
42249b5e25fSSatish Balay 
4231ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
4241ebc52fbSHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
425efee365bSSatish Balay   ierr = PetscLogFlops(32*(a->nz*2 - A->m) - A->m);CHKERRQ(ierr);
42649b5e25fSSatish Balay   PetscFunctionReturn(0);
42749b5e25fSSatish Balay }
42849b5e25fSSatish Balay 
4294a2ae208SSatish Balay #undef __FUNCT__
4304a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqSBAIJ_5"
431dfbe8321SBarry Smith PetscErrorCode MatMult_SeqSBAIJ_5(Mat A,Vec xx,Vec zz)
43249b5e25fSSatish Balay {
43349b5e25fSSatish Balay   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
43487828ca2SBarry Smith   PetscScalar    *x,*z,*xb,x1,x2,x3,x4,x5,zero=0.0;
43549b5e25fSSatish Balay   MatScalar      *v;
4366849ba73SBarry Smith   PetscErrorCode ierr;
43713f74950SBarry Smith   PetscInt       mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,*ib,cval,j,jmin;
43849b5e25fSSatish Balay 
43949b5e25fSSatish Balay   PetscFunctionBegin;
4402dcb1b2aSMatthew Knepley   ierr = VecSet(zz,zero);CHKERRQ(ierr);
4411ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
4421ebc52fbSHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
44349b5e25fSSatish Balay 
44449b5e25fSSatish Balay   v     = a->a;
44549b5e25fSSatish Balay   xb = x;
44649b5e25fSSatish Balay 
44749b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
44849b5e25fSSatish Balay     n  = ai[1] - ai[0]; /* length of i_th block row of A */
44949b5e25fSSatish Balay     x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5=xb[4];
45049b5e25fSSatish Balay     ib = aj + *ai;
451831a3094SHong Zhang     jmin = 0;
4527fbae186SHong Zhang     if (*ib == i){      /* (diag of A)*x */
45349b5e25fSSatish Balay       z[5*i]   += v[0]*x1  + v[5]*x2 + v[10]*x3 + v[15]*x4+ v[20]*x5;
45449b5e25fSSatish Balay       z[5*i+1] += v[5]*x1  + v[6]*x2 + v[11]*x3 + v[16]*x4+ v[21]*x5;
45549b5e25fSSatish Balay       z[5*i+2] += v[10]*x1 +v[11]*x2 + v[12]*x3 + v[17]*x4+ v[22]*x5;
45649b5e25fSSatish Balay       z[5*i+3] += v[15]*x1 +v[16]*x2 + v[17]*x3 + v[18]*x4+ v[23]*x5;
45749b5e25fSSatish Balay       z[5*i+4] += v[20]*x1 +v[21]*x2 + v[22]*x3 + v[23]*x4+ v[24]*x5;
458831a3094SHong Zhang       v += 25; jmin++;
4597fbae186SHong Zhang     }
460831a3094SHong Zhang     for (j=jmin; j<n; j++) {
46149b5e25fSSatish Balay       /* (strict lower triangular part of A)*x  */
46249b5e25fSSatish Balay       cval       = ib[j]*5;
46349b5e25fSSatish Balay       z[cval]     += v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4 + v[4]*x5;
46449b5e25fSSatish Balay       z[cval+1]   += v[5]*x1 + v[6]*x2 + v[7]*x3 + v[8]*x4 + v[9]*x5;
46549b5e25fSSatish Balay       z[cval+2]   += v[10]*x1 + v[11]*x2 + v[12]*x3 + v[13]*x4+ v[14]*x5;
46649b5e25fSSatish Balay       z[cval+3]   += v[15]*x1 + v[16]*x2 + v[17]*x3 + v[18]*x4+ v[19]*x5;
46749b5e25fSSatish Balay       z[cval+4]   += v[20]*x1 + v[21]*x2 + v[22]*x3 + v[23]*x4+ v[24]*x5;
46849b5e25fSSatish Balay       /* (strict upper triangular part of A)*x  */
46949b5e25fSSatish 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];
47049b5e25fSSatish 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];
47149b5e25fSSatish 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];
47249b5e25fSSatish 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];
47349b5e25fSSatish 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];
47449b5e25fSSatish Balay       v  += 25;
47549b5e25fSSatish Balay     }
47649b5e25fSSatish Balay     xb +=5; ai++;
47749b5e25fSSatish Balay   }
47849b5e25fSSatish Balay 
4791ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
4801ebc52fbSHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
481efee365bSSatish Balay   ierr = PetscLogFlops(50*(a->nz*2 - A->m) - A->m);CHKERRQ(ierr);
48249b5e25fSSatish Balay   PetscFunctionReturn(0);
48349b5e25fSSatish Balay }
48449b5e25fSSatish Balay 
48549b5e25fSSatish Balay 
4864a2ae208SSatish Balay #undef __FUNCT__
4874a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqSBAIJ_6"
488dfbe8321SBarry Smith PetscErrorCode MatMult_SeqSBAIJ_6(Mat A,Vec xx,Vec zz)
48949b5e25fSSatish Balay {
49049b5e25fSSatish Balay   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
49187828ca2SBarry Smith   PetscScalar    *x,*z,*xb,x1,x2,x3,x4,x5,x6,zero=0.0;
49249b5e25fSSatish Balay   MatScalar      *v;
4936849ba73SBarry Smith   PetscErrorCode ierr;
49413f74950SBarry Smith   PetscInt       mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,*ib,cval,j,jmin;
49549b5e25fSSatish Balay 
49649b5e25fSSatish Balay   PetscFunctionBegin;
4972dcb1b2aSMatthew Knepley   ierr = VecSet(zz,zero);CHKERRQ(ierr);
4981ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
4991ebc52fbSHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
50049b5e25fSSatish Balay 
50149b5e25fSSatish Balay   v     = a->a;
50249b5e25fSSatish Balay   xb = x;
50349b5e25fSSatish Balay 
50449b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
50549b5e25fSSatish Balay     n  = ai[1] - ai[0]; /* length of i_th block row of A */
50649b5e25fSSatish Balay     x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5=xb[4]; x6=xb[5];
50749b5e25fSSatish Balay     ib = aj + *ai;
508831a3094SHong Zhang     jmin = 0;
5097fbae186SHong Zhang     if (*ib == i){      /* (diag of A)*x */
51049b5e25fSSatish Balay       z[6*i]   += v[0]*x1  + v[6]*x2 + v[12]*x3 + v[18]*x4+ v[24]*x5 + v[30]*x6;
51149b5e25fSSatish Balay       z[6*i+1] += v[6]*x1  + v[7]*x2 + v[13]*x3 + v[19]*x4+ v[25]*x5 + v[31]*x6;
51249b5e25fSSatish Balay       z[6*i+2] += v[12]*x1 +v[13]*x2 + v[14]*x3 + v[20]*x4+ v[26]*x5 + v[32]*x6;
51349b5e25fSSatish Balay       z[6*i+3] += v[18]*x1 +v[19]*x2 + v[20]*x3 + v[21]*x4+ v[27]*x5 + v[33]*x6;
51449b5e25fSSatish Balay       z[6*i+4] += v[24]*x1 +v[25]*x2 + v[26]*x3 + v[27]*x4+ v[28]*x5 + v[34]*x6;
51549b5e25fSSatish Balay       z[6*i+5] += v[30]*x1 +v[31]*x2 + v[32]*x3 + v[33]*x4+ v[34]*x5 + v[35]*x6;
516831a3094SHong Zhang       v += 36; jmin++;
5177fbae186SHong Zhang     }
518831a3094SHong Zhang     for (j=jmin; j<n; j++) {
51949b5e25fSSatish Balay       /* (strict lower triangular part of A)*x  */
52049b5e25fSSatish Balay       cval       = ib[j]*6;
52149b5e25fSSatish Balay       z[cval]   += v[0]*x1  + v[1]*x2 + v[2]*x3 + v[3]*x4+ v[4]*x5 + v[5]*x6;
52249b5e25fSSatish Balay       z[cval+1] += v[6]*x1  + v[7]*x2 + v[8]*x3 + v[9]*x4+ v[10]*x5 + v[11]*x6;
52349b5e25fSSatish Balay       z[cval+2] += v[12]*x1  + v[13]*x2 + v[14]*x3 + v[15]*x4+ v[16]*x5 + v[17]*x6;
52449b5e25fSSatish Balay       z[cval+3] += v[18]*x1  + v[19]*x2 + v[20]*x3 + v[21]*x4+ v[22]*x5 + v[23]*x6;
52549b5e25fSSatish Balay       z[cval+4] += v[24]*x1  + v[25]*x2 + v[26]*x3 + v[27]*x4+ v[28]*x5 + v[29]*x6;
52649b5e25fSSatish Balay       z[cval+5] += v[30]*x1  + v[31]*x2 + v[32]*x3 + v[33]*x4+ v[34]*x5 + v[35]*x6;
52749b5e25fSSatish Balay       /* (strict upper triangular part of A)*x  */
52849b5e25fSSatish 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];
52949b5e25fSSatish 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];
53049b5e25fSSatish 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];
53149b5e25fSSatish 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];
53249b5e25fSSatish 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];
53349b5e25fSSatish 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];
53449b5e25fSSatish Balay       v  += 36;
53549b5e25fSSatish Balay     }
53649b5e25fSSatish Balay     xb +=6; ai++;
53749b5e25fSSatish Balay   }
53849b5e25fSSatish Balay 
5391ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
5401ebc52fbSHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
541efee365bSSatish Balay   ierr = PetscLogFlops(72*(a->nz*2 - A->m) - A->m);CHKERRQ(ierr);
54249b5e25fSSatish Balay   PetscFunctionReturn(0);
54349b5e25fSSatish Balay }
5444a2ae208SSatish Balay #undef __FUNCT__
5454a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqSBAIJ_7"
546dfbe8321SBarry Smith PetscErrorCode MatMult_SeqSBAIJ_7(Mat A,Vec xx,Vec zz)
54749b5e25fSSatish Balay {
54849b5e25fSSatish Balay   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
54987828ca2SBarry Smith   PetscScalar    *x,*z,*xb,x1,x2,x3,x4,x5,x6,x7,zero=0.0;
55049b5e25fSSatish Balay   MatScalar      *v;
5516849ba73SBarry Smith   PetscErrorCode ierr;
55213f74950SBarry Smith   PetscInt       mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,*ib,cval,j,jmin;
55349b5e25fSSatish Balay 
55449b5e25fSSatish Balay   PetscFunctionBegin;
5552dcb1b2aSMatthew Knepley   ierr = VecSet(zz,zero);CHKERRQ(ierr);
5561ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
5571ebc52fbSHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
55849b5e25fSSatish Balay 
55949b5e25fSSatish Balay   v     = a->a;
56049b5e25fSSatish Balay   xb = x;
56149b5e25fSSatish Balay 
56249b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
56349b5e25fSSatish Balay     n  = ai[1] - ai[0]; /* length of i_th block row of A */
56449b5e25fSSatish Balay     x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5=xb[4]; x6=xb[5]; x7=xb[6];
56549b5e25fSSatish Balay     ib = aj + *ai;
566831a3094SHong Zhang     jmin = 0;
5677fbae186SHong Zhang     if (*ib == i){      /* (diag of A)*x */
56849b5e25fSSatish 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;
56949b5e25fSSatish 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;
57049b5e25fSSatish 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;
57149b5e25fSSatish 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;
57249b5e25fSSatish 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;
57349b5e25fSSatish 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;
57449b5e25fSSatish 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;
575831a3094SHong Zhang       v += 49; jmin++;
5767fbae186SHong Zhang     }
577831a3094SHong Zhang     for (j=jmin; j<n; j++) {
57849b5e25fSSatish Balay       /* (strict lower triangular part of A)*x  */
57949b5e25fSSatish Balay       cval       = ib[j]*7;
58049b5e25fSSatish Balay       z[cval]   += v[0]*x1  + v[1]*x2 + v[2]*x3 + v[3]*x4+ v[4]*x5 + v[5]*x6+ v[6]*x7;
58149b5e25fSSatish 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;
58249b5e25fSSatish 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;
58349b5e25fSSatish 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;
58449b5e25fSSatish 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;
58549b5e25fSSatish 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;
58649b5e25fSSatish 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;
58749b5e25fSSatish Balay       /* (strict upper triangular part of A)*x  */
58849b5e25fSSatish 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];
58949b5e25fSSatish 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];
59049b5e25fSSatish 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];
59149b5e25fSSatish 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];
59249b5e25fSSatish 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];
59349b5e25fSSatish 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];
59449b5e25fSSatish 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];
59549b5e25fSSatish Balay       v  += 49;
59649b5e25fSSatish Balay     }
59749b5e25fSSatish Balay     xb +=7; ai++;
59849b5e25fSSatish Balay   }
5991ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
6001ebc52fbSHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
601efee365bSSatish Balay   ierr = PetscLogFlops(98*(a->nz*2 - A->m) - A->m);CHKERRQ(ierr);
60249b5e25fSSatish Balay   PetscFunctionReturn(0);
60349b5e25fSSatish Balay }
60449b5e25fSSatish Balay 
60549b5e25fSSatish Balay /*
60649b5e25fSSatish Balay     This will not work with MatScalar == float because it calls the BLAS
60749b5e25fSSatish Balay */
6084a2ae208SSatish Balay #undef __FUNCT__
6094a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqSBAIJ_N"
610dfbe8321SBarry Smith PetscErrorCode MatMult_SeqSBAIJ_N(Mat A,Vec xx,Vec zz)
61149b5e25fSSatish Balay {
61249b5e25fSSatish Balay   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
61387828ca2SBarry Smith   PetscScalar    *x,*x_ptr,*z,*z_ptr,*xb,*zb,*work,*workt,zero=0.0;
6140b60a74dSHong Zhang   MatScalar      *v;
615dfbe8321SBarry Smith   PetscErrorCode ierr;
616521d7252SBarry Smith   PetscInt       mbs=a->mbs,i,*idx,*aj,*ii,bs=A->bs,j,n,bs2=a->bs2,ncols,k;
61749b5e25fSSatish Balay 
61849b5e25fSSatish Balay   PetscFunctionBegin;
6192dcb1b2aSMatthew Knepley   ierr = VecSet(zz,zero);CHKERRQ(ierr);
6201ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr); x_ptr=x;
6211ebc52fbSHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr); z_ptr=z;
62249b5e25fSSatish Balay 
62349b5e25fSSatish Balay   aj   = a->j;
62449b5e25fSSatish Balay   v    = a->a;
62549b5e25fSSatish Balay   ii   = a->i;
62649b5e25fSSatish Balay 
62749b5e25fSSatish Balay   if (!a->mult_work) {
62887828ca2SBarry Smith     ierr = PetscMalloc((A->m+1)*sizeof(PetscScalar),&a->mult_work);CHKERRQ(ierr);
62949b5e25fSSatish Balay   }
63049b5e25fSSatish Balay   work = a->mult_work;
63149b5e25fSSatish Balay 
63249b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
63349b5e25fSSatish Balay     n     = ii[1] - ii[0]; ncols = n*bs;
63449b5e25fSSatish Balay     workt = work; idx=aj+ii[0];
63549b5e25fSSatish Balay 
63649b5e25fSSatish Balay     /* upper triangular part */
63749b5e25fSSatish Balay     for (j=0; j<n; j++) {
63849b5e25fSSatish Balay       xb = x_ptr + bs*(*idx++);
63949b5e25fSSatish Balay       for (k=0; k<bs; k++) workt[k] = xb[k];
64049b5e25fSSatish Balay       workt += bs;
64149b5e25fSSatish Balay     }
64249b5e25fSSatish Balay     /* z(i*bs:(i+1)*bs-1) += A(i,:)*x */
64349b5e25fSSatish Balay     Kernel_w_gets_w_plus_Ar_times_v(bs,ncols,work,v,z);
64449b5e25fSSatish Balay 
64549b5e25fSSatish Balay     /* strict lower triangular part */
646831a3094SHong Zhang     idx = aj+ii[0];
647831a3094SHong Zhang     if (*idx == i){
64896b9376eSHong Zhang       ncols -= bs; v += bs2; idx++; n--;
649831a3094SHong Zhang     }
65096b9376eSHong Zhang 
65149b5e25fSSatish Balay     if (ncols > 0){
65249b5e25fSSatish Balay       workt = work;
65387828ca2SBarry Smith       ierr  = PetscMemzero(workt,ncols*sizeof(PetscScalar));CHKERRQ(ierr);
654831a3094SHong Zhang       Kernel_w_gets_w_plus_trans_Ar_times_v(bs,ncols,x,v,workt);
655831a3094SHong Zhang       for (j=0; j<n; j++) {
656831a3094SHong Zhang         zb = z_ptr + bs*(*idx++);
65749b5e25fSSatish Balay         for (k=0; k<bs; k++) zb[k] += workt[k] ;
65849b5e25fSSatish Balay         workt += bs;
65949b5e25fSSatish Balay       }
66049b5e25fSSatish Balay     }
66149b5e25fSSatish Balay     x += bs; v += n*bs2; z += bs; ii++;
66249b5e25fSSatish Balay   }
66349b5e25fSSatish Balay 
6641ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
6651ebc52fbSHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
666efee365bSSatish Balay   ierr = PetscLogFlops(2*(a->nz*2 - A->m)*bs2 - A->m);CHKERRQ(ierr);
66749b5e25fSSatish Balay   PetscFunctionReturn(0);
66849b5e25fSSatish Balay }
66949b5e25fSSatish Balay 
6704a2ae208SSatish Balay #undef __FUNCT__
6714a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqSBAIJ_1"
672dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqSBAIJ_1(Mat A,Vec xx,Vec yy,Vec zz)
67349b5e25fSSatish Balay {
67449b5e25fSSatish Balay   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
67587828ca2SBarry Smith   PetscScalar    *x,*y,*z,*xb,x1;
67649b5e25fSSatish Balay   MatScalar      *v;
6776849ba73SBarry Smith   PetscErrorCode ierr;
67813f74950SBarry Smith   PetscInt       mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,*ib,cval,j,jmin;
67949b5e25fSSatish Balay 
68049b5e25fSSatish Balay   PetscFunctionBegin;
6811ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
68249b5e25fSSatish Balay   if (yy != xx) {
6831ebc52fbSHong Zhang     ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
68449b5e25fSSatish Balay   } else {
68549b5e25fSSatish Balay     y = x;
68649b5e25fSSatish Balay   }
68749b5e25fSSatish Balay   if (zz != yy) {
688a9d4b620SHong Zhang     /* ierr = VecCopy(yy,zz);CHKERRQ(ierr); */
6891ebc52fbSHong Zhang     ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
690a9d4b620SHong Zhang     ierr = PetscMemcpy(z,y,yy->n*sizeof(PetscScalar));CHKERRQ(ierr);
69149b5e25fSSatish Balay   } else {
69249b5e25fSSatish Balay     z = y;
69349b5e25fSSatish Balay   }
69449b5e25fSSatish Balay 
69549b5e25fSSatish Balay   v  = a->a;
69649b5e25fSSatish Balay   xb = x;
69749b5e25fSSatish Balay 
69849b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
69949b5e25fSSatish Balay     n  = ai[1] - ai[0];  /* length of i_th row of A */
70049b5e25fSSatish Balay     x1 = xb[0];
70149b5e25fSSatish Balay     ib = aj + *ai;
702831a3094SHong Zhang     jmin = 0;
703831a3094SHong Zhang     if (*ib == i) {            /* (diag of A)*x */
704831a3094SHong Zhang       z[i] += *v++ * x[*ib++]; jmin++;
705831a3094SHong Zhang     }
706831a3094SHong Zhang     for (j=jmin; j<n; j++) {
70749b5e25fSSatish Balay       cval    = *ib;
70849b5e25fSSatish Balay       z[cval] += *v * x1;      /* (strict lower triangular part of A)*x  */
70949b5e25fSSatish Balay       z[i] += *v++ * x[*ib++]; /* (strict upper triangular part of A)*x  */
71049b5e25fSSatish Balay     }
71149b5e25fSSatish Balay     xb++; ai++;
71249b5e25fSSatish Balay   }
71349b5e25fSSatish Balay 
7141ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
7151ebc52fbSHong Zhang   if (yy != xx) ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
7161ebc52fbSHong Zhang   if (zz != yy) ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
71749b5e25fSSatish Balay 
718efee365bSSatish Balay   ierr = PetscLogFlops(2*(a->nz*2 - A->m));CHKERRQ(ierr);
71949b5e25fSSatish Balay   PetscFunctionReturn(0);
72049b5e25fSSatish Balay }
72149b5e25fSSatish Balay 
7224a2ae208SSatish Balay #undef __FUNCT__
7234a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqSBAIJ_2"
724dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqSBAIJ_2(Mat A,Vec xx,Vec yy,Vec zz)
72549b5e25fSSatish Balay {
72649b5e25fSSatish Balay   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
72787828ca2SBarry Smith   PetscScalar    *x,*y,*z,*xb,x1,x2;
72849b5e25fSSatish Balay   MatScalar      *v;
7296849ba73SBarry Smith   PetscErrorCode ierr;
73013f74950SBarry Smith   PetscInt       mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,*ib,cval,j,jmin;
73149b5e25fSSatish Balay 
73249b5e25fSSatish Balay   PetscFunctionBegin;
7331ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
73449b5e25fSSatish Balay   if (yy != xx) {
7351ebc52fbSHong Zhang     ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
73649b5e25fSSatish Balay   } else {
73749b5e25fSSatish Balay     y = x;
73849b5e25fSSatish Balay   }
73949b5e25fSSatish Balay   if (zz != yy) {
740a9d4b620SHong Zhang     /* ierr = VecCopy(yy,zz);CHKERRQ(ierr); */
7411ebc52fbSHong Zhang     ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
742a9d4b620SHong Zhang     ierr = PetscMemcpy(z,y,yy->n*sizeof(PetscScalar));CHKERRQ(ierr);
74349b5e25fSSatish Balay   } else {
74449b5e25fSSatish Balay     z = y;
74549b5e25fSSatish Balay   }
74649b5e25fSSatish Balay 
74749b5e25fSSatish Balay   v     = a->a;
74849b5e25fSSatish Balay   xb = x;
74949b5e25fSSatish Balay 
75049b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
75149b5e25fSSatish Balay     n  = ai[1] - ai[0]; /* length of i_th block row of A */
75249b5e25fSSatish Balay     x1 = xb[0]; x2 = xb[1];
75349b5e25fSSatish Balay     ib = aj + *ai;
754831a3094SHong Zhang     jmin = 0;
7557fbae186SHong Zhang     if (*ib == i){      /* (diag of A)*x */
75649b5e25fSSatish Balay       z[2*i]   += v[0]*x1 + v[2]*x2;
75749b5e25fSSatish Balay       z[2*i+1] += v[2]*x1 + v[3]*x2;
758831a3094SHong Zhang       v += 4; jmin++;
7597fbae186SHong Zhang     }
760831a3094SHong Zhang     for (j=jmin; j<n; j++) {
76149b5e25fSSatish Balay       /* (strict lower triangular part of A)*x  */
76249b5e25fSSatish Balay       cval       = ib[j]*2;
76349b5e25fSSatish Balay       z[cval]     += v[0]*x1 + v[1]*x2;
76449b5e25fSSatish Balay       z[cval+1]   += v[2]*x1 + v[3]*x2;
76549b5e25fSSatish Balay       /* (strict upper triangular part of A)*x  */
76649b5e25fSSatish Balay       z[2*i]   += v[0]*x[cval] + v[2]*x[cval+1];
76749b5e25fSSatish Balay       z[2*i+1] += v[1]*x[cval] + v[3]*x[cval+1];
76849b5e25fSSatish Balay       v  += 4;
76949b5e25fSSatish Balay     }
77049b5e25fSSatish Balay     xb +=2; ai++;
77149b5e25fSSatish Balay   }
77249b5e25fSSatish Balay 
7731ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
7741ebc52fbSHong Zhang   if (yy != xx) ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
7751ebc52fbSHong Zhang   if (zz != yy) ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
77649b5e25fSSatish Balay 
777efee365bSSatish Balay   ierr = PetscLogFlops(4*(a->nz*2 - A->m));CHKERRQ(ierr);
77849b5e25fSSatish Balay   PetscFunctionReturn(0);
77949b5e25fSSatish Balay }
78049b5e25fSSatish Balay 
7814a2ae208SSatish Balay #undef __FUNCT__
7824a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqSBAIJ_3"
783dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqSBAIJ_3(Mat A,Vec xx,Vec yy,Vec zz)
78449b5e25fSSatish Balay {
78549b5e25fSSatish Balay   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
78687828ca2SBarry Smith   PetscScalar    *x,*y,*z,*xb,x1,x2,x3;
78749b5e25fSSatish Balay   MatScalar      *v;
7886849ba73SBarry Smith   PetscErrorCode ierr;
78913f74950SBarry Smith   PetscInt       mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,*ib,cval,j,jmin;
79049b5e25fSSatish Balay 
79149b5e25fSSatish Balay   PetscFunctionBegin;
7921ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
79349b5e25fSSatish Balay   if (yy != xx) {
7941ebc52fbSHong Zhang     ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
79549b5e25fSSatish Balay   } else {
79649b5e25fSSatish Balay     y = x;
79749b5e25fSSatish Balay   }
79849b5e25fSSatish Balay   if (zz != yy) {
799a9d4b620SHong Zhang     /* ierr = VecCopy(yy,zz);CHKERRQ(ierr); */
8001ebc52fbSHong Zhang     ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
801a9d4b620SHong Zhang     ierr = PetscMemcpy(z,y,yy->n*sizeof(PetscScalar));CHKERRQ(ierr);
80249b5e25fSSatish Balay   } else {
80349b5e25fSSatish Balay     z = y;
80449b5e25fSSatish Balay   }
80549b5e25fSSatish Balay 
80649b5e25fSSatish Balay   v     = a->a;
80749b5e25fSSatish Balay   xb = x;
80849b5e25fSSatish Balay 
80949b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
81049b5e25fSSatish Balay     n  = ai[1] - ai[0]; /* length of i_th block row of A */
81149b5e25fSSatish Balay     x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
81249b5e25fSSatish Balay     ib = aj + *ai;
813831a3094SHong Zhang     jmin = 0;
8147fbae186SHong Zhang     if (*ib == i){     /* (diag of A)*x */
81549b5e25fSSatish Balay      z[3*i]   += v[0]*x1 + v[3]*x2 + v[6]*x3;
81649b5e25fSSatish Balay      z[3*i+1] += v[3]*x1 + v[4]*x2 + v[7]*x3;
81749b5e25fSSatish Balay      z[3*i+2] += v[6]*x1 + v[7]*x2 + v[8]*x3;
818831a3094SHong Zhang      v += 9; jmin++;
8197fbae186SHong Zhang     }
820831a3094SHong Zhang     for (j=jmin; j<n; j++) {
82149b5e25fSSatish Balay       /* (strict lower triangular part of A)*x  */
82249b5e25fSSatish Balay       cval       = ib[j]*3;
82349b5e25fSSatish Balay       z[cval]     += v[0]*x1 + v[1]*x2 + v[2]*x3;
82449b5e25fSSatish Balay       z[cval+1]   += v[3]*x1 + v[4]*x2 + v[5]*x3;
82549b5e25fSSatish Balay       z[cval+2]   += v[6]*x1 + v[7]*x2 + v[8]*x3;
82649b5e25fSSatish Balay       /* (strict upper triangular part of A)*x  */
82749b5e25fSSatish Balay       z[3*i]   += v[0]*x[cval] + v[3]*x[cval+1]+ v[6]*x[cval+2];
82849b5e25fSSatish Balay       z[3*i+1] += v[1]*x[cval] + v[4]*x[cval+1]+ v[7]*x[cval+2];
82949b5e25fSSatish Balay       z[3*i+2] += v[2]*x[cval] + v[5]*x[cval+1]+ v[8]*x[cval+2];
83049b5e25fSSatish Balay       v  += 9;
83149b5e25fSSatish Balay     }
83249b5e25fSSatish Balay     xb +=3; ai++;
83349b5e25fSSatish Balay   }
83449b5e25fSSatish Balay 
8351ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
8361ebc52fbSHong Zhang   if (yy != xx) ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
8371ebc52fbSHong Zhang   if (zz != yy) ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
83849b5e25fSSatish Balay 
839efee365bSSatish Balay   ierr = PetscLogFlops(18*(a->nz*2 - A->m));CHKERRQ(ierr);
84049b5e25fSSatish Balay   PetscFunctionReturn(0);
84149b5e25fSSatish Balay }
84249b5e25fSSatish Balay 
8434a2ae208SSatish Balay #undef __FUNCT__
8444a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqSBAIJ_4"
845dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqSBAIJ_4(Mat A,Vec xx,Vec yy,Vec zz)
84649b5e25fSSatish Balay {
84749b5e25fSSatish Balay   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
84887828ca2SBarry Smith   PetscScalar    *x,*y,*z,*xb,x1,x2,x3,x4;
84949b5e25fSSatish Balay   MatScalar      *v;
8506849ba73SBarry Smith   PetscErrorCode ierr;
85113f74950SBarry Smith   PetscInt       mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,*ib,cval,j,jmin;
85249b5e25fSSatish Balay 
85349b5e25fSSatish Balay   PetscFunctionBegin;
8541ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
85549b5e25fSSatish Balay   if (yy != xx) {
8561ebc52fbSHong Zhang     ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
85749b5e25fSSatish Balay   } else {
85849b5e25fSSatish Balay     y = x;
85949b5e25fSSatish Balay   }
86049b5e25fSSatish Balay   if (zz != yy) {
861a9d4b620SHong Zhang     /* ierr = VecCopy(yy,zz);CHKERRQ(ierr); */
8621ebc52fbSHong Zhang     ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
863a9d4b620SHong Zhang     ierr = PetscMemcpy(z,y,yy->n*sizeof(PetscScalar));CHKERRQ(ierr);
86449b5e25fSSatish Balay   } else {
86549b5e25fSSatish Balay     z = y;
86649b5e25fSSatish Balay   }
86749b5e25fSSatish Balay 
86849b5e25fSSatish Balay   v     = a->a;
86949b5e25fSSatish Balay   xb = x;
87049b5e25fSSatish Balay 
87149b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
87249b5e25fSSatish Balay     n  = ai[1] - ai[0]; /* length of i_th block row of A */
87349b5e25fSSatish Balay     x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
87449b5e25fSSatish Balay     ib = aj + *ai;
875831a3094SHong Zhang     jmin = 0;
8767fbae186SHong Zhang     if (*ib == i){      /* (diag of A)*x */
87749b5e25fSSatish Balay       z[4*i]   += v[0]*x1 + v[4]*x2 +  v[8]*x3 + v[12]*x4;
87849b5e25fSSatish Balay       z[4*i+1] += v[4]*x1 + v[5]*x2 +  v[9]*x3 + v[13]*x4;
87949b5e25fSSatish Balay       z[4*i+2] += v[8]*x1 + v[9]*x2 + v[10]*x3 + v[14]*x4;
88049b5e25fSSatish Balay       z[4*i+3] += v[12]*x1+ v[13]*x2+ v[14]*x3 + v[15]*x4;
881831a3094SHong Zhang       v += 16; jmin++;
8827fbae186SHong Zhang     }
883831a3094SHong Zhang     for (j=jmin; j<n; j++) {
88449b5e25fSSatish Balay       /* (strict lower triangular part of A)*x  */
88549b5e25fSSatish Balay       cval       = ib[j]*4;
88649b5e25fSSatish Balay       z[cval]     += v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4;
88749b5e25fSSatish Balay       z[cval+1]   += v[4]*x1 + v[5]*x2 + v[6]*x3 + v[7]*x4;
88849b5e25fSSatish Balay       z[cval+2]   += v[8]*x1 + v[9]*x2 + v[10]*x3 + v[11]*x4;
88949b5e25fSSatish Balay       z[cval+3]   += v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4;
89049b5e25fSSatish Balay       /* (strict upper triangular part of A)*x  */
89149b5e25fSSatish Balay       z[4*i]   += v[0]*x[cval] + v[4]*x[cval+1]+ v[8]*x[cval+2] + v[12]*x[cval+3];
89249b5e25fSSatish Balay       z[4*i+1] += v[1]*x[cval] + v[5]*x[cval+1]+ v[9]*x[cval+2] + v[13]*x[cval+3];
89349b5e25fSSatish Balay       z[4*i+2] += v[2]*x[cval] + v[6]*x[cval+1]+ v[10]*x[cval+2]+ v[14]*x[cval+3];
89449b5e25fSSatish Balay       z[4*i+3] += v[3]*x[cval] + v[7]*x[cval+1]+ v[11]*x[cval+2]+ v[15]*x[cval+3];
89549b5e25fSSatish Balay       v  += 16;
89649b5e25fSSatish Balay     }
89749b5e25fSSatish Balay     xb +=4; ai++;
89849b5e25fSSatish Balay   }
89949b5e25fSSatish Balay 
9001ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
9011ebc52fbSHong Zhang   if (yy != xx) ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
9021ebc52fbSHong Zhang   if (zz != yy) ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
90349b5e25fSSatish Balay 
904efee365bSSatish Balay   ierr = PetscLogFlops(32*(a->nz*2 - A->m));CHKERRQ(ierr);
90549b5e25fSSatish Balay   PetscFunctionReturn(0);
90649b5e25fSSatish Balay }
90749b5e25fSSatish Balay 
9084a2ae208SSatish Balay #undef __FUNCT__
9094a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqSBAIJ_5"
910dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqSBAIJ_5(Mat A,Vec xx,Vec yy,Vec zz)
91149b5e25fSSatish Balay {
91249b5e25fSSatish Balay   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
91387828ca2SBarry Smith   PetscScalar    *x,*y,*z,*xb,x1,x2,x3,x4,x5;
91449b5e25fSSatish Balay   MatScalar      *v;
9156849ba73SBarry Smith   PetscErrorCode ierr;
91613f74950SBarry Smith   PetscInt       mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,*ib,cval,j,jmin;
91749b5e25fSSatish Balay 
91849b5e25fSSatish Balay   PetscFunctionBegin;
9191ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
92049b5e25fSSatish Balay   if (yy != xx) {
9211ebc52fbSHong Zhang     ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
92249b5e25fSSatish Balay   } else {
92349b5e25fSSatish Balay     y = x;
92449b5e25fSSatish Balay   }
92549b5e25fSSatish Balay   if (zz != yy) {
926a9d4b620SHong Zhang     /* ierr = VecCopy(yy,zz);CHKERRQ(ierr); */
9271ebc52fbSHong Zhang     ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
928a9d4b620SHong Zhang     ierr = PetscMemcpy(z,y,yy->n*sizeof(PetscScalar));CHKERRQ(ierr);
92949b5e25fSSatish Balay   } else {
93049b5e25fSSatish Balay     z = y;
93149b5e25fSSatish Balay   }
93249b5e25fSSatish Balay 
93349b5e25fSSatish Balay   v     = a->a;
93449b5e25fSSatish Balay   xb = x;
93549b5e25fSSatish Balay 
93649b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
93749b5e25fSSatish Balay     n  = ai[1] - ai[0]; /* length of i_th block row of A */
93849b5e25fSSatish Balay     x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5=xb[4];
93949b5e25fSSatish Balay     ib = aj + *ai;
940831a3094SHong Zhang     jmin = 0;
9417fbae186SHong Zhang     if (*ib == i){      /* (diag of A)*x */
94249b5e25fSSatish Balay       z[5*i]   += v[0]*x1  + v[5]*x2 + v[10]*x3 + v[15]*x4+ v[20]*x5;
94349b5e25fSSatish Balay       z[5*i+1] += v[5]*x1  + v[6]*x2 + v[11]*x3 + v[16]*x4+ v[21]*x5;
94449b5e25fSSatish Balay       z[5*i+2] += v[10]*x1 +v[11]*x2 + v[12]*x3 + v[17]*x4+ v[22]*x5;
94549b5e25fSSatish Balay       z[5*i+3] += v[15]*x1 +v[16]*x2 + v[17]*x3 + v[18]*x4+ v[23]*x5;
94649b5e25fSSatish Balay       z[5*i+4] += v[20]*x1 +v[21]*x2 + v[22]*x3 + v[23]*x4+ v[24]*x5;
947831a3094SHong Zhang       v += 25; jmin++;
9487fbae186SHong Zhang     }
949831a3094SHong Zhang     for (j=jmin; j<n; j++) {
95049b5e25fSSatish Balay       /* (strict lower triangular part of A)*x  */
95149b5e25fSSatish Balay       cval       = ib[j]*5;
95249b5e25fSSatish Balay       z[cval]     += v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4 + v[4]*x5;
95349b5e25fSSatish Balay       z[cval+1]   += v[5]*x1 + v[6]*x2 + v[7]*x3 + v[8]*x4 + v[9]*x5;
95449b5e25fSSatish Balay       z[cval+2]   += v[10]*x1 + v[11]*x2 + v[12]*x3 + v[13]*x4+ v[14]*x5;
95549b5e25fSSatish Balay       z[cval+3]   += v[15]*x1 + v[16]*x2 + v[17]*x3 + v[18]*x4+ v[19]*x5;
95649b5e25fSSatish Balay       z[cval+4]   += v[20]*x1 + v[21]*x2 + v[22]*x3 + v[23]*x4+ v[24]*x5;
95749b5e25fSSatish Balay       /* (strict upper triangular part of A)*x  */
95849b5e25fSSatish 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];
95949b5e25fSSatish 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];
96049b5e25fSSatish 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];
96149b5e25fSSatish 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];
96249b5e25fSSatish 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];
96349b5e25fSSatish Balay       v  += 25;
96449b5e25fSSatish Balay     }
96549b5e25fSSatish Balay     xb +=5; ai++;
96649b5e25fSSatish Balay   }
96749b5e25fSSatish Balay 
9681ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
9691ebc52fbSHong Zhang   if (yy != xx) ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
9701ebc52fbSHong Zhang   if (zz != yy) ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
97149b5e25fSSatish Balay 
972efee365bSSatish Balay   ierr = PetscLogFlops(50*(a->nz*2 - A->m));CHKERRQ(ierr);
97349b5e25fSSatish Balay   PetscFunctionReturn(0);
97449b5e25fSSatish Balay }
9754a2ae208SSatish Balay #undef __FUNCT__
9764a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqSBAIJ_6"
977dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqSBAIJ_6(Mat A,Vec xx,Vec yy,Vec zz)
97849b5e25fSSatish Balay {
97949b5e25fSSatish Balay   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
98087828ca2SBarry Smith   PetscScalar    *x,*y,*z,*xb,x1,x2,x3,x4,x5,x6;
98149b5e25fSSatish Balay   MatScalar      *v;
9826849ba73SBarry Smith   PetscErrorCode ierr;
98313f74950SBarry Smith   PetscInt       mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,*ib,cval,j,jmin;
98449b5e25fSSatish Balay 
98549b5e25fSSatish Balay   PetscFunctionBegin;
9861ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
98749b5e25fSSatish Balay   if (yy != xx) {
9881ebc52fbSHong Zhang     ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
98949b5e25fSSatish Balay   } else {
99049b5e25fSSatish Balay     y = x;
99149b5e25fSSatish Balay   }
99249b5e25fSSatish Balay   if (zz != yy) {
9931ebc52fbSHong Zhang     ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
994a9d4b620SHong Zhang     ierr = PetscMemcpy(z,y,yy->n*sizeof(PetscScalar));CHKERRQ(ierr);
99549b5e25fSSatish Balay   } else {
99649b5e25fSSatish Balay     z = y;
99749b5e25fSSatish Balay   }
99849b5e25fSSatish Balay 
99949b5e25fSSatish Balay   v     = a->a;
100049b5e25fSSatish Balay   xb = x;
100149b5e25fSSatish Balay 
100249b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
100349b5e25fSSatish Balay     n  = ai[1] - ai[0]; /* length of i_th block row of A */
100449b5e25fSSatish Balay     x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5=xb[4]; x6=xb[5];
100549b5e25fSSatish Balay     ib = aj + *ai;
1006831a3094SHong Zhang     jmin = 0;
10077fbae186SHong Zhang     if (*ib == i){     /* (diag of A)*x */
100849b5e25fSSatish Balay       z[6*i]   += v[0]*x1  + v[6]*x2 + v[12]*x3 + v[18]*x4+ v[24]*x5 + v[30]*x6;
100949b5e25fSSatish Balay       z[6*i+1] += v[6]*x1  + v[7]*x2 + v[13]*x3 + v[19]*x4+ v[25]*x5 + v[31]*x6;
101049b5e25fSSatish Balay       z[6*i+2] += v[12]*x1 +v[13]*x2 + v[14]*x3 + v[20]*x4+ v[26]*x5 + v[32]*x6;
101149b5e25fSSatish Balay       z[6*i+3] += v[18]*x1 +v[19]*x2 + v[20]*x3 + v[21]*x4+ v[27]*x5 + v[33]*x6;
101249b5e25fSSatish Balay       z[6*i+4] += v[24]*x1 +v[25]*x2 + v[26]*x3 + v[27]*x4+ v[28]*x5 + v[34]*x6;
101349b5e25fSSatish Balay       z[6*i+5] += v[30]*x1 +v[31]*x2 + v[32]*x3 + v[33]*x4+ v[34]*x5 + v[35]*x6;
1014831a3094SHong Zhang       v += 36; jmin++;
10157fbae186SHong Zhang     }
1016831a3094SHong Zhang     for (j=jmin; j<n; j++) {
101749b5e25fSSatish Balay       /* (strict lower triangular part of A)*x  */
101849b5e25fSSatish Balay       cval       = ib[j]*6;
101949b5e25fSSatish Balay       z[cval]   += v[0]*x1  + v[1]*x2 + v[2]*x3 + v[3]*x4+ v[4]*x5 + v[5]*x6;
102049b5e25fSSatish Balay       z[cval+1] += v[6]*x1  + v[7]*x2 + v[8]*x3 + v[9]*x4+ v[10]*x5 + v[11]*x6;
102149b5e25fSSatish Balay       z[cval+2] += v[12]*x1  + v[13]*x2 + v[14]*x3 + v[15]*x4+ v[16]*x5 + v[17]*x6;
102249b5e25fSSatish Balay       z[cval+3] += v[18]*x1  + v[19]*x2 + v[20]*x3 + v[21]*x4+ v[22]*x5 + v[23]*x6;
102349b5e25fSSatish Balay       z[cval+4] += v[24]*x1  + v[25]*x2 + v[26]*x3 + v[27]*x4+ v[28]*x5 + v[29]*x6;
102449b5e25fSSatish Balay       z[cval+5] += v[30]*x1  + v[31]*x2 + v[32]*x3 + v[33]*x4+ v[34]*x5 + v[35]*x6;
102549b5e25fSSatish Balay       /* (strict upper triangular part of A)*x  */
102649b5e25fSSatish 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];
102749b5e25fSSatish 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];
102849b5e25fSSatish 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];
102949b5e25fSSatish 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];
103049b5e25fSSatish 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];
103149b5e25fSSatish 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];
103249b5e25fSSatish Balay       v  += 36;
103349b5e25fSSatish Balay     }
103449b5e25fSSatish Balay     xb +=6; ai++;
103549b5e25fSSatish Balay   }
103649b5e25fSSatish Balay 
10371ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
10381ebc52fbSHong Zhang   if (yy != xx) ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
10391ebc52fbSHong Zhang   if (zz != yy) ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
104049b5e25fSSatish Balay 
1041efee365bSSatish Balay   ierr = PetscLogFlops(72*(a->nz*2 - A->m));CHKERRQ(ierr);
104249b5e25fSSatish Balay   PetscFunctionReturn(0);
104349b5e25fSSatish Balay }
104449b5e25fSSatish Balay 
10454a2ae208SSatish Balay #undef __FUNCT__
10464a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqSBAIJ_7"
1047dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqSBAIJ_7(Mat A,Vec xx,Vec yy,Vec zz)
104849b5e25fSSatish Balay {
104949b5e25fSSatish Balay   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
105087828ca2SBarry Smith   PetscScalar    *x,*y,*z,*xb,x1,x2,x3,x4,x5,x6,x7;
105149b5e25fSSatish Balay   MatScalar      *v;
10526849ba73SBarry Smith   PetscErrorCode ierr;
105313f74950SBarry Smith   PetscInt       mbs=a->mbs,i,*aj=a->j,*ai=a->i,n,*ib,cval,j,jmin;
105449b5e25fSSatish Balay 
105549b5e25fSSatish Balay   PetscFunctionBegin;
10561ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
105749b5e25fSSatish Balay   if (yy != xx) {
10581ebc52fbSHong Zhang     ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
105949b5e25fSSatish Balay   } else {
106049b5e25fSSatish Balay     y = x;
106149b5e25fSSatish Balay   }
106249b5e25fSSatish Balay   if (zz != yy) {
10631ebc52fbSHong Zhang     ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
1064a9d4b620SHong Zhang     ierr = PetscMemcpy(z,y,yy->n*sizeof(PetscScalar));CHKERRQ(ierr);
106549b5e25fSSatish Balay   } else {
106649b5e25fSSatish Balay     z = y;
106749b5e25fSSatish Balay   }
106849b5e25fSSatish Balay 
106949b5e25fSSatish Balay   v     = a->a;
107049b5e25fSSatish Balay   xb = x;
107149b5e25fSSatish Balay 
107249b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
107349b5e25fSSatish Balay     n  = ai[1] - ai[0]; /* length of i_th block row of A */
107449b5e25fSSatish Balay     x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5=xb[4]; x6=xb[5]; x7=xb[6];
107549b5e25fSSatish Balay     ib = aj + *ai;
1076831a3094SHong Zhang     jmin = 0;
10777fbae186SHong Zhang     if (*ib == i){     /* (diag of A)*x */
107849b5e25fSSatish 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;
107949b5e25fSSatish 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;
108049b5e25fSSatish 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;
108149b5e25fSSatish 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;
108249b5e25fSSatish 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;
108349b5e25fSSatish 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;
108449b5e25fSSatish 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;
1085831a3094SHong Zhang       v += 49; jmin++;
10867fbae186SHong Zhang     }
1087831a3094SHong Zhang     for (j=jmin; j<n; j++) {
108849b5e25fSSatish Balay       /* (strict lower triangular part of A)*x  */
108949b5e25fSSatish Balay       cval       = ib[j]*7;
109049b5e25fSSatish Balay       z[cval]   += v[0]*x1  + v[1]*x2 + v[2]*x3 + v[3]*x4+ v[4]*x5 + v[5]*x6+ v[6]*x7;
109149b5e25fSSatish 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;
109249b5e25fSSatish 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;
109349b5e25fSSatish 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;
109449b5e25fSSatish 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;
109549b5e25fSSatish 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;
109649b5e25fSSatish 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;
109749b5e25fSSatish Balay       /* (strict upper triangular part of A)*x  */
109849b5e25fSSatish 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];
109949b5e25fSSatish 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];
110049b5e25fSSatish 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];
110149b5e25fSSatish 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];
110249b5e25fSSatish 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];
110349b5e25fSSatish 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];
110449b5e25fSSatish 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];
110549b5e25fSSatish Balay       v  += 49;
110649b5e25fSSatish Balay     }
110749b5e25fSSatish Balay     xb +=7; ai++;
110849b5e25fSSatish Balay   }
110949b5e25fSSatish Balay 
11101ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
11111ebc52fbSHong Zhang   if (yy != xx) ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
11121ebc52fbSHong Zhang   if (zz != yy) ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
111349b5e25fSSatish Balay 
1114efee365bSSatish Balay   ierr = PetscLogFlops(98*(a->nz*2 - A->m));CHKERRQ(ierr);
111549b5e25fSSatish Balay   PetscFunctionReturn(0);
111649b5e25fSSatish Balay }
111749b5e25fSSatish Balay 
11184a2ae208SSatish Balay #undef __FUNCT__
11194a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqSBAIJ_N"
1120dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqSBAIJ_N(Mat A,Vec xx,Vec yy,Vec zz)
112149b5e25fSSatish Balay {
112249b5e25fSSatish Balay   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
112387828ca2SBarry Smith   PetscScalar    *x,*x_ptr,*y,*z,*z_ptr=0,*xb,*zb,*work,*workt;
1124066653e3SSatish Balay   MatScalar      *v;
1125dfbe8321SBarry Smith   PetscErrorCode ierr;
1126521d7252SBarry Smith   PetscInt       mbs=a->mbs,i,*idx,*aj,*ii,bs=A->bs,j,n,bs2=a->bs2,ncols,k;
112749b5e25fSSatish Balay 
112849b5e25fSSatish Balay   PetscFunctionBegin;
11291ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr); x_ptr=x;
113049b5e25fSSatish Balay   if (yy != xx) {
11311ebc52fbSHong Zhang     ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
113249b5e25fSSatish Balay   } else {
113349b5e25fSSatish Balay     y = x;
113449b5e25fSSatish Balay   }
113549b5e25fSSatish Balay   if (zz != yy) {
11361ebc52fbSHong Zhang     ierr = VecGetArray(zz,&z);CHKERRQ(ierr); z_ptr=z;
1137a9d4b620SHong Zhang     ierr = PetscMemcpy(z,y,yy->n*sizeof(PetscScalar));CHKERRQ(ierr);
113849b5e25fSSatish Balay   } else {
113949b5e25fSSatish Balay     z = y;
114049b5e25fSSatish Balay   }
114149b5e25fSSatish Balay 
114249b5e25fSSatish Balay   aj   = a->j;
114349b5e25fSSatish Balay   v    = a->a;
114449b5e25fSSatish Balay   ii   = a->i;
114549b5e25fSSatish Balay 
114649b5e25fSSatish Balay   if (!a->mult_work) {
114787828ca2SBarry Smith     ierr = PetscMalloc((A->m+1)*sizeof(PetscScalar),&a->mult_work);CHKERRQ(ierr);
114849b5e25fSSatish Balay   }
114949b5e25fSSatish Balay   work = a->mult_work;
115049b5e25fSSatish Balay 
115149b5e25fSSatish Balay 
115249b5e25fSSatish Balay   for (i=0; i<mbs; i++) {
115349b5e25fSSatish Balay     n     = ii[1] - ii[0]; ncols = n*bs;
115449b5e25fSSatish Balay     workt = work; idx=aj+ii[0];
115549b5e25fSSatish Balay 
115649b5e25fSSatish Balay     /* upper triangular part */
115749b5e25fSSatish Balay     for (j=0; j<n; j++) {
115849b5e25fSSatish Balay       xb = x_ptr + bs*(*idx++);
115949b5e25fSSatish Balay       for (k=0; k<bs; k++) workt[k] = xb[k];
116049b5e25fSSatish Balay       workt += bs;
116149b5e25fSSatish Balay     }
116249b5e25fSSatish Balay     /* z(i*bs:(i+1)*bs-1) += A(i,:)*x */
116349b5e25fSSatish Balay     Kernel_w_gets_w_plus_Ar_times_v(bs,ncols,work,v,z);
116449b5e25fSSatish Balay 
116549b5e25fSSatish Balay     /* strict lower triangular part */
1166831a3094SHong Zhang     idx = aj+ii[0];
1167831a3094SHong Zhang     if (*idx == i){
116896b9376eSHong Zhang       ncols -= bs; v += bs2; idx++; n--;
1169831a3094SHong Zhang     }
117049b5e25fSSatish Balay     if (ncols > 0){
117149b5e25fSSatish Balay       workt = work;
117287828ca2SBarry Smith       ierr  = PetscMemzero(workt,ncols*sizeof(PetscScalar));CHKERRQ(ierr);
1173831a3094SHong Zhang       Kernel_w_gets_w_plus_trans_Ar_times_v(bs,ncols,x,v,workt);
1174831a3094SHong Zhang       for (j=0; j<n; j++) {
1175831a3094SHong Zhang         zb = z_ptr + bs*(*idx++);
117649b5e25fSSatish Balay         for (k=0; k<bs; k++) zb[k] += workt[k] ;
117749b5e25fSSatish Balay         workt += bs;
117849b5e25fSSatish Balay       }
117949b5e25fSSatish Balay     }
118049b5e25fSSatish Balay 
118149b5e25fSSatish Balay     x += bs; v += n*bs2; z += bs; ii++;
118249b5e25fSSatish Balay   }
118349b5e25fSSatish Balay 
11841ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
11851ebc52fbSHong Zhang   if (yy != xx) ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
11861ebc52fbSHong Zhang   if (zz != yy) ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
118749b5e25fSSatish Balay 
1188efee365bSSatish Balay   ierr = PetscLogFlops(2*(a->nz*2 - A->m));CHKERRQ(ierr);
118949b5e25fSSatish Balay   PetscFunctionReturn(0);
119049b5e25fSSatish Balay }
119149b5e25fSSatish Balay 
11924a2ae208SSatish Balay #undef __FUNCT__
11934a2ae208SSatish Balay #define __FUNCT__ "MatScale_SeqSBAIJ"
1194dfbe8321SBarry Smith PetscErrorCode MatScale_SeqSBAIJ(const PetscScalar *alpha,Mat inA)
119549b5e25fSSatish Balay {
119649b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)inA->data;
11974ce68768SBarry Smith   PetscBLASInt one = 1,totalnz = (PetscBLASInt)a->bs2*a->nz;
1198efee365bSSatish Balay   PetscErrorCode ierr;
119949b5e25fSSatish Balay 
120049b5e25fSSatish Balay   PetscFunctionBegin;
120171044d3cSBarry Smith   BLASscal_(&totalnz,(PetscScalar*)alpha,a->a,&one);
1202efee365bSSatish Balay   ierr = PetscLogFlops(totalnz);CHKERRQ(ierr);
120349b5e25fSSatish Balay   PetscFunctionReturn(0);
120449b5e25fSSatish Balay }
120549b5e25fSSatish Balay 
12064a2ae208SSatish Balay #undef __FUNCT__
12074a2ae208SSatish Balay #define __FUNCT__ "MatNorm_SeqSBAIJ"
1208dfbe8321SBarry Smith PetscErrorCode MatNorm_SeqSBAIJ(Mat A,NormType type,PetscReal *norm)
120949b5e25fSSatish Balay {
121049b5e25fSSatish Balay   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
121149b5e25fSSatish Balay   MatScalar      *v = a->a;
121249b5e25fSSatish Balay   PetscReal      sum_diag = 0.0, sum_off = 0.0, *sum;
121313f74950SBarry Smith   PetscInt       i,j,k,bs = A->bs,bs2=a->bs2,k1,mbs=a->mbs,*aj=a->j;
12146849ba73SBarry Smith   PetscErrorCode ierr;
121513f74950SBarry Smith   PetscInt       *jl,*il,jmin,jmax,nexti,ik,*col;
121649b5e25fSSatish Balay 
121749b5e25fSSatish Balay   PetscFunctionBegin;
121849b5e25fSSatish Balay   if (type == NORM_FROBENIUS) {
121949b5e25fSSatish Balay     for (k=0; k<mbs; k++){
122049b5e25fSSatish Balay       jmin = a->i[k]; jmax = a->i[k+1];
1221831a3094SHong Zhang       col  = aj + jmin;
1222831a3094SHong Zhang       if (*col == k){         /* diagonal block */
122349b5e25fSSatish Balay         for (i=0; i<bs2; i++){
122449b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX)
122549b5e25fSSatish Balay           sum_diag += PetscRealPart(PetscConj(*v)*(*v)); v++;
122649b5e25fSSatish Balay #else
122749b5e25fSSatish Balay           sum_diag += (*v)*(*v); v++;
122849b5e25fSSatish Balay #endif
122949b5e25fSSatish Balay         }
1230831a3094SHong Zhang         jmin++;
1231831a3094SHong Zhang       }
1232831a3094SHong Zhang       for (j=jmin; j<jmax; j++){  /* off-diagonal blocks */
123349b5e25fSSatish Balay         for (i=0; i<bs2; i++){
123449b5e25fSSatish Balay #if defined(PETSC_USE_COMPLEX)
123549b5e25fSSatish Balay           sum_off += PetscRealPart(PetscConj(*v)*(*v)); v++;
123649b5e25fSSatish Balay #else
123749b5e25fSSatish Balay           sum_off += (*v)*(*v); v++;
123849b5e25fSSatish Balay #endif
123949b5e25fSSatish Balay         }
124049b5e25fSSatish Balay       }
124149b5e25fSSatish Balay     }
124249b5e25fSSatish Balay     *norm = sqrt(sum_diag + 2*sum_off);
12430b8dc8d2SHong Zhang   }  else if (type == NORM_INFINITY || type == NORM_1) { /* maximum row/column sum */
12440b8dc8d2SHong Zhang     ierr = PetscMalloc((2*mbs+1)*sizeof(PetscInt)+bs*sizeof(PetscReal),&il);CHKERRQ(ierr);
12450b8dc8d2SHong Zhang     jl   = il + mbs;
12460b8dc8d2SHong Zhang     sum  = (PetscReal*)(jl + mbs);
12470b8dc8d2SHong Zhang     for (i=0; i<mbs; i++) jl[i] = mbs;
12480b8dc8d2SHong Zhang     il[0] = 0;
124949b5e25fSSatish Balay 
125049b5e25fSSatish Balay     *norm = 0.0;
125149b5e25fSSatish Balay     for (k=0; k<mbs; k++) { /* k_th block row */
125249b5e25fSSatish Balay       for (j=0; j<bs; j++) sum[j]=0.0;
125349b5e25fSSatish Balay       /*-- col sum --*/
125449b5e25fSSatish Balay       i = jl[k]; /* first |A(i,k)| to be added */
1255831a3094SHong Zhang       /* jl[k]=i: first nozero element in row i for submatrix A(1:k,k:n) (active window)
1256831a3094SHong Zhang                   at step k */
125749b5e25fSSatish Balay       while (i<mbs){
125849b5e25fSSatish Balay         nexti = jl[i];  /* next block row to be added */
125949b5e25fSSatish Balay         ik    = il[i];  /* block index of A(i,k) in the array a */
126049b5e25fSSatish Balay         for (j=0; j<bs; j++){
126149b5e25fSSatish Balay           v = a->a + ik*bs2 + j*bs;
126249b5e25fSSatish Balay           for (k1=0; k1<bs; k1++) {
126349b5e25fSSatish Balay             sum[j] += PetscAbsScalar(*v); v++;
126449b5e25fSSatish Balay           }
126549b5e25fSSatish Balay         }
126649b5e25fSSatish Balay         /* update il, jl */
1267831a3094SHong Zhang         jmin = ik + 1; /* block index of array a: points to the next nonzero of A in row i */
1268831a3094SHong Zhang         jmax = a->i[i+1];
126949b5e25fSSatish Balay         if (jmin < jmax){
127049b5e25fSSatish Balay           il[i] = jmin;
127149b5e25fSSatish Balay           j   = a->j[jmin];
127249b5e25fSSatish Balay           jl[i] = jl[j]; jl[j]=i;
127349b5e25fSSatish Balay         }
127449b5e25fSSatish Balay         i = nexti;
127549b5e25fSSatish Balay       }
127649b5e25fSSatish Balay       /*-- row sum --*/
127749b5e25fSSatish Balay       jmin = a->i[k]; jmax = a->i[k+1];
127849b5e25fSSatish Balay       for (i=jmin; i<jmax; i++) {
127949b5e25fSSatish Balay         for (j=0; j<bs; j++){
128049b5e25fSSatish Balay           v = a->a + i*bs2 + j;
128149b5e25fSSatish Balay           for (k1=0; k1<bs; k1++){
12820b8dc8d2SHong Zhang             sum[j] += PetscAbsScalar(*v); v += bs;
128349b5e25fSSatish Balay           }
128449b5e25fSSatish Balay         }
128549b5e25fSSatish Balay       }
128649b5e25fSSatish Balay       /* add k_th block row to il, jl */
1287831a3094SHong Zhang       col = aj+jmin;
1288831a3094SHong Zhang       if (*col == k) jmin++;
128949b5e25fSSatish Balay       if (jmin < jmax){
129049b5e25fSSatish Balay         il[k] = jmin;
12910b8dc8d2SHong Zhang         j = a->j[jmin]; jl[k] = jl[j]; jl[j] = k;
129249b5e25fSSatish Balay       }
129349b5e25fSSatish Balay       for (j=0; j<bs; j++){
129449b5e25fSSatish Balay         if (sum[j] > *norm) *norm = sum[j];
129549b5e25fSSatish Balay       }
129649b5e25fSSatish Balay     }
129749b5e25fSSatish Balay     ierr = PetscFree(il);CHKERRQ(ierr);
129849b5e25fSSatish Balay   } else {
1299347d480fSBarry Smith     SETERRQ(PETSC_ERR_SUP,"No support for this norm yet");
130049b5e25fSSatish Balay   }
130149b5e25fSSatish Balay   PetscFunctionReturn(0);
130249b5e25fSSatish Balay }
130349b5e25fSSatish Balay 
13044a2ae208SSatish Balay #undef __FUNCT__
13054a2ae208SSatish Balay #define __FUNCT__ "MatEqual_SeqSBAIJ"
1306dfbe8321SBarry Smith PetscErrorCode MatEqual_SeqSBAIJ(Mat A,Mat B,PetscTruth* flg)
130749b5e25fSSatish Balay {
130849b5e25fSSatish Balay   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ *)A->data,*b = (Mat_SeqSBAIJ *)B->data;
1309dfbe8321SBarry Smith   PetscErrorCode ierr;
131049b5e25fSSatish Balay 
131149b5e25fSSatish Balay   PetscFunctionBegin;
131249b5e25fSSatish Balay 
131349b5e25fSSatish Balay   /* If the  matrix/block dimensions are not equal, or no of nonzeros or shift */
1314521d7252SBarry Smith   if ((A->m != B->m) || (A->n != B->n) || (A->bs != B->bs)|| (a->nz != b->nz)) {
1315ef511fbeSHong Zhang     *flg = PETSC_FALSE;
1316ef511fbeSHong Zhang     PetscFunctionReturn(0);
131749b5e25fSSatish Balay   }
131849b5e25fSSatish Balay 
131949b5e25fSSatish Balay   /* if the a->i are the same */
132013f74950SBarry Smith   ierr = PetscMemcmp(a->i,b->i,(a->mbs+1)*sizeof(PetscInt),flg);CHKERRQ(ierr);
1321abc0a331SBarry Smith   if (!*flg) {
132249b5e25fSSatish Balay     PetscFunctionReturn(0);
132349b5e25fSSatish Balay   }
132449b5e25fSSatish Balay 
132549b5e25fSSatish Balay   /* if a->j are the same */
132613f74950SBarry Smith   ierr = PetscMemcmp(a->j,b->j,(a->nz)*sizeof(PetscInt),flg);CHKERRQ(ierr);
1327abc0a331SBarry Smith   if (!*flg) {
132849b5e25fSSatish Balay     PetscFunctionReturn(0);
132949b5e25fSSatish Balay   }
133049b5e25fSSatish Balay   /* if a->a are the same */
1331521d7252SBarry Smith   ierr = PetscMemcmp(a->a,b->a,(a->nz)*(A->bs)*(A->bs)*sizeof(PetscScalar),flg);CHKERRQ(ierr);
1332935af2e7SHong Zhang   PetscFunctionReturn(0);
133349b5e25fSSatish Balay }
133449b5e25fSSatish Balay 
13354a2ae208SSatish Balay #undef __FUNCT__
13364a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_SeqSBAIJ"
1337dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_SeqSBAIJ(Mat A,Vec v)
133849b5e25fSSatish Balay {
133949b5e25fSSatish Balay   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
1340dfbe8321SBarry Smith   PetscErrorCode ierr;
134113f74950SBarry Smith   PetscInt       i,j,k,n,row,bs,*ai,*aj,ambs,bs2;
134287828ca2SBarry Smith   PetscScalar    *x,zero = 0.0;
134349b5e25fSSatish Balay   MatScalar      *aa,*aa_j;
134449b5e25fSSatish Balay 
134549b5e25fSSatish Balay   PetscFunctionBegin;
1346521d7252SBarry Smith   bs   = A->bs;
134782799104SHong Zhang   if (A->factor && bs>1) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix with bs>1");
134882799104SHong Zhang 
134949b5e25fSSatish Balay   aa   = a->a;
135049b5e25fSSatish Balay   ai   = a->i;
135149b5e25fSSatish Balay   aj   = a->j;
135249b5e25fSSatish Balay   ambs = a->mbs;
135349b5e25fSSatish Balay   bs2  = a->bs2;
135449b5e25fSSatish Balay 
13552dcb1b2aSMatthew Knepley   ierr = VecSet(v,zero);CHKERRQ(ierr);
13561ebc52fbSHong Zhang   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
135749b5e25fSSatish Balay   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
1358ef511fbeSHong Zhang   if (n != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
135949b5e25fSSatish Balay   for (i=0; i<ambs; i++) {
136049b5e25fSSatish Balay     j=ai[i];
136149b5e25fSSatish Balay     if (aj[j] == i) {             /* if this is a diagonal element */
136249b5e25fSSatish Balay       row  = i*bs;
136349b5e25fSSatish Balay       aa_j = aa + j*bs2;
136482799104SHong Zhang       if (A->factor && bs==1){
136582799104SHong Zhang         for (k=0; k<bs2; k+=(bs+1),row++) x[row] = 1.0/aa_j[k];
136682799104SHong Zhang       } else {
136749b5e25fSSatish Balay         for (k=0; k<bs2; k+=(bs+1),row++) x[row] = aa_j[k];
136849b5e25fSSatish Balay       }
136949b5e25fSSatish Balay     }
137082799104SHong Zhang   }
137182799104SHong Zhang 
13721ebc52fbSHong Zhang   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
137349b5e25fSSatish Balay   PetscFunctionReturn(0);
137449b5e25fSSatish Balay }
137549b5e25fSSatish Balay 
13764a2ae208SSatish Balay #undef __FUNCT__
13774a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_SeqSBAIJ"
1378dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_SeqSBAIJ(Mat A,Vec ll,Vec rr)
137949b5e25fSSatish Balay {
138049b5e25fSSatish Balay   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
13815e90f9d9SHong Zhang   PetscScalar    *l,x,*li,*ri;
138249b5e25fSSatish Balay   MatScalar      *aa,*v;
1383dfbe8321SBarry Smith   PetscErrorCode ierr;
13845e90f9d9SHong Zhang   PetscInt       i,j,k,lm,M,m,*ai,*aj,mbs,tmp,bs,bs2;
1385b3bf805bSHong Zhang   PetscTruth     flg;
138649b5e25fSSatish Balay 
138749b5e25fSSatish Balay   PetscFunctionBegin;
1388b3bf805bSHong Zhang   if (ll != rr){
1389b3bf805bSHong Zhang     ierr = VecEqual(ll,rr,&flg);CHKERRQ(ierr);
1390b3bf805bSHong Zhang     if (!flg)
1391b3bf805bSHong Zhang       SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"For symmetric format, left and right scaling vectors must be same\n");
1392b3bf805bSHong Zhang   }
1393b3bf805bSHong Zhang   if (!ll) PetscFunctionReturn(0);
139449b5e25fSSatish Balay   ai  = a->i;
139549b5e25fSSatish Balay   aj  = a->j;
139649b5e25fSSatish Balay   aa  = a->a;
1397ef511fbeSHong Zhang   m   = A->m;
1398521d7252SBarry Smith   bs  = A->bs;
139949b5e25fSSatish Balay   mbs = a->mbs;
140049b5e25fSSatish Balay   bs2 = a->bs2;
140149b5e25fSSatish Balay 
14021ebc52fbSHong Zhang   ierr = VecGetArray(ll,&l);CHKERRQ(ierr);
140349b5e25fSSatish Balay   ierr = VecGetLocalSize(ll,&lm);CHKERRQ(ierr);
1404347d480fSBarry Smith   if (lm != m) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vector wrong length");
140549b5e25fSSatish Balay   for (i=0; i<mbs; i++) { /* for each block row */
140649b5e25fSSatish Balay     M  = ai[i+1] - ai[i];
140749b5e25fSSatish Balay     li = l + i*bs;
140849b5e25fSSatish Balay     v  = aa + bs2*ai[i];
140949b5e25fSSatish Balay     for (j=0; j<M; j++) { /* for each block */
141049b5e25fSSatish Balay       ri = l + bs*aj[ai[i]+j];
14115e90f9d9SHong Zhang       for (k=0; k<bs; k++) {
141249b5e25fSSatish Balay         x = ri[k];
141349b5e25fSSatish Balay         for (tmp=0; tmp<bs; tmp++) (*v++) *= li[tmp]*x;
141449b5e25fSSatish Balay       }
141549b5e25fSSatish Balay     }
141649b5e25fSSatish Balay   }
14171ebc52fbSHong Zhang   ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr);
1418efee365bSSatish Balay   ierr = PetscLogFlops(2*a->nz);CHKERRQ(ierr);
141949b5e25fSSatish Balay   PetscFunctionReturn(0);
142049b5e25fSSatish Balay }
142149b5e25fSSatish Balay 
14224a2ae208SSatish Balay #undef __FUNCT__
14234a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_SeqSBAIJ"
1424dfbe8321SBarry Smith PetscErrorCode MatGetInfo_SeqSBAIJ(Mat A,MatInfoType flag,MatInfo *info)
142549b5e25fSSatish Balay {
142649b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
142749b5e25fSSatish Balay 
142849b5e25fSSatish Balay   PetscFunctionBegin;
1429ef511fbeSHong Zhang   info->rows_global    = (double)A->m;
1430ef511fbeSHong Zhang   info->columns_global = (double)A->m;
1431ef511fbeSHong Zhang   info->rows_local     = (double)A->m;
1432ef511fbeSHong Zhang   info->columns_local  = (double)A->m;
143349b5e25fSSatish Balay   info->block_size     = a->bs2;
14346c6c5352SBarry Smith   info->nz_allocated   = a->maxnz; /*num. of nonzeros in upper triangular part */
14356c6c5352SBarry Smith   info->nz_used        = a->bs2*a->nz; /*num. of nonzeros in upper triangular part */
143649b5e25fSSatish Balay   info->nz_unneeded    = (double)(info->nz_allocated - info->nz_used);
143749b5e25fSSatish Balay   info->assemblies   = A->num_ass;
143849b5e25fSSatish Balay   info->mallocs      = a->reallocs;
143949b5e25fSSatish Balay   info->memory       = A->mem;
144049b5e25fSSatish Balay   if (A->factor) {
144149b5e25fSSatish Balay     info->fill_ratio_given  = A->info.fill_ratio_given;
144249b5e25fSSatish Balay     info->fill_ratio_needed = A->info.fill_ratio_needed;
144349b5e25fSSatish Balay     info->factor_mallocs    = A->info.factor_mallocs;
144449b5e25fSSatish Balay   } else {
144549b5e25fSSatish Balay     info->fill_ratio_given  = 0;
144649b5e25fSSatish Balay     info->fill_ratio_needed = 0;
144749b5e25fSSatish Balay     info->factor_mallocs    = 0;
144849b5e25fSSatish Balay   }
144949b5e25fSSatish Balay   PetscFunctionReturn(0);
145049b5e25fSSatish Balay }
145149b5e25fSSatish Balay 
145249b5e25fSSatish Balay 
14534a2ae208SSatish Balay #undef __FUNCT__
14544a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_SeqSBAIJ"
1455dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_SeqSBAIJ(Mat A)
145649b5e25fSSatish Balay {
145749b5e25fSSatish Balay   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
1458dfbe8321SBarry Smith   PetscErrorCode ierr;
145949b5e25fSSatish Balay 
146049b5e25fSSatish Balay   PetscFunctionBegin;
146149b5e25fSSatish Balay   ierr = PetscMemzero(a->a,a->bs2*a->i[a->mbs]*sizeof(MatScalar));CHKERRQ(ierr);
146249b5e25fSSatish Balay   PetscFunctionReturn(0);
146349b5e25fSSatish Balay }
1464dc354874SHong Zhang 
14654a2ae208SSatish Balay #undef __FUNCT__
14664a2ae208SSatish Balay #define __FUNCT__ "MatGetRowMax_SeqSBAIJ"
1467dfbe8321SBarry Smith PetscErrorCode MatGetRowMax_SeqSBAIJ(Mat A,Vec v)
1468dc354874SHong Zhang {
1469dc354874SHong Zhang   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
1470dfbe8321SBarry Smith   PetscErrorCode ierr;
147113f74950SBarry Smith   PetscInt i,j,n,row,col,bs,*ai,*aj,mbs;
1472c3fca9a7SHong Zhang   PetscReal    atmp;
1473273d9f13SBarry Smith   MatScalar    *aa;
147487828ca2SBarry Smith   PetscScalar  zero = 0.0,*x;
147513f74950SBarry Smith   PetscInt          ncols,brow,bcol,krow,kcol;
1476dc354874SHong Zhang 
1477dc354874SHong Zhang   PetscFunctionBegin;
1478dc354874SHong Zhang   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1479521d7252SBarry Smith   bs   = A->bs;
1480dc354874SHong Zhang   aa   = a->a;
1481dc354874SHong Zhang   ai   = a->i;
1482dc354874SHong Zhang   aj   = a->j;
148344117c81SHong Zhang   mbs = a->mbs;
1484dc354874SHong Zhang 
14852dcb1b2aSMatthew Knepley   ierr = VecSet(v,zero);CHKERRQ(ierr);
14861ebc52fbSHong Zhang   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1487dc354874SHong Zhang   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
1488ef511fbeSHong Zhang   if (n != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
148944117c81SHong Zhang   for (i=0; i<mbs; i++) {
1490d0f6400bSHong Zhang     ncols = ai[1] - ai[0]; ai++;
1491d0f6400bSHong Zhang     brow  = bs*i;
149244117c81SHong Zhang     for (j=0; j<ncols; j++){
1493d0f6400bSHong Zhang       bcol = bs*(*aj);
149444117c81SHong Zhang       for (kcol=0; kcol<bs; kcol++){
1495d0f6400bSHong Zhang         col = bcol + kcol;      /* col index */
149644117c81SHong Zhang         for (krow=0; krow<bs; krow++){
1497d0f6400bSHong Zhang           atmp = PetscAbsScalar(*aa); aa++;
1498d0f6400bSHong Zhang           row = brow + krow;    /* row index */
1499c3fca9a7SHong Zhang           if (PetscRealPart(x[row]) < atmp) x[row] = atmp;
1500c3fca9a7SHong Zhang           if (*aj > i && PetscRealPart(x[col]) < atmp) x[col] = atmp;
150144117c81SHong Zhang         }
150244117c81SHong Zhang       }
1503d0f6400bSHong Zhang       aj++;
1504dc354874SHong Zhang     }
1505dc354874SHong Zhang   }
15061ebc52fbSHong Zhang   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1507dc354874SHong Zhang   PetscFunctionReturn(0);
1508dc354874SHong Zhang }
1509