xref: /petsc/src/mat/impls/baij/seq/baij2.c (revision 5d0c19d75c660d4fec594a5399ec8d8ba29c54a8)
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;
14*5d0c19d7SBarry Smith   PetscInt       row,i,j,k,l,m,n,*nidx,isz,val,ival;
15*5d0c19d7SBarry Smith   const PetscInt *idx;
16690b6cddSBarry Smith   PetscInt       start,end,*ai,*aj,bs,*nidx2;
17f1af5d2fSBarry Smith   PetscBT        table;
18a3192f15SSatish Balay 
193a40ed3dSBarry Smith   PetscFunctionBegin;
20a3192f15SSatish Balay   m     = a->mbs;
21a3192f15SSatish Balay   ai    = a->i;
22a3192f15SSatish Balay   aj    = a->j;
23d0f46423SBarry Smith   bs    = A->rmap->bs;
24a3192f15SSatish Balay 
2529bbc08cSBarry Smith   if (ov < 0)  SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap specified");
26a3192f15SSatish Balay 
276831982aSBarry Smith   ierr = PetscBTCreate(m,table);CHKERRQ(ierr);
28690b6cddSBarry Smith   ierr = PetscMalloc((m+1)*sizeof(PetscInt),&nidx);CHKERRQ(ierr);
29d0f46423SBarry Smith   ierr = PetscMalloc((A->rmap->N+1)*sizeof(PetscInt),&nidx2);CHKERRQ(ierr);
30a3192f15SSatish Balay 
31a3192f15SSatish Balay   for (i=0; i<is_max; i++) {
32a3192f15SSatish Balay     /* Initialise the two local arrays */
33a3192f15SSatish Balay     isz  = 0;
346831982aSBarry Smith     ierr = PetscBTMemzero(m,table);CHKERRQ(ierr);
35a3192f15SSatish Balay 
36a3192f15SSatish Balay     /* Extract the indices, assume there can be duplicate entries */
37a3192f15SSatish Balay     ierr = ISGetIndices(is[i],&idx);CHKERRQ(ierr);
38b9b97703SBarry Smith     ierr = ISGetLocalSize(is[i],&n);CHKERRQ(ierr);
39a3192f15SSatish Balay 
40a3192f15SSatish Balay     /* Enter these into the temp arrays i.e mark table[row], enter row into new index */
41a3192f15SSatish Balay     for (j=0; j<n ; ++j){
42218c64b6SSatish Balay       ival = idx[j]/bs; /* convert the indices into block indices */
435eee224dSHong Zhang       if (ival>=m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"index greater than mat-dim");
44f1af5d2fSBarry Smith       if(!PetscBTLookupSet(table,ival)) { nidx[isz++] = ival;}
45a3192f15SSatish Balay     }
46a3192f15SSatish Balay     ierr = ISRestoreIndices(is[i],&idx);CHKERRQ(ierr);
47a3192f15SSatish Balay     ierr = ISDestroy(is[i]);CHKERRQ(ierr);
48a3192f15SSatish Balay 
49a3192f15SSatish Balay     k = 0;
50a3192f15SSatish Balay     for (j=0; j<ov; j++){ /* for each overlap*/
51a3192f15SSatish Balay       n = isz;
52a3192f15SSatish Balay       for (; k<n ; k++){ /* do only those rows in nidx[k], which are not done yet */
53a3192f15SSatish Balay         row   = nidx[k];
54a3192f15SSatish Balay         start = ai[row];
55a3192f15SSatish Balay         end   = ai[row+1];
56a3192f15SSatish Balay         for (l = start; l<end ; l++){
57a3192f15SSatish Balay           val = aj[l];
58f1af5d2fSBarry Smith           if (!PetscBTLookupSet(table,val)) {nidx[isz++] = val;}
59a3192f15SSatish Balay         }
60a3192f15SSatish Balay       }
61a3192f15SSatish Balay     }
62218c64b6SSatish Balay     /* expand the Index Set */
63218c64b6SSatish Balay     for (j=0; j<isz; j++) {
64218c64b6SSatish Balay       for (k=0; k<bs; k++)
65218c64b6SSatish Balay         nidx2[j*bs+k] = nidx[j]*bs+k;
66218c64b6SSatish Balay     }
67329f5518SBarry Smith     ierr = ISCreateGeneral(PETSC_COMM_SELF,isz*bs,nidx2,is+i);CHKERRQ(ierr);
68a3192f15SSatish Balay   }
696831982aSBarry Smith   ierr = PetscBTDestroy(table);CHKERRQ(ierr);
70606d414cSSatish Balay   ierr = PetscFree(nidx);CHKERRQ(ierr);
71606d414cSSatish Balay   ierr = PetscFree(nidx2);CHKERRQ(ierr);
723a40ed3dSBarry Smith   PetscFunctionReturn(0);
73a3192f15SSatish Balay }
741c351548SSatish Balay 
754a2ae208SSatish Balay #undef __FUNCT__
764a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqBAIJ_Private"
77690b6cddSBarry Smith PetscErrorCode MatGetSubMatrix_SeqBAIJ_Private(Mat A,IS isrow,IS iscol,PetscInt cs,MatReuse scall,Mat *B)
78736121d4SSatish Balay {
79736121d4SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*c;
806849ba73SBarry Smith   PetscErrorCode ierr;
81690b6cddSBarry Smith   PetscInt       *smap,i,k,kstart,kend,oldcols = a->nbs,*lens;
82690b6cddSBarry Smith   PetscInt       row,mat_i,*mat_j,tcol,*mat_ilen;
83*5d0c19d7SBarry Smith   const PetscInt *irow,*icol;
84*5d0c19d7SBarry Smith   PetscInt       nrows,ncols,*ssmap,bs=A->rmap->bs,bs2=a->bs2;
85690b6cddSBarry Smith   PetscInt       *aj = a->j,*ai = a->i;
863f1db9ecSBarry Smith   MatScalar      *mat_a;
87736121d4SSatish Balay   Mat            C;
880f5bd95cSBarry Smith   PetscTruth     flag;
89736121d4SSatish Balay 
903a40ed3dSBarry Smith   PetscFunctionBegin;
91888f2ed8SSatish Balay   ierr = ISSorted(iscol,(PetscTruth*)&i);CHKERRQ(ierr);
9229bbc08cSBarry Smith   if (!i) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"IS is not sorted");
93736121d4SSatish Balay 
94736121d4SSatish Balay   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
95218c64b6SSatish Balay   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
96b9b97703SBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
97b9b97703SBarry Smith   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
98736121d4SSatish Balay 
99690b6cddSBarry Smith   ierr = PetscMalloc((1+oldcols)*sizeof(PetscInt),&smap);CHKERRQ(ierr);
100736121d4SSatish Balay   ssmap = smap;
101690b6cddSBarry Smith   ierr = PetscMalloc((1+nrows)*sizeof(PetscInt),&lens);CHKERRQ(ierr);
102690b6cddSBarry Smith   ierr  = PetscMemzero(smap,oldcols*sizeof(PetscInt));CHKERRQ(ierr);
103736121d4SSatish Balay   for (i=0; i<ncols; i++) smap[icol[i]] = i+1;
104736121d4SSatish Balay   /* determine lens of each row */
105736121d4SSatish Balay   for (i=0; i<nrows; i++) {
106736121d4SSatish Balay     kstart  = ai[irow[i]];
107736121d4SSatish Balay     kend    = kstart + a->ilen[irow[i]];
108736121d4SSatish Balay     lens[i] = 0;
109736121d4SSatish Balay       for (k=kstart; k<kend; k++) {
110736121d4SSatish Balay         if (ssmap[aj[k]]) {
111736121d4SSatish Balay           lens[i]++;
112736121d4SSatish Balay         }
113736121d4SSatish Balay       }
114736121d4SSatish Balay     }
115736121d4SSatish Balay   /* Create and fill new matrix */
116736121d4SSatish Balay   if (scall == MAT_REUSE_MATRIX) {
117736121d4SSatish Balay     c = (Mat_SeqBAIJ *)((*B)->data);
118736121d4SSatish Balay 
119d0f46423SBarry Smith     if (c->mbs!=nrows || c->nbs!=ncols || (*B)->rmap->bs!=bs) SETERRQ(PETSC_ERR_ARG_SIZ,"Submatrix wrong size");
120690b6cddSBarry Smith     ierr = PetscMemcmp(c->ilen,lens,c->mbs *sizeof(PetscInt),&flag);CHKERRQ(ierr);
121abc0a331SBarry Smith     if (!flag) {
12229bbc08cSBarry Smith       SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong no of nonzeros");
123736121d4SSatish Balay     }
124690b6cddSBarry Smith     ierr = PetscMemzero(c->ilen,c->mbs*sizeof(PetscInt));CHKERRQ(ierr);
125736121d4SSatish Balay     C = *B;
1263a40ed3dSBarry Smith   } else {
1277adad957SLisandro Dalcin     ierr = MatCreate(((PetscObject)A)->comm,&C);CHKERRQ(ierr);
128f69a0ea3SMatthew Knepley     ierr = MatSetSizes(C,nrows*bs,ncols*bs,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1297adad957SLisandro Dalcin     ierr = MatSetType(C,((PetscObject)A)->type_name);CHKERRQ(ierr);
130ab93d7beSBarry Smith     ierr = MatSeqBAIJSetPreallocation_SeqBAIJ(C,bs,0,lens);CHKERRQ(ierr);
131736121d4SSatish Balay   }
132736121d4SSatish Balay   c = (Mat_SeqBAIJ *)(C->data);
133736121d4SSatish Balay   for (i=0; i<nrows; i++) {
134736121d4SSatish Balay     row    = irow[i];
135736121d4SSatish Balay     kstart = ai[row];
136736121d4SSatish Balay     kend   = kstart + a->ilen[row];
137736121d4SSatish Balay     mat_i  = c->i[i];
138736121d4SSatish Balay     mat_j  = c->j + mat_i;
139218c64b6SSatish Balay     mat_a  = c->a + mat_i*bs2;
140736121d4SSatish Balay     mat_ilen = c->ilen + i;
141736121d4SSatish Balay     for (k=kstart; k<kend; k++) {
142736121d4SSatish Balay       if ((tcol=ssmap[a->j[k]])) {
143736121d4SSatish Balay         *mat_j++ = tcol - 1;
144549d3d68SSatish Balay         ierr     = PetscMemcpy(mat_a,a->a+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
145549d3d68SSatish Balay         mat_a   += bs2;
146736121d4SSatish Balay         (*mat_ilen)++;
147736121d4SSatish Balay       }
148736121d4SSatish Balay     }
149736121d4SSatish Balay   }
150218c64b6SSatish Balay 
151736121d4SSatish Balay   /* Free work space */
152736121d4SSatish Balay   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
153606d414cSSatish Balay   ierr = PetscFree(smap);CHKERRQ(ierr);
154606d414cSSatish Balay   ierr = PetscFree(lens);CHKERRQ(ierr);
1556d4a8577SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1566d4a8577SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
157736121d4SSatish Balay 
158736121d4SSatish Balay   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
159736121d4SSatish Balay   *B = C;
1603a40ed3dSBarry Smith   PetscFunctionReturn(0);
161736121d4SSatish Balay }
162736121d4SSatish Balay 
1634a2ae208SSatish Balay #undef __FUNCT__
1644a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqBAIJ"
165690b6cddSBarry Smith PetscErrorCode MatGetSubMatrix_SeqBAIJ(Mat A,IS isrow,IS iscol,PetscInt cs,MatReuse scall,Mat *B)
166218c64b6SSatish Balay {
167218c64b6SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
168218c64b6SSatish Balay   IS             is1,is2;
1696849ba73SBarry Smith   PetscErrorCode ierr;
170*5d0c19d7SBarry Smith   PetscInt       *vary,*iary,nrows,ncols,i,bs=A->rmap->bs,count;
171*5d0c19d7SBarry Smith   const PetscInt *irow,*icol;
172218c64b6SSatish Balay 
1733a40ed3dSBarry Smith   PetscFunctionBegin;
174218c64b6SSatish Balay   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
175218c64b6SSatish Balay   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
176b9b97703SBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
177b9b97703SBarry Smith   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
178218c64b6SSatish Balay 
179218c64b6SSatish Balay   /* Verify if the indices corespond to each element in a block
180218c64b6SSatish Balay    and form the IS with compressed IS */
181690b6cddSBarry Smith   ierr = PetscMalloc(2*(a->mbs+1)*sizeof(PetscInt),&vary);CHKERRQ(ierr);
182218c64b6SSatish Balay   iary = vary + a->mbs;
183690b6cddSBarry Smith   ierr = PetscMemzero(vary,(a->mbs)*sizeof(PetscInt));CHKERRQ(ierr);
184218c64b6SSatish Balay   for (i=0; i<nrows; i++) vary[irow[i]/bs]++;
185218c64b6SSatish Balay   count = 0;
186218c64b6SSatish Balay   for (i=0; i<a->mbs; i++) {
187634064b4SBarry Smith     if (vary[i]!=0 && vary[i]!=bs) SETERRQ(PETSC_ERR_ARG_SIZ,"Index set does not match blocks");
188218c64b6SSatish Balay     if (vary[i]==bs) iary[count++] = i;
189218c64b6SSatish Balay   }
190029af93fSBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,count,iary,&is1);CHKERRQ(ierr);
191218c64b6SSatish Balay 
192690b6cddSBarry Smith   ierr = PetscMemzero(vary,(a->mbs)*sizeof(PetscInt));CHKERRQ(ierr);
193218c64b6SSatish Balay   for (i=0; i<ncols; i++) vary[icol[i]/bs]++;
194218c64b6SSatish Balay   count = 0;
195218c64b6SSatish Balay   for (i=0; i<a->mbs; i++) {
196634064b4SBarry Smith     if (vary[i]!=0 && vary[i]!=bs) SETERRQ(PETSC_ERR_PLIB,"Internal error in PETSc");
197218c64b6SSatish Balay     if (vary[i]==bs) iary[count++] = i;
198218c64b6SSatish Balay   }
199029af93fSBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,count,iary,&is2);CHKERRQ(ierr);
200218c64b6SSatish Balay   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
201218c64b6SSatish Balay   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
202606d414cSSatish Balay   ierr = PetscFree(vary);CHKERRQ(ierr);
203218c64b6SSatish Balay 
2046a6a5d1dSBarry Smith   ierr = MatGetSubMatrix_SeqBAIJ_Private(A,is1,is2,cs,scall,B);CHKERRQ(ierr);
205218c64b6SSatish Balay   ISDestroy(is1);
206218c64b6SSatish Balay   ISDestroy(is2);
2073a40ed3dSBarry Smith   PetscFunctionReturn(0);
208218c64b6SSatish Balay }
209218c64b6SSatish Balay 
2104a2ae208SSatish Balay #undef __FUNCT__
2114a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_SeqBAIJ"
212690b6cddSBarry Smith PetscErrorCode MatGetSubMatrices_SeqBAIJ(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
213736121d4SSatish Balay {
2146849ba73SBarry Smith   PetscErrorCode ierr;
215690b6cddSBarry Smith   PetscInt       i;
216736121d4SSatish Balay 
2173a40ed3dSBarry Smith   PetscFunctionBegin;
218736121d4SSatish Balay   if (scall == MAT_INITIAL_MATRIX) {
21982502324SSatish Balay     ierr = PetscMalloc((n+1)*sizeof(Mat),B);CHKERRQ(ierr);
220736121d4SSatish Balay   }
221736121d4SSatish Balay 
222736121d4SSatish Balay   for (i=0; i<n; i++) {
2236a6a5d1dSBarry Smith     ierr = MatGetSubMatrix_SeqBAIJ(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr);
224736121d4SSatish Balay   }
2253a40ed3dSBarry Smith   PetscFunctionReturn(0);
226736121d4SSatish Balay }
227218c64b6SSatish Balay 
228218c64b6SSatish Balay 
2292d61bbb3SSatish Balay /* -------------------------------------------------------*/
2302d61bbb3SSatish Balay /* Should check that shapes of vectors and matrices match */
2312d61bbb3SSatish Balay /* -------------------------------------------------------*/
232d9eff348SSatish Balay #include "petscblaslapack.h"
2332d61bbb3SSatish Balay 
2344a2ae208SSatish Balay #undef __FUNCT__
2354a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_1"
236dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_1(Mat A,Vec xx,Vec zz)
2372d61bbb3SSatish Balay {
2382d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
239d9fead3dSBarry Smith   PetscScalar       *z,sum;
240d9fead3dSBarry Smith   const PetscScalar *x;
241d9fead3dSBarry Smith   const MatScalar   *v;
2426849ba73SBarry Smith   PetscErrorCode    ierr;
24398c9bda7SSatish Balay   PetscInt          mbs,i,*idx,*ii,n,*ridx=PETSC_NULL,nonzerorow=0;
24426e093fcSHong Zhang   PetscTruth        usecprow=a->compressedrow.use;
2452d61bbb3SSatish Balay 
2462d61bbb3SSatish Balay   PetscFunctionBegin;
247d9fead3dSBarry Smith   ierr = VecGetArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
2481ebc52fbSHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
2492d61bbb3SSatish Balay 
2502d61bbb3SSatish Balay   idx = a->j;
2512d61bbb3SSatish Balay   v   = a->a;
25226e093fcSHong Zhang   if (usecprow){
25326e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
25426e093fcSHong Zhang     ii   = a->compressedrow.i;
2557b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
25626e093fcSHong Zhang   } else {
25726e093fcSHong Zhang     mbs = a->mbs;
2582d61bbb3SSatish Balay     ii  = a->i;
25926e093fcSHong Zhang   }
2602d61bbb3SSatish Balay 
2612d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
2622d61bbb3SSatish Balay     n    = ii[1] - ii[0]; ii++;
2632d61bbb3SSatish Balay     sum  = 0.0;
26498c9bda7SSatish Balay     nonzerorow += (n>0);
2652d61bbb3SSatish Balay     while (n--) sum += *v++ * x[*idx++];
26626e093fcSHong Zhang     if (usecprow){
2677b2bb3b9SHong Zhang       z[ridx[i]] = sum;
26826e093fcSHong Zhang     } else {
2692d61bbb3SSatish Balay       z[i] = sum;
2702d61bbb3SSatish Balay     }
27126e093fcSHong Zhang   }
272d9fead3dSBarry Smith   ierr = VecRestoreArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
2731ebc52fbSHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
27498c9bda7SSatish Balay   ierr = PetscLogFlops(2*a->nz - nonzerorow);CHKERRQ(ierr);
2752d61bbb3SSatish Balay   PetscFunctionReturn(0);
2762d61bbb3SSatish Balay }
2772d61bbb3SSatish Balay 
2784a2ae208SSatish Balay #undef __FUNCT__
2794a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_2"
280dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_2(Mat A,Vec xx,Vec zz)
2812d61bbb3SSatish Balay {
2822d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
283d9fead3dSBarry Smith   PetscScalar       *z = 0,sum1,sum2,*zarray;
284d9fead3dSBarry Smith   const PetscScalar *x,*xb;
28587828ca2SBarry Smith   PetscScalar       x1,x2;
286d9fead3dSBarry Smith   const MatScalar   *v;
287dfbe8321SBarry Smith   PetscErrorCode    ierr;
28898c9bda7SSatish Balay   PetscInt          mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL,nonzerorow=0;
28926e093fcSHong Zhang   PetscTruth        usecprow=a->compressedrow.use;
2902d61bbb3SSatish Balay 
2912d61bbb3SSatish Balay   PetscFunctionBegin;
292d9fead3dSBarry Smith   ierr = VecGetArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
29326e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
2942d61bbb3SSatish Balay 
2952d61bbb3SSatish Balay   idx = a->j;
2962d61bbb3SSatish Balay   v   = a->a;
29726e093fcSHong Zhang   if (usecprow){
29826e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
29926e093fcSHong Zhang     ii   = a->compressedrow.i;
3007b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
30126e093fcSHong Zhang   } else {
30226e093fcSHong Zhang     mbs = a->mbs;
3032d61bbb3SSatish Balay     ii  = a->i;
30426e093fcSHong Zhang     z   = zarray;
30526e093fcSHong Zhang   }
3062d61bbb3SSatish Balay 
3072d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
3082d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
3092d61bbb3SSatish Balay     sum1 = 0.0; sum2 = 0.0;
31098c9bda7SSatish Balay     nonzerorow += (n>0);
3112d61bbb3SSatish Balay     for (j=0; j<n; j++) {
3122d61bbb3SSatish Balay       xb = x + 2*(*idx++); x1 = xb[0]; x2 = xb[1];
3132d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[2]*x2;
3142d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[3]*x2;
3152d61bbb3SSatish Balay       v += 4;
3162d61bbb3SSatish Balay     }
3177b2bb3b9SHong Zhang     if (usecprow) z = zarray + 2*ridx[i];
3182d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2;
31926e093fcSHong Zhang     if (!usecprow) z += 2;
3202d61bbb3SSatish Balay   }
321d9fead3dSBarry Smith   ierr = VecRestoreArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
32226e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
32398c9bda7SSatish Balay   ierr = PetscLogFlops(8*a->nz - 2*nonzerorow);CHKERRQ(ierr);
3242d61bbb3SSatish Balay   PetscFunctionReturn(0);
3252d61bbb3SSatish Balay }
3262d61bbb3SSatish Balay 
3274a2ae208SSatish Balay #undef __FUNCT__
3284a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_3"
329dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_3(Mat A,Vec xx,Vec zz)
3302d61bbb3SSatish Balay {
3312d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
332d9fead3dSBarry Smith   PetscScalar       *z = 0,sum1,sum2,sum3,x1,x2,x3,*zarray;
333d9fead3dSBarry Smith   const PetscScalar *x,*xb;
334d9fead3dSBarry Smith   const MatScalar   *v;
335dfbe8321SBarry Smith   PetscErrorCode    ierr;
336028cd4eaSSatish Balay   PetscInt          mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL,nonzerorow=0;
337028cd4eaSSatish Balay   PetscTruth        usecprow=a->compressedrow.use;
33826e093fcSHong Zhang 
3392d61bbb3SSatish Balay 
340b6410449SSatish Balay #if defined(PETSC_HAVE_PRAGMA_DISJOINT)
341fee21e36SBarry Smith #pragma disjoint(*v,*z,*xb)
342fee21e36SBarry Smith #endif
343fee21e36SBarry Smith 
3442d61bbb3SSatish Balay   PetscFunctionBegin;
345d9fead3dSBarry Smith   ierr = VecGetArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
34626e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
3472d61bbb3SSatish Balay 
3482d61bbb3SSatish Balay   idx = a->j;
3492d61bbb3SSatish Balay   v   = a->a;
35026e093fcSHong Zhang   if (usecprow){
35126e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
35226e093fcSHong Zhang     ii   = a->compressedrow.i;
3537b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
35426e093fcSHong Zhang   } else {
35526e093fcSHong Zhang     mbs = a->mbs;
3562d61bbb3SSatish Balay     ii  = a->i;
35726e093fcSHong Zhang     z   = zarray;
35826e093fcSHong Zhang   }
3592d61bbb3SSatish Balay 
3602d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
3612d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
3622d61bbb3SSatish Balay     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0;
36398c9bda7SSatish Balay     nonzerorow += (n>0);
3642d61bbb3SSatish Balay     for (j=0; j<n; j++) {
3652d61bbb3SSatish Balay       xb = x + 3*(*idx++); x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
3662d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[3]*x2 + v[6]*x3;
3672d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[4]*x2 + v[7]*x3;
3682d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[5]*x2 + v[8]*x3;
3692d61bbb3SSatish Balay       v += 9;
3702d61bbb3SSatish Balay     }
3717b2bb3b9SHong Zhang     if (usecprow) z = zarray + 3*ridx[i];
3722d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3;
37326e093fcSHong Zhang     if (!usecprow) z += 3;
3742d61bbb3SSatish Balay   }
375d9fead3dSBarry Smith   ierr = VecRestoreArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
37626e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
37798c9bda7SSatish Balay   ierr = PetscLogFlops(18*a->nz - 3*nonzerorow);CHKERRQ(ierr);
3782d61bbb3SSatish Balay   PetscFunctionReturn(0);
3792d61bbb3SSatish Balay }
3802d61bbb3SSatish Balay 
3814a2ae208SSatish Balay #undef __FUNCT__
3824a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_4"
383dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_4(Mat A,Vec xx,Vec zz)
3842d61bbb3SSatish Balay {
3852d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
386d9fead3dSBarry Smith   PetscScalar       *z = 0,sum1,sum2,sum3,sum4,x1,x2,x3,x4,*zarray;
387d9fead3dSBarry Smith   const PetscScalar *x,*xb;
388d9fead3dSBarry Smith   const MatScalar   *v;
389dfbe8321SBarry Smith   PetscErrorCode    ierr;
39098c9bda7SSatish Balay   PetscInt          mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL,nonzerorow=0;
39126e093fcSHong Zhang   PetscTruth        usecprow=a->compressedrow.use;
3922d61bbb3SSatish Balay 
3932d61bbb3SSatish Balay   PetscFunctionBegin;
394d9fead3dSBarry Smith   ierr = VecGetArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
39526e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
3962d61bbb3SSatish Balay 
3972d61bbb3SSatish Balay   idx = a->j;
3982d61bbb3SSatish Balay   v   = a->a;
39926e093fcSHong Zhang   if (usecprow){
40026e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
40126e093fcSHong Zhang     ii   = a->compressedrow.i;
4027b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
40326e093fcSHong Zhang   } else {
40426e093fcSHong Zhang     mbs = a->mbs;
4052d61bbb3SSatish Balay     ii  = a->i;
40626e093fcSHong Zhang     z   = zarray;
40726e093fcSHong Zhang   }
4082d61bbb3SSatish Balay 
4092d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
4102d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
4112d61bbb3SSatish Balay     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0;
41298c9bda7SSatish Balay     nonzerorow += (n>0);
4132d61bbb3SSatish Balay     for (j=0; j<n; j++) {
4142d61bbb3SSatish Balay       xb = x + 4*(*idx++);
4152d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
4162d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[4]*x2 + v[8]*x3   + v[12]*x4;
4172d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[5]*x2 + v[9]*x3   + v[13]*x4;
4182d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[6]*x2 + v[10]*x3  + v[14]*x4;
4192d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[7]*x2 + v[11]*x3  + v[15]*x4;
4202d61bbb3SSatish Balay       v += 16;
4212d61bbb3SSatish Balay     }
4227b2bb3b9SHong Zhang     if (usecprow) z = zarray + 4*ridx[i];
4232d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4;
42426e093fcSHong Zhang     if (!usecprow) z += 4;
4252d61bbb3SSatish Balay   }
426d9fead3dSBarry Smith   ierr = VecRestoreArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
42726e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
42898c9bda7SSatish Balay   ierr = PetscLogFlops(32*a->nz - 4*nonzerorow);CHKERRQ(ierr);
4292d61bbb3SSatish Balay   PetscFunctionReturn(0);
4302d61bbb3SSatish Balay }
4312d61bbb3SSatish Balay 
4324a2ae208SSatish Balay #undef __FUNCT__
4334a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_5"
434dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_5(Mat A,Vec xx,Vec zz)
4352d61bbb3SSatish Balay {
4362d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
437d9fead3dSBarry Smith   PetscScalar       sum1,sum2,sum3,sum4,sum5,x1,x2,x3,x4,x5,*z = 0,*zarray;
438d9fead3dSBarry Smith   const PetscScalar *xb,*x;
439d9fead3dSBarry Smith   const MatScalar   *v;
440dfbe8321SBarry Smith   PetscErrorCode    ierr;
44198c9bda7SSatish Balay   PetscInt          mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL,nonzerorow=0;
44226e093fcSHong Zhang   PetscTruth        usecprow=a->compressedrow.use;
4432d61bbb3SSatish Balay 
444433994e6SBarry Smith   PetscFunctionBegin;
445d9fead3dSBarry Smith   ierr = VecGetArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
44626e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
4472d61bbb3SSatish Balay 
4482d61bbb3SSatish Balay   idx = a->j;
4492d61bbb3SSatish Balay   v   = a->a;
45026e093fcSHong Zhang   if (usecprow){
45126e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
45226e093fcSHong Zhang     ii   = a->compressedrow.i;
4537b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
45426e093fcSHong Zhang   } else {
45526e093fcSHong Zhang     mbs = a->mbs;
4562d61bbb3SSatish Balay     ii  = a->i;
45726e093fcSHong Zhang     z   = zarray;
45826e093fcSHong Zhang   }
4592d61bbb3SSatish Balay 
4602d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
4612d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
4622d61bbb3SSatish Balay     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0;
46398c9bda7SSatish Balay     nonzerorow += (n>0);
4642d61bbb3SSatish Balay     for (j=0; j<n; j++) {
4652d61bbb3SSatish Balay       xb = x + 5*(*idx++);
4662d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4];
4672d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[5]*x2 + v[10]*x3  + v[15]*x4 + v[20]*x5;
4682d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[6]*x2 + v[11]*x3  + v[16]*x4 + v[21]*x5;
4692d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[7]*x2 + v[12]*x3  + v[17]*x4 + v[22]*x5;
4702d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[8]*x2 + v[13]*x3  + v[18]*x4 + v[23]*x5;
4712d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[9]*x2 + v[14]*x3  + v[19]*x4 + v[24]*x5;
4722d61bbb3SSatish Balay       v += 25;
4732d61bbb3SSatish Balay     }
4747b2bb3b9SHong Zhang     if (usecprow) z = zarray + 5*ridx[i];
4752d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5;
47626e093fcSHong Zhang     if (!usecprow) z += 5;
4772d61bbb3SSatish Balay   }
478d9fead3dSBarry Smith   ierr = VecRestoreArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
47926e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
48098c9bda7SSatish Balay   ierr = PetscLogFlops(50*a->nz - 5*nonzerorow);CHKERRQ(ierr);
4812d61bbb3SSatish Balay   PetscFunctionReturn(0);
4822d61bbb3SSatish Balay }
4832d61bbb3SSatish Balay 
48415091d37SBarry Smith 
4854a2ae208SSatish Balay #undef __FUNCT__
4864a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_6"
487dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_6(Mat A,Vec xx,Vec zz)
48815091d37SBarry Smith {
48915091d37SBarry Smith   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
490d9fead3dSBarry Smith   PetscScalar       *z = 0,sum1,sum2,sum3,sum4,sum5,sum6;
491d9fead3dSBarry Smith   const PetscScalar *x,*xb;
49226e093fcSHong Zhang   PetscScalar       x1,x2,x3,x4,x5,x6,*zarray;
493d9fead3dSBarry Smith   const MatScalar   *v;
494dfbe8321SBarry Smith   PetscErrorCode    ierr;
49598c9bda7SSatish Balay   PetscInt          mbs=a->mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL,nonzerorow=0;
49626e093fcSHong Zhang   PetscTruth        usecprow=a->compressedrow.use;
49715091d37SBarry Smith 
498433994e6SBarry Smith   PetscFunctionBegin;
499d9fead3dSBarry Smith   ierr = VecGetArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
50026e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
50115091d37SBarry Smith 
50215091d37SBarry Smith   idx = a->j;
50315091d37SBarry Smith   v   = a->a;
50426e093fcSHong Zhang   if (usecprow){
50526e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
50626e093fcSHong Zhang     ii   = a->compressedrow.i;
5077b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
50826e093fcSHong Zhang   } else {
50926e093fcSHong Zhang     mbs = a->mbs;
51015091d37SBarry Smith     ii  = a->i;
51126e093fcSHong Zhang     z   = zarray;
51226e093fcSHong Zhang   }
51315091d37SBarry Smith 
51415091d37SBarry Smith   for (i=0; i<mbs; i++) {
51515091d37SBarry Smith     n  = ii[1] - ii[0]; ii++;
51615091d37SBarry Smith     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0;
51798c9bda7SSatish Balay     nonzerorow += (n>0);
51815091d37SBarry Smith     for (j=0; j<n; j++) {
51915091d37SBarry Smith       xb = x + 6*(*idx++);
52015091d37SBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5];
52115091d37SBarry Smith       sum1 += v[0]*x1 + v[6]*x2  + v[12]*x3  + v[18]*x4 + v[24]*x5 + v[30]*x6;
52215091d37SBarry Smith       sum2 += v[1]*x1 + v[7]*x2  + v[13]*x3  + v[19]*x4 + v[25]*x5 + v[31]*x6;
52315091d37SBarry Smith       sum3 += v[2]*x1 + v[8]*x2  + v[14]*x3  + v[20]*x4 + v[26]*x5 + v[32]*x6;
52415091d37SBarry Smith       sum4 += v[3]*x1 + v[9]*x2  + v[15]*x3  + v[21]*x4 + v[27]*x5 + v[33]*x6;
52515091d37SBarry Smith       sum5 += v[4]*x1 + v[10]*x2 + v[16]*x3  + v[22]*x4 + v[28]*x5 + v[34]*x6;
52615091d37SBarry Smith       sum6 += v[5]*x1 + v[11]*x2 + v[17]*x3  + v[23]*x4 + v[29]*x5 + v[35]*x6;
52715091d37SBarry Smith       v += 36;
52815091d37SBarry Smith     }
5297b2bb3b9SHong Zhang     if (usecprow) z = zarray + 6*ridx[i];
53015091d37SBarry Smith     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6;
53126e093fcSHong Zhang     if (!usecprow) z += 6;
53215091d37SBarry Smith   }
53315091d37SBarry Smith 
534d9fead3dSBarry Smith   ierr = VecRestoreArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
53526e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
53698c9bda7SSatish Balay   ierr = PetscLogFlops(72*a->nz - 6*nonzerorow);CHKERRQ(ierr);
53715091d37SBarry Smith   PetscFunctionReturn(0);
53815091d37SBarry Smith }
5394a2ae208SSatish Balay #undef __FUNCT__
5404a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_7"
541dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_7(Mat A,Vec xx,Vec zz)
5422d61bbb3SSatish Balay {
5432d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
544d9fead3dSBarry Smith   PetscScalar       *z = 0,sum1,sum2,sum3,sum4,sum5,sum6,sum7;
545d9fead3dSBarry Smith   const PetscScalar *x,*xb;
54626e093fcSHong Zhang   PetscScalar       x1,x2,x3,x4,x5,x6,x7,*zarray;
547d9fead3dSBarry Smith   const MatScalar   *v;
548dfbe8321SBarry Smith   PetscErrorCode    ierr;
54998c9bda7SSatish Balay   PetscInt          mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL,nonzerorow=0;
55026e093fcSHong Zhang   PetscTruth        usecprow=a->compressedrow.use;
5512d61bbb3SSatish Balay 
552433994e6SBarry Smith   PetscFunctionBegin;
553d9fead3dSBarry Smith   ierr = VecGetArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
55426e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
5552d61bbb3SSatish Balay 
5562d61bbb3SSatish Balay   idx = a->j;
5572d61bbb3SSatish Balay   v   = a->a;
55826e093fcSHong Zhang   if (usecprow){
55926e093fcSHong Zhang     mbs    = a->compressedrow.nrows;
56026e093fcSHong Zhang     ii     = a->compressedrow.i;
5617b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
56226e093fcSHong Zhang   } else {
56326e093fcSHong Zhang     mbs = a->mbs;
5642d61bbb3SSatish Balay     ii  = a->i;
56526e093fcSHong Zhang     z   = zarray;
56626e093fcSHong Zhang   }
5672d61bbb3SSatish Balay 
5682d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
5692d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
5702d61bbb3SSatish Balay     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0; sum7 = 0.0;
57198c9bda7SSatish Balay     nonzerorow += (n>0);
5722d61bbb3SSatish Balay     for (j=0; j<n; j++) {
5732d61bbb3SSatish Balay       xb = x + 7*(*idx++);
5742d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6];
5752d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[7]*x2  + v[14]*x3  + v[21]*x4 + v[28]*x5 + v[35]*x6 + v[42]*x7;
5762d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[8]*x2  + v[15]*x3  + v[22]*x4 + v[29]*x5 + v[36]*x6 + v[43]*x7;
5772d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[9]*x2  + v[16]*x3  + v[23]*x4 + v[30]*x5 + v[37]*x6 + v[44]*x7;
5782d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[10]*x2 + v[17]*x3  + v[24]*x4 + v[31]*x5 + v[38]*x6 + v[45]*x7;
5792d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[11]*x2 + v[18]*x3  + v[25]*x4 + v[32]*x5 + v[39]*x6 + v[46]*x7;
5802d61bbb3SSatish Balay       sum6 += v[5]*x1 + v[12]*x2 + v[19]*x3  + v[26]*x4 + v[33]*x5 + v[40]*x6 + v[47]*x7;
5812d61bbb3SSatish Balay       sum7 += v[6]*x1 + v[13]*x2 + v[20]*x3  + v[27]*x4 + v[34]*x5 + v[41]*x6 + v[48]*x7;
5822d61bbb3SSatish Balay       v += 49;
5832d61bbb3SSatish Balay     }
5847b2bb3b9SHong Zhang     if (usecprow) z = zarray + 7*ridx[i];
5852d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
58626e093fcSHong Zhang     if (!usecprow) z += 7;
5872d61bbb3SSatish Balay   }
5882d61bbb3SSatish Balay 
589d9fead3dSBarry Smith   ierr = VecRestoreArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
59026e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
59198c9bda7SSatish Balay   ierr = PetscLogFlops(98*a->nz - 7*nonzerorow);CHKERRQ(ierr);
5922d61bbb3SSatish Balay   PetscFunctionReturn(0);
5932d61bbb3SSatish Balay }
5942d61bbb3SSatish Balay 
5953f1db9ecSBarry Smith /*
5963f1db9ecSBarry Smith     This will not work with MatScalar == float because it calls the BLAS
5973f1db9ecSBarry Smith */
5984a2ae208SSatish Balay #undef __FUNCT__
5994a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_N"
600dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_N(Mat A,Vec xx,Vec zz)
6012d61bbb3SSatish Balay {
6022d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
603dcf5cc72SBarry Smith   PetscScalar    *x,*z = 0,*xb,*work,*workt,*zarray;
6043f1db9ecSBarry Smith   MatScalar      *v;
605dfbe8321SBarry Smith   PetscErrorCode ierr;
606d0f46423SBarry Smith   PetscInt       mbs=a->mbs,i,*idx,*ii,bs=A->rmap->bs,j,n,bs2=a->bs2;
60798c9bda7SSatish Balay   PetscInt       ncols,k,*ridx=PETSC_NULL,nonzerorow=0;
60826e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
6092d61bbb3SSatish Balay 
6102d61bbb3SSatish Balay   PetscFunctionBegin;
6111ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
61226e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
6132d61bbb3SSatish Balay 
6142d61bbb3SSatish Balay   idx = a->j;
6152d61bbb3SSatish Balay   v   = a->a;
61626e093fcSHong Zhang   if (usecprow){
61726e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
61826e093fcSHong Zhang     ii   = a->compressedrow.i;
6197b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
62026e093fcSHong Zhang   } else {
62126e093fcSHong Zhang     mbs = a->mbs;
6222d61bbb3SSatish Balay     ii  = a->i;
62326e093fcSHong Zhang     z   = zarray;
62426e093fcSHong Zhang   }
625218c64b6SSatish Balay 
6262d61bbb3SSatish Balay   if (!a->mult_work) {
627d0f46423SBarry Smith     k    = PetscMax(A->rmap->n,A->cmap->n);
62887828ca2SBarry Smith     ierr = PetscMalloc((k+1)*sizeof(PetscScalar),&a->mult_work);CHKERRQ(ierr);
6292d61bbb3SSatish Balay   }
6302d61bbb3SSatish Balay   work = a->mult_work;
6312d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
6322d61bbb3SSatish Balay     n     = ii[1] - ii[0]; ii++;
6332d61bbb3SSatish Balay     ncols = n*bs;
6342d61bbb3SSatish Balay     workt = work;
63598c9bda7SSatish Balay     nonzerorow += (n>0);
6362d61bbb3SSatish Balay     for (j=0; j<n; j++) {
6372d61bbb3SSatish Balay       xb = x + bs*(*idx++);
6382d61bbb3SSatish Balay       for (k=0; k<bs; k++) workt[k] = xb[k];
6392d61bbb3SSatish Balay       workt += bs;
6402d61bbb3SSatish Balay     }
6417b2bb3b9SHong Zhang     if (usecprow) z = zarray + bs*ridx[i];
642737d121aSSatish Balay     Kernel_w_gets_Ar_times_v(bs,ncols,work,v,z);
64371044d3cSBarry Smith     /* BLASgemv_("N",&bs,&ncols,&_DOne,v,&bs,work,&_One,&_DZero,z,&_One); */
6442d61bbb3SSatish Balay     v += n*bs2;
64526e093fcSHong Zhang     if (!usecprow) z += bs;
6462d61bbb3SSatish Balay   }
6471ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
64826e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
64998c9bda7SSatish Balay   ierr = PetscLogFlops(2*a->nz*bs2 - bs*nonzerorow);CHKERRQ(ierr);
6502d61bbb3SSatish Balay   PetscFunctionReturn(0);
6512d61bbb3SSatish Balay }
6522d61bbb3SSatish Balay 
6536a65c766SHong Zhang extern PetscErrorCode VecCopy_Seq(Vec,Vec);
6544a2ae208SSatish Balay #undef __FUNCT__
6554a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_1"
656dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_1(Mat A,Vec xx,Vec yy,Vec zz)
6572d61bbb3SSatish Balay {
6582d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
659dcf5cc72SBarry Smith   PetscScalar    *x,*y = 0,*z = 0,sum,*yarray,*zarray;
6603f1db9ecSBarry Smith   MatScalar      *v;
661dfbe8321SBarry Smith   PetscErrorCode ierr;
6624eb6d288SHong Zhang   PetscInt       mbs=a->mbs,i,*idx,*ii,n,*ridx=PETSC_NULL;
66326e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
6642d61bbb3SSatish Balay 
6652d61bbb3SSatish Balay   PetscFunctionBegin;
6661ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
66726e093fcSHong Zhang   ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr);
6682e8a6d31SBarry Smith   if (zz != yy) {
66926e093fcSHong Zhang     ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
6702e8a6d31SBarry Smith   } else {
67126e093fcSHong Zhang     zarray = yarray;
6722e8a6d31SBarry Smith   }
6732d61bbb3SSatish Balay 
6742d61bbb3SSatish Balay   idx = a->j;
6752d61bbb3SSatish Balay   v   = a->a;
67626e093fcSHong Zhang   if (usecprow){
6774eb6d288SHong Zhang     if (zz != yy){
6784eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,mbs*sizeof(PetscScalar));CHKERRQ(ierr);
6794eb6d288SHong Zhang     }
68026e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
68126e093fcSHong Zhang     ii   = a->compressedrow.i;
6827b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
68326e093fcSHong Zhang   } else {
6842d61bbb3SSatish Balay     ii  = a->i;
68526e093fcSHong Zhang     y   = yarray;
68626e093fcSHong Zhang     z   = zarray;
68726e093fcSHong Zhang   }
6882d61bbb3SSatish Balay 
6892d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
6902d61bbb3SSatish Balay     n    = ii[1] - ii[0]; ii++;
69126e093fcSHong Zhang     if (usecprow){
6927b2bb3b9SHong Zhang       z = zarray + ridx[i];
6937b2bb3b9SHong Zhang       y = yarray + ridx[i];
69426e093fcSHong Zhang     }
69526e093fcSHong Zhang     sum = y[0];
6962d61bbb3SSatish Balay     while (n--) sum += *v++ * x[*idx++];
69726e093fcSHong Zhang     z[0] = sum;
69826e093fcSHong Zhang     if (!usecprow){
69926e093fcSHong Zhang       z++; y++;
70026e093fcSHong Zhang     }
7012d61bbb3SSatish Balay   }
7021ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
70326e093fcSHong Zhang   ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr);
7042e8a6d31SBarry Smith   if (zz != yy) {
70526e093fcSHong Zhang     ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
7062e8a6d31SBarry Smith   }
707efee365bSSatish Balay   ierr = PetscLogFlops(2*a->nz);CHKERRQ(ierr);
7082d61bbb3SSatish Balay   PetscFunctionReturn(0);
7092d61bbb3SSatish Balay }
7102d61bbb3SSatish Balay 
7114a2ae208SSatish Balay #undef __FUNCT__
7124a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_2"
713dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_2(Mat A,Vec xx,Vec yy,Vec zz)
7142d61bbb3SSatish Balay {
7152d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
716dcf5cc72SBarry Smith   PetscScalar    *x,*y = 0,*z = 0,*xb,sum1,sum2;
71726e093fcSHong Zhang   PetscScalar    x1,x2,*yarray,*zarray;
7183f1db9ecSBarry Smith   MatScalar      *v;
719dfbe8321SBarry Smith   PetscErrorCode ierr;
7204eb6d288SHong Zhang   PetscInt       mbs=a->mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL;
72126e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
7222d61bbb3SSatish Balay 
7232d61bbb3SSatish Balay   PetscFunctionBegin;
7241ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
72526e093fcSHong Zhang   ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr);
7262e8a6d31SBarry Smith   if (zz != yy) {
72726e093fcSHong Zhang     ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
7282e8a6d31SBarry Smith   } else {
72926e093fcSHong Zhang     zarray = yarray;
7302e8a6d31SBarry Smith   }
7312d61bbb3SSatish Balay 
7322d61bbb3SSatish Balay   idx = a->j;
7332d61bbb3SSatish Balay   v   = a->a;
73426e093fcSHong Zhang   if (usecprow){
7354eb6d288SHong Zhang     if (zz != yy){
7364eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,2*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
7374eb6d288SHong Zhang     }
73826e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
73926e093fcSHong Zhang     ii   = a->compressedrow.i;
7407b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
7414eb6d288SHong Zhang     if (zz != yy){
7424eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,a->mbs*sizeof(PetscScalar));CHKERRQ(ierr);
7434eb6d288SHong Zhang     }
74426e093fcSHong Zhang   } else {
7452d61bbb3SSatish Balay     ii  = a->i;
74626e093fcSHong Zhang     y   = yarray;
74726e093fcSHong Zhang     z   = zarray;
74826e093fcSHong Zhang   }
7492d61bbb3SSatish Balay 
7502d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
7512d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
75226e093fcSHong Zhang     if (usecprow){
7537b2bb3b9SHong Zhang       z = zarray + 2*ridx[i];
7547b2bb3b9SHong Zhang       y = yarray + 2*ridx[i];
75526e093fcSHong Zhang     }
7562d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1];
7572d61bbb3SSatish Balay     for (j=0; j<n; j++) {
7582d61bbb3SSatish Balay       xb = x + 2*(*idx++); x1 = xb[0]; x2 = xb[1];
7592d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[2]*x2;
7602d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[3]*x2;
7612d61bbb3SSatish Balay       v += 4;
7622d61bbb3SSatish Balay     }
7632d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2;
76426e093fcSHong Zhang     if (!usecprow){
7652d61bbb3SSatish Balay       z += 2; y += 2;
7662d61bbb3SSatish Balay     }
76726e093fcSHong Zhang   }
7681ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
76926e093fcSHong Zhang   ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr);
7702e8a6d31SBarry Smith   if (zz != yy) {
77126e093fcSHong Zhang     ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
7722e8a6d31SBarry Smith   }
773efee365bSSatish Balay   ierr = PetscLogFlops(4*a->nz);CHKERRQ(ierr);
7742d61bbb3SSatish Balay   PetscFunctionReturn(0);
7752d61bbb3SSatish Balay }
7762d61bbb3SSatish Balay 
7774a2ae208SSatish Balay #undef __FUNCT__
7784a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_3"
779dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_3(Mat A,Vec xx,Vec yy,Vec zz)
7802d61bbb3SSatish Balay {
7812d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
782dcf5cc72SBarry Smith   PetscScalar    *x,*y = 0,*z = 0,*xb,sum1,sum2,sum3,x1,x2,x3,*yarray,*zarray;
7833f1db9ecSBarry Smith   MatScalar      *v;
784dfbe8321SBarry Smith   PetscErrorCode ierr;
7854eb6d288SHong Zhang   PetscInt       mbs=a->mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL;
78626e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
7872d61bbb3SSatish Balay 
7882d61bbb3SSatish Balay   PetscFunctionBegin;
7891ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
79026e093fcSHong Zhang   ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr);
7912e8a6d31SBarry Smith   if (zz != yy) {
79226e093fcSHong Zhang     ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
7932e8a6d31SBarry Smith   } else {
79426e093fcSHong Zhang     zarray = yarray;
7952e8a6d31SBarry Smith   }
7962d61bbb3SSatish Balay 
7972d61bbb3SSatish Balay   idx = a->j;
7982d61bbb3SSatish Balay   v   = a->a;
79926e093fcSHong Zhang   if (usecprow){
8004eb6d288SHong Zhang     if (zz != yy){
8014eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,3*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
8024eb6d288SHong Zhang     }
80326e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
80426e093fcSHong Zhang     ii   = a->compressedrow.i;
8057b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
80626e093fcSHong Zhang   } else {
8072d61bbb3SSatish Balay     ii  = a->i;
80826e093fcSHong Zhang     y   = yarray;
80926e093fcSHong Zhang     z   = zarray;
81026e093fcSHong Zhang   }
8112d61bbb3SSatish Balay 
8122d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
8132d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
81426e093fcSHong Zhang     if (usecprow){
8157b2bb3b9SHong Zhang       z = zarray + 3*ridx[i];
8167b2bb3b9SHong Zhang       y = yarray + 3*ridx[i];
81726e093fcSHong Zhang     }
8182d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2];
8192d61bbb3SSatish Balay     for (j=0; j<n; j++) {
8202d61bbb3SSatish Balay       xb = x + 3*(*idx++); x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
8212d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[3]*x2 + v[6]*x3;
8222d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[4]*x2 + v[7]*x3;
8232d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[5]*x2 + v[8]*x3;
8242d61bbb3SSatish Balay       v += 9;
8252d61bbb3SSatish Balay     }
8262d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3;
82726e093fcSHong Zhang     if (!usecprow){
8282d61bbb3SSatish Balay       z += 3; y += 3;
8292d61bbb3SSatish Balay     }
83026e093fcSHong Zhang   }
8311ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
83226e093fcSHong Zhang   ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr);
8332e8a6d31SBarry Smith   if (zz != yy) {
83426e093fcSHong Zhang     ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
8352e8a6d31SBarry Smith   }
836efee365bSSatish Balay   ierr = PetscLogFlops(18*a->nz);CHKERRQ(ierr);
8372d61bbb3SSatish Balay   PetscFunctionReturn(0);
8382d61bbb3SSatish Balay }
8392d61bbb3SSatish Balay 
8404a2ae208SSatish Balay #undef __FUNCT__
8414a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_4"
842dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_4(Mat A,Vec xx,Vec yy,Vec zz)
8432d61bbb3SSatish Balay {
8442d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
845dcf5cc72SBarry Smith   PetscScalar    *x,*y = 0,*z = 0,*xb,sum1,sum2,sum3,sum4,x1,x2,x3,x4,*yarray,*zarray;
8463f1db9ecSBarry Smith   MatScalar      *v;
847dfbe8321SBarry Smith   PetscErrorCode ierr;
8484eb6d288SHong Zhang   PetscInt       mbs=a->mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL;
84926e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
8502d61bbb3SSatish Balay 
8512d61bbb3SSatish Balay   PetscFunctionBegin;
8521ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
85326e093fcSHong Zhang   ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr);
8542e8a6d31SBarry Smith   if (zz != yy) {
85526e093fcSHong Zhang     ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
8562e8a6d31SBarry Smith   } else {
85726e093fcSHong Zhang     zarray = yarray;
8582e8a6d31SBarry Smith   }
8592d61bbb3SSatish Balay 
8602d61bbb3SSatish Balay   idx   = a->j;
8612d61bbb3SSatish Balay   v     = a->a;
86226e093fcSHong Zhang   if (usecprow){
8634eb6d288SHong Zhang     if (zz != yy){
8644eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,4*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
8654eb6d288SHong Zhang     }
86626e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
86726e093fcSHong Zhang     ii   = a->compressedrow.i;
8687b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
86926e093fcSHong Zhang   } else {
8702d61bbb3SSatish Balay     ii  = a->i;
87126e093fcSHong Zhang     y   = yarray;
87226e093fcSHong Zhang     z   = zarray;
87326e093fcSHong Zhang   }
8742d61bbb3SSatish Balay 
8752d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
8762d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
87726e093fcSHong Zhang     if (usecprow){
8787b2bb3b9SHong Zhang       z = zarray + 4*ridx[i];
8797b2bb3b9SHong Zhang       y = yarray + 4*ridx[i];
88026e093fcSHong Zhang     }
8812d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3];
8822d61bbb3SSatish Balay     for (j=0; j<n; j++) {
8832d61bbb3SSatish Balay       xb = x + 4*(*idx++);
8842d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
8852d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[4]*x2 + v[8]*x3   + v[12]*x4;
8862d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[5]*x2 + v[9]*x3   + v[13]*x4;
8872d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[6]*x2 + v[10]*x3  + v[14]*x4;
8882d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[7]*x2 + v[11]*x3  + v[15]*x4;
8892d61bbb3SSatish Balay       v += 16;
8902d61bbb3SSatish Balay     }
8912d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4;
89226e093fcSHong Zhang     if (!usecprow){
8932d61bbb3SSatish Balay       z += 4; y += 4;
8942d61bbb3SSatish Balay     }
89526e093fcSHong Zhang   }
8961ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
89726e093fcSHong Zhang   ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr);
8982e8a6d31SBarry Smith   if (zz != yy) {
89926e093fcSHong Zhang     ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
9002e8a6d31SBarry Smith   }
901efee365bSSatish Balay   ierr = PetscLogFlops(32*a->nz);CHKERRQ(ierr);
9022d61bbb3SSatish Balay   PetscFunctionReturn(0);
9032d61bbb3SSatish Balay }
9042d61bbb3SSatish Balay 
9054a2ae208SSatish Balay #undef __FUNCT__
9064a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_5"
907dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_5(Mat A,Vec xx,Vec yy,Vec zz)
9082d61bbb3SSatish Balay {
9092d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
910dcf5cc72SBarry Smith   PetscScalar    *x,*y = 0,*z = 0,*xb,sum1,sum2,sum3,sum4,sum5,x1,x2,x3,x4,x5;
91126e093fcSHong Zhang   PetscScalar    *yarray,*zarray;
9123f1db9ecSBarry Smith   MatScalar      *v;
913dfbe8321SBarry Smith   PetscErrorCode ierr;
9144eb6d288SHong Zhang   PetscInt       mbs=a->mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL;
91526e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
9162d61bbb3SSatish Balay 
9172d61bbb3SSatish Balay   PetscFunctionBegin;
9181ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
91926e093fcSHong Zhang   ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr);
9202e8a6d31SBarry Smith   if (zz != yy) {
92126e093fcSHong Zhang     ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
9222e8a6d31SBarry Smith   } else {
92326e093fcSHong Zhang     zarray = yarray;
9242e8a6d31SBarry Smith   }
9252d61bbb3SSatish Balay 
9262d61bbb3SSatish Balay   idx = a->j;
9272d61bbb3SSatish Balay   v   = a->a;
92826e093fcSHong Zhang   if (usecprow){
9294eb6d288SHong Zhang     if (zz != yy){
9304eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,5*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
9314eb6d288SHong Zhang     }
93226e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
93326e093fcSHong Zhang     ii   = a->compressedrow.i;
9347b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
93526e093fcSHong Zhang   } else {
9362d61bbb3SSatish Balay     ii  = a->i;
93726e093fcSHong Zhang     y   = yarray;
93826e093fcSHong Zhang     z   = zarray;
93926e093fcSHong Zhang   }
9402d61bbb3SSatish Balay 
9412d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
9422d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
94326e093fcSHong Zhang     if (usecprow){
9447b2bb3b9SHong Zhang       z = zarray + 5*ridx[i];
9457b2bb3b9SHong Zhang       y = yarray + 5*ridx[i];
94626e093fcSHong Zhang     }
9472d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4];
9482d61bbb3SSatish Balay     for (j=0; j<n; j++) {
9492d61bbb3SSatish Balay       xb = x + 5*(*idx++);
9502d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4];
9512d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[5]*x2 + v[10]*x3  + v[15]*x4 + v[20]*x5;
9522d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[6]*x2 + v[11]*x3  + v[16]*x4 + v[21]*x5;
9532d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[7]*x2 + v[12]*x3  + v[17]*x4 + v[22]*x5;
9542d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[8]*x2 + v[13]*x3  + v[18]*x4 + v[23]*x5;
9552d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[9]*x2 + v[14]*x3  + v[19]*x4 + v[24]*x5;
9562d61bbb3SSatish Balay       v += 25;
9572d61bbb3SSatish Balay     }
9582d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5;
95926e093fcSHong Zhang     if (!usecprow){
9602d61bbb3SSatish Balay       z += 5; y += 5;
9612d61bbb3SSatish Balay     }
96226e093fcSHong Zhang   }
9631ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
96426e093fcSHong Zhang   ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr);
9652e8a6d31SBarry Smith   if (zz != yy) {
96626e093fcSHong Zhang     ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
9672e8a6d31SBarry Smith   }
968efee365bSSatish Balay   ierr = PetscLogFlops(50*a->nz);CHKERRQ(ierr);
9692d61bbb3SSatish Balay   PetscFunctionReturn(0);
9702d61bbb3SSatish Balay }
9714a2ae208SSatish Balay #undef __FUNCT__
9724a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_6"
973dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_6(Mat A,Vec xx,Vec yy,Vec zz)
97415091d37SBarry Smith {
97515091d37SBarry Smith   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
976dcf5cc72SBarry Smith   PetscScalar    *x,*y = 0,*z = 0,*xb,sum1,sum2,sum3,sum4,sum5,sum6;
97726e093fcSHong Zhang   PetscScalar    x1,x2,x3,x4,x5,x6,*yarray,*zarray;
97815091d37SBarry Smith   MatScalar      *v;
979dfbe8321SBarry Smith   PetscErrorCode ierr;
9804eb6d288SHong Zhang   PetscInt       mbs=a->mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL;
98126e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
98215091d37SBarry Smith 
98315091d37SBarry Smith   PetscFunctionBegin;
9841ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
98526e093fcSHong Zhang   ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr);
98615091d37SBarry Smith   if (zz != yy) {
98726e093fcSHong Zhang     ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
98815091d37SBarry Smith   } else {
98926e093fcSHong Zhang     zarray = yarray;
99015091d37SBarry Smith   }
99115091d37SBarry Smith 
99215091d37SBarry Smith   idx = a->j;
99315091d37SBarry Smith   v   = a->a;
99426e093fcSHong Zhang   if (usecprow){
9954eb6d288SHong Zhang     if (zz != yy){
9964eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,6*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
9974eb6d288SHong Zhang     }
99826e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
99926e093fcSHong Zhang     ii   = a->compressedrow.i;
10007b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
100126e093fcSHong Zhang   } else {
100215091d37SBarry Smith     ii  = a->i;
100326e093fcSHong Zhang     y   = yarray;
100426e093fcSHong Zhang     z   = zarray;
100526e093fcSHong Zhang   }
100615091d37SBarry Smith 
100715091d37SBarry Smith   for (i=0; i<mbs; i++) {
100815091d37SBarry Smith     n  = ii[1] - ii[0]; ii++;
100926e093fcSHong Zhang     if (usecprow){
10107b2bb3b9SHong Zhang       z = zarray + 6*ridx[i];
10117b2bb3b9SHong Zhang       y = yarray + 6*ridx[i];
101226e093fcSHong Zhang     }
101315091d37SBarry Smith     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4]; sum6 = y[5];
101415091d37SBarry Smith     for (j=0; j<n; j++) {
10153b95cb0eSSatish Balay       xb = x + 6*(*idx++);
101615091d37SBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5];
101715091d37SBarry Smith       sum1 += v[0]*x1 + v[6]*x2  + v[12]*x3  + v[18]*x4 + v[24]*x5 + v[30]*x6;
101815091d37SBarry Smith       sum2 += v[1]*x1 + v[7]*x2  + v[13]*x3  + v[19]*x4 + v[25]*x5 + v[31]*x6;
101915091d37SBarry Smith       sum3 += v[2]*x1 + v[8]*x2  + v[14]*x3  + v[20]*x4 + v[26]*x5 + v[32]*x6;
102015091d37SBarry Smith       sum4 += v[3]*x1 + v[9]*x2  + v[15]*x3  + v[21]*x4 + v[27]*x5 + v[33]*x6;
102115091d37SBarry Smith       sum5 += v[4]*x1 + v[10]*x2 + v[16]*x3  + v[22]*x4 + v[28]*x5 + v[34]*x6;
102215091d37SBarry Smith       sum6 += v[5]*x1 + v[11]*x2 + v[17]*x3  + v[23]*x4 + v[29]*x5 + v[35]*x6;
102315091d37SBarry Smith       v += 36;
102415091d37SBarry Smith     }
102515091d37SBarry Smith     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6;
102626e093fcSHong Zhang     if (!usecprow){
102715091d37SBarry Smith       z += 6; y += 6;
102815091d37SBarry Smith     }
102926e093fcSHong Zhang   }
10301ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
103126e093fcSHong Zhang   ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr);
103215091d37SBarry Smith   if (zz != yy) {
103326e093fcSHong Zhang     ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
103415091d37SBarry Smith   }
1035efee365bSSatish Balay   ierr = PetscLogFlops(72*a->nz);CHKERRQ(ierr);
103615091d37SBarry Smith   PetscFunctionReturn(0);
103715091d37SBarry Smith }
10382d61bbb3SSatish Balay 
10394a2ae208SSatish Balay #undef __FUNCT__
10404a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_7"
1041dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_7(Mat A,Vec xx,Vec yy,Vec zz)
10422d61bbb3SSatish Balay {
10432d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
1044dcf5cc72SBarry Smith   PetscScalar    *x,*y = 0,*z = 0,*xb,sum1,sum2,sum3,sum4,sum5,sum6,sum7;
104526e093fcSHong Zhang   PetscScalar    x1,x2,x3,x4,x5,x6,x7,*yarray,*zarray;
10463f1db9ecSBarry Smith   MatScalar      *v;
1047dfbe8321SBarry Smith   PetscErrorCode ierr;
10487b2bb3b9SHong Zhang   PetscInt       mbs=a->mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL;
104926e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
10502d61bbb3SSatish Balay 
10512d61bbb3SSatish Balay   PetscFunctionBegin;
10521ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
105326e093fcSHong Zhang   ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr);
10542e8a6d31SBarry Smith   if (zz != yy) {
105526e093fcSHong Zhang     ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
10562e8a6d31SBarry Smith   } else {
105726e093fcSHong Zhang     zarray = yarray;
10582e8a6d31SBarry Smith   }
10592d61bbb3SSatish Balay 
10602d61bbb3SSatish Balay   idx = a->j;
10612d61bbb3SSatish Balay   v   = a->a;
106226e093fcSHong Zhang   if (usecprow){
10634eb6d288SHong Zhang     if (zz != yy){
10644eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,7*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
10654eb6d288SHong Zhang     }
106626e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
106726e093fcSHong Zhang     ii   = a->compressedrow.i;
10687b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
106926e093fcSHong Zhang   } else {
10702d61bbb3SSatish Balay     ii  = a->i;
107126e093fcSHong Zhang     y   = yarray;
107226e093fcSHong Zhang     z   = zarray;
107326e093fcSHong Zhang   }
10742d61bbb3SSatish Balay 
10752d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
10762d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
107726e093fcSHong Zhang     if (usecprow){
10787b2bb3b9SHong Zhang       z = zarray + 7*ridx[i];
10797b2bb3b9SHong Zhang       y = yarray + 7*ridx[i];
108026e093fcSHong Zhang     }
10812d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4]; sum6 = y[5]; sum7 = y[6];
10822d61bbb3SSatish Balay     for (j=0; j<n; j++) {
10832d61bbb3SSatish Balay       xb = x + 7*(*idx++);
10842d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6];
10852d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[7]*x2  + v[14]*x3  + v[21]*x4 + v[28]*x5 + v[35]*x6 + v[42]*x7;
10862d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[8]*x2  + v[15]*x3  + v[22]*x4 + v[29]*x5 + v[36]*x6 + v[43]*x7;
10872d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[9]*x2  + v[16]*x3  + v[23]*x4 + v[30]*x5 + v[37]*x6 + v[44]*x7;
10882d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[10]*x2 + v[17]*x3  + v[24]*x4 + v[31]*x5 + v[38]*x6 + v[45]*x7;
10892d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[11]*x2 + v[18]*x3  + v[25]*x4 + v[32]*x5 + v[39]*x6 + v[46]*x7;
10902d61bbb3SSatish Balay       sum6 += v[5]*x1 + v[12]*x2 + v[19]*x3  + v[26]*x4 + v[33]*x5 + v[40]*x6 + v[47]*x7;
10912d61bbb3SSatish Balay       sum7 += v[6]*x1 + v[13]*x2 + v[20]*x3  + v[27]*x4 + v[34]*x5 + v[41]*x6 + v[48]*x7;
10922d61bbb3SSatish Balay       v += 49;
10932d61bbb3SSatish Balay     }
10942d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
109526e093fcSHong Zhang     if (!usecprow){
10962d61bbb3SSatish Balay       z += 7; y += 7;
10972d61bbb3SSatish Balay     }
109826e093fcSHong Zhang   }
10991ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
110026e093fcSHong Zhang   ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr);
11012e8a6d31SBarry Smith   if (zz != yy) {
110226e093fcSHong Zhang     ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
11032e8a6d31SBarry Smith   }
1104efee365bSSatish Balay   ierr = PetscLogFlops(98*a->nz);CHKERRQ(ierr);
11052d61bbb3SSatish Balay   PetscFunctionReturn(0);
11062d61bbb3SSatish Balay }
1107218c64b6SSatish Balay 
11084a2ae208SSatish Balay #undef __FUNCT__
11094a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_N"
1110dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_N(Mat A,Vec xx,Vec yy,Vec zz)
11112d61bbb3SSatish Balay {
11122d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
1113bba805e6SBarry Smith   PetscScalar    *x,*z = 0,*xb,*work,*workt,*zarray;
11143f1db9ecSBarry Smith   MatScalar      *v;
11156849ba73SBarry Smith   PetscErrorCode ierr;
1116d0f46423SBarry Smith   PetscInt       mbs,i,*idx,*ii,bs=A->rmap->bs,j,n,bs2=a->bs2;
11177b2bb3b9SHong Zhang   PetscInt       ncols,k,*ridx=PETSC_NULL;
111826e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
1119218c64b6SSatish Balay 
11202d61bbb3SSatish Balay   PetscFunctionBegin;
11216a65c766SHong Zhang   ierr = VecCopy_Seq(yy,zz);CHKERRQ(ierr);
11221ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
112326e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
11242d61bbb3SSatish Balay 
11252d61bbb3SSatish Balay   idx = a->j;
11262d61bbb3SSatish Balay   v   = a->a;
112726e093fcSHong Zhang   if (usecprow){
112826e093fcSHong Zhang     mbs    = a->compressedrow.nrows;
112926e093fcSHong Zhang     ii     = a->compressedrow.i;
11307b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
113126e093fcSHong Zhang   } else {
113226e093fcSHong Zhang     mbs = a->mbs;
11332d61bbb3SSatish Balay     ii  = a->i;
113426e093fcSHong Zhang     z   = zarray;
113526e093fcSHong Zhang   }
11362d61bbb3SSatish Balay 
11372d61bbb3SSatish Balay   if (!a->mult_work) {
1138d0f46423SBarry Smith     k    = PetscMax(A->rmap->n,A->cmap->n);
113987828ca2SBarry Smith     ierr = PetscMalloc((k+1)*sizeof(PetscScalar),&a->mult_work);CHKERRQ(ierr);
11402d61bbb3SSatish Balay   }
11412d61bbb3SSatish Balay   work = a->mult_work;
11422d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
11432d61bbb3SSatish Balay     n     = ii[1] - ii[0]; ii++;
11442d61bbb3SSatish Balay     ncols = n*bs;
11452d61bbb3SSatish Balay     workt = work;
11462d61bbb3SSatish Balay     for (j=0; j<n; j++) {
11472d61bbb3SSatish Balay       xb = x + bs*(*idx++);
11482d61bbb3SSatish Balay       for (k=0; k<bs; k++) workt[k] = xb[k];
11492d61bbb3SSatish Balay       workt += bs;
11502d61bbb3SSatish Balay     }
11517b2bb3b9SHong Zhang     if (usecprow) z = zarray + bs*ridx[i];
11523f1db9ecSBarry Smith     Kernel_w_gets_w_plus_Ar_times_v(bs,ncols,work,v,z);
115371044d3cSBarry Smith     /* BLASgemv_("N",&bs,&ncols,&_DOne,v,&bs,work,&_One,&_DOne,z,&_One); */
11542d61bbb3SSatish Balay     v += n*bs2;
115526e093fcSHong Zhang     if (!usecprow){
11562d61bbb3SSatish Balay       z += bs;
11572d61bbb3SSatish Balay     }
115826e093fcSHong Zhang   }
11591ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
116026e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
1161efee365bSSatish Balay   ierr = PetscLogFlops(2*a->nz*bs2);CHKERRQ(ierr);
11622d61bbb3SSatish Balay   PetscFunctionReturn(0);
11632d61bbb3SSatish Balay }
11642d61bbb3SSatish Balay 
11654a2ae208SSatish Balay #undef __FUNCT__
11664a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_SeqBAIJ"
1167dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_SeqBAIJ(Mat A,Vec xx,Vec zz)
11682d61bbb3SSatish Balay {
11693447b6efSHong Zhang   PetscScalar    zero = 0.0;
11706849ba73SBarry Smith   PetscErrorCode ierr;
11712d61bbb3SSatish Balay 
11722d61bbb3SSatish Balay   PetscFunctionBegin;
11732dcb1b2aSMatthew Knepley   ierr = VecSet(zz,zero);CHKERRQ(ierr);
11743447b6efSHong Zhang   ierr = MatMultTransposeAdd_SeqBAIJ(A,xx,zz,zz);CHKERRQ(ierr);
11752d61bbb3SSatish Balay   PetscFunctionReturn(0);
11762d61bbb3SSatish Balay }
11772d61bbb3SSatish Balay 
11784a2ae208SSatish Balay #undef __FUNCT__
11794a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_SeqBAIJ"
1180dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_SeqBAIJ(Mat A,Vec xx,Vec yy,Vec zz)
11812d61bbb3SSatish Balay 
11822d61bbb3SSatish Balay {
11832d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
1184dcf5cc72SBarry Smith   PetscScalar       *zb,*x,*z,*xb = 0,x1,x2,x3,x4,x5;
11853f1db9ecSBarry Smith   MatScalar         *v;
11866849ba73SBarry Smith   PetscErrorCode    ierr;
1187d0f46423SBarry Smith   PetscInt          mbs,i,*idx,*ii,rval,bs=A->rmap->bs,j,n,bs2=a->bs2,*ib,*ridx=PETSC_NULL;
11883447b6efSHong Zhang   Mat_CompressedRow cprow = a->compressedrow;
11893447b6efSHong Zhang   PetscTruth        usecprow=cprow.use;
11902d61bbb3SSatish Balay 
11912d61bbb3SSatish Balay   PetscFunctionBegin;
11926a65c766SHong Zhang   if (yy != zz) { ierr = VecCopy_Seq(yy,zz);CHKERRQ(ierr); }
11933447b6efSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
11943447b6efSHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
11952d61bbb3SSatish Balay 
11962d61bbb3SSatish Balay   idx = a->j;
11972d61bbb3SSatish Balay   v   = a->a;
11983447b6efSHong Zhang   if (usecprow){
11993447b6efSHong Zhang     mbs  = cprow.nrows;
12003447b6efSHong Zhang     ii   = cprow.i;
12017b2bb3b9SHong Zhang     ridx = cprow.rindex;
12023447b6efSHong Zhang   } else {
12033447b6efSHong Zhang     mbs=a->mbs;
12042d61bbb3SSatish Balay     ii = a->i;
1205f1af5d2fSBarry Smith     xb = x;
12063447b6efSHong Zhang   }
12072d61bbb3SSatish Balay 
12082d61bbb3SSatish Balay   switch (bs) {
12092d61bbb3SSatish Balay   case 1:
12102d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
12117b2bb3b9SHong Zhang       if (usecprow) xb = x + ridx[i];
1212f1af5d2fSBarry Smith       x1 = xb[0];
12133447b6efSHong Zhang       ib = idx + ii[0];
12143447b6efSHong Zhang       n  = ii[1] - ii[0]; ii++;
12152d61bbb3SSatish Balay       for (j=0; j<n; j++) {
12162d61bbb3SSatish Balay         rval    = ib[j];
1217f1af5d2fSBarry Smith         z[rval] += *v * x1;
1218f1af5d2fSBarry Smith         v++;
12192d61bbb3SSatish Balay       }
12203447b6efSHong Zhang       if (!usecprow) xb++;
12212d61bbb3SSatish Balay     }
12222d61bbb3SSatish Balay     break;
12232d61bbb3SSatish Balay   case 2:
12242d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
12257b2bb3b9SHong Zhang       if (usecprow) xb = x + 2*ridx[i];
1226f1af5d2fSBarry Smith       x1 = xb[0]; x2 = xb[1];
12273447b6efSHong Zhang       ib = idx + ii[0];
12283447b6efSHong Zhang       n  = ii[1] - ii[0]; ii++;
12292d61bbb3SSatish Balay       for (j=0; j<n; j++) {
12302d61bbb3SSatish Balay         rval      = ib[j]*2;
12312d61bbb3SSatish Balay         z[rval++] += v[0]*x1 + v[1]*x2;
12322d61bbb3SSatish Balay         z[rval++] += v[2]*x1 + v[3]*x2;
12332d61bbb3SSatish Balay         v  += 4;
12342d61bbb3SSatish Balay       }
12353447b6efSHong Zhang       if (!usecprow) xb += 2;
12362d61bbb3SSatish Balay     }
12372d61bbb3SSatish Balay     break;
12382d61bbb3SSatish Balay   case 3:
12392d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
12407b2bb3b9SHong Zhang       if (usecprow) xb = x + 3*ridx[i];
1241f1af5d2fSBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
12423447b6efSHong Zhang       ib = idx + ii[0];
12433447b6efSHong Zhang       n  = ii[1] - ii[0]; ii++;
12442d61bbb3SSatish Balay       for (j=0; j<n; j++) {
12452d61bbb3SSatish Balay         rval      = ib[j]*3;
12462d61bbb3SSatish Balay         z[rval++] += v[0]*x1 + v[1]*x2 + v[2]*x3;
12472d61bbb3SSatish Balay         z[rval++] += v[3]*x1 + v[4]*x2 + v[5]*x3;
12482d61bbb3SSatish Balay         z[rval++] += v[6]*x1 + v[7]*x2 + v[8]*x3;
12492d61bbb3SSatish Balay         v  += 9;
12502d61bbb3SSatish Balay       }
12513447b6efSHong Zhang       if (!usecprow) xb += 3;
12522d61bbb3SSatish Balay     }
12532d61bbb3SSatish Balay     break;
12542d61bbb3SSatish Balay   case 4:
12552d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
12567b2bb3b9SHong Zhang       if (usecprow) xb = x + 4*ridx[i];
1257f1af5d2fSBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
12583447b6efSHong Zhang       ib = idx + ii[0];
12593447b6efSHong Zhang       n  = ii[1] - ii[0]; ii++;
12602d61bbb3SSatish Balay       for (j=0; j<n; j++) {
12612d61bbb3SSatish Balay         rval      = ib[j]*4;
12622d61bbb3SSatish Balay         z[rval++] +=  v[0]*x1 +  v[1]*x2 +  v[2]*x3 +  v[3]*x4;
12632d61bbb3SSatish Balay         z[rval++] +=  v[4]*x1 +  v[5]*x2 +  v[6]*x3 +  v[7]*x4;
12642d61bbb3SSatish Balay         z[rval++] +=  v[8]*x1 +  v[9]*x2 + v[10]*x3 + v[11]*x4;
12652d61bbb3SSatish Balay         z[rval++] += v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4;
12662d61bbb3SSatish Balay         v  += 16;
12672d61bbb3SSatish Balay       }
12683447b6efSHong Zhang       if (!usecprow) xb += 4;
12692d61bbb3SSatish Balay     }
12702d61bbb3SSatish Balay     break;
12712d61bbb3SSatish Balay   case 5:
12722d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
12737b2bb3b9SHong Zhang       if (usecprow) xb = x + 5*ridx[i];
1274f1af5d2fSBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
12752d61bbb3SSatish Balay       x4 = xb[3]; x5 = xb[4];
12763447b6efSHong Zhang       ib = idx + ii[0];
12773447b6efSHong Zhang       n  = ii[1] - ii[0]; ii++;
12782d61bbb3SSatish Balay       for (j=0; j<n; j++) {
12792d61bbb3SSatish Balay         rval      = ib[j]*5;
12802d61bbb3SSatish Balay         z[rval++] +=  v[0]*x1 +  v[1]*x2 +  v[2]*x3 +  v[3]*x4 +  v[4]*x5;
12812d61bbb3SSatish Balay         z[rval++] +=  v[5]*x1 +  v[6]*x2 +  v[7]*x3 +  v[8]*x4 +  v[9]*x5;
12822d61bbb3SSatish Balay         z[rval++] += v[10]*x1 + v[11]*x2 + v[12]*x3 + v[13]*x4 + v[14]*x5;
12832d61bbb3SSatish Balay         z[rval++] += v[15]*x1 + v[16]*x2 + v[17]*x3 + v[18]*x4 + v[19]*x5;
12842d61bbb3SSatish Balay         z[rval++] += v[20]*x1 + v[21]*x2 + v[22]*x3 + v[23]*x4 + v[24]*x5;
12852d61bbb3SSatish Balay         v  += 25;
12862d61bbb3SSatish Balay       }
12873447b6efSHong Zhang       if (!usecprow) xb += 5;
12882d61bbb3SSatish Balay     }
12892d61bbb3SSatish Balay     break;
1290f1af5d2fSBarry Smith   default: {      /* block sizes larger then 5 by 5 are handled by BLAS */
1291690b6cddSBarry Smith       PetscInt     ncols,k;
12923447b6efSHong Zhang       PetscScalar  *work,*workt,*xtmp;
12933f1db9ecSBarry Smith 
12942d61bbb3SSatish Balay       if (!a->mult_work) {
1295d0f46423SBarry Smith         k = PetscMax(A->rmap->n,A->cmap->n);
129687828ca2SBarry Smith         ierr = PetscMalloc((k+1)*sizeof(PetscScalar),&a->mult_work);CHKERRQ(ierr);
12972d61bbb3SSatish Balay       }
12982d61bbb3SSatish Balay       work = a->mult_work;
12993447b6efSHong Zhang       xtmp = x;
13002d61bbb3SSatish Balay       for (i=0; i<mbs; i++) {
13012d61bbb3SSatish Balay         n     = ii[1] - ii[0]; ii++;
13022d61bbb3SSatish Balay         ncols = n*bs;
130387828ca2SBarry Smith         ierr  = PetscMemzero(work,ncols*sizeof(PetscScalar));CHKERRQ(ierr);
13043447b6efSHong Zhang         if (usecprow) {
13057b2bb3b9SHong Zhang           xtmp = x + bs*ridx[i];
13063447b6efSHong Zhang         }
13073447b6efSHong Zhang         Kernel_w_gets_w_plus_trans_Ar_times_v(bs,ncols,xtmp,v,work);
130871044d3cSBarry Smith         /* BLASgemv_("T",&bs,&ncols,&_DOne,v,&bs,xtmp,&_One,&_DOne,work,&_One); */
13092d61bbb3SSatish Balay         v += n*bs2;
13103447b6efSHong Zhang         if (!usecprow) xtmp += bs;
13112d61bbb3SSatish Balay         workt = work;
13122d61bbb3SSatish Balay         for (j=0; j<n; j++) {
13132d61bbb3SSatish Balay           zb = z + bs*(*idx++);
13142d61bbb3SSatish Balay           for (k=0; k<bs; k++) zb[k] += workt[k] ;
13152d61bbb3SSatish Balay           workt += bs;
13162d61bbb3SSatish Balay         }
13172d61bbb3SSatish Balay       }
13182d61bbb3SSatish Balay     }
13192d61bbb3SSatish Balay   }
13203447b6efSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
13213447b6efSHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
1322efee365bSSatish Balay   ierr = PetscLogFlops(2*a->nz*a->bs2);CHKERRQ(ierr);
13232d61bbb3SSatish Balay   PetscFunctionReturn(0);
13242d61bbb3SSatish Balay }
13252d61bbb3SSatish Balay 
13264a2ae208SSatish Balay #undef __FUNCT__
13274a2ae208SSatish Balay #define __FUNCT__ "MatScale_SeqBAIJ"
1328f4df32b1SMatthew Knepley PetscErrorCode MatScale_SeqBAIJ(Mat inA,PetscScalar alpha)
13292d61bbb3SSatish Balay {
13302d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)inA->data;
1331690b6cddSBarry Smith   PetscInt       totalnz = a->bs2*a->nz;
1332f4df32b1SMatthew Knepley   PetscScalar    oalpha = alpha;
1333efee365bSSatish Balay   PetscErrorCode ierr;
13340805154bSBarry Smith   PetscBLASInt   one = 1,tnz = PetscBLASIntCast(totalnz);
13352d61bbb3SSatish Balay 
13362d61bbb3SSatish Balay   PetscFunctionBegin;
1337f4df32b1SMatthew Knepley   BLASscal_(&tnz,&oalpha,a->a,&one);
1338efee365bSSatish Balay   ierr = PetscLogFlops(totalnz);CHKERRQ(ierr);
13392d61bbb3SSatish Balay   PetscFunctionReturn(0);
13402d61bbb3SSatish Balay }
13412d61bbb3SSatish Balay 
13424a2ae208SSatish Balay #undef __FUNCT__
13434a2ae208SSatish Balay #define __FUNCT__ "MatNorm_SeqBAIJ"
1344dfbe8321SBarry Smith PetscErrorCode MatNorm_SeqBAIJ(Mat A,NormType type,PetscReal *norm)
13452d61bbb3SSatish Balay {
13468a62d963SHong Zhang   PetscErrorCode ierr;
13472d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
13483f1db9ecSBarry Smith   MatScalar      *v = a->a;
1349329f5518SBarry Smith   PetscReal      sum = 0.0;
1350d0f46423SBarry Smith   PetscInt       i,j,k,bs=A->rmap->bs,nz=a->nz,bs2=a->bs2,k1;
13512d61bbb3SSatish Balay 
13522d61bbb3SSatish Balay   PetscFunctionBegin;
13532d61bbb3SSatish Balay   if (type == NORM_FROBENIUS) {
13542d61bbb3SSatish Balay     for (i=0; i< bs2*nz; i++) {
1355aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
1356329f5518SBarry Smith       sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
13572d61bbb3SSatish Balay #else
13582d61bbb3SSatish Balay       sum += (*v)*(*v); v++;
13592d61bbb3SSatish Balay #endif
13602d61bbb3SSatish Balay     }
13612d61bbb3SSatish Balay     *norm = sqrt(sum);
13628a62d963SHong Zhang   } else if (type == NORM_1) { /* maximum column sum */
13638a62d963SHong Zhang     PetscReal *tmp;
13648a62d963SHong Zhang     PetscInt  *bcol = a->j;
1365d0f46423SBarry Smith     ierr = PetscMalloc((A->cmap->n+1)*sizeof(PetscReal),&tmp);CHKERRQ(ierr);
1366d0f46423SBarry Smith     ierr = PetscMemzero(tmp,A->cmap->n*sizeof(PetscReal));CHKERRQ(ierr);
13678a62d963SHong Zhang     for (i=0; i<nz; i++){
13688a62d963SHong Zhang       for (j=0; j<bs; j++){
13698a62d963SHong Zhang         k1 = bs*(*bcol) + j; /* column index */
13708a62d963SHong Zhang         for (k=0; k<bs; k++){
13718a62d963SHong Zhang           tmp[k1] += PetscAbsScalar(*v); v++;
13728a62d963SHong Zhang         }
13738a62d963SHong Zhang       }
13748a62d963SHong Zhang       bcol++;
13758a62d963SHong Zhang     }
13768a62d963SHong Zhang     *norm = 0.0;
1377d0f46423SBarry Smith     for (j=0; j<A->cmap->n; j++) {
13788a62d963SHong Zhang       if (tmp[j] > *norm) *norm = tmp[j];
13798a62d963SHong Zhang     }
13808a62d963SHong Zhang     ierr = PetscFree(tmp);CHKERRQ(ierr);
1381596552b5SBarry Smith   } else if (type == NORM_INFINITY) { /* maximum row sum */
1382596552b5SBarry Smith     *norm = 0.0;
1383596552b5SBarry Smith     for (k=0; k<bs; k++) {
138474f84c7bSSatish Balay       for (j=0; j<a->mbs; j++) {
1385596552b5SBarry Smith         v = a->a + bs2*a->i[j] + k;
1386596552b5SBarry Smith         sum = 0.0;
1387596552b5SBarry Smith         for (i=0; i<a->i[j+1]-a->i[j]; i++) {
13880e90e235SBarry Smith           for (k1=0; k1<bs; k1++){
1389596552b5SBarry Smith             sum += PetscAbsScalar(*v);
1390596552b5SBarry Smith             v   += bs;
13912d61bbb3SSatish Balay           }
13920e90e235SBarry Smith         }
1393596552b5SBarry Smith         if (sum > *norm) *norm = sum;
1394596552b5SBarry Smith       }
1395596552b5SBarry Smith     }
1396596552b5SBarry Smith   } else {
139729bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"No support for this norm yet");
13982d61bbb3SSatish Balay   }
13992d61bbb3SSatish Balay   PetscFunctionReturn(0);
14002d61bbb3SSatish Balay }
14012d61bbb3SSatish Balay 
14022d61bbb3SSatish Balay 
14034a2ae208SSatish Balay #undef __FUNCT__
14044a2ae208SSatish Balay #define __FUNCT__ "MatEqual_SeqBAIJ"
1405dfbe8321SBarry Smith PetscErrorCode MatEqual_SeqBAIJ(Mat A,Mat B,PetscTruth* flg)
14062d61bbb3SSatish Balay {
14072d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ *)A->data,*b = (Mat_SeqBAIJ *)B->data;
1408dfbe8321SBarry Smith   PetscErrorCode ierr;
14092d61bbb3SSatish Balay 
14102d61bbb3SSatish Balay   PetscFunctionBegin;
14112d61bbb3SSatish Balay   /* If the  matrix/block dimensions are not equal, or no of nonzeros or shift */
1412d0f46423SBarry Smith   if ((A->rmap->N != B->rmap->N) || (A->cmap->n != B->cmap->n) || (A->rmap->bs != B->rmap->bs)|| (a->nz != b->nz)) {
1413273d9f13SBarry Smith     *flg = PETSC_FALSE;
1414273d9f13SBarry Smith     PetscFunctionReturn(0);
14152d61bbb3SSatish Balay   }
14162d61bbb3SSatish Balay 
14172d61bbb3SSatish Balay   /* if the a->i are the same */
1418690b6cddSBarry Smith   ierr = PetscMemcmp(a->i,b->i,(a->mbs+1)*sizeof(PetscInt),flg);CHKERRQ(ierr);
1419abc0a331SBarry Smith   if (!*flg) {
14200f5bd95cSBarry Smith     PetscFunctionReturn(0);
14212d61bbb3SSatish Balay   }
14222d61bbb3SSatish Balay 
14232d61bbb3SSatish Balay   /* if a->j are the same */
1424690b6cddSBarry Smith   ierr = PetscMemcmp(a->j,b->j,(a->nz)*sizeof(PetscInt),flg);CHKERRQ(ierr);
1425abc0a331SBarry Smith   if (!*flg) {
14260f5bd95cSBarry Smith     PetscFunctionReturn(0);
14272d61bbb3SSatish Balay   }
14282d61bbb3SSatish Balay   /* if a->a are the same */
1429d0f46423SBarry Smith   ierr = PetscMemcmp(a->a,b->a,(a->nz)*(A->rmap->bs)*(B->rmap->bs)*sizeof(PetscScalar),flg);CHKERRQ(ierr);
14302d61bbb3SSatish Balay   PetscFunctionReturn(0);
14312d61bbb3SSatish Balay 
14322d61bbb3SSatish Balay }
14332d61bbb3SSatish Balay 
14344a2ae208SSatish Balay #undef __FUNCT__
14354a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_SeqBAIJ"
1436dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_SeqBAIJ(Mat A,Vec v)
14372d61bbb3SSatish Balay {
14382d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
1439dfbe8321SBarry Smith   PetscErrorCode ierr;
1440690b6cddSBarry Smith   PetscInt       i,j,k,n,row,bs,*ai,*aj,ambs,bs2;
144187828ca2SBarry Smith   PetscScalar    *x,zero = 0.0;
14423f1db9ecSBarry Smith   MatScalar      *aa,*aa_j;
14432d61bbb3SSatish Balay 
14442d61bbb3SSatish Balay   PetscFunctionBegin;
144529bbc08cSBarry Smith   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1446d0f46423SBarry Smith   bs   = A->rmap->bs;
14472d61bbb3SSatish Balay   aa   = a->a;
14482d61bbb3SSatish Balay   ai   = a->i;
14492d61bbb3SSatish Balay   aj   = a->j;
14502d61bbb3SSatish Balay   ambs = a->mbs;
14512d61bbb3SSatish Balay   bs2  = a->bs2;
14522d61bbb3SSatish Balay 
14532dcb1b2aSMatthew Knepley   ierr = VecSet(v,zero);CHKERRQ(ierr);
14541ebc52fbSHong Zhang   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1455e1311b90SBarry Smith   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
1456d0f46423SBarry Smith   if (n != A->rmap->N) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
14572d61bbb3SSatish Balay   for (i=0; i<ambs; i++) {
14582d61bbb3SSatish Balay     for (j=ai[i]; j<ai[i+1]; j++) {
14592d61bbb3SSatish Balay       if (aj[j] == i) {
14602d61bbb3SSatish Balay         row  = i*bs;
14612d61bbb3SSatish Balay         aa_j = aa+j*bs2;
14622d61bbb3SSatish Balay         for (k=0; k<bs2; k+=(bs+1),row++) x[row] = aa_j[k];
14632d61bbb3SSatish Balay         break;
14642d61bbb3SSatish Balay       }
14652d61bbb3SSatish Balay     }
14662d61bbb3SSatish Balay   }
14671ebc52fbSHong Zhang   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
14682d61bbb3SSatish Balay   PetscFunctionReturn(0);
14692d61bbb3SSatish Balay }
14702d61bbb3SSatish Balay 
14714a2ae208SSatish Balay #undef __FUNCT__
14724a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_SeqBAIJ"
1473dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_SeqBAIJ(Mat A,Vec ll,Vec rr)
14742d61bbb3SSatish Balay {
14752d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
147687828ca2SBarry Smith   PetscScalar    *l,*r,x,*li,*ri;
14773f1db9ecSBarry Smith   MatScalar      *aa,*v;
1478dfbe8321SBarry Smith   PetscErrorCode ierr;
1479690b6cddSBarry Smith   PetscInt       i,j,k,lm,rn,M,m,n,*ai,*aj,mbs,tmp,bs,bs2;
14802d61bbb3SSatish Balay 
14812d61bbb3SSatish Balay   PetscFunctionBegin;
14822d61bbb3SSatish Balay   ai  = a->i;
14832d61bbb3SSatish Balay   aj  = a->j;
14842d61bbb3SSatish Balay   aa  = a->a;
1485d0f46423SBarry Smith   m   = A->rmap->n;
1486d0f46423SBarry Smith   n   = A->cmap->n;
1487d0f46423SBarry Smith   bs  = A->rmap->bs;
14882d61bbb3SSatish Balay   mbs = a->mbs;
14892d61bbb3SSatish Balay   bs2 = a->bs2;
14902d61bbb3SSatish Balay   if (ll) {
14911ebc52fbSHong Zhang     ierr = VecGetArray(ll,&l);CHKERRQ(ierr);
1492e1311b90SBarry Smith     ierr = VecGetLocalSize(ll,&lm);CHKERRQ(ierr);
149329bbc08cSBarry Smith     if (lm != m) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vector wrong length");
14942d61bbb3SSatish Balay     for (i=0; i<mbs; i++) { /* for each block row */
14952d61bbb3SSatish Balay       M  = ai[i+1] - ai[i];
14962d61bbb3SSatish Balay       li = l + i*bs;
14972d61bbb3SSatish Balay       v  = aa + bs2*ai[i];
14982d61bbb3SSatish Balay       for (j=0; j<M; j++) { /* for each block */
14992d61bbb3SSatish Balay         for (k=0; k<bs2; k++) {
15002d61bbb3SSatish Balay           (*v++) *= li[k%bs];
15012d61bbb3SSatish Balay         }
15022d61bbb3SSatish Balay       }
15032d61bbb3SSatish Balay     }
15041ebc52fbSHong Zhang     ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr);
1505efee365bSSatish Balay     ierr = PetscLogFlops(a->nz);CHKERRQ(ierr);
15062d61bbb3SSatish Balay   }
15072d61bbb3SSatish Balay 
15082d61bbb3SSatish Balay   if (rr) {
15091ebc52fbSHong Zhang     ierr = VecGetArray(rr,&r);CHKERRQ(ierr);
1510e1311b90SBarry Smith     ierr = VecGetLocalSize(rr,&rn);CHKERRQ(ierr);
151129bbc08cSBarry Smith     if (rn != n) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vector wrong length");
15122d61bbb3SSatish Balay     for (i=0; i<mbs; i++) { /* for each block row */
15132d61bbb3SSatish Balay       M  = ai[i+1] - ai[i];
15142d61bbb3SSatish Balay       v  = aa + bs2*ai[i];
15152d61bbb3SSatish Balay       for (j=0; j<M; j++) { /* for each block */
15162d61bbb3SSatish Balay         ri = r + bs*aj[ai[i]+j];
15172d61bbb3SSatish Balay         for (k=0; k<bs; k++) {
15182d61bbb3SSatish Balay           x = ri[k];
15192d61bbb3SSatish Balay           for (tmp=0; tmp<bs; tmp++) (*v++) *= x;
15202d61bbb3SSatish Balay         }
15212d61bbb3SSatish Balay       }
15222d61bbb3SSatish Balay     }
15231ebc52fbSHong Zhang     ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr);
1524efee365bSSatish Balay     ierr = PetscLogFlops(a->nz);CHKERRQ(ierr);
15252d61bbb3SSatish Balay   }
15262d61bbb3SSatish Balay   PetscFunctionReturn(0);
15272d61bbb3SSatish Balay }
15282d61bbb3SSatish Balay 
15292d61bbb3SSatish Balay 
15304a2ae208SSatish Balay #undef __FUNCT__
15314a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_SeqBAIJ"
1532dfbe8321SBarry Smith PetscErrorCode MatGetInfo_SeqBAIJ(Mat A,MatInfoType flag,MatInfo *info)
15332d61bbb3SSatish Balay {
15342d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
15352d61bbb3SSatish Balay 
15362d61bbb3SSatish Balay   PetscFunctionBegin;
1537d0f46423SBarry Smith   info->rows_global    = (double)A->rmap->N;
1538d0f46423SBarry Smith   info->columns_global = (double)A->cmap->n;
1539d0f46423SBarry Smith   info->rows_local     = (double)A->rmap->n;
1540d0f46423SBarry Smith   info->columns_local  = (double)A->cmap->n;
15412d61bbb3SSatish Balay   info->block_size     = a->bs2;
15422d61bbb3SSatish Balay   info->nz_allocated   = a->maxnz;
15432d61bbb3SSatish Balay   info->nz_used        = a->bs2*a->nz;
15442d61bbb3SSatish Balay   info->nz_unneeded    = (double)(info->nz_allocated - info->nz_used);
15452d61bbb3SSatish Balay   info->assemblies   = A->num_ass;
15462d61bbb3SSatish Balay   info->mallocs      = a->reallocs;
15477adad957SLisandro Dalcin   info->memory       = ((PetscObject)A)->mem;
15482d61bbb3SSatish Balay   if (A->factor) {
15492d61bbb3SSatish Balay     info->fill_ratio_given  = A->info.fill_ratio_given;
15502d61bbb3SSatish Balay     info->fill_ratio_needed = A->info.fill_ratio_needed;
15512d61bbb3SSatish Balay     info->factor_mallocs    = A->info.factor_mallocs;
15522d61bbb3SSatish Balay   } else {
15532d61bbb3SSatish Balay     info->fill_ratio_given  = 0;
15542d61bbb3SSatish Balay     info->fill_ratio_needed = 0;
15552d61bbb3SSatish Balay     info->factor_mallocs    = 0;
15562d61bbb3SSatish Balay   }
15572d61bbb3SSatish Balay   PetscFunctionReturn(0);
15582d61bbb3SSatish Balay }
15592d61bbb3SSatish Balay 
15602d61bbb3SSatish Balay 
15614a2ae208SSatish Balay #undef __FUNCT__
15624a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_SeqBAIJ"
1563dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_SeqBAIJ(Mat A)
15642d61bbb3SSatish Balay {
15652d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
1566dfbe8321SBarry Smith   PetscErrorCode ierr;
15672d61bbb3SSatish Balay 
15682d61bbb3SSatish Balay   PetscFunctionBegin;
1569549d3d68SSatish Balay   ierr = PetscMemzero(a->a,a->bs2*a->i[a->mbs]*sizeof(MatScalar));CHKERRQ(ierr);
15702d61bbb3SSatish Balay   PetscFunctionReturn(0);
15712d61bbb3SSatish Balay }
1572