xref: /petsc/src/mat/impls/baij/seq/baij2.c (revision 7adad9577d7f34aa90d02da3975546c0571d04a7)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
2cac129eeSSatish Balay 
370f55243SBarry Smith #include "src/mat/impls/baij/seq/baij.h"
42d61bbb3SSatish Balay #include "src/inline/spops.h"
53f1db9ecSBarry Smith #include "src/inline/ilu.h"
60a835dfdSSatish Balay #include "petscbt.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;
14690b6cddSBarry Smith   PetscInt       row,i,j,k,l,m,n,*idx,*nidx,isz,val,ival;
15690b6cddSBarry Smith   PetscInt       start,end,*ai,*aj,bs,*nidx2;
16f1af5d2fSBarry Smith   PetscBT        table;
17a3192f15SSatish Balay 
183a40ed3dSBarry Smith   PetscFunctionBegin;
19a3192f15SSatish Balay   m     = a->mbs;
20a3192f15SSatish Balay   ai    = a->i;
21a3192f15SSatish Balay   aj    = a->j;
22899cda47SBarry Smith   bs    = A->rmap.bs;
23a3192f15SSatish Balay 
2429bbc08cSBarry Smith   if (ov < 0)  SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap specified");
25a3192f15SSatish Balay 
266831982aSBarry Smith   ierr = PetscBTCreate(m,table);CHKERRQ(ierr);
27690b6cddSBarry Smith   ierr = PetscMalloc((m+1)*sizeof(PetscInt),&nidx);CHKERRQ(ierr);
28899cda47SBarry Smith   ierr = PetscMalloc((A->rmap.N+1)*sizeof(PetscInt),&nidx2);CHKERRQ(ierr);
29a3192f15SSatish Balay 
30a3192f15SSatish Balay   for (i=0; i<is_max; i++) {
31a3192f15SSatish Balay     /* Initialise the two local arrays */
32a3192f15SSatish Balay     isz  = 0;
336831982aSBarry Smith     ierr = PetscBTMemzero(m,table);CHKERRQ(ierr);
34a3192f15SSatish Balay 
35a3192f15SSatish Balay     /* Extract the indices, assume there can be duplicate entries */
36a3192f15SSatish Balay     ierr = ISGetIndices(is[i],&idx);CHKERRQ(ierr);
37b9b97703SBarry Smith     ierr = ISGetLocalSize(is[i],&n);CHKERRQ(ierr);
38a3192f15SSatish Balay 
39a3192f15SSatish Balay     /* Enter these into the temp arrays i.e mark table[row], enter row into new index */
40a3192f15SSatish Balay     for (j=0; j<n ; ++j){
41218c64b6SSatish Balay       ival = idx[j]/bs; /* convert the indices into block indices */
425eee224dSHong Zhang       if (ival>=m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"index greater than mat-dim");
43f1af5d2fSBarry Smith       if(!PetscBTLookupSet(table,ival)) { nidx[isz++] = ival;}
44a3192f15SSatish Balay     }
45a3192f15SSatish Balay     ierr = ISRestoreIndices(is[i],&idx);CHKERRQ(ierr);
46a3192f15SSatish Balay     ierr = ISDestroy(is[i]);CHKERRQ(ierr);
47a3192f15SSatish Balay 
48a3192f15SSatish Balay     k = 0;
49a3192f15SSatish Balay     for (j=0; j<ov; j++){ /* for each overlap*/
50a3192f15SSatish Balay       n = isz;
51a3192f15SSatish Balay       for (; k<n ; k++){ /* do only those rows in nidx[k], which are not done yet */
52a3192f15SSatish Balay         row   = nidx[k];
53a3192f15SSatish Balay         start = ai[row];
54a3192f15SSatish Balay         end   = ai[row+1];
55a3192f15SSatish Balay         for (l = start; l<end ; l++){
56a3192f15SSatish Balay           val = aj[l];
57f1af5d2fSBarry Smith           if (!PetscBTLookupSet(table,val)) {nidx[isz++] = val;}
58a3192f15SSatish Balay         }
59a3192f15SSatish Balay       }
60a3192f15SSatish Balay     }
61218c64b6SSatish Balay     /* expand the Index Set */
62218c64b6SSatish Balay     for (j=0; j<isz; j++) {
63218c64b6SSatish Balay       for (k=0; k<bs; k++)
64218c64b6SSatish Balay         nidx2[j*bs+k] = nidx[j]*bs+k;
65218c64b6SSatish Balay     }
66329f5518SBarry Smith     ierr = ISCreateGeneral(PETSC_COMM_SELF,isz*bs,nidx2,is+i);CHKERRQ(ierr);
67a3192f15SSatish Balay   }
686831982aSBarry Smith   ierr = PetscBTDestroy(table);CHKERRQ(ierr);
69606d414cSSatish Balay   ierr = PetscFree(nidx);CHKERRQ(ierr);
70606d414cSSatish Balay   ierr = PetscFree(nidx2);CHKERRQ(ierr);
713a40ed3dSBarry Smith   PetscFunctionReturn(0);
72a3192f15SSatish Balay }
731c351548SSatish Balay 
744a2ae208SSatish Balay #undef __FUNCT__
754a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqBAIJ_Private"
76690b6cddSBarry Smith PetscErrorCode MatGetSubMatrix_SeqBAIJ_Private(Mat A,IS isrow,IS iscol,PetscInt cs,MatReuse scall,Mat *B)
77736121d4SSatish Balay {
78736121d4SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*c;
796849ba73SBarry Smith   PetscErrorCode ierr;
80690b6cddSBarry Smith   PetscInt       *smap,i,k,kstart,kend,oldcols = a->nbs,*lens;
81690b6cddSBarry Smith   PetscInt       row,mat_i,*mat_j,tcol,*mat_ilen;
82899cda47SBarry Smith   PetscInt       *irow,*icol,nrows,ncols,*ssmap,bs=A->rmap.bs,bs2=a->bs2;
83690b6cddSBarry Smith   PetscInt       *aj = a->j,*ai = a->i;
843f1db9ecSBarry Smith   MatScalar      *mat_a;
85736121d4SSatish Balay   Mat            C;
860f5bd95cSBarry Smith   PetscTruth     flag;
87736121d4SSatish Balay 
883a40ed3dSBarry Smith   PetscFunctionBegin;
89888f2ed8SSatish Balay   ierr = ISSorted(iscol,(PetscTruth*)&i);CHKERRQ(ierr);
9029bbc08cSBarry Smith   if (!i) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"IS is not sorted");
91736121d4SSatish Balay 
92736121d4SSatish Balay   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
93218c64b6SSatish Balay   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
94b9b97703SBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
95b9b97703SBarry Smith   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
96736121d4SSatish Balay 
97690b6cddSBarry Smith   ierr = PetscMalloc((1+oldcols)*sizeof(PetscInt),&smap);CHKERRQ(ierr);
98736121d4SSatish Balay   ssmap = smap;
99690b6cddSBarry Smith   ierr = PetscMalloc((1+nrows)*sizeof(PetscInt),&lens);CHKERRQ(ierr);
100690b6cddSBarry Smith   ierr  = PetscMemzero(smap,oldcols*sizeof(PetscInt));CHKERRQ(ierr);
101736121d4SSatish Balay   for (i=0; i<ncols; i++) smap[icol[i]] = i+1;
102736121d4SSatish Balay   /* determine lens of each row */
103736121d4SSatish Balay   for (i=0; i<nrows; i++) {
104736121d4SSatish Balay     kstart  = ai[irow[i]];
105736121d4SSatish Balay     kend    = kstart + a->ilen[irow[i]];
106736121d4SSatish Balay     lens[i] = 0;
107736121d4SSatish Balay       for (k=kstart; k<kend; k++) {
108736121d4SSatish Balay         if (ssmap[aj[k]]) {
109736121d4SSatish Balay           lens[i]++;
110736121d4SSatish Balay         }
111736121d4SSatish Balay       }
112736121d4SSatish Balay     }
113736121d4SSatish Balay   /* Create and fill new matrix */
114736121d4SSatish Balay   if (scall == MAT_REUSE_MATRIX) {
115736121d4SSatish Balay     c = (Mat_SeqBAIJ *)((*B)->data);
116736121d4SSatish Balay 
117899cda47SBarry Smith     if (c->mbs!=nrows || c->nbs!=ncols || (*B)->rmap.bs!=bs) SETERRQ(PETSC_ERR_ARG_SIZ,"Submatrix wrong size");
118690b6cddSBarry Smith     ierr = PetscMemcmp(c->ilen,lens,c->mbs *sizeof(PetscInt),&flag);CHKERRQ(ierr);
119abc0a331SBarry Smith     if (!flag) {
12029bbc08cSBarry Smith       SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong no of nonzeros");
121736121d4SSatish Balay     }
122690b6cddSBarry Smith     ierr = PetscMemzero(c->ilen,c->mbs*sizeof(PetscInt));CHKERRQ(ierr);
123736121d4SSatish Balay     C = *B;
1243a40ed3dSBarry Smith   } else {
125*7adad957SLisandro Dalcin     ierr = MatCreate(((PetscObject)A)->comm,&C);CHKERRQ(ierr);
126f69a0ea3SMatthew Knepley     ierr = MatSetSizes(C,nrows*bs,ncols*bs,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
127*7adad957SLisandro Dalcin     ierr = MatSetType(C,((PetscObject)A)->type_name);CHKERRQ(ierr);
128ab93d7beSBarry Smith     ierr = MatSeqBAIJSetPreallocation_SeqBAIJ(C,bs,0,lens);CHKERRQ(ierr);
129736121d4SSatish Balay   }
130736121d4SSatish Balay   c = (Mat_SeqBAIJ *)(C->data);
131736121d4SSatish Balay   for (i=0; i<nrows; i++) {
132736121d4SSatish Balay     row    = irow[i];
133736121d4SSatish Balay     kstart = ai[row];
134736121d4SSatish Balay     kend   = kstart + a->ilen[row];
135736121d4SSatish Balay     mat_i  = c->i[i];
136736121d4SSatish Balay     mat_j  = c->j + mat_i;
137218c64b6SSatish Balay     mat_a  = c->a + mat_i*bs2;
138736121d4SSatish Balay     mat_ilen = c->ilen + i;
139736121d4SSatish Balay     for (k=kstart; k<kend; k++) {
140736121d4SSatish Balay       if ((tcol=ssmap[a->j[k]])) {
141736121d4SSatish Balay         *mat_j++ = tcol - 1;
142549d3d68SSatish Balay         ierr     = PetscMemcpy(mat_a,a->a+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
143549d3d68SSatish Balay         mat_a   += bs2;
144736121d4SSatish Balay         (*mat_ilen)++;
145736121d4SSatish Balay       }
146736121d4SSatish Balay     }
147736121d4SSatish Balay   }
148218c64b6SSatish Balay 
149736121d4SSatish Balay   /* Free work space */
150736121d4SSatish Balay   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
151606d414cSSatish Balay   ierr = PetscFree(smap);CHKERRQ(ierr);
152606d414cSSatish Balay   ierr = PetscFree(lens);CHKERRQ(ierr);
1536d4a8577SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1546d4a8577SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
155736121d4SSatish Balay 
156736121d4SSatish Balay   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
157736121d4SSatish Balay   *B = C;
1583a40ed3dSBarry Smith   PetscFunctionReturn(0);
159736121d4SSatish Balay }
160736121d4SSatish Balay 
1614a2ae208SSatish Balay #undef __FUNCT__
1624a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqBAIJ"
163690b6cddSBarry Smith PetscErrorCode MatGetSubMatrix_SeqBAIJ(Mat A,IS isrow,IS iscol,PetscInt cs,MatReuse scall,Mat *B)
164218c64b6SSatish Balay {
165218c64b6SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
166218c64b6SSatish Balay   IS             is1,is2;
1676849ba73SBarry Smith   PetscErrorCode ierr;
168899cda47SBarry Smith   PetscInt       *vary,*iary,*irow,*icol,nrows,ncols,i,bs=A->rmap.bs,count;
169218c64b6SSatish Balay 
1703a40ed3dSBarry Smith   PetscFunctionBegin;
171218c64b6SSatish Balay   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
172218c64b6SSatish Balay   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
173b9b97703SBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
174b9b97703SBarry Smith   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
175218c64b6SSatish Balay 
176218c64b6SSatish Balay   /* Verify if the indices corespond to each element in a block
177218c64b6SSatish Balay    and form the IS with compressed IS */
178690b6cddSBarry Smith   ierr = PetscMalloc(2*(a->mbs+1)*sizeof(PetscInt),&vary);CHKERRQ(ierr);
179218c64b6SSatish Balay   iary = vary + a->mbs;
180690b6cddSBarry Smith   ierr = PetscMemzero(vary,(a->mbs)*sizeof(PetscInt));CHKERRQ(ierr);
181218c64b6SSatish Balay   for (i=0; i<nrows; i++) vary[irow[i]/bs]++;
182218c64b6SSatish Balay   count = 0;
183218c64b6SSatish Balay   for (i=0; i<a->mbs; i++) {
184634064b4SBarry Smith     if (vary[i]!=0 && vary[i]!=bs) SETERRQ(PETSC_ERR_ARG_SIZ,"Index set does not match blocks");
185218c64b6SSatish Balay     if (vary[i]==bs) iary[count++] = i;
186218c64b6SSatish Balay   }
187029af93fSBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,count,iary,&is1);CHKERRQ(ierr);
188218c64b6SSatish Balay 
189690b6cddSBarry Smith   ierr = PetscMemzero(vary,(a->mbs)*sizeof(PetscInt));CHKERRQ(ierr);
190218c64b6SSatish Balay   for (i=0; i<ncols; i++) vary[icol[i]/bs]++;
191218c64b6SSatish Balay   count = 0;
192218c64b6SSatish Balay   for (i=0; i<a->mbs; i++) {
193634064b4SBarry Smith     if (vary[i]!=0 && vary[i]!=bs) SETERRQ(PETSC_ERR_PLIB,"Internal error in PETSc");
194218c64b6SSatish Balay     if (vary[i]==bs) iary[count++] = i;
195218c64b6SSatish Balay   }
196029af93fSBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,count,iary,&is2);CHKERRQ(ierr);
197218c64b6SSatish Balay   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
198218c64b6SSatish Balay   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
199606d414cSSatish Balay   ierr = PetscFree(vary);CHKERRQ(ierr);
200218c64b6SSatish Balay 
2016a6a5d1dSBarry Smith   ierr = MatGetSubMatrix_SeqBAIJ_Private(A,is1,is2,cs,scall,B);CHKERRQ(ierr);
202218c64b6SSatish Balay   ISDestroy(is1);
203218c64b6SSatish Balay   ISDestroy(is2);
2043a40ed3dSBarry Smith   PetscFunctionReturn(0);
205218c64b6SSatish Balay }
206218c64b6SSatish Balay 
2074a2ae208SSatish Balay #undef __FUNCT__
2084a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_SeqBAIJ"
209690b6cddSBarry Smith PetscErrorCode MatGetSubMatrices_SeqBAIJ(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
210736121d4SSatish Balay {
2116849ba73SBarry Smith   PetscErrorCode ierr;
212690b6cddSBarry Smith   PetscInt       i;
213736121d4SSatish Balay 
2143a40ed3dSBarry Smith   PetscFunctionBegin;
215736121d4SSatish Balay   if (scall == MAT_INITIAL_MATRIX) {
21682502324SSatish Balay     ierr = PetscMalloc((n+1)*sizeof(Mat),B);CHKERRQ(ierr);
217736121d4SSatish Balay   }
218736121d4SSatish Balay 
219736121d4SSatish Balay   for (i=0; i<n; i++) {
2206a6a5d1dSBarry Smith     ierr = MatGetSubMatrix_SeqBAIJ(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr);
221736121d4SSatish Balay   }
2223a40ed3dSBarry Smith   PetscFunctionReturn(0);
223736121d4SSatish Balay }
224218c64b6SSatish Balay 
225218c64b6SSatish Balay 
2262d61bbb3SSatish Balay /* -------------------------------------------------------*/
2272d61bbb3SSatish Balay /* Should check that shapes of vectors and matrices match */
2282d61bbb3SSatish Balay /* -------------------------------------------------------*/
229d9eff348SSatish Balay #include "petscblaslapack.h"
2302d61bbb3SSatish Balay 
2314a2ae208SSatish Balay #undef __FUNCT__
2324a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_1"
233dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_1(Mat A,Vec xx,Vec zz)
2342d61bbb3SSatish Balay {
2352d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
23687828ca2SBarry Smith   PetscScalar    *x,*z,sum;
2373f1db9ecSBarry Smith   MatScalar      *v;
2386849ba73SBarry Smith   PetscErrorCode ierr;
2397b2bb3b9SHong Zhang   PetscInt       mbs,i,*idx,*ii,n,*ridx=PETSC_NULL;
24026e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
2412d61bbb3SSatish Balay 
2422d61bbb3SSatish Balay   PetscFunctionBegin;
2431ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2441ebc52fbSHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
2452d61bbb3SSatish Balay 
2462d61bbb3SSatish Balay   idx = a->j;
2472d61bbb3SSatish Balay   v   = a->a;
24826e093fcSHong Zhang   if (usecprow){
24926e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
25026e093fcSHong Zhang     ii   = a->compressedrow.i;
2517b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
25226e093fcSHong Zhang   } else {
25326e093fcSHong Zhang     mbs = a->mbs;
2542d61bbb3SSatish Balay     ii  = a->i;
25526e093fcSHong Zhang   }
2562d61bbb3SSatish Balay 
2572d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
2582d61bbb3SSatish Balay     n    = ii[1] - ii[0]; ii++;
2592d61bbb3SSatish Balay     sum  = 0.0;
2602d61bbb3SSatish Balay     while (n--) sum += *v++ * x[*idx++];
26126e093fcSHong Zhang     if (usecprow){
2627b2bb3b9SHong Zhang       z[ridx[i]] = sum;
26326e093fcSHong Zhang     } else {
2642d61bbb3SSatish Balay       z[i] = sum;
2652d61bbb3SSatish Balay     }
26626e093fcSHong Zhang   }
2671ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2681ebc52fbSHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
269efee365bSSatish Balay   ierr = PetscLogFlops(2*a->nz - mbs);CHKERRQ(ierr);
2702d61bbb3SSatish Balay   PetscFunctionReturn(0);
2712d61bbb3SSatish Balay }
2722d61bbb3SSatish Balay 
2734a2ae208SSatish Balay #undef __FUNCT__
2744a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_2"
275dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_2(Mat A,Vec xx,Vec zz)
2762d61bbb3SSatish Balay {
2772d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
278dcf5cc72SBarry Smith   PetscScalar    *x,*z = 0,*xb,sum1,sum2,*zarray;
27987828ca2SBarry Smith   PetscScalar    x1,x2;
2803f1db9ecSBarry Smith   MatScalar      *v;
281dfbe8321SBarry Smith   PetscErrorCode ierr;
2827b2bb3b9SHong Zhang   PetscInt       mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL;
28326e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
2842d61bbb3SSatish Balay 
2852d61bbb3SSatish Balay   PetscFunctionBegin;
2861ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
28726e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
2882d61bbb3SSatish Balay 
2892d61bbb3SSatish Balay   idx = a->j;
2902d61bbb3SSatish Balay   v   = a->a;
29126e093fcSHong Zhang   if (usecprow){
29226e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
29326e093fcSHong Zhang     ii   = a->compressedrow.i;
2947b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
29526e093fcSHong Zhang   } else {
29626e093fcSHong Zhang     mbs = a->mbs;
2972d61bbb3SSatish Balay     ii  = a->i;
29826e093fcSHong Zhang     z   = zarray;
29926e093fcSHong Zhang   }
3002d61bbb3SSatish Balay 
3012d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
3022d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
3032d61bbb3SSatish Balay     sum1 = 0.0; sum2 = 0.0;
3042d61bbb3SSatish Balay     for (j=0; j<n; j++) {
3052d61bbb3SSatish Balay       xb = x + 2*(*idx++); x1 = xb[0]; x2 = xb[1];
3062d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[2]*x2;
3072d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[3]*x2;
3082d61bbb3SSatish Balay       v += 4;
3092d61bbb3SSatish Balay     }
3107b2bb3b9SHong Zhang     if (usecprow) z = zarray + 2*ridx[i];
3112d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2;
31226e093fcSHong Zhang     if (!usecprow) z += 2;
3132d61bbb3SSatish Balay   }
3141ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
31526e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
316efee365bSSatish Balay   ierr = PetscLogFlops(8*a->nz - 2*mbs);CHKERRQ(ierr);
3172d61bbb3SSatish Balay   PetscFunctionReturn(0);
3182d61bbb3SSatish Balay }
3192d61bbb3SSatish Balay 
3204a2ae208SSatish Balay #undef __FUNCT__
3214a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_3"
322dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_3(Mat A,Vec xx,Vec zz)
3232d61bbb3SSatish Balay {
3242d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
325dcf5cc72SBarry Smith   PetscScalar    *x,*z = 0,*xb,sum1,sum2,sum3,x1,x2,x3,*zarray;
3263f1db9ecSBarry Smith   MatScalar      *v;
327dfbe8321SBarry Smith   PetscErrorCode ierr;
3287b2bb3b9SHong Zhang   PetscInt       mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL;
32926e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
33026e093fcSHong Zhang 
3312d61bbb3SSatish Balay 
332b6410449SSatish Balay #if defined(PETSC_HAVE_PRAGMA_DISJOINT)
333fee21e36SBarry Smith #pragma disjoint(*v,*z,*xb)
334fee21e36SBarry Smith #endif
335fee21e36SBarry Smith 
3362d61bbb3SSatish Balay   PetscFunctionBegin;
3371ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
33826e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
3392d61bbb3SSatish Balay 
3402d61bbb3SSatish Balay   idx = a->j;
3412d61bbb3SSatish Balay   v   = a->a;
34226e093fcSHong Zhang   if (usecprow){
34326e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
34426e093fcSHong Zhang     ii   = a->compressedrow.i;
3457b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
34626e093fcSHong Zhang   } else {
34726e093fcSHong Zhang     mbs = a->mbs;
3482d61bbb3SSatish Balay     ii  = a->i;
34926e093fcSHong Zhang     z   = zarray;
35026e093fcSHong Zhang   }
3512d61bbb3SSatish Balay 
3522d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
3532d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
3542d61bbb3SSatish Balay     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0;
3552d61bbb3SSatish Balay     for (j=0; j<n; j++) {
3562d61bbb3SSatish Balay       xb = x + 3*(*idx++); x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
3572d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[3]*x2 + v[6]*x3;
3582d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[4]*x2 + v[7]*x3;
3592d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[5]*x2 + v[8]*x3;
3602d61bbb3SSatish Balay       v += 9;
3612d61bbb3SSatish Balay     }
3627b2bb3b9SHong Zhang     if (usecprow) z = zarray + 3*ridx[i];
3632d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3;
36426e093fcSHong Zhang     if (!usecprow) z += 3;
3652d61bbb3SSatish Balay   }
3661ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
36726e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
368efee365bSSatish Balay   ierr = PetscLogFlops(18*a->nz - 3*mbs);CHKERRQ(ierr);
3692d61bbb3SSatish Balay   PetscFunctionReturn(0);
3702d61bbb3SSatish Balay }
3712d61bbb3SSatish Balay 
3724a2ae208SSatish Balay #undef __FUNCT__
3734a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_4"
374dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_4(Mat A,Vec xx,Vec zz)
3752d61bbb3SSatish Balay {
3762d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
377dcf5cc72SBarry Smith   PetscScalar    *x,*z = 0,*xb,sum1,sum2,sum3,sum4,x1,x2,x3,x4,*zarray;
3783f1db9ecSBarry Smith   MatScalar      *v;
379dfbe8321SBarry Smith   PetscErrorCode ierr;
3807b2bb3b9SHong Zhang   PetscInt       mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL;
38126e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
3822d61bbb3SSatish Balay 
3832d61bbb3SSatish Balay   PetscFunctionBegin;
3841ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
38526e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
3862d61bbb3SSatish Balay 
3872d61bbb3SSatish Balay   idx = a->j;
3882d61bbb3SSatish Balay   v   = a->a;
38926e093fcSHong Zhang   if (usecprow){
39026e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
39126e093fcSHong Zhang     ii   = a->compressedrow.i;
3927b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
39326e093fcSHong Zhang   } else {
39426e093fcSHong Zhang     mbs = a->mbs;
3952d61bbb3SSatish Balay     ii  = a->i;
39626e093fcSHong Zhang     z   = zarray;
39726e093fcSHong Zhang   }
3982d61bbb3SSatish Balay 
3992d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
4002d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
4012d61bbb3SSatish Balay     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0;
4022d61bbb3SSatish Balay     for (j=0; j<n; j++) {
4032d61bbb3SSatish Balay       xb = x + 4*(*idx++);
4042d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
4052d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[4]*x2 + v[8]*x3   + v[12]*x4;
4062d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[5]*x2 + v[9]*x3   + v[13]*x4;
4072d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[6]*x2 + v[10]*x3  + v[14]*x4;
4082d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[7]*x2 + v[11]*x3  + v[15]*x4;
4092d61bbb3SSatish Balay       v += 16;
4102d61bbb3SSatish Balay     }
4117b2bb3b9SHong Zhang     if (usecprow) z = zarray + 4*ridx[i];
4122d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4;
41326e093fcSHong Zhang     if (!usecprow) z += 4;
4142d61bbb3SSatish Balay   }
4151ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
41626e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
417efee365bSSatish Balay   ierr = PetscLogFlops(32*a->nz - 4*mbs);CHKERRQ(ierr);
4182d61bbb3SSatish Balay   PetscFunctionReturn(0);
4192d61bbb3SSatish Balay }
4202d61bbb3SSatish Balay 
4214a2ae208SSatish Balay #undef __FUNCT__
4224a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_5"
423dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_5(Mat A,Vec xx,Vec zz)
4242d61bbb3SSatish Balay {
4252d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
426dcf5cc72SBarry Smith   PetscScalar    sum1,sum2,sum3,sum4,sum5,x1,x2,x3,x4,x5,*xb,*z = 0,*x,*zarray;
4273f1db9ecSBarry Smith   MatScalar      *v;
428dfbe8321SBarry Smith   PetscErrorCode ierr;
4297b2bb3b9SHong Zhang   PetscInt       mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL;
43026e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
4312d61bbb3SSatish Balay 
432433994e6SBarry Smith   PetscFunctionBegin;
4331ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
43426e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
4352d61bbb3SSatish Balay 
4362d61bbb3SSatish Balay   idx = a->j;
4372d61bbb3SSatish Balay   v   = a->a;
43826e093fcSHong Zhang   if (usecprow){
43926e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
44026e093fcSHong Zhang     ii   = a->compressedrow.i;
4417b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
44226e093fcSHong Zhang   } else {
44326e093fcSHong Zhang     mbs = a->mbs;
4442d61bbb3SSatish Balay     ii  = a->i;
44526e093fcSHong Zhang     z   = zarray;
44626e093fcSHong Zhang   }
4472d61bbb3SSatish Balay 
4482d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
4492d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
4502d61bbb3SSatish Balay     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0;
4512d61bbb3SSatish Balay     for (j=0; j<n; j++) {
4522d61bbb3SSatish Balay       xb = x + 5*(*idx++);
4532d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4];
4542d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[5]*x2 + v[10]*x3  + v[15]*x4 + v[20]*x5;
4552d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[6]*x2 + v[11]*x3  + v[16]*x4 + v[21]*x5;
4562d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[7]*x2 + v[12]*x3  + v[17]*x4 + v[22]*x5;
4572d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[8]*x2 + v[13]*x3  + v[18]*x4 + v[23]*x5;
4582d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[9]*x2 + v[14]*x3  + v[19]*x4 + v[24]*x5;
4592d61bbb3SSatish Balay       v += 25;
4602d61bbb3SSatish Balay     }
4617b2bb3b9SHong Zhang     if (usecprow) z = zarray + 5*ridx[i];
4622d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5;
46326e093fcSHong Zhang     if (!usecprow) z += 5;
4642d61bbb3SSatish Balay   }
4651ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
46626e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
467efee365bSSatish Balay   ierr = PetscLogFlops(50*a->nz - 5*mbs);CHKERRQ(ierr);
4682d61bbb3SSatish Balay   PetscFunctionReturn(0);
4692d61bbb3SSatish Balay }
4702d61bbb3SSatish Balay 
47115091d37SBarry Smith 
4724a2ae208SSatish Balay #undef __FUNCT__
4734a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_6"
474dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_6(Mat A,Vec xx,Vec zz)
47515091d37SBarry Smith {
47615091d37SBarry Smith   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
477dcf5cc72SBarry Smith   PetscScalar    *x,*z = 0,*xb,sum1,sum2,sum3,sum4,sum5,sum6;
47826e093fcSHong Zhang   PetscScalar    x1,x2,x3,x4,x5,x6,*zarray;
47915091d37SBarry Smith   MatScalar      *v;
480dfbe8321SBarry Smith   PetscErrorCode ierr;
4817b2bb3b9SHong Zhang   PetscInt       mbs=a->mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL;
48226e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
48315091d37SBarry Smith 
484433994e6SBarry Smith   PetscFunctionBegin;
4851ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
48626e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
48715091d37SBarry Smith 
48815091d37SBarry Smith   idx = a->j;
48915091d37SBarry Smith   v   = a->a;
49026e093fcSHong Zhang   if (usecprow){
49126e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
49226e093fcSHong Zhang     ii   = a->compressedrow.i;
4937b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
49426e093fcSHong Zhang   } else {
49526e093fcSHong Zhang     mbs = a->mbs;
49615091d37SBarry Smith     ii  = a->i;
49726e093fcSHong Zhang     z   = zarray;
49826e093fcSHong Zhang   }
49915091d37SBarry Smith 
50015091d37SBarry Smith   for (i=0; i<mbs; i++) {
50115091d37SBarry Smith     n  = ii[1] - ii[0]; ii++;
50215091d37SBarry Smith     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0;
50315091d37SBarry Smith     for (j=0; j<n; j++) {
50415091d37SBarry Smith       xb = x + 6*(*idx++);
50515091d37SBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5];
50615091d37SBarry Smith       sum1 += v[0]*x1 + v[6]*x2  + v[12]*x3  + v[18]*x4 + v[24]*x5 + v[30]*x6;
50715091d37SBarry Smith       sum2 += v[1]*x1 + v[7]*x2  + v[13]*x3  + v[19]*x4 + v[25]*x5 + v[31]*x6;
50815091d37SBarry Smith       sum3 += v[2]*x1 + v[8]*x2  + v[14]*x3  + v[20]*x4 + v[26]*x5 + v[32]*x6;
50915091d37SBarry Smith       sum4 += v[3]*x1 + v[9]*x2  + v[15]*x3  + v[21]*x4 + v[27]*x5 + v[33]*x6;
51015091d37SBarry Smith       sum5 += v[4]*x1 + v[10]*x2 + v[16]*x3  + v[22]*x4 + v[28]*x5 + v[34]*x6;
51115091d37SBarry Smith       sum6 += v[5]*x1 + v[11]*x2 + v[17]*x3  + v[23]*x4 + v[29]*x5 + v[35]*x6;
51215091d37SBarry Smith       v += 36;
51315091d37SBarry Smith     }
5147b2bb3b9SHong Zhang     if (usecprow) z = zarray + 6*ridx[i];
51515091d37SBarry Smith     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6;
51626e093fcSHong Zhang     if (!usecprow) z += 6;
51715091d37SBarry Smith   }
51815091d37SBarry Smith 
5191ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
52026e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
521efee365bSSatish Balay   ierr = PetscLogFlops(72*a->nz - 6*mbs);CHKERRQ(ierr);
52215091d37SBarry Smith   PetscFunctionReturn(0);
52315091d37SBarry Smith }
5244a2ae208SSatish Balay #undef __FUNCT__
5254a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_7"
526dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_7(Mat A,Vec xx,Vec zz)
5272d61bbb3SSatish Balay {
5282d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
529dcf5cc72SBarry Smith   PetscScalar    *x,*z = 0,*xb,sum1,sum2,sum3,sum4,sum5,sum6,sum7;
53026e093fcSHong Zhang   PetscScalar    x1,x2,x3,x4,x5,x6,x7,*zarray;
5313f1db9ecSBarry Smith   MatScalar      *v;
532dfbe8321SBarry Smith   PetscErrorCode ierr;
5337b2bb3b9SHong Zhang   PetscInt       mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL;
53426e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
5352d61bbb3SSatish Balay 
536433994e6SBarry Smith   PetscFunctionBegin;
5371ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
53826e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
5392d61bbb3SSatish Balay 
5402d61bbb3SSatish Balay   idx = a->j;
5412d61bbb3SSatish Balay   v   = a->a;
54226e093fcSHong Zhang   if (usecprow){
54326e093fcSHong Zhang     mbs    = a->compressedrow.nrows;
54426e093fcSHong Zhang     ii     = a->compressedrow.i;
5457b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
54626e093fcSHong Zhang   } else {
54726e093fcSHong Zhang     mbs = a->mbs;
5482d61bbb3SSatish Balay     ii  = a->i;
54926e093fcSHong Zhang     z   = zarray;
55026e093fcSHong Zhang   }
5512d61bbb3SSatish Balay 
5522d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
5532d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
5542d61bbb3SSatish Balay     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0; sum7 = 0.0;
5552d61bbb3SSatish Balay     for (j=0; j<n; j++) {
5562d61bbb3SSatish Balay       xb = x + 7*(*idx++);
5572d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6];
5582d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[7]*x2  + v[14]*x3  + v[21]*x4 + v[28]*x5 + v[35]*x6 + v[42]*x7;
5592d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[8]*x2  + v[15]*x3  + v[22]*x4 + v[29]*x5 + v[36]*x6 + v[43]*x7;
5602d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[9]*x2  + v[16]*x3  + v[23]*x4 + v[30]*x5 + v[37]*x6 + v[44]*x7;
5612d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[10]*x2 + v[17]*x3  + v[24]*x4 + v[31]*x5 + v[38]*x6 + v[45]*x7;
5622d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[11]*x2 + v[18]*x3  + v[25]*x4 + v[32]*x5 + v[39]*x6 + v[46]*x7;
5632d61bbb3SSatish Balay       sum6 += v[5]*x1 + v[12]*x2 + v[19]*x3  + v[26]*x4 + v[33]*x5 + v[40]*x6 + v[47]*x7;
5642d61bbb3SSatish Balay       sum7 += v[6]*x1 + v[13]*x2 + v[20]*x3  + v[27]*x4 + v[34]*x5 + v[41]*x6 + v[48]*x7;
5652d61bbb3SSatish Balay       v += 49;
5662d61bbb3SSatish Balay     }
5677b2bb3b9SHong Zhang     if (usecprow) z = zarray + 7*ridx[i];
5682d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
56926e093fcSHong Zhang     if (!usecprow) z += 7;
5702d61bbb3SSatish Balay   }
5712d61bbb3SSatish Balay 
5721ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
57326e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
574efee365bSSatish Balay   ierr = PetscLogFlops(98*a->nz - 7*mbs);CHKERRQ(ierr);
5752d61bbb3SSatish Balay   PetscFunctionReturn(0);
5762d61bbb3SSatish Balay }
5772d61bbb3SSatish Balay 
5783f1db9ecSBarry Smith /*
5793f1db9ecSBarry Smith     This will not work with MatScalar == float because it calls the BLAS
5803f1db9ecSBarry Smith */
5814a2ae208SSatish Balay #undef __FUNCT__
5824a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_N"
583dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_N(Mat A,Vec xx,Vec zz)
5842d61bbb3SSatish Balay {
5852d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
586dcf5cc72SBarry Smith   PetscScalar    *x,*z = 0,*xb,*work,*workt,*zarray;
5873f1db9ecSBarry Smith   MatScalar      *v;
588dfbe8321SBarry Smith   PetscErrorCode ierr;
589899cda47SBarry Smith   PetscInt       mbs=a->mbs,i,*idx,*ii,bs=A->rmap.bs,j,n,bs2=a->bs2;
5907b2bb3b9SHong Zhang   PetscInt       ncols,k,*ridx=PETSC_NULL;
59126e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
5922d61bbb3SSatish Balay 
5932d61bbb3SSatish Balay   PetscFunctionBegin;
5941ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
59526e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
5962d61bbb3SSatish Balay 
5972d61bbb3SSatish Balay   idx = a->j;
5982d61bbb3SSatish Balay   v   = a->a;
59926e093fcSHong Zhang   if (usecprow){
60026e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
60126e093fcSHong Zhang     ii   = a->compressedrow.i;
6027b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
60326e093fcSHong Zhang   } else {
60426e093fcSHong Zhang     mbs = a->mbs;
6052d61bbb3SSatish Balay     ii  = a->i;
60626e093fcSHong Zhang     z   = zarray;
60726e093fcSHong Zhang   }
608218c64b6SSatish Balay 
6092d61bbb3SSatish Balay   if (!a->mult_work) {
610899cda47SBarry Smith     k    = PetscMax(A->rmap.n,A->cmap.n);
61187828ca2SBarry Smith     ierr = PetscMalloc((k+1)*sizeof(PetscScalar),&a->mult_work);CHKERRQ(ierr);
6122d61bbb3SSatish Balay   }
6132d61bbb3SSatish Balay   work = a->mult_work;
6142d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
6152d61bbb3SSatish Balay     n     = ii[1] - ii[0]; ii++;
6162d61bbb3SSatish Balay     ncols = n*bs;
6172d61bbb3SSatish Balay     workt = work;
6182d61bbb3SSatish Balay     for (j=0; j<n; j++) {
6192d61bbb3SSatish Balay       xb = x + bs*(*idx++);
6202d61bbb3SSatish Balay       for (k=0; k<bs; k++) workt[k] = xb[k];
6212d61bbb3SSatish Balay       workt += bs;
6222d61bbb3SSatish Balay     }
6237b2bb3b9SHong Zhang     if (usecprow) z = zarray + bs*ridx[i];
624737d121aSSatish Balay     Kernel_w_gets_Ar_times_v(bs,ncols,work,v,z);
62571044d3cSBarry Smith     /* BLASgemv_("N",&bs,&ncols,&_DOne,v,&bs,work,&_One,&_DZero,z,&_One); */
6262d61bbb3SSatish Balay     v += n*bs2;
62726e093fcSHong Zhang     if (!usecprow) z += bs;
6282d61bbb3SSatish Balay   }
6291ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
63026e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
631efee365bSSatish Balay   ierr = PetscLogFlops(2*a->nz*bs2 - bs*mbs);CHKERRQ(ierr);
6322d61bbb3SSatish Balay   PetscFunctionReturn(0);
6332d61bbb3SSatish Balay }
6342d61bbb3SSatish Balay 
6356a65c766SHong Zhang extern PetscErrorCode VecCopy_Seq(Vec,Vec);
6364a2ae208SSatish Balay #undef __FUNCT__
6374a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_1"
638dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_1(Mat A,Vec xx,Vec yy,Vec zz)
6392d61bbb3SSatish Balay {
6402d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
641dcf5cc72SBarry Smith   PetscScalar    *x,*y = 0,*z = 0,sum,*yarray,*zarray;
6423f1db9ecSBarry Smith   MatScalar      *v;
643dfbe8321SBarry Smith   PetscErrorCode ierr;
6444eb6d288SHong Zhang   PetscInt       mbs=a->mbs,i,*idx,*ii,n,*ridx=PETSC_NULL;
64526e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
6462d61bbb3SSatish Balay 
6472d61bbb3SSatish Balay   PetscFunctionBegin;
6481ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
64926e093fcSHong Zhang   ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr);
6502e8a6d31SBarry Smith   if (zz != yy) {
65126e093fcSHong Zhang     ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
6522e8a6d31SBarry Smith   } else {
65326e093fcSHong Zhang     zarray = yarray;
6542e8a6d31SBarry Smith   }
6552d61bbb3SSatish Balay 
6562d61bbb3SSatish Balay   idx = a->j;
6572d61bbb3SSatish Balay   v   = a->a;
65826e093fcSHong Zhang   if (usecprow){
6594eb6d288SHong Zhang     if (zz != yy){
6604eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,mbs*sizeof(PetscScalar));CHKERRQ(ierr);
6614eb6d288SHong Zhang     }
66226e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
66326e093fcSHong Zhang     ii   = a->compressedrow.i;
6647b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
66526e093fcSHong Zhang   } else {
6662d61bbb3SSatish Balay     ii  = a->i;
66726e093fcSHong Zhang     y   = yarray;
66826e093fcSHong Zhang     z   = zarray;
66926e093fcSHong Zhang   }
6702d61bbb3SSatish Balay 
6712d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
6722d61bbb3SSatish Balay     n    = ii[1] - ii[0]; ii++;
67326e093fcSHong Zhang     if (usecprow){
6747b2bb3b9SHong Zhang       z = zarray + ridx[i];
6757b2bb3b9SHong Zhang       y = yarray + ridx[i];
67626e093fcSHong Zhang     }
67726e093fcSHong Zhang     sum = y[0];
6782d61bbb3SSatish Balay     while (n--) sum += *v++ * x[*idx++];
67926e093fcSHong Zhang     z[0] = sum;
68026e093fcSHong Zhang     if (!usecprow){
68126e093fcSHong Zhang       z++; y++;
68226e093fcSHong Zhang     }
6832d61bbb3SSatish Balay   }
6841ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
68526e093fcSHong Zhang   ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr);
6862e8a6d31SBarry Smith   if (zz != yy) {
68726e093fcSHong Zhang     ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
6882e8a6d31SBarry Smith   }
689efee365bSSatish Balay   ierr = PetscLogFlops(2*a->nz);CHKERRQ(ierr);
6902d61bbb3SSatish Balay   PetscFunctionReturn(0);
6912d61bbb3SSatish Balay }
6922d61bbb3SSatish Balay 
6934a2ae208SSatish Balay #undef __FUNCT__
6944a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_2"
695dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_2(Mat A,Vec xx,Vec yy,Vec zz)
6962d61bbb3SSatish Balay {
6972d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
698dcf5cc72SBarry Smith   PetscScalar    *x,*y = 0,*z = 0,*xb,sum1,sum2;
69926e093fcSHong Zhang   PetscScalar    x1,x2,*yarray,*zarray;
7003f1db9ecSBarry Smith   MatScalar      *v;
701dfbe8321SBarry Smith   PetscErrorCode ierr;
7024eb6d288SHong Zhang   PetscInt       mbs=a->mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL;
70326e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
7042d61bbb3SSatish Balay 
7052d61bbb3SSatish Balay   PetscFunctionBegin;
7061ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
70726e093fcSHong Zhang   ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr);
7082e8a6d31SBarry Smith   if (zz != yy) {
70926e093fcSHong Zhang     ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
7102e8a6d31SBarry Smith   } else {
71126e093fcSHong Zhang     zarray = yarray;
7122e8a6d31SBarry Smith   }
7132d61bbb3SSatish Balay 
7142d61bbb3SSatish Balay   idx = a->j;
7152d61bbb3SSatish Balay   v   = a->a;
71626e093fcSHong Zhang   if (usecprow){
7174eb6d288SHong Zhang     if (zz != yy){
7184eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,2*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
7194eb6d288SHong Zhang     }
72026e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
72126e093fcSHong Zhang     ii   = a->compressedrow.i;
7227b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
7234eb6d288SHong Zhang     if (zz != yy){
7244eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,a->mbs*sizeof(PetscScalar));CHKERRQ(ierr);
7254eb6d288SHong Zhang     }
72626e093fcSHong Zhang   } else {
7272d61bbb3SSatish Balay     ii  = a->i;
72826e093fcSHong Zhang     y   = yarray;
72926e093fcSHong Zhang     z   = zarray;
73026e093fcSHong Zhang   }
7312d61bbb3SSatish Balay 
7322d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
7332d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
73426e093fcSHong Zhang     if (usecprow){
7357b2bb3b9SHong Zhang       z = zarray + 2*ridx[i];
7367b2bb3b9SHong Zhang       y = yarray + 2*ridx[i];
73726e093fcSHong Zhang     }
7382d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1];
7392d61bbb3SSatish Balay     for (j=0; j<n; j++) {
7402d61bbb3SSatish Balay       xb = x + 2*(*idx++); x1 = xb[0]; x2 = xb[1];
7412d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[2]*x2;
7422d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[3]*x2;
7432d61bbb3SSatish Balay       v += 4;
7442d61bbb3SSatish Balay     }
7452d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2;
74626e093fcSHong Zhang     if (!usecprow){
7472d61bbb3SSatish Balay       z += 2; y += 2;
7482d61bbb3SSatish Balay     }
74926e093fcSHong Zhang   }
7501ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
75126e093fcSHong Zhang   ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr);
7522e8a6d31SBarry Smith   if (zz != yy) {
75326e093fcSHong Zhang     ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
7542e8a6d31SBarry Smith   }
755efee365bSSatish Balay   ierr = PetscLogFlops(4*a->nz);CHKERRQ(ierr);
7562d61bbb3SSatish Balay   PetscFunctionReturn(0);
7572d61bbb3SSatish Balay }
7582d61bbb3SSatish Balay 
7594a2ae208SSatish Balay #undef __FUNCT__
7604a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_3"
761dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_3(Mat A,Vec xx,Vec yy,Vec zz)
7622d61bbb3SSatish Balay {
7632d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
764dcf5cc72SBarry Smith   PetscScalar    *x,*y = 0,*z = 0,*xb,sum1,sum2,sum3,x1,x2,x3,*yarray,*zarray;
7653f1db9ecSBarry Smith   MatScalar      *v;
766dfbe8321SBarry Smith   PetscErrorCode ierr;
7674eb6d288SHong Zhang   PetscInt       mbs=a->mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL;
76826e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
7692d61bbb3SSatish Balay 
7702d61bbb3SSatish Balay   PetscFunctionBegin;
7711ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
77226e093fcSHong Zhang   ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr);
7732e8a6d31SBarry Smith   if (zz != yy) {
77426e093fcSHong Zhang     ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
7752e8a6d31SBarry Smith   } else {
77626e093fcSHong Zhang     zarray = yarray;
7772e8a6d31SBarry Smith   }
7782d61bbb3SSatish Balay 
7792d61bbb3SSatish Balay   idx = a->j;
7802d61bbb3SSatish Balay   v   = a->a;
78126e093fcSHong Zhang   if (usecprow){
7824eb6d288SHong Zhang     if (zz != yy){
7834eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,3*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
7844eb6d288SHong Zhang     }
78526e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
78626e093fcSHong Zhang     ii   = a->compressedrow.i;
7877b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
78826e093fcSHong Zhang   } else {
7892d61bbb3SSatish Balay     ii  = a->i;
79026e093fcSHong Zhang     y   = yarray;
79126e093fcSHong Zhang     z   = zarray;
79226e093fcSHong Zhang   }
7932d61bbb3SSatish Balay 
7942d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
7952d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
79626e093fcSHong Zhang     if (usecprow){
7977b2bb3b9SHong Zhang       z = zarray + 3*ridx[i];
7987b2bb3b9SHong Zhang       y = yarray + 3*ridx[i];
79926e093fcSHong Zhang     }
8002d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2];
8012d61bbb3SSatish Balay     for (j=0; j<n; j++) {
8022d61bbb3SSatish Balay       xb = x + 3*(*idx++); x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
8032d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[3]*x2 + v[6]*x3;
8042d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[4]*x2 + v[7]*x3;
8052d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[5]*x2 + v[8]*x3;
8062d61bbb3SSatish Balay       v += 9;
8072d61bbb3SSatish Balay     }
8082d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3;
80926e093fcSHong Zhang     if (!usecprow){
8102d61bbb3SSatish Balay       z += 3; y += 3;
8112d61bbb3SSatish Balay     }
81226e093fcSHong Zhang   }
8131ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
81426e093fcSHong Zhang   ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr);
8152e8a6d31SBarry Smith   if (zz != yy) {
81626e093fcSHong Zhang     ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
8172e8a6d31SBarry Smith   }
818efee365bSSatish Balay   ierr = PetscLogFlops(18*a->nz);CHKERRQ(ierr);
8192d61bbb3SSatish Balay   PetscFunctionReturn(0);
8202d61bbb3SSatish Balay }
8212d61bbb3SSatish Balay 
8224a2ae208SSatish Balay #undef __FUNCT__
8234a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_4"
824dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_4(Mat A,Vec xx,Vec yy,Vec zz)
8252d61bbb3SSatish Balay {
8262d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
827dcf5cc72SBarry Smith   PetscScalar    *x,*y = 0,*z = 0,*xb,sum1,sum2,sum3,sum4,x1,x2,x3,x4,*yarray,*zarray;
8283f1db9ecSBarry Smith   MatScalar      *v;
829dfbe8321SBarry Smith   PetscErrorCode ierr;
8304eb6d288SHong Zhang   PetscInt       mbs=a->mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL;
83126e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
8322d61bbb3SSatish Balay 
8332d61bbb3SSatish Balay   PetscFunctionBegin;
8341ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
83526e093fcSHong Zhang   ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr);
8362e8a6d31SBarry Smith   if (zz != yy) {
83726e093fcSHong Zhang     ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
8382e8a6d31SBarry Smith   } else {
83926e093fcSHong Zhang     zarray = yarray;
8402e8a6d31SBarry Smith   }
8412d61bbb3SSatish Balay 
8422d61bbb3SSatish Balay   idx   = a->j;
8432d61bbb3SSatish Balay   v     = a->a;
84426e093fcSHong Zhang   if (usecprow){
8454eb6d288SHong Zhang     if (zz != yy){
8464eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,4*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
8474eb6d288SHong Zhang     }
84826e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
84926e093fcSHong Zhang     ii   = a->compressedrow.i;
8507b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
85126e093fcSHong Zhang   } else {
8522d61bbb3SSatish Balay     ii  = a->i;
85326e093fcSHong Zhang     y   = yarray;
85426e093fcSHong Zhang     z   = zarray;
85526e093fcSHong Zhang   }
8562d61bbb3SSatish Balay 
8572d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
8582d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
85926e093fcSHong Zhang     if (usecprow){
8607b2bb3b9SHong Zhang       z = zarray + 4*ridx[i];
8617b2bb3b9SHong Zhang       y = yarray + 4*ridx[i];
86226e093fcSHong Zhang     }
8632d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3];
8642d61bbb3SSatish Balay     for (j=0; j<n; j++) {
8652d61bbb3SSatish Balay       xb = x + 4*(*idx++);
8662d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
8672d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[4]*x2 + v[8]*x3   + v[12]*x4;
8682d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[5]*x2 + v[9]*x3   + v[13]*x4;
8692d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[6]*x2 + v[10]*x3  + v[14]*x4;
8702d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[7]*x2 + v[11]*x3  + v[15]*x4;
8712d61bbb3SSatish Balay       v += 16;
8722d61bbb3SSatish Balay     }
8732d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4;
87426e093fcSHong Zhang     if (!usecprow){
8752d61bbb3SSatish Balay       z += 4; y += 4;
8762d61bbb3SSatish Balay     }
87726e093fcSHong Zhang   }
8781ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
87926e093fcSHong Zhang   ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr);
8802e8a6d31SBarry Smith   if (zz != yy) {
88126e093fcSHong Zhang     ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
8822e8a6d31SBarry Smith   }
883efee365bSSatish Balay   ierr = PetscLogFlops(32*a->nz);CHKERRQ(ierr);
8842d61bbb3SSatish Balay   PetscFunctionReturn(0);
8852d61bbb3SSatish Balay }
8862d61bbb3SSatish Balay 
8874a2ae208SSatish Balay #undef __FUNCT__
8884a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_5"
889dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_5(Mat A,Vec xx,Vec yy,Vec zz)
8902d61bbb3SSatish Balay {
8912d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
892dcf5cc72SBarry Smith   PetscScalar    *x,*y = 0,*z = 0,*xb,sum1,sum2,sum3,sum4,sum5,x1,x2,x3,x4,x5;
89326e093fcSHong Zhang   PetscScalar    *yarray,*zarray;
8943f1db9ecSBarry Smith   MatScalar      *v;
895dfbe8321SBarry Smith   PetscErrorCode ierr;
8964eb6d288SHong Zhang   PetscInt       mbs=a->mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL;
89726e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
8982d61bbb3SSatish Balay 
8992d61bbb3SSatish Balay   PetscFunctionBegin;
9001ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
90126e093fcSHong Zhang   ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr);
9022e8a6d31SBarry Smith   if (zz != yy) {
90326e093fcSHong Zhang     ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
9042e8a6d31SBarry Smith   } else {
90526e093fcSHong Zhang     zarray = yarray;
9062e8a6d31SBarry Smith   }
9072d61bbb3SSatish Balay 
9082d61bbb3SSatish Balay   idx = a->j;
9092d61bbb3SSatish Balay   v   = a->a;
91026e093fcSHong Zhang   if (usecprow){
9114eb6d288SHong Zhang     if (zz != yy){
9124eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,5*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
9134eb6d288SHong Zhang     }
91426e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
91526e093fcSHong Zhang     ii   = a->compressedrow.i;
9167b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
91726e093fcSHong Zhang   } else {
9182d61bbb3SSatish Balay     ii  = a->i;
91926e093fcSHong Zhang     y   = yarray;
92026e093fcSHong Zhang     z   = zarray;
92126e093fcSHong Zhang   }
9222d61bbb3SSatish Balay 
9232d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
9242d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
92526e093fcSHong Zhang     if (usecprow){
9267b2bb3b9SHong Zhang       z = zarray + 5*ridx[i];
9277b2bb3b9SHong Zhang       y = yarray + 5*ridx[i];
92826e093fcSHong Zhang     }
9292d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4];
9302d61bbb3SSatish Balay     for (j=0; j<n; j++) {
9312d61bbb3SSatish Balay       xb = x + 5*(*idx++);
9322d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4];
9332d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[5]*x2 + v[10]*x3  + v[15]*x4 + v[20]*x5;
9342d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[6]*x2 + v[11]*x3  + v[16]*x4 + v[21]*x5;
9352d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[7]*x2 + v[12]*x3  + v[17]*x4 + v[22]*x5;
9362d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[8]*x2 + v[13]*x3  + v[18]*x4 + v[23]*x5;
9372d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[9]*x2 + v[14]*x3  + v[19]*x4 + v[24]*x5;
9382d61bbb3SSatish Balay       v += 25;
9392d61bbb3SSatish Balay     }
9402d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5;
94126e093fcSHong Zhang     if (!usecprow){
9422d61bbb3SSatish Balay       z += 5; y += 5;
9432d61bbb3SSatish Balay     }
94426e093fcSHong Zhang   }
9451ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
94626e093fcSHong Zhang   ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr);
9472e8a6d31SBarry Smith   if (zz != yy) {
94826e093fcSHong Zhang     ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
9492e8a6d31SBarry Smith   }
950efee365bSSatish Balay   ierr = PetscLogFlops(50*a->nz);CHKERRQ(ierr);
9512d61bbb3SSatish Balay   PetscFunctionReturn(0);
9522d61bbb3SSatish Balay }
9534a2ae208SSatish Balay #undef __FUNCT__
9544a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_6"
955dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_6(Mat A,Vec xx,Vec yy,Vec zz)
95615091d37SBarry Smith {
95715091d37SBarry Smith   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
958dcf5cc72SBarry Smith   PetscScalar    *x,*y = 0,*z = 0,*xb,sum1,sum2,sum3,sum4,sum5,sum6;
95926e093fcSHong Zhang   PetscScalar    x1,x2,x3,x4,x5,x6,*yarray,*zarray;
96015091d37SBarry Smith   MatScalar      *v;
961dfbe8321SBarry Smith   PetscErrorCode ierr;
9624eb6d288SHong Zhang   PetscInt       mbs=a->mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL;
96326e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
96415091d37SBarry Smith 
96515091d37SBarry Smith   PetscFunctionBegin;
9661ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
96726e093fcSHong Zhang   ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr);
96815091d37SBarry Smith   if (zz != yy) {
96926e093fcSHong Zhang     ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
97015091d37SBarry Smith   } else {
97126e093fcSHong Zhang     zarray = yarray;
97215091d37SBarry Smith   }
97315091d37SBarry Smith 
97415091d37SBarry Smith   idx = a->j;
97515091d37SBarry Smith   v   = a->a;
97626e093fcSHong Zhang   if (usecprow){
9774eb6d288SHong Zhang     if (zz != yy){
9784eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,6*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
9794eb6d288SHong Zhang     }
98026e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
98126e093fcSHong Zhang     ii   = a->compressedrow.i;
9827b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
98326e093fcSHong Zhang   } else {
98415091d37SBarry Smith     ii  = a->i;
98526e093fcSHong Zhang     y   = yarray;
98626e093fcSHong Zhang     z   = zarray;
98726e093fcSHong Zhang   }
98815091d37SBarry Smith 
98915091d37SBarry Smith   for (i=0; i<mbs; i++) {
99015091d37SBarry Smith     n  = ii[1] - ii[0]; ii++;
99126e093fcSHong Zhang     if (usecprow){
9927b2bb3b9SHong Zhang       z = zarray + 6*ridx[i];
9937b2bb3b9SHong Zhang       y = yarray + 6*ridx[i];
99426e093fcSHong Zhang     }
99515091d37SBarry Smith     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4]; sum6 = y[5];
99615091d37SBarry Smith     for (j=0; j<n; j++) {
9973b95cb0eSSatish Balay       xb = x + 6*(*idx++);
99815091d37SBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5];
99915091d37SBarry Smith       sum1 += v[0]*x1 + v[6]*x2  + v[12]*x3  + v[18]*x4 + v[24]*x5 + v[30]*x6;
100015091d37SBarry Smith       sum2 += v[1]*x1 + v[7]*x2  + v[13]*x3  + v[19]*x4 + v[25]*x5 + v[31]*x6;
100115091d37SBarry Smith       sum3 += v[2]*x1 + v[8]*x2  + v[14]*x3  + v[20]*x4 + v[26]*x5 + v[32]*x6;
100215091d37SBarry Smith       sum4 += v[3]*x1 + v[9]*x2  + v[15]*x3  + v[21]*x4 + v[27]*x5 + v[33]*x6;
100315091d37SBarry Smith       sum5 += v[4]*x1 + v[10]*x2 + v[16]*x3  + v[22]*x4 + v[28]*x5 + v[34]*x6;
100415091d37SBarry Smith       sum6 += v[5]*x1 + v[11]*x2 + v[17]*x3  + v[23]*x4 + v[29]*x5 + v[35]*x6;
100515091d37SBarry Smith       v += 36;
100615091d37SBarry Smith     }
100715091d37SBarry Smith     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6;
100826e093fcSHong Zhang     if (!usecprow){
100915091d37SBarry Smith       z += 6; y += 6;
101015091d37SBarry Smith     }
101126e093fcSHong Zhang   }
10121ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
101326e093fcSHong Zhang   ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr);
101415091d37SBarry Smith   if (zz != yy) {
101526e093fcSHong Zhang     ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
101615091d37SBarry Smith   }
1017efee365bSSatish Balay   ierr = PetscLogFlops(72*a->nz);CHKERRQ(ierr);
101815091d37SBarry Smith   PetscFunctionReturn(0);
101915091d37SBarry Smith }
10202d61bbb3SSatish Balay 
10214a2ae208SSatish Balay #undef __FUNCT__
10224a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_7"
1023dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_7(Mat A,Vec xx,Vec yy,Vec zz)
10242d61bbb3SSatish Balay {
10252d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
1026dcf5cc72SBarry Smith   PetscScalar    *x,*y = 0,*z = 0,*xb,sum1,sum2,sum3,sum4,sum5,sum6,sum7;
102726e093fcSHong Zhang   PetscScalar    x1,x2,x3,x4,x5,x6,x7,*yarray,*zarray;
10283f1db9ecSBarry Smith   MatScalar      *v;
1029dfbe8321SBarry Smith   PetscErrorCode ierr;
10307b2bb3b9SHong Zhang   PetscInt       mbs=a->mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL;
103126e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
10322d61bbb3SSatish Balay 
10332d61bbb3SSatish Balay   PetscFunctionBegin;
10341ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
103526e093fcSHong Zhang   ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr);
10362e8a6d31SBarry Smith   if (zz != yy) {
103726e093fcSHong Zhang     ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
10382e8a6d31SBarry Smith   } else {
103926e093fcSHong Zhang     zarray = yarray;
10402e8a6d31SBarry Smith   }
10412d61bbb3SSatish Balay 
10422d61bbb3SSatish Balay   idx = a->j;
10432d61bbb3SSatish Balay   v   = a->a;
104426e093fcSHong Zhang   if (usecprow){
10454eb6d288SHong Zhang     if (zz != yy){
10464eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,7*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
10474eb6d288SHong Zhang     }
104826e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
104926e093fcSHong Zhang     ii   = a->compressedrow.i;
10507b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
105126e093fcSHong Zhang   } else {
10522d61bbb3SSatish Balay     ii  = a->i;
105326e093fcSHong Zhang     y   = yarray;
105426e093fcSHong Zhang     z   = zarray;
105526e093fcSHong Zhang   }
10562d61bbb3SSatish Balay 
10572d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
10582d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
105926e093fcSHong Zhang     if (usecprow){
10607b2bb3b9SHong Zhang       z = zarray + 7*ridx[i];
10617b2bb3b9SHong Zhang       y = yarray + 7*ridx[i];
106226e093fcSHong Zhang     }
10632d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4]; sum6 = y[5]; sum7 = y[6];
10642d61bbb3SSatish Balay     for (j=0; j<n; j++) {
10652d61bbb3SSatish Balay       xb = x + 7*(*idx++);
10662d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6];
10672d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[7]*x2  + v[14]*x3  + v[21]*x4 + v[28]*x5 + v[35]*x6 + v[42]*x7;
10682d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[8]*x2  + v[15]*x3  + v[22]*x4 + v[29]*x5 + v[36]*x6 + v[43]*x7;
10692d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[9]*x2  + v[16]*x3  + v[23]*x4 + v[30]*x5 + v[37]*x6 + v[44]*x7;
10702d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[10]*x2 + v[17]*x3  + v[24]*x4 + v[31]*x5 + v[38]*x6 + v[45]*x7;
10712d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[11]*x2 + v[18]*x3  + v[25]*x4 + v[32]*x5 + v[39]*x6 + v[46]*x7;
10722d61bbb3SSatish Balay       sum6 += v[5]*x1 + v[12]*x2 + v[19]*x3  + v[26]*x4 + v[33]*x5 + v[40]*x6 + v[47]*x7;
10732d61bbb3SSatish Balay       sum7 += v[6]*x1 + v[13]*x2 + v[20]*x3  + v[27]*x4 + v[34]*x5 + v[41]*x6 + v[48]*x7;
10742d61bbb3SSatish Balay       v += 49;
10752d61bbb3SSatish Balay     }
10762d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
107726e093fcSHong Zhang     if (!usecprow){
10782d61bbb3SSatish Balay       z += 7; y += 7;
10792d61bbb3SSatish Balay     }
108026e093fcSHong Zhang   }
10811ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
108226e093fcSHong Zhang   ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr);
10832e8a6d31SBarry Smith   if (zz != yy) {
108426e093fcSHong Zhang     ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
10852e8a6d31SBarry Smith   }
1086efee365bSSatish Balay   ierr = PetscLogFlops(98*a->nz);CHKERRQ(ierr);
10872d61bbb3SSatish Balay   PetscFunctionReturn(0);
10882d61bbb3SSatish Balay }
1089218c64b6SSatish Balay 
10904a2ae208SSatish Balay #undef __FUNCT__
10914a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_N"
1092dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_N(Mat A,Vec xx,Vec yy,Vec zz)
10932d61bbb3SSatish Balay {
10942d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
1095bba805e6SBarry Smith   PetscScalar    *x,*z = 0,*xb,*work,*workt,*zarray;
10963f1db9ecSBarry Smith   MatScalar      *v;
10976849ba73SBarry Smith   PetscErrorCode ierr;
1098899cda47SBarry Smith   PetscInt       mbs,i,*idx,*ii,bs=A->rmap.bs,j,n,bs2=a->bs2;
10997b2bb3b9SHong Zhang   PetscInt       ncols,k,*ridx=PETSC_NULL;
110026e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
1101218c64b6SSatish Balay 
11022d61bbb3SSatish Balay   PetscFunctionBegin;
11036a65c766SHong Zhang   ierr = VecCopy_Seq(yy,zz);CHKERRQ(ierr);
11041ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
110526e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
11062d61bbb3SSatish Balay 
11072d61bbb3SSatish Balay   idx = a->j;
11082d61bbb3SSatish Balay   v   = a->a;
110926e093fcSHong Zhang   if (usecprow){
111026e093fcSHong Zhang     mbs    = a->compressedrow.nrows;
111126e093fcSHong Zhang     ii     = a->compressedrow.i;
11127b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
111326e093fcSHong Zhang   } else {
111426e093fcSHong Zhang     mbs = a->mbs;
11152d61bbb3SSatish Balay     ii  = a->i;
111626e093fcSHong Zhang     z   = zarray;
111726e093fcSHong Zhang   }
11182d61bbb3SSatish Balay 
11192d61bbb3SSatish Balay   if (!a->mult_work) {
1120899cda47SBarry Smith     k    = PetscMax(A->rmap.n,A->cmap.n);
112187828ca2SBarry Smith     ierr = PetscMalloc((k+1)*sizeof(PetscScalar),&a->mult_work);CHKERRQ(ierr);
11222d61bbb3SSatish Balay   }
11232d61bbb3SSatish Balay   work = a->mult_work;
11242d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
11252d61bbb3SSatish Balay     n     = ii[1] - ii[0]; ii++;
11262d61bbb3SSatish Balay     ncols = n*bs;
11272d61bbb3SSatish Balay     workt = work;
11282d61bbb3SSatish Balay     for (j=0; j<n; j++) {
11292d61bbb3SSatish Balay       xb = x + bs*(*idx++);
11302d61bbb3SSatish Balay       for (k=0; k<bs; k++) workt[k] = xb[k];
11312d61bbb3SSatish Balay       workt += bs;
11322d61bbb3SSatish Balay     }
11337b2bb3b9SHong Zhang     if (usecprow) z = zarray + bs*ridx[i];
11343f1db9ecSBarry Smith     Kernel_w_gets_w_plus_Ar_times_v(bs,ncols,work,v,z);
113571044d3cSBarry Smith     /* BLASgemv_("N",&bs,&ncols,&_DOne,v,&bs,work,&_One,&_DOne,z,&_One); */
11362d61bbb3SSatish Balay     v += n*bs2;
113726e093fcSHong Zhang     if (!usecprow){
11382d61bbb3SSatish Balay       z += bs;
11392d61bbb3SSatish Balay     }
114026e093fcSHong Zhang   }
11411ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
114226e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
1143efee365bSSatish Balay   ierr = PetscLogFlops(2*a->nz*bs2);CHKERRQ(ierr);
11442d61bbb3SSatish Balay   PetscFunctionReturn(0);
11452d61bbb3SSatish Balay }
11462d61bbb3SSatish Balay 
11474a2ae208SSatish Balay #undef __FUNCT__
11484a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_SeqBAIJ"
1149dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_SeqBAIJ(Mat A,Vec xx,Vec zz)
11502d61bbb3SSatish Balay {
11513447b6efSHong Zhang   PetscScalar    zero = 0.0;
11526849ba73SBarry Smith   PetscErrorCode ierr;
1153e8e7597cSSatish Balay #if defined (PETSC_USE_LOG)
1154e8e7597cSSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
1155e8e7597cSSatish Balay #endif
11562d61bbb3SSatish Balay 
11572d61bbb3SSatish Balay   PetscFunctionBegin;
11582dcb1b2aSMatthew Knepley   ierr = VecSet(zz,zero);CHKERRQ(ierr);
11593447b6efSHong Zhang   ierr = MatMultTransposeAdd_SeqBAIJ(A,xx,zz,zz);CHKERRQ(ierr);
11602d61bbb3SSatish Balay 
1161899cda47SBarry Smith   ierr = PetscLogFlops(2*a->nz*a->bs2 - A->cmap.n);CHKERRQ(ierr);
11622d61bbb3SSatish Balay   PetscFunctionReturn(0);
11632d61bbb3SSatish Balay }
11642d61bbb3SSatish Balay 
11654a2ae208SSatish Balay #undef __FUNCT__
11664a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_SeqBAIJ"
1167dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_SeqBAIJ(Mat A,Vec xx,Vec yy,Vec zz)
11682d61bbb3SSatish Balay 
11692d61bbb3SSatish Balay {
11702d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
1171dcf5cc72SBarry Smith   PetscScalar    *zb,*x,*z,*xb = 0,x1,x2,x3,x4,x5;
11723f1db9ecSBarry Smith   MatScalar      *v;
11736849ba73SBarry Smith   PetscErrorCode ierr;
1174899cda47SBarry Smith   PetscInt       mbs,i,*idx,*ii,rval,bs=A->rmap.bs,j,n,bs2=a->bs2,*ib,*ridx=PETSC_NULL;
11753447b6efSHong Zhang   Mat_CompressedRow cprow = a->compressedrow;
11763447b6efSHong Zhang   PetscTruth        usecprow=cprow.use;
11772d61bbb3SSatish Balay 
11782d61bbb3SSatish Balay   PetscFunctionBegin;
11796a65c766SHong Zhang   if (yy != zz) { ierr = VecCopy_Seq(yy,zz);CHKERRQ(ierr); }
11803447b6efSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
11813447b6efSHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
11822d61bbb3SSatish Balay 
11832d61bbb3SSatish Balay   idx = a->j;
11842d61bbb3SSatish Balay   v   = a->a;
11853447b6efSHong Zhang   if (usecprow){
11863447b6efSHong Zhang     mbs  = cprow.nrows;
11873447b6efSHong Zhang     ii   = cprow.i;
11887b2bb3b9SHong Zhang     ridx = cprow.rindex;
11893447b6efSHong Zhang   } else {
11903447b6efSHong Zhang     mbs=a->mbs;
11912d61bbb3SSatish Balay     ii = a->i;
1192f1af5d2fSBarry Smith     xb = x;
11933447b6efSHong Zhang   }
11942d61bbb3SSatish Balay 
11952d61bbb3SSatish Balay   switch (bs) {
11962d61bbb3SSatish Balay   case 1:
11972d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
11987b2bb3b9SHong Zhang       if (usecprow) xb = x + ridx[i];
1199f1af5d2fSBarry Smith       x1 = xb[0];
12003447b6efSHong Zhang       ib = idx + ii[0];
12013447b6efSHong Zhang       n  = ii[1] - ii[0]; ii++;
12022d61bbb3SSatish Balay       for (j=0; j<n; j++) {
12032d61bbb3SSatish Balay         rval    = ib[j];
1204f1af5d2fSBarry Smith         z[rval] += *v * x1;
1205f1af5d2fSBarry Smith         v++;
12062d61bbb3SSatish Balay       }
12073447b6efSHong Zhang       if (!usecprow) xb++;
12082d61bbb3SSatish Balay     }
12092d61bbb3SSatish Balay     break;
12102d61bbb3SSatish Balay   case 2:
12112d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
12127b2bb3b9SHong Zhang       if (usecprow) xb = x + 2*ridx[i];
1213f1af5d2fSBarry Smith       x1 = xb[0]; x2 = xb[1];
12143447b6efSHong Zhang       ib = idx + ii[0];
12153447b6efSHong Zhang       n  = ii[1] - ii[0]; ii++;
12162d61bbb3SSatish Balay       for (j=0; j<n; j++) {
12172d61bbb3SSatish Balay         rval      = ib[j]*2;
12182d61bbb3SSatish Balay         z[rval++] += v[0]*x1 + v[1]*x2;
12192d61bbb3SSatish Balay         z[rval++] += v[2]*x1 + v[3]*x2;
12202d61bbb3SSatish Balay         v  += 4;
12212d61bbb3SSatish Balay       }
12223447b6efSHong Zhang       if (!usecprow) xb += 2;
12232d61bbb3SSatish Balay     }
12242d61bbb3SSatish Balay     break;
12252d61bbb3SSatish Balay   case 3:
12262d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
12277b2bb3b9SHong Zhang       if (usecprow) xb = x + 3*ridx[i];
1228f1af5d2fSBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
12293447b6efSHong Zhang       ib = idx + ii[0];
12303447b6efSHong Zhang       n  = ii[1] - ii[0]; ii++;
12312d61bbb3SSatish Balay       for (j=0; j<n; j++) {
12322d61bbb3SSatish Balay         rval      = ib[j]*3;
12332d61bbb3SSatish Balay         z[rval++] += v[0]*x1 + v[1]*x2 + v[2]*x3;
12342d61bbb3SSatish Balay         z[rval++] += v[3]*x1 + v[4]*x2 + v[5]*x3;
12352d61bbb3SSatish Balay         z[rval++] += v[6]*x1 + v[7]*x2 + v[8]*x3;
12362d61bbb3SSatish Balay         v  += 9;
12372d61bbb3SSatish Balay       }
12383447b6efSHong Zhang       if (!usecprow) xb += 3;
12392d61bbb3SSatish Balay     }
12402d61bbb3SSatish Balay     break;
12412d61bbb3SSatish Balay   case 4:
12422d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
12437b2bb3b9SHong Zhang       if (usecprow) xb = x + 4*ridx[i];
1244f1af5d2fSBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
12453447b6efSHong Zhang       ib = idx + ii[0];
12463447b6efSHong Zhang       n  = ii[1] - ii[0]; ii++;
12472d61bbb3SSatish Balay       for (j=0; j<n; j++) {
12482d61bbb3SSatish Balay         rval      = ib[j]*4;
12492d61bbb3SSatish Balay         z[rval++] +=  v[0]*x1 +  v[1]*x2 +  v[2]*x3 +  v[3]*x4;
12502d61bbb3SSatish Balay         z[rval++] +=  v[4]*x1 +  v[5]*x2 +  v[6]*x3 +  v[7]*x4;
12512d61bbb3SSatish Balay         z[rval++] +=  v[8]*x1 +  v[9]*x2 + v[10]*x3 + v[11]*x4;
12522d61bbb3SSatish Balay         z[rval++] += v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4;
12532d61bbb3SSatish Balay         v  += 16;
12542d61bbb3SSatish Balay       }
12553447b6efSHong Zhang       if (!usecprow) xb += 4;
12562d61bbb3SSatish Balay     }
12572d61bbb3SSatish Balay     break;
12582d61bbb3SSatish Balay   case 5:
12592d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
12607b2bb3b9SHong Zhang       if (usecprow) xb = x + 5*ridx[i];
1261f1af5d2fSBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
12622d61bbb3SSatish Balay       x4 = xb[3]; x5 = xb[4];
12633447b6efSHong Zhang       ib = idx + ii[0];
12643447b6efSHong Zhang       n  = ii[1] - ii[0]; ii++;
12652d61bbb3SSatish Balay       for (j=0; j<n; j++) {
12662d61bbb3SSatish Balay         rval      = ib[j]*5;
12672d61bbb3SSatish Balay         z[rval++] +=  v[0]*x1 +  v[1]*x2 +  v[2]*x3 +  v[3]*x4 +  v[4]*x5;
12682d61bbb3SSatish Balay         z[rval++] +=  v[5]*x1 +  v[6]*x2 +  v[7]*x3 +  v[8]*x4 +  v[9]*x5;
12692d61bbb3SSatish Balay         z[rval++] += v[10]*x1 + v[11]*x2 + v[12]*x3 + v[13]*x4 + v[14]*x5;
12702d61bbb3SSatish Balay         z[rval++] += v[15]*x1 + v[16]*x2 + v[17]*x3 + v[18]*x4 + v[19]*x5;
12712d61bbb3SSatish Balay         z[rval++] += v[20]*x1 + v[21]*x2 + v[22]*x3 + v[23]*x4 + v[24]*x5;
12722d61bbb3SSatish Balay         v  += 25;
12732d61bbb3SSatish Balay       }
12743447b6efSHong Zhang       if (!usecprow) xb += 5;
12752d61bbb3SSatish Balay     }
12762d61bbb3SSatish Balay     break;
1277f1af5d2fSBarry Smith   default: {      /* block sizes larger then 5 by 5 are handled by BLAS */
1278690b6cddSBarry Smith       PetscInt     ncols,k;
12793447b6efSHong Zhang       PetscScalar  *work,*workt,*xtmp;
12803f1db9ecSBarry Smith 
12812d61bbb3SSatish Balay       if (!a->mult_work) {
1282899cda47SBarry Smith         k = PetscMax(A->rmap.n,A->cmap.n);
128387828ca2SBarry Smith         ierr = PetscMalloc((k+1)*sizeof(PetscScalar),&a->mult_work);CHKERRQ(ierr);
12842d61bbb3SSatish Balay       }
12852d61bbb3SSatish Balay       work = a->mult_work;
12863447b6efSHong Zhang       xtmp = x;
12872d61bbb3SSatish Balay       for (i=0; i<mbs; i++) {
12882d61bbb3SSatish Balay         n     = ii[1] - ii[0]; ii++;
12892d61bbb3SSatish Balay         ncols = n*bs;
129087828ca2SBarry Smith         ierr  = PetscMemzero(work,ncols*sizeof(PetscScalar));CHKERRQ(ierr);
12913447b6efSHong Zhang         if (usecprow) {
12927b2bb3b9SHong Zhang           xtmp = x + bs*ridx[i];
12933447b6efSHong Zhang         }
12943447b6efSHong Zhang         Kernel_w_gets_w_plus_trans_Ar_times_v(bs,ncols,xtmp,v,work);
129571044d3cSBarry Smith         /* BLASgemv_("T",&bs,&ncols,&_DOne,v,&bs,xtmp,&_One,&_DOne,work,&_One); */
12962d61bbb3SSatish Balay         v += n*bs2;
12973447b6efSHong Zhang         if (!usecprow) xtmp += bs;
12982d61bbb3SSatish Balay         workt = work;
12992d61bbb3SSatish Balay         for (j=0; j<n; j++) {
13002d61bbb3SSatish Balay           zb = z + bs*(*idx++);
13012d61bbb3SSatish Balay           for (k=0; k<bs; k++) zb[k] += workt[k] ;
13022d61bbb3SSatish Balay           workt += bs;
13032d61bbb3SSatish Balay         }
13042d61bbb3SSatish Balay       }
13052d61bbb3SSatish Balay     }
13062d61bbb3SSatish Balay   }
13073447b6efSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
13083447b6efSHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
1309efee365bSSatish Balay   ierr = PetscLogFlops(2*a->nz*a->bs2);CHKERRQ(ierr);
13102d61bbb3SSatish Balay   PetscFunctionReturn(0);
13112d61bbb3SSatish Balay }
13122d61bbb3SSatish Balay 
13134a2ae208SSatish Balay #undef __FUNCT__
13144a2ae208SSatish Balay #define __FUNCT__ "MatScale_SeqBAIJ"
1315f4df32b1SMatthew Knepley PetscErrorCode MatScale_SeqBAIJ(Mat inA,PetscScalar alpha)
13162d61bbb3SSatish Balay {
13172d61bbb3SSatish Balay   Mat_SeqBAIJ  *a = (Mat_SeqBAIJ*)inA->data;
1318690b6cddSBarry Smith   PetscInt     totalnz = a->bs2*a->nz;
1319f4df32b1SMatthew Knepley   PetscScalar  oalpha = alpha;
1320efee365bSSatish Balay   PetscErrorCode ierr;
13213eda8832SBarry Smith #if defined(PETSC_USE_MAT_SINGLE)
1322690b6cddSBarry Smith   PetscInt     i;
132387052302SDinesh Kaushik #else
13244ce68768SBarry Smith   PetscBLASInt tnz = (PetscBLASInt) totalnz,one = 1;
13253eda8832SBarry Smith #endif
13262d61bbb3SSatish Balay 
13272d61bbb3SSatish Balay   PetscFunctionBegin;
13283eda8832SBarry Smith #if defined(PETSC_USE_MAT_SINGLE)
1329f4df32b1SMatthew Knepley   for (i=0; i<totalnz; i++) a->a[i] *= oalpha;
13303eda8832SBarry Smith #else
1331f4df32b1SMatthew Knepley   BLASscal_(&tnz,&oalpha,a->a,&one);
13323eda8832SBarry Smith #endif
1333efee365bSSatish Balay   ierr = PetscLogFlops(totalnz);CHKERRQ(ierr);
13342d61bbb3SSatish Balay   PetscFunctionReturn(0);
13352d61bbb3SSatish Balay }
13362d61bbb3SSatish Balay 
13374a2ae208SSatish Balay #undef __FUNCT__
13384a2ae208SSatish Balay #define __FUNCT__ "MatNorm_SeqBAIJ"
1339dfbe8321SBarry Smith PetscErrorCode MatNorm_SeqBAIJ(Mat A,NormType type,PetscReal *norm)
13402d61bbb3SSatish Balay {
13418a62d963SHong Zhang   PetscErrorCode ierr;
13422d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
13433f1db9ecSBarry Smith   MatScalar      *v = a->a;
1344329f5518SBarry Smith   PetscReal      sum = 0.0;
1345899cda47SBarry Smith   PetscInt       i,j,k,bs=A->rmap.bs,nz=a->nz,bs2=a->bs2,k1;
13462d61bbb3SSatish Balay 
13472d61bbb3SSatish Balay   PetscFunctionBegin;
13482d61bbb3SSatish Balay   if (type == NORM_FROBENIUS) {
13492d61bbb3SSatish Balay     for (i=0; i< bs2*nz; i++) {
1350aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
1351329f5518SBarry Smith       sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
13522d61bbb3SSatish Balay #else
13532d61bbb3SSatish Balay       sum += (*v)*(*v); v++;
13542d61bbb3SSatish Balay #endif
13552d61bbb3SSatish Balay     }
13562d61bbb3SSatish Balay     *norm = sqrt(sum);
13578a62d963SHong Zhang   } else if (type == NORM_1) { /* maximum column sum */
13588a62d963SHong Zhang     PetscReal *tmp;
13598a62d963SHong Zhang     PetscInt  *bcol = a->j;
1360899cda47SBarry Smith     ierr = PetscMalloc((A->cmap.n+1)*sizeof(PetscReal),&tmp);CHKERRQ(ierr);
1361899cda47SBarry Smith     ierr = PetscMemzero(tmp,A->cmap.n*sizeof(PetscReal));CHKERRQ(ierr);
13628a62d963SHong Zhang     for (i=0; i<nz; i++){
13638a62d963SHong Zhang       for (j=0; j<bs; j++){
13648a62d963SHong Zhang         k1 = bs*(*bcol) + j; /* column index */
13658a62d963SHong Zhang         for (k=0; k<bs; k++){
13668a62d963SHong Zhang           tmp[k1] += PetscAbsScalar(*v); v++;
13678a62d963SHong Zhang         }
13688a62d963SHong Zhang       }
13698a62d963SHong Zhang       bcol++;
13708a62d963SHong Zhang     }
13718a62d963SHong Zhang     *norm = 0.0;
1372899cda47SBarry Smith     for (j=0; j<A->cmap.n; j++) {
13738a62d963SHong Zhang       if (tmp[j] > *norm) *norm = tmp[j];
13748a62d963SHong Zhang     }
13758a62d963SHong Zhang     ierr = PetscFree(tmp);CHKERRQ(ierr);
1376596552b5SBarry Smith   } else if (type == NORM_INFINITY) { /* maximum row sum */
1377596552b5SBarry Smith     *norm = 0.0;
1378596552b5SBarry Smith     for (k=0; k<bs; k++) {
137974f84c7bSSatish Balay       for (j=0; j<a->mbs; j++) {
1380596552b5SBarry Smith         v = a->a + bs2*a->i[j] + k;
1381596552b5SBarry Smith         sum = 0.0;
1382596552b5SBarry Smith         for (i=0; i<a->i[j+1]-a->i[j]; i++) {
13830e90e235SBarry Smith           for (k1=0; k1<bs; k1++){
1384596552b5SBarry Smith             sum += PetscAbsScalar(*v);
1385596552b5SBarry Smith             v   += bs;
13862d61bbb3SSatish Balay           }
13870e90e235SBarry Smith         }
1388596552b5SBarry Smith         if (sum > *norm) *norm = sum;
1389596552b5SBarry Smith       }
1390596552b5SBarry Smith     }
1391596552b5SBarry Smith   } else {
139229bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"No support for this norm yet");
13932d61bbb3SSatish Balay   }
13942d61bbb3SSatish Balay   PetscFunctionReturn(0);
13952d61bbb3SSatish Balay }
13962d61bbb3SSatish Balay 
13972d61bbb3SSatish Balay 
13984a2ae208SSatish Balay #undef __FUNCT__
13994a2ae208SSatish Balay #define __FUNCT__ "MatEqual_SeqBAIJ"
1400dfbe8321SBarry Smith PetscErrorCode MatEqual_SeqBAIJ(Mat A,Mat B,PetscTruth* flg)
14012d61bbb3SSatish Balay {
14022d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ *)A->data,*b = (Mat_SeqBAIJ *)B->data;
1403dfbe8321SBarry Smith   PetscErrorCode ierr;
14042d61bbb3SSatish Balay 
14052d61bbb3SSatish Balay   PetscFunctionBegin;
14062d61bbb3SSatish Balay   /* If the  matrix/block dimensions are not equal, or no of nonzeros or shift */
1407899cda47SBarry Smith   if ((A->rmap.N != B->rmap.N) || (A->cmap.n != B->cmap.n) || (A->rmap.bs != B->rmap.bs)|| (a->nz != b->nz)) {
1408273d9f13SBarry Smith     *flg = PETSC_FALSE;
1409273d9f13SBarry Smith     PetscFunctionReturn(0);
14102d61bbb3SSatish Balay   }
14112d61bbb3SSatish Balay 
14122d61bbb3SSatish Balay   /* if the a->i are the same */
1413690b6cddSBarry Smith   ierr = PetscMemcmp(a->i,b->i,(a->mbs+1)*sizeof(PetscInt),flg);CHKERRQ(ierr);
1414abc0a331SBarry Smith   if (!*flg) {
14150f5bd95cSBarry Smith     PetscFunctionReturn(0);
14162d61bbb3SSatish Balay   }
14172d61bbb3SSatish Balay 
14182d61bbb3SSatish Balay   /* if a->j are the same */
1419690b6cddSBarry Smith   ierr = PetscMemcmp(a->j,b->j,(a->nz)*sizeof(PetscInt),flg);CHKERRQ(ierr);
1420abc0a331SBarry Smith   if (!*flg) {
14210f5bd95cSBarry Smith     PetscFunctionReturn(0);
14222d61bbb3SSatish Balay   }
14232d61bbb3SSatish Balay   /* if a->a are the same */
1424899cda47SBarry Smith   ierr = PetscMemcmp(a->a,b->a,(a->nz)*(A->rmap.bs)*(B->rmap.bs)*sizeof(PetscScalar),flg);CHKERRQ(ierr);
14252d61bbb3SSatish Balay   PetscFunctionReturn(0);
14262d61bbb3SSatish Balay 
14272d61bbb3SSatish Balay }
14282d61bbb3SSatish Balay 
14294a2ae208SSatish Balay #undef __FUNCT__
14304a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_SeqBAIJ"
1431dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_SeqBAIJ(Mat A,Vec v)
14322d61bbb3SSatish Balay {
14332d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
1434dfbe8321SBarry Smith   PetscErrorCode ierr;
1435690b6cddSBarry Smith   PetscInt       i,j,k,n,row,bs,*ai,*aj,ambs,bs2;
143687828ca2SBarry Smith   PetscScalar    *x,zero = 0.0;
14373f1db9ecSBarry Smith   MatScalar      *aa,*aa_j;
14382d61bbb3SSatish Balay 
14392d61bbb3SSatish Balay   PetscFunctionBegin;
144029bbc08cSBarry Smith   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1441899cda47SBarry Smith   bs   = A->rmap.bs;
14422d61bbb3SSatish Balay   aa   = a->a;
14432d61bbb3SSatish Balay   ai   = a->i;
14442d61bbb3SSatish Balay   aj   = a->j;
14452d61bbb3SSatish Balay   ambs = a->mbs;
14462d61bbb3SSatish Balay   bs2  = a->bs2;
14472d61bbb3SSatish Balay 
14482dcb1b2aSMatthew Knepley   ierr = VecSet(v,zero);CHKERRQ(ierr);
14491ebc52fbSHong Zhang   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1450e1311b90SBarry Smith   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
1451899cda47SBarry Smith   if (n != A->rmap.N) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
14522d61bbb3SSatish Balay   for (i=0; i<ambs; i++) {
14532d61bbb3SSatish Balay     for (j=ai[i]; j<ai[i+1]; j++) {
14542d61bbb3SSatish Balay       if (aj[j] == i) {
14552d61bbb3SSatish Balay         row  = i*bs;
14562d61bbb3SSatish Balay         aa_j = aa+j*bs2;
14572d61bbb3SSatish Balay         for (k=0; k<bs2; k+=(bs+1),row++) x[row] = aa_j[k];
14582d61bbb3SSatish Balay         break;
14592d61bbb3SSatish Balay       }
14602d61bbb3SSatish Balay     }
14612d61bbb3SSatish Balay   }
14621ebc52fbSHong Zhang   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
14632d61bbb3SSatish Balay   PetscFunctionReturn(0);
14642d61bbb3SSatish Balay }
14652d61bbb3SSatish Balay 
14664a2ae208SSatish Balay #undef __FUNCT__
14674a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_SeqBAIJ"
1468dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_SeqBAIJ(Mat A,Vec ll,Vec rr)
14692d61bbb3SSatish Balay {
14702d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
147187828ca2SBarry Smith   PetscScalar    *l,*r,x,*li,*ri;
14723f1db9ecSBarry Smith   MatScalar      *aa,*v;
1473dfbe8321SBarry Smith   PetscErrorCode ierr;
1474690b6cddSBarry Smith   PetscInt       i,j,k,lm,rn,M,m,n,*ai,*aj,mbs,tmp,bs,bs2;
14752d61bbb3SSatish Balay 
14762d61bbb3SSatish Balay   PetscFunctionBegin;
14772d61bbb3SSatish Balay   ai  = a->i;
14782d61bbb3SSatish Balay   aj  = a->j;
14792d61bbb3SSatish Balay   aa  = a->a;
1480899cda47SBarry Smith   m   = A->rmap.n;
1481899cda47SBarry Smith   n   = A->cmap.n;
1482899cda47SBarry Smith   bs  = A->rmap.bs;
14832d61bbb3SSatish Balay   mbs = a->mbs;
14842d61bbb3SSatish Balay   bs2 = a->bs2;
14852d61bbb3SSatish Balay   if (ll) {
14861ebc52fbSHong Zhang     ierr = VecGetArray(ll,&l);CHKERRQ(ierr);
1487e1311b90SBarry Smith     ierr = VecGetLocalSize(ll,&lm);CHKERRQ(ierr);
148829bbc08cSBarry Smith     if (lm != m) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vector wrong length");
14892d61bbb3SSatish Balay     for (i=0; i<mbs; i++) { /* for each block row */
14902d61bbb3SSatish Balay       M  = ai[i+1] - ai[i];
14912d61bbb3SSatish Balay       li = l + i*bs;
14922d61bbb3SSatish Balay       v  = aa + bs2*ai[i];
14932d61bbb3SSatish Balay       for (j=0; j<M; j++) { /* for each block */
14942d61bbb3SSatish Balay         for (k=0; k<bs2; k++) {
14952d61bbb3SSatish Balay           (*v++) *= li[k%bs];
14962d61bbb3SSatish Balay         }
14972d61bbb3SSatish Balay       }
14982d61bbb3SSatish Balay     }
14991ebc52fbSHong Zhang     ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr);
1500efee365bSSatish Balay     ierr = PetscLogFlops(a->nz);CHKERRQ(ierr);
15012d61bbb3SSatish Balay   }
15022d61bbb3SSatish Balay 
15032d61bbb3SSatish Balay   if (rr) {
15041ebc52fbSHong Zhang     ierr = VecGetArray(rr,&r);CHKERRQ(ierr);
1505e1311b90SBarry Smith     ierr = VecGetLocalSize(rr,&rn);CHKERRQ(ierr);
150629bbc08cSBarry Smith     if (rn != n) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vector wrong length");
15072d61bbb3SSatish Balay     for (i=0; i<mbs; i++) { /* for each block row */
15082d61bbb3SSatish Balay       M  = ai[i+1] - ai[i];
15092d61bbb3SSatish Balay       v  = aa + bs2*ai[i];
15102d61bbb3SSatish Balay       for (j=0; j<M; j++) { /* for each block */
15112d61bbb3SSatish Balay         ri = r + bs*aj[ai[i]+j];
15122d61bbb3SSatish Balay         for (k=0; k<bs; k++) {
15132d61bbb3SSatish Balay           x = ri[k];
15142d61bbb3SSatish Balay           for (tmp=0; tmp<bs; tmp++) (*v++) *= x;
15152d61bbb3SSatish Balay         }
15162d61bbb3SSatish Balay       }
15172d61bbb3SSatish Balay     }
15181ebc52fbSHong Zhang     ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr);
1519efee365bSSatish Balay     ierr = PetscLogFlops(a->nz);CHKERRQ(ierr);
15202d61bbb3SSatish Balay   }
15212d61bbb3SSatish Balay   PetscFunctionReturn(0);
15222d61bbb3SSatish Balay }
15232d61bbb3SSatish Balay 
15242d61bbb3SSatish Balay 
15254a2ae208SSatish Balay #undef __FUNCT__
15264a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_SeqBAIJ"
1527dfbe8321SBarry Smith PetscErrorCode MatGetInfo_SeqBAIJ(Mat A,MatInfoType flag,MatInfo *info)
15282d61bbb3SSatish Balay {
15292d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
15302d61bbb3SSatish Balay 
15312d61bbb3SSatish Balay   PetscFunctionBegin;
1532899cda47SBarry Smith   info->rows_global    = (double)A->rmap.N;
1533899cda47SBarry Smith   info->columns_global = (double)A->cmap.n;
1534899cda47SBarry Smith   info->rows_local     = (double)A->rmap.n;
1535899cda47SBarry Smith   info->columns_local  = (double)A->cmap.n;
15362d61bbb3SSatish Balay   info->block_size     = a->bs2;
15372d61bbb3SSatish Balay   info->nz_allocated   = a->maxnz;
15382d61bbb3SSatish Balay   info->nz_used        = a->bs2*a->nz;
15392d61bbb3SSatish Balay   info->nz_unneeded    = (double)(info->nz_allocated - info->nz_used);
15402d61bbb3SSatish Balay   info->assemblies   = A->num_ass;
15412d61bbb3SSatish Balay   info->mallocs      = a->reallocs;
1542*7adad957SLisandro Dalcin   info->memory       = ((PetscObject)A)->mem;
15432d61bbb3SSatish Balay   if (A->factor) {
15442d61bbb3SSatish Balay     info->fill_ratio_given  = A->info.fill_ratio_given;
15452d61bbb3SSatish Balay     info->fill_ratio_needed = A->info.fill_ratio_needed;
15462d61bbb3SSatish Balay     info->factor_mallocs    = A->info.factor_mallocs;
15472d61bbb3SSatish Balay   } else {
15482d61bbb3SSatish Balay     info->fill_ratio_given  = 0;
15492d61bbb3SSatish Balay     info->fill_ratio_needed = 0;
15502d61bbb3SSatish Balay     info->factor_mallocs    = 0;
15512d61bbb3SSatish Balay   }
15522d61bbb3SSatish Balay   PetscFunctionReturn(0);
15532d61bbb3SSatish Balay }
15542d61bbb3SSatish Balay 
15552d61bbb3SSatish Balay 
15564a2ae208SSatish Balay #undef __FUNCT__
15574a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_SeqBAIJ"
1558dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_SeqBAIJ(Mat A)
15592d61bbb3SSatish Balay {
15602d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
1561dfbe8321SBarry Smith   PetscErrorCode ierr;
15622d61bbb3SSatish Balay 
15632d61bbb3SSatish Balay   PetscFunctionBegin;
1564549d3d68SSatish Balay   ierr = PetscMemzero(a->a,a->bs2*a->i[a->mbs]*sizeof(MatScalar));CHKERRQ(ierr);
15652d61bbb3SSatish Balay   PetscFunctionReturn(0);
15662d61bbb3SSatish Balay }
1567