xref: /petsc/src/mat/impls/baij/seq/baij2.c (revision c60f02097ec1bad1d486ea4e6635560ad6843df9)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
2cac129eeSSatish Balay 
37c4f633dSBarry Smith #include "../src/mat/impls/baij/seq/baij.h"
4*c60f0209SBarry Smith #include "../src/mat/blockinvert.h"
50a835dfdSSatish Balay #include "petscbt.h"
6cac129eeSSatish Balay 
74a2ae208SSatish Balay #undef __FUNCT__
84a2ae208SSatish Balay #define __FUNCT__ "MatIncreaseOverlap_SeqBAIJ"
9690b6cddSBarry Smith PetscErrorCode MatIncreaseOverlap_SeqBAIJ(Mat A,PetscInt is_max,IS is[],PetscInt ov)
10a3192f15SSatish Balay {
11a3192f15SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
126849ba73SBarry Smith   PetscErrorCode ierr;
135d0c19d7SBarry Smith   PetscInt       row,i,j,k,l,m,n,*nidx,isz,val,ival;
145d0c19d7SBarry Smith   const PetscInt *idx;
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;
22d0f46423SBarry 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);
28d0f46423SBarry 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"
764aa3045dSJed Brown PetscErrorCode MatGetSubMatrix_SeqBAIJ_Private(Mat A,IS isrow,IS iscol,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;
825d0c19d7SBarry Smith   const PetscInt *irow,*icol;
835d0c19d7SBarry Smith   PetscInt       nrows,ncols,*ssmap,bs=A->rmap->bs,bs2=a->bs2;
84690b6cddSBarry Smith   PetscInt       *aj = a->j,*ai = a->i;
853f1db9ecSBarry Smith   MatScalar      *mat_a;
86736121d4SSatish Balay   Mat            C;
8714ca34e6SBarry Smith   PetscTruth     flag,sorted;
88736121d4SSatish Balay 
893a40ed3dSBarry Smith   PetscFunctionBegin;
9014ca34e6SBarry Smith   ierr = ISSorted(iscol,&sorted);CHKERRQ(ierr);
9114ca34e6SBarry Smith   if (!sorted) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"IS is not sorted");
92736121d4SSatish Balay 
93736121d4SSatish Balay   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
94218c64b6SSatish Balay   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
95b9b97703SBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
96b9b97703SBarry Smith   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
97736121d4SSatish Balay 
98690b6cddSBarry Smith   ierr = PetscMalloc((1+oldcols)*sizeof(PetscInt),&smap);CHKERRQ(ierr);
99736121d4SSatish Balay   ssmap = smap;
100690b6cddSBarry Smith   ierr = PetscMalloc((1+nrows)*sizeof(PetscInt),&lens);CHKERRQ(ierr);
101690b6cddSBarry Smith   ierr  = PetscMemzero(smap,oldcols*sizeof(PetscInt));CHKERRQ(ierr);
102736121d4SSatish Balay   for (i=0; i<ncols; i++) smap[icol[i]] = i+1;
103736121d4SSatish Balay   /* determine lens of each row */
104736121d4SSatish Balay   for (i=0; i<nrows; i++) {
105736121d4SSatish Balay     kstart  = ai[irow[i]];
106736121d4SSatish Balay     kend    = kstart + a->ilen[irow[i]];
107736121d4SSatish Balay     lens[i] = 0;
108736121d4SSatish Balay       for (k=kstart; k<kend; k++) {
109736121d4SSatish Balay         if (ssmap[aj[k]]) {
110736121d4SSatish Balay           lens[i]++;
111736121d4SSatish Balay         }
112736121d4SSatish Balay       }
113736121d4SSatish Balay     }
114736121d4SSatish Balay   /* Create and fill new matrix */
115736121d4SSatish Balay   if (scall == MAT_REUSE_MATRIX) {
116736121d4SSatish Balay     c = (Mat_SeqBAIJ *)((*B)->data);
117736121d4SSatish Balay 
118d0f46423SBarry Smith     if (c->mbs!=nrows || c->nbs!=ncols || (*B)->rmap->bs!=bs) SETERRQ(PETSC_ERR_ARG_SIZ,"Submatrix wrong size");
119690b6cddSBarry Smith     ierr = PetscMemcmp(c->ilen,lens,c->mbs *sizeof(PetscInt),&flag);CHKERRQ(ierr);
120abc0a331SBarry Smith     if (!flag) {
12129bbc08cSBarry Smith       SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong no of nonzeros");
122736121d4SSatish Balay     }
123690b6cddSBarry Smith     ierr = PetscMemzero(c->ilen,c->mbs*sizeof(PetscInt));CHKERRQ(ierr);
124736121d4SSatish Balay     C = *B;
1253a40ed3dSBarry Smith   } else {
1267adad957SLisandro Dalcin     ierr = MatCreate(((PetscObject)A)->comm,&C);CHKERRQ(ierr);
127f69a0ea3SMatthew Knepley     ierr = MatSetSizes(C,nrows*bs,ncols*bs,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1287adad957SLisandro Dalcin     ierr = MatSetType(C,((PetscObject)A)->type_name);CHKERRQ(ierr);
129ab93d7beSBarry Smith     ierr = MatSeqBAIJSetPreallocation_SeqBAIJ(C,bs,0,lens);CHKERRQ(ierr);
130736121d4SSatish Balay   }
131736121d4SSatish Balay   c = (Mat_SeqBAIJ *)(C->data);
132736121d4SSatish Balay   for (i=0; i<nrows; i++) {
133736121d4SSatish Balay     row    = irow[i];
134736121d4SSatish Balay     kstart = ai[row];
135736121d4SSatish Balay     kend   = kstart + a->ilen[row];
136736121d4SSatish Balay     mat_i  = c->i[i];
137736121d4SSatish Balay     mat_j  = c->j + mat_i;
138218c64b6SSatish Balay     mat_a  = c->a + mat_i*bs2;
139736121d4SSatish Balay     mat_ilen = c->ilen + i;
140736121d4SSatish Balay     for (k=kstart; k<kend; k++) {
141736121d4SSatish Balay       if ((tcol=ssmap[a->j[k]])) {
142736121d4SSatish Balay         *mat_j++ = tcol - 1;
143549d3d68SSatish Balay         ierr     = PetscMemcpy(mat_a,a->a+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
144549d3d68SSatish Balay         mat_a   += bs2;
145736121d4SSatish Balay         (*mat_ilen)++;
146736121d4SSatish Balay       }
147736121d4SSatish Balay     }
148736121d4SSatish Balay   }
149218c64b6SSatish Balay 
150736121d4SSatish Balay   /* Free work space */
151736121d4SSatish Balay   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
152606d414cSSatish Balay   ierr = PetscFree(smap);CHKERRQ(ierr);
153606d414cSSatish Balay   ierr = PetscFree(lens);CHKERRQ(ierr);
1546d4a8577SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1556d4a8577SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
156736121d4SSatish Balay 
157736121d4SSatish Balay   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
158736121d4SSatish Balay   *B = C;
1593a40ed3dSBarry Smith   PetscFunctionReturn(0);
160736121d4SSatish Balay }
161736121d4SSatish Balay 
1624a2ae208SSatish Balay #undef __FUNCT__
1634a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqBAIJ"
1644aa3045dSJed Brown PetscErrorCode MatGetSubMatrix_SeqBAIJ(Mat A,IS isrow,IS iscol,MatReuse scall,Mat *B)
165218c64b6SSatish Balay {
166218c64b6SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
167218c64b6SSatish Balay   IS             is1,is2;
1686849ba73SBarry Smith   PetscErrorCode ierr;
1695d0c19d7SBarry Smith   PetscInt       *vary,*iary,nrows,ncols,i,bs=A->rmap->bs,count;
1705d0c19d7SBarry Smith   const PetscInt *irow,*icol;
171218c64b6SSatish Balay 
1723a40ed3dSBarry Smith   PetscFunctionBegin;
173218c64b6SSatish Balay   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
174218c64b6SSatish Balay   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
175b9b97703SBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
176b9b97703SBarry Smith   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
177218c64b6SSatish Balay 
178218c64b6SSatish Balay   /* Verify if the indices corespond to each element in a block
179218c64b6SSatish Balay    and form the IS with compressed IS */
180690b6cddSBarry Smith   ierr = PetscMalloc(2*(a->mbs+1)*sizeof(PetscInt),&vary);CHKERRQ(ierr);
181218c64b6SSatish Balay   iary = vary + a->mbs;
182690b6cddSBarry 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++) {
186634064b4SBarry Smith     if (vary[i]!=0 && vary[i]!=bs) SETERRQ(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++) {
195634064b4SBarry Smith     if (vary[i]!=0 && vary[i]!=bs) SETERRQ(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);
201606d414cSSatish Balay   ierr = PetscFree(vary);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 /* -------------------------------------------------------*/
231d9eff348SSatish Balay #include "petscblaslapack.h"
2322d61bbb3SSatish Balay 
2334a2ae208SSatish Balay #undef __FUNCT__
2344a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_1"
235dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_1(Mat A,Vec xx,Vec zz)
2362d61bbb3SSatish Balay {
2372d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
238d9fead3dSBarry Smith   PetscScalar       *z,sum;
239d9fead3dSBarry Smith   const PetscScalar *x;
240d9fead3dSBarry Smith   const MatScalar   *v;
2416849ba73SBarry Smith   PetscErrorCode    ierr;
24298c9bda7SSatish Balay   PetscInt          mbs,i,*idx,*ii,n,*ridx=PETSC_NULL,nonzerorow=0;
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 
2492d61bbb3SSatish Balay   idx = a->j;
2502d61bbb3SSatish Balay   v   = a->a;
25126e093fcSHong Zhang   if (usecprow){
25226e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
25326e093fcSHong Zhang     ii   = a->compressedrow.i;
2547b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
25526e093fcSHong Zhang   } else {
25626e093fcSHong Zhang     mbs = a->mbs;
2572d61bbb3SSatish Balay     ii  = a->i;
25826e093fcSHong Zhang   }
2592d61bbb3SSatish Balay 
2602d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
2612d61bbb3SSatish Balay     n    = ii[1] - ii[0]; ii++;
2622d61bbb3SSatish Balay     sum  = 0.0;
26398c9bda7SSatish Balay     nonzerorow += (n>0);
2642d61bbb3SSatish Balay     while (n--) sum += *v++ * x[*idx++];
26526e093fcSHong Zhang     if (usecprow){
2667b2bb3b9SHong Zhang       z[ridx[i]] = sum;
26726e093fcSHong Zhang     } else {
2682d61bbb3SSatish Balay       z[i] = sum;
2692d61bbb3SSatish Balay     }
27026e093fcSHong Zhang   }
271d9fead3dSBarry Smith   ierr = VecRestoreArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
2721ebc52fbSHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
273dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*a->nz - nonzerorow);CHKERRQ(ierr);
2742d61bbb3SSatish Balay   PetscFunctionReturn(0);
2752d61bbb3SSatish Balay }
2762d61bbb3SSatish Balay 
2774a2ae208SSatish Balay #undef __FUNCT__
2784a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_2"
279dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_2(Mat A,Vec xx,Vec zz)
2802d61bbb3SSatish Balay {
2812d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
282d9fead3dSBarry Smith   PetscScalar       *z = 0,sum1,sum2,*zarray;
283d9fead3dSBarry Smith   const PetscScalar *x,*xb;
28487828ca2SBarry Smith   PetscScalar       x1,x2;
285d9fead3dSBarry Smith   const MatScalar   *v;
286dfbe8321SBarry Smith   PetscErrorCode    ierr;
28798c9bda7SSatish Balay   PetscInt          mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL,nonzerorow=0;
28826e093fcSHong Zhang   PetscTruth        usecprow=a->compressedrow.use;
2892d61bbb3SSatish Balay 
2902d61bbb3SSatish Balay   PetscFunctionBegin;
291d9fead3dSBarry Smith   ierr = VecGetArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
29226e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
2932d61bbb3SSatish Balay 
2942d61bbb3SSatish Balay   idx = a->j;
2952d61bbb3SSatish Balay   v   = a->a;
29626e093fcSHong Zhang   if (usecprow){
29726e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
29826e093fcSHong Zhang     ii   = a->compressedrow.i;
2997b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
30026e093fcSHong Zhang   } else {
30126e093fcSHong Zhang     mbs = a->mbs;
3022d61bbb3SSatish Balay     ii  = a->i;
30326e093fcSHong Zhang     z   = zarray;
30426e093fcSHong Zhang   }
3052d61bbb3SSatish Balay 
3062d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
3072d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
3082d61bbb3SSatish Balay     sum1 = 0.0; sum2 = 0.0;
30998c9bda7SSatish Balay     nonzerorow += (n>0);
3102d61bbb3SSatish Balay     for (j=0; j<n; j++) {
3112d61bbb3SSatish Balay       xb = x + 2*(*idx++); x1 = xb[0]; x2 = xb[1];
3122d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[2]*x2;
3132d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[3]*x2;
3142d61bbb3SSatish Balay       v += 4;
3152d61bbb3SSatish Balay     }
3167b2bb3b9SHong Zhang     if (usecprow) z = zarray + 2*ridx[i];
3172d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2;
31826e093fcSHong Zhang     if (!usecprow) z += 2;
3192d61bbb3SSatish Balay   }
320d9fead3dSBarry Smith   ierr = VecRestoreArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
32126e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
322dc0b31edSSatish Balay   ierr = PetscLogFlops(8.0*a->nz - 2.0*nonzerorow);CHKERRQ(ierr);
3232d61bbb3SSatish Balay   PetscFunctionReturn(0);
3242d61bbb3SSatish Balay }
3252d61bbb3SSatish Balay 
3264a2ae208SSatish Balay #undef __FUNCT__
3274a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_3"
328dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_3(Mat A,Vec xx,Vec zz)
3292d61bbb3SSatish Balay {
3302d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
331d9fead3dSBarry Smith   PetscScalar       *z = 0,sum1,sum2,sum3,x1,x2,x3,*zarray;
332d9fead3dSBarry Smith   const PetscScalar *x,*xb;
333d9fead3dSBarry Smith   const MatScalar   *v;
334dfbe8321SBarry Smith   PetscErrorCode    ierr;
335028cd4eaSSatish Balay   PetscInt          mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL,nonzerorow=0;
336028cd4eaSSatish Balay   PetscTruth        usecprow=a->compressedrow.use;
33726e093fcSHong Zhang 
3382d61bbb3SSatish Balay 
339b6410449SSatish Balay #if defined(PETSC_HAVE_PRAGMA_DISJOINT)
340fee21e36SBarry Smith #pragma disjoint(*v,*z,*xb)
341fee21e36SBarry Smith #endif
342fee21e36SBarry Smith 
3432d61bbb3SSatish Balay   PetscFunctionBegin;
344d9fead3dSBarry Smith   ierr = VecGetArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
34526e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
3462d61bbb3SSatish Balay 
3472d61bbb3SSatish Balay   idx = a->j;
3482d61bbb3SSatish Balay   v   = a->a;
34926e093fcSHong Zhang   if (usecprow){
35026e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
35126e093fcSHong Zhang     ii   = a->compressedrow.i;
3527b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
35326e093fcSHong Zhang   } else {
35426e093fcSHong Zhang     mbs = a->mbs;
3552d61bbb3SSatish Balay     ii  = a->i;
35626e093fcSHong Zhang     z   = zarray;
35726e093fcSHong Zhang   }
3582d61bbb3SSatish Balay 
3592d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
3602d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
3612d61bbb3SSatish Balay     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0;
36298c9bda7SSatish Balay     nonzerorow += (n>0);
3632d61bbb3SSatish Balay     for (j=0; j<n; j++) {
3642d61bbb3SSatish Balay       xb = x + 3*(*idx++); x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
3652d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[3]*x2 + v[6]*x3;
3662d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[4]*x2 + v[7]*x3;
3672d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[5]*x2 + v[8]*x3;
3682d61bbb3SSatish Balay       v += 9;
3692d61bbb3SSatish Balay     }
3707b2bb3b9SHong Zhang     if (usecprow) z = zarray + 3*ridx[i];
3712d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3;
37226e093fcSHong Zhang     if (!usecprow) z += 3;
3732d61bbb3SSatish Balay   }
374d9fead3dSBarry Smith   ierr = VecRestoreArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
37526e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
376dc0b31edSSatish Balay   ierr = PetscLogFlops(18.0*a->nz - 3.0*nonzerorow);CHKERRQ(ierr);
3772d61bbb3SSatish Balay   PetscFunctionReturn(0);
3782d61bbb3SSatish Balay }
3792d61bbb3SSatish Balay 
3804a2ae208SSatish Balay #undef __FUNCT__
3814a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_4"
382dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_4(Mat A,Vec xx,Vec zz)
3832d61bbb3SSatish Balay {
3842d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
385d9fead3dSBarry Smith   PetscScalar       *z = 0,sum1,sum2,sum3,sum4,x1,x2,x3,x4,*zarray;
386d9fead3dSBarry Smith   const PetscScalar *x,*xb;
387d9fead3dSBarry Smith   const MatScalar   *v;
388dfbe8321SBarry Smith   PetscErrorCode    ierr;
38998c9bda7SSatish Balay   PetscInt          mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL,nonzerorow=0;
39026e093fcSHong Zhang   PetscTruth        usecprow=a->compressedrow.use;
3912d61bbb3SSatish Balay 
3922d61bbb3SSatish Balay   PetscFunctionBegin;
393d9fead3dSBarry Smith   ierr = VecGetArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
39426e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
3952d61bbb3SSatish Balay 
3962d61bbb3SSatish Balay   idx = a->j;
3972d61bbb3SSatish Balay   v   = a->a;
39826e093fcSHong Zhang   if (usecprow){
39926e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
40026e093fcSHong Zhang     ii   = a->compressedrow.i;
4017b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
40226e093fcSHong Zhang   } else {
40326e093fcSHong Zhang     mbs = a->mbs;
4042d61bbb3SSatish Balay     ii  = a->i;
40526e093fcSHong Zhang     z   = zarray;
40626e093fcSHong Zhang   }
4072d61bbb3SSatish Balay 
4082d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
4092d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
4102d61bbb3SSatish Balay     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0;
41198c9bda7SSatish Balay     nonzerorow += (n>0);
4122d61bbb3SSatish Balay     for (j=0; j<n; j++) {
4132d61bbb3SSatish Balay       xb = x + 4*(*idx++);
4142d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
4152d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[4]*x2 + v[8]*x3   + v[12]*x4;
4162d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[5]*x2 + v[9]*x3   + v[13]*x4;
4172d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[6]*x2 + v[10]*x3  + v[14]*x4;
4182d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[7]*x2 + v[11]*x3  + v[15]*x4;
4192d61bbb3SSatish Balay       v += 16;
4202d61bbb3SSatish Balay     }
4217b2bb3b9SHong Zhang     if (usecprow) z = zarray + 4*ridx[i];
4222d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4;
42326e093fcSHong Zhang     if (!usecprow) z += 4;
4242d61bbb3SSatish Balay   }
425d9fead3dSBarry Smith   ierr = VecRestoreArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
42626e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
427dc0b31edSSatish Balay   ierr = PetscLogFlops(32.0*a->nz - 4.0*nonzerorow);CHKERRQ(ierr);
4282d61bbb3SSatish Balay   PetscFunctionReturn(0);
4292d61bbb3SSatish Balay }
4302d61bbb3SSatish Balay 
4314a2ae208SSatish Balay #undef __FUNCT__
4324a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_5"
433dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_5(Mat A,Vec xx,Vec zz)
4342d61bbb3SSatish Balay {
4352d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
436d9fead3dSBarry Smith   PetscScalar       sum1,sum2,sum3,sum4,sum5,x1,x2,x3,x4,x5,*z = 0,*zarray;
437d9fead3dSBarry Smith   const PetscScalar *xb,*x;
438d9fead3dSBarry Smith   const MatScalar   *v;
439dfbe8321SBarry Smith   PetscErrorCode    ierr;
44098c9bda7SSatish Balay   PetscInt          mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL,nonzerorow=0;
44126e093fcSHong Zhang   PetscTruth        usecprow=a->compressedrow.use;
4422d61bbb3SSatish Balay 
443433994e6SBarry Smith   PetscFunctionBegin;
444d9fead3dSBarry Smith   ierr = VecGetArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
44526e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
4462d61bbb3SSatish Balay 
4472d61bbb3SSatish Balay   idx = a->j;
4482d61bbb3SSatish Balay   v   = a->a;
44926e093fcSHong Zhang   if (usecprow){
45026e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
45126e093fcSHong Zhang     ii   = a->compressedrow.i;
4527b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
45326e093fcSHong Zhang   } else {
45426e093fcSHong Zhang     mbs = a->mbs;
4552d61bbb3SSatish Balay     ii  = a->i;
45626e093fcSHong Zhang     z   = zarray;
45726e093fcSHong Zhang   }
4582d61bbb3SSatish Balay 
4592d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
4602d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
4612d61bbb3SSatish Balay     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0;
46298c9bda7SSatish Balay     nonzerorow += (n>0);
4632d61bbb3SSatish Balay     for (j=0; j<n; j++) {
4642d61bbb3SSatish Balay       xb = x + 5*(*idx++);
4652d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4];
4662d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[5]*x2 + v[10]*x3  + v[15]*x4 + v[20]*x5;
4672d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[6]*x2 + v[11]*x3  + v[16]*x4 + v[21]*x5;
4682d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[7]*x2 + v[12]*x3  + v[17]*x4 + v[22]*x5;
4692d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[8]*x2 + v[13]*x3  + v[18]*x4 + v[23]*x5;
4702d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[9]*x2 + v[14]*x3  + v[19]*x4 + v[24]*x5;
4712d61bbb3SSatish Balay       v += 25;
4722d61bbb3SSatish Balay     }
4737b2bb3b9SHong Zhang     if (usecprow) z = zarray + 5*ridx[i];
4742d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5;
47526e093fcSHong Zhang     if (!usecprow) z += 5;
4762d61bbb3SSatish Balay   }
477d9fead3dSBarry Smith   ierr = VecRestoreArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
47826e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
479dc0b31edSSatish Balay   ierr = PetscLogFlops(50.0*a->nz - 5.0*nonzerorow);CHKERRQ(ierr);
4802d61bbb3SSatish Balay   PetscFunctionReturn(0);
4812d61bbb3SSatish Balay }
4822d61bbb3SSatish Balay 
48315091d37SBarry Smith 
4844a2ae208SSatish Balay #undef __FUNCT__
4854a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_6"
486dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_6(Mat A,Vec xx,Vec zz)
48715091d37SBarry Smith {
48815091d37SBarry Smith   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
489d9fead3dSBarry Smith   PetscScalar       *z = 0,sum1,sum2,sum3,sum4,sum5,sum6;
490d9fead3dSBarry Smith   const PetscScalar *x,*xb;
49126e093fcSHong Zhang   PetscScalar       x1,x2,x3,x4,x5,x6,*zarray;
492d9fead3dSBarry Smith   const MatScalar   *v;
493dfbe8321SBarry Smith   PetscErrorCode    ierr;
49498c9bda7SSatish Balay   PetscInt          mbs=a->mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL,nonzerorow=0;
49526e093fcSHong Zhang   PetscTruth        usecprow=a->compressedrow.use;
49615091d37SBarry Smith 
497433994e6SBarry Smith   PetscFunctionBegin;
498d9fead3dSBarry Smith   ierr = VecGetArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
49926e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
50015091d37SBarry Smith 
50115091d37SBarry Smith   idx = a->j;
50215091d37SBarry Smith   v   = a->a;
50326e093fcSHong Zhang   if (usecprow){
50426e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
50526e093fcSHong Zhang     ii   = a->compressedrow.i;
5067b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
50726e093fcSHong Zhang   } else {
50826e093fcSHong Zhang     mbs = a->mbs;
50915091d37SBarry Smith     ii  = a->i;
51026e093fcSHong Zhang     z   = zarray;
51126e093fcSHong Zhang   }
51215091d37SBarry Smith 
51315091d37SBarry Smith   for (i=0; i<mbs; i++) {
51415091d37SBarry Smith     n  = ii[1] - ii[0]; ii++;
51515091d37SBarry Smith     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0;
51698c9bda7SSatish Balay     nonzerorow += (n>0);
51715091d37SBarry Smith     for (j=0; j<n; j++) {
51815091d37SBarry Smith       xb = x + 6*(*idx++);
51915091d37SBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5];
52015091d37SBarry Smith       sum1 += v[0]*x1 + v[6]*x2  + v[12]*x3  + v[18]*x4 + v[24]*x5 + v[30]*x6;
52115091d37SBarry Smith       sum2 += v[1]*x1 + v[7]*x2  + v[13]*x3  + v[19]*x4 + v[25]*x5 + v[31]*x6;
52215091d37SBarry Smith       sum3 += v[2]*x1 + v[8]*x2  + v[14]*x3  + v[20]*x4 + v[26]*x5 + v[32]*x6;
52315091d37SBarry Smith       sum4 += v[3]*x1 + v[9]*x2  + v[15]*x3  + v[21]*x4 + v[27]*x5 + v[33]*x6;
52415091d37SBarry Smith       sum5 += v[4]*x1 + v[10]*x2 + v[16]*x3  + v[22]*x4 + v[28]*x5 + v[34]*x6;
52515091d37SBarry Smith       sum6 += v[5]*x1 + v[11]*x2 + v[17]*x3  + v[23]*x4 + v[29]*x5 + v[35]*x6;
52615091d37SBarry Smith       v += 36;
52715091d37SBarry Smith     }
5287b2bb3b9SHong Zhang     if (usecprow) z = zarray + 6*ridx[i];
52915091d37SBarry Smith     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6;
53026e093fcSHong Zhang     if (!usecprow) z += 6;
53115091d37SBarry Smith   }
53215091d37SBarry Smith 
533d9fead3dSBarry Smith   ierr = VecRestoreArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
53426e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
535dc0b31edSSatish Balay   ierr = PetscLogFlops(72.0*a->nz - 6.0*nonzerorow);CHKERRQ(ierr);
53615091d37SBarry Smith   PetscFunctionReturn(0);
53715091d37SBarry Smith }
5384a2ae208SSatish Balay #undef __FUNCT__
5394a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_7"
540dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_7(Mat A,Vec xx,Vec zz)
5412d61bbb3SSatish Balay {
5422d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
543d9fead3dSBarry Smith   PetscScalar       *z = 0,sum1,sum2,sum3,sum4,sum5,sum6,sum7;
544d9fead3dSBarry Smith   const PetscScalar *x,*xb;
54526e093fcSHong Zhang   PetscScalar       x1,x2,x3,x4,x5,x6,x7,*zarray;
546d9fead3dSBarry Smith   const MatScalar   *v;
547dfbe8321SBarry Smith   PetscErrorCode    ierr;
54898c9bda7SSatish Balay   PetscInt          mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL,nonzerorow=0;
54926e093fcSHong Zhang   PetscTruth        usecprow=a->compressedrow.use;
5502d61bbb3SSatish Balay 
551433994e6SBarry Smith   PetscFunctionBegin;
552d9fead3dSBarry Smith   ierr = VecGetArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
55326e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
5542d61bbb3SSatish Balay 
5552d61bbb3SSatish Balay   idx = a->j;
5562d61bbb3SSatish Balay   v   = a->a;
55726e093fcSHong Zhang   if (usecprow){
55826e093fcSHong Zhang     mbs    = a->compressedrow.nrows;
55926e093fcSHong Zhang     ii     = a->compressedrow.i;
5607b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
56126e093fcSHong Zhang   } else {
56226e093fcSHong Zhang     mbs = a->mbs;
5632d61bbb3SSatish Balay     ii  = a->i;
56426e093fcSHong Zhang     z   = zarray;
56526e093fcSHong Zhang   }
5662d61bbb3SSatish Balay 
5672d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
5682d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
5692d61bbb3SSatish Balay     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0; sum7 = 0.0;
57098c9bda7SSatish Balay     nonzerorow += (n>0);
5712d61bbb3SSatish Balay     for (j=0; j<n; j++) {
5722d61bbb3SSatish Balay       xb = x + 7*(*idx++);
5732d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6];
5742d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[7]*x2  + v[14]*x3  + v[21]*x4 + v[28]*x5 + v[35]*x6 + v[42]*x7;
5752d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[8]*x2  + v[15]*x3  + v[22]*x4 + v[29]*x5 + v[36]*x6 + v[43]*x7;
5762d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[9]*x2  + v[16]*x3  + v[23]*x4 + v[30]*x5 + v[37]*x6 + v[44]*x7;
5772d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[10]*x2 + v[17]*x3  + v[24]*x4 + v[31]*x5 + v[38]*x6 + v[45]*x7;
5782d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[11]*x2 + v[18]*x3  + v[25]*x4 + v[32]*x5 + v[39]*x6 + v[46]*x7;
5792d61bbb3SSatish Balay       sum6 += v[5]*x1 + v[12]*x2 + v[19]*x3  + v[26]*x4 + v[33]*x5 + v[40]*x6 + v[47]*x7;
5802d61bbb3SSatish Balay       sum7 += v[6]*x1 + v[13]*x2 + v[20]*x3  + v[27]*x4 + v[34]*x5 + v[41]*x6 + v[48]*x7;
5812d61bbb3SSatish Balay       v += 49;
5822d61bbb3SSatish Balay     }
5837b2bb3b9SHong Zhang     if (usecprow) z = zarray + 7*ridx[i];
5842d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
58526e093fcSHong Zhang     if (!usecprow) z += 7;
5862d61bbb3SSatish Balay   }
5872d61bbb3SSatish Balay 
588d9fead3dSBarry Smith   ierr = VecRestoreArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
58926e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
590dc0b31edSSatish Balay   ierr = PetscLogFlops(98.0*a->nz - 7.0*nonzerorow);CHKERRQ(ierr);
5912d61bbb3SSatish Balay   PetscFunctionReturn(0);
5922d61bbb3SSatish Balay }
5932d61bbb3SSatish Balay 
5943f1db9ecSBarry Smith /*
5953f1db9ecSBarry Smith     This will not work with MatScalar == float because it calls the BLAS
5963f1db9ecSBarry Smith */
5974a2ae208SSatish Balay #undef __FUNCT__
5984a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_N"
599dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_N(Mat A,Vec xx,Vec zz)
6002d61bbb3SSatish Balay {
6012d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
602dcf5cc72SBarry Smith   PetscScalar    *x,*z = 0,*xb,*work,*workt,*zarray;
6033f1db9ecSBarry Smith   MatScalar      *v;
604dfbe8321SBarry Smith   PetscErrorCode ierr;
605d0f46423SBarry Smith   PetscInt       mbs=a->mbs,i,*idx,*ii,bs=A->rmap->bs,j,n,bs2=a->bs2;
60698c9bda7SSatish Balay   PetscInt       ncols,k,*ridx=PETSC_NULL,nonzerorow=0;
60726e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
6082d61bbb3SSatish Balay 
6092d61bbb3SSatish Balay   PetscFunctionBegin;
6101ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
61126e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
6122d61bbb3SSatish Balay 
6132d61bbb3SSatish Balay   idx = a->j;
6142d61bbb3SSatish Balay   v   = a->a;
61526e093fcSHong Zhang   if (usecprow){
61626e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
61726e093fcSHong Zhang     ii   = a->compressedrow.i;
6187b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
61926e093fcSHong Zhang   } else {
62026e093fcSHong Zhang     mbs = a->mbs;
6212d61bbb3SSatish Balay     ii  = a->i;
62226e093fcSHong Zhang     z   = zarray;
62326e093fcSHong Zhang   }
624218c64b6SSatish Balay 
6252d61bbb3SSatish Balay   if (!a->mult_work) {
626d0f46423SBarry Smith     k    = PetscMax(A->rmap->n,A->cmap->n);
62787828ca2SBarry Smith     ierr = PetscMalloc((k+1)*sizeof(PetscScalar),&a->mult_work);CHKERRQ(ierr);
6282d61bbb3SSatish Balay   }
6292d61bbb3SSatish Balay   work = a->mult_work;
6302d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
6312d61bbb3SSatish Balay     n     = ii[1] - ii[0]; ii++;
6322d61bbb3SSatish Balay     ncols = n*bs;
6332d61bbb3SSatish Balay     workt = work;
63498c9bda7SSatish Balay     nonzerorow += (n>0);
6352d61bbb3SSatish Balay     for (j=0; j<n; j++) {
6362d61bbb3SSatish Balay       xb = x + bs*(*idx++);
6372d61bbb3SSatish Balay       for (k=0; k<bs; k++) workt[k] = xb[k];
6382d61bbb3SSatish Balay       workt += bs;
6392d61bbb3SSatish Balay     }
6407b2bb3b9SHong Zhang     if (usecprow) z = zarray + bs*ridx[i];
641737d121aSSatish Balay     Kernel_w_gets_Ar_times_v(bs,ncols,work,v,z);
64271044d3cSBarry Smith     /* BLASgemv_("N",&bs,&ncols,&_DOne,v,&bs,work,&_One,&_DZero,z,&_One); */
6432d61bbb3SSatish Balay     v += n*bs2;
64426e093fcSHong Zhang     if (!usecprow) z += bs;
6452d61bbb3SSatish Balay   }
6461ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
64726e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
648dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*a->nz*bs2 - bs*nonzerorow);CHKERRQ(ierr);
6492d61bbb3SSatish Balay   PetscFunctionReturn(0);
6502d61bbb3SSatish Balay }
6512d61bbb3SSatish Balay 
6526a65c766SHong Zhang extern PetscErrorCode VecCopy_Seq(Vec,Vec);
6534a2ae208SSatish Balay #undef __FUNCT__
6544a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_1"
655dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_1(Mat A,Vec xx,Vec yy,Vec zz)
6562d61bbb3SSatish Balay {
6572d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
658dcf5cc72SBarry Smith   PetscScalar    *x,*y = 0,*z = 0,sum,*yarray,*zarray;
6593f1db9ecSBarry Smith   MatScalar      *v;
660dfbe8321SBarry Smith   PetscErrorCode ierr;
6614eb6d288SHong Zhang   PetscInt       mbs=a->mbs,i,*idx,*ii,n,*ridx=PETSC_NULL;
66226e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
6632d61bbb3SSatish Balay 
6642d61bbb3SSatish Balay   PetscFunctionBegin;
6651ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
66626e093fcSHong Zhang   ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr);
6672e8a6d31SBarry Smith   if (zz != yy) {
66826e093fcSHong Zhang     ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
6692e8a6d31SBarry Smith   } else {
67026e093fcSHong Zhang     zarray = yarray;
6712e8a6d31SBarry Smith   }
6722d61bbb3SSatish Balay 
6732d61bbb3SSatish Balay   idx = a->j;
6742d61bbb3SSatish Balay   v   = a->a;
67526e093fcSHong Zhang   if (usecprow){
6764eb6d288SHong Zhang     if (zz != yy){
6774eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,mbs*sizeof(PetscScalar));CHKERRQ(ierr);
6784eb6d288SHong Zhang     }
67926e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
68026e093fcSHong Zhang     ii   = a->compressedrow.i;
6817b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
68226e093fcSHong Zhang   } else {
6832d61bbb3SSatish Balay     ii  = a->i;
68426e093fcSHong Zhang     y   = yarray;
68526e093fcSHong Zhang     z   = zarray;
68626e093fcSHong Zhang   }
6872d61bbb3SSatish Balay 
6882d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
6892d61bbb3SSatish Balay     n    = ii[1] - ii[0]; ii++;
69026e093fcSHong Zhang     if (usecprow){
6917b2bb3b9SHong Zhang       z = zarray + ridx[i];
6927b2bb3b9SHong Zhang       y = yarray + ridx[i];
69326e093fcSHong Zhang     }
69426e093fcSHong Zhang     sum = y[0];
6952d61bbb3SSatish Balay     while (n--) sum += *v++ * x[*idx++];
69626e093fcSHong Zhang     z[0] = sum;
69726e093fcSHong Zhang     if (!usecprow){
69826e093fcSHong Zhang       z++; y++;
69926e093fcSHong Zhang     }
7002d61bbb3SSatish Balay   }
7011ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
70226e093fcSHong Zhang   ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr);
7032e8a6d31SBarry Smith   if (zz != yy) {
70426e093fcSHong Zhang     ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
7052e8a6d31SBarry Smith   }
706dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*a->nz);CHKERRQ(ierr);
7072d61bbb3SSatish Balay   PetscFunctionReturn(0);
7082d61bbb3SSatish Balay }
7092d61bbb3SSatish Balay 
7104a2ae208SSatish Balay #undef __FUNCT__
7114a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_2"
712dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_2(Mat A,Vec xx,Vec yy,Vec zz)
7132d61bbb3SSatish Balay {
7142d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
715dcf5cc72SBarry Smith   PetscScalar    *x,*y = 0,*z = 0,*xb,sum1,sum2;
71626e093fcSHong Zhang   PetscScalar    x1,x2,*yarray,*zarray;
7173f1db9ecSBarry Smith   MatScalar      *v;
718dfbe8321SBarry Smith   PetscErrorCode ierr;
7194eb6d288SHong Zhang   PetscInt       mbs=a->mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL;
72026e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
7212d61bbb3SSatish Balay 
7222d61bbb3SSatish Balay   PetscFunctionBegin;
7231ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
72426e093fcSHong Zhang   ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr);
7252e8a6d31SBarry Smith   if (zz != yy) {
72626e093fcSHong Zhang     ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
7272e8a6d31SBarry Smith   } else {
72826e093fcSHong Zhang     zarray = yarray;
7292e8a6d31SBarry Smith   }
7302d61bbb3SSatish Balay 
7312d61bbb3SSatish Balay   idx = a->j;
7322d61bbb3SSatish Balay   v   = a->a;
73326e093fcSHong Zhang   if (usecprow){
7344eb6d288SHong Zhang     if (zz != yy){
7354eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,2*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
7364eb6d288SHong Zhang     }
73726e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
73826e093fcSHong Zhang     ii   = a->compressedrow.i;
7397b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
7404eb6d288SHong Zhang     if (zz != yy){
7414eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,a->mbs*sizeof(PetscScalar));CHKERRQ(ierr);
7424eb6d288SHong Zhang     }
74326e093fcSHong Zhang   } else {
7442d61bbb3SSatish Balay     ii  = a->i;
74526e093fcSHong Zhang     y   = yarray;
74626e093fcSHong Zhang     z   = zarray;
74726e093fcSHong Zhang   }
7482d61bbb3SSatish Balay 
7492d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
7502d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
75126e093fcSHong Zhang     if (usecprow){
7527b2bb3b9SHong Zhang       z = zarray + 2*ridx[i];
7537b2bb3b9SHong Zhang       y = yarray + 2*ridx[i];
75426e093fcSHong Zhang     }
7552d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1];
7562d61bbb3SSatish Balay     for (j=0; j<n; j++) {
7572d61bbb3SSatish Balay       xb = x + 2*(*idx++); x1 = xb[0]; x2 = xb[1];
7582d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[2]*x2;
7592d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[3]*x2;
7602d61bbb3SSatish Balay       v += 4;
7612d61bbb3SSatish Balay     }
7622d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2;
76326e093fcSHong Zhang     if (!usecprow){
7642d61bbb3SSatish Balay       z += 2; y += 2;
7652d61bbb3SSatish Balay     }
76626e093fcSHong Zhang   }
7671ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
76826e093fcSHong Zhang   ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr);
7692e8a6d31SBarry Smith   if (zz != yy) {
77026e093fcSHong Zhang     ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
7712e8a6d31SBarry Smith   }
772dc0b31edSSatish Balay   ierr = PetscLogFlops(4.0*a->nz);CHKERRQ(ierr);
7732d61bbb3SSatish Balay   PetscFunctionReturn(0);
7742d61bbb3SSatish Balay }
7752d61bbb3SSatish Balay 
7764a2ae208SSatish Balay #undef __FUNCT__
7774a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_3"
778dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_3(Mat A,Vec xx,Vec yy,Vec zz)
7792d61bbb3SSatish Balay {
7802d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
781dcf5cc72SBarry Smith   PetscScalar    *x,*y = 0,*z = 0,*xb,sum1,sum2,sum3,x1,x2,x3,*yarray,*zarray;
7823f1db9ecSBarry Smith   MatScalar      *v;
783dfbe8321SBarry Smith   PetscErrorCode ierr;
7844eb6d288SHong Zhang   PetscInt       mbs=a->mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL;
78526e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
7862d61bbb3SSatish Balay 
7872d61bbb3SSatish Balay   PetscFunctionBegin;
7881ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
78926e093fcSHong Zhang   ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr);
7902e8a6d31SBarry Smith   if (zz != yy) {
79126e093fcSHong Zhang     ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
7922e8a6d31SBarry Smith   } else {
79326e093fcSHong Zhang     zarray = yarray;
7942e8a6d31SBarry Smith   }
7952d61bbb3SSatish Balay 
7962d61bbb3SSatish Balay   idx = a->j;
7972d61bbb3SSatish Balay   v   = a->a;
79826e093fcSHong Zhang   if (usecprow){
7994eb6d288SHong Zhang     if (zz != yy){
8004eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,3*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
8014eb6d288SHong Zhang     }
80226e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
80326e093fcSHong Zhang     ii   = a->compressedrow.i;
8047b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
80526e093fcSHong Zhang   } else {
8062d61bbb3SSatish Balay     ii  = a->i;
80726e093fcSHong Zhang     y   = yarray;
80826e093fcSHong Zhang     z   = zarray;
80926e093fcSHong Zhang   }
8102d61bbb3SSatish Balay 
8112d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
8122d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
81326e093fcSHong Zhang     if (usecprow){
8147b2bb3b9SHong Zhang       z = zarray + 3*ridx[i];
8157b2bb3b9SHong Zhang       y = yarray + 3*ridx[i];
81626e093fcSHong Zhang     }
8172d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2];
8182d61bbb3SSatish Balay     for (j=0; j<n; j++) {
8192d61bbb3SSatish Balay       xb = x + 3*(*idx++); x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
8202d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[3]*x2 + v[6]*x3;
8212d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[4]*x2 + v[7]*x3;
8222d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[5]*x2 + v[8]*x3;
8232d61bbb3SSatish Balay       v += 9;
8242d61bbb3SSatish Balay     }
8252d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3;
82626e093fcSHong Zhang     if (!usecprow){
8272d61bbb3SSatish Balay       z += 3; y += 3;
8282d61bbb3SSatish Balay     }
82926e093fcSHong Zhang   }
8301ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
83126e093fcSHong Zhang   ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr);
8322e8a6d31SBarry Smith   if (zz != yy) {
83326e093fcSHong Zhang     ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
8342e8a6d31SBarry Smith   }
835dc0b31edSSatish Balay   ierr = PetscLogFlops(18.0*a->nz);CHKERRQ(ierr);
8362d61bbb3SSatish Balay   PetscFunctionReturn(0);
8372d61bbb3SSatish Balay }
8382d61bbb3SSatish Balay 
8394a2ae208SSatish Balay #undef __FUNCT__
8404a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_4"
841dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_4(Mat A,Vec xx,Vec yy,Vec zz)
8422d61bbb3SSatish Balay {
8432d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
844dcf5cc72SBarry Smith   PetscScalar    *x,*y = 0,*z = 0,*xb,sum1,sum2,sum3,sum4,x1,x2,x3,x4,*yarray,*zarray;
8453f1db9ecSBarry Smith   MatScalar      *v;
846dfbe8321SBarry Smith   PetscErrorCode ierr;
8474eb6d288SHong Zhang   PetscInt       mbs=a->mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL;
84826e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
8492d61bbb3SSatish Balay 
8502d61bbb3SSatish Balay   PetscFunctionBegin;
8511ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
85226e093fcSHong Zhang   ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr);
8532e8a6d31SBarry Smith   if (zz != yy) {
85426e093fcSHong Zhang     ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
8552e8a6d31SBarry Smith   } else {
85626e093fcSHong Zhang     zarray = yarray;
8572e8a6d31SBarry Smith   }
8582d61bbb3SSatish Balay 
8592d61bbb3SSatish Balay   idx   = a->j;
8602d61bbb3SSatish Balay   v     = a->a;
86126e093fcSHong Zhang   if (usecprow){
8624eb6d288SHong Zhang     if (zz != yy){
8634eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,4*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
8644eb6d288SHong Zhang     }
86526e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
86626e093fcSHong Zhang     ii   = a->compressedrow.i;
8677b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
86826e093fcSHong Zhang   } else {
8692d61bbb3SSatish Balay     ii  = a->i;
87026e093fcSHong Zhang     y   = yarray;
87126e093fcSHong Zhang     z   = zarray;
87226e093fcSHong Zhang   }
8732d61bbb3SSatish Balay 
8742d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
8752d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
87626e093fcSHong Zhang     if (usecprow){
8777b2bb3b9SHong Zhang       z = zarray + 4*ridx[i];
8787b2bb3b9SHong Zhang       y = yarray + 4*ridx[i];
87926e093fcSHong Zhang     }
8802d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3];
8812d61bbb3SSatish Balay     for (j=0; j<n; j++) {
8822d61bbb3SSatish Balay       xb = x + 4*(*idx++);
8832d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
8842d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[4]*x2 + v[8]*x3   + v[12]*x4;
8852d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[5]*x2 + v[9]*x3   + v[13]*x4;
8862d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[6]*x2 + v[10]*x3  + v[14]*x4;
8872d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[7]*x2 + v[11]*x3  + v[15]*x4;
8882d61bbb3SSatish Balay       v += 16;
8892d61bbb3SSatish Balay     }
8902d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4;
89126e093fcSHong Zhang     if (!usecprow){
8922d61bbb3SSatish Balay       z += 4; y += 4;
8932d61bbb3SSatish Balay     }
89426e093fcSHong Zhang   }
8951ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
89626e093fcSHong Zhang   ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr);
8972e8a6d31SBarry Smith   if (zz != yy) {
89826e093fcSHong Zhang     ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
8992e8a6d31SBarry Smith   }
900dc0b31edSSatish Balay   ierr = PetscLogFlops(32.0*a->nz);CHKERRQ(ierr);
9012d61bbb3SSatish Balay   PetscFunctionReturn(0);
9022d61bbb3SSatish Balay }
9032d61bbb3SSatish Balay 
9044a2ae208SSatish Balay #undef __FUNCT__
9054a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_5"
906dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_5(Mat A,Vec xx,Vec yy,Vec zz)
9072d61bbb3SSatish Balay {
9082d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
909dcf5cc72SBarry Smith   PetscScalar    *x,*y = 0,*z = 0,*xb,sum1,sum2,sum3,sum4,sum5,x1,x2,x3,x4,x5;
91026e093fcSHong Zhang   PetscScalar    *yarray,*zarray;
9113f1db9ecSBarry Smith   MatScalar      *v;
912dfbe8321SBarry Smith   PetscErrorCode ierr;
9134eb6d288SHong Zhang   PetscInt       mbs=a->mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL;
91426e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
9152d61bbb3SSatish Balay 
9162d61bbb3SSatish Balay   PetscFunctionBegin;
9171ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
91826e093fcSHong Zhang   ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr);
9192e8a6d31SBarry Smith   if (zz != yy) {
92026e093fcSHong Zhang     ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
9212e8a6d31SBarry Smith   } else {
92226e093fcSHong Zhang     zarray = yarray;
9232e8a6d31SBarry Smith   }
9242d61bbb3SSatish Balay 
9252d61bbb3SSatish Balay   idx = a->j;
9262d61bbb3SSatish Balay   v   = a->a;
92726e093fcSHong Zhang   if (usecprow){
9284eb6d288SHong Zhang     if (zz != yy){
9294eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,5*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
9304eb6d288SHong Zhang     }
93126e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
93226e093fcSHong Zhang     ii   = a->compressedrow.i;
9337b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
93426e093fcSHong Zhang   } else {
9352d61bbb3SSatish Balay     ii  = a->i;
93626e093fcSHong Zhang     y   = yarray;
93726e093fcSHong Zhang     z   = zarray;
93826e093fcSHong Zhang   }
9392d61bbb3SSatish Balay 
9402d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
9412d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
94226e093fcSHong Zhang     if (usecprow){
9437b2bb3b9SHong Zhang       z = zarray + 5*ridx[i];
9447b2bb3b9SHong Zhang       y = yarray + 5*ridx[i];
94526e093fcSHong Zhang     }
9462d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4];
9472d61bbb3SSatish Balay     for (j=0; j<n; j++) {
9482d61bbb3SSatish Balay       xb = x + 5*(*idx++);
9492d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4];
9502d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[5]*x2 + v[10]*x3  + v[15]*x4 + v[20]*x5;
9512d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[6]*x2 + v[11]*x3  + v[16]*x4 + v[21]*x5;
9522d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[7]*x2 + v[12]*x3  + v[17]*x4 + v[22]*x5;
9532d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[8]*x2 + v[13]*x3  + v[18]*x4 + v[23]*x5;
9542d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[9]*x2 + v[14]*x3  + v[19]*x4 + v[24]*x5;
9552d61bbb3SSatish Balay       v += 25;
9562d61bbb3SSatish Balay     }
9572d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5;
95826e093fcSHong Zhang     if (!usecprow){
9592d61bbb3SSatish Balay       z += 5; y += 5;
9602d61bbb3SSatish Balay     }
96126e093fcSHong Zhang   }
9621ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
96326e093fcSHong Zhang   ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr);
9642e8a6d31SBarry Smith   if (zz != yy) {
96526e093fcSHong Zhang     ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
9662e8a6d31SBarry Smith   }
967dc0b31edSSatish Balay   ierr = PetscLogFlops(50.0*a->nz);CHKERRQ(ierr);
9682d61bbb3SSatish Balay   PetscFunctionReturn(0);
9692d61bbb3SSatish Balay }
9704a2ae208SSatish Balay #undef __FUNCT__
9714a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_6"
972dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_6(Mat A,Vec xx,Vec yy,Vec zz)
97315091d37SBarry Smith {
97415091d37SBarry Smith   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
975dcf5cc72SBarry Smith   PetscScalar    *x,*y = 0,*z = 0,*xb,sum1,sum2,sum3,sum4,sum5,sum6;
97626e093fcSHong Zhang   PetscScalar    x1,x2,x3,x4,x5,x6,*yarray,*zarray;
97715091d37SBarry Smith   MatScalar      *v;
978dfbe8321SBarry Smith   PetscErrorCode ierr;
9794eb6d288SHong Zhang   PetscInt       mbs=a->mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL;
98026e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
98115091d37SBarry Smith 
98215091d37SBarry Smith   PetscFunctionBegin;
9831ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
98426e093fcSHong Zhang   ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr);
98515091d37SBarry Smith   if (zz != yy) {
98626e093fcSHong Zhang     ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
98715091d37SBarry Smith   } else {
98826e093fcSHong Zhang     zarray = yarray;
98915091d37SBarry Smith   }
99015091d37SBarry Smith 
99115091d37SBarry Smith   idx = a->j;
99215091d37SBarry Smith   v   = a->a;
99326e093fcSHong Zhang   if (usecprow){
9944eb6d288SHong Zhang     if (zz != yy){
9954eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,6*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
9964eb6d288SHong Zhang     }
99726e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
99826e093fcSHong Zhang     ii   = a->compressedrow.i;
9997b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
100026e093fcSHong Zhang   } else {
100115091d37SBarry Smith     ii  = a->i;
100226e093fcSHong Zhang     y   = yarray;
100326e093fcSHong Zhang     z   = zarray;
100426e093fcSHong Zhang   }
100515091d37SBarry Smith 
100615091d37SBarry Smith   for (i=0; i<mbs; i++) {
100715091d37SBarry Smith     n  = ii[1] - ii[0]; ii++;
100826e093fcSHong Zhang     if (usecprow){
10097b2bb3b9SHong Zhang       z = zarray + 6*ridx[i];
10107b2bb3b9SHong Zhang       y = yarray + 6*ridx[i];
101126e093fcSHong Zhang     }
101215091d37SBarry Smith     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4]; sum6 = y[5];
101315091d37SBarry Smith     for (j=0; j<n; j++) {
10143b95cb0eSSatish Balay       xb = x + 6*(*idx++);
101515091d37SBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5];
101615091d37SBarry Smith       sum1 += v[0]*x1 + v[6]*x2  + v[12]*x3  + v[18]*x4 + v[24]*x5 + v[30]*x6;
101715091d37SBarry Smith       sum2 += v[1]*x1 + v[7]*x2  + v[13]*x3  + v[19]*x4 + v[25]*x5 + v[31]*x6;
101815091d37SBarry Smith       sum3 += v[2]*x1 + v[8]*x2  + v[14]*x3  + v[20]*x4 + v[26]*x5 + v[32]*x6;
101915091d37SBarry Smith       sum4 += v[3]*x1 + v[9]*x2  + v[15]*x3  + v[21]*x4 + v[27]*x5 + v[33]*x6;
102015091d37SBarry Smith       sum5 += v[4]*x1 + v[10]*x2 + v[16]*x3  + v[22]*x4 + v[28]*x5 + v[34]*x6;
102115091d37SBarry Smith       sum6 += v[5]*x1 + v[11]*x2 + v[17]*x3  + v[23]*x4 + v[29]*x5 + v[35]*x6;
102215091d37SBarry Smith       v += 36;
102315091d37SBarry Smith     }
102415091d37SBarry Smith     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6;
102526e093fcSHong Zhang     if (!usecprow){
102615091d37SBarry Smith       z += 6; y += 6;
102715091d37SBarry Smith     }
102826e093fcSHong Zhang   }
10291ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
103026e093fcSHong Zhang   ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr);
103115091d37SBarry Smith   if (zz != yy) {
103226e093fcSHong Zhang     ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
103315091d37SBarry Smith   }
1034dc0b31edSSatish Balay   ierr = PetscLogFlops(72.0*a->nz);CHKERRQ(ierr);
103515091d37SBarry Smith   PetscFunctionReturn(0);
103615091d37SBarry Smith }
10372d61bbb3SSatish Balay 
10384a2ae208SSatish Balay #undef __FUNCT__
10394a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_7"
1040dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_7(Mat A,Vec xx,Vec yy,Vec zz)
10412d61bbb3SSatish Balay {
10422d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
1043dcf5cc72SBarry Smith   PetscScalar    *x,*y = 0,*z = 0,*xb,sum1,sum2,sum3,sum4,sum5,sum6,sum7;
104426e093fcSHong Zhang   PetscScalar    x1,x2,x3,x4,x5,x6,x7,*yarray,*zarray;
10453f1db9ecSBarry Smith   MatScalar      *v;
1046dfbe8321SBarry Smith   PetscErrorCode ierr;
10477b2bb3b9SHong Zhang   PetscInt       mbs=a->mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL;
104826e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
10492d61bbb3SSatish Balay 
10502d61bbb3SSatish Balay   PetscFunctionBegin;
10511ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
105226e093fcSHong Zhang   ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr);
10532e8a6d31SBarry Smith   if (zz != yy) {
105426e093fcSHong Zhang     ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
10552e8a6d31SBarry Smith   } else {
105626e093fcSHong Zhang     zarray = yarray;
10572e8a6d31SBarry Smith   }
10582d61bbb3SSatish Balay 
10592d61bbb3SSatish Balay   idx = a->j;
10602d61bbb3SSatish Balay   v   = a->a;
106126e093fcSHong Zhang   if (usecprow){
10624eb6d288SHong Zhang     if (zz != yy){
10634eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,7*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
10644eb6d288SHong Zhang     }
106526e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
106626e093fcSHong Zhang     ii   = a->compressedrow.i;
10677b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
106826e093fcSHong Zhang   } else {
10692d61bbb3SSatish Balay     ii  = a->i;
107026e093fcSHong Zhang     y   = yarray;
107126e093fcSHong Zhang     z   = zarray;
107226e093fcSHong Zhang   }
10732d61bbb3SSatish Balay 
10742d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
10752d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
107626e093fcSHong Zhang     if (usecprow){
10777b2bb3b9SHong Zhang       z = zarray + 7*ridx[i];
10787b2bb3b9SHong Zhang       y = yarray + 7*ridx[i];
107926e093fcSHong Zhang     }
10802d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4]; sum6 = y[5]; sum7 = y[6];
10812d61bbb3SSatish Balay     for (j=0; j<n; j++) {
10822d61bbb3SSatish Balay       xb = x + 7*(*idx++);
10832d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6];
10842d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[7]*x2  + v[14]*x3  + v[21]*x4 + v[28]*x5 + v[35]*x6 + v[42]*x7;
10852d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[8]*x2  + v[15]*x3  + v[22]*x4 + v[29]*x5 + v[36]*x6 + v[43]*x7;
10862d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[9]*x2  + v[16]*x3  + v[23]*x4 + v[30]*x5 + v[37]*x6 + v[44]*x7;
10872d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[10]*x2 + v[17]*x3  + v[24]*x4 + v[31]*x5 + v[38]*x6 + v[45]*x7;
10882d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[11]*x2 + v[18]*x3  + v[25]*x4 + v[32]*x5 + v[39]*x6 + v[46]*x7;
10892d61bbb3SSatish Balay       sum6 += v[5]*x1 + v[12]*x2 + v[19]*x3  + v[26]*x4 + v[33]*x5 + v[40]*x6 + v[47]*x7;
10902d61bbb3SSatish Balay       sum7 += v[6]*x1 + v[13]*x2 + v[20]*x3  + v[27]*x4 + v[34]*x5 + v[41]*x6 + v[48]*x7;
10912d61bbb3SSatish Balay       v += 49;
10922d61bbb3SSatish Balay     }
10932d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
109426e093fcSHong Zhang     if (!usecprow){
10952d61bbb3SSatish Balay       z += 7; y += 7;
10962d61bbb3SSatish Balay     }
109726e093fcSHong Zhang   }
10981ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
109926e093fcSHong Zhang   ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr);
11002e8a6d31SBarry Smith   if (zz != yy) {
110126e093fcSHong Zhang     ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
11022e8a6d31SBarry Smith   }
1103dc0b31edSSatish Balay   ierr = PetscLogFlops(98.0*a->nz);CHKERRQ(ierr);
11042d61bbb3SSatish Balay   PetscFunctionReturn(0);
11052d61bbb3SSatish Balay }
1106218c64b6SSatish Balay 
11074a2ae208SSatish Balay #undef __FUNCT__
11084a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_N"
1109dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_N(Mat A,Vec xx,Vec yy,Vec zz)
11102d61bbb3SSatish Balay {
11112d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
1112bba805e6SBarry Smith   PetscScalar    *x,*z = 0,*xb,*work,*workt,*zarray;
11133f1db9ecSBarry Smith   MatScalar      *v;
11146849ba73SBarry Smith   PetscErrorCode ierr;
1115d0f46423SBarry Smith   PetscInt       mbs,i,*idx,*ii,bs=A->rmap->bs,j,n,bs2=a->bs2;
11167b2bb3b9SHong Zhang   PetscInt       ncols,k,*ridx=PETSC_NULL;
111726e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
1118218c64b6SSatish Balay 
11192d61bbb3SSatish Balay   PetscFunctionBegin;
11206a65c766SHong Zhang   ierr = VecCopy_Seq(yy,zz);CHKERRQ(ierr);
11211ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
112226e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
11232d61bbb3SSatish Balay 
11242d61bbb3SSatish Balay   idx = a->j;
11252d61bbb3SSatish Balay   v   = a->a;
112626e093fcSHong Zhang   if (usecprow){
112726e093fcSHong Zhang     mbs    = a->compressedrow.nrows;
112826e093fcSHong Zhang     ii     = a->compressedrow.i;
11297b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
113026e093fcSHong Zhang   } else {
113126e093fcSHong Zhang     mbs = a->mbs;
11322d61bbb3SSatish Balay     ii  = a->i;
113326e093fcSHong Zhang     z   = zarray;
113426e093fcSHong Zhang   }
11352d61bbb3SSatish Balay 
11362d61bbb3SSatish Balay   if (!a->mult_work) {
1137d0f46423SBarry Smith     k    = PetscMax(A->rmap->n,A->cmap->n);
113887828ca2SBarry Smith     ierr = PetscMalloc((k+1)*sizeof(PetscScalar),&a->mult_work);CHKERRQ(ierr);
11392d61bbb3SSatish Balay   }
11402d61bbb3SSatish Balay   work = a->mult_work;
11412d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
11422d61bbb3SSatish Balay     n     = ii[1] - ii[0]; ii++;
11432d61bbb3SSatish Balay     ncols = n*bs;
11442d61bbb3SSatish Balay     workt = work;
11452d61bbb3SSatish Balay     for (j=0; j<n; j++) {
11462d61bbb3SSatish Balay       xb = x + bs*(*idx++);
11472d61bbb3SSatish Balay       for (k=0; k<bs; k++) workt[k] = xb[k];
11482d61bbb3SSatish Balay       workt += bs;
11492d61bbb3SSatish Balay     }
11507b2bb3b9SHong Zhang     if (usecprow) z = zarray + bs*ridx[i];
11513f1db9ecSBarry Smith     Kernel_w_gets_w_plus_Ar_times_v(bs,ncols,work,v,z);
115271044d3cSBarry Smith     /* BLASgemv_("N",&bs,&ncols,&_DOne,v,&bs,work,&_One,&_DOne,z,&_One); */
11532d61bbb3SSatish Balay     v += n*bs2;
115426e093fcSHong Zhang     if (!usecprow){
11552d61bbb3SSatish Balay       z += bs;
11562d61bbb3SSatish Balay     }
115726e093fcSHong Zhang   }
11581ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
115926e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
1160dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*a->nz*bs2);CHKERRQ(ierr);
11612d61bbb3SSatish Balay   PetscFunctionReturn(0);
11622d61bbb3SSatish Balay }
11632d61bbb3SSatish Balay 
11644a2ae208SSatish Balay #undef __FUNCT__
11654a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_SeqBAIJ"
1166dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_SeqBAIJ(Mat A,Vec xx,Vec zz)
11672d61bbb3SSatish Balay {
11683447b6efSHong Zhang   PetscScalar    zero = 0.0;
11696849ba73SBarry Smith   PetscErrorCode ierr;
11702d61bbb3SSatish Balay 
11712d61bbb3SSatish Balay   PetscFunctionBegin;
11722dcb1b2aSMatthew Knepley   ierr = VecSet(zz,zero);CHKERRQ(ierr);
11733447b6efSHong Zhang   ierr = MatMultTransposeAdd_SeqBAIJ(A,xx,zz,zz);CHKERRQ(ierr);
11742d61bbb3SSatish Balay   PetscFunctionReturn(0);
11752d61bbb3SSatish Balay }
11762d61bbb3SSatish Balay 
11774a2ae208SSatish Balay #undef __FUNCT__
11784a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_SeqBAIJ"
1179dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_SeqBAIJ(Mat A,Vec xx,Vec yy,Vec zz)
11802d61bbb3SSatish Balay 
11812d61bbb3SSatish Balay {
11822d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
1183dcf5cc72SBarry Smith   PetscScalar       *zb,*x,*z,*xb = 0,x1,x2,x3,x4,x5;
11843f1db9ecSBarry Smith   MatScalar         *v;
11856849ba73SBarry Smith   PetscErrorCode    ierr;
1186d0f46423SBarry Smith   PetscInt          mbs,i,*idx,*ii,rval,bs=A->rmap->bs,j,n,bs2=a->bs2,*ib,*ridx=PETSC_NULL;
11873447b6efSHong Zhang   Mat_CompressedRow cprow = a->compressedrow;
11883447b6efSHong Zhang   PetscTruth        usecprow=cprow.use;
11892d61bbb3SSatish Balay 
11902d61bbb3SSatish Balay   PetscFunctionBegin;
11916a65c766SHong Zhang   if (yy != zz) { ierr = VecCopy_Seq(yy,zz);CHKERRQ(ierr); }
11923447b6efSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
11933447b6efSHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
11942d61bbb3SSatish Balay 
11952d61bbb3SSatish Balay   idx = a->j;
11962d61bbb3SSatish Balay   v   = a->a;
11973447b6efSHong Zhang   if (usecprow){
11983447b6efSHong Zhang     mbs  = cprow.nrows;
11993447b6efSHong Zhang     ii   = cprow.i;
12007b2bb3b9SHong Zhang     ridx = cprow.rindex;
12013447b6efSHong Zhang   } else {
12023447b6efSHong Zhang     mbs=a->mbs;
12032d61bbb3SSatish Balay     ii = a->i;
1204f1af5d2fSBarry Smith     xb = x;
12053447b6efSHong Zhang   }
12062d61bbb3SSatish Balay 
12072d61bbb3SSatish Balay   switch (bs) {
12082d61bbb3SSatish Balay   case 1:
12092d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
12107b2bb3b9SHong Zhang       if (usecprow) xb = x + ridx[i];
1211f1af5d2fSBarry Smith       x1 = xb[0];
12123447b6efSHong Zhang       ib = idx + ii[0];
12133447b6efSHong Zhang       n  = ii[1] - ii[0]; ii++;
12142d61bbb3SSatish Balay       for (j=0; j<n; j++) {
12152d61bbb3SSatish Balay         rval    = ib[j];
1216f1af5d2fSBarry Smith         z[rval] += *v * x1;
1217f1af5d2fSBarry Smith         v++;
12182d61bbb3SSatish Balay       }
12193447b6efSHong Zhang       if (!usecprow) xb++;
12202d61bbb3SSatish Balay     }
12212d61bbb3SSatish Balay     break;
12222d61bbb3SSatish Balay   case 2:
12232d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
12247b2bb3b9SHong Zhang       if (usecprow) xb = x + 2*ridx[i];
1225f1af5d2fSBarry Smith       x1 = xb[0]; x2 = xb[1];
12263447b6efSHong Zhang       ib = idx + ii[0];
12273447b6efSHong Zhang       n  = ii[1] - ii[0]; ii++;
12282d61bbb3SSatish Balay       for (j=0; j<n; j++) {
12292d61bbb3SSatish Balay         rval      = ib[j]*2;
12302d61bbb3SSatish Balay         z[rval++] += v[0]*x1 + v[1]*x2;
12312d61bbb3SSatish Balay         z[rval++] += v[2]*x1 + v[3]*x2;
12322d61bbb3SSatish Balay         v  += 4;
12332d61bbb3SSatish Balay       }
12343447b6efSHong Zhang       if (!usecprow) xb += 2;
12352d61bbb3SSatish Balay     }
12362d61bbb3SSatish Balay     break;
12372d61bbb3SSatish Balay   case 3:
12382d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
12397b2bb3b9SHong Zhang       if (usecprow) xb = x + 3*ridx[i];
1240f1af5d2fSBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
12413447b6efSHong Zhang       ib = idx + ii[0];
12423447b6efSHong Zhang       n  = ii[1] - ii[0]; ii++;
12432d61bbb3SSatish Balay       for (j=0; j<n; j++) {
12442d61bbb3SSatish Balay         rval      = ib[j]*3;
12452d61bbb3SSatish Balay         z[rval++] += v[0]*x1 + v[1]*x2 + v[2]*x3;
12462d61bbb3SSatish Balay         z[rval++] += v[3]*x1 + v[4]*x2 + v[5]*x3;
12472d61bbb3SSatish Balay         z[rval++] += v[6]*x1 + v[7]*x2 + v[8]*x3;
12482d61bbb3SSatish Balay         v  += 9;
12492d61bbb3SSatish Balay       }
12503447b6efSHong Zhang       if (!usecprow) xb += 3;
12512d61bbb3SSatish Balay     }
12522d61bbb3SSatish Balay     break;
12532d61bbb3SSatish Balay   case 4:
12542d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
12557b2bb3b9SHong Zhang       if (usecprow) xb = x + 4*ridx[i];
1256f1af5d2fSBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
12573447b6efSHong Zhang       ib = idx + ii[0];
12583447b6efSHong Zhang       n  = ii[1] - ii[0]; ii++;
12592d61bbb3SSatish Balay       for (j=0; j<n; j++) {
12602d61bbb3SSatish Balay         rval      = ib[j]*4;
12612d61bbb3SSatish Balay         z[rval++] +=  v[0]*x1 +  v[1]*x2 +  v[2]*x3 +  v[3]*x4;
12622d61bbb3SSatish Balay         z[rval++] +=  v[4]*x1 +  v[5]*x2 +  v[6]*x3 +  v[7]*x4;
12632d61bbb3SSatish Balay         z[rval++] +=  v[8]*x1 +  v[9]*x2 + v[10]*x3 + v[11]*x4;
12642d61bbb3SSatish Balay         z[rval++] += v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4;
12652d61bbb3SSatish Balay         v  += 16;
12662d61bbb3SSatish Balay       }
12673447b6efSHong Zhang       if (!usecprow) xb += 4;
12682d61bbb3SSatish Balay     }
12692d61bbb3SSatish Balay     break;
12702d61bbb3SSatish Balay   case 5:
12712d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
12727b2bb3b9SHong Zhang       if (usecprow) xb = x + 5*ridx[i];
1273f1af5d2fSBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
12742d61bbb3SSatish Balay       x4 = xb[3]; x5 = xb[4];
12753447b6efSHong Zhang       ib = idx + ii[0];
12763447b6efSHong Zhang       n  = ii[1] - ii[0]; ii++;
12772d61bbb3SSatish Balay       for (j=0; j<n; j++) {
12782d61bbb3SSatish Balay         rval      = ib[j]*5;
12792d61bbb3SSatish Balay         z[rval++] +=  v[0]*x1 +  v[1]*x2 +  v[2]*x3 +  v[3]*x4 +  v[4]*x5;
12802d61bbb3SSatish Balay         z[rval++] +=  v[5]*x1 +  v[6]*x2 +  v[7]*x3 +  v[8]*x4 +  v[9]*x5;
12812d61bbb3SSatish Balay         z[rval++] += v[10]*x1 + v[11]*x2 + v[12]*x3 + v[13]*x4 + v[14]*x5;
12822d61bbb3SSatish Balay         z[rval++] += v[15]*x1 + v[16]*x2 + v[17]*x3 + v[18]*x4 + v[19]*x5;
12832d61bbb3SSatish Balay         z[rval++] += v[20]*x1 + v[21]*x2 + v[22]*x3 + v[23]*x4 + v[24]*x5;
12842d61bbb3SSatish Balay         v  += 25;
12852d61bbb3SSatish Balay       }
12863447b6efSHong Zhang       if (!usecprow) xb += 5;
12872d61bbb3SSatish Balay     }
12882d61bbb3SSatish Balay     break;
1289f1af5d2fSBarry Smith   default: {      /* block sizes larger then 5 by 5 are handled by BLAS */
1290690b6cddSBarry Smith       PetscInt     ncols,k;
12913447b6efSHong Zhang       PetscScalar  *work,*workt,*xtmp;
12923f1db9ecSBarry Smith 
12932d61bbb3SSatish Balay       if (!a->mult_work) {
1294d0f46423SBarry Smith         k = PetscMax(A->rmap->n,A->cmap->n);
129587828ca2SBarry Smith         ierr = PetscMalloc((k+1)*sizeof(PetscScalar),&a->mult_work);CHKERRQ(ierr);
12962d61bbb3SSatish Balay       }
12972d61bbb3SSatish Balay       work = a->mult_work;
12983447b6efSHong Zhang       xtmp = x;
12992d61bbb3SSatish Balay       for (i=0; i<mbs; i++) {
13002d61bbb3SSatish Balay         n     = ii[1] - ii[0]; ii++;
13012d61bbb3SSatish Balay         ncols = n*bs;
130287828ca2SBarry Smith         ierr  = PetscMemzero(work,ncols*sizeof(PetscScalar));CHKERRQ(ierr);
13033447b6efSHong Zhang         if (usecprow) {
13047b2bb3b9SHong Zhang           xtmp = x + bs*ridx[i];
13053447b6efSHong Zhang         }
13063447b6efSHong Zhang         Kernel_w_gets_w_plus_trans_Ar_times_v(bs,ncols,xtmp,v,work);
130771044d3cSBarry Smith         /* BLASgemv_("T",&bs,&ncols,&_DOne,v,&bs,xtmp,&_One,&_DOne,work,&_One); */
13082d61bbb3SSatish Balay         v += n*bs2;
13093447b6efSHong Zhang         if (!usecprow) xtmp += bs;
13102d61bbb3SSatish Balay         workt = work;
13112d61bbb3SSatish Balay         for (j=0; j<n; j++) {
13122d61bbb3SSatish Balay           zb = z + bs*(*idx++);
13132d61bbb3SSatish Balay           for (k=0; k<bs; k++) zb[k] += workt[k] ;
13142d61bbb3SSatish Balay           workt += bs;
13152d61bbb3SSatish Balay         }
13162d61bbb3SSatish Balay       }
13172d61bbb3SSatish Balay     }
13182d61bbb3SSatish Balay   }
13193447b6efSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
13203447b6efSHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
1321dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*a->nz*a->bs2);CHKERRQ(ierr);
13222d61bbb3SSatish Balay   PetscFunctionReturn(0);
13232d61bbb3SSatish Balay }
13242d61bbb3SSatish Balay 
13254a2ae208SSatish Balay #undef __FUNCT__
13264a2ae208SSatish Balay #define __FUNCT__ "MatScale_SeqBAIJ"
1327f4df32b1SMatthew Knepley PetscErrorCode MatScale_SeqBAIJ(Mat inA,PetscScalar alpha)
13282d61bbb3SSatish Balay {
13292d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)inA->data;
1330690b6cddSBarry Smith   PetscInt       totalnz = a->bs2*a->nz;
1331f4df32b1SMatthew Knepley   PetscScalar    oalpha = alpha;
1332efee365bSSatish Balay   PetscErrorCode ierr;
13330805154bSBarry Smith   PetscBLASInt   one = 1,tnz = PetscBLASIntCast(totalnz);
13342d61bbb3SSatish Balay 
13352d61bbb3SSatish Balay   PetscFunctionBegin;
1336f4df32b1SMatthew Knepley   BLASscal_(&tnz,&oalpha,a->a,&one);
1337efee365bSSatish Balay   ierr = PetscLogFlops(totalnz);CHKERRQ(ierr);
13382d61bbb3SSatish Balay   PetscFunctionReturn(0);
13392d61bbb3SSatish Balay }
13402d61bbb3SSatish Balay 
13414a2ae208SSatish Balay #undef __FUNCT__
13424a2ae208SSatish Balay #define __FUNCT__ "MatNorm_SeqBAIJ"
1343dfbe8321SBarry Smith PetscErrorCode MatNorm_SeqBAIJ(Mat A,NormType type,PetscReal *norm)
13442d61bbb3SSatish Balay {
13458a62d963SHong Zhang   PetscErrorCode ierr;
13462d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
13473f1db9ecSBarry Smith   MatScalar      *v = a->a;
1348329f5518SBarry Smith   PetscReal      sum = 0.0;
1349d0f46423SBarry Smith   PetscInt       i,j,k,bs=A->rmap->bs,nz=a->nz,bs2=a->bs2,k1;
13502d61bbb3SSatish Balay 
13512d61bbb3SSatish Balay   PetscFunctionBegin;
13522d61bbb3SSatish Balay   if (type == NORM_FROBENIUS) {
13532d61bbb3SSatish Balay     for (i=0; i< bs2*nz; i++) {
1354aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
1355329f5518SBarry Smith       sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
13562d61bbb3SSatish Balay #else
13572d61bbb3SSatish Balay       sum += (*v)*(*v); v++;
13582d61bbb3SSatish Balay #endif
13592d61bbb3SSatish Balay     }
13602d61bbb3SSatish Balay     *norm = sqrt(sum);
13618a62d963SHong Zhang   } else if (type == NORM_1) { /* maximum column sum */
13628a62d963SHong Zhang     PetscReal *tmp;
13638a62d963SHong Zhang     PetscInt  *bcol = a->j;
1364d0f46423SBarry Smith     ierr = PetscMalloc((A->cmap->n+1)*sizeof(PetscReal),&tmp);CHKERRQ(ierr);
1365d0f46423SBarry Smith     ierr = PetscMemzero(tmp,A->cmap->n*sizeof(PetscReal));CHKERRQ(ierr);
13668a62d963SHong Zhang     for (i=0; i<nz; i++){
13678a62d963SHong Zhang       for (j=0; j<bs; j++){
13688a62d963SHong Zhang         k1 = bs*(*bcol) + j; /* column index */
13698a62d963SHong Zhang         for (k=0; k<bs; k++){
13708a62d963SHong Zhang           tmp[k1] += PetscAbsScalar(*v); v++;
13718a62d963SHong Zhang         }
13728a62d963SHong Zhang       }
13738a62d963SHong Zhang       bcol++;
13748a62d963SHong Zhang     }
13758a62d963SHong Zhang     *norm = 0.0;
1376d0f46423SBarry Smith     for (j=0; j<A->cmap->n; j++) {
13778a62d963SHong Zhang       if (tmp[j] > *norm) *norm = tmp[j];
13788a62d963SHong Zhang     }
13798a62d963SHong Zhang     ierr = PetscFree(tmp);CHKERRQ(ierr);
1380596552b5SBarry Smith   } else if (type == NORM_INFINITY) { /* maximum row sum */
1381596552b5SBarry Smith     *norm = 0.0;
1382596552b5SBarry Smith     for (k=0; k<bs; k++) {
138374f84c7bSSatish Balay       for (j=0; j<a->mbs; j++) {
1384596552b5SBarry Smith         v = a->a + bs2*a->i[j] + k;
1385596552b5SBarry Smith         sum = 0.0;
1386596552b5SBarry Smith         for (i=0; i<a->i[j+1]-a->i[j]; i++) {
13870e90e235SBarry Smith           for (k1=0; k1<bs; k1++){
1388596552b5SBarry Smith             sum += PetscAbsScalar(*v);
1389596552b5SBarry Smith             v   += bs;
13902d61bbb3SSatish Balay           }
13910e90e235SBarry Smith         }
1392596552b5SBarry Smith         if (sum > *norm) *norm = sum;
1393596552b5SBarry Smith       }
1394596552b5SBarry Smith     }
1395596552b5SBarry Smith   } else {
139629bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"No support for this norm yet");
13972d61bbb3SSatish Balay   }
13982d61bbb3SSatish Balay   PetscFunctionReturn(0);
13992d61bbb3SSatish Balay }
14002d61bbb3SSatish Balay 
14012d61bbb3SSatish Balay 
14024a2ae208SSatish Balay #undef __FUNCT__
14034a2ae208SSatish Balay #define __FUNCT__ "MatEqual_SeqBAIJ"
1404dfbe8321SBarry Smith PetscErrorCode MatEqual_SeqBAIJ(Mat A,Mat B,PetscTruth* flg)
14052d61bbb3SSatish Balay {
14062d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ *)A->data,*b = (Mat_SeqBAIJ *)B->data;
1407dfbe8321SBarry Smith   PetscErrorCode ierr;
14082d61bbb3SSatish Balay 
14092d61bbb3SSatish Balay   PetscFunctionBegin;
14102d61bbb3SSatish Balay   /* If the  matrix/block dimensions are not equal, or no of nonzeros or shift */
1411d0f46423SBarry Smith   if ((A->rmap->N != B->rmap->N) || (A->cmap->n != B->cmap->n) || (A->rmap->bs != B->rmap->bs)|| (a->nz != b->nz)) {
1412273d9f13SBarry Smith     *flg = PETSC_FALSE;
1413273d9f13SBarry Smith     PetscFunctionReturn(0);
14142d61bbb3SSatish Balay   }
14152d61bbb3SSatish Balay 
14162d61bbb3SSatish Balay   /* if the a->i are the same */
1417690b6cddSBarry Smith   ierr = PetscMemcmp(a->i,b->i,(a->mbs+1)*sizeof(PetscInt),flg);CHKERRQ(ierr);
1418abc0a331SBarry Smith   if (!*flg) {
14190f5bd95cSBarry Smith     PetscFunctionReturn(0);
14202d61bbb3SSatish Balay   }
14212d61bbb3SSatish Balay 
14222d61bbb3SSatish Balay   /* if a->j are the same */
1423690b6cddSBarry Smith   ierr = PetscMemcmp(a->j,b->j,(a->nz)*sizeof(PetscInt),flg);CHKERRQ(ierr);
1424abc0a331SBarry Smith   if (!*flg) {
14250f5bd95cSBarry Smith     PetscFunctionReturn(0);
14262d61bbb3SSatish Balay   }
14272d61bbb3SSatish Balay   /* if a->a are the same */
1428d0f46423SBarry Smith   ierr = PetscMemcmp(a->a,b->a,(a->nz)*(A->rmap->bs)*(B->rmap->bs)*sizeof(PetscScalar),flg);CHKERRQ(ierr);
14292d61bbb3SSatish Balay   PetscFunctionReturn(0);
14302d61bbb3SSatish Balay 
14312d61bbb3SSatish Balay }
14322d61bbb3SSatish Balay 
14334a2ae208SSatish Balay #undef __FUNCT__
14344a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_SeqBAIJ"
1435dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_SeqBAIJ(Mat A,Vec v)
14362d61bbb3SSatish Balay {
14372d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
1438dfbe8321SBarry Smith   PetscErrorCode ierr;
1439690b6cddSBarry Smith   PetscInt       i,j,k,n,row,bs,*ai,*aj,ambs,bs2;
144087828ca2SBarry Smith   PetscScalar    *x,zero = 0.0;
14413f1db9ecSBarry Smith   MatScalar      *aa,*aa_j;
14422d61bbb3SSatish Balay 
14432d61bbb3SSatish Balay   PetscFunctionBegin;
144429bbc08cSBarry Smith   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1445d0f46423SBarry Smith   bs   = A->rmap->bs;
14462d61bbb3SSatish Balay   aa   = a->a;
14472d61bbb3SSatish Balay   ai   = a->i;
14482d61bbb3SSatish Balay   aj   = a->j;
14492d61bbb3SSatish Balay   ambs = a->mbs;
14502d61bbb3SSatish Balay   bs2  = a->bs2;
14512d61bbb3SSatish Balay 
14522dcb1b2aSMatthew Knepley   ierr = VecSet(v,zero);CHKERRQ(ierr);
14531ebc52fbSHong Zhang   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1454e1311b90SBarry Smith   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
1455d0f46423SBarry Smith   if (n != A->rmap->N) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
14562d61bbb3SSatish Balay   for (i=0; i<ambs; i++) {
14572d61bbb3SSatish Balay     for (j=ai[i]; j<ai[i+1]; j++) {
14582d61bbb3SSatish Balay       if (aj[j] == i) {
14592d61bbb3SSatish Balay         row  = i*bs;
14602d61bbb3SSatish Balay         aa_j = aa+j*bs2;
14612d61bbb3SSatish Balay         for (k=0; k<bs2; k+=(bs+1),row++) x[row] = aa_j[k];
14622d61bbb3SSatish Balay         break;
14632d61bbb3SSatish Balay       }
14642d61bbb3SSatish Balay     }
14652d61bbb3SSatish Balay   }
14661ebc52fbSHong Zhang   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
14672d61bbb3SSatish Balay   PetscFunctionReturn(0);
14682d61bbb3SSatish Balay }
14692d61bbb3SSatish Balay 
14704a2ae208SSatish Balay #undef __FUNCT__
14714a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_SeqBAIJ"
1472dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_SeqBAIJ(Mat A,Vec ll,Vec rr)
14732d61bbb3SSatish Balay {
14742d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
147587828ca2SBarry Smith   PetscScalar    *l,*r,x,*li,*ri;
14763f1db9ecSBarry Smith   MatScalar      *aa,*v;
1477dfbe8321SBarry Smith   PetscErrorCode ierr;
1478690b6cddSBarry Smith   PetscInt       i,j,k,lm,rn,M,m,n,*ai,*aj,mbs,tmp,bs,bs2;
14792d61bbb3SSatish Balay 
14802d61bbb3SSatish Balay   PetscFunctionBegin;
14812d61bbb3SSatish Balay   ai  = a->i;
14822d61bbb3SSatish Balay   aj  = a->j;
14832d61bbb3SSatish Balay   aa  = a->a;
1484d0f46423SBarry Smith   m   = A->rmap->n;
1485d0f46423SBarry Smith   n   = A->cmap->n;
1486d0f46423SBarry Smith   bs  = A->rmap->bs;
14872d61bbb3SSatish Balay   mbs = a->mbs;
14882d61bbb3SSatish Balay   bs2 = a->bs2;
14892d61bbb3SSatish Balay   if (ll) {
14901ebc52fbSHong Zhang     ierr = VecGetArray(ll,&l);CHKERRQ(ierr);
1491e1311b90SBarry Smith     ierr = VecGetLocalSize(ll,&lm);CHKERRQ(ierr);
149229bbc08cSBarry Smith     if (lm != m) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vector wrong length");
14932d61bbb3SSatish Balay     for (i=0; i<mbs; i++) { /* for each block row */
14942d61bbb3SSatish Balay       M  = ai[i+1] - ai[i];
14952d61bbb3SSatish Balay       li = l + i*bs;
14962d61bbb3SSatish Balay       v  = aa + bs2*ai[i];
14972d61bbb3SSatish Balay       for (j=0; j<M; j++) { /* for each block */
14982d61bbb3SSatish Balay         for (k=0; k<bs2; k++) {
14992d61bbb3SSatish Balay           (*v++) *= li[k%bs];
15002d61bbb3SSatish Balay         }
15012d61bbb3SSatish Balay       }
15022d61bbb3SSatish Balay     }
15031ebc52fbSHong Zhang     ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr);
1504efee365bSSatish Balay     ierr = PetscLogFlops(a->nz);CHKERRQ(ierr);
15052d61bbb3SSatish Balay   }
15062d61bbb3SSatish Balay 
15072d61bbb3SSatish Balay   if (rr) {
15081ebc52fbSHong Zhang     ierr = VecGetArray(rr,&r);CHKERRQ(ierr);
1509e1311b90SBarry Smith     ierr = VecGetLocalSize(rr,&rn);CHKERRQ(ierr);
151029bbc08cSBarry Smith     if (rn != n) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vector wrong length");
15112d61bbb3SSatish Balay     for (i=0; i<mbs; i++) { /* for each block row */
15122d61bbb3SSatish Balay       M  = ai[i+1] - ai[i];
15132d61bbb3SSatish Balay       v  = aa + bs2*ai[i];
15142d61bbb3SSatish Balay       for (j=0; j<M; j++) { /* for each block */
15152d61bbb3SSatish Balay         ri = r + bs*aj[ai[i]+j];
15162d61bbb3SSatish Balay         for (k=0; k<bs; k++) {
15172d61bbb3SSatish Balay           x = ri[k];
15182d61bbb3SSatish Balay           for (tmp=0; tmp<bs; tmp++) (*v++) *= x;
15192d61bbb3SSatish Balay         }
15202d61bbb3SSatish Balay       }
15212d61bbb3SSatish Balay     }
15221ebc52fbSHong Zhang     ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr);
1523efee365bSSatish Balay     ierr = PetscLogFlops(a->nz);CHKERRQ(ierr);
15242d61bbb3SSatish Balay   }
15252d61bbb3SSatish Balay   PetscFunctionReturn(0);
15262d61bbb3SSatish Balay }
15272d61bbb3SSatish Balay 
15282d61bbb3SSatish Balay 
15294a2ae208SSatish Balay #undef __FUNCT__
15304a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_SeqBAIJ"
1531dfbe8321SBarry Smith PetscErrorCode MatGetInfo_SeqBAIJ(Mat A,MatInfoType flag,MatInfo *info)
15322d61bbb3SSatish Balay {
15332d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
15342d61bbb3SSatish Balay 
15352d61bbb3SSatish Balay   PetscFunctionBegin;
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;
15427adad957SLisandro 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