xref: /petsc/src/mat/impls/baij/seq/baij2.c (revision e32f2f54e699d0aa6e733466c00da7e34666fe5e)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
2cac129eeSSatish Balay 
37c4f633dSBarry Smith #include "../src/mat/impls/baij/seq/baij.h"
4c60f0209SBarry Smith #include "../src/mat/blockinvert.h"
50a835dfdSSatish Balay #include "petscbt.h"
6f3da1532SBarry Smith #include "petscblaslapack.h"
7cac129eeSSatish Balay 
84a2ae208SSatish Balay #undef __FUNCT__
94a2ae208SSatish Balay #define __FUNCT__ "MatIncreaseOverlap_SeqBAIJ"
10690b6cddSBarry Smith PetscErrorCode MatIncreaseOverlap_SeqBAIJ(Mat A,PetscInt is_max,IS is[],PetscInt ov)
11a3192f15SSatish Balay {
12a3192f15SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
136849ba73SBarry Smith   PetscErrorCode ierr;
145d0c19d7SBarry Smith   PetscInt       row,i,j,k,l,m,n,*nidx,isz,val,ival;
155d0c19d7SBarry Smith   const PetscInt *idx;
16690b6cddSBarry Smith   PetscInt       start,end,*ai,*aj,bs,*nidx2;
17f1af5d2fSBarry Smith   PetscBT        table;
18a3192f15SSatish Balay 
193a40ed3dSBarry Smith   PetscFunctionBegin;
20a3192f15SSatish Balay   m     = a->mbs;
21a3192f15SSatish Balay   ai    = a->i;
22a3192f15SSatish Balay   aj    = a->j;
23d0f46423SBarry Smith   bs    = A->rmap->bs;
24a3192f15SSatish Balay 
25*e32f2f54SBarry Smith   if (ov < 0)  SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap specified");
26a3192f15SSatish Balay 
276831982aSBarry Smith   ierr = PetscBTCreate(m,table);CHKERRQ(ierr);
28690b6cddSBarry Smith   ierr = PetscMalloc((m+1)*sizeof(PetscInt),&nidx);CHKERRQ(ierr);
29d0f46423SBarry Smith   ierr = PetscMalloc((A->rmap->N+1)*sizeof(PetscInt),&nidx2);CHKERRQ(ierr);
30a3192f15SSatish Balay 
31a3192f15SSatish Balay   for (i=0; i<is_max; i++) {
32a3192f15SSatish Balay     /* Initialise the two local arrays */
33a3192f15SSatish Balay     isz  = 0;
346831982aSBarry Smith     ierr = PetscBTMemzero(m,table);CHKERRQ(ierr);
35a3192f15SSatish Balay 
36a3192f15SSatish Balay     /* Extract the indices, assume there can be duplicate entries */
37a3192f15SSatish Balay     ierr = ISGetIndices(is[i],&idx);CHKERRQ(ierr);
38b9b97703SBarry Smith     ierr = ISGetLocalSize(is[i],&n);CHKERRQ(ierr);
39a3192f15SSatish Balay 
40a3192f15SSatish Balay     /* Enter these into the temp arrays i.e mark table[row], enter row into new index */
41a3192f15SSatish Balay     for (j=0; j<n ; ++j){
42218c64b6SSatish Balay       ival = idx[j]/bs; /* convert the indices into block indices */
43*e32f2f54SBarry Smith       if (ival>=m) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"index greater than mat-dim");
44f1af5d2fSBarry Smith       if(!PetscBTLookupSet(table,ival)) { nidx[isz++] = ival;}
45a3192f15SSatish Balay     }
46a3192f15SSatish Balay     ierr = ISRestoreIndices(is[i],&idx);CHKERRQ(ierr);
47a3192f15SSatish Balay     ierr = ISDestroy(is[i]);CHKERRQ(ierr);
48a3192f15SSatish Balay 
49a3192f15SSatish Balay     k = 0;
50a3192f15SSatish Balay     for (j=0; j<ov; j++){ /* for each overlap*/
51a3192f15SSatish Balay       n = isz;
52a3192f15SSatish Balay       for (; k<n ; k++){ /* do only those rows in nidx[k], which are not done yet */
53a3192f15SSatish Balay         row   = nidx[k];
54a3192f15SSatish Balay         start = ai[row];
55a3192f15SSatish Balay         end   = ai[row+1];
56a3192f15SSatish Balay         for (l = start; l<end ; l++){
57a3192f15SSatish Balay           val = aj[l];
58f1af5d2fSBarry Smith           if (!PetscBTLookupSet(table,val)) {nidx[isz++] = val;}
59a3192f15SSatish Balay         }
60a3192f15SSatish Balay       }
61a3192f15SSatish Balay     }
62218c64b6SSatish Balay     /* expand the Index Set */
63218c64b6SSatish Balay     for (j=0; j<isz; j++) {
64218c64b6SSatish Balay       for (k=0; k<bs; k++)
65218c64b6SSatish Balay         nidx2[j*bs+k] = nidx[j]*bs+k;
66218c64b6SSatish Balay     }
67329f5518SBarry Smith     ierr = ISCreateGeneral(PETSC_COMM_SELF,isz*bs,nidx2,is+i);CHKERRQ(ierr);
68a3192f15SSatish Balay   }
696831982aSBarry Smith   ierr = PetscBTDestroy(table);CHKERRQ(ierr);
70606d414cSSatish Balay   ierr = PetscFree(nidx);CHKERRQ(ierr);
71606d414cSSatish Balay   ierr = PetscFree(nidx2);CHKERRQ(ierr);
723a40ed3dSBarry Smith   PetscFunctionReturn(0);
73a3192f15SSatish Balay }
741c351548SSatish Balay 
754a2ae208SSatish Balay #undef __FUNCT__
764a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqBAIJ_Private"
774aa3045dSJed Brown PetscErrorCode MatGetSubMatrix_SeqBAIJ_Private(Mat A,IS isrow,IS iscol,MatReuse scall,Mat *B)
78736121d4SSatish Balay {
79736121d4SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*c;
806849ba73SBarry Smith   PetscErrorCode ierr;
81690b6cddSBarry Smith   PetscInt       *smap,i,k,kstart,kend,oldcols = a->nbs,*lens;
82690b6cddSBarry Smith   PetscInt       row,mat_i,*mat_j,tcol,*mat_ilen;
835d0c19d7SBarry Smith   const PetscInt *irow,*icol;
845d0c19d7SBarry Smith   PetscInt       nrows,ncols,*ssmap,bs=A->rmap->bs,bs2=a->bs2;
85690b6cddSBarry Smith   PetscInt       *aj = a->j,*ai = a->i;
863f1db9ecSBarry Smith   MatScalar      *mat_a;
87736121d4SSatish Balay   Mat            C;
8814ca34e6SBarry Smith   PetscTruth     flag,sorted;
89736121d4SSatish Balay 
903a40ed3dSBarry Smith   PetscFunctionBegin;
9114ca34e6SBarry Smith   ierr = ISSorted(iscol,&sorted);CHKERRQ(ierr);
92*e32f2f54SBarry Smith   if (!sorted) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"IS is not sorted");
93736121d4SSatish Balay 
94736121d4SSatish Balay   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
95218c64b6SSatish Balay   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
96b9b97703SBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
97b9b97703SBarry Smith   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
98736121d4SSatish Balay 
99690b6cddSBarry Smith   ierr = PetscMalloc((1+oldcols)*sizeof(PetscInt),&smap);CHKERRQ(ierr);
100736121d4SSatish Balay   ssmap = smap;
101690b6cddSBarry Smith   ierr = PetscMalloc((1+nrows)*sizeof(PetscInt),&lens);CHKERRQ(ierr);
102690b6cddSBarry Smith   ierr  = PetscMemzero(smap,oldcols*sizeof(PetscInt));CHKERRQ(ierr);
103736121d4SSatish Balay   for (i=0; i<ncols; i++) smap[icol[i]] = i+1;
104736121d4SSatish Balay   /* determine lens of each row */
105736121d4SSatish Balay   for (i=0; i<nrows; i++) {
106736121d4SSatish Balay     kstart  = ai[irow[i]];
107736121d4SSatish Balay     kend    = kstart + a->ilen[irow[i]];
108736121d4SSatish Balay     lens[i] = 0;
109736121d4SSatish Balay       for (k=kstart; k<kend; k++) {
110736121d4SSatish Balay         if (ssmap[aj[k]]) {
111736121d4SSatish Balay           lens[i]++;
112736121d4SSatish Balay         }
113736121d4SSatish Balay       }
114736121d4SSatish Balay     }
115736121d4SSatish Balay   /* Create and fill new matrix */
116736121d4SSatish Balay   if (scall == MAT_REUSE_MATRIX) {
117736121d4SSatish Balay     c = (Mat_SeqBAIJ *)((*B)->data);
118736121d4SSatish Balay 
119*e32f2f54SBarry Smith     if (c->mbs!=nrows || c->nbs!=ncols || (*B)->rmap->bs!=bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Submatrix wrong size");
120690b6cddSBarry Smith     ierr = PetscMemcmp(c->ilen,lens,c->mbs *sizeof(PetscInt),&flag);CHKERRQ(ierr);
121abc0a331SBarry Smith     if (!flag) {
122*e32f2f54SBarry Smith       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong no of nonzeros");
123736121d4SSatish Balay     }
124690b6cddSBarry Smith     ierr = PetscMemzero(c->ilen,c->mbs*sizeof(PetscInt));CHKERRQ(ierr);
125736121d4SSatish Balay     C = *B;
1263a40ed3dSBarry Smith   } else {
1277adad957SLisandro Dalcin     ierr = MatCreate(((PetscObject)A)->comm,&C);CHKERRQ(ierr);
128f69a0ea3SMatthew Knepley     ierr = MatSetSizes(C,nrows*bs,ncols*bs,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1297adad957SLisandro Dalcin     ierr = MatSetType(C,((PetscObject)A)->type_name);CHKERRQ(ierr);
130ab93d7beSBarry Smith     ierr = MatSeqBAIJSetPreallocation_SeqBAIJ(C,bs,0,lens);CHKERRQ(ierr);
131736121d4SSatish Balay   }
132736121d4SSatish Balay   c = (Mat_SeqBAIJ *)(C->data);
133736121d4SSatish Balay   for (i=0; i<nrows; i++) {
134736121d4SSatish Balay     row    = irow[i];
135736121d4SSatish Balay     kstart = ai[row];
136736121d4SSatish Balay     kend   = kstart + a->ilen[row];
137736121d4SSatish Balay     mat_i  = c->i[i];
138736121d4SSatish Balay     mat_j  = c->j + mat_i;
139218c64b6SSatish Balay     mat_a  = c->a + mat_i*bs2;
140736121d4SSatish Balay     mat_ilen = c->ilen + i;
141736121d4SSatish Balay     for (k=kstart; k<kend; k++) {
142736121d4SSatish Balay       if ((tcol=ssmap[a->j[k]])) {
143736121d4SSatish Balay         *mat_j++ = tcol - 1;
144549d3d68SSatish Balay         ierr     = PetscMemcpy(mat_a,a->a+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
145549d3d68SSatish Balay         mat_a   += bs2;
146736121d4SSatish Balay         (*mat_ilen)++;
147736121d4SSatish Balay       }
148736121d4SSatish Balay     }
149736121d4SSatish Balay   }
150218c64b6SSatish Balay 
151736121d4SSatish Balay   /* Free work space */
152736121d4SSatish Balay   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
153606d414cSSatish Balay   ierr = PetscFree(smap);CHKERRQ(ierr);
154606d414cSSatish Balay   ierr = PetscFree(lens);CHKERRQ(ierr);
1556d4a8577SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1566d4a8577SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
157736121d4SSatish Balay 
158736121d4SSatish Balay   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
159736121d4SSatish Balay   *B = C;
1603a40ed3dSBarry Smith   PetscFunctionReturn(0);
161736121d4SSatish Balay }
162736121d4SSatish Balay 
1634a2ae208SSatish Balay #undef __FUNCT__
1644a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqBAIJ"
1654aa3045dSJed Brown PetscErrorCode MatGetSubMatrix_SeqBAIJ(Mat A,IS isrow,IS iscol,MatReuse scall,Mat *B)
166218c64b6SSatish Balay {
167218c64b6SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
168218c64b6SSatish Balay   IS             is1,is2;
1696849ba73SBarry Smith   PetscErrorCode ierr;
1705d0c19d7SBarry Smith   PetscInt       *vary,*iary,nrows,ncols,i,bs=A->rmap->bs,count;
1715d0c19d7SBarry Smith   const PetscInt *irow,*icol;
172218c64b6SSatish Balay 
1733a40ed3dSBarry Smith   PetscFunctionBegin;
174218c64b6SSatish Balay   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
175218c64b6SSatish Balay   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
176b9b97703SBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
177b9b97703SBarry Smith   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
178218c64b6SSatish Balay 
179218c64b6SSatish Balay   /* Verify if the indices corespond to each element in a block
180218c64b6SSatish Balay    and form the IS with compressed IS */
181fca92195SBarry Smith   ierr = PetscMalloc2(a->mbs,PetscInt,&vary,a->mbs,PetscInt,&iary);CHKERRQ(ierr);
182fca92195SBarry Smith   ierr = PetscMemzero(vary,a->mbs*sizeof(PetscInt));CHKERRQ(ierr);
183218c64b6SSatish Balay   for (i=0; i<nrows; i++) vary[irow[i]/bs]++;
184218c64b6SSatish Balay   count = 0;
185218c64b6SSatish Balay   for (i=0; i<a->mbs; i++) {
186*e32f2f54SBarry Smith     if (vary[i]!=0 && vary[i]!=bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Index set does not match blocks");
187218c64b6SSatish Balay     if (vary[i]==bs) iary[count++] = i;
188218c64b6SSatish Balay   }
189029af93fSBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,count,iary,&is1);CHKERRQ(ierr);
190218c64b6SSatish Balay 
191690b6cddSBarry Smith   ierr = PetscMemzero(vary,(a->mbs)*sizeof(PetscInt));CHKERRQ(ierr);
192218c64b6SSatish Balay   for (i=0; i<ncols; i++) vary[icol[i]/bs]++;
193218c64b6SSatish Balay   count = 0;
194218c64b6SSatish Balay   for (i=0; i<a->mbs; i++) {
195*e32f2f54SBarry Smith     if (vary[i]!=0 && vary[i]!=bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Internal error in PETSc");
196218c64b6SSatish Balay     if (vary[i]==bs) iary[count++] = i;
197218c64b6SSatish Balay   }
198029af93fSBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,count,iary,&is2);CHKERRQ(ierr);
199218c64b6SSatish Balay   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
200218c64b6SSatish Balay   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
201fca92195SBarry Smith   ierr = PetscFree2(vary,iary);CHKERRQ(ierr);
202218c64b6SSatish Balay 
2034aa3045dSJed Brown   ierr = MatGetSubMatrix_SeqBAIJ_Private(A,is1,is2,scall,B);CHKERRQ(ierr);
204218c64b6SSatish Balay   ISDestroy(is1);
205218c64b6SSatish Balay   ISDestroy(is2);
2063a40ed3dSBarry Smith   PetscFunctionReturn(0);
207218c64b6SSatish Balay }
208218c64b6SSatish Balay 
2094a2ae208SSatish Balay #undef __FUNCT__
2104a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_SeqBAIJ"
211690b6cddSBarry Smith PetscErrorCode MatGetSubMatrices_SeqBAIJ(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
212736121d4SSatish Balay {
2136849ba73SBarry Smith   PetscErrorCode ierr;
214690b6cddSBarry Smith   PetscInt       i;
215736121d4SSatish Balay 
2163a40ed3dSBarry Smith   PetscFunctionBegin;
217736121d4SSatish Balay   if (scall == MAT_INITIAL_MATRIX) {
21882502324SSatish Balay     ierr = PetscMalloc((n+1)*sizeof(Mat),B);CHKERRQ(ierr);
219736121d4SSatish Balay   }
220736121d4SSatish Balay 
221736121d4SSatish Balay   for (i=0; i<n; i++) {
2224aa3045dSJed Brown     ierr = MatGetSubMatrix_SeqBAIJ(A,irow[i],icol[i],scall,&(*B)[i]);CHKERRQ(ierr);
223736121d4SSatish Balay   }
2243a40ed3dSBarry Smith   PetscFunctionReturn(0);
225736121d4SSatish Balay }
226218c64b6SSatish Balay 
227218c64b6SSatish Balay 
2282d61bbb3SSatish Balay /* -------------------------------------------------------*/
2292d61bbb3SSatish Balay /* Should check that shapes of vectors and matrices match */
2302d61bbb3SSatish Balay /* -------------------------------------------------------*/
2312d61bbb3SSatish Balay 
2324a2ae208SSatish Balay #undef __FUNCT__
2334a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_1"
234dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_1(Mat A,Vec xx,Vec zz)
2352d61bbb3SSatish Balay {
2362d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
237d9fead3dSBarry Smith   PetscScalar       *z,sum;
238d9fead3dSBarry Smith   const PetscScalar *x;
239d9fead3dSBarry Smith   const MatScalar   *v;
2406849ba73SBarry Smith   PetscErrorCode    ierr;
2412162cab8SBarry Smith   PetscInt          mbs,i,n,nonzerorow=0;
2422162cab8SBarry Smith   const PetscInt    *idx,*ii,*ridx=PETSC_NULL;
24326e093fcSHong Zhang   PetscTruth        usecprow=a->compressedrow.use;
2442d61bbb3SSatish Balay 
2452d61bbb3SSatish Balay   PetscFunctionBegin;
246d9fead3dSBarry Smith   ierr = VecGetArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
2471ebc52fbSHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
2482d61bbb3SSatish Balay 
24926e093fcSHong Zhang   if (usecprow){
25026e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
25126e093fcSHong Zhang     ii   = a->compressedrow.i;
2527b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
253a1c3900fSBarry Smith     ierr = PetscMemzero(z,mbs*sizeof(PetscScalar));CHKERRQ(ierr);
25426e093fcSHong Zhang   } else {
25526e093fcSHong Zhang     mbs = a->mbs;
2562d61bbb3SSatish Balay     ii  = a->i;
25726e093fcSHong Zhang   }
2582d61bbb3SSatish Balay 
2592d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
260ee54c7eeSHong Zhang     n    = ii[1] - ii[0];
261ee54c7eeSHong Zhang     v    = a->a + ii[0];
262ee54c7eeSHong Zhang     idx  = a->j + ii[0];
263ee54c7eeSHong Zhang     ii++;
2642d61bbb3SSatish Balay     sum  = 0.0;
2652162cab8SBarry Smith     PetscSparseDensePlusDot(sum,x,v,idx,n);
26626e093fcSHong Zhang     if (usecprow){
2677b2bb3b9SHong Zhang       z[ridx[i]] = sum;
26826e093fcSHong Zhang     } else {
269122f12eaSBarry Smith       nonzerorow += (n>0);
2702d61bbb3SSatish Balay       z[i] = sum;
2712d61bbb3SSatish Balay     }
27226e093fcSHong Zhang   }
273d9fead3dSBarry Smith   ierr = VecRestoreArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
2741ebc52fbSHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
275dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*a->nz - nonzerorow);CHKERRQ(ierr);
2762d61bbb3SSatish Balay   PetscFunctionReturn(0);
2772d61bbb3SSatish Balay }
2782d61bbb3SSatish Balay 
2794a2ae208SSatish Balay #undef __FUNCT__
2804a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_2"
281dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_2(Mat A,Vec xx,Vec zz)
2822d61bbb3SSatish Balay {
2832d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
284d9fead3dSBarry Smith   PetscScalar       *z = 0,sum1,sum2,*zarray;
285d9fead3dSBarry Smith   const PetscScalar *x,*xb;
28687828ca2SBarry Smith   PetscScalar       x1,x2;
287d9fead3dSBarry Smith   const MatScalar   *v;
288dfbe8321SBarry Smith   PetscErrorCode    ierr;
28998c9bda7SSatish Balay   PetscInt          mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL,nonzerorow=0;
29026e093fcSHong Zhang   PetscTruth        usecprow=a->compressedrow.use;
2912d61bbb3SSatish Balay 
2922d61bbb3SSatish Balay   PetscFunctionBegin;
293d9fead3dSBarry Smith   ierr = VecGetArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
29426e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
2952d61bbb3SSatish Balay 
2962d61bbb3SSatish Balay   idx = a->j;
2972d61bbb3SSatish Balay   v   = a->a;
29826e093fcSHong Zhang   if (usecprow){
29926e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
30026e093fcSHong Zhang     ii   = a->compressedrow.i;
3017b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
30226e093fcSHong Zhang   } else {
30326e093fcSHong Zhang     mbs = a->mbs;
3042d61bbb3SSatish Balay     ii  = a->i;
30526e093fcSHong Zhang     z   = zarray;
30626e093fcSHong Zhang   }
3072d61bbb3SSatish Balay 
3082d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
3092d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
3102d61bbb3SSatish Balay     sum1 = 0.0; sum2 = 0.0;
31198c9bda7SSatish Balay     nonzerorow += (n>0);
3122d61bbb3SSatish Balay     for (j=0; j<n; j++) {
3132d61bbb3SSatish Balay       xb = x + 2*(*idx++); x1 = xb[0]; x2 = xb[1];
3142d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[2]*x2;
3152d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[3]*x2;
3162d61bbb3SSatish Balay       v += 4;
3172d61bbb3SSatish Balay     }
3187b2bb3b9SHong Zhang     if (usecprow) z = zarray + 2*ridx[i];
3192d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2;
32026e093fcSHong Zhang     if (!usecprow) z += 2;
3212d61bbb3SSatish Balay   }
322d9fead3dSBarry Smith   ierr = VecRestoreArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
32326e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
324dc0b31edSSatish Balay   ierr = PetscLogFlops(8.0*a->nz - 2.0*nonzerorow);CHKERRQ(ierr);
3252d61bbb3SSatish Balay   PetscFunctionReturn(0);
3262d61bbb3SSatish Balay }
3272d61bbb3SSatish Balay 
3284a2ae208SSatish Balay #undef __FUNCT__
3294a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_3"
330dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_3(Mat A,Vec xx,Vec zz)
3312d61bbb3SSatish Balay {
3322d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
333d9fead3dSBarry Smith   PetscScalar       *z = 0,sum1,sum2,sum3,x1,x2,x3,*zarray;
334d9fead3dSBarry Smith   const PetscScalar *x,*xb;
335d9fead3dSBarry Smith   const MatScalar   *v;
336dfbe8321SBarry Smith   PetscErrorCode    ierr;
337028cd4eaSSatish Balay   PetscInt          mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL,nonzerorow=0;
338028cd4eaSSatish Balay   PetscTruth        usecprow=a->compressedrow.use;
33926e093fcSHong Zhang 
3402d61bbb3SSatish Balay 
341b6410449SSatish Balay #if defined(PETSC_HAVE_PRAGMA_DISJOINT)
342fee21e36SBarry Smith #pragma disjoint(*v,*z,*xb)
343fee21e36SBarry Smith #endif
344fee21e36SBarry Smith 
3452d61bbb3SSatish Balay   PetscFunctionBegin;
346d9fead3dSBarry Smith   ierr = VecGetArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
34726e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
3482d61bbb3SSatish Balay 
3492d61bbb3SSatish Balay   idx = a->j;
3502d61bbb3SSatish Balay   v   = a->a;
35126e093fcSHong Zhang   if (usecprow){
35226e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
35326e093fcSHong Zhang     ii   = a->compressedrow.i;
3547b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
35526e093fcSHong Zhang   } else {
35626e093fcSHong Zhang     mbs = a->mbs;
3572d61bbb3SSatish Balay     ii  = a->i;
35826e093fcSHong Zhang     z   = zarray;
35926e093fcSHong Zhang   }
3602d61bbb3SSatish Balay 
3612d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
3622d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
3632d61bbb3SSatish Balay     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0;
36498c9bda7SSatish Balay     nonzerorow += (n>0);
3652d61bbb3SSatish Balay     for (j=0; j<n; j++) {
3662d61bbb3SSatish Balay       xb = x + 3*(*idx++); x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
3672d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[3]*x2 + v[6]*x3;
3682d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[4]*x2 + v[7]*x3;
3692d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[5]*x2 + v[8]*x3;
3702d61bbb3SSatish Balay       v += 9;
3712d61bbb3SSatish Balay     }
3727b2bb3b9SHong Zhang     if (usecprow) z = zarray + 3*ridx[i];
3732d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3;
37426e093fcSHong Zhang     if (!usecprow) z += 3;
3752d61bbb3SSatish Balay   }
376d9fead3dSBarry Smith   ierr = VecRestoreArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
37726e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
378dc0b31edSSatish Balay   ierr = PetscLogFlops(18.0*a->nz - 3.0*nonzerorow);CHKERRQ(ierr);
3792d61bbb3SSatish Balay   PetscFunctionReturn(0);
3802d61bbb3SSatish Balay }
3812d61bbb3SSatish Balay 
3824a2ae208SSatish Balay #undef __FUNCT__
3834a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_4"
384dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_4(Mat A,Vec xx,Vec zz)
3852d61bbb3SSatish Balay {
3862d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
387d9fead3dSBarry Smith   PetscScalar       *z = 0,sum1,sum2,sum3,sum4,x1,x2,x3,x4,*zarray;
388d9fead3dSBarry Smith   const PetscScalar *x,*xb;
389d9fead3dSBarry Smith   const MatScalar   *v;
390dfbe8321SBarry Smith   PetscErrorCode    ierr;
39198c9bda7SSatish Balay   PetscInt          mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL,nonzerorow=0;
39226e093fcSHong Zhang   PetscTruth        usecprow=a->compressedrow.use;
3932d61bbb3SSatish Balay 
3942d61bbb3SSatish Balay   PetscFunctionBegin;
395d9fead3dSBarry Smith   ierr = VecGetArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
39626e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
3972d61bbb3SSatish Balay 
3982d61bbb3SSatish Balay   idx = a->j;
3992d61bbb3SSatish Balay   v   = a->a;
40026e093fcSHong Zhang   if (usecprow){
40126e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
40226e093fcSHong Zhang     ii   = a->compressedrow.i;
4037b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
40426e093fcSHong Zhang   } else {
40526e093fcSHong Zhang     mbs = a->mbs;
4062d61bbb3SSatish Balay     ii  = a->i;
40726e093fcSHong Zhang     z   = zarray;
40826e093fcSHong Zhang   }
4092d61bbb3SSatish Balay 
4102d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
4112d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
4122d61bbb3SSatish Balay     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0;
41398c9bda7SSatish Balay     nonzerorow += (n>0);
4142d61bbb3SSatish Balay     for (j=0; j<n; j++) {
4152d61bbb3SSatish Balay       xb = x + 4*(*idx++);
4162d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
4172d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[4]*x2 + v[8]*x3   + v[12]*x4;
4182d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[5]*x2 + v[9]*x3   + v[13]*x4;
4192d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[6]*x2 + v[10]*x3  + v[14]*x4;
4202d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[7]*x2 + v[11]*x3  + v[15]*x4;
4212d61bbb3SSatish Balay       v += 16;
4222d61bbb3SSatish Balay     }
4237b2bb3b9SHong Zhang     if (usecprow) z = zarray + 4*ridx[i];
4242d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4;
42526e093fcSHong Zhang     if (!usecprow) z += 4;
4262d61bbb3SSatish Balay   }
427d9fead3dSBarry Smith   ierr = VecRestoreArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
42826e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
429dc0b31edSSatish Balay   ierr = PetscLogFlops(32.0*a->nz - 4.0*nonzerorow);CHKERRQ(ierr);
4302d61bbb3SSatish Balay   PetscFunctionReturn(0);
4312d61bbb3SSatish Balay }
4322d61bbb3SSatish Balay 
4334a2ae208SSatish Balay #undef __FUNCT__
4344a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_5"
435dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_5(Mat A,Vec xx,Vec zz)
4362d61bbb3SSatish Balay {
4372d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
438d9fead3dSBarry Smith   PetscScalar       sum1,sum2,sum3,sum4,sum5,x1,x2,x3,x4,x5,*z = 0,*zarray;
439d9fead3dSBarry Smith   const PetscScalar *xb,*x;
440d9fead3dSBarry Smith   const MatScalar   *v;
441dfbe8321SBarry Smith   PetscErrorCode    ierr;
442766f9fbaSBarry Smith   const PetscInt    *idx,*ii,*ridx=PETSC_NULL;
443766f9fbaSBarry Smith   PetscInt          mbs,i,j,n,nonzerorow=0;
44426e093fcSHong Zhang   PetscTruth        usecprow=a->compressedrow.use;
4452d61bbb3SSatish Balay 
446433994e6SBarry Smith   PetscFunctionBegin;
447d9fead3dSBarry Smith   ierr = VecGetArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
44826e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
4492d61bbb3SSatish Balay 
4502d61bbb3SSatish Balay   idx = a->j;
4512d61bbb3SSatish Balay   v   = a->a;
45226e093fcSHong Zhang   if (usecprow){
45326e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
45426e093fcSHong Zhang     ii   = a->compressedrow.i;
4557b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
45626e093fcSHong Zhang   } else {
45726e093fcSHong Zhang     mbs = a->mbs;
4582d61bbb3SSatish Balay     ii  = a->i;
45926e093fcSHong Zhang     z   = zarray;
46026e093fcSHong Zhang   }
4612d61bbb3SSatish Balay 
4622d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
4632d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
4642d61bbb3SSatish Balay     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0;
46598c9bda7SSatish Balay     nonzerorow += (n>0);
4662d61bbb3SSatish Balay     for (j=0; j<n; j++) {
4672d61bbb3SSatish Balay       xb = x + 5*(*idx++);
4682d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4];
4692d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[5]*x2 + v[10]*x3  + v[15]*x4 + v[20]*x5;
4702d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[6]*x2 + v[11]*x3  + v[16]*x4 + v[21]*x5;
4712d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[7]*x2 + v[12]*x3  + v[17]*x4 + v[22]*x5;
4722d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[8]*x2 + v[13]*x3  + v[18]*x4 + v[23]*x5;
4732d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[9]*x2 + v[14]*x3  + v[19]*x4 + v[24]*x5;
4742d61bbb3SSatish Balay       v += 25;
4752d61bbb3SSatish Balay     }
4767b2bb3b9SHong Zhang     if (usecprow) z = zarray + 5*ridx[i];
4772d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5;
47826e093fcSHong Zhang     if (!usecprow) z += 5;
4792d61bbb3SSatish Balay   }
480d9fead3dSBarry Smith   ierr = VecRestoreArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
48126e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
482dc0b31edSSatish Balay   ierr = PetscLogFlops(50.0*a->nz - 5.0*nonzerorow);CHKERRQ(ierr);
4832d61bbb3SSatish Balay   PetscFunctionReturn(0);
4842d61bbb3SSatish Balay }
4852d61bbb3SSatish Balay 
48615091d37SBarry Smith 
4874a2ae208SSatish Balay #undef __FUNCT__
4884a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_6"
489dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_6(Mat A,Vec xx,Vec zz)
49015091d37SBarry Smith {
49115091d37SBarry Smith   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
492d9fead3dSBarry Smith   PetscScalar       *z = 0,sum1,sum2,sum3,sum4,sum5,sum6;
493d9fead3dSBarry Smith   const PetscScalar *x,*xb;
49426e093fcSHong Zhang   PetscScalar       x1,x2,x3,x4,x5,x6,*zarray;
495d9fead3dSBarry Smith   const MatScalar   *v;
496dfbe8321SBarry Smith   PetscErrorCode    ierr;
49798c9bda7SSatish Balay   PetscInt          mbs=a->mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL,nonzerorow=0;
49826e093fcSHong Zhang   PetscTruth        usecprow=a->compressedrow.use;
49915091d37SBarry Smith 
500433994e6SBarry Smith   PetscFunctionBegin;
501d9fead3dSBarry Smith   ierr = VecGetArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
50226e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
50315091d37SBarry Smith 
50415091d37SBarry Smith   idx = a->j;
50515091d37SBarry Smith   v   = a->a;
50626e093fcSHong Zhang   if (usecprow){
50726e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
50826e093fcSHong Zhang     ii   = a->compressedrow.i;
5097b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
51026e093fcSHong Zhang   } else {
51126e093fcSHong Zhang     mbs = a->mbs;
51215091d37SBarry Smith     ii  = a->i;
51326e093fcSHong Zhang     z   = zarray;
51426e093fcSHong Zhang   }
51515091d37SBarry Smith 
51615091d37SBarry Smith   for (i=0; i<mbs; i++) {
51715091d37SBarry Smith     n  = ii[1] - ii[0]; ii++;
51815091d37SBarry Smith     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0;
51998c9bda7SSatish Balay     nonzerorow += (n>0);
52015091d37SBarry Smith     for (j=0; j<n; j++) {
52115091d37SBarry Smith       xb = x + 6*(*idx++);
52215091d37SBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5];
52315091d37SBarry Smith       sum1 += v[0]*x1 + v[6]*x2  + v[12]*x3  + v[18]*x4 + v[24]*x5 + v[30]*x6;
52415091d37SBarry Smith       sum2 += v[1]*x1 + v[7]*x2  + v[13]*x3  + v[19]*x4 + v[25]*x5 + v[31]*x6;
52515091d37SBarry Smith       sum3 += v[2]*x1 + v[8]*x2  + v[14]*x3  + v[20]*x4 + v[26]*x5 + v[32]*x6;
52615091d37SBarry Smith       sum4 += v[3]*x1 + v[9]*x2  + v[15]*x3  + v[21]*x4 + v[27]*x5 + v[33]*x6;
52715091d37SBarry Smith       sum5 += v[4]*x1 + v[10]*x2 + v[16]*x3  + v[22]*x4 + v[28]*x5 + v[34]*x6;
52815091d37SBarry Smith       sum6 += v[5]*x1 + v[11]*x2 + v[17]*x3  + v[23]*x4 + v[29]*x5 + v[35]*x6;
52915091d37SBarry Smith       v += 36;
53015091d37SBarry Smith     }
5317b2bb3b9SHong Zhang     if (usecprow) z = zarray + 6*ridx[i];
53215091d37SBarry Smith     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6;
53326e093fcSHong Zhang     if (!usecprow) z += 6;
53415091d37SBarry Smith   }
53515091d37SBarry Smith 
536d9fead3dSBarry Smith   ierr = VecRestoreArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
53726e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
538dc0b31edSSatish Balay   ierr = PetscLogFlops(72.0*a->nz - 6.0*nonzerorow);CHKERRQ(ierr);
53915091d37SBarry Smith   PetscFunctionReturn(0);
54015091d37SBarry Smith }
5418ab949d8SShri Abhyankar 
5424a2ae208SSatish Balay #undef __FUNCT__
5434a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_7"
544dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_7(Mat A,Vec xx,Vec zz)
5452d61bbb3SSatish Balay {
5462d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
547d9fead3dSBarry Smith   PetscScalar       *z = 0,sum1,sum2,sum3,sum4,sum5,sum6,sum7;
548d9fead3dSBarry Smith   const PetscScalar *x,*xb;
54926e093fcSHong Zhang   PetscScalar       x1,x2,x3,x4,x5,x6,x7,*zarray;
550d9fead3dSBarry Smith   const MatScalar   *v;
551dfbe8321SBarry Smith   PetscErrorCode    ierr;
55298c9bda7SSatish Balay   PetscInt          mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL,nonzerorow=0;
55326e093fcSHong Zhang   PetscTruth        usecprow=a->compressedrow.use;
5542d61bbb3SSatish Balay 
555433994e6SBarry Smith   PetscFunctionBegin;
556d9fead3dSBarry Smith   ierr = VecGetArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
55726e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
5582d61bbb3SSatish Balay 
5592d61bbb3SSatish Balay   idx = a->j;
5602d61bbb3SSatish Balay   v   = a->a;
56126e093fcSHong Zhang   if (usecprow){
56226e093fcSHong Zhang     mbs    = a->compressedrow.nrows;
56326e093fcSHong Zhang     ii     = a->compressedrow.i;
5647b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
56526e093fcSHong Zhang   } else {
56626e093fcSHong Zhang     mbs = a->mbs;
5672d61bbb3SSatish Balay     ii  = a->i;
56826e093fcSHong Zhang     z   = zarray;
56926e093fcSHong Zhang   }
5702d61bbb3SSatish Balay 
5712d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
5722d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
5732d61bbb3SSatish Balay     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0; sum7 = 0.0;
57498c9bda7SSatish Balay     nonzerorow += (n>0);
5752d61bbb3SSatish Balay     for (j=0; j<n; j++) {
5762d61bbb3SSatish Balay       xb = x + 7*(*idx++);
5772d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6];
5782d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[7]*x2  + v[14]*x3  + v[21]*x4 + v[28]*x5 + v[35]*x6 + v[42]*x7;
5792d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[8]*x2  + v[15]*x3  + v[22]*x4 + v[29]*x5 + v[36]*x6 + v[43]*x7;
5802d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[9]*x2  + v[16]*x3  + v[23]*x4 + v[30]*x5 + v[37]*x6 + v[44]*x7;
5812d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[10]*x2 + v[17]*x3  + v[24]*x4 + v[31]*x5 + v[38]*x6 + v[45]*x7;
5822d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[11]*x2 + v[18]*x3  + v[25]*x4 + v[32]*x5 + v[39]*x6 + v[46]*x7;
5832d61bbb3SSatish Balay       sum6 += v[5]*x1 + v[12]*x2 + v[19]*x3  + v[26]*x4 + v[33]*x5 + v[40]*x6 + v[47]*x7;
5842d61bbb3SSatish Balay       sum7 += v[6]*x1 + v[13]*x2 + v[20]*x3  + v[27]*x4 + v[34]*x5 + v[41]*x6 + v[48]*x7;
5852d61bbb3SSatish Balay       v += 49;
5862d61bbb3SSatish Balay     }
5877b2bb3b9SHong Zhang     if (usecprow) z = zarray + 7*ridx[i];
5882d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
58926e093fcSHong Zhang     if (!usecprow) z += 7;
5902d61bbb3SSatish Balay   }
5912d61bbb3SSatish Balay 
592d9fead3dSBarry Smith   ierr = VecRestoreArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
59326e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
594dc0b31edSSatish Balay   ierr = PetscLogFlops(98.0*a->nz - 7.0*nonzerorow);CHKERRQ(ierr);
5952d61bbb3SSatish Balay   PetscFunctionReturn(0);
5962d61bbb3SSatish Balay }
5972d61bbb3SSatish Balay 
5988ab949d8SShri Abhyankar /* MatMult_SeqBAIJ_15 version 1: Columns in the block are accessed one at a time */
599832cc040SShri Abhyankar /* Default MatMult for block size 15 */
6008ab949d8SShri Abhyankar 
6018ab949d8SShri Abhyankar #undef __FUNCT__
6028ab949d8SShri Abhyankar #define __FUNCT__ "MatMult_SeqBAIJ_15_ver1"
6038ab949d8SShri Abhyankar PetscErrorCode MatMult_SeqBAIJ_15_ver1(Mat A,Vec xx,Vec zz)
6048ab949d8SShri Abhyankar {
6058ab949d8SShri Abhyankar   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
6068ab949d8SShri Abhyankar   PetscScalar       *z = 0,sum1,sum2,sum3,sum4,sum5,sum6,sum7,sum8,sum9,sum10,sum11,sum12,sum13,sum14,sum15;
6078ab949d8SShri Abhyankar   const PetscScalar *x,*xb;
6088ab949d8SShri Abhyankar   PetscScalar        x1,x2,x3,x4,x5,x6,x7,x8,x9,x10,x11,x12,x13,x14,x15,*zarray;
6098ab949d8SShri Abhyankar   const MatScalar   *v;
6108ab949d8SShri Abhyankar   PetscErrorCode    ierr;
6118ab949d8SShri Abhyankar   const PetscInt    *ii,*ij=a->j,*idx;
6128ab949d8SShri Abhyankar   PetscInt          mbs,i,j,k,n,*ridx=PETSC_NULL,nonzerorow=0;
6138ab949d8SShri Abhyankar   PetscTruth        usecprow=a->compressedrow.use;
6148ab949d8SShri Abhyankar 
6158ab949d8SShri Abhyankar   PetscFunctionBegin;
6168ab949d8SShri Abhyankar   ierr = VecGetArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
6178ab949d8SShri Abhyankar   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
6188ab949d8SShri Abhyankar 
6198ab949d8SShri Abhyankar   v   = a->a;
6208ab949d8SShri Abhyankar   if (usecprow){
6218ab949d8SShri Abhyankar     mbs    = a->compressedrow.nrows;
6228ab949d8SShri Abhyankar     ii     = a->compressedrow.i;
6238ab949d8SShri Abhyankar     ridx = a->compressedrow.rindex;
6248ab949d8SShri Abhyankar   } else {
6258ab949d8SShri Abhyankar     mbs = a->mbs;
6268ab949d8SShri Abhyankar     ii  = a->i;
6278ab949d8SShri Abhyankar     z   = zarray;
6288ab949d8SShri Abhyankar   }
6298ab949d8SShri Abhyankar 
6308ab949d8SShri Abhyankar   for (i=0; i<mbs; i++) {
6318ab949d8SShri Abhyankar     n  = ii[i+1] - ii[i];
6328ab949d8SShri Abhyankar     idx = ij + ii[i];
6338ab949d8SShri Abhyankar     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0; sum7 = 0.0;
6348ab949d8SShri Abhyankar     sum8 = 0.0; sum9 = 0.0; sum10 = 0.0; sum11 = 0.0; sum12 = 0.0; sum13 = 0.0; sum14 = 0.0;sum15 = 0.0;
6358ab949d8SShri Abhyankar 
6368ab949d8SShri Abhyankar     nonzerorow += (n>0);
6378ab949d8SShri Abhyankar     for (j=0; j<n; j++) {
6388ab949d8SShri Abhyankar       xb = x + 15*(idx[j]);
6398ab949d8SShri Abhyankar       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6];
6408ab949d8SShri Abhyankar       x8 = xb[7]; x9 = xb[8]; x10 = xb[9]; x11 = xb[10]; x12 = xb[11]; x13 = xb[12]; x14 = xb[13];x15 = xb[14];
6418ab949d8SShri Abhyankar 
6428ab949d8SShri Abhyankar       for(k=0;k<15;k++){
6438ab949d8SShri Abhyankar 	sum1  += v[0]*xb[k];
6448ab949d8SShri Abhyankar         sum2  += v[1]*xb[k];
6458ab949d8SShri Abhyankar 	sum3  += v[2]*xb[k];
6468ab949d8SShri Abhyankar 	sum4  += v[3]*xb[k];
6478ab949d8SShri Abhyankar 	sum5  += v[4]*xb[k];
6488ab949d8SShri Abhyankar         sum6  += v[5]*xb[k];
6498ab949d8SShri Abhyankar 	sum7  += v[6]*xb[k];
6508ab949d8SShri Abhyankar 	sum8  += v[7]*xb[k];
6518ab949d8SShri Abhyankar         sum9  += v[8]*xb[k];
6528ab949d8SShri Abhyankar         sum10 += v[9]*xb[k];
6538ab949d8SShri Abhyankar 	sum11 += v[10]*xb[k];
6548ab949d8SShri Abhyankar 	sum12 += v[11]*xb[k];
6558ab949d8SShri Abhyankar 	sum13 += v[12]*xb[k];
6568ab949d8SShri Abhyankar         sum14 += v[13]*xb[k];
6578ab949d8SShri Abhyankar 	sum15 += v[14]*xb[k];
6588ab949d8SShri Abhyankar 	v += 15;
6598ab949d8SShri Abhyankar       }
6608ab949d8SShri Abhyankar     }
6618ab949d8SShri Abhyankar     if (usecprow) z = zarray + 15*ridx[i];
6628ab949d8SShri Abhyankar     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
6638ab949d8SShri Abhyankar     z[7] = sum8; z[8] = sum9; z[9] = sum10; z[10] = sum11; z[11] = sum12; z[12] = sum13; z[13] = sum14;z[14] = sum15;
6648ab949d8SShri Abhyankar 
6658ab949d8SShri Abhyankar     if (!usecprow) z += 15;
6668ab949d8SShri Abhyankar   }
6678ab949d8SShri Abhyankar 
6688ab949d8SShri Abhyankar   ierr = VecRestoreArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
6698ab949d8SShri Abhyankar   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
6708ab949d8SShri Abhyankar   ierr = PetscLogFlops(450.0*a->nz - 15.0*nonzerorow);CHKERRQ(ierr);
6718ab949d8SShri Abhyankar   PetscFunctionReturn(0);
6728ab949d8SShri Abhyankar }
6738ab949d8SShri Abhyankar 
6748ab949d8SShri Abhyankar /* MatMult_SeqBAIJ_15_ver2 : Columns in the block are accessed in sets of 4,4,4,3 */
6758ab949d8SShri Abhyankar #undef __FUNCT__
6768ab949d8SShri Abhyankar #define __FUNCT__ "MatMult_SeqBAIJ_15_ver2"
6778ab949d8SShri Abhyankar PetscErrorCode MatMult_SeqBAIJ_15_ver2(Mat A,Vec xx,Vec zz)
6788ab949d8SShri Abhyankar {
6798ab949d8SShri Abhyankar   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
6808ab949d8SShri Abhyankar   PetscScalar       *z = 0,sum1,sum2,sum3,sum4,sum5,sum6,sum7,sum8,sum9,sum10,sum11,sum12,sum13,sum14,sum15;
6818ab949d8SShri Abhyankar   const PetscScalar *x,*xb;
6820b8f6341SShri Abhyankar   PetscScalar        x1,x2,x3,x4,*zarray;
6838ab949d8SShri Abhyankar   const MatScalar   *v;
6848ab949d8SShri Abhyankar   PetscErrorCode    ierr;
6858ab949d8SShri Abhyankar   const PetscInt    *ii,*ij=a->j,*idx;
6868ab949d8SShri Abhyankar   PetscInt          mbs,i,j,n,*ridx=PETSC_NULL,nonzerorow=0;
6878ab949d8SShri Abhyankar   PetscTruth        usecprow=a->compressedrow.use;
6888ab949d8SShri Abhyankar 
6898ab949d8SShri Abhyankar   PetscFunctionBegin;
6908ab949d8SShri Abhyankar   ierr = VecGetArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
6918ab949d8SShri Abhyankar   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
6928ab949d8SShri Abhyankar 
6938ab949d8SShri Abhyankar   v   = a->a;
6948ab949d8SShri Abhyankar   if (usecprow){
6958ab949d8SShri Abhyankar     mbs    = a->compressedrow.nrows;
6968ab949d8SShri Abhyankar     ii     = a->compressedrow.i;
6978ab949d8SShri Abhyankar     ridx = a->compressedrow.rindex;
6988ab949d8SShri Abhyankar   } else {
6998ab949d8SShri Abhyankar     mbs = a->mbs;
7008ab949d8SShri Abhyankar     ii  = a->i;
7018ab949d8SShri Abhyankar     z   = zarray;
7028ab949d8SShri Abhyankar   }
7038ab949d8SShri Abhyankar 
7048ab949d8SShri Abhyankar   for (i=0; i<mbs; i++) {
7058ab949d8SShri Abhyankar     n  = ii[i+1] - ii[i];
7068ab949d8SShri Abhyankar     idx = ij + ii[i];
7078ab949d8SShri Abhyankar     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0; sum7 = 0.0;
7088ab949d8SShri Abhyankar     sum8 = 0.0; sum9 = 0.0; sum10 = 0.0; sum11 = 0.0; sum12 = 0.0; sum13 = 0.0; sum14 = 0.0;sum15 = 0.0;
7098ab949d8SShri Abhyankar 
7108ab949d8SShri Abhyankar     nonzerorow += (n>0);
7118ab949d8SShri Abhyankar     for (j=0; j<n; j++) {
7128ab949d8SShri Abhyankar       xb = x + 15*(idx[j]);
7130b8f6341SShri Abhyankar       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
7148ab949d8SShri Abhyankar 
7158ab949d8SShri Abhyankar       sum1  += v[0]*x1 + v[15]*x2 + v[30]*x3   + v[45]*x4;
7168ab949d8SShri Abhyankar       sum2  += v[1]*x1 + v[16]*x2 + v[31]*x3   + v[46]*x4;
7178ab949d8SShri Abhyankar       sum3  += v[2]*x1 + v[17]*x2 + v[32]*x3  + v[47]*x4;
7188ab949d8SShri Abhyankar       sum4  += v[3]*x1 + v[18]*x2 + v[33]*x3  + v[48]*x4;
7198ab949d8SShri Abhyankar       sum5  += v[4]*x1 + v[19]*x2 + v[34]*x3   + v[49]*x4;
7208ab949d8SShri Abhyankar       sum6  += v[5]*x1 + v[20]*x2 + v[35]*x3   + v[50]*x4;
7218ab949d8SShri Abhyankar       sum7  += v[6]*x1 + v[21]*x2 + v[36]*x3  + v[51]*x4;
7228ab949d8SShri Abhyankar       sum8  += v[7]*x1 + v[22]*x2 + v[37]*x3  + v[52]*x4;
7238ab949d8SShri Abhyankar       sum9  += v[8]*x1 + v[23]*x2 + v[38]*x3   + v[53]*x4;
7248ab949d8SShri Abhyankar       sum10 += v[9]*x1 + v[24]*x2 + v[39]*x3   + v[54]*x4;
7258ab949d8SShri Abhyankar       sum11 += v[10]*x1 + v[25]*x2 + v[40]*x3  + v[55]*x4;
7268ab949d8SShri Abhyankar       sum12 += v[11]*x1 + v[26]*x2 + v[41]*x3  + v[56]*x4;
7278ab949d8SShri Abhyankar       sum13 += v[12]*x1 + v[27]*x2 + v[42]*x3   + v[57]*x4;
7288ab949d8SShri Abhyankar       sum14 += v[13]*x1 + v[28]*x2 + v[43]*x3   + v[58]*x4;
7298ab949d8SShri Abhyankar       sum15 += v[14]*x1 + v[29]*x2 + v[44]*x3  + v[59]*x4;
7308ab949d8SShri Abhyankar 
7318ab949d8SShri Abhyankar       v += 60;
7328ab949d8SShri Abhyankar 
7330b8f6341SShri Abhyankar       x1 = xb[4]; x2 = xb[5]; x3 = xb[6]; x4 = xb[7];
7348ab949d8SShri Abhyankar 
7358ab949d8SShri Abhyankar       sum1  += v[0]*x1 + v[15]*x2 + v[30]*x3   + v[45]*x4;
7368ab949d8SShri Abhyankar       sum2  += v[1]*x1 + v[16]*x2 + v[31]*x3   + v[46]*x4;
7378ab949d8SShri Abhyankar       sum3  += v[2]*x1 + v[17]*x2 + v[32]*x3  + v[47]*x4;
7388ab949d8SShri Abhyankar       sum4  += v[3]*x1 + v[18]*x2 + v[33]*x3  + v[48]*x4;
7398ab949d8SShri Abhyankar       sum5  += v[4]*x1 + v[19]*x2 + v[34]*x3   + v[49]*x4;
7408ab949d8SShri Abhyankar       sum6  += v[5]*x1 + v[20]*x2 + v[35]*x3   + v[50]*x4;
7418ab949d8SShri Abhyankar       sum7  += v[6]*x1 + v[21]*x2 + v[36]*x3  + v[51]*x4;
7428ab949d8SShri Abhyankar       sum8  += v[7]*x1 + v[22]*x2 + v[37]*x3  + v[52]*x4;
7438ab949d8SShri Abhyankar       sum9  += v[8]*x1 + v[23]*x2 + v[38]*x3   + v[53]*x4;
7448ab949d8SShri Abhyankar       sum10 += v[9]*x1 + v[24]*x2 + v[39]*x3   + v[54]*x4;
7458ab949d8SShri Abhyankar       sum11 += v[10]*x1 + v[25]*x2 + v[40]*x3  + v[55]*x4;
7468ab949d8SShri Abhyankar       sum12 += v[11]*x1 + v[26]*x2 + v[41]*x3  + v[56]*x4;
7478ab949d8SShri Abhyankar       sum13 += v[12]*x1 + v[27]*x2 + v[42]*x3   + v[57]*x4;
7488ab949d8SShri Abhyankar       sum14 += v[13]*x1 + v[28]*x2 + v[43]*x3   + v[58]*x4;
7498ab949d8SShri Abhyankar       sum15 += v[14]*x1 + v[29]*x2 + v[44]*x3  + v[59]*x4;
7508ab949d8SShri Abhyankar       v += 60;
7518ab949d8SShri Abhyankar 
7520b8f6341SShri Abhyankar       x1 = xb[8]; x2 = xb[9]; x3 = xb[10]; x4 = xb[11];
7530b8f6341SShri Abhyankar       sum1  += v[0]*x1 + v[15]*x2 + v[30]*x3   + v[45]*x4;
7540b8f6341SShri Abhyankar       sum2  += v[1]*x1 + v[16]*x2 + v[31]*x3   + v[46]*x4;
7550b8f6341SShri Abhyankar       sum3  += v[2]*x1 + v[17]*x2 + v[32]*x3  + v[47]*x4;
7560b8f6341SShri Abhyankar       sum4  += v[3]*x1 + v[18]*x2 + v[33]*x3  + v[48]*x4;
7570b8f6341SShri Abhyankar       sum5  += v[4]*x1 + v[19]*x2 + v[34]*x3   + v[49]*x4;
7580b8f6341SShri Abhyankar       sum6  += v[5]*x1 + v[20]*x2 + v[35]*x3   + v[50]*x4;
7590b8f6341SShri Abhyankar       sum7  += v[6]*x1 + v[21]*x2 + v[36]*x3  + v[51]*x4;
7600b8f6341SShri Abhyankar       sum8  += v[7]*x1 + v[22]*x2 + v[37]*x3  + v[52]*x4;
7610b8f6341SShri Abhyankar       sum9  += v[8]*x1 + v[23]*x2 + v[38]*x3   + v[53]*x4;
7620b8f6341SShri Abhyankar       sum10 += v[9]*x1 + v[24]*x2 + v[39]*x3   + v[54]*x4;
7630b8f6341SShri Abhyankar       sum11 += v[10]*x1 + v[25]*x2 + v[40]*x3  + v[55]*x4;
7640b8f6341SShri Abhyankar       sum12 += v[11]*x1 + v[26]*x2 + v[41]*x3  + v[56]*x4;
7650b8f6341SShri Abhyankar       sum13 += v[12]*x1 + v[27]*x2 + v[42]*x3   + v[57]*x4;
7660b8f6341SShri Abhyankar       sum14 += v[13]*x1 + v[28]*x2 + v[43]*x3   + v[58]*x4;
7670b8f6341SShri Abhyankar       sum15 += v[14]*x1 + v[29]*x2 + v[44]*x3  + v[59]*x4;
7680b8f6341SShri Abhyankar       v  += 60;
7690b8f6341SShri Abhyankar 
7700b8f6341SShri Abhyankar       x1 = xb[12]; x2 = xb[13]; x3 = xb[14];
7718ab949d8SShri Abhyankar       sum1  += v[0]*x1 + v[15]*x2 + v[30]*x3;
7728ab949d8SShri Abhyankar       sum2  += v[1]*x1 + v[16]*x2 + v[31]*x3;
7738ab949d8SShri Abhyankar       sum3  += v[2]*x1 + v[17]*x2 + v[32]*x3;
7748ab949d8SShri Abhyankar       sum4  += v[3]*x1 + v[18]*x2 + v[33]*x3;
7758ab949d8SShri Abhyankar       sum5  += v[4]*x1 + v[19]*x2 + v[34]*x3;
7768ab949d8SShri Abhyankar       sum6  += v[5]*x1 + v[20]*x2 + v[35]*x3;
7778ab949d8SShri Abhyankar       sum7  += v[6]*x1 + v[21]*x2 + v[36]*x3;
7788ab949d8SShri Abhyankar       sum8  += v[7]*x1 + v[22]*x2 + v[37]*x3;
7798ab949d8SShri Abhyankar       sum9  += v[8]*x1 + v[23]*x2 + v[38]*x3;
7808ab949d8SShri Abhyankar       sum10 += v[9]*x1 + v[24]*x2 + v[39]*x3;
7818ab949d8SShri Abhyankar       sum11 += v[10]*x1 + v[25]*x2 + v[40]*x3;
7828ab949d8SShri Abhyankar       sum12 += v[11]*x1 + v[26]*x2 + v[41]*x3;
7838ab949d8SShri Abhyankar       sum13 += v[12]*x1 + v[27]*x2 + v[42]*x3;
7848ab949d8SShri Abhyankar       sum14 += v[13]*x1 + v[28]*x2 + v[43]*x3;
7858ab949d8SShri Abhyankar       sum15 += v[14]*x1 + v[29]*x2 + v[44]*x3;
7868ab949d8SShri Abhyankar       v += 45;
7878ab949d8SShri Abhyankar     }
7888ab949d8SShri Abhyankar     if (usecprow) z = zarray + 15*ridx[i];
7898ab949d8SShri Abhyankar     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
7908ab949d8SShri Abhyankar     z[7] = sum8; z[8] = sum9; z[9] = sum10; z[10] = sum11; z[11] = sum12; z[12] = sum13; z[13] = sum14;z[14] = sum15;
7918ab949d8SShri Abhyankar 
7928ab949d8SShri Abhyankar     if (!usecprow) z += 15;
7938ab949d8SShri Abhyankar   }
7948ab949d8SShri Abhyankar 
7958ab949d8SShri Abhyankar   ierr = VecRestoreArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
7968ab949d8SShri Abhyankar   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
7978ab949d8SShri Abhyankar   ierr = PetscLogFlops(450.0*a->nz - 15.0*nonzerorow);CHKERRQ(ierr);
7988ab949d8SShri Abhyankar   PetscFunctionReturn(0);
7998ab949d8SShri Abhyankar }
8008ab949d8SShri Abhyankar 
8018ab949d8SShri Abhyankar /* MatMult_SeqBAIJ_15_ver3 : Columns in the block are accessed in sets of 8,7 */
8028ab949d8SShri Abhyankar #undef __FUNCT__
8038ab949d8SShri Abhyankar #define __FUNCT__ "MatMult_SeqBAIJ_15_ver3"
8048ab949d8SShri Abhyankar PetscErrorCode MatMult_SeqBAIJ_15_ver3(Mat A,Vec xx,Vec zz)
8058ab949d8SShri Abhyankar {
8068ab949d8SShri Abhyankar   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
8078ab949d8SShri Abhyankar   PetscScalar       *z = 0,sum1,sum2,sum3,sum4,sum5,sum6,sum7,sum8,sum9,sum10,sum11,sum12,sum13,sum14,sum15;
8088ab949d8SShri Abhyankar   const PetscScalar *x,*xb;
8090b8f6341SShri Abhyankar   PetscScalar        x1,x2,x3,x4,x5,x6,x7,x8,*zarray;
8108ab949d8SShri Abhyankar   const MatScalar   *v;
8118ab949d8SShri Abhyankar   PetscErrorCode    ierr;
8128ab949d8SShri Abhyankar   const PetscInt    *ii,*ij=a->j,*idx;
8138ab949d8SShri Abhyankar   PetscInt          mbs,i,j,n,*ridx=PETSC_NULL,nonzerorow=0;
8148ab949d8SShri Abhyankar   PetscTruth        usecprow=a->compressedrow.use;
8158ab949d8SShri Abhyankar 
8168ab949d8SShri Abhyankar   PetscFunctionBegin;
8178ab949d8SShri Abhyankar   ierr = VecGetArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
8188ab949d8SShri Abhyankar   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
8198ab949d8SShri Abhyankar 
8208ab949d8SShri Abhyankar   v   = a->a;
8218ab949d8SShri Abhyankar   if (usecprow){
8228ab949d8SShri Abhyankar     mbs    = a->compressedrow.nrows;
8238ab949d8SShri Abhyankar     ii     = a->compressedrow.i;
8248ab949d8SShri Abhyankar     ridx = a->compressedrow.rindex;
8258ab949d8SShri Abhyankar   } else {
8268ab949d8SShri Abhyankar     mbs = a->mbs;
8278ab949d8SShri Abhyankar     ii  = a->i;
8288ab949d8SShri Abhyankar     z   = zarray;
8298ab949d8SShri Abhyankar   }
8308ab949d8SShri Abhyankar 
8318ab949d8SShri Abhyankar   for (i=0; i<mbs; i++) {
8328ab949d8SShri Abhyankar     n  = ii[i+1] - ii[i];
8338ab949d8SShri Abhyankar     idx = ij + ii[i];
8348ab949d8SShri Abhyankar     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0; sum7 = 0.0;
8358ab949d8SShri Abhyankar     sum8 = 0.0; sum9 = 0.0; sum10 = 0.0; sum11 = 0.0; sum12 = 0.0; sum13 = 0.0; sum14 = 0.0;sum15 = 0.0;
8368ab949d8SShri Abhyankar 
8378ab949d8SShri Abhyankar     nonzerorow += (n>0);
8388ab949d8SShri Abhyankar     for (j=0; j<n; j++) {
8398ab949d8SShri Abhyankar       xb = x + 15*(idx[j]);
8408ab949d8SShri Abhyankar       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6];
8410b8f6341SShri Abhyankar       x8 = xb[7];
8428ab949d8SShri Abhyankar 
8438ab949d8SShri Abhyankar       sum1  += v[0]*x1 + v[15]*x2  + v[30]*x3  + v[45]*x4 + v[60]*x5 + v[75]*x6 + v[90]*x7 + v[105]*x8;
8448ab949d8SShri Abhyankar       sum2  += v[1]*x1 + v[16]*x2  + v[31]*x3  + v[46]*x4 + v[61]*x5 + v[76]*x6 + v[91]*x7 + v[106]*x8;
8458ab949d8SShri Abhyankar       sum3  += v[2]*x1 + v[17]*x2  + v[32]*x3  + v[47]*x4 + v[62]*x5 + v[77]*x6 + v[92]*x7 + v[107]*x8;
8468ab949d8SShri Abhyankar       sum4  += v[3]*x1 + v[18]*x2 + v[33]*x3  + v[48]*x4 + v[63]*x5 + v[78]*x6 + v[93]*x7 + v[108]*x8;
8478ab949d8SShri Abhyankar       sum5  += v[4]*x1 + v[19]*x2 + v[34]*x3  + v[49]*x4 + v[64]*x5 + v[79]*x6 + v[94]*x7 + v[109]*x8;
8488ab949d8SShri Abhyankar       sum6  += v[5]*x1 + v[20]*x2 + v[35]*x3  + v[50]*x4 + v[65]*x5 + v[80]*x6 + v[95]*x7 + v[110]*x8;
8498ab949d8SShri Abhyankar       sum7  += v[6]*x1 + v[21]*x2 + v[36]*x3  + v[51]*x4 + v[66]*x5 + v[81]*x6 + v[96]*x7 + v[111]*x8;
8508ab949d8SShri Abhyankar       sum8  += v[7]*x1 + v[22]*x2  + v[37]*x3  + v[52]*x4 + v[67]*x5 + v[82]*x6 + v[97]*x7 + v[112]*x8;
8518ab949d8SShri Abhyankar       sum9  += v[8]*x1 + v[23]*x2  + v[38]*x3  + v[53]*x4 + v[68]*x5 + v[83]*x6 + v[98]*x7 + v[113]*x8;
8528ab949d8SShri Abhyankar       sum10 += v[9]*x1 + v[24]*x2  + v[39]*x3  + v[54]*x4 + v[69]*x5 + v[84]*x6 + v[99]*x7 + v[114]*x8;
8538ab949d8SShri Abhyankar       sum11 += v[10]*x1 + v[25]*x2 + v[40]*x3  + v[55]*x4 + v[70]*x5 + v[85]*x6 + v[100]*x7 + v[115]*x8;
8548ab949d8SShri Abhyankar       sum12 += v[11]*x1 + v[26]*x2 + v[41]*x3  + v[56]*x4 + v[71]*x5 + v[86]*x6 + v[101]*x7 + v[116]*x8;
8558ab949d8SShri Abhyankar       sum13 += v[12]*x1 + v[27]*x2 + v[42]*x3  + v[57]*x4 + v[72]*x5 + v[87]*x6 + v[102]*x7 + v[117]*x8;
8568ab949d8SShri Abhyankar       sum14 += v[13]*x1 + v[28]*x2 + v[43]*x3  + v[58]*x4 + v[73]*x5 + v[88]*x6 + v[103]*x7 + v[118]*x8;
8578ab949d8SShri Abhyankar       sum15 += v[14]*x1 + v[29]*x2 + v[44]*x3  + v[59]*x4 + v[74]*x5 + v[89]*x6 + v[104]*x7 + v[119]*x8;
8588ab949d8SShri Abhyankar       v += 120;
8598ab949d8SShri Abhyankar 
8600b8f6341SShri Abhyankar       x1 = xb[8]; x2 = xb[9]; x3 = xb[10]; x4 = xb[11]; x5 = xb[12]; x6 = xb[13]; x7 = xb[14];
8610b8f6341SShri Abhyankar 
8628ab949d8SShri Abhyankar       sum1  += v[0]*x1 + v[15]*x2  + v[30]*x3  + v[45]*x4 + v[60]*x5 + v[75]*x6 + v[90]*x7;
8638ab949d8SShri Abhyankar       sum2  += v[1]*x1 + v[16]*x2  + v[31]*x3  + v[46]*x4 + v[61]*x5 + v[76]*x6 + v[91]*x7;
8648ab949d8SShri Abhyankar       sum3  += v[2]*x1 + v[17]*x2  + v[32]*x3  + v[47]*x4 + v[62]*x5 + v[77]*x6 + v[92]*x7;
8658ab949d8SShri Abhyankar       sum4  += v[3]*x1 + v[18]*x2 + v[33]*x3  + v[48]*x4 + v[63]*x5 + v[78]*x6 + v[93]*x7;
8668ab949d8SShri Abhyankar       sum5  += v[4]*x1 + v[19]*x2 + v[34]*x3  + v[49]*x4 + v[64]*x5 + v[79]*x6 + v[94]*x7;
8678ab949d8SShri Abhyankar       sum6  += v[5]*x1 + v[20]*x2 + v[35]*x3  + v[50]*x4 + v[65]*x5 + v[80]*x6 + v[95]*x7;
8688ab949d8SShri Abhyankar       sum7  += v[6]*x1 + v[21]*x2 + v[36]*x3  + v[51]*x4 + v[66]*x5 + v[81]*x6 + v[96]*x7;
8698ab949d8SShri Abhyankar       sum8  += v[7]*x1 + v[22]*x2  + v[37]*x3  + v[52]*x4 + v[67]*x5 + v[82]*x6 + v[97]*x7;
8708ab949d8SShri Abhyankar       sum9  += v[8]*x1 + v[23]*x2  + v[38]*x3  + v[53]*x4 + v[68]*x5 + v[83]*x6 + v[98]*x7;
8718ab949d8SShri Abhyankar       sum10 += v[9]*x1 + v[24]*x2  + v[39]*x3  + v[54]*x4 + v[69]*x5 + v[84]*x6 + v[99]*x7;
8728ab949d8SShri Abhyankar       sum11 += v[10]*x1 + v[25]*x2 + v[40]*x3  + v[55]*x4 + v[70]*x5 + v[85]*x6 + v[100]*x7;
8738ab949d8SShri Abhyankar       sum12 += v[11]*x1 + v[26]*x2 + v[41]*x3  + v[56]*x4 + v[71]*x5 + v[86]*x6 + v[101]*x7;
8748ab949d8SShri Abhyankar       sum13 += v[12]*x1 + v[27]*x2 + v[42]*x3  + v[57]*x4 + v[72]*x5 + v[87]*x6 + v[102]*x7;
8758ab949d8SShri Abhyankar       sum14 += v[13]*x1 + v[28]*x2 + v[43]*x3  + v[58]*x4 + v[73]*x5 + v[88]*x6 + v[103]*x7;
8768ab949d8SShri Abhyankar       sum15 += v[14]*x1 + v[29]*x2 + v[44]*x3  + v[59]*x4 + v[74]*x5 + v[89]*x6 + v[104]*x7;
8778ab949d8SShri Abhyankar       v += 105;
8788ab949d8SShri Abhyankar     }
8798ab949d8SShri Abhyankar     if (usecprow) z = zarray + 15*ridx[i];
8808ab949d8SShri Abhyankar     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
8818ab949d8SShri Abhyankar     z[7] = sum8; z[8] = sum9; z[9] = sum10; z[10] = sum11; z[11] = sum12; z[12] = sum13; z[13] = sum14;z[14] = sum15;
8828ab949d8SShri Abhyankar 
8838ab949d8SShri Abhyankar     if (!usecprow) z += 15;
8848ab949d8SShri Abhyankar   }
8858ab949d8SShri Abhyankar 
8868ab949d8SShri Abhyankar   ierr = VecRestoreArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
8878ab949d8SShri Abhyankar   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
8888ab949d8SShri Abhyankar   ierr = PetscLogFlops(450.0*a->nz - 15.0*nonzerorow);CHKERRQ(ierr);
8898ab949d8SShri Abhyankar   PetscFunctionReturn(0);
8908ab949d8SShri Abhyankar }
8918ab949d8SShri Abhyankar 
8928ab949d8SShri Abhyankar /* MatMult_SeqBAIJ_15_ver4 : All columns in the block are accessed at once */
8938ab949d8SShri Abhyankar 
8948ab949d8SShri Abhyankar #undef __FUNCT__
8958ab949d8SShri Abhyankar #define __FUNCT__ "MatMult_SeqBAIJ_15_ver4"
8968ab949d8SShri Abhyankar PetscErrorCode MatMult_SeqBAIJ_15_ver4(Mat A,Vec xx,Vec zz)
8978ab949d8SShri Abhyankar {
8988ab949d8SShri Abhyankar   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
8998ab949d8SShri Abhyankar   PetscScalar       *z = 0,sum1,sum2,sum3,sum4,sum5,sum6,sum7,sum8,sum9,sum10,sum11,sum12,sum13,sum14,sum15;
9008ab949d8SShri Abhyankar   const PetscScalar *x,*xb;
9018ab949d8SShri Abhyankar   PetscScalar        x1,x2,x3,x4,x5,x6,x7,x8,x9,x10,x11,x12,x13,x14,x15,*zarray;
9028ab949d8SShri Abhyankar   const MatScalar   *v;
9038ab949d8SShri Abhyankar   PetscErrorCode    ierr;
9048ab949d8SShri Abhyankar   const PetscInt    *ii,*ij=a->j,*idx;
9058ab949d8SShri Abhyankar   PetscInt          mbs,i,j,n,*ridx=PETSC_NULL,nonzerorow=0;
9068ab949d8SShri Abhyankar   PetscTruth        usecprow=a->compressedrow.use;
9078ab949d8SShri Abhyankar 
9088ab949d8SShri Abhyankar   PetscFunctionBegin;
9098ab949d8SShri Abhyankar   ierr = VecGetArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
9108ab949d8SShri Abhyankar   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
9118ab949d8SShri Abhyankar 
9128ab949d8SShri Abhyankar   v   = a->a;
9138ab949d8SShri Abhyankar   if (usecprow){
9148ab949d8SShri Abhyankar     mbs    = a->compressedrow.nrows;
9158ab949d8SShri Abhyankar     ii     = a->compressedrow.i;
9168ab949d8SShri Abhyankar     ridx = a->compressedrow.rindex;
9178ab949d8SShri Abhyankar   } else {
9188ab949d8SShri Abhyankar     mbs = a->mbs;
9198ab949d8SShri Abhyankar     ii  = a->i;
9208ab949d8SShri Abhyankar     z   = zarray;
9218ab949d8SShri Abhyankar   }
9228ab949d8SShri Abhyankar 
9238ab949d8SShri Abhyankar   for (i=0; i<mbs; i++) {
9248ab949d8SShri Abhyankar     n  = ii[i+1] - ii[i];
9258ab949d8SShri Abhyankar     idx = ij + ii[i];
9268ab949d8SShri Abhyankar     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0; sum7 = 0.0;
9278ab949d8SShri Abhyankar     sum8 = 0.0; sum9 = 0.0; sum10 = 0.0; sum11 = 0.0; sum12 = 0.0; sum13 = 0.0; sum14 = 0.0;sum15 = 0.0;
9288ab949d8SShri Abhyankar 
9298ab949d8SShri Abhyankar     nonzerorow += (n>0);
9308ab949d8SShri Abhyankar     for (j=0; j<n; j++) {
9318ab949d8SShri Abhyankar       xb = x + 15*(idx[j]);
9328ab949d8SShri Abhyankar       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6];
9338ab949d8SShri Abhyankar       x8 = xb[7]; x9 = xb[8]; x10 = xb[9]; x11 = xb[10]; x12 = xb[11]; x13 = xb[12]; x14 = xb[13];x15 = xb[14];
9348ab949d8SShri Abhyankar 
9358ab949d8SShri Abhyankar       sum1 +=  v[0]*x1  + v[15]*x2 + v[30]*x3 + v[45]*x4 + v[60]*x5 + v[75]*x6 + v[90]*x7  + v[105]*x8 + v[120]*x9 + v[135]*x10 + v[150]*x11 + v[165]*x12 + v[180]*x13 + v[195]*x14 + v[210]*x15;
9368ab949d8SShri Abhyankar       sum2 +=  v[1]*x1  + v[16]*x2 + v[31]*x3 + v[46]*x4 + v[61]*x5 + v[76]*x6 + v[91]*x7  + v[106]*x8 + v[121]*x9 + v[136]*x10 + v[151]*x11 + v[166]*x12 + v[181]*x13 + v[196]*x14 + v[211]*x15;
9378ab949d8SShri Abhyankar       sum3 +=  v[2]*x1  + v[17]*x2 + v[32]*x3 + v[47]*x4 + v[62]*x5 + v[77]*x6 + v[92]*x7  + v[107]*x8 + v[122]*x9 + v[137]*x10 + v[152]*x11 + v[167]*x12 + v[182]*x13 + v[197]*x14 + v[212]*x15;
9388ab949d8SShri Abhyankar       sum4 +=  v[3]*x1  + v[18]*x2 + v[33]*x3 + v[48]*x4 + v[63]*x5 + v[78]*x6 + v[93]*x7  + v[108]*x8 + v[123]*x9 + v[138]*x10 + v[153]*x11 + v[168]*x12 + v[183]*x13 + v[198]*x14 + v[213]*x15;
9398ab949d8SShri Abhyankar       sum5  += v[4]*x1  + v[19]*x2 + v[34]*x3 + v[49]*x4 + v[64]*x5 + v[79]*x6 + v[94]*x7  + v[109]*x8 + v[124]*x9 + v[139]*x10 + v[154]*x11 + v[169]*x12 + v[184]*x13 + v[199]*x14 + v[214]*x15;
9408ab949d8SShri Abhyankar       sum6  += v[5]*x1  + v[20]*x2 + v[35]*x3 + v[50]*x4 + v[65]*x5 + v[80]*x6 + v[95]*x7  + v[110]*x8 + v[125]*x9 + v[140]*x10 + v[155]*x11 + v[170]*x12 + v[185]*x13 + v[200]*x14 + v[215]*x15;
9418ab949d8SShri Abhyankar       sum7  += v[6]*x1  + v[21]*x2 + v[36]*x3 + v[51]*x4 + v[66]*x5 + v[81]*x6 + v[96]*x7  + v[111]*x8 + v[126]*x9 + v[141]*x10 + v[156]*x11 + v[171]*x12 + v[186]*x13 + v[201]*x14 + v[216]*x15;
9428ab949d8SShri Abhyankar       sum8  += v[7]*x1  + v[22]*x2 + v[37]*x3 + v[52]*x4 + v[67]*x5 + v[82]*x6 + v[97]*x7  + v[112]*x8 + v[127]*x9 + v[142]*x10 + v[157]*x11 + v[172]*x12 + v[187]*x13 + v[202]*x14 + v[217]*x15;
9438ab949d8SShri Abhyankar       sum9  += v[8]*x1  + v[23]*x2 + v[38]*x3 + v[53]*x4 + v[68]*x5 + v[83]*x6 + v[98]*x7  + v[113]*x8 + v[128]*x9 + v[143]*x10 + v[158]*x11 + v[173]*x12 + v[188]*x13 + v[203]*x14 + v[218]*x15;
9448ab949d8SShri Abhyankar       sum10 += v[9]*x1  + v[24]*x2 + v[39]*x3 + v[54]*x4 + v[69]*x5 + v[84]*x6 + v[99]*x7  + v[114]*x8 + v[129]*x9 + v[144]*x10 + v[159]*x11 + v[174]*x12 + v[189]*x13 + v[204]*x14 + v[219]*x15;
9458ab949d8SShri Abhyankar       sum11 += v[10]*x1 + v[25]*x2 + v[40]*x3 + v[55]*x4 + v[70]*x5 + v[85]*x6 + v[100]*x7 + v[115]*x8 + v[130]*x9 + v[145]*x10 + v[160]*x11 + v[175]*x12 + v[190]*x13 + v[205]*x14 + v[220]*x15;
9468ab949d8SShri Abhyankar       sum12 += v[11]*x1 + v[26]*x2 + v[41]*x3 + v[56]*x4 + v[71]*x5 + v[86]*x6 + v[101]*x7 + v[116]*x8 + v[131]*x9 + v[146]*x10 + v[161]*x11 + v[176]*x12 + v[191]*x13 + v[206]*x14 + v[221]*x15;
9478ab949d8SShri Abhyankar       sum13 += v[12]*x1 + v[27]*x2 + v[42]*x3 + v[57]*x4 + v[72]*x5 + v[87]*x6 + v[102]*x7 + v[117]*x8 + v[132]*x9 + v[147]*x10 + v[162]*x11 + v[177]*x12 + v[192]*x13 + v[207]*x14 + v[222]*x15;
9488ab949d8SShri Abhyankar       sum14 += v[13]*x1 + v[28]*x2 + v[43]*x3 + v[58]*x4 + v[73]*x5 + v[88]*x6 + v[103]*x7 + v[118]*x8 + v[133]*x9 + v[148]*x10 + v[163]*x11 + v[178]*x12 + v[193]*x13 + v[208]*x14 + v[223]*x15;
9498ab949d8SShri Abhyankar       sum15 += v[14]*x1 + v[29]*x2 + v[44]*x3 + v[59]*x4 + v[74]*x5 + v[89]*x6 + v[104]*x7 + v[119]*x8 + v[134]*x9 + v[149]*x10 + v[164]*x11 + v[179]*x12 + v[194]*x13 + v[209]*x14 + v[224]*x15;
9508ab949d8SShri Abhyankar       v += 225;
9518ab949d8SShri Abhyankar     }
9528ab949d8SShri Abhyankar     if (usecprow) z = zarray + 15*ridx[i];
9538ab949d8SShri Abhyankar     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
9548ab949d8SShri Abhyankar     z[7] = sum8; z[8] = sum9; z[9] = sum10; z[10] = sum11; z[11] = sum12; z[12] = sum13; z[13] = sum14;z[14] = sum15;
9558ab949d8SShri Abhyankar 
9568ab949d8SShri Abhyankar     if (!usecprow) z += 15;
9578ab949d8SShri Abhyankar   }
9588ab949d8SShri Abhyankar 
9598ab949d8SShri Abhyankar   ierr = VecRestoreArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
9608ab949d8SShri Abhyankar   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
9618ab949d8SShri Abhyankar   ierr = PetscLogFlops(450.0*a->nz - 15.0*nonzerorow);CHKERRQ(ierr);
9628ab949d8SShri Abhyankar   PetscFunctionReturn(0);
9638ab949d8SShri Abhyankar }
9648ab949d8SShri Abhyankar 
9658ab949d8SShri Abhyankar 
9663f1db9ecSBarry Smith /*
9673f1db9ecSBarry Smith     This will not work with MatScalar == float because it calls the BLAS
9683f1db9ecSBarry Smith */
9694a2ae208SSatish Balay #undef __FUNCT__
9704a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_N"
971dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_N(Mat A,Vec xx,Vec zz)
9722d61bbb3SSatish Balay {
9732d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
974dcf5cc72SBarry Smith   PetscScalar    *x,*z = 0,*xb,*work,*workt,*zarray;
9753f1db9ecSBarry Smith   MatScalar      *v;
976dfbe8321SBarry Smith   PetscErrorCode ierr;
977d0f46423SBarry Smith   PetscInt       mbs=a->mbs,i,*idx,*ii,bs=A->rmap->bs,j,n,bs2=a->bs2;
97898c9bda7SSatish Balay   PetscInt       ncols,k,*ridx=PETSC_NULL,nonzerorow=0;
97926e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
9802d61bbb3SSatish Balay 
9812d61bbb3SSatish Balay   PetscFunctionBegin;
9821ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
98326e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
9842d61bbb3SSatish Balay 
9852d61bbb3SSatish Balay   idx = a->j;
9862d61bbb3SSatish Balay   v   = a->a;
98726e093fcSHong Zhang   if (usecprow){
98826e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
98926e093fcSHong Zhang     ii   = a->compressedrow.i;
9907b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
99126e093fcSHong Zhang   } else {
99226e093fcSHong Zhang     mbs = a->mbs;
9932d61bbb3SSatish Balay     ii  = a->i;
99426e093fcSHong Zhang     z   = zarray;
99526e093fcSHong Zhang   }
996218c64b6SSatish Balay 
9972d61bbb3SSatish Balay   if (!a->mult_work) {
998d0f46423SBarry Smith     k    = PetscMax(A->rmap->n,A->cmap->n);
99987828ca2SBarry Smith     ierr = PetscMalloc((k+1)*sizeof(PetscScalar),&a->mult_work);CHKERRQ(ierr);
10002d61bbb3SSatish Balay   }
10012d61bbb3SSatish Balay   work = a->mult_work;
10022d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
10032d61bbb3SSatish Balay     n     = ii[1] - ii[0]; ii++;
10042d61bbb3SSatish Balay     ncols = n*bs;
10052d61bbb3SSatish Balay     workt = work;
100698c9bda7SSatish Balay     nonzerorow += (n>0);
10072d61bbb3SSatish Balay     for (j=0; j<n; j++) {
10082d61bbb3SSatish Balay       xb = x + bs*(*idx++);
10092d61bbb3SSatish Balay       for (k=0; k<bs; k++) workt[k] = xb[k];
10102d61bbb3SSatish Balay       workt += bs;
10112d61bbb3SSatish Balay     }
10127b2bb3b9SHong Zhang     if (usecprow) z = zarray + bs*ridx[i];
1013737d121aSSatish Balay     Kernel_w_gets_Ar_times_v(bs,ncols,work,v,z);
101471044d3cSBarry Smith     /* BLASgemv_("N",&bs,&ncols,&_DOne,v,&bs,work,&_One,&_DZero,z,&_One); */
10152d61bbb3SSatish Balay     v += n*bs2;
101626e093fcSHong Zhang     if (!usecprow) z += bs;
10172d61bbb3SSatish Balay   }
10181ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
101926e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
1020dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*a->nz*bs2 - bs*nonzerorow);CHKERRQ(ierr);
10212d61bbb3SSatish Balay   PetscFunctionReturn(0);
10222d61bbb3SSatish Balay }
10232d61bbb3SSatish Balay 
10246a65c766SHong Zhang extern PetscErrorCode VecCopy_Seq(Vec,Vec);
10254a2ae208SSatish Balay #undef __FUNCT__
10264a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_1"
1027dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_1(Mat A,Vec xx,Vec yy,Vec zz)
10282d61bbb3SSatish Balay {
10292d61bbb3SSatish Balay   Mat_SeqBAIJ        *a = (Mat_SeqBAIJ*)A->data;
1030122f12eaSBarry Smith   const PetscScalar  *x;
1031122f12eaSBarry Smith   PetscScalar        *y,*z,sum;
1032122f12eaSBarry Smith   const MatScalar    *v;
1033dfbe8321SBarry Smith   PetscErrorCode     ierr;
1034122f12eaSBarry Smith   PetscInt           mbs=a->mbs,i,n,*ridx=PETSC_NULL,nonzerorow=0;
1035122f12eaSBarry Smith   const PetscInt     *idx,*ii;
103626e093fcSHong Zhang   PetscTruth         usecprow=a->compressedrow.use;
10372d61bbb3SSatish Balay 
10382d61bbb3SSatish Balay   PetscFunctionBegin;
1039122f12eaSBarry Smith   ierr = VecGetArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
1040122f12eaSBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
10412e8a6d31SBarry Smith   if (zz != yy) {
1042122f12eaSBarry Smith     ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
10432e8a6d31SBarry Smith   } else {
1044122f12eaSBarry Smith     z = y;
10452e8a6d31SBarry Smith   }
10462d61bbb3SSatish Balay 
10472d61bbb3SSatish Balay   idx = a->j;
10482d61bbb3SSatish Balay   v   = a->a;
104926e093fcSHong Zhang   if (usecprow){
10504eb6d288SHong Zhang     if (zz != yy){
1051122f12eaSBarry Smith       ierr = PetscMemcpy(z,y,mbs*sizeof(PetscScalar));CHKERRQ(ierr);
10524eb6d288SHong Zhang     }
105326e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
105426e093fcSHong Zhang     ii   = a->compressedrow.i;
10557b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
105626e093fcSHong Zhang   } else {
10572d61bbb3SSatish Balay     ii  = a->i;
105826e093fcSHong Zhang   }
10592d61bbb3SSatish Balay 
10602d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
1061122f12eaSBarry Smith     n    = ii[1] - ii[0];
1062122f12eaSBarry Smith     ii++;
106326e093fcSHong Zhang     if (!usecprow){
1064122f12eaSBarry Smith       nonzerorow += (n>0);
1065122f12eaSBarry Smith       sum = y[i];
1066122f12eaSBarry Smith     } else {
1067122f12eaSBarry Smith       sum = y[ridx[i]];
1068122f12eaSBarry Smith     }
1069122f12eaSBarry Smith     PetscSparseDensePlusDot(sum,x,v,idx,n);
1070122f12eaSBarry Smith     v += n;
1071122f12eaSBarry Smith     idx += n;
1072122f12eaSBarry Smith     if (usecprow){
1073122f12eaSBarry Smith       z[ridx[i]] = sum;
1074122f12eaSBarry Smith     } else {
1075122f12eaSBarry Smith       z[i] = sum;
107626e093fcSHong Zhang     }
10772d61bbb3SSatish Balay   }
1078122f12eaSBarry Smith   ierr = VecRestoreArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
1079122f12eaSBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
10802e8a6d31SBarry Smith   if (zz != yy) {
1081122f12eaSBarry Smith     ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
10822e8a6d31SBarry Smith   }
1083122f12eaSBarry Smith   ierr = PetscLogFlops(2.0*a->nz - nonzerorow);CHKERRQ(ierr);
10842d61bbb3SSatish Balay   PetscFunctionReturn(0);
10852d61bbb3SSatish Balay }
10862d61bbb3SSatish Balay 
10874a2ae208SSatish Balay #undef __FUNCT__
10884a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_2"
1089dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_2(Mat A,Vec xx,Vec yy,Vec zz)
10902d61bbb3SSatish Balay {
10912d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
1092dcf5cc72SBarry Smith   PetscScalar    *x,*y = 0,*z = 0,*xb,sum1,sum2;
109326e093fcSHong Zhang   PetscScalar    x1,x2,*yarray,*zarray;
10943f1db9ecSBarry Smith   MatScalar      *v;
1095dfbe8321SBarry Smith   PetscErrorCode ierr;
10964eb6d288SHong Zhang   PetscInt       mbs=a->mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL;
109726e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
10982d61bbb3SSatish Balay 
10992d61bbb3SSatish Balay   PetscFunctionBegin;
11001ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
110126e093fcSHong Zhang   ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr);
11022e8a6d31SBarry Smith   if (zz != yy) {
110326e093fcSHong Zhang     ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
11042e8a6d31SBarry Smith   } else {
110526e093fcSHong Zhang     zarray = yarray;
11062e8a6d31SBarry Smith   }
11072d61bbb3SSatish Balay 
11082d61bbb3SSatish Balay   idx = a->j;
11092d61bbb3SSatish Balay   v   = a->a;
111026e093fcSHong Zhang   if (usecprow){
11114eb6d288SHong Zhang     if (zz != yy){
11124eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,2*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
11134eb6d288SHong Zhang     }
111426e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
111526e093fcSHong Zhang     ii   = a->compressedrow.i;
11167b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
11174eb6d288SHong Zhang     if (zz != yy){
11184eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,a->mbs*sizeof(PetscScalar));CHKERRQ(ierr);
11194eb6d288SHong Zhang     }
112026e093fcSHong Zhang   } else {
11212d61bbb3SSatish Balay     ii  = a->i;
112226e093fcSHong Zhang     y   = yarray;
112326e093fcSHong Zhang     z   = zarray;
112426e093fcSHong Zhang   }
11252d61bbb3SSatish Balay 
11262d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
11272d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
112826e093fcSHong Zhang     if (usecprow){
11297b2bb3b9SHong Zhang       z = zarray + 2*ridx[i];
11307b2bb3b9SHong Zhang       y = yarray + 2*ridx[i];
113126e093fcSHong Zhang     }
11322d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1];
11332d61bbb3SSatish Balay     for (j=0; j<n; j++) {
11342d61bbb3SSatish Balay       xb = x + 2*(*idx++); x1 = xb[0]; x2 = xb[1];
11352d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[2]*x2;
11362d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[3]*x2;
11372d61bbb3SSatish Balay       v += 4;
11382d61bbb3SSatish Balay     }
11392d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2;
114026e093fcSHong Zhang     if (!usecprow){
11412d61bbb3SSatish Balay       z += 2; y += 2;
11422d61bbb3SSatish Balay     }
114326e093fcSHong Zhang   }
11441ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
114526e093fcSHong Zhang   ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr);
11462e8a6d31SBarry Smith   if (zz != yy) {
114726e093fcSHong Zhang     ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
11482e8a6d31SBarry Smith   }
1149dc0b31edSSatish Balay   ierr = PetscLogFlops(4.0*a->nz);CHKERRQ(ierr);
11502d61bbb3SSatish Balay   PetscFunctionReturn(0);
11512d61bbb3SSatish Balay }
11522d61bbb3SSatish Balay 
11534a2ae208SSatish Balay #undef __FUNCT__
11544a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_3"
1155dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_3(Mat A,Vec xx,Vec yy,Vec zz)
11562d61bbb3SSatish Balay {
11572d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
1158dcf5cc72SBarry Smith   PetscScalar    *x,*y = 0,*z = 0,*xb,sum1,sum2,sum3,x1,x2,x3,*yarray,*zarray;
11593f1db9ecSBarry Smith   MatScalar      *v;
1160dfbe8321SBarry Smith   PetscErrorCode ierr;
11614eb6d288SHong Zhang   PetscInt       mbs=a->mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL;
116226e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
11632d61bbb3SSatish Balay 
11642d61bbb3SSatish Balay   PetscFunctionBegin;
11651ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
116626e093fcSHong Zhang   ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr);
11672e8a6d31SBarry Smith   if (zz != yy) {
116826e093fcSHong Zhang     ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
11692e8a6d31SBarry Smith   } else {
117026e093fcSHong Zhang     zarray = yarray;
11712e8a6d31SBarry Smith   }
11722d61bbb3SSatish Balay 
11732d61bbb3SSatish Balay   idx = a->j;
11742d61bbb3SSatish Balay   v   = a->a;
117526e093fcSHong Zhang   if (usecprow){
11764eb6d288SHong Zhang     if (zz != yy){
11774eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,3*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
11784eb6d288SHong Zhang     }
117926e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
118026e093fcSHong Zhang     ii   = a->compressedrow.i;
11817b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
118226e093fcSHong Zhang   } else {
11832d61bbb3SSatish Balay     ii  = a->i;
118426e093fcSHong Zhang     y   = yarray;
118526e093fcSHong Zhang     z   = zarray;
118626e093fcSHong Zhang   }
11872d61bbb3SSatish Balay 
11882d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
11892d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
119026e093fcSHong Zhang     if (usecprow){
11917b2bb3b9SHong Zhang       z = zarray + 3*ridx[i];
11927b2bb3b9SHong Zhang       y = yarray + 3*ridx[i];
119326e093fcSHong Zhang     }
11942d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2];
11952d61bbb3SSatish Balay     for (j=0; j<n; j++) {
11962d61bbb3SSatish Balay       xb = x + 3*(*idx++); x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
11972d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[3]*x2 + v[6]*x3;
11982d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[4]*x2 + v[7]*x3;
11992d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[5]*x2 + v[8]*x3;
12002d61bbb3SSatish Balay       v += 9;
12012d61bbb3SSatish Balay     }
12022d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3;
120326e093fcSHong Zhang     if (!usecprow){
12042d61bbb3SSatish Balay       z += 3; y += 3;
12052d61bbb3SSatish Balay     }
120626e093fcSHong Zhang   }
12071ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
120826e093fcSHong Zhang   ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr);
12092e8a6d31SBarry Smith   if (zz != yy) {
121026e093fcSHong Zhang     ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
12112e8a6d31SBarry Smith   }
1212dc0b31edSSatish Balay   ierr = PetscLogFlops(18.0*a->nz);CHKERRQ(ierr);
12132d61bbb3SSatish Balay   PetscFunctionReturn(0);
12142d61bbb3SSatish Balay }
12152d61bbb3SSatish Balay 
12164a2ae208SSatish Balay #undef __FUNCT__
12174a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_4"
1218dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_4(Mat A,Vec xx,Vec yy,Vec zz)
12192d61bbb3SSatish Balay {
12202d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
1221dcf5cc72SBarry Smith   PetscScalar    *x,*y = 0,*z = 0,*xb,sum1,sum2,sum3,sum4,x1,x2,x3,x4,*yarray,*zarray;
12223f1db9ecSBarry Smith   MatScalar      *v;
1223dfbe8321SBarry Smith   PetscErrorCode ierr;
12244eb6d288SHong Zhang   PetscInt       mbs=a->mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL;
122526e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
12262d61bbb3SSatish Balay 
12272d61bbb3SSatish Balay   PetscFunctionBegin;
12281ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
122926e093fcSHong Zhang   ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr);
12302e8a6d31SBarry Smith   if (zz != yy) {
123126e093fcSHong Zhang     ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
12322e8a6d31SBarry Smith   } else {
123326e093fcSHong Zhang     zarray = yarray;
12342e8a6d31SBarry Smith   }
12352d61bbb3SSatish Balay 
12362d61bbb3SSatish Balay   idx   = a->j;
12372d61bbb3SSatish Balay   v     = a->a;
123826e093fcSHong Zhang   if (usecprow){
12394eb6d288SHong Zhang     if (zz != yy){
12404eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,4*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
12414eb6d288SHong Zhang     }
124226e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
124326e093fcSHong Zhang     ii   = a->compressedrow.i;
12447b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
124526e093fcSHong Zhang   } else {
12462d61bbb3SSatish Balay     ii  = a->i;
124726e093fcSHong Zhang     y   = yarray;
124826e093fcSHong Zhang     z   = zarray;
124926e093fcSHong Zhang   }
12502d61bbb3SSatish Balay 
12512d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
12522d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
125326e093fcSHong Zhang     if (usecprow){
12547b2bb3b9SHong Zhang       z = zarray + 4*ridx[i];
12557b2bb3b9SHong Zhang       y = yarray + 4*ridx[i];
125626e093fcSHong Zhang     }
12572d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3];
12582d61bbb3SSatish Balay     for (j=0; j<n; j++) {
12592d61bbb3SSatish Balay       xb = x + 4*(*idx++);
12602d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
12612d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[4]*x2 + v[8]*x3   + v[12]*x4;
12622d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[5]*x2 + v[9]*x3   + v[13]*x4;
12632d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[6]*x2 + v[10]*x3  + v[14]*x4;
12642d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[7]*x2 + v[11]*x3  + v[15]*x4;
12652d61bbb3SSatish Balay       v += 16;
12662d61bbb3SSatish Balay     }
12672d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4;
126826e093fcSHong Zhang     if (!usecprow){
12692d61bbb3SSatish Balay       z += 4; y += 4;
12702d61bbb3SSatish Balay     }
127126e093fcSHong Zhang   }
12721ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
127326e093fcSHong Zhang   ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr);
12742e8a6d31SBarry Smith   if (zz != yy) {
127526e093fcSHong Zhang     ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
12762e8a6d31SBarry Smith   }
1277dc0b31edSSatish Balay   ierr = PetscLogFlops(32.0*a->nz);CHKERRQ(ierr);
12782d61bbb3SSatish Balay   PetscFunctionReturn(0);
12792d61bbb3SSatish Balay }
12802d61bbb3SSatish Balay 
12814a2ae208SSatish Balay #undef __FUNCT__
12824a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_5"
1283dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_5(Mat A,Vec xx,Vec yy,Vec zz)
12842d61bbb3SSatish Balay {
12852d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
1286dcf5cc72SBarry Smith   PetscScalar    *x,*y = 0,*z = 0,*xb,sum1,sum2,sum3,sum4,sum5,x1,x2,x3,x4,x5;
128726e093fcSHong Zhang   PetscScalar    *yarray,*zarray;
12883f1db9ecSBarry Smith   MatScalar      *v;
1289dfbe8321SBarry Smith   PetscErrorCode ierr;
12904eb6d288SHong Zhang   PetscInt       mbs=a->mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL;
129126e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
12922d61bbb3SSatish Balay 
12932d61bbb3SSatish Balay   PetscFunctionBegin;
12941ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
129526e093fcSHong Zhang   ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr);
12962e8a6d31SBarry Smith   if (zz != yy) {
129726e093fcSHong Zhang     ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
12982e8a6d31SBarry Smith   } else {
129926e093fcSHong Zhang     zarray = yarray;
13002e8a6d31SBarry Smith   }
13012d61bbb3SSatish Balay 
13022d61bbb3SSatish Balay   idx = a->j;
13032d61bbb3SSatish Balay   v   = a->a;
130426e093fcSHong Zhang   if (usecprow){
13054eb6d288SHong Zhang     if (zz != yy){
13064eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,5*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
13074eb6d288SHong Zhang     }
130826e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
130926e093fcSHong Zhang     ii   = a->compressedrow.i;
13107b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
131126e093fcSHong Zhang   } else {
13122d61bbb3SSatish Balay     ii  = a->i;
131326e093fcSHong Zhang     y   = yarray;
131426e093fcSHong Zhang     z   = zarray;
131526e093fcSHong Zhang   }
13162d61bbb3SSatish Balay 
13172d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
13182d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
131926e093fcSHong Zhang     if (usecprow){
13207b2bb3b9SHong Zhang       z = zarray + 5*ridx[i];
13217b2bb3b9SHong Zhang       y = yarray + 5*ridx[i];
132226e093fcSHong Zhang     }
13232d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4];
13242d61bbb3SSatish Balay     for (j=0; j<n; j++) {
13252d61bbb3SSatish Balay       xb = x + 5*(*idx++);
13262d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4];
13272d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[5]*x2 + v[10]*x3  + v[15]*x4 + v[20]*x5;
13282d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[6]*x2 + v[11]*x3  + v[16]*x4 + v[21]*x5;
13292d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[7]*x2 + v[12]*x3  + v[17]*x4 + v[22]*x5;
13302d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[8]*x2 + v[13]*x3  + v[18]*x4 + v[23]*x5;
13312d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[9]*x2 + v[14]*x3  + v[19]*x4 + v[24]*x5;
13322d61bbb3SSatish Balay       v += 25;
13332d61bbb3SSatish Balay     }
13342d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5;
133526e093fcSHong Zhang     if (!usecprow){
13362d61bbb3SSatish Balay       z += 5; y += 5;
13372d61bbb3SSatish Balay     }
133826e093fcSHong Zhang   }
13391ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
134026e093fcSHong Zhang   ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr);
13412e8a6d31SBarry Smith   if (zz != yy) {
134226e093fcSHong Zhang     ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
13432e8a6d31SBarry Smith   }
1344dc0b31edSSatish Balay   ierr = PetscLogFlops(50.0*a->nz);CHKERRQ(ierr);
13452d61bbb3SSatish Balay   PetscFunctionReturn(0);
13462d61bbb3SSatish Balay }
13474a2ae208SSatish Balay #undef __FUNCT__
13484a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_6"
1349dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_6(Mat A,Vec xx,Vec yy,Vec zz)
135015091d37SBarry Smith {
135115091d37SBarry Smith   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
1352dcf5cc72SBarry Smith   PetscScalar    *x,*y = 0,*z = 0,*xb,sum1,sum2,sum3,sum4,sum5,sum6;
135326e093fcSHong Zhang   PetscScalar    x1,x2,x3,x4,x5,x6,*yarray,*zarray;
135415091d37SBarry Smith   MatScalar      *v;
1355dfbe8321SBarry Smith   PetscErrorCode ierr;
13564eb6d288SHong Zhang   PetscInt       mbs=a->mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL;
135726e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
135815091d37SBarry Smith 
135915091d37SBarry Smith   PetscFunctionBegin;
13601ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
136126e093fcSHong Zhang   ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr);
136215091d37SBarry Smith   if (zz != yy) {
136326e093fcSHong Zhang     ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
136415091d37SBarry Smith   } else {
136526e093fcSHong Zhang     zarray = yarray;
136615091d37SBarry Smith   }
136715091d37SBarry Smith 
136815091d37SBarry Smith   idx = a->j;
136915091d37SBarry Smith   v   = a->a;
137026e093fcSHong Zhang   if (usecprow){
13714eb6d288SHong Zhang     if (zz != yy){
13724eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,6*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
13734eb6d288SHong Zhang     }
137426e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
137526e093fcSHong Zhang     ii   = a->compressedrow.i;
13767b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
137726e093fcSHong Zhang   } else {
137815091d37SBarry Smith     ii  = a->i;
137926e093fcSHong Zhang     y   = yarray;
138026e093fcSHong Zhang     z   = zarray;
138126e093fcSHong Zhang   }
138215091d37SBarry Smith 
138315091d37SBarry Smith   for (i=0; i<mbs; i++) {
138415091d37SBarry Smith     n  = ii[1] - ii[0]; ii++;
138526e093fcSHong Zhang     if (usecprow){
13867b2bb3b9SHong Zhang       z = zarray + 6*ridx[i];
13877b2bb3b9SHong Zhang       y = yarray + 6*ridx[i];
138826e093fcSHong Zhang     }
138915091d37SBarry Smith     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4]; sum6 = y[5];
139015091d37SBarry Smith     for (j=0; j<n; j++) {
13913b95cb0eSSatish Balay       xb = x + 6*(*idx++);
139215091d37SBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5];
139315091d37SBarry Smith       sum1 += v[0]*x1 + v[6]*x2  + v[12]*x3  + v[18]*x4 + v[24]*x5 + v[30]*x6;
139415091d37SBarry Smith       sum2 += v[1]*x1 + v[7]*x2  + v[13]*x3  + v[19]*x4 + v[25]*x5 + v[31]*x6;
139515091d37SBarry Smith       sum3 += v[2]*x1 + v[8]*x2  + v[14]*x3  + v[20]*x4 + v[26]*x5 + v[32]*x6;
139615091d37SBarry Smith       sum4 += v[3]*x1 + v[9]*x2  + v[15]*x3  + v[21]*x4 + v[27]*x5 + v[33]*x6;
139715091d37SBarry Smith       sum5 += v[4]*x1 + v[10]*x2 + v[16]*x3  + v[22]*x4 + v[28]*x5 + v[34]*x6;
139815091d37SBarry Smith       sum6 += v[5]*x1 + v[11]*x2 + v[17]*x3  + v[23]*x4 + v[29]*x5 + v[35]*x6;
139915091d37SBarry Smith       v += 36;
140015091d37SBarry Smith     }
140115091d37SBarry Smith     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6;
140226e093fcSHong Zhang     if (!usecprow){
140315091d37SBarry Smith       z += 6; y += 6;
140415091d37SBarry Smith     }
140526e093fcSHong Zhang   }
14061ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
140726e093fcSHong Zhang   ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr);
140815091d37SBarry Smith   if (zz != yy) {
140926e093fcSHong Zhang     ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
141015091d37SBarry Smith   }
1411dc0b31edSSatish Balay   ierr = PetscLogFlops(72.0*a->nz);CHKERRQ(ierr);
141215091d37SBarry Smith   PetscFunctionReturn(0);
141315091d37SBarry Smith }
14142d61bbb3SSatish Balay 
14154a2ae208SSatish Balay #undef __FUNCT__
14164a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_7"
1417dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_7(Mat A,Vec xx,Vec yy,Vec zz)
14182d61bbb3SSatish Balay {
14192d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
1420dcf5cc72SBarry Smith   PetscScalar    *x,*y = 0,*z = 0,*xb,sum1,sum2,sum3,sum4,sum5,sum6,sum7;
142126e093fcSHong Zhang   PetscScalar    x1,x2,x3,x4,x5,x6,x7,*yarray,*zarray;
14223f1db9ecSBarry Smith   MatScalar      *v;
1423dfbe8321SBarry Smith   PetscErrorCode ierr;
14247b2bb3b9SHong Zhang   PetscInt       mbs=a->mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL;
142526e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
14262d61bbb3SSatish Balay 
14272d61bbb3SSatish Balay   PetscFunctionBegin;
14281ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
142926e093fcSHong Zhang   ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr);
14302e8a6d31SBarry Smith   if (zz != yy) {
143126e093fcSHong Zhang     ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
14322e8a6d31SBarry Smith   } else {
143326e093fcSHong Zhang     zarray = yarray;
14342e8a6d31SBarry Smith   }
14352d61bbb3SSatish Balay 
14362d61bbb3SSatish Balay   idx = a->j;
14372d61bbb3SSatish Balay   v   = a->a;
143826e093fcSHong Zhang   if (usecprow){
14394eb6d288SHong Zhang     if (zz != yy){
14404eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,7*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
14414eb6d288SHong Zhang     }
144226e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
144326e093fcSHong Zhang     ii   = a->compressedrow.i;
14447b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
144526e093fcSHong Zhang   } else {
14462d61bbb3SSatish Balay     ii  = a->i;
144726e093fcSHong Zhang     y   = yarray;
144826e093fcSHong Zhang     z   = zarray;
144926e093fcSHong Zhang   }
14502d61bbb3SSatish Balay 
14512d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
14522d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
145326e093fcSHong Zhang     if (usecprow){
14547b2bb3b9SHong Zhang       z = zarray + 7*ridx[i];
14557b2bb3b9SHong Zhang       y = yarray + 7*ridx[i];
145626e093fcSHong Zhang     }
14572d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4]; sum6 = y[5]; sum7 = y[6];
14582d61bbb3SSatish Balay     for (j=0; j<n; j++) {
14592d61bbb3SSatish Balay       xb = x + 7*(*idx++);
14602d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6];
14612d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[7]*x2  + v[14]*x3  + v[21]*x4 + v[28]*x5 + v[35]*x6 + v[42]*x7;
14622d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[8]*x2  + v[15]*x3  + v[22]*x4 + v[29]*x5 + v[36]*x6 + v[43]*x7;
14632d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[9]*x2  + v[16]*x3  + v[23]*x4 + v[30]*x5 + v[37]*x6 + v[44]*x7;
14642d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[10]*x2 + v[17]*x3  + v[24]*x4 + v[31]*x5 + v[38]*x6 + v[45]*x7;
14652d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[11]*x2 + v[18]*x3  + v[25]*x4 + v[32]*x5 + v[39]*x6 + v[46]*x7;
14662d61bbb3SSatish Balay       sum6 += v[5]*x1 + v[12]*x2 + v[19]*x3  + v[26]*x4 + v[33]*x5 + v[40]*x6 + v[47]*x7;
14672d61bbb3SSatish Balay       sum7 += v[6]*x1 + v[13]*x2 + v[20]*x3  + v[27]*x4 + v[34]*x5 + v[41]*x6 + v[48]*x7;
14682d61bbb3SSatish Balay       v += 49;
14692d61bbb3SSatish Balay     }
14702d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
147126e093fcSHong Zhang     if (!usecprow){
14722d61bbb3SSatish Balay       z += 7; y += 7;
14732d61bbb3SSatish Balay     }
147426e093fcSHong Zhang   }
14751ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
147626e093fcSHong Zhang   ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr);
14772e8a6d31SBarry Smith   if (zz != yy) {
147826e093fcSHong Zhang     ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
14792e8a6d31SBarry Smith   }
1480dc0b31edSSatish Balay   ierr = PetscLogFlops(98.0*a->nz);CHKERRQ(ierr);
14812d61bbb3SSatish Balay   PetscFunctionReturn(0);
14822d61bbb3SSatish Balay }
1483218c64b6SSatish Balay 
14844a2ae208SSatish Balay #undef __FUNCT__
14854a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_N"
1486dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_N(Mat A,Vec xx,Vec yy,Vec zz)
14872d61bbb3SSatish Balay {
14882d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
1489bba805e6SBarry Smith   PetscScalar    *x,*z = 0,*xb,*work,*workt,*zarray;
14903f1db9ecSBarry Smith   MatScalar      *v;
14916849ba73SBarry Smith   PetscErrorCode ierr;
1492d0f46423SBarry Smith   PetscInt       mbs,i,*idx,*ii,bs=A->rmap->bs,j,n,bs2=a->bs2;
14937b2bb3b9SHong Zhang   PetscInt       ncols,k,*ridx=PETSC_NULL;
149426e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
1495218c64b6SSatish Balay 
14962d61bbb3SSatish Balay   PetscFunctionBegin;
14976a65c766SHong Zhang   ierr = VecCopy_Seq(yy,zz);CHKERRQ(ierr);
14981ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
149926e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
15002d61bbb3SSatish Balay 
15012d61bbb3SSatish Balay   idx = a->j;
15022d61bbb3SSatish Balay   v   = a->a;
150326e093fcSHong Zhang   if (usecprow){
150426e093fcSHong Zhang     mbs    = a->compressedrow.nrows;
150526e093fcSHong Zhang     ii     = a->compressedrow.i;
15067b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
150726e093fcSHong Zhang   } else {
150826e093fcSHong Zhang     mbs = a->mbs;
15092d61bbb3SSatish Balay     ii  = a->i;
151026e093fcSHong Zhang     z   = zarray;
151126e093fcSHong Zhang   }
15122d61bbb3SSatish Balay 
15132d61bbb3SSatish Balay   if (!a->mult_work) {
1514d0f46423SBarry Smith     k    = PetscMax(A->rmap->n,A->cmap->n);
151587828ca2SBarry Smith     ierr = PetscMalloc((k+1)*sizeof(PetscScalar),&a->mult_work);CHKERRQ(ierr);
15162d61bbb3SSatish Balay   }
15172d61bbb3SSatish Balay   work = a->mult_work;
15182d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
15192d61bbb3SSatish Balay     n     = ii[1] - ii[0]; ii++;
15202d61bbb3SSatish Balay     ncols = n*bs;
15212d61bbb3SSatish Balay     workt = work;
15222d61bbb3SSatish Balay     for (j=0; j<n; j++) {
15232d61bbb3SSatish Balay       xb = x + bs*(*idx++);
15242d61bbb3SSatish Balay       for (k=0; k<bs; k++) workt[k] = xb[k];
15252d61bbb3SSatish Balay       workt += bs;
15262d61bbb3SSatish Balay     }
15277b2bb3b9SHong Zhang     if (usecprow) z = zarray + bs*ridx[i];
15283f1db9ecSBarry Smith     Kernel_w_gets_w_plus_Ar_times_v(bs,ncols,work,v,z);
152971044d3cSBarry Smith     /* BLASgemv_("N",&bs,&ncols,&_DOne,v,&bs,work,&_One,&_DOne,z,&_One); */
15302d61bbb3SSatish Balay     v += n*bs2;
153126e093fcSHong Zhang     if (!usecprow){
15322d61bbb3SSatish Balay       z += bs;
15332d61bbb3SSatish Balay     }
153426e093fcSHong Zhang   }
15351ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
153626e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
1537dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*a->nz*bs2);CHKERRQ(ierr);
15382d61bbb3SSatish Balay   PetscFunctionReturn(0);
15392d61bbb3SSatish Balay }
15402d61bbb3SSatish Balay 
15414a2ae208SSatish Balay #undef __FUNCT__
1542547795f9SHong Zhang #define __FUNCT__ "MatMultHermitianTranspose_SeqBAIJ"
1543547795f9SHong Zhang PetscErrorCode MatMultHermitianTranspose_SeqBAIJ(Mat A,Vec xx,Vec zz)
1544547795f9SHong Zhang {
1545547795f9SHong Zhang   PetscScalar    zero = 0.0;
1546547795f9SHong Zhang   PetscErrorCode ierr;
1547547795f9SHong Zhang 
1548547795f9SHong Zhang   PetscFunctionBegin;
1549547795f9SHong Zhang   ierr = VecSet(zz,zero);CHKERRQ(ierr);
1550547795f9SHong Zhang   ierr = MatMultHermitianTransposeAdd_SeqBAIJ(A,xx,zz,zz);CHKERRQ(ierr);
1551547795f9SHong Zhang   PetscFunctionReturn(0);
1552547795f9SHong Zhang }
1553547795f9SHong Zhang 
1554547795f9SHong Zhang #undef __FUNCT__
15554a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_SeqBAIJ"
1556dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_SeqBAIJ(Mat A,Vec xx,Vec zz)
15572d61bbb3SSatish Balay {
15583447b6efSHong Zhang   PetscScalar    zero = 0.0;
15596849ba73SBarry Smith   PetscErrorCode ierr;
15602d61bbb3SSatish Balay 
15612d61bbb3SSatish Balay   PetscFunctionBegin;
15622dcb1b2aSMatthew Knepley   ierr = VecSet(zz,zero);CHKERRQ(ierr);
15633447b6efSHong Zhang   ierr = MatMultTransposeAdd_SeqBAIJ(A,xx,zz,zz);CHKERRQ(ierr);
15642d61bbb3SSatish Balay   PetscFunctionReturn(0);
15652d61bbb3SSatish Balay }
15662d61bbb3SSatish Balay 
15674a2ae208SSatish Balay #undef __FUNCT__
1568547795f9SHong Zhang #define __FUNCT__ "MatMultHermitianTransposeAdd_SeqBAIJ"
1569547795f9SHong Zhang PetscErrorCode MatMultHermitianTransposeAdd_SeqBAIJ(Mat A,Vec xx,Vec yy,Vec zz)
1570547795f9SHong Zhang 
1571547795f9SHong Zhang {
1572547795f9SHong Zhang   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
1573547795f9SHong Zhang   PetscScalar       *zb,*x,*z,*xb = 0,x1,x2,x3,x4,x5;
1574547795f9SHong Zhang   MatScalar         *v;
1575547795f9SHong Zhang   PetscErrorCode    ierr;
1576547795f9SHong Zhang   PetscInt          mbs,i,*idx,*ii,rval,bs=A->rmap->bs,j,n,bs2=a->bs2,*ib,*ridx=PETSC_NULL;
1577547795f9SHong Zhang   Mat_CompressedRow cprow = a->compressedrow;
1578547795f9SHong Zhang   PetscTruth        usecprow=cprow.use;
1579547795f9SHong Zhang 
1580547795f9SHong Zhang   PetscFunctionBegin;
1581547795f9SHong Zhang   if (yy != zz) { ierr = VecCopy_Seq(yy,zz);CHKERRQ(ierr); }
1582547795f9SHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1583547795f9SHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
1584547795f9SHong Zhang 
1585547795f9SHong Zhang   idx = a->j;
1586547795f9SHong Zhang   v   = a->a;
1587547795f9SHong Zhang   if (usecprow){
1588547795f9SHong Zhang     mbs  = cprow.nrows;
1589547795f9SHong Zhang     ii   = cprow.i;
1590547795f9SHong Zhang     ridx = cprow.rindex;
1591547795f9SHong Zhang   } else {
1592547795f9SHong Zhang     mbs=a->mbs;
1593547795f9SHong Zhang     ii = a->i;
1594547795f9SHong Zhang     xb = x;
1595547795f9SHong Zhang   }
1596547795f9SHong Zhang 
1597547795f9SHong Zhang   switch (bs) {
1598547795f9SHong Zhang   case 1:
1599547795f9SHong Zhang     for (i=0; i<mbs; i++) {
1600547795f9SHong Zhang       if (usecprow) xb = x + ridx[i];
1601547795f9SHong Zhang       x1 = xb[0];
1602547795f9SHong Zhang       ib = idx + ii[0];
1603547795f9SHong Zhang       n  = ii[1] - ii[0]; ii++;
1604547795f9SHong Zhang       for (j=0; j<n; j++) {
1605547795f9SHong Zhang         rval    = ib[j];
1606547795f9SHong Zhang         z[rval] += PetscConj(*v) * x1;
1607547795f9SHong Zhang         v++;
1608547795f9SHong Zhang       }
1609547795f9SHong Zhang       if (!usecprow) xb++;
1610547795f9SHong Zhang     }
1611547795f9SHong Zhang     break;
1612547795f9SHong Zhang   case 2:
1613547795f9SHong Zhang     for (i=0; i<mbs; i++) {
1614547795f9SHong Zhang       if (usecprow) xb = x + 2*ridx[i];
1615547795f9SHong Zhang       x1 = xb[0]; x2 = xb[1];
1616547795f9SHong Zhang       ib = idx + ii[0];
1617547795f9SHong Zhang       n  = ii[1] - ii[0]; ii++;
1618547795f9SHong Zhang       for (j=0; j<n; j++) {
1619547795f9SHong Zhang         rval      = ib[j]*2;
1620547795f9SHong Zhang         z[rval++] += PetscConj(v[0])*x1 + PetscConj(v[1])*x2;
1621547795f9SHong Zhang         z[rval++] += PetscConj(v[2])*x1 + PetscConj(v[3])*x2;
1622547795f9SHong Zhang         v  += 4;
1623547795f9SHong Zhang       }
1624547795f9SHong Zhang       if (!usecprow) xb += 2;
1625547795f9SHong Zhang     }
1626547795f9SHong Zhang     break;
1627547795f9SHong Zhang   case 3:
1628547795f9SHong Zhang     for (i=0; i<mbs; i++) {
1629547795f9SHong Zhang       if (usecprow) xb = x + 3*ridx[i];
1630547795f9SHong Zhang       x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
1631547795f9SHong Zhang       ib = idx + ii[0];
1632547795f9SHong Zhang       n  = ii[1] - ii[0]; ii++;
1633547795f9SHong Zhang       for (j=0; j<n; j++) {
1634547795f9SHong Zhang         rval      = ib[j]*3;
1635547795f9SHong Zhang         z[rval++] += PetscConj(v[0])*x1 + PetscConj(v[1])*x2 + PetscConj(v[2])*x3;
1636547795f9SHong Zhang         z[rval++] += PetscConj(v[3])*x1 + PetscConj(v[4])*x2 + PetscConj(v[5])*x3;
1637547795f9SHong Zhang         z[rval++] += PetscConj(v[6])*x1 + PetscConj(v[7])*x2 + PetscConj(v[8])*x3;
1638547795f9SHong Zhang         v  += 9;
1639547795f9SHong Zhang       }
1640547795f9SHong Zhang       if (!usecprow) xb += 3;
1641547795f9SHong Zhang     }
1642547795f9SHong Zhang     break;
1643547795f9SHong Zhang   case 4:
1644547795f9SHong Zhang     for (i=0; i<mbs; i++) {
1645547795f9SHong Zhang       if (usecprow) xb = x + 4*ridx[i];
1646547795f9SHong Zhang       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
1647547795f9SHong Zhang       ib = idx + ii[0];
1648547795f9SHong Zhang       n  = ii[1] - ii[0]; ii++;
1649547795f9SHong Zhang       for (j=0; j<n; j++) {
1650547795f9SHong Zhang         rval      = ib[j]*4;
1651547795f9SHong Zhang         z[rval++] +=  PetscConj(v[0])*x1 + PetscConj(v[1])*x2  + PetscConj(v[2])*x3  + PetscConj(v[3])*x4;
1652547795f9SHong Zhang         z[rval++] +=  PetscConj(v[4])*x1 + PetscConj(v[5])*x2  + PetscConj(v[6])*x3  + PetscConj(v[7])*x4;
1653547795f9SHong Zhang         z[rval++] +=  PetscConj(v[8])*x1 + PetscConj(v[9])*x2  + PetscConj(v[10])*x3 + PetscConj(v[11])*x4;
1654547795f9SHong Zhang         z[rval++] += PetscConj(v[12])*x1 + PetscConj(v[13])*x2 + PetscConj(v[14])*x3 + PetscConj(v[15])*x4;
1655547795f9SHong Zhang         v  += 16;
1656547795f9SHong Zhang       }
1657547795f9SHong Zhang       if (!usecprow) xb += 4;
1658547795f9SHong Zhang     }
1659547795f9SHong Zhang     break;
1660547795f9SHong Zhang   case 5:
1661547795f9SHong Zhang     for (i=0; i<mbs; i++) {
1662547795f9SHong Zhang       if (usecprow) xb = x + 5*ridx[i];
1663547795f9SHong Zhang       x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
1664547795f9SHong Zhang       x4 = xb[3]; x5 = xb[4];
1665547795f9SHong Zhang       ib = idx + ii[0];
1666547795f9SHong Zhang       n  = ii[1] - ii[0]; ii++;
1667547795f9SHong Zhang       for (j=0; j<n; j++) {
1668547795f9SHong Zhang         rval      = ib[j]*5;
1669547795f9SHong Zhang         z[rval++] +=  PetscConj(v[0])*x1 +  PetscConj(v[1])*x2 +  PetscConj(v[2])*x3 +  PetscConj(v[3])*x4 +  PetscConj(v[4])*x5;
1670547795f9SHong Zhang         z[rval++] +=  PetscConj(v[5])*x1 +  PetscConj(v[6])*x2 +  PetscConj(v[7])*x3 +  PetscConj(v[8])*x4 +  PetscConj(v[9])*x5;
1671547795f9SHong Zhang         z[rval++] += PetscConj(v[10])*x1 + PetscConj(v[11])*x2 + PetscConj(v[12])*x3 + PetscConj(v[13])*x4 + PetscConj(v[14])*x5;
1672547795f9SHong Zhang         z[rval++] += PetscConj(v[15])*x1 + PetscConj(v[16])*x2 + PetscConj(v[17])*x3 + PetscConj(v[18])*x4 + PetscConj(v[19])*x5;
1673547795f9SHong Zhang         z[rval++] += PetscConj(v[20])*x1 + PetscConj(v[21])*x2 + PetscConj(v[22])*x3 + PetscConj(v[23])*x4 + PetscConj(v[24])*x5;
1674547795f9SHong Zhang         v  += 25;
1675547795f9SHong Zhang       }
1676547795f9SHong Zhang       if (!usecprow) xb += 5;
1677547795f9SHong Zhang     }
1678547795f9SHong Zhang     break;
1679547795f9SHong Zhang   default: {      /* block sizes larger than 5 by 5 are handled by BLAS */
1680547795f9SHong Zhang       PetscInt     ncols,k;
1681547795f9SHong Zhang       PetscScalar  *work,*workt,*xtmp;
1682547795f9SHong Zhang 
1683*e32f2f54SBarry Smith       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"block size larger than 5 is not supported yet");
1684547795f9SHong Zhang       if (!a->mult_work) {
1685547795f9SHong Zhang         k = PetscMax(A->rmap->n,A->cmap->n);
1686547795f9SHong Zhang         ierr = PetscMalloc((k+1)*sizeof(PetscScalar),&a->mult_work);CHKERRQ(ierr);
1687547795f9SHong Zhang       }
1688547795f9SHong Zhang       work = a->mult_work;
1689547795f9SHong Zhang       xtmp = x;
1690547795f9SHong Zhang       for (i=0; i<mbs; i++) {
1691547795f9SHong Zhang         n     = ii[1] - ii[0]; ii++;
1692547795f9SHong Zhang         ncols = n*bs;
1693547795f9SHong Zhang         ierr  = PetscMemzero(work,ncols*sizeof(PetscScalar));CHKERRQ(ierr);
1694547795f9SHong Zhang         if (usecprow) {
1695547795f9SHong Zhang           xtmp = x + bs*ridx[i];
1696547795f9SHong Zhang         }
1697547795f9SHong Zhang         Kernel_w_gets_w_plus_trans_Ar_times_v(bs,ncols,xtmp,v,work);
1698547795f9SHong Zhang         /* BLASgemv_("T",&bs,&ncols,&_DOne,v,&bs,xtmp,&_One,&_DOne,work,&_One); */
1699547795f9SHong Zhang         v += n*bs2;
1700547795f9SHong Zhang         if (!usecprow) xtmp += bs;
1701547795f9SHong Zhang         workt = work;
1702547795f9SHong Zhang         for (j=0; j<n; j++) {
1703547795f9SHong Zhang           zb = z + bs*(*idx++);
1704547795f9SHong Zhang           for (k=0; k<bs; k++) zb[k] += workt[k] ;
1705547795f9SHong Zhang           workt += bs;
1706547795f9SHong Zhang         }
1707547795f9SHong Zhang       }
1708547795f9SHong Zhang     }
1709547795f9SHong Zhang   }
1710547795f9SHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
1711547795f9SHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
1712547795f9SHong Zhang   ierr = PetscLogFlops(2.0*a->nz*a->bs2);CHKERRQ(ierr);
1713547795f9SHong Zhang   PetscFunctionReturn(0);
1714547795f9SHong Zhang }
1715547795f9SHong Zhang 
1716547795f9SHong Zhang #undef __FUNCT__
17174a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_SeqBAIJ"
1718dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_SeqBAIJ(Mat A,Vec xx,Vec yy,Vec zz)
17192d61bbb3SSatish Balay {
17202d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
1721dcf5cc72SBarry Smith   PetscScalar       *zb,*x,*z,*xb = 0,x1,x2,x3,x4,x5;
17223f1db9ecSBarry Smith   MatScalar         *v;
17236849ba73SBarry Smith   PetscErrorCode    ierr;
1724d0f46423SBarry Smith   PetscInt          mbs,i,*idx,*ii,rval,bs=A->rmap->bs,j,n,bs2=a->bs2,*ib,*ridx=PETSC_NULL;
17253447b6efSHong Zhang   Mat_CompressedRow cprow = a->compressedrow;
17263447b6efSHong Zhang   PetscTruth        usecprow=cprow.use;
17272d61bbb3SSatish Balay 
17282d61bbb3SSatish Balay   PetscFunctionBegin;
17296a65c766SHong Zhang   if (yy != zz) { ierr = VecCopy_Seq(yy,zz);CHKERRQ(ierr); }
17303447b6efSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
17313447b6efSHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
17322d61bbb3SSatish Balay 
17332d61bbb3SSatish Balay   idx = a->j;
17342d61bbb3SSatish Balay   v   = a->a;
17353447b6efSHong Zhang   if (usecprow){
17363447b6efSHong Zhang     mbs  = cprow.nrows;
17373447b6efSHong Zhang     ii   = cprow.i;
17387b2bb3b9SHong Zhang     ridx = cprow.rindex;
17393447b6efSHong Zhang   } else {
17403447b6efSHong Zhang     mbs=a->mbs;
17412d61bbb3SSatish Balay     ii = a->i;
1742f1af5d2fSBarry Smith     xb = x;
17433447b6efSHong Zhang   }
17442d61bbb3SSatish Balay 
17452d61bbb3SSatish Balay   switch (bs) {
17462d61bbb3SSatish Balay   case 1:
17472d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
17487b2bb3b9SHong Zhang       if (usecprow) xb = x + ridx[i];
1749f1af5d2fSBarry Smith       x1 = xb[0];
17503447b6efSHong Zhang       ib = idx + ii[0];
17513447b6efSHong Zhang       n  = ii[1] - ii[0]; ii++;
17522d61bbb3SSatish Balay       for (j=0; j<n; j++) {
17532d61bbb3SSatish Balay         rval    = ib[j];
1754f1af5d2fSBarry Smith         z[rval] += *v * x1;
1755f1af5d2fSBarry Smith         v++;
17562d61bbb3SSatish Balay       }
17573447b6efSHong Zhang       if (!usecprow) xb++;
17582d61bbb3SSatish Balay     }
17592d61bbb3SSatish Balay     break;
17602d61bbb3SSatish Balay   case 2:
17612d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
17627b2bb3b9SHong Zhang       if (usecprow) xb = x + 2*ridx[i];
1763f1af5d2fSBarry Smith       x1 = xb[0]; x2 = xb[1];
17643447b6efSHong Zhang       ib = idx + ii[0];
17653447b6efSHong Zhang       n  = ii[1] - ii[0]; ii++;
17662d61bbb3SSatish Balay       for (j=0; j<n; j++) {
17672d61bbb3SSatish Balay         rval      = ib[j]*2;
17682d61bbb3SSatish Balay         z[rval++] += v[0]*x1 + v[1]*x2;
17692d61bbb3SSatish Balay         z[rval++] += v[2]*x1 + v[3]*x2;
17702d61bbb3SSatish Balay         v  += 4;
17712d61bbb3SSatish Balay       }
17723447b6efSHong Zhang       if (!usecprow) xb += 2;
17732d61bbb3SSatish Balay     }
17742d61bbb3SSatish Balay     break;
17752d61bbb3SSatish Balay   case 3:
17762d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
17777b2bb3b9SHong Zhang       if (usecprow) xb = x + 3*ridx[i];
1778f1af5d2fSBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
17793447b6efSHong Zhang       ib = idx + ii[0];
17803447b6efSHong Zhang       n  = ii[1] - ii[0]; ii++;
17812d61bbb3SSatish Balay       for (j=0; j<n; j++) {
17822d61bbb3SSatish Balay         rval      = ib[j]*3;
17832d61bbb3SSatish Balay         z[rval++] += v[0]*x1 + v[1]*x2 + v[2]*x3;
17842d61bbb3SSatish Balay         z[rval++] += v[3]*x1 + v[4]*x2 + v[5]*x3;
17852d61bbb3SSatish Balay         z[rval++] += v[6]*x1 + v[7]*x2 + v[8]*x3;
17862d61bbb3SSatish Balay         v  += 9;
17872d61bbb3SSatish Balay       }
17883447b6efSHong Zhang       if (!usecprow) xb += 3;
17892d61bbb3SSatish Balay     }
17902d61bbb3SSatish Balay     break;
17912d61bbb3SSatish Balay   case 4:
17922d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
17937b2bb3b9SHong Zhang       if (usecprow) xb = x + 4*ridx[i];
1794f1af5d2fSBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
17953447b6efSHong Zhang       ib = idx + ii[0];
17963447b6efSHong Zhang       n  = ii[1] - ii[0]; ii++;
17972d61bbb3SSatish Balay       for (j=0; j<n; j++) {
17982d61bbb3SSatish Balay         rval      = ib[j]*4;
17992d61bbb3SSatish Balay         z[rval++] +=  v[0]*x1 +  v[1]*x2 +  v[2]*x3 +  v[3]*x4;
18002d61bbb3SSatish Balay         z[rval++] +=  v[4]*x1 +  v[5]*x2 +  v[6]*x3 +  v[7]*x4;
18012d61bbb3SSatish Balay         z[rval++] +=  v[8]*x1 +  v[9]*x2 + v[10]*x3 + v[11]*x4;
18022d61bbb3SSatish Balay         z[rval++] += v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4;
18032d61bbb3SSatish Balay         v  += 16;
18042d61bbb3SSatish Balay       }
18053447b6efSHong Zhang       if (!usecprow) xb += 4;
18062d61bbb3SSatish Balay     }
18072d61bbb3SSatish Balay     break;
18082d61bbb3SSatish Balay   case 5:
18092d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
18107b2bb3b9SHong Zhang       if (usecprow) xb = x + 5*ridx[i];
1811f1af5d2fSBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
18122d61bbb3SSatish Balay       x4 = xb[3]; x5 = xb[4];
18133447b6efSHong Zhang       ib = idx + ii[0];
18143447b6efSHong Zhang       n  = ii[1] - ii[0]; ii++;
18152d61bbb3SSatish Balay       for (j=0; j<n; j++) {
18162d61bbb3SSatish Balay         rval      = ib[j]*5;
18172d61bbb3SSatish Balay         z[rval++] +=  v[0]*x1 +  v[1]*x2 +  v[2]*x3 +  v[3]*x4 +  v[4]*x5;
18182d61bbb3SSatish Balay         z[rval++] +=  v[5]*x1 +  v[6]*x2 +  v[7]*x3 +  v[8]*x4 +  v[9]*x5;
18192d61bbb3SSatish Balay         z[rval++] += v[10]*x1 + v[11]*x2 + v[12]*x3 + v[13]*x4 + v[14]*x5;
18202d61bbb3SSatish Balay         z[rval++] += v[15]*x1 + v[16]*x2 + v[17]*x3 + v[18]*x4 + v[19]*x5;
18212d61bbb3SSatish Balay         z[rval++] += v[20]*x1 + v[21]*x2 + v[22]*x3 + v[23]*x4 + v[24]*x5;
18222d61bbb3SSatish Balay         v  += 25;
18232d61bbb3SSatish Balay       }
18243447b6efSHong Zhang       if (!usecprow) xb += 5;
18252d61bbb3SSatish Balay     }
18262d61bbb3SSatish Balay     break;
1827f1af5d2fSBarry Smith   default: {      /* block sizes larger then 5 by 5 are handled by BLAS */
1828690b6cddSBarry Smith       PetscInt     ncols,k;
18293447b6efSHong Zhang       PetscScalar  *work,*workt,*xtmp;
18303f1db9ecSBarry Smith 
18312d61bbb3SSatish Balay       if (!a->mult_work) {
1832d0f46423SBarry Smith         k = PetscMax(A->rmap->n,A->cmap->n);
183387828ca2SBarry Smith         ierr = PetscMalloc((k+1)*sizeof(PetscScalar),&a->mult_work);CHKERRQ(ierr);
18342d61bbb3SSatish Balay       }
18352d61bbb3SSatish Balay       work = a->mult_work;
18363447b6efSHong Zhang       xtmp = x;
18372d61bbb3SSatish Balay       for (i=0; i<mbs; i++) {
18382d61bbb3SSatish Balay         n     = ii[1] - ii[0]; ii++;
18392d61bbb3SSatish Balay         ncols = n*bs;
184087828ca2SBarry Smith         ierr  = PetscMemzero(work,ncols*sizeof(PetscScalar));CHKERRQ(ierr);
18413447b6efSHong Zhang         if (usecprow) {
18427b2bb3b9SHong Zhang           xtmp = x + bs*ridx[i];
18433447b6efSHong Zhang         }
18443447b6efSHong Zhang         Kernel_w_gets_w_plus_trans_Ar_times_v(bs,ncols,xtmp,v,work);
184571044d3cSBarry Smith         /* BLASgemv_("T",&bs,&ncols,&_DOne,v,&bs,xtmp,&_One,&_DOne,work,&_One); */
18462d61bbb3SSatish Balay         v += n*bs2;
18473447b6efSHong Zhang         if (!usecprow) xtmp += bs;
18482d61bbb3SSatish Balay         workt = work;
18492d61bbb3SSatish Balay         for (j=0; j<n; j++) {
18502d61bbb3SSatish Balay           zb = z + bs*(*idx++);
18512d61bbb3SSatish Balay           for (k=0; k<bs; k++) zb[k] += workt[k] ;
18522d61bbb3SSatish Balay           workt += bs;
18532d61bbb3SSatish Balay         }
18542d61bbb3SSatish Balay       }
18552d61bbb3SSatish Balay     }
18562d61bbb3SSatish Balay   }
18573447b6efSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
18583447b6efSHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
1859dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*a->nz*a->bs2);CHKERRQ(ierr);
18602d61bbb3SSatish Balay   PetscFunctionReturn(0);
18612d61bbb3SSatish Balay }
18622d61bbb3SSatish Balay 
18634a2ae208SSatish Balay #undef __FUNCT__
18644a2ae208SSatish Balay #define __FUNCT__ "MatScale_SeqBAIJ"
1865f4df32b1SMatthew Knepley PetscErrorCode MatScale_SeqBAIJ(Mat inA,PetscScalar alpha)
18662d61bbb3SSatish Balay {
18672d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)inA->data;
1868690b6cddSBarry Smith   PetscInt       totalnz = a->bs2*a->nz;
1869f4df32b1SMatthew Knepley   PetscScalar    oalpha = alpha;
1870efee365bSSatish Balay   PetscErrorCode ierr;
18710805154bSBarry Smith   PetscBLASInt   one = 1,tnz = PetscBLASIntCast(totalnz);
18722d61bbb3SSatish Balay 
18732d61bbb3SSatish Balay   PetscFunctionBegin;
1874f4df32b1SMatthew Knepley   BLASscal_(&tnz,&oalpha,a->a,&one);
1875efee365bSSatish Balay   ierr = PetscLogFlops(totalnz);CHKERRQ(ierr);
18762d61bbb3SSatish Balay   PetscFunctionReturn(0);
18772d61bbb3SSatish Balay }
18782d61bbb3SSatish Balay 
18794a2ae208SSatish Balay #undef __FUNCT__
18804a2ae208SSatish Balay #define __FUNCT__ "MatNorm_SeqBAIJ"
1881dfbe8321SBarry Smith PetscErrorCode MatNorm_SeqBAIJ(Mat A,NormType type,PetscReal *norm)
18822d61bbb3SSatish Balay {
18838a62d963SHong Zhang   PetscErrorCode ierr;
18842d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
18853f1db9ecSBarry Smith   MatScalar      *v = a->a;
1886329f5518SBarry Smith   PetscReal      sum = 0.0;
1887d0f46423SBarry Smith   PetscInt       i,j,k,bs=A->rmap->bs,nz=a->nz,bs2=a->bs2,k1;
18882d61bbb3SSatish Balay 
18892d61bbb3SSatish Balay   PetscFunctionBegin;
18902d61bbb3SSatish Balay   if (type == NORM_FROBENIUS) {
18912d61bbb3SSatish Balay     for (i=0; i< bs2*nz; i++) {
1892aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
1893329f5518SBarry Smith       sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
18942d61bbb3SSatish Balay #else
18952d61bbb3SSatish Balay       sum += (*v)*(*v); v++;
18962d61bbb3SSatish Balay #endif
18972d61bbb3SSatish Balay     }
18982d61bbb3SSatish Balay     *norm = sqrt(sum);
18998a62d963SHong Zhang   } else if (type == NORM_1) { /* maximum column sum */
19008a62d963SHong Zhang     PetscReal *tmp;
19018a62d963SHong Zhang     PetscInt  *bcol = a->j;
1902d0f46423SBarry Smith     ierr = PetscMalloc((A->cmap->n+1)*sizeof(PetscReal),&tmp);CHKERRQ(ierr);
1903d0f46423SBarry Smith     ierr = PetscMemzero(tmp,A->cmap->n*sizeof(PetscReal));CHKERRQ(ierr);
19048a62d963SHong Zhang     for (i=0; i<nz; i++){
19058a62d963SHong Zhang       for (j=0; j<bs; j++){
19068a62d963SHong Zhang         k1 = bs*(*bcol) + j; /* column index */
19078a62d963SHong Zhang         for (k=0; k<bs; k++){
19088a62d963SHong Zhang           tmp[k1] += PetscAbsScalar(*v); v++;
19098a62d963SHong Zhang         }
19108a62d963SHong Zhang       }
19118a62d963SHong Zhang       bcol++;
19128a62d963SHong Zhang     }
19138a62d963SHong Zhang     *norm = 0.0;
1914d0f46423SBarry Smith     for (j=0; j<A->cmap->n; j++) {
19158a62d963SHong Zhang       if (tmp[j] > *norm) *norm = tmp[j];
19168a62d963SHong Zhang     }
19178a62d963SHong Zhang     ierr = PetscFree(tmp);CHKERRQ(ierr);
1918596552b5SBarry Smith   } else if (type == NORM_INFINITY) { /* maximum row sum */
1919596552b5SBarry Smith     *norm = 0.0;
1920596552b5SBarry Smith     for (k=0; k<bs; k++) {
192174f84c7bSSatish Balay       for (j=0; j<a->mbs; j++) {
1922596552b5SBarry Smith         v = a->a + bs2*a->i[j] + k;
1923596552b5SBarry Smith         sum = 0.0;
1924596552b5SBarry Smith         for (i=0; i<a->i[j+1]-a->i[j]; i++) {
19250e90e235SBarry Smith           for (k1=0; k1<bs; k1++){
1926596552b5SBarry Smith             sum += PetscAbsScalar(*v);
1927596552b5SBarry Smith             v   += bs;
19282d61bbb3SSatish Balay           }
19290e90e235SBarry Smith         }
1930596552b5SBarry Smith         if (sum > *norm) *norm = sum;
1931596552b5SBarry Smith       }
1932596552b5SBarry Smith     }
1933596552b5SBarry Smith   } else {
1934*e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for this norm yet");
19352d61bbb3SSatish Balay   }
19362d61bbb3SSatish Balay   PetscFunctionReturn(0);
19372d61bbb3SSatish Balay }
19382d61bbb3SSatish Balay 
19392d61bbb3SSatish Balay 
19404a2ae208SSatish Balay #undef __FUNCT__
19414a2ae208SSatish Balay #define __FUNCT__ "MatEqual_SeqBAIJ"
1942dfbe8321SBarry Smith PetscErrorCode MatEqual_SeqBAIJ(Mat A,Mat B,PetscTruth* flg)
19432d61bbb3SSatish Balay {
19442d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ *)A->data,*b = (Mat_SeqBAIJ *)B->data;
1945dfbe8321SBarry Smith   PetscErrorCode ierr;
19462d61bbb3SSatish Balay 
19472d61bbb3SSatish Balay   PetscFunctionBegin;
19482d61bbb3SSatish Balay   /* If the  matrix/block dimensions are not equal, or no of nonzeros or shift */
1949d0f46423SBarry Smith   if ((A->rmap->N != B->rmap->N) || (A->cmap->n != B->cmap->n) || (A->rmap->bs != B->rmap->bs)|| (a->nz != b->nz)) {
1950273d9f13SBarry Smith     *flg = PETSC_FALSE;
1951273d9f13SBarry Smith     PetscFunctionReturn(0);
19522d61bbb3SSatish Balay   }
19532d61bbb3SSatish Balay 
19542d61bbb3SSatish Balay   /* if the a->i are the same */
1955690b6cddSBarry Smith   ierr = PetscMemcmp(a->i,b->i,(a->mbs+1)*sizeof(PetscInt),flg);CHKERRQ(ierr);
1956abc0a331SBarry Smith   if (!*flg) {
19570f5bd95cSBarry Smith     PetscFunctionReturn(0);
19582d61bbb3SSatish Balay   }
19592d61bbb3SSatish Balay 
19602d61bbb3SSatish Balay   /* if a->j are the same */
1961690b6cddSBarry Smith   ierr = PetscMemcmp(a->j,b->j,(a->nz)*sizeof(PetscInt),flg);CHKERRQ(ierr);
1962abc0a331SBarry Smith   if (!*flg) {
19630f5bd95cSBarry Smith     PetscFunctionReturn(0);
19642d61bbb3SSatish Balay   }
19652d61bbb3SSatish Balay   /* if a->a are the same */
1966d0f46423SBarry Smith   ierr = PetscMemcmp(a->a,b->a,(a->nz)*(A->rmap->bs)*(B->rmap->bs)*sizeof(PetscScalar),flg);CHKERRQ(ierr);
19672d61bbb3SSatish Balay   PetscFunctionReturn(0);
19682d61bbb3SSatish Balay 
19692d61bbb3SSatish Balay }
19702d61bbb3SSatish Balay 
19714a2ae208SSatish Balay #undef __FUNCT__
19724a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_SeqBAIJ"
1973dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_SeqBAIJ(Mat A,Vec v)
19742d61bbb3SSatish Balay {
19752d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
1976dfbe8321SBarry Smith   PetscErrorCode ierr;
1977690b6cddSBarry Smith   PetscInt       i,j,k,n,row,bs,*ai,*aj,ambs,bs2;
197887828ca2SBarry Smith   PetscScalar    *x,zero = 0.0;
19793f1db9ecSBarry Smith   MatScalar      *aa,*aa_j;
19802d61bbb3SSatish Balay 
19812d61bbb3SSatish Balay   PetscFunctionBegin;
1982*e32f2f54SBarry Smith   if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1983d0f46423SBarry Smith   bs   = A->rmap->bs;
19842d61bbb3SSatish Balay   aa   = a->a;
19852d61bbb3SSatish Balay   ai   = a->i;
19862d61bbb3SSatish Balay   aj   = a->j;
19872d61bbb3SSatish Balay   ambs = a->mbs;
19882d61bbb3SSatish Balay   bs2  = a->bs2;
19892d61bbb3SSatish Balay 
19902dcb1b2aSMatthew Knepley   ierr = VecSet(v,zero);CHKERRQ(ierr);
19911ebc52fbSHong Zhang   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1992e1311b90SBarry Smith   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
1993*e32f2f54SBarry Smith   if (n != A->rmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
19942d61bbb3SSatish Balay   for (i=0; i<ambs; i++) {
19952d61bbb3SSatish Balay     for (j=ai[i]; j<ai[i+1]; j++) {
19962d61bbb3SSatish Balay       if (aj[j] == i) {
19972d61bbb3SSatish Balay         row  = i*bs;
19982d61bbb3SSatish Balay         aa_j = aa+j*bs2;
19992d61bbb3SSatish Balay         for (k=0; k<bs2; k+=(bs+1),row++) x[row] = aa_j[k];
20002d61bbb3SSatish Balay         break;
20012d61bbb3SSatish Balay       }
20022d61bbb3SSatish Balay     }
20032d61bbb3SSatish Balay   }
20041ebc52fbSHong Zhang   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
20052d61bbb3SSatish Balay   PetscFunctionReturn(0);
20062d61bbb3SSatish Balay }
20072d61bbb3SSatish Balay 
20084a2ae208SSatish Balay #undef __FUNCT__
20094a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_SeqBAIJ"
2010dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_SeqBAIJ(Mat A,Vec ll,Vec rr)
20112d61bbb3SSatish Balay {
20122d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
201387828ca2SBarry Smith   PetscScalar    *l,*r,x,*li,*ri;
20143f1db9ecSBarry Smith   MatScalar      *aa,*v;
2015dfbe8321SBarry Smith   PetscErrorCode ierr;
2016690b6cddSBarry Smith   PetscInt       i,j,k,lm,rn,M,m,n,*ai,*aj,mbs,tmp,bs,bs2;
20172d61bbb3SSatish Balay 
20182d61bbb3SSatish Balay   PetscFunctionBegin;
20192d61bbb3SSatish Balay   ai  = a->i;
20202d61bbb3SSatish Balay   aj  = a->j;
20212d61bbb3SSatish Balay   aa  = a->a;
2022d0f46423SBarry Smith   m   = A->rmap->n;
2023d0f46423SBarry Smith   n   = A->cmap->n;
2024d0f46423SBarry Smith   bs  = A->rmap->bs;
20252d61bbb3SSatish Balay   mbs = a->mbs;
20262d61bbb3SSatish Balay   bs2 = a->bs2;
20272d61bbb3SSatish Balay   if (ll) {
20281ebc52fbSHong Zhang     ierr = VecGetArray(ll,&l);CHKERRQ(ierr);
2029e1311b90SBarry Smith     ierr = VecGetLocalSize(ll,&lm);CHKERRQ(ierr);
2030*e32f2f54SBarry Smith     if (lm != m) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Left scaling vector wrong length");
20312d61bbb3SSatish Balay     for (i=0; i<mbs; i++) { /* for each block row */
20322d61bbb3SSatish Balay       M  = ai[i+1] - ai[i];
20332d61bbb3SSatish Balay       li = l + i*bs;
20342d61bbb3SSatish Balay       v  = aa + bs2*ai[i];
20352d61bbb3SSatish Balay       for (j=0; j<M; j++) { /* for each block */
20362d61bbb3SSatish Balay         for (k=0; k<bs2; k++) {
20372d61bbb3SSatish Balay           (*v++) *= li[k%bs];
20382d61bbb3SSatish Balay         }
20392d61bbb3SSatish Balay       }
20402d61bbb3SSatish Balay     }
20411ebc52fbSHong Zhang     ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr);
2042efee365bSSatish Balay     ierr = PetscLogFlops(a->nz);CHKERRQ(ierr);
20432d61bbb3SSatish Balay   }
20442d61bbb3SSatish Balay 
20452d61bbb3SSatish Balay   if (rr) {
20461ebc52fbSHong Zhang     ierr = VecGetArray(rr,&r);CHKERRQ(ierr);
2047e1311b90SBarry Smith     ierr = VecGetLocalSize(rr,&rn);CHKERRQ(ierr);
2048*e32f2f54SBarry Smith     if (rn != n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Right scaling vector wrong length");
20492d61bbb3SSatish Balay     for (i=0; i<mbs; i++) { /* for each block row */
20502d61bbb3SSatish Balay       M  = ai[i+1] - ai[i];
20512d61bbb3SSatish Balay       v  = aa + bs2*ai[i];
20522d61bbb3SSatish Balay       for (j=0; j<M; j++) { /* for each block */
20532d61bbb3SSatish Balay         ri = r + bs*aj[ai[i]+j];
20542d61bbb3SSatish Balay         for (k=0; k<bs; k++) {
20552d61bbb3SSatish Balay           x = ri[k];
20562d61bbb3SSatish Balay           for (tmp=0; tmp<bs; tmp++) (*v++) *= x;
20572d61bbb3SSatish Balay         }
20582d61bbb3SSatish Balay       }
20592d61bbb3SSatish Balay     }
20601ebc52fbSHong Zhang     ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr);
2061efee365bSSatish Balay     ierr = PetscLogFlops(a->nz);CHKERRQ(ierr);
20622d61bbb3SSatish Balay   }
20632d61bbb3SSatish Balay   PetscFunctionReturn(0);
20642d61bbb3SSatish Balay }
20652d61bbb3SSatish Balay 
20662d61bbb3SSatish Balay 
20674a2ae208SSatish Balay #undef __FUNCT__
20684a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_SeqBAIJ"
2069dfbe8321SBarry Smith PetscErrorCode MatGetInfo_SeqBAIJ(Mat A,MatInfoType flag,MatInfo *info)
20702d61bbb3SSatish Balay {
20712d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
20722d61bbb3SSatish Balay 
20732d61bbb3SSatish Balay   PetscFunctionBegin;
20742d61bbb3SSatish Balay   info->block_size     = a->bs2;
20752d61bbb3SSatish Balay   info->nz_allocated   = a->maxnz;
20762d61bbb3SSatish Balay   info->nz_used        = a->bs2*a->nz;
20772d61bbb3SSatish Balay   info->nz_unneeded    = (double)(info->nz_allocated - info->nz_used);
20782d61bbb3SSatish Balay   info->assemblies   = A->num_ass;
20792d61bbb3SSatish Balay   info->mallocs      = a->reallocs;
20807adad957SLisandro Dalcin   info->memory       = ((PetscObject)A)->mem;
2081d5f3da31SBarry Smith   if (A->factortype) {
20822d61bbb3SSatish Balay     info->fill_ratio_given  = A->info.fill_ratio_given;
20832d61bbb3SSatish Balay     info->fill_ratio_needed = A->info.fill_ratio_needed;
20842d61bbb3SSatish Balay     info->factor_mallocs    = A->info.factor_mallocs;
20852d61bbb3SSatish Balay   } else {
20862d61bbb3SSatish Balay     info->fill_ratio_given  = 0;
20872d61bbb3SSatish Balay     info->fill_ratio_needed = 0;
20882d61bbb3SSatish Balay     info->factor_mallocs    = 0;
20892d61bbb3SSatish Balay   }
20902d61bbb3SSatish Balay   PetscFunctionReturn(0);
20912d61bbb3SSatish Balay }
20922d61bbb3SSatish Balay 
20932d61bbb3SSatish Balay 
20944a2ae208SSatish Balay #undef __FUNCT__
20954a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_SeqBAIJ"
2096dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_SeqBAIJ(Mat A)
20972d61bbb3SSatish Balay {
20982d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
2099dfbe8321SBarry Smith   PetscErrorCode ierr;
21002d61bbb3SSatish Balay 
21012d61bbb3SSatish Balay   PetscFunctionBegin;
2102549d3d68SSatish Balay   ierr = PetscMemzero(a->a,a->bs2*a->i[a->mbs]*sizeof(MatScalar));CHKERRQ(ierr);
21032d61bbb3SSatish Balay   PetscFunctionReturn(0);
21042d61bbb3SSatish Balay }
2105