xref: /petsc/src/mat/impls/baij/seq/baij2.c (revision 3649974ff2fe2c90f327f29aac9f68ec978ade9e)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
2cac129eeSSatish Balay 
37c4f633dSBarry Smith #include "../src/mat/impls/baij/seq/baij.h"
4c60f0209SBarry Smith #include "../src/mat/blockinvert.h"
50a835dfdSSatish Balay #include "petscbt.h"
6f3da1532SBarry Smith #include "petscblaslapack.h"
7cac129eeSSatish Balay 
84a2ae208SSatish Balay #undef __FUNCT__
94a2ae208SSatish Balay #define __FUNCT__ "MatIncreaseOverlap_SeqBAIJ"
10690b6cddSBarry Smith PetscErrorCode MatIncreaseOverlap_SeqBAIJ(Mat A,PetscInt is_max,IS is[],PetscInt ov)
11a3192f15SSatish Balay {
12a3192f15SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
136849ba73SBarry Smith   PetscErrorCode ierr;
145d0c19d7SBarry Smith   PetscInt       row,i,j,k,l,m,n,*nidx,isz,val,ival;
155d0c19d7SBarry Smith   const PetscInt *idx;
16690b6cddSBarry Smith   PetscInt       start,end,*ai,*aj,bs,*nidx2;
17f1af5d2fSBarry Smith   PetscBT        table;
18a3192f15SSatish Balay 
193a40ed3dSBarry Smith   PetscFunctionBegin;
20a3192f15SSatish Balay   m     = a->mbs;
21a3192f15SSatish Balay   ai    = a->i;
22a3192f15SSatish Balay   aj    = a->j;
23d0f46423SBarry Smith   bs    = A->rmap->bs;
24a3192f15SSatish Balay 
25e32f2f54SBarry Smith   if (ov < 0)  SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap specified");
26a3192f15SSatish Balay 
276831982aSBarry Smith   ierr = PetscBTCreate(m,table);CHKERRQ(ierr);
28690b6cddSBarry Smith   ierr = PetscMalloc((m+1)*sizeof(PetscInt),&nidx);CHKERRQ(ierr);
29d0f46423SBarry Smith   ierr = PetscMalloc((A->rmap->N+1)*sizeof(PetscInt),&nidx2);CHKERRQ(ierr);
30a3192f15SSatish Balay 
31a3192f15SSatish Balay   for (i=0; i<is_max; i++) {
32a3192f15SSatish Balay     /* Initialise the two local arrays */
33a3192f15SSatish Balay     isz  = 0;
346831982aSBarry Smith     ierr = PetscBTMemzero(m,table);CHKERRQ(ierr);
35a3192f15SSatish Balay 
36a3192f15SSatish Balay     /* Extract the indices, assume there can be duplicate entries */
37a3192f15SSatish Balay     ierr = ISGetIndices(is[i],&idx);CHKERRQ(ierr);
38b9b97703SBarry Smith     ierr = ISGetLocalSize(is[i],&n);CHKERRQ(ierr);
39a3192f15SSatish Balay 
40a3192f15SSatish Balay     /* Enter these into the temp arrays i.e mark table[row], enter row into new index */
41a3192f15SSatish Balay     for (j=0; j<n ; ++j){
42218c64b6SSatish Balay       ival = idx[j]/bs; /* convert the indices into block indices */
43e32f2f54SBarry Smith       if (ival>=m) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"index greater than mat-dim");
44f1af5d2fSBarry Smith       if(!PetscBTLookupSet(table,ival)) { nidx[isz++] = ival;}
45a3192f15SSatish Balay     }
46a3192f15SSatish Balay     ierr = ISRestoreIndices(is[i],&idx);CHKERRQ(ierr);
47a3192f15SSatish Balay     ierr = ISDestroy(is[i]);CHKERRQ(ierr);
48a3192f15SSatish Balay 
49a3192f15SSatish Balay     k = 0;
50a3192f15SSatish Balay     for (j=0; j<ov; j++){ /* for each overlap*/
51a3192f15SSatish Balay       n = isz;
52a3192f15SSatish Balay       for (; k<n ; k++){ /* do only those rows in nidx[k], which are not done yet */
53a3192f15SSatish Balay         row   = nidx[k];
54a3192f15SSatish Balay         start = ai[row];
55a3192f15SSatish Balay         end   = ai[row+1];
56a3192f15SSatish Balay         for (l = start; l<end ; l++){
57a3192f15SSatish Balay           val = aj[l];
58f1af5d2fSBarry Smith           if (!PetscBTLookupSet(table,val)) {nidx[isz++] = val;}
59a3192f15SSatish Balay         }
60a3192f15SSatish Balay       }
61a3192f15SSatish Balay     }
62218c64b6SSatish Balay     /* expand the Index Set */
63218c64b6SSatish Balay     for (j=0; j<isz; j++) {
64218c64b6SSatish Balay       for (k=0; k<bs; k++)
65218c64b6SSatish Balay         nidx2[j*bs+k] = nidx[j]*bs+k;
66218c64b6SSatish Balay     }
67329f5518SBarry Smith     ierr = ISCreateGeneral(PETSC_COMM_SELF,isz*bs,nidx2,is+i);CHKERRQ(ierr);
68a3192f15SSatish Balay   }
696831982aSBarry Smith   ierr = PetscBTDestroy(table);CHKERRQ(ierr);
70606d414cSSatish Balay   ierr = PetscFree(nidx);CHKERRQ(ierr);
71606d414cSSatish Balay   ierr = PetscFree(nidx2);CHKERRQ(ierr);
723a40ed3dSBarry Smith   PetscFunctionReturn(0);
73a3192f15SSatish Balay }
741c351548SSatish Balay 
754a2ae208SSatish Balay #undef __FUNCT__
764a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqBAIJ_Private"
774aa3045dSJed Brown PetscErrorCode MatGetSubMatrix_SeqBAIJ_Private(Mat A,IS isrow,IS iscol,MatReuse scall,Mat *B)
78736121d4SSatish Balay {
79736121d4SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*c;
806849ba73SBarry Smith   PetscErrorCode ierr;
81690b6cddSBarry Smith   PetscInt       *smap,i,k,kstart,kend,oldcols = a->nbs,*lens;
82690b6cddSBarry Smith   PetscInt       row,mat_i,*mat_j,tcol,*mat_ilen;
835d0c19d7SBarry Smith   const PetscInt *irow,*icol;
845d0c19d7SBarry Smith   PetscInt       nrows,ncols,*ssmap,bs=A->rmap->bs,bs2=a->bs2;
85690b6cddSBarry Smith   PetscInt       *aj = a->j,*ai = a->i;
863f1db9ecSBarry Smith   MatScalar      *mat_a;
87736121d4SSatish Balay   Mat            C;
8814ca34e6SBarry Smith   PetscTruth     flag,sorted;
89736121d4SSatish Balay 
903a40ed3dSBarry Smith   PetscFunctionBegin;
9114ca34e6SBarry Smith   ierr = ISSorted(iscol,&sorted);CHKERRQ(ierr);
92e32f2f54SBarry Smith   if (!sorted) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"IS is not sorted");
93736121d4SSatish Balay 
94736121d4SSatish Balay   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
95218c64b6SSatish Balay   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
96b9b97703SBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
97b9b97703SBarry Smith   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
98736121d4SSatish Balay 
99690b6cddSBarry Smith   ierr = PetscMalloc((1+oldcols)*sizeof(PetscInt),&smap);CHKERRQ(ierr);
100736121d4SSatish Balay   ssmap = smap;
101690b6cddSBarry Smith   ierr = PetscMalloc((1+nrows)*sizeof(PetscInt),&lens);CHKERRQ(ierr);
102690b6cddSBarry Smith   ierr  = PetscMemzero(smap,oldcols*sizeof(PetscInt));CHKERRQ(ierr);
103736121d4SSatish Balay   for (i=0; i<ncols; i++) smap[icol[i]] = i+1;
104736121d4SSatish Balay   /* determine lens of each row */
105736121d4SSatish Balay   for (i=0; i<nrows; i++) {
106736121d4SSatish Balay     kstart  = ai[irow[i]];
107736121d4SSatish Balay     kend    = kstart + a->ilen[irow[i]];
108736121d4SSatish Balay     lens[i] = 0;
109736121d4SSatish Balay       for (k=kstart; k<kend; k++) {
110736121d4SSatish Balay         if (ssmap[aj[k]]) {
111736121d4SSatish Balay           lens[i]++;
112736121d4SSatish Balay         }
113736121d4SSatish Balay       }
114736121d4SSatish Balay     }
115736121d4SSatish Balay   /* Create and fill new matrix */
116736121d4SSatish Balay   if (scall == MAT_REUSE_MATRIX) {
117736121d4SSatish Balay     c = (Mat_SeqBAIJ *)((*B)->data);
118736121d4SSatish Balay 
119e32f2f54SBarry Smith     if (c->mbs!=nrows || c->nbs!=ncols || (*B)->rmap->bs!=bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Submatrix wrong size");
120690b6cddSBarry Smith     ierr = PetscMemcmp(c->ilen,lens,c->mbs *sizeof(PetscInt),&flag);CHKERRQ(ierr);
121e7e72b3dSBarry Smith     if (!flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong no of nonzeros");
122690b6cddSBarry Smith     ierr = PetscMemzero(c->ilen,c->mbs*sizeof(PetscInt));CHKERRQ(ierr);
123736121d4SSatish Balay     C = *B;
1243a40ed3dSBarry Smith   } else {
1257adad957SLisandro Dalcin     ierr = MatCreate(((PetscObject)A)->comm,&C);CHKERRQ(ierr);
126f69a0ea3SMatthew Knepley     ierr = MatSetSizes(C,nrows*bs,ncols*bs,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1277adad957SLisandro Dalcin     ierr = MatSetType(C,((PetscObject)A)->type_name);CHKERRQ(ierr);
128ab93d7beSBarry Smith     ierr = MatSeqBAIJSetPreallocation_SeqBAIJ(C,bs,0,lens);CHKERRQ(ierr);
129736121d4SSatish Balay   }
130736121d4SSatish Balay   c = (Mat_SeqBAIJ *)(C->data);
131736121d4SSatish Balay   for (i=0; i<nrows; i++) {
132736121d4SSatish Balay     row    = irow[i];
133736121d4SSatish Balay     kstart = ai[row];
134736121d4SSatish Balay     kend   = kstart + a->ilen[row];
135736121d4SSatish Balay     mat_i  = c->i[i];
136736121d4SSatish Balay     mat_j  = c->j + mat_i;
137218c64b6SSatish Balay     mat_a  = c->a + mat_i*bs2;
138736121d4SSatish Balay     mat_ilen = c->ilen + i;
139736121d4SSatish Balay     for (k=kstart; k<kend; k++) {
140736121d4SSatish Balay       if ((tcol=ssmap[a->j[k]])) {
141736121d4SSatish Balay         *mat_j++ = tcol - 1;
142549d3d68SSatish Balay         ierr     = PetscMemcpy(mat_a,a->a+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
143549d3d68SSatish Balay         mat_a   += bs2;
144736121d4SSatish Balay         (*mat_ilen)++;
145736121d4SSatish Balay       }
146736121d4SSatish Balay     }
147736121d4SSatish Balay   }
148218c64b6SSatish Balay 
149736121d4SSatish Balay   /* Free work space */
150736121d4SSatish Balay   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
151606d414cSSatish Balay   ierr = PetscFree(smap);CHKERRQ(ierr);
152606d414cSSatish Balay   ierr = PetscFree(lens);CHKERRQ(ierr);
1536d4a8577SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1546d4a8577SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
155736121d4SSatish Balay 
156736121d4SSatish Balay   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
157736121d4SSatish Balay   *B = C;
1583a40ed3dSBarry Smith   PetscFunctionReturn(0);
159736121d4SSatish Balay }
160736121d4SSatish Balay 
1614a2ae208SSatish Balay #undef __FUNCT__
1624a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqBAIJ"
1634aa3045dSJed Brown PetscErrorCode MatGetSubMatrix_SeqBAIJ(Mat A,IS isrow,IS iscol,MatReuse scall,Mat *B)
164218c64b6SSatish Balay {
165218c64b6SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
166218c64b6SSatish Balay   IS             is1,is2;
1676849ba73SBarry Smith   PetscErrorCode ierr;
1685d0c19d7SBarry Smith   PetscInt       *vary,*iary,nrows,ncols,i,bs=A->rmap->bs,count;
1695d0c19d7SBarry Smith   const PetscInt *irow,*icol;
170218c64b6SSatish Balay 
1713a40ed3dSBarry Smith   PetscFunctionBegin;
172218c64b6SSatish Balay   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
173218c64b6SSatish Balay   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
174b9b97703SBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
175b9b97703SBarry Smith   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
176218c64b6SSatish Balay 
177218c64b6SSatish Balay   /* Verify if the indices corespond to each element in a block
178218c64b6SSatish Balay    and form the IS with compressed IS */
179fca92195SBarry Smith   ierr = PetscMalloc2(a->mbs,PetscInt,&vary,a->mbs,PetscInt,&iary);CHKERRQ(ierr);
180fca92195SBarry Smith   ierr = PetscMemzero(vary,a->mbs*sizeof(PetscInt));CHKERRQ(ierr);
181218c64b6SSatish Balay   for (i=0; i<nrows; i++) vary[irow[i]/bs]++;
182218c64b6SSatish Balay   count = 0;
183218c64b6SSatish Balay   for (i=0; i<a->mbs; i++) {
184e32f2f54SBarry Smith     if (vary[i]!=0 && vary[i]!=bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Index set does not match blocks");
185218c64b6SSatish Balay     if (vary[i]==bs) iary[count++] = i;
186218c64b6SSatish Balay   }
187029af93fSBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,count,iary,&is1);CHKERRQ(ierr);
188218c64b6SSatish Balay 
189690b6cddSBarry Smith   ierr = PetscMemzero(vary,(a->mbs)*sizeof(PetscInt));CHKERRQ(ierr);
190218c64b6SSatish Balay   for (i=0; i<ncols; i++) vary[icol[i]/bs]++;
191218c64b6SSatish Balay   count = 0;
192218c64b6SSatish Balay   for (i=0; i<a->mbs; i++) {
193e32f2f54SBarry Smith     if (vary[i]!=0 && vary[i]!=bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Internal error in PETSc");
194218c64b6SSatish Balay     if (vary[i]==bs) iary[count++] = i;
195218c64b6SSatish Balay   }
196029af93fSBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,count,iary,&is2);CHKERRQ(ierr);
197218c64b6SSatish Balay   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
198218c64b6SSatish Balay   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
199fca92195SBarry Smith   ierr = PetscFree2(vary,iary);CHKERRQ(ierr);
200218c64b6SSatish Balay 
2014aa3045dSJed Brown   ierr = MatGetSubMatrix_SeqBAIJ_Private(A,is1,is2,scall,B);CHKERRQ(ierr);
202218c64b6SSatish Balay   ISDestroy(is1);
203218c64b6SSatish Balay   ISDestroy(is2);
2043a40ed3dSBarry Smith   PetscFunctionReturn(0);
205218c64b6SSatish Balay }
206218c64b6SSatish Balay 
2074a2ae208SSatish Balay #undef __FUNCT__
2084a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_SeqBAIJ"
209690b6cddSBarry Smith PetscErrorCode MatGetSubMatrices_SeqBAIJ(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
210736121d4SSatish Balay {
2116849ba73SBarry Smith   PetscErrorCode ierr;
212690b6cddSBarry Smith   PetscInt       i;
213736121d4SSatish Balay 
2143a40ed3dSBarry Smith   PetscFunctionBegin;
215736121d4SSatish Balay   if (scall == MAT_INITIAL_MATRIX) {
21682502324SSatish Balay     ierr = PetscMalloc((n+1)*sizeof(Mat),B);CHKERRQ(ierr);
217736121d4SSatish Balay   }
218736121d4SSatish Balay 
219736121d4SSatish Balay   for (i=0; i<n; i++) {
2204aa3045dSJed Brown     ierr = MatGetSubMatrix_SeqBAIJ(A,irow[i],icol[i],scall,&(*B)[i]);CHKERRQ(ierr);
221736121d4SSatish Balay   }
2223a40ed3dSBarry Smith   PetscFunctionReturn(0);
223736121d4SSatish Balay }
224218c64b6SSatish Balay 
225218c64b6SSatish Balay 
2262d61bbb3SSatish Balay /* -------------------------------------------------------*/
2272d61bbb3SSatish Balay /* Should check that shapes of vectors and matrices match */
2282d61bbb3SSatish Balay /* -------------------------------------------------------*/
2292d61bbb3SSatish Balay 
2304a2ae208SSatish Balay #undef __FUNCT__
2314a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_1"
232dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_1(Mat A,Vec xx,Vec zz)
2332d61bbb3SSatish Balay {
2342d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
235d9fead3dSBarry Smith   PetscScalar       *z,sum;
236d9fead3dSBarry Smith   const PetscScalar *x;
237d9fead3dSBarry Smith   const MatScalar   *v;
2386849ba73SBarry Smith   PetscErrorCode    ierr;
2392162cab8SBarry Smith   PetscInt          mbs,i,n,nonzerorow=0;
2402162cab8SBarry Smith   const PetscInt    *idx,*ii,*ridx=PETSC_NULL;
24126e093fcSHong Zhang   PetscTruth        usecprow=a->compressedrow.use;
2422d61bbb3SSatish Balay 
2432d61bbb3SSatish Balay   PetscFunctionBegin;
244*3649974fSBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
2451ebc52fbSHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
2462d61bbb3SSatish Balay 
24726e093fcSHong Zhang   if (usecprow){
24826e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
24926e093fcSHong Zhang     ii   = a->compressedrow.i;
2507b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
251a1c3900fSBarry Smith     ierr = PetscMemzero(z,mbs*sizeof(PetscScalar));CHKERRQ(ierr);
25226e093fcSHong Zhang   } else {
25326e093fcSHong Zhang     mbs = a->mbs;
2542d61bbb3SSatish Balay     ii  = a->i;
25526e093fcSHong Zhang   }
2562d61bbb3SSatish Balay 
2572d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
258ee54c7eeSHong Zhang     n    = ii[1] - ii[0];
259ee54c7eeSHong Zhang     v    = a->a + ii[0];
260ee54c7eeSHong Zhang     idx  = a->j + ii[0];
261ee54c7eeSHong Zhang     ii++;
2622d61bbb3SSatish Balay     sum  = 0.0;
2632162cab8SBarry Smith     PetscSparseDensePlusDot(sum,x,v,idx,n);
26426e093fcSHong Zhang     if (usecprow){
2657b2bb3b9SHong Zhang       z[ridx[i]] = sum;
26626e093fcSHong Zhang     } else {
267122f12eaSBarry Smith       nonzerorow += (n>0);
2682d61bbb3SSatish Balay       z[i] = sum;
2692d61bbb3SSatish Balay     }
27026e093fcSHong Zhang   }
271*3649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&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;
291*3649974fSBarry Smith   ierr = VecGetArrayRead(xx,&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   }
320*3649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&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;
344*3649974fSBarry Smith   ierr = VecGetArrayRead(xx,&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   }
374*3649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&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;
393*3649974fSBarry Smith   ierr = VecGetArrayRead(xx,&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   }
425*3649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&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;
440766f9fbaSBarry Smith   const PetscInt    *idx,*ii,*ridx=PETSC_NULL;
441766f9fbaSBarry Smith   PetscInt          mbs,i,j,n,nonzerorow=0;
44226e093fcSHong Zhang   PetscTruth        usecprow=a->compressedrow.use;
4432d61bbb3SSatish Balay 
444433994e6SBarry Smith   PetscFunctionBegin;
445*3649974fSBarry Smith   ierr = VecGetArrayRead(xx,&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   }
478*3649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
47926e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
480dc0b31edSSatish Balay   ierr = PetscLogFlops(50.0*a->nz - 5.0*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;
499*3649974fSBarry Smith   ierr = VecGetArrayRead(xx,&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 
534*3649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
53526e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
536dc0b31edSSatish Balay   ierr = PetscLogFlops(72.0*a->nz - 6.0*nonzerorow);CHKERRQ(ierr);
53715091d37SBarry Smith   PetscFunctionReturn(0);
53815091d37SBarry Smith }
5398ab949d8SShri Abhyankar 
5404a2ae208SSatish Balay #undef __FUNCT__
5414a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_7"
542dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_7(Mat A,Vec xx,Vec zz)
5432d61bbb3SSatish Balay {
5442d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
545d9fead3dSBarry Smith   PetscScalar       *z = 0,sum1,sum2,sum3,sum4,sum5,sum6,sum7;
546d9fead3dSBarry Smith   const PetscScalar *x,*xb;
54726e093fcSHong Zhang   PetscScalar       x1,x2,x3,x4,x5,x6,x7,*zarray;
548d9fead3dSBarry Smith   const MatScalar   *v;
549dfbe8321SBarry Smith   PetscErrorCode    ierr;
55098c9bda7SSatish Balay   PetscInt          mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL,nonzerorow=0;
55126e093fcSHong Zhang   PetscTruth        usecprow=a->compressedrow.use;
5522d61bbb3SSatish Balay 
553433994e6SBarry Smith   PetscFunctionBegin;
554*3649974fSBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
55526e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
5562d61bbb3SSatish Balay 
5572d61bbb3SSatish Balay   idx = a->j;
5582d61bbb3SSatish Balay   v   = a->a;
55926e093fcSHong Zhang   if (usecprow){
56026e093fcSHong Zhang     mbs    = a->compressedrow.nrows;
56126e093fcSHong Zhang     ii     = a->compressedrow.i;
5627b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
56326e093fcSHong Zhang   } else {
56426e093fcSHong Zhang     mbs = a->mbs;
5652d61bbb3SSatish Balay     ii  = a->i;
56626e093fcSHong Zhang     z   = zarray;
56726e093fcSHong Zhang   }
5682d61bbb3SSatish Balay 
5692d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
5702d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
5712d61bbb3SSatish Balay     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0; sum7 = 0.0;
57298c9bda7SSatish Balay     nonzerorow += (n>0);
5732d61bbb3SSatish Balay     for (j=0; j<n; j++) {
5742d61bbb3SSatish Balay       xb = x + 7*(*idx++);
5752d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6];
5762d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[7]*x2  + v[14]*x3  + v[21]*x4 + v[28]*x5 + v[35]*x6 + v[42]*x7;
5772d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[8]*x2  + v[15]*x3  + v[22]*x4 + v[29]*x5 + v[36]*x6 + v[43]*x7;
5782d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[9]*x2  + v[16]*x3  + v[23]*x4 + v[30]*x5 + v[37]*x6 + v[44]*x7;
5792d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[10]*x2 + v[17]*x3  + v[24]*x4 + v[31]*x5 + v[38]*x6 + v[45]*x7;
5802d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[11]*x2 + v[18]*x3  + v[25]*x4 + v[32]*x5 + v[39]*x6 + v[46]*x7;
5812d61bbb3SSatish Balay       sum6 += v[5]*x1 + v[12]*x2 + v[19]*x3  + v[26]*x4 + v[33]*x5 + v[40]*x6 + v[47]*x7;
5822d61bbb3SSatish Balay       sum7 += v[6]*x1 + v[13]*x2 + v[20]*x3  + v[27]*x4 + v[34]*x5 + v[41]*x6 + v[48]*x7;
5832d61bbb3SSatish Balay       v += 49;
5842d61bbb3SSatish Balay     }
5857b2bb3b9SHong Zhang     if (usecprow) z = zarray + 7*ridx[i];
5862d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
58726e093fcSHong Zhang     if (!usecprow) z += 7;
5882d61bbb3SSatish Balay   }
5892d61bbb3SSatish Balay 
590*3649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
59126e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
592dc0b31edSSatish Balay   ierr = PetscLogFlops(98.0*a->nz - 7.0*nonzerorow);CHKERRQ(ierr);
5932d61bbb3SSatish Balay   PetscFunctionReturn(0);
5942d61bbb3SSatish Balay }
5952d61bbb3SSatish Balay 
5968ab949d8SShri Abhyankar /* MatMult_SeqBAIJ_15 version 1: Columns in the block are accessed one at a time */
597832cc040SShri Abhyankar /* Default MatMult for block size 15 */
5988ab949d8SShri Abhyankar 
5998ab949d8SShri Abhyankar #undef __FUNCT__
6008ab949d8SShri Abhyankar #define __FUNCT__ "MatMult_SeqBAIJ_15_ver1"
6018ab949d8SShri Abhyankar PetscErrorCode MatMult_SeqBAIJ_15_ver1(Mat A,Vec xx,Vec zz)
6028ab949d8SShri Abhyankar {
6038ab949d8SShri Abhyankar   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
6048ab949d8SShri Abhyankar   PetscScalar       *z = 0,sum1,sum2,sum3,sum4,sum5,sum6,sum7,sum8,sum9,sum10,sum11,sum12,sum13,sum14,sum15;
6058ab949d8SShri Abhyankar   const PetscScalar *x,*xb;
6068ab949d8SShri Abhyankar   PetscScalar        x1,x2,x3,x4,x5,x6,x7,x8,x9,x10,x11,x12,x13,x14,x15,*zarray;
6078ab949d8SShri Abhyankar   const MatScalar   *v;
6088ab949d8SShri Abhyankar   PetscErrorCode    ierr;
6098ab949d8SShri Abhyankar   const PetscInt    *ii,*ij=a->j,*idx;
6108ab949d8SShri Abhyankar   PetscInt          mbs,i,j,k,n,*ridx=PETSC_NULL,nonzerorow=0;
6118ab949d8SShri Abhyankar   PetscTruth        usecprow=a->compressedrow.use;
6128ab949d8SShri Abhyankar 
6138ab949d8SShri Abhyankar   PetscFunctionBegin;
614*3649974fSBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
6158ab949d8SShri Abhyankar   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
6168ab949d8SShri Abhyankar 
6178ab949d8SShri Abhyankar   v   = a->a;
6188ab949d8SShri Abhyankar   if (usecprow){
6198ab949d8SShri Abhyankar     mbs    = a->compressedrow.nrows;
6208ab949d8SShri Abhyankar     ii     = a->compressedrow.i;
6218ab949d8SShri Abhyankar     ridx = a->compressedrow.rindex;
6228ab949d8SShri Abhyankar   } else {
6238ab949d8SShri Abhyankar     mbs = a->mbs;
6248ab949d8SShri Abhyankar     ii  = a->i;
6258ab949d8SShri Abhyankar     z   = zarray;
6268ab949d8SShri Abhyankar   }
6278ab949d8SShri Abhyankar 
6288ab949d8SShri Abhyankar   for (i=0; i<mbs; i++) {
6298ab949d8SShri Abhyankar     n  = ii[i+1] - ii[i];
6308ab949d8SShri Abhyankar     idx = ij + ii[i];
6318ab949d8SShri Abhyankar     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0; sum7 = 0.0;
6328ab949d8SShri Abhyankar     sum8 = 0.0; sum9 = 0.0; sum10 = 0.0; sum11 = 0.0; sum12 = 0.0; sum13 = 0.0; sum14 = 0.0;sum15 = 0.0;
6338ab949d8SShri Abhyankar 
6348ab949d8SShri Abhyankar     nonzerorow += (n>0);
6358ab949d8SShri Abhyankar     for (j=0; j<n; j++) {
6368ab949d8SShri Abhyankar       xb = x + 15*(idx[j]);
6378ab949d8SShri Abhyankar       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6];
6388ab949d8SShri Abhyankar       x8 = xb[7]; x9 = xb[8]; x10 = xb[9]; x11 = xb[10]; x12 = xb[11]; x13 = xb[12]; x14 = xb[13];x15 = xb[14];
6398ab949d8SShri Abhyankar 
6408ab949d8SShri Abhyankar       for(k=0;k<15;k++){
6418ab949d8SShri Abhyankar 	sum1  += v[0]*xb[k];
6428ab949d8SShri Abhyankar         sum2  += v[1]*xb[k];
6438ab949d8SShri Abhyankar 	sum3  += v[2]*xb[k];
6448ab949d8SShri Abhyankar 	sum4  += v[3]*xb[k];
6458ab949d8SShri Abhyankar 	sum5  += v[4]*xb[k];
6468ab949d8SShri Abhyankar         sum6  += v[5]*xb[k];
6478ab949d8SShri Abhyankar 	sum7  += v[6]*xb[k];
6488ab949d8SShri Abhyankar 	sum8  += v[7]*xb[k];
6498ab949d8SShri Abhyankar         sum9  += v[8]*xb[k];
6508ab949d8SShri Abhyankar         sum10 += v[9]*xb[k];
6518ab949d8SShri Abhyankar 	sum11 += v[10]*xb[k];
6528ab949d8SShri Abhyankar 	sum12 += v[11]*xb[k];
6538ab949d8SShri Abhyankar 	sum13 += v[12]*xb[k];
6548ab949d8SShri Abhyankar         sum14 += v[13]*xb[k];
6558ab949d8SShri Abhyankar 	sum15 += v[14]*xb[k];
6568ab949d8SShri Abhyankar 	v += 15;
6578ab949d8SShri Abhyankar       }
6588ab949d8SShri Abhyankar     }
6598ab949d8SShri Abhyankar     if (usecprow) z = zarray + 15*ridx[i];
6608ab949d8SShri Abhyankar     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
6618ab949d8SShri Abhyankar     z[7] = sum8; z[8] = sum9; z[9] = sum10; z[10] = sum11; z[11] = sum12; z[12] = sum13; z[13] = sum14;z[14] = sum15;
6628ab949d8SShri Abhyankar 
6638ab949d8SShri Abhyankar     if (!usecprow) z += 15;
6648ab949d8SShri Abhyankar   }
6658ab949d8SShri Abhyankar 
666*3649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
6678ab949d8SShri Abhyankar   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
6688ab949d8SShri Abhyankar   ierr = PetscLogFlops(450.0*a->nz - 15.0*nonzerorow);CHKERRQ(ierr);
6698ab949d8SShri Abhyankar   PetscFunctionReturn(0);
6708ab949d8SShri Abhyankar }
6718ab949d8SShri Abhyankar 
6728ab949d8SShri Abhyankar /* MatMult_SeqBAIJ_15_ver2 : Columns in the block are accessed in sets of 4,4,4,3 */
6738ab949d8SShri Abhyankar #undef __FUNCT__
6748ab949d8SShri Abhyankar #define __FUNCT__ "MatMult_SeqBAIJ_15_ver2"
6758ab949d8SShri Abhyankar PetscErrorCode MatMult_SeqBAIJ_15_ver2(Mat A,Vec xx,Vec zz)
6768ab949d8SShri Abhyankar {
6778ab949d8SShri Abhyankar   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
6788ab949d8SShri Abhyankar   PetscScalar       *z = 0,sum1,sum2,sum3,sum4,sum5,sum6,sum7,sum8,sum9,sum10,sum11,sum12,sum13,sum14,sum15;
6798ab949d8SShri Abhyankar   const PetscScalar *x,*xb;
6800b8f6341SShri Abhyankar   PetscScalar        x1,x2,x3,x4,*zarray;
6818ab949d8SShri Abhyankar   const MatScalar   *v;
6828ab949d8SShri Abhyankar   PetscErrorCode    ierr;
6838ab949d8SShri Abhyankar   const PetscInt    *ii,*ij=a->j,*idx;
6848ab949d8SShri Abhyankar   PetscInt          mbs,i,j,n,*ridx=PETSC_NULL,nonzerorow=0;
6858ab949d8SShri Abhyankar   PetscTruth        usecprow=a->compressedrow.use;
6868ab949d8SShri Abhyankar 
6878ab949d8SShri Abhyankar   PetscFunctionBegin;
688*3649974fSBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
6898ab949d8SShri Abhyankar   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
6908ab949d8SShri Abhyankar 
6918ab949d8SShri Abhyankar   v   = a->a;
6928ab949d8SShri Abhyankar   if (usecprow){
6938ab949d8SShri Abhyankar     mbs    = a->compressedrow.nrows;
6948ab949d8SShri Abhyankar     ii     = a->compressedrow.i;
6958ab949d8SShri Abhyankar     ridx = a->compressedrow.rindex;
6968ab949d8SShri Abhyankar   } else {
6978ab949d8SShri Abhyankar     mbs = a->mbs;
6988ab949d8SShri Abhyankar     ii  = a->i;
6998ab949d8SShri Abhyankar     z   = zarray;
7008ab949d8SShri Abhyankar   }
7018ab949d8SShri Abhyankar 
7028ab949d8SShri Abhyankar   for (i=0; i<mbs; i++) {
7038ab949d8SShri Abhyankar     n  = ii[i+1] - ii[i];
7048ab949d8SShri Abhyankar     idx = ij + ii[i];
7058ab949d8SShri Abhyankar     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0; sum7 = 0.0;
7068ab949d8SShri Abhyankar     sum8 = 0.0; sum9 = 0.0; sum10 = 0.0; sum11 = 0.0; sum12 = 0.0; sum13 = 0.0; sum14 = 0.0;sum15 = 0.0;
7078ab949d8SShri Abhyankar 
7088ab949d8SShri Abhyankar     nonzerorow += (n>0);
7098ab949d8SShri Abhyankar     for (j=0; j<n; j++) {
7108ab949d8SShri Abhyankar       xb = x + 15*(idx[j]);
7110b8f6341SShri Abhyankar       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
7128ab949d8SShri Abhyankar 
7138ab949d8SShri Abhyankar       sum1  += v[0]*x1 + v[15]*x2 + v[30]*x3   + v[45]*x4;
7148ab949d8SShri Abhyankar       sum2  += v[1]*x1 + v[16]*x2 + v[31]*x3   + v[46]*x4;
7158ab949d8SShri Abhyankar       sum3  += v[2]*x1 + v[17]*x2 + v[32]*x3  + v[47]*x4;
7168ab949d8SShri Abhyankar       sum4  += v[3]*x1 + v[18]*x2 + v[33]*x3  + v[48]*x4;
7178ab949d8SShri Abhyankar       sum5  += v[4]*x1 + v[19]*x2 + v[34]*x3   + v[49]*x4;
7188ab949d8SShri Abhyankar       sum6  += v[5]*x1 + v[20]*x2 + v[35]*x3   + v[50]*x4;
7198ab949d8SShri Abhyankar       sum7  += v[6]*x1 + v[21]*x2 + v[36]*x3  + v[51]*x4;
7208ab949d8SShri Abhyankar       sum8  += v[7]*x1 + v[22]*x2 + v[37]*x3  + v[52]*x4;
7218ab949d8SShri Abhyankar       sum9  += v[8]*x1 + v[23]*x2 + v[38]*x3   + v[53]*x4;
7228ab949d8SShri Abhyankar       sum10 += v[9]*x1 + v[24]*x2 + v[39]*x3   + v[54]*x4;
7238ab949d8SShri Abhyankar       sum11 += v[10]*x1 + v[25]*x2 + v[40]*x3  + v[55]*x4;
7248ab949d8SShri Abhyankar       sum12 += v[11]*x1 + v[26]*x2 + v[41]*x3  + v[56]*x4;
7258ab949d8SShri Abhyankar       sum13 += v[12]*x1 + v[27]*x2 + v[42]*x3   + v[57]*x4;
7268ab949d8SShri Abhyankar       sum14 += v[13]*x1 + v[28]*x2 + v[43]*x3   + v[58]*x4;
7278ab949d8SShri Abhyankar       sum15 += v[14]*x1 + v[29]*x2 + v[44]*x3  + v[59]*x4;
7288ab949d8SShri Abhyankar 
7298ab949d8SShri Abhyankar       v += 60;
7308ab949d8SShri Abhyankar 
7310b8f6341SShri Abhyankar       x1 = xb[4]; x2 = xb[5]; x3 = xb[6]; x4 = xb[7];
7328ab949d8SShri Abhyankar 
7338ab949d8SShri Abhyankar       sum1  += v[0]*x1 + v[15]*x2 + v[30]*x3   + v[45]*x4;
7348ab949d8SShri Abhyankar       sum2  += v[1]*x1 + v[16]*x2 + v[31]*x3   + v[46]*x4;
7358ab949d8SShri Abhyankar       sum3  += v[2]*x1 + v[17]*x2 + v[32]*x3  + v[47]*x4;
7368ab949d8SShri Abhyankar       sum4  += v[3]*x1 + v[18]*x2 + v[33]*x3  + v[48]*x4;
7378ab949d8SShri Abhyankar       sum5  += v[4]*x1 + v[19]*x2 + v[34]*x3   + v[49]*x4;
7388ab949d8SShri Abhyankar       sum6  += v[5]*x1 + v[20]*x2 + v[35]*x3   + v[50]*x4;
7398ab949d8SShri Abhyankar       sum7  += v[6]*x1 + v[21]*x2 + v[36]*x3  + v[51]*x4;
7408ab949d8SShri Abhyankar       sum8  += v[7]*x1 + v[22]*x2 + v[37]*x3  + v[52]*x4;
7418ab949d8SShri Abhyankar       sum9  += v[8]*x1 + v[23]*x2 + v[38]*x3   + v[53]*x4;
7428ab949d8SShri Abhyankar       sum10 += v[9]*x1 + v[24]*x2 + v[39]*x3   + v[54]*x4;
7438ab949d8SShri Abhyankar       sum11 += v[10]*x1 + v[25]*x2 + v[40]*x3  + v[55]*x4;
7448ab949d8SShri Abhyankar       sum12 += v[11]*x1 + v[26]*x2 + v[41]*x3  + v[56]*x4;
7458ab949d8SShri Abhyankar       sum13 += v[12]*x1 + v[27]*x2 + v[42]*x3   + v[57]*x4;
7468ab949d8SShri Abhyankar       sum14 += v[13]*x1 + v[28]*x2 + v[43]*x3   + v[58]*x4;
7478ab949d8SShri Abhyankar       sum15 += v[14]*x1 + v[29]*x2 + v[44]*x3  + v[59]*x4;
7488ab949d8SShri Abhyankar       v += 60;
7498ab949d8SShri Abhyankar 
7500b8f6341SShri Abhyankar       x1 = xb[8]; x2 = xb[9]; x3 = xb[10]; x4 = xb[11];
7510b8f6341SShri Abhyankar       sum1  += v[0]*x1 + v[15]*x2 + v[30]*x3   + v[45]*x4;
7520b8f6341SShri Abhyankar       sum2  += v[1]*x1 + v[16]*x2 + v[31]*x3   + v[46]*x4;
7530b8f6341SShri Abhyankar       sum3  += v[2]*x1 + v[17]*x2 + v[32]*x3  + v[47]*x4;
7540b8f6341SShri Abhyankar       sum4  += v[3]*x1 + v[18]*x2 + v[33]*x3  + v[48]*x4;
7550b8f6341SShri Abhyankar       sum5  += v[4]*x1 + v[19]*x2 + v[34]*x3   + v[49]*x4;
7560b8f6341SShri Abhyankar       sum6  += v[5]*x1 + v[20]*x2 + v[35]*x3   + v[50]*x4;
7570b8f6341SShri Abhyankar       sum7  += v[6]*x1 + v[21]*x2 + v[36]*x3  + v[51]*x4;
7580b8f6341SShri Abhyankar       sum8  += v[7]*x1 + v[22]*x2 + v[37]*x3  + v[52]*x4;
7590b8f6341SShri Abhyankar       sum9  += v[8]*x1 + v[23]*x2 + v[38]*x3   + v[53]*x4;
7600b8f6341SShri Abhyankar       sum10 += v[9]*x1 + v[24]*x2 + v[39]*x3   + v[54]*x4;
7610b8f6341SShri Abhyankar       sum11 += v[10]*x1 + v[25]*x2 + v[40]*x3  + v[55]*x4;
7620b8f6341SShri Abhyankar       sum12 += v[11]*x1 + v[26]*x2 + v[41]*x3  + v[56]*x4;
7630b8f6341SShri Abhyankar       sum13 += v[12]*x1 + v[27]*x2 + v[42]*x3   + v[57]*x4;
7640b8f6341SShri Abhyankar       sum14 += v[13]*x1 + v[28]*x2 + v[43]*x3   + v[58]*x4;
7650b8f6341SShri Abhyankar       sum15 += v[14]*x1 + v[29]*x2 + v[44]*x3  + v[59]*x4;
7660b8f6341SShri Abhyankar       v  += 60;
7670b8f6341SShri Abhyankar 
7680b8f6341SShri Abhyankar       x1 = xb[12]; x2 = xb[13]; x3 = xb[14];
7698ab949d8SShri Abhyankar       sum1  += v[0]*x1 + v[15]*x2 + v[30]*x3;
7708ab949d8SShri Abhyankar       sum2  += v[1]*x1 + v[16]*x2 + v[31]*x3;
7718ab949d8SShri Abhyankar       sum3  += v[2]*x1 + v[17]*x2 + v[32]*x3;
7728ab949d8SShri Abhyankar       sum4  += v[3]*x1 + v[18]*x2 + v[33]*x3;
7738ab949d8SShri Abhyankar       sum5  += v[4]*x1 + v[19]*x2 + v[34]*x3;
7748ab949d8SShri Abhyankar       sum6  += v[5]*x1 + v[20]*x2 + v[35]*x3;
7758ab949d8SShri Abhyankar       sum7  += v[6]*x1 + v[21]*x2 + v[36]*x3;
7768ab949d8SShri Abhyankar       sum8  += v[7]*x1 + v[22]*x2 + v[37]*x3;
7778ab949d8SShri Abhyankar       sum9  += v[8]*x1 + v[23]*x2 + v[38]*x3;
7788ab949d8SShri Abhyankar       sum10 += v[9]*x1 + v[24]*x2 + v[39]*x3;
7798ab949d8SShri Abhyankar       sum11 += v[10]*x1 + v[25]*x2 + v[40]*x3;
7808ab949d8SShri Abhyankar       sum12 += v[11]*x1 + v[26]*x2 + v[41]*x3;
7818ab949d8SShri Abhyankar       sum13 += v[12]*x1 + v[27]*x2 + v[42]*x3;
7828ab949d8SShri Abhyankar       sum14 += v[13]*x1 + v[28]*x2 + v[43]*x3;
7838ab949d8SShri Abhyankar       sum15 += v[14]*x1 + v[29]*x2 + v[44]*x3;
7848ab949d8SShri Abhyankar       v += 45;
7858ab949d8SShri Abhyankar     }
7868ab949d8SShri Abhyankar     if (usecprow) z = zarray + 15*ridx[i];
7878ab949d8SShri Abhyankar     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
7888ab949d8SShri Abhyankar     z[7] = sum8; z[8] = sum9; z[9] = sum10; z[10] = sum11; z[11] = sum12; z[12] = sum13; z[13] = sum14;z[14] = sum15;
7898ab949d8SShri Abhyankar 
7908ab949d8SShri Abhyankar     if (!usecprow) z += 15;
7918ab949d8SShri Abhyankar   }
7928ab949d8SShri Abhyankar 
793*3649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
7948ab949d8SShri Abhyankar   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
7958ab949d8SShri Abhyankar   ierr = PetscLogFlops(450.0*a->nz - 15.0*nonzerorow);CHKERRQ(ierr);
7968ab949d8SShri Abhyankar   PetscFunctionReturn(0);
7978ab949d8SShri Abhyankar }
7988ab949d8SShri Abhyankar 
7998ab949d8SShri Abhyankar /* MatMult_SeqBAIJ_15_ver3 : Columns in the block are accessed in sets of 8,7 */
8008ab949d8SShri Abhyankar #undef __FUNCT__
8018ab949d8SShri Abhyankar #define __FUNCT__ "MatMult_SeqBAIJ_15_ver3"
8028ab949d8SShri Abhyankar PetscErrorCode MatMult_SeqBAIJ_15_ver3(Mat A,Vec xx,Vec zz)
8038ab949d8SShri Abhyankar {
8048ab949d8SShri Abhyankar   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
8058ab949d8SShri Abhyankar   PetscScalar       *z = 0,sum1,sum2,sum3,sum4,sum5,sum6,sum7,sum8,sum9,sum10,sum11,sum12,sum13,sum14,sum15;
8068ab949d8SShri Abhyankar   const PetscScalar *x,*xb;
8070b8f6341SShri Abhyankar   PetscScalar        x1,x2,x3,x4,x5,x6,x7,x8,*zarray;
8088ab949d8SShri Abhyankar   const MatScalar   *v;
8098ab949d8SShri Abhyankar   PetscErrorCode    ierr;
8108ab949d8SShri Abhyankar   const PetscInt    *ii,*ij=a->j,*idx;
8118ab949d8SShri Abhyankar   PetscInt          mbs,i,j,n,*ridx=PETSC_NULL,nonzerorow=0;
8128ab949d8SShri Abhyankar   PetscTruth        usecprow=a->compressedrow.use;
8138ab949d8SShri Abhyankar 
8148ab949d8SShri Abhyankar   PetscFunctionBegin;
815*3649974fSBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
8168ab949d8SShri Abhyankar   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
8178ab949d8SShri Abhyankar 
8188ab949d8SShri Abhyankar   v   = a->a;
8198ab949d8SShri Abhyankar   if (usecprow){
8208ab949d8SShri Abhyankar     mbs    = a->compressedrow.nrows;
8218ab949d8SShri Abhyankar     ii     = a->compressedrow.i;
8228ab949d8SShri Abhyankar     ridx = a->compressedrow.rindex;
8238ab949d8SShri Abhyankar   } else {
8248ab949d8SShri Abhyankar     mbs = a->mbs;
8258ab949d8SShri Abhyankar     ii  = a->i;
8268ab949d8SShri Abhyankar     z   = zarray;
8278ab949d8SShri Abhyankar   }
8288ab949d8SShri Abhyankar 
8298ab949d8SShri Abhyankar   for (i=0; i<mbs; i++) {
8308ab949d8SShri Abhyankar     n  = ii[i+1] - ii[i];
8318ab949d8SShri Abhyankar     idx = ij + ii[i];
8328ab949d8SShri Abhyankar     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0; sum7 = 0.0;
8338ab949d8SShri Abhyankar     sum8 = 0.0; sum9 = 0.0; sum10 = 0.0; sum11 = 0.0; sum12 = 0.0; sum13 = 0.0; sum14 = 0.0;sum15 = 0.0;
8348ab949d8SShri Abhyankar 
8358ab949d8SShri Abhyankar     nonzerorow += (n>0);
8368ab949d8SShri Abhyankar     for (j=0; j<n; j++) {
8378ab949d8SShri Abhyankar       xb = x + 15*(idx[j]);
8388ab949d8SShri Abhyankar       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6];
8390b8f6341SShri Abhyankar       x8 = xb[7];
8408ab949d8SShri Abhyankar 
8418ab949d8SShri Abhyankar       sum1  += v[0]*x1 + v[15]*x2  + v[30]*x3  + v[45]*x4 + v[60]*x5 + v[75]*x6 + v[90]*x7 + v[105]*x8;
8428ab949d8SShri Abhyankar       sum2  += v[1]*x1 + v[16]*x2  + v[31]*x3  + v[46]*x4 + v[61]*x5 + v[76]*x6 + v[91]*x7 + v[106]*x8;
8438ab949d8SShri Abhyankar       sum3  += v[2]*x1 + v[17]*x2  + v[32]*x3  + v[47]*x4 + v[62]*x5 + v[77]*x6 + v[92]*x7 + v[107]*x8;
8448ab949d8SShri Abhyankar       sum4  += v[3]*x1 + v[18]*x2 + v[33]*x3  + v[48]*x4 + v[63]*x5 + v[78]*x6 + v[93]*x7 + v[108]*x8;
8458ab949d8SShri Abhyankar       sum5  += v[4]*x1 + v[19]*x2 + v[34]*x3  + v[49]*x4 + v[64]*x5 + v[79]*x6 + v[94]*x7 + v[109]*x8;
8468ab949d8SShri Abhyankar       sum6  += v[5]*x1 + v[20]*x2 + v[35]*x3  + v[50]*x4 + v[65]*x5 + v[80]*x6 + v[95]*x7 + v[110]*x8;
8478ab949d8SShri Abhyankar       sum7  += v[6]*x1 + v[21]*x2 + v[36]*x3  + v[51]*x4 + v[66]*x5 + v[81]*x6 + v[96]*x7 + v[111]*x8;
8488ab949d8SShri Abhyankar       sum8  += v[7]*x1 + v[22]*x2  + v[37]*x3  + v[52]*x4 + v[67]*x5 + v[82]*x6 + v[97]*x7 + v[112]*x8;
8498ab949d8SShri Abhyankar       sum9  += v[8]*x1 + v[23]*x2  + v[38]*x3  + v[53]*x4 + v[68]*x5 + v[83]*x6 + v[98]*x7 + v[113]*x8;
8508ab949d8SShri Abhyankar       sum10 += v[9]*x1 + v[24]*x2  + v[39]*x3  + v[54]*x4 + v[69]*x5 + v[84]*x6 + v[99]*x7 + v[114]*x8;
8518ab949d8SShri Abhyankar       sum11 += v[10]*x1 + v[25]*x2 + v[40]*x3  + v[55]*x4 + v[70]*x5 + v[85]*x6 + v[100]*x7 + v[115]*x8;
8528ab949d8SShri Abhyankar       sum12 += v[11]*x1 + v[26]*x2 + v[41]*x3  + v[56]*x4 + v[71]*x5 + v[86]*x6 + v[101]*x7 + v[116]*x8;
8538ab949d8SShri Abhyankar       sum13 += v[12]*x1 + v[27]*x2 + v[42]*x3  + v[57]*x4 + v[72]*x5 + v[87]*x6 + v[102]*x7 + v[117]*x8;
8548ab949d8SShri Abhyankar       sum14 += v[13]*x1 + v[28]*x2 + v[43]*x3  + v[58]*x4 + v[73]*x5 + v[88]*x6 + v[103]*x7 + v[118]*x8;
8558ab949d8SShri Abhyankar       sum15 += v[14]*x1 + v[29]*x2 + v[44]*x3  + v[59]*x4 + v[74]*x5 + v[89]*x6 + v[104]*x7 + v[119]*x8;
8568ab949d8SShri Abhyankar       v += 120;
8578ab949d8SShri Abhyankar 
8580b8f6341SShri Abhyankar       x1 = xb[8]; x2 = xb[9]; x3 = xb[10]; x4 = xb[11]; x5 = xb[12]; x6 = xb[13]; x7 = xb[14];
8590b8f6341SShri Abhyankar 
8608ab949d8SShri Abhyankar       sum1  += v[0]*x1 + v[15]*x2  + v[30]*x3  + v[45]*x4 + v[60]*x5 + v[75]*x6 + v[90]*x7;
8618ab949d8SShri Abhyankar       sum2  += v[1]*x1 + v[16]*x2  + v[31]*x3  + v[46]*x4 + v[61]*x5 + v[76]*x6 + v[91]*x7;
8628ab949d8SShri Abhyankar       sum3  += v[2]*x1 + v[17]*x2  + v[32]*x3  + v[47]*x4 + v[62]*x5 + v[77]*x6 + v[92]*x7;
8638ab949d8SShri Abhyankar       sum4  += v[3]*x1 + v[18]*x2 + v[33]*x3  + v[48]*x4 + v[63]*x5 + v[78]*x6 + v[93]*x7;
8648ab949d8SShri Abhyankar       sum5  += v[4]*x1 + v[19]*x2 + v[34]*x3  + v[49]*x4 + v[64]*x5 + v[79]*x6 + v[94]*x7;
8658ab949d8SShri Abhyankar       sum6  += v[5]*x1 + v[20]*x2 + v[35]*x3  + v[50]*x4 + v[65]*x5 + v[80]*x6 + v[95]*x7;
8668ab949d8SShri Abhyankar       sum7  += v[6]*x1 + v[21]*x2 + v[36]*x3  + v[51]*x4 + v[66]*x5 + v[81]*x6 + v[96]*x7;
8678ab949d8SShri Abhyankar       sum8  += v[7]*x1 + v[22]*x2  + v[37]*x3  + v[52]*x4 + v[67]*x5 + v[82]*x6 + v[97]*x7;
8688ab949d8SShri Abhyankar       sum9  += v[8]*x1 + v[23]*x2  + v[38]*x3  + v[53]*x4 + v[68]*x5 + v[83]*x6 + v[98]*x7;
8698ab949d8SShri Abhyankar       sum10 += v[9]*x1 + v[24]*x2  + v[39]*x3  + v[54]*x4 + v[69]*x5 + v[84]*x6 + v[99]*x7;
8708ab949d8SShri Abhyankar       sum11 += v[10]*x1 + v[25]*x2 + v[40]*x3  + v[55]*x4 + v[70]*x5 + v[85]*x6 + v[100]*x7;
8718ab949d8SShri Abhyankar       sum12 += v[11]*x1 + v[26]*x2 + v[41]*x3  + v[56]*x4 + v[71]*x5 + v[86]*x6 + v[101]*x7;
8728ab949d8SShri Abhyankar       sum13 += v[12]*x1 + v[27]*x2 + v[42]*x3  + v[57]*x4 + v[72]*x5 + v[87]*x6 + v[102]*x7;
8738ab949d8SShri Abhyankar       sum14 += v[13]*x1 + v[28]*x2 + v[43]*x3  + v[58]*x4 + v[73]*x5 + v[88]*x6 + v[103]*x7;
8748ab949d8SShri Abhyankar       sum15 += v[14]*x1 + v[29]*x2 + v[44]*x3  + v[59]*x4 + v[74]*x5 + v[89]*x6 + v[104]*x7;
8758ab949d8SShri Abhyankar       v += 105;
8768ab949d8SShri Abhyankar     }
8778ab949d8SShri Abhyankar     if (usecprow) z = zarray + 15*ridx[i];
8788ab949d8SShri Abhyankar     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
8798ab949d8SShri Abhyankar     z[7] = sum8; z[8] = sum9; z[9] = sum10; z[10] = sum11; z[11] = sum12; z[12] = sum13; z[13] = sum14;z[14] = sum15;
8808ab949d8SShri Abhyankar 
8818ab949d8SShri Abhyankar     if (!usecprow) z += 15;
8828ab949d8SShri Abhyankar   }
8838ab949d8SShri Abhyankar 
884*3649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
8858ab949d8SShri Abhyankar   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
8868ab949d8SShri Abhyankar   ierr = PetscLogFlops(450.0*a->nz - 15.0*nonzerorow);CHKERRQ(ierr);
8878ab949d8SShri Abhyankar   PetscFunctionReturn(0);
8888ab949d8SShri Abhyankar }
8898ab949d8SShri Abhyankar 
8908ab949d8SShri Abhyankar /* MatMult_SeqBAIJ_15_ver4 : All columns in the block are accessed at once */
8918ab949d8SShri Abhyankar 
8928ab949d8SShri Abhyankar #undef __FUNCT__
8938ab949d8SShri Abhyankar #define __FUNCT__ "MatMult_SeqBAIJ_15_ver4"
8948ab949d8SShri Abhyankar PetscErrorCode MatMult_SeqBAIJ_15_ver4(Mat A,Vec xx,Vec zz)
8958ab949d8SShri Abhyankar {
8968ab949d8SShri Abhyankar   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
8978ab949d8SShri Abhyankar   PetscScalar       *z = 0,sum1,sum2,sum3,sum4,sum5,sum6,sum7,sum8,sum9,sum10,sum11,sum12,sum13,sum14,sum15;
8988ab949d8SShri Abhyankar   const PetscScalar *x,*xb;
8998ab949d8SShri Abhyankar   PetscScalar        x1,x2,x3,x4,x5,x6,x7,x8,x9,x10,x11,x12,x13,x14,x15,*zarray;
9008ab949d8SShri Abhyankar   const MatScalar   *v;
9018ab949d8SShri Abhyankar   PetscErrorCode    ierr;
9028ab949d8SShri Abhyankar   const PetscInt    *ii,*ij=a->j,*idx;
9038ab949d8SShri Abhyankar   PetscInt          mbs,i,j,n,*ridx=PETSC_NULL,nonzerorow=0;
9048ab949d8SShri Abhyankar   PetscTruth        usecprow=a->compressedrow.use;
9058ab949d8SShri Abhyankar 
9068ab949d8SShri Abhyankar   PetscFunctionBegin;
907*3649974fSBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
9088ab949d8SShri Abhyankar   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
9098ab949d8SShri Abhyankar 
9108ab949d8SShri Abhyankar   v   = a->a;
9118ab949d8SShri Abhyankar   if (usecprow){
9128ab949d8SShri Abhyankar     mbs    = a->compressedrow.nrows;
9138ab949d8SShri Abhyankar     ii     = a->compressedrow.i;
9148ab949d8SShri Abhyankar     ridx = a->compressedrow.rindex;
9158ab949d8SShri Abhyankar   } else {
9168ab949d8SShri Abhyankar     mbs = a->mbs;
9178ab949d8SShri Abhyankar     ii  = a->i;
9188ab949d8SShri Abhyankar     z   = zarray;
9198ab949d8SShri Abhyankar   }
9208ab949d8SShri Abhyankar 
9218ab949d8SShri Abhyankar   for (i=0; i<mbs; i++) {
9228ab949d8SShri Abhyankar     n  = ii[i+1] - ii[i];
9238ab949d8SShri Abhyankar     idx = ij + ii[i];
9248ab949d8SShri Abhyankar     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0; sum7 = 0.0;
9258ab949d8SShri Abhyankar     sum8 = 0.0; sum9 = 0.0; sum10 = 0.0; sum11 = 0.0; sum12 = 0.0; sum13 = 0.0; sum14 = 0.0;sum15 = 0.0;
9268ab949d8SShri Abhyankar 
9278ab949d8SShri Abhyankar     nonzerorow += (n>0);
9288ab949d8SShri Abhyankar     for (j=0; j<n; j++) {
9298ab949d8SShri Abhyankar       xb = x + 15*(idx[j]);
9308ab949d8SShri Abhyankar       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6];
9318ab949d8SShri Abhyankar       x8 = xb[7]; x9 = xb[8]; x10 = xb[9]; x11 = xb[10]; x12 = xb[11]; x13 = xb[12]; x14 = xb[13];x15 = xb[14];
9328ab949d8SShri Abhyankar 
9338ab949d8SShri Abhyankar       sum1 +=  v[0]*x1  + v[15]*x2 + v[30]*x3 + v[45]*x4 + v[60]*x5 + v[75]*x6 + v[90]*x7  + v[105]*x8 + v[120]*x9 + v[135]*x10 + v[150]*x11 + v[165]*x12 + v[180]*x13 + v[195]*x14 + v[210]*x15;
9348ab949d8SShri Abhyankar       sum2 +=  v[1]*x1  + v[16]*x2 + v[31]*x3 + v[46]*x4 + v[61]*x5 + v[76]*x6 + v[91]*x7  + v[106]*x8 + v[121]*x9 + v[136]*x10 + v[151]*x11 + v[166]*x12 + v[181]*x13 + v[196]*x14 + v[211]*x15;
9358ab949d8SShri Abhyankar       sum3 +=  v[2]*x1  + v[17]*x2 + v[32]*x3 + v[47]*x4 + v[62]*x5 + v[77]*x6 + v[92]*x7  + v[107]*x8 + v[122]*x9 + v[137]*x10 + v[152]*x11 + v[167]*x12 + v[182]*x13 + v[197]*x14 + v[212]*x15;
9368ab949d8SShri Abhyankar       sum4 +=  v[3]*x1  + v[18]*x2 + v[33]*x3 + v[48]*x4 + v[63]*x5 + v[78]*x6 + v[93]*x7  + v[108]*x8 + v[123]*x9 + v[138]*x10 + v[153]*x11 + v[168]*x12 + v[183]*x13 + v[198]*x14 + v[213]*x15;
9378ab949d8SShri Abhyankar       sum5  += v[4]*x1  + v[19]*x2 + v[34]*x3 + v[49]*x4 + v[64]*x5 + v[79]*x6 + v[94]*x7  + v[109]*x8 + v[124]*x9 + v[139]*x10 + v[154]*x11 + v[169]*x12 + v[184]*x13 + v[199]*x14 + v[214]*x15;
9388ab949d8SShri Abhyankar       sum6  += v[5]*x1  + v[20]*x2 + v[35]*x3 + v[50]*x4 + v[65]*x5 + v[80]*x6 + v[95]*x7  + v[110]*x8 + v[125]*x9 + v[140]*x10 + v[155]*x11 + v[170]*x12 + v[185]*x13 + v[200]*x14 + v[215]*x15;
9398ab949d8SShri Abhyankar       sum7  += v[6]*x1  + v[21]*x2 + v[36]*x3 + v[51]*x4 + v[66]*x5 + v[81]*x6 + v[96]*x7  + v[111]*x8 + v[126]*x9 + v[141]*x10 + v[156]*x11 + v[171]*x12 + v[186]*x13 + v[201]*x14 + v[216]*x15;
9408ab949d8SShri Abhyankar       sum8  += v[7]*x1  + v[22]*x2 + v[37]*x3 + v[52]*x4 + v[67]*x5 + v[82]*x6 + v[97]*x7  + v[112]*x8 + v[127]*x9 + v[142]*x10 + v[157]*x11 + v[172]*x12 + v[187]*x13 + v[202]*x14 + v[217]*x15;
9418ab949d8SShri Abhyankar       sum9  += v[8]*x1  + v[23]*x2 + v[38]*x3 + v[53]*x4 + v[68]*x5 + v[83]*x6 + v[98]*x7  + v[113]*x8 + v[128]*x9 + v[143]*x10 + v[158]*x11 + v[173]*x12 + v[188]*x13 + v[203]*x14 + v[218]*x15;
9428ab949d8SShri Abhyankar       sum10 += v[9]*x1  + v[24]*x2 + v[39]*x3 + v[54]*x4 + v[69]*x5 + v[84]*x6 + v[99]*x7  + v[114]*x8 + v[129]*x9 + v[144]*x10 + v[159]*x11 + v[174]*x12 + v[189]*x13 + v[204]*x14 + v[219]*x15;
9438ab949d8SShri Abhyankar       sum11 += v[10]*x1 + v[25]*x2 + v[40]*x3 + v[55]*x4 + v[70]*x5 + v[85]*x6 + v[100]*x7 + v[115]*x8 + v[130]*x9 + v[145]*x10 + v[160]*x11 + v[175]*x12 + v[190]*x13 + v[205]*x14 + v[220]*x15;
9448ab949d8SShri Abhyankar       sum12 += v[11]*x1 + v[26]*x2 + v[41]*x3 + v[56]*x4 + v[71]*x5 + v[86]*x6 + v[101]*x7 + v[116]*x8 + v[131]*x9 + v[146]*x10 + v[161]*x11 + v[176]*x12 + v[191]*x13 + v[206]*x14 + v[221]*x15;
9458ab949d8SShri Abhyankar       sum13 += v[12]*x1 + v[27]*x2 + v[42]*x3 + v[57]*x4 + v[72]*x5 + v[87]*x6 + v[102]*x7 + v[117]*x8 + v[132]*x9 + v[147]*x10 + v[162]*x11 + v[177]*x12 + v[192]*x13 + v[207]*x14 + v[222]*x15;
9468ab949d8SShri Abhyankar       sum14 += v[13]*x1 + v[28]*x2 + v[43]*x3 + v[58]*x4 + v[73]*x5 + v[88]*x6 + v[103]*x7 + v[118]*x8 + v[133]*x9 + v[148]*x10 + v[163]*x11 + v[178]*x12 + v[193]*x13 + v[208]*x14 + v[223]*x15;
9478ab949d8SShri Abhyankar       sum15 += v[14]*x1 + v[29]*x2 + v[44]*x3 + v[59]*x4 + v[74]*x5 + v[89]*x6 + v[104]*x7 + v[119]*x8 + v[134]*x9 + v[149]*x10 + v[164]*x11 + v[179]*x12 + v[194]*x13 + v[209]*x14 + v[224]*x15;
9488ab949d8SShri Abhyankar       v += 225;
9498ab949d8SShri Abhyankar     }
9508ab949d8SShri Abhyankar     if (usecprow) z = zarray + 15*ridx[i];
9518ab949d8SShri Abhyankar     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
9528ab949d8SShri Abhyankar     z[7] = sum8; z[8] = sum9; z[9] = sum10; z[10] = sum11; z[11] = sum12; z[12] = sum13; z[13] = sum14;z[14] = sum15;
9538ab949d8SShri Abhyankar 
9548ab949d8SShri Abhyankar     if (!usecprow) z += 15;
9558ab949d8SShri Abhyankar   }
9568ab949d8SShri Abhyankar 
957*3649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
9588ab949d8SShri Abhyankar   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
9598ab949d8SShri Abhyankar   ierr = PetscLogFlops(450.0*a->nz - 15.0*nonzerorow);CHKERRQ(ierr);
9608ab949d8SShri Abhyankar   PetscFunctionReturn(0);
9618ab949d8SShri Abhyankar }
9628ab949d8SShri Abhyankar 
9638ab949d8SShri Abhyankar 
9643f1db9ecSBarry Smith /*
9653f1db9ecSBarry Smith     This will not work with MatScalar == float because it calls the BLAS
9663f1db9ecSBarry Smith */
9674a2ae208SSatish Balay #undef __FUNCT__
9684a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_N"
969dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_N(Mat A,Vec xx,Vec zz)
9702d61bbb3SSatish Balay {
9712d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
972dcf5cc72SBarry Smith   PetscScalar    *x,*z = 0,*xb,*work,*workt,*zarray;
9733f1db9ecSBarry Smith   MatScalar      *v;
974dfbe8321SBarry Smith   PetscErrorCode ierr;
975d0f46423SBarry Smith   PetscInt       mbs=a->mbs,i,*idx,*ii,bs=A->rmap->bs,j,n,bs2=a->bs2;
97698c9bda7SSatish Balay   PetscInt       ncols,k,*ridx=PETSC_NULL,nonzerorow=0;
97726e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
9782d61bbb3SSatish Balay 
9792d61bbb3SSatish Balay   PetscFunctionBegin;
9801ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
98126e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
9822d61bbb3SSatish Balay 
9832d61bbb3SSatish Balay   idx = a->j;
9842d61bbb3SSatish Balay   v   = a->a;
98526e093fcSHong Zhang   if (usecprow){
98626e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
98726e093fcSHong Zhang     ii   = a->compressedrow.i;
9887b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
98926e093fcSHong Zhang   } else {
99026e093fcSHong Zhang     mbs = a->mbs;
9912d61bbb3SSatish Balay     ii  = a->i;
99226e093fcSHong Zhang     z   = zarray;
99326e093fcSHong Zhang   }
994218c64b6SSatish Balay 
9952d61bbb3SSatish Balay   if (!a->mult_work) {
996d0f46423SBarry Smith     k    = PetscMax(A->rmap->n,A->cmap->n);
99787828ca2SBarry Smith     ierr = PetscMalloc((k+1)*sizeof(PetscScalar),&a->mult_work);CHKERRQ(ierr);
9982d61bbb3SSatish Balay   }
9992d61bbb3SSatish Balay   work = a->mult_work;
10002d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
10012d61bbb3SSatish Balay     n     = ii[1] - ii[0]; ii++;
10022d61bbb3SSatish Balay     ncols = n*bs;
10032d61bbb3SSatish Balay     workt = work;
100498c9bda7SSatish Balay     nonzerorow += (n>0);
10052d61bbb3SSatish Balay     for (j=0; j<n; j++) {
10062d61bbb3SSatish Balay       xb = x + bs*(*idx++);
10072d61bbb3SSatish Balay       for (k=0; k<bs; k++) workt[k] = xb[k];
10082d61bbb3SSatish Balay       workt += bs;
10092d61bbb3SSatish Balay     }
10107b2bb3b9SHong Zhang     if (usecprow) z = zarray + bs*ridx[i];
1011737d121aSSatish Balay     Kernel_w_gets_Ar_times_v(bs,ncols,work,v,z);
101271044d3cSBarry Smith     /* BLASgemv_("N",&bs,&ncols,&_DOne,v,&bs,work,&_One,&_DZero,z,&_One); */
10132d61bbb3SSatish Balay     v += n*bs2;
101426e093fcSHong Zhang     if (!usecprow) z += bs;
10152d61bbb3SSatish Balay   }
10161ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
101726e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
1018dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*a->nz*bs2 - bs*nonzerorow);CHKERRQ(ierr);
10192d61bbb3SSatish Balay   PetscFunctionReturn(0);
10202d61bbb3SSatish Balay }
10212d61bbb3SSatish Balay 
10226a65c766SHong Zhang extern PetscErrorCode VecCopy_Seq(Vec,Vec);
10234a2ae208SSatish Balay #undef __FUNCT__
10244a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_1"
1025dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_1(Mat A,Vec xx,Vec yy,Vec zz)
10262d61bbb3SSatish Balay {
10272d61bbb3SSatish Balay   Mat_SeqBAIJ        *a = (Mat_SeqBAIJ*)A->data;
1028122f12eaSBarry Smith   const PetscScalar  *x;
1029122f12eaSBarry Smith   PetscScalar        *y,*z,sum;
1030122f12eaSBarry Smith   const MatScalar    *v;
1031dfbe8321SBarry Smith   PetscErrorCode     ierr;
1032122f12eaSBarry Smith   PetscInt           mbs=a->mbs,i,n,*ridx=PETSC_NULL,nonzerorow=0;
1033122f12eaSBarry Smith   const PetscInt     *idx,*ii;
103426e093fcSHong Zhang   PetscTruth         usecprow=a->compressedrow.use;
10352d61bbb3SSatish Balay 
10362d61bbb3SSatish Balay   PetscFunctionBegin;
1037*3649974fSBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
1038122f12eaSBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
10392e8a6d31SBarry Smith   if (zz != yy) {
1040122f12eaSBarry Smith     ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
10412e8a6d31SBarry Smith   } else {
1042122f12eaSBarry Smith     z = y;
10432e8a6d31SBarry Smith   }
10442d61bbb3SSatish Balay 
10452d61bbb3SSatish Balay   idx = a->j;
10462d61bbb3SSatish Balay   v   = a->a;
104726e093fcSHong Zhang   if (usecprow){
10484eb6d288SHong Zhang     if (zz != yy){
1049122f12eaSBarry Smith       ierr = PetscMemcpy(z,y,mbs*sizeof(PetscScalar));CHKERRQ(ierr);
10504eb6d288SHong Zhang     }
105126e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
105226e093fcSHong Zhang     ii   = a->compressedrow.i;
10537b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
105426e093fcSHong Zhang   } else {
10552d61bbb3SSatish Balay     ii  = a->i;
105626e093fcSHong Zhang   }
10572d61bbb3SSatish Balay 
10582d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
1059122f12eaSBarry Smith     n    = ii[1] - ii[0];
1060122f12eaSBarry Smith     ii++;
106126e093fcSHong Zhang     if (!usecprow){
1062122f12eaSBarry Smith       nonzerorow += (n>0);
1063122f12eaSBarry Smith       sum = y[i];
1064122f12eaSBarry Smith     } else {
1065122f12eaSBarry Smith       sum = y[ridx[i]];
1066122f12eaSBarry Smith     }
1067122f12eaSBarry Smith     PetscSparseDensePlusDot(sum,x,v,idx,n);
1068122f12eaSBarry Smith     v += n;
1069122f12eaSBarry Smith     idx += n;
1070122f12eaSBarry Smith     if (usecprow){
1071122f12eaSBarry Smith       z[ridx[i]] = sum;
1072122f12eaSBarry Smith     } else {
1073122f12eaSBarry Smith       z[i] = sum;
107426e093fcSHong Zhang     }
10752d61bbb3SSatish Balay   }
1076*3649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
1077122f12eaSBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
10782e8a6d31SBarry Smith   if (zz != yy) {
1079122f12eaSBarry Smith     ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
10802e8a6d31SBarry Smith   }
1081122f12eaSBarry Smith   ierr = PetscLogFlops(2.0*a->nz - nonzerorow);CHKERRQ(ierr);
10822d61bbb3SSatish Balay   PetscFunctionReturn(0);
10832d61bbb3SSatish Balay }
10842d61bbb3SSatish Balay 
10854a2ae208SSatish Balay #undef __FUNCT__
10864a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_2"
1087dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_2(Mat A,Vec xx,Vec yy,Vec zz)
10882d61bbb3SSatish Balay {
10892d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
1090dcf5cc72SBarry Smith   PetscScalar    *x,*y = 0,*z = 0,*xb,sum1,sum2;
109126e093fcSHong Zhang   PetscScalar    x1,x2,*yarray,*zarray;
10923f1db9ecSBarry Smith   MatScalar      *v;
1093dfbe8321SBarry Smith   PetscErrorCode ierr;
10944eb6d288SHong Zhang   PetscInt       mbs=a->mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL;
109526e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
10962d61bbb3SSatish Balay 
10972d61bbb3SSatish Balay   PetscFunctionBegin;
10981ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
109926e093fcSHong Zhang   ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr);
11002e8a6d31SBarry Smith   if (zz != yy) {
110126e093fcSHong Zhang     ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
11022e8a6d31SBarry Smith   } else {
110326e093fcSHong Zhang     zarray = yarray;
11042e8a6d31SBarry Smith   }
11052d61bbb3SSatish Balay 
11062d61bbb3SSatish Balay   idx = a->j;
11072d61bbb3SSatish Balay   v   = a->a;
110826e093fcSHong Zhang   if (usecprow){
11094eb6d288SHong Zhang     if (zz != yy){
11104eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,2*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
11114eb6d288SHong Zhang     }
111226e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
111326e093fcSHong Zhang     ii   = a->compressedrow.i;
11147b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
11154eb6d288SHong Zhang     if (zz != yy){
11164eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,a->mbs*sizeof(PetscScalar));CHKERRQ(ierr);
11174eb6d288SHong Zhang     }
111826e093fcSHong Zhang   } else {
11192d61bbb3SSatish Balay     ii  = a->i;
112026e093fcSHong Zhang     y   = yarray;
112126e093fcSHong Zhang     z   = zarray;
112226e093fcSHong Zhang   }
11232d61bbb3SSatish Balay 
11242d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
11252d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
112626e093fcSHong Zhang     if (usecprow){
11277b2bb3b9SHong Zhang       z = zarray + 2*ridx[i];
11287b2bb3b9SHong Zhang       y = yarray + 2*ridx[i];
112926e093fcSHong Zhang     }
11302d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1];
11312d61bbb3SSatish Balay     for (j=0; j<n; j++) {
11322d61bbb3SSatish Balay       xb = x + 2*(*idx++); x1 = xb[0]; x2 = xb[1];
11332d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[2]*x2;
11342d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[3]*x2;
11352d61bbb3SSatish Balay       v += 4;
11362d61bbb3SSatish Balay     }
11372d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2;
113826e093fcSHong Zhang     if (!usecprow){
11392d61bbb3SSatish Balay       z += 2; y += 2;
11402d61bbb3SSatish Balay     }
114126e093fcSHong Zhang   }
11421ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
114326e093fcSHong Zhang   ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr);
11442e8a6d31SBarry Smith   if (zz != yy) {
114526e093fcSHong Zhang     ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
11462e8a6d31SBarry Smith   }
1147dc0b31edSSatish Balay   ierr = PetscLogFlops(4.0*a->nz);CHKERRQ(ierr);
11482d61bbb3SSatish Balay   PetscFunctionReturn(0);
11492d61bbb3SSatish Balay }
11502d61bbb3SSatish Balay 
11514a2ae208SSatish Balay #undef __FUNCT__
11524a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_3"
1153dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_3(Mat A,Vec xx,Vec yy,Vec zz)
11542d61bbb3SSatish Balay {
11552d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
1156dcf5cc72SBarry Smith   PetscScalar    *x,*y = 0,*z = 0,*xb,sum1,sum2,sum3,x1,x2,x3,*yarray,*zarray;
11573f1db9ecSBarry Smith   MatScalar      *v;
1158dfbe8321SBarry Smith   PetscErrorCode ierr;
11594eb6d288SHong Zhang   PetscInt       mbs=a->mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL;
116026e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
11612d61bbb3SSatish Balay 
11622d61bbb3SSatish Balay   PetscFunctionBegin;
11631ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
116426e093fcSHong Zhang   ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr);
11652e8a6d31SBarry Smith   if (zz != yy) {
116626e093fcSHong Zhang     ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
11672e8a6d31SBarry Smith   } else {
116826e093fcSHong Zhang     zarray = yarray;
11692e8a6d31SBarry Smith   }
11702d61bbb3SSatish Balay 
11712d61bbb3SSatish Balay   idx = a->j;
11722d61bbb3SSatish Balay   v   = a->a;
117326e093fcSHong Zhang   if (usecprow){
11744eb6d288SHong Zhang     if (zz != yy){
11754eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,3*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
11764eb6d288SHong Zhang     }
117726e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
117826e093fcSHong Zhang     ii   = a->compressedrow.i;
11797b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
118026e093fcSHong Zhang   } else {
11812d61bbb3SSatish Balay     ii  = a->i;
118226e093fcSHong Zhang     y   = yarray;
118326e093fcSHong Zhang     z   = zarray;
118426e093fcSHong Zhang   }
11852d61bbb3SSatish Balay 
11862d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
11872d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
118826e093fcSHong Zhang     if (usecprow){
11897b2bb3b9SHong Zhang       z = zarray + 3*ridx[i];
11907b2bb3b9SHong Zhang       y = yarray + 3*ridx[i];
119126e093fcSHong Zhang     }
11922d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2];
11932d61bbb3SSatish Balay     for (j=0; j<n; j++) {
11942d61bbb3SSatish Balay       xb = x + 3*(*idx++); x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
11952d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[3]*x2 + v[6]*x3;
11962d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[4]*x2 + v[7]*x3;
11972d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[5]*x2 + v[8]*x3;
11982d61bbb3SSatish Balay       v += 9;
11992d61bbb3SSatish Balay     }
12002d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3;
120126e093fcSHong Zhang     if (!usecprow){
12022d61bbb3SSatish Balay       z += 3; y += 3;
12032d61bbb3SSatish Balay     }
120426e093fcSHong Zhang   }
12051ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
120626e093fcSHong Zhang   ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr);
12072e8a6d31SBarry Smith   if (zz != yy) {
120826e093fcSHong Zhang     ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
12092e8a6d31SBarry Smith   }
1210dc0b31edSSatish Balay   ierr = PetscLogFlops(18.0*a->nz);CHKERRQ(ierr);
12112d61bbb3SSatish Balay   PetscFunctionReturn(0);
12122d61bbb3SSatish Balay }
12132d61bbb3SSatish Balay 
12144a2ae208SSatish Balay #undef __FUNCT__
12154a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_4"
1216dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_4(Mat A,Vec xx,Vec yy,Vec zz)
12172d61bbb3SSatish Balay {
12182d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
1219dcf5cc72SBarry Smith   PetscScalar    *x,*y = 0,*z = 0,*xb,sum1,sum2,sum3,sum4,x1,x2,x3,x4,*yarray,*zarray;
12203f1db9ecSBarry Smith   MatScalar      *v;
1221dfbe8321SBarry Smith   PetscErrorCode ierr;
12224eb6d288SHong Zhang   PetscInt       mbs=a->mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL;
122326e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
12242d61bbb3SSatish Balay 
12252d61bbb3SSatish Balay   PetscFunctionBegin;
12261ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
122726e093fcSHong Zhang   ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr);
12282e8a6d31SBarry Smith   if (zz != yy) {
122926e093fcSHong Zhang     ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
12302e8a6d31SBarry Smith   } else {
123126e093fcSHong Zhang     zarray = yarray;
12322e8a6d31SBarry Smith   }
12332d61bbb3SSatish Balay 
12342d61bbb3SSatish Balay   idx   = a->j;
12352d61bbb3SSatish Balay   v     = a->a;
123626e093fcSHong Zhang   if (usecprow){
12374eb6d288SHong Zhang     if (zz != yy){
12384eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,4*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
12394eb6d288SHong Zhang     }
124026e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
124126e093fcSHong Zhang     ii   = a->compressedrow.i;
12427b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
124326e093fcSHong Zhang   } else {
12442d61bbb3SSatish Balay     ii  = a->i;
124526e093fcSHong Zhang     y   = yarray;
124626e093fcSHong Zhang     z   = zarray;
124726e093fcSHong Zhang   }
12482d61bbb3SSatish Balay 
12492d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
12502d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
125126e093fcSHong Zhang     if (usecprow){
12527b2bb3b9SHong Zhang       z = zarray + 4*ridx[i];
12537b2bb3b9SHong Zhang       y = yarray + 4*ridx[i];
125426e093fcSHong Zhang     }
12552d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3];
12562d61bbb3SSatish Balay     for (j=0; j<n; j++) {
12572d61bbb3SSatish Balay       xb = x + 4*(*idx++);
12582d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
12592d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[4]*x2 + v[8]*x3   + v[12]*x4;
12602d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[5]*x2 + v[9]*x3   + v[13]*x4;
12612d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[6]*x2 + v[10]*x3  + v[14]*x4;
12622d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[7]*x2 + v[11]*x3  + v[15]*x4;
12632d61bbb3SSatish Balay       v += 16;
12642d61bbb3SSatish Balay     }
12652d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4;
126626e093fcSHong Zhang     if (!usecprow){
12672d61bbb3SSatish Balay       z += 4; y += 4;
12682d61bbb3SSatish Balay     }
126926e093fcSHong Zhang   }
12701ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
127126e093fcSHong Zhang   ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr);
12722e8a6d31SBarry Smith   if (zz != yy) {
127326e093fcSHong Zhang     ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
12742e8a6d31SBarry Smith   }
1275dc0b31edSSatish Balay   ierr = PetscLogFlops(32.0*a->nz);CHKERRQ(ierr);
12762d61bbb3SSatish Balay   PetscFunctionReturn(0);
12772d61bbb3SSatish Balay }
12782d61bbb3SSatish Balay 
12794a2ae208SSatish Balay #undef __FUNCT__
12804a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_5"
1281dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_5(Mat A,Vec xx,Vec yy,Vec zz)
12822d61bbb3SSatish Balay {
12832d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
1284dcf5cc72SBarry Smith   PetscScalar    *x,*y = 0,*z = 0,*xb,sum1,sum2,sum3,sum4,sum5,x1,x2,x3,x4,x5;
128526e093fcSHong Zhang   PetscScalar    *yarray,*zarray;
12863f1db9ecSBarry Smith   MatScalar      *v;
1287dfbe8321SBarry Smith   PetscErrorCode ierr;
12884eb6d288SHong Zhang   PetscInt       mbs=a->mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL;
128926e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
12902d61bbb3SSatish Balay 
12912d61bbb3SSatish Balay   PetscFunctionBegin;
12921ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
129326e093fcSHong Zhang   ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr);
12942e8a6d31SBarry Smith   if (zz != yy) {
129526e093fcSHong Zhang     ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
12962e8a6d31SBarry Smith   } else {
129726e093fcSHong Zhang     zarray = yarray;
12982e8a6d31SBarry Smith   }
12992d61bbb3SSatish Balay 
13002d61bbb3SSatish Balay   idx = a->j;
13012d61bbb3SSatish Balay   v   = a->a;
130226e093fcSHong Zhang   if (usecprow){
13034eb6d288SHong Zhang     if (zz != yy){
13044eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,5*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
13054eb6d288SHong Zhang     }
130626e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
130726e093fcSHong Zhang     ii   = a->compressedrow.i;
13087b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
130926e093fcSHong Zhang   } else {
13102d61bbb3SSatish Balay     ii  = a->i;
131126e093fcSHong Zhang     y   = yarray;
131226e093fcSHong Zhang     z   = zarray;
131326e093fcSHong Zhang   }
13142d61bbb3SSatish Balay 
13152d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
13162d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
131726e093fcSHong Zhang     if (usecprow){
13187b2bb3b9SHong Zhang       z = zarray + 5*ridx[i];
13197b2bb3b9SHong Zhang       y = yarray + 5*ridx[i];
132026e093fcSHong Zhang     }
13212d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4];
13222d61bbb3SSatish Balay     for (j=0; j<n; j++) {
13232d61bbb3SSatish Balay       xb = x + 5*(*idx++);
13242d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4];
13252d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[5]*x2 + v[10]*x3  + v[15]*x4 + v[20]*x5;
13262d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[6]*x2 + v[11]*x3  + v[16]*x4 + v[21]*x5;
13272d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[7]*x2 + v[12]*x3  + v[17]*x4 + v[22]*x5;
13282d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[8]*x2 + v[13]*x3  + v[18]*x4 + v[23]*x5;
13292d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[9]*x2 + v[14]*x3  + v[19]*x4 + v[24]*x5;
13302d61bbb3SSatish Balay       v += 25;
13312d61bbb3SSatish Balay     }
13322d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5;
133326e093fcSHong Zhang     if (!usecprow){
13342d61bbb3SSatish Balay       z += 5; y += 5;
13352d61bbb3SSatish Balay     }
133626e093fcSHong Zhang   }
13371ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
133826e093fcSHong Zhang   ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr);
13392e8a6d31SBarry Smith   if (zz != yy) {
134026e093fcSHong Zhang     ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
13412e8a6d31SBarry Smith   }
1342dc0b31edSSatish Balay   ierr = PetscLogFlops(50.0*a->nz);CHKERRQ(ierr);
13432d61bbb3SSatish Balay   PetscFunctionReturn(0);
13442d61bbb3SSatish Balay }
13454a2ae208SSatish Balay #undef __FUNCT__
13464a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_6"
1347dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_6(Mat A,Vec xx,Vec yy,Vec zz)
134815091d37SBarry Smith {
134915091d37SBarry Smith   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
1350dcf5cc72SBarry Smith   PetscScalar    *x,*y = 0,*z = 0,*xb,sum1,sum2,sum3,sum4,sum5,sum6;
135126e093fcSHong Zhang   PetscScalar    x1,x2,x3,x4,x5,x6,*yarray,*zarray;
135215091d37SBarry Smith   MatScalar      *v;
1353dfbe8321SBarry Smith   PetscErrorCode ierr;
13544eb6d288SHong Zhang   PetscInt       mbs=a->mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL;
135526e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
135615091d37SBarry Smith 
135715091d37SBarry Smith   PetscFunctionBegin;
13581ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
135926e093fcSHong Zhang   ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr);
136015091d37SBarry Smith   if (zz != yy) {
136126e093fcSHong Zhang     ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
136215091d37SBarry Smith   } else {
136326e093fcSHong Zhang     zarray = yarray;
136415091d37SBarry Smith   }
136515091d37SBarry Smith 
136615091d37SBarry Smith   idx = a->j;
136715091d37SBarry Smith   v   = a->a;
136826e093fcSHong Zhang   if (usecprow){
13694eb6d288SHong Zhang     if (zz != yy){
13704eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,6*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
13714eb6d288SHong Zhang     }
137226e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
137326e093fcSHong Zhang     ii   = a->compressedrow.i;
13747b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
137526e093fcSHong Zhang   } else {
137615091d37SBarry Smith     ii  = a->i;
137726e093fcSHong Zhang     y   = yarray;
137826e093fcSHong Zhang     z   = zarray;
137926e093fcSHong Zhang   }
138015091d37SBarry Smith 
138115091d37SBarry Smith   for (i=0; i<mbs; i++) {
138215091d37SBarry Smith     n  = ii[1] - ii[0]; ii++;
138326e093fcSHong Zhang     if (usecprow){
13847b2bb3b9SHong Zhang       z = zarray + 6*ridx[i];
13857b2bb3b9SHong Zhang       y = yarray + 6*ridx[i];
138626e093fcSHong Zhang     }
138715091d37SBarry Smith     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4]; sum6 = y[5];
138815091d37SBarry Smith     for (j=0; j<n; j++) {
13893b95cb0eSSatish Balay       xb = x + 6*(*idx++);
139015091d37SBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5];
139115091d37SBarry Smith       sum1 += v[0]*x1 + v[6]*x2  + v[12]*x3  + v[18]*x4 + v[24]*x5 + v[30]*x6;
139215091d37SBarry Smith       sum2 += v[1]*x1 + v[7]*x2  + v[13]*x3  + v[19]*x4 + v[25]*x5 + v[31]*x6;
139315091d37SBarry Smith       sum3 += v[2]*x1 + v[8]*x2  + v[14]*x3  + v[20]*x4 + v[26]*x5 + v[32]*x6;
139415091d37SBarry Smith       sum4 += v[3]*x1 + v[9]*x2  + v[15]*x3  + v[21]*x4 + v[27]*x5 + v[33]*x6;
139515091d37SBarry Smith       sum5 += v[4]*x1 + v[10]*x2 + v[16]*x3  + v[22]*x4 + v[28]*x5 + v[34]*x6;
139615091d37SBarry Smith       sum6 += v[5]*x1 + v[11]*x2 + v[17]*x3  + v[23]*x4 + v[29]*x5 + v[35]*x6;
139715091d37SBarry Smith       v += 36;
139815091d37SBarry Smith     }
139915091d37SBarry Smith     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6;
140026e093fcSHong Zhang     if (!usecprow){
140115091d37SBarry Smith       z += 6; y += 6;
140215091d37SBarry Smith     }
140326e093fcSHong Zhang   }
14041ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
140526e093fcSHong Zhang   ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr);
140615091d37SBarry Smith   if (zz != yy) {
140726e093fcSHong Zhang     ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
140815091d37SBarry Smith   }
1409dc0b31edSSatish Balay   ierr = PetscLogFlops(72.0*a->nz);CHKERRQ(ierr);
141015091d37SBarry Smith   PetscFunctionReturn(0);
141115091d37SBarry Smith }
14122d61bbb3SSatish Balay 
14134a2ae208SSatish Balay #undef __FUNCT__
14144a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_7"
1415dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_7(Mat A,Vec xx,Vec yy,Vec zz)
14162d61bbb3SSatish Balay {
14172d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
1418dcf5cc72SBarry Smith   PetscScalar    *x,*y = 0,*z = 0,*xb,sum1,sum2,sum3,sum4,sum5,sum6,sum7;
141926e093fcSHong Zhang   PetscScalar    x1,x2,x3,x4,x5,x6,x7,*yarray,*zarray;
14203f1db9ecSBarry Smith   MatScalar      *v;
1421dfbe8321SBarry Smith   PetscErrorCode ierr;
14227b2bb3b9SHong Zhang   PetscInt       mbs=a->mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL;
142326e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
14242d61bbb3SSatish Balay 
14252d61bbb3SSatish Balay   PetscFunctionBegin;
14261ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
142726e093fcSHong Zhang   ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr);
14282e8a6d31SBarry Smith   if (zz != yy) {
142926e093fcSHong Zhang     ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
14302e8a6d31SBarry Smith   } else {
143126e093fcSHong Zhang     zarray = yarray;
14322e8a6d31SBarry Smith   }
14332d61bbb3SSatish Balay 
14342d61bbb3SSatish Balay   idx = a->j;
14352d61bbb3SSatish Balay   v   = a->a;
143626e093fcSHong Zhang   if (usecprow){
14374eb6d288SHong Zhang     if (zz != yy){
14384eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,7*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
14394eb6d288SHong Zhang     }
144026e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
144126e093fcSHong Zhang     ii   = a->compressedrow.i;
14427b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
144326e093fcSHong Zhang   } else {
14442d61bbb3SSatish Balay     ii  = a->i;
144526e093fcSHong Zhang     y   = yarray;
144626e093fcSHong Zhang     z   = zarray;
144726e093fcSHong Zhang   }
14482d61bbb3SSatish Balay 
14492d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
14502d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
145126e093fcSHong Zhang     if (usecprow){
14527b2bb3b9SHong Zhang       z = zarray + 7*ridx[i];
14537b2bb3b9SHong Zhang       y = yarray + 7*ridx[i];
145426e093fcSHong Zhang     }
14552d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4]; sum6 = y[5]; sum7 = y[6];
14562d61bbb3SSatish Balay     for (j=0; j<n; j++) {
14572d61bbb3SSatish Balay       xb = x + 7*(*idx++);
14582d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6];
14592d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[7]*x2  + v[14]*x3  + v[21]*x4 + v[28]*x5 + v[35]*x6 + v[42]*x7;
14602d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[8]*x2  + v[15]*x3  + v[22]*x4 + v[29]*x5 + v[36]*x6 + v[43]*x7;
14612d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[9]*x2  + v[16]*x3  + v[23]*x4 + v[30]*x5 + v[37]*x6 + v[44]*x7;
14622d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[10]*x2 + v[17]*x3  + v[24]*x4 + v[31]*x5 + v[38]*x6 + v[45]*x7;
14632d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[11]*x2 + v[18]*x3  + v[25]*x4 + v[32]*x5 + v[39]*x6 + v[46]*x7;
14642d61bbb3SSatish Balay       sum6 += v[5]*x1 + v[12]*x2 + v[19]*x3  + v[26]*x4 + v[33]*x5 + v[40]*x6 + v[47]*x7;
14652d61bbb3SSatish Balay       sum7 += v[6]*x1 + v[13]*x2 + v[20]*x3  + v[27]*x4 + v[34]*x5 + v[41]*x6 + v[48]*x7;
14662d61bbb3SSatish Balay       v += 49;
14672d61bbb3SSatish Balay     }
14682d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
146926e093fcSHong Zhang     if (!usecprow){
14702d61bbb3SSatish Balay       z += 7; y += 7;
14712d61bbb3SSatish Balay     }
147226e093fcSHong Zhang   }
14731ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
147426e093fcSHong Zhang   ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr);
14752e8a6d31SBarry Smith   if (zz != yy) {
147626e093fcSHong Zhang     ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
14772e8a6d31SBarry Smith   }
1478dc0b31edSSatish Balay   ierr = PetscLogFlops(98.0*a->nz);CHKERRQ(ierr);
14792d61bbb3SSatish Balay   PetscFunctionReturn(0);
14802d61bbb3SSatish Balay }
1481218c64b6SSatish Balay 
14824a2ae208SSatish Balay #undef __FUNCT__
14834a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_N"
1484dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_N(Mat A,Vec xx,Vec yy,Vec zz)
14852d61bbb3SSatish Balay {
14862d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
1487bba805e6SBarry Smith   PetscScalar    *x,*z = 0,*xb,*work,*workt,*zarray;
14883f1db9ecSBarry Smith   MatScalar      *v;
14896849ba73SBarry Smith   PetscErrorCode ierr;
1490d0f46423SBarry Smith   PetscInt       mbs,i,*idx,*ii,bs=A->rmap->bs,j,n,bs2=a->bs2;
14917b2bb3b9SHong Zhang   PetscInt       ncols,k,*ridx=PETSC_NULL;
149226e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
1493218c64b6SSatish Balay 
14942d61bbb3SSatish Balay   PetscFunctionBegin;
14956a65c766SHong Zhang   ierr = VecCopy_Seq(yy,zz);CHKERRQ(ierr);
14961ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
149726e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
14982d61bbb3SSatish Balay 
14992d61bbb3SSatish Balay   idx = a->j;
15002d61bbb3SSatish Balay   v   = a->a;
150126e093fcSHong Zhang   if (usecprow){
150226e093fcSHong Zhang     mbs    = a->compressedrow.nrows;
150326e093fcSHong Zhang     ii     = a->compressedrow.i;
15047b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
150526e093fcSHong Zhang   } else {
150626e093fcSHong Zhang     mbs = a->mbs;
15072d61bbb3SSatish Balay     ii  = a->i;
150826e093fcSHong Zhang     z   = zarray;
150926e093fcSHong Zhang   }
15102d61bbb3SSatish Balay 
15112d61bbb3SSatish Balay   if (!a->mult_work) {
1512d0f46423SBarry Smith     k    = PetscMax(A->rmap->n,A->cmap->n);
151387828ca2SBarry Smith     ierr = PetscMalloc((k+1)*sizeof(PetscScalar),&a->mult_work);CHKERRQ(ierr);
15142d61bbb3SSatish Balay   }
15152d61bbb3SSatish Balay   work = a->mult_work;
15162d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
15172d61bbb3SSatish Balay     n     = ii[1] - ii[0]; ii++;
15182d61bbb3SSatish Balay     ncols = n*bs;
15192d61bbb3SSatish Balay     workt = work;
15202d61bbb3SSatish Balay     for (j=0; j<n; j++) {
15212d61bbb3SSatish Balay       xb = x + bs*(*idx++);
15222d61bbb3SSatish Balay       for (k=0; k<bs; k++) workt[k] = xb[k];
15232d61bbb3SSatish Balay       workt += bs;
15242d61bbb3SSatish Balay     }
15257b2bb3b9SHong Zhang     if (usecprow) z = zarray + bs*ridx[i];
15263f1db9ecSBarry Smith     Kernel_w_gets_w_plus_Ar_times_v(bs,ncols,work,v,z);
152771044d3cSBarry Smith     /* BLASgemv_("N",&bs,&ncols,&_DOne,v,&bs,work,&_One,&_DOne,z,&_One); */
15282d61bbb3SSatish Balay     v += n*bs2;
152926e093fcSHong Zhang     if (!usecprow){
15302d61bbb3SSatish Balay       z += bs;
15312d61bbb3SSatish Balay     }
153226e093fcSHong Zhang   }
15331ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
153426e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
1535dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*a->nz*bs2);CHKERRQ(ierr);
15362d61bbb3SSatish Balay   PetscFunctionReturn(0);
15372d61bbb3SSatish Balay }
15382d61bbb3SSatish Balay 
15394a2ae208SSatish Balay #undef __FUNCT__
1540547795f9SHong Zhang #define __FUNCT__ "MatMultHermitianTranspose_SeqBAIJ"
1541547795f9SHong Zhang PetscErrorCode MatMultHermitianTranspose_SeqBAIJ(Mat A,Vec xx,Vec zz)
1542547795f9SHong Zhang {
1543547795f9SHong Zhang   PetscScalar    zero = 0.0;
1544547795f9SHong Zhang   PetscErrorCode ierr;
1545547795f9SHong Zhang 
1546547795f9SHong Zhang   PetscFunctionBegin;
1547547795f9SHong Zhang   ierr = VecSet(zz,zero);CHKERRQ(ierr);
1548547795f9SHong Zhang   ierr = MatMultHermitianTransposeAdd_SeqBAIJ(A,xx,zz,zz);CHKERRQ(ierr);
1549547795f9SHong Zhang   PetscFunctionReturn(0);
1550547795f9SHong Zhang }
1551547795f9SHong Zhang 
1552547795f9SHong Zhang #undef __FUNCT__
15534a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_SeqBAIJ"
1554dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_SeqBAIJ(Mat A,Vec xx,Vec zz)
15552d61bbb3SSatish Balay {
15563447b6efSHong Zhang   PetscScalar    zero = 0.0;
15576849ba73SBarry Smith   PetscErrorCode ierr;
15582d61bbb3SSatish Balay 
15592d61bbb3SSatish Balay   PetscFunctionBegin;
15602dcb1b2aSMatthew Knepley   ierr = VecSet(zz,zero);CHKERRQ(ierr);
15613447b6efSHong Zhang   ierr = MatMultTransposeAdd_SeqBAIJ(A,xx,zz,zz);CHKERRQ(ierr);
15622d61bbb3SSatish Balay   PetscFunctionReturn(0);
15632d61bbb3SSatish Balay }
15642d61bbb3SSatish Balay 
15654a2ae208SSatish Balay #undef __FUNCT__
1566547795f9SHong Zhang #define __FUNCT__ "MatMultHermitianTransposeAdd_SeqBAIJ"
1567547795f9SHong Zhang PetscErrorCode MatMultHermitianTransposeAdd_SeqBAIJ(Mat A,Vec xx,Vec yy,Vec zz)
1568547795f9SHong Zhang 
1569547795f9SHong Zhang {
1570547795f9SHong Zhang   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
1571547795f9SHong Zhang   PetscScalar       *zb,*x,*z,*xb = 0,x1,x2,x3,x4,x5;
1572547795f9SHong Zhang   MatScalar         *v;
1573547795f9SHong Zhang   PetscErrorCode    ierr;
1574547795f9SHong Zhang   PetscInt          mbs,i,*idx,*ii,rval,bs=A->rmap->bs,j,n,bs2=a->bs2,*ib,*ridx=PETSC_NULL;
1575547795f9SHong Zhang   Mat_CompressedRow cprow = a->compressedrow;
1576547795f9SHong Zhang   PetscTruth        usecprow=cprow.use;
1577547795f9SHong Zhang 
1578547795f9SHong Zhang   PetscFunctionBegin;
1579547795f9SHong Zhang   if (yy != zz) { ierr = VecCopy_Seq(yy,zz);CHKERRQ(ierr); }
1580547795f9SHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1581547795f9SHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
1582547795f9SHong Zhang 
1583547795f9SHong Zhang   idx = a->j;
1584547795f9SHong Zhang   v   = a->a;
1585547795f9SHong Zhang   if (usecprow){
1586547795f9SHong Zhang     mbs  = cprow.nrows;
1587547795f9SHong Zhang     ii   = cprow.i;
1588547795f9SHong Zhang     ridx = cprow.rindex;
1589547795f9SHong Zhang   } else {
1590547795f9SHong Zhang     mbs=a->mbs;
1591547795f9SHong Zhang     ii = a->i;
1592547795f9SHong Zhang     xb = x;
1593547795f9SHong Zhang   }
1594547795f9SHong Zhang 
1595547795f9SHong Zhang   switch (bs) {
1596547795f9SHong Zhang   case 1:
1597547795f9SHong Zhang     for (i=0; i<mbs; i++) {
1598547795f9SHong Zhang       if (usecprow) xb = x + ridx[i];
1599547795f9SHong Zhang       x1 = xb[0];
1600547795f9SHong Zhang       ib = idx + ii[0];
1601547795f9SHong Zhang       n  = ii[1] - ii[0]; ii++;
1602547795f9SHong Zhang       for (j=0; j<n; j++) {
1603547795f9SHong Zhang         rval    = ib[j];
1604547795f9SHong Zhang         z[rval] += PetscConj(*v) * x1;
1605547795f9SHong Zhang         v++;
1606547795f9SHong Zhang       }
1607547795f9SHong Zhang       if (!usecprow) xb++;
1608547795f9SHong Zhang     }
1609547795f9SHong Zhang     break;
1610547795f9SHong Zhang   case 2:
1611547795f9SHong Zhang     for (i=0; i<mbs; i++) {
1612547795f9SHong Zhang       if (usecprow) xb = x + 2*ridx[i];
1613547795f9SHong Zhang       x1 = xb[0]; x2 = xb[1];
1614547795f9SHong Zhang       ib = idx + ii[0];
1615547795f9SHong Zhang       n  = ii[1] - ii[0]; ii++;
1616547795f9SHong Zhang       for (j=0; j<n; j++) {
1617547795f9SHong Zhang         rval      = ib[j]*2;
1618547795f9SHong Zhang         z[rval++] += PetscConj(v[0])*x1 + PetscConj(v[1])*x2;
1619547795f9SHong Zhang         z[rval++] += PetscConj(v[2])*x1 + PetscConj(v[3])*x2;
1620547795f9SHong Zhang         v  += 4;
1621547795f9SHong Zhang       }
1622547795f9SHong Zhang       if (!usecprow) xb += 2;
1623547795f9SHong Zhang     }
1624547795f9SHong Zhang     break;
1625547795f9SHong Zhang   case 3:
1626547795f9SHong Zhang     for (i=0; i<mbs; i++) {
1627547795f9SHong Zhang       if (usecprow) xb = x + 3*ridx[i];
1628547795f9SHong Zhang       x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
1629547795f9SHong Zhang       ib = idx + ii[0];
1630547795f9SHong Zhang       n  = ii[1] - ii[0]; ii++;
1631547795f9SHong Zhang       for (j=0; j<n; j++) {
1632547795f9SHong Zhang         rval      = ib[j]*3;
1633547795f9SHong Zhang         z[rval++] += PetscConj(v[0])*x1 + PetscConj(v[1])*x2 + PetscConj(v[2])*x3;
1634547795f9SHong Zhang         z[rval++] += PetscConj(v[3])*x1 + PetscConj(v[4])*x2 + PetscConj(v[5])*x3;
1635547795f9SHong Zhang         z[rval++] += PetscConj(v[6])*x1 + PetscConj(v[7])*x2 + PetscConj(v[8])*x3;
1636547795f9SHong Zhang         v  += 9;
1637547795f9SHong Zhang       }
1638547795f9SHong Zhang       if (!usecprow) xb += 3;
1639547795f9SHong Zhang     }
1640547795f9SHong Zhang     break;
1641547795f9SHong Zhang   case 4:
1642547795f9SHong Zhang     for (i=0; i<mbs; i++) {
1643547795f9SHong Zhang       if (usecprow) xb = x + 4*ridx[i];
1644547795f9SHong Zhang       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
1645547795f9SHong Zhang       ib = idx + ii[0];
1646547795f9SHong Zhang       n  = ii[1] - ii[0]; ii++;
1647547795f9SHong Zhang       for (j=0; j<n; j++) {
1648547795f9SHong Zhang         rval      = ib[j]*4;
1649547795f9SHong Zhang         z[rval++] +=  PetscConj(v[0])*x1 + PetscConj(v[1])*x2  + PetscConj(v[2])*x3  + PetscConj(v[3])*x4;
1650547795f9SHong Zhang         z[rval++] +=  PetscConj(v[4])*x1 + PetscConj(v[5])*x2  + PetscConj(v[6])*x3  + PetscConj(v[7])*x4;
1651547795f9SHong Zhang         z[rval++] +=  PetscConj(v[8])*x1 + PetscConj(v[9])*x2  + PetscConj(v[10])*x3 + PetscConj(v[11])*x4;
1652547795f9SHong Zhang         z[rval++] += PetscConj(v[12])*x1 + PetscConj(v[13])*x2 + PetscConj(v[14])*x3 + PetscConj(v[15])*x4;
1653547795f9SHong Zhang         v  += 16;
1654547795f9SHong Zhang       }
1655547795f9SHong Zhang       if (!usecprow) xb += 4;
1656547795f9SHong Zhang     }
1657547795f9SHong Zhang     break;
1658547795f9SHong Zhang   case 5:
1659547795f9SHong Zhang     for (i=0; i<mbs; i++) {
1660547795f9SHong Zhang       if (usecprow) xb = x + 5*ridx[i];
1661547795f9SHong Zhang       x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
1662547795f9SHong Zhang       x4 = xb[3]; x5 = xb[4];
1663547795f9SHong Zhang       ib = idx + ii[0];
1664547795f9SHong Zhang       n  = ii[1] - ii[0]; ii++;
1665547795f9SHong Zhang       for (j=0; j<n; j++) {
1666547795f9SHong Zhang         rval      = ib[j]*5;
1667547795f9SHong Zhang         z[rval++] +=  PetscConj(v[0])*x1 +  PetscConj(v[1])*x2 +  PetscConj(v[2])*x3 +  PetscConj(v[3])*x4 +  PetscConj(v[4])*x5;
1668547795f9SHong Zhang         z[rval++] +=  PetscConj(v[5])*x1 +  PetscConj(v[6])*x2 +  PetscConj(v[7])*x3 +  PetscConj(v[8])*x4 +  PetscConj(v[9])*x5;
1669547795f9SHong Zhang         z[rval++] += PetscConj(v[10])*x1 + PetscConj(v[11])*x2 + PetscConj(v[12])*x3 + PetscConj(v[13])*x4 + PetscConj(v[14])*x5;
1670547795f9SHong Zhang         z[rval++] += PetscConj(v[15])*x1 + PetscConj(v[16])*x2 + PetscConj(v[17])*x3 + PetscConj(v[18])*x4 + PetscConj(v[19])*x5;
1671547795f9SHong Zhang         z[rval++] += PetscConj(v[20])*x1 + PetscConj(v[21])*x2 + PetscConj(v[22])*x3 + PetscConj(v[23])*x4 + PetscConj(v[24])*x5;
1672547795f9SHong Zhang         v  += 25;
1673547795f9SHong Zhang       }
1674547795f9SHong Zhang       if (!usecprow) xb += 5;
1675547795f9SHong Zhang     }
1676547795f9SHong Zhang     break;
1677547795f9SHong Zhang   default: {      /* block sizes larger than 5 by 5 are handled by BLAS */
1678547795f9SHong Zhang       PetscInt     ncols,k;
1679547795f9SHong Zhang       PetscScalar  *work,*workt,*xtmp;
1680547795f9SHong Zhang 
1681e32f2f54SBarry Smith       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"block size larger than 5 is not supported yet");
1682547795f9SHong Zhang       if (!a->mult_work) {
1683547795f9SHong Zhang         k = PetscMax(A->rmap->n,A->cmap->n);
1684547795f9SHong Zhang         ierr = PetscMalloc((k+1)*sizeof(PetscScalar),&a->mult_work);CHKERRQ(ierr);
1685547795f9SHong Zhang       }
1686547795f9SHong Zhang       work = a->mult_work;
1687547795f9SHong Zhang       xtmp = x;
1688547795f9SHong Zhang       for (i=0; i<mbs; i++) {
1689547795f9SHong Zhang         n     = ii[1] - ii[0]; ii++;
1690547795f9SHong Zhang         ncols = n*bs;
1691547795f9SHong Zhang         ierr  = PetscMemzero(work,ncols*sizeof(PetscScalar));CHKERRQ(ierr);
1692547795f9SHong Zhang         if (usecprow) {
1693547795f9SHong Zhang           xtmp = x + bs*ridx[i];
1694547795f9SHong Zhang         }
1695547795f9SHong Zhang         Kernel_w_gets_w_plus_trans_Ar_times_v(bs,ncols,xtmp,v,work);
1696547795f9SHong Zhang         /* BLASgemv_("T",&bs,&ncols,&_DOne,v,&bs,xtmp,&_One,&_DOne,work,&_One); */
1697547795f9SHong Zhang         v += n*bs2;
1698547795f9SHong Zhang         if (!usecprow) xtmp += bs;
1699547795f9SHong Zhang         workt = work;
1700547795f9SHong Zhang         for (j=0; j<n; j++) {
1701547795f9SHong Zhang           zb = z + bs*(*idx++);
1702547795f9SHong Zhang           for (k=0; k<bs; k++) zb[k] += workt[k] ;
1703547795f9SHong Zhang           workt += bs;
1704547795f9SHong Zhang         }
1705547795f9SHong Zhang       }
1706547795f9SHong Zhang     }
1707547795f9SHong Zhang   }
1708547795f9SHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
1709547795f9SHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
1710547795f9SHong Zhang   ierr = PetscLogFlops(2.0*a->nz*a->bs2);CHKERRQ(ierr);
1711547795f9SHong Zhang   PetscFunctionReturn(0);
1712547795f9SHong Zhang }
1713547795f9SHong Zhang 
1714547795f9SHong Zhang #undef __FUNCT__
17154a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_SeqBAIJ"
1716dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_SeqBAIJ(Mat A,Vec xx,Vec yy,Vec zz)
17172d61bbb3SSatish Balay {
17182d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
1719dcf5cc72SBarry Smith   PetscScalar       *zb,*x,*z,*xb = 0,x1,x2,x3,x4,x5;
17203f1db9ecSBarry Smith   MatScalar         *v;
17216849ba73SBarry Smith   PetscErrorCode    ierr;
1722d0f46423SBarry Smith   PetscInt          mbs,i,*idx,*ii,rval,bs=A->rmap->bs,j,n,bs2=a->bs2,*ib,*ridx=PETSC_NULL;
17233447b6efSHong Zhang   Mat_CompressedRow cprow = a->compressedrow;
17243447b6efSHong Zhang   PetscTruth        usecprow=cprow.use;
17252d61bbb3SSatish Balay 
17262d61bbb3SSatish Balay   PetscFunctionBegin;
17276a65c766SHong Zhang   if (yy != zz) { ierr = VecCopy_Seq(yy,zz);CHKERRQ(ierr); }
17283447b6efSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
17293447b6efSHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
17302d61bbb3SSatish Balay 
17312d61bbb3SSatish Balay   idx = a->j;
17322d61bbb3SSatish Balay   v   = a->a;
17333447b6efSHong Zhang   if (usecprow){
17343447b6efSHong Zhang     mbs  = cprow.nrows;
17353447b6efSHong Zhang     ii   = cprow.i;
17367b2bb3b9SHong Zhang     ridx = cprow.rindex;
17373447b6efSHong Zhang   } else {
17383447b6efSHong Zhang     mbs=a->mbs;
17392d61bbb3SSatish Balay     ii = a->i;
1740f1af5d2fSBarry Smith     xb = x;
17413447b6efSHong Zhang   }
17422d61bbb3SSatish Balay 
17432d61bbb3SSatish Balay   switch (bs) {
17442d61bbb3SSatish Balay   case 1:
17452d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
17467b2bb3b9SHong Zhang       if (usecprow) xb = x + ridx[i];
1747f1af5d2fSBarry Smith       x1 = xb[0];
17483447b6efSHong Zhang       ib = idx + ii[0];
17493447b6efSHong Zhang       n  = ii[1] - ii[0]; ii++;
17502d61bbb3SSatish Balay       for (j=0; j<n; j++) {
17512d61bbb3SSatish Balay         rval    = ib[j];
1752f1af5d2fSBarry Smith         z[rval] += *v * x1;
1753f1af5d2fSBarry Smith         v++;
17542d61bbb3SSatish Balay       }
17553447b6efSHong Zhang       if (!usecprow) xb++;
17562d61bbb3SSatish Balay     }
17572d61bbb3SSatish Balay     break;
17582d61bbb3SSatish Balay   case 2:
17592d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
17607b2bb3b9SHong Zhang       if (usecprow) xb = x + 2*ridx[i];
1761f1af5d2fSBarry Smith       x1 = xb[0]; x2 = xb[1];
17623447b6efSHong Zhang       ib = idx + ii[0];
17633447b6efSHong Zhang       n  = ii[1] - ii[0]; ii++;
17642d61bbb3SSatish Balay       for (j=0; j<n; j++) {
17652d61bbb3SSatish Balay         rval      = ib[j]*2;
17662d61bbb3SSatish Balay         z[rval++] += v[0]*x1 + v[1]*x2;
17672d61bbb3SSatish Balay         z[rval++] += v[2]*x1 + v[3]*x2;
17682d61bbb3SSatish Balay         v  += 4;
17692d61bbb3SSatish Balay       }
17703447b6efSHong Zhang       if (!usecprow) xb += 2;
17712d61bbb3SSatish Balay     }
17722d61bbb3SSatish Balay     break;
17732d61bbb3SSatish Balay   case 3:
17742d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
17757b2bb3b9SHong Zhang       if (usecprow) xb = x + 3*ridx[i];
1776f1af5d2fSBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
17773447b6efSHong Zhang       ib = idx + ii[0];
17783447b6efSHong Zhang       n  = ii[1] - ii[0]; ii++;
17792d61bbb3SSatish Balay       for (j=0; j<n; j++) {
17802d61bbb3SSatish Balay         rval      = ib[j]*3;
17812d61bbb3SSatish Balay         z[rval++] += v[0]*x1 + v[1]*x2 + v[2]*x3;
17822d61bbb3SSatish Balay         z[rval++] += v[3]*x1 + v[4]*x2 + v[5]*x3;
17832d61bbb3SSatish Balay         z[rval++] += v[6]*x1 + v[7]*x2 + v[8]*x3;
17842d61bbb3SSatish Balay         v  += 9;
17852d61bbb3SSatish Balay       }
17863447b6efSHong Zhang       if (!usecprow) xb += 3;
17872d61bbb3SSatish Balay     }
17882d61bbb3SSatish Balay     break;
17892d61bbb3SSatish Balay   case 4:
17902d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
17917b2bb3b9SHong Zhang       if (usecprow) xb = x + 4*ridx[i];
1792f1af5d2fSBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
17933447b6efSHong Zhang       ib = idx + ii[0];
17943447b6efSHong Zhang       n  = ii[1] - ii[0]; ii++;
17952d61bbb3SSatish Balay       for (j=0; j<n; j++) {
17962d61bbb3SSatish Balay         rval      = ib[j]*4;
17972d61bbb3SSatish Balay         z[rval++] +=  v[0]*x1 +  v[1]*x2 +  v[2]*x3 +  v[3]*x4;
17982d61bbb3SSatish Balay         z[rval++] +=  v[4]*x1 +  v[5]*x2 +  v[6]*x3 +  v[7]*x4;
17992d61bbb3SSatish Balay         z[rval++] +=  v[8]*x1 +  v[9]*x2 + v[10]*x3 + v[11]*x4;
18002d61bbb3SSatish Balay         z[rval++] += v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4;
18012d61bbb3SSatish Balay         v  += 16;
18022d61bbb3SSatish Balay       }
18033447b6efSHong Zhang       if (!usecprow) xb += 4;
18042d61bbb3SSatish Balay     }
18052d61bbb3SSatish Balay     break;
18062d61bbb3SSatish Balay   case 5:
18072d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
18087b2bb3b9SHong Zhang       if (usecprow) xb = x + 5*ridx[i];
1809f1af5d2fSBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
18102d61bbb3SSatish Balay       x4 = xb[3]; x5 = xb[4];
18113447b6efSHong Zhang       ib = idx + ii[0];
18123447b6efSHong Zhang       n  = ii[1] - ii[0]; ii++;
18132d61bbb3SSatish Balay       for (j=0; j<n; j++) {
18142d61bbb3SSatish Balay         rval      = ib[j]*5;
18152d61bbb3SSatish Balay         z[rval++] +=  v[0]*x1 +  v[1]*x2 +  v[2]*x3 +  v[3]*x4 +  v[4]*x5;
18162d61bbb3SSatish Balay         z[rval++] +=  v[5]*x1 +  v[6]*x2 +  v[7]*x3 +  v[8]*x4 +  v[9]*x5;
18172d61bbb3SSatish Balay         z[rval++] += v[10]*x1 + v[11]*x2 + v[12]*x3 + v[13]*x4 + v[14]*x5;
18182d61bbb3SSatish Balay         z[rval++] += v[15]*x1 + v[16]*x2 + v[17]*x3 + v[18]*x4 + v[19]*x5;
18192d61bbb3SSatish Balay         z[rval++] += v[20]*x1 + v[21]*x2 + v[22]*x3 + v[23]*x4 + v[24]*x5;
18202d61bbb3SSatish Balay         v  += 25;
18212d61bbb3SSatish Balay       }
18223447b6efSHong Zhang       if (!usecprow) xb += 5;
18232d61bbb3SSatish Balay     }
18242d61bbb3SSatish Balay     break;
1825f1af5d2fSBarry Smith   default: {      /* block sizes larger then 5 by 5 are handled by BLAS */
1826690b6cddSBarry Smith       PetscInt     ncols,k;
18273447b6efSHong Zhang       PetscScalar  *work,*workt,*xtmp;
18283f1db9ecSBarry Smith 
18292d61bbb3SSatish Balay       if (!a->mult_work) {
1830d0f46423SBarry Smith         k = PetscMax(A->rmap->n,A->cmap->n);
183187828ca2SBarry Smith         ierr = PetscMalloc((k+1)*sizeof(PetscScalar),&a->mult_work);CHKERRQ(ierr);
18322d61bbb3SSatish Balay       }
18332d61bbb3SSatish Balay       work = a->mult_work;
18343447b6efSHong Zhang       xtmp = x;
18352d61bbb3SSatish Balay       for (i=0; i<mbs; i++) {
18362d61bbb3SSatish Balay         n     = ii[1] - ii[0]; ii++;
18372d61bbb3SSatish Balay         ncols = n*bs;
183887828ca2SBarry Smith         ierr  = PetscMemzero(work,ncols*sizeof(PetscScalar));CHKERRQ(ierr);
18393447b6efSHong Zhang         if (usecprow) {
18407b2bb3b9SHong Zhang           xtmp = x + bs*ridx[i];
18413447b6efSHong Zhang         }
18423447b6efSHong Zhang         Kernel_w_gets_w_plus_trans_Ar_times_v(bs,ncols,xtmp,v,work);
184371044d3cSBarry Smith         /* BLASgemv_("T",&bs,&ncols,&_DOne,v,&bs,xtmp,&_One,&_DOne,work,&_One); */
18442d61bbb3SSatish Balay         v += n*bs2;
18453447b6efSHong Zhang         if (!usecprow) xtmp += bs;
18462d61bbb3SSatish Balay         workt = work;
18472d61bbb3SSatish Balay         for (j=0; j<n; j++) {
18482d61bbb3SSatish Balay           zb = z + bs*(*idx++);
18492d61bbb3SSatish Balay           for (k=0; k<bs; k++) zb[k] += workt[k] ;
18502d61bbb3SSatish Balay           workt += bs;
18512d61bbb3SSatish Balay         }
18522d61bbb3SSatish Balay       }
18532d61bbb3SSatish Balay     }
18542d61bbb3SSatish Balay   }
18553447b6efSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
18563447b6efSHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
1857dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*a->nz*a->bs2);CHKERRQ(ierr);
18582d61bbb3SSatish Balay   PetscFunctionReturn(0);
18592d61bbb3SSatish Balay }
18602d61bbb3SSatish Balay 
18614a2ae208SSatish Balay #undef __FUNCT__
18624a2ae208SSatish Balay #define __FUNCT__ "MatScale_SeqBAIJ"
1863f4df32b1SMatthew Knepley PetscErrorCode MatScale_SeqBAIJ(Mat inA,PetscScalar alpha)
18642d61bbb3SSatish Balay {
18652d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)inA->data;
1866690b6cddSBarry Smith   PetscInt       totalnz = a->bs2*a->nz;
1867f4df32b1SMatthew Knepley   PetscScalar    oalpha = alpha;
1868efee365bSSatish Balay   PetscErrorCode ierr;
18690805154bSBarry Smith   PetscBLASInt   one = 1,tnz = PetscBLASIntCast(totalnz);
18702d61bbb3SSatish Balay 
18712d61bbb3SSatish Balay   PetscFunctionBegin;
1872f4df32b1SMatthew Knepley   BLASscal_(&tnz,&oalpha,a->a,&one);
1873efee365bSSatish Balay   ierr = PetscLogFlops(totalnz);CHKERRQ(ierr);
18742d61bbb3SSatish Balay   PetscFunctionReturn(0);
18752d61bbb3SSatish Balay }
18762d61bbb3SSatish Balay 
18774a2ae208SSatish Balay #undef __FUNCT__
18784a2ae208SSatish Balay #define __FUNCT__ "MatNorm_SeqBAIJ"
1879dfbe8321SBarry Smith PetscErrorCode MatNorm_SeqBAIJ(Mat A,NormType type,PetscReal *norm)
18802d61bbb3SSatish Balay {
18818a62d963SHong Zhang   PetscErrorCode ierr;
18822d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
18833f1db9ecSBarry Smith   MatScalar      *v = a->a;
1884329f5518SBarry Smith   PetscReal      sum = 0.0;
1885d0f46423SBarry Smith   PetscInt       i,j,k,bs=A->rmap->bs,nz=a->nz,bs2=a->bs2,k1;
18862d61bbb3SSatish Balay 
18872d61bbb3SSatish Balay   PetscFunctionBegin;
18882d61bbb3SSatish Balay   if (type == NORM_FROBENIUS) {
18892d61bbb3SSatish Balay     for (i=0; i< bs2*nz; i++) {
1890aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
1891329f5518SBarry Smith       sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
18922d61bbb3SSatish Balay #else
18932d61bbb3SSatish Balay       sum += (*v)*(*v); v++;
18942d61bbb3SSatish Balay #endif
18952d61bbb3SSatish Balay     }
18962d61bbb3SSatish Balay     *norm = sqrt(sum);
18978a62d963SHong Zhang   } else if (type == NORM_1) { /* maximum column sum */
18988a62d963SHong Zhang     PetscReal *tmp;
18998a62d963SHong Zhang     PetscInt  *bcol = a->j;
1900d0f46423SBarry Smith     ierr = PetscMalloc((A->cmap->n+1)*sizeof(PetscReal),&tmp);CHKERRQ(ierr);
1901d0f46423SBarry Smith     ierr = PetscMemzero(tmp,A->cmap->n*sizeof(PetscReal));CHKERRQ(ierr);
19028a62d963SHong Zhang     for (i=0; i<nz; i++){
19038a62d963SHong Zhang       for (j=0; j<bs; j++){
19048a62d963SHong Zhang         k1 = bs*(*bcol) + j; /* column index */
19058a62d963SHong Zhang         for (k=0; k<bs; k++){
19068a62d963SHong Zhang           tmp[k1] += PetscAbsScalar(*v); v++;
19078a62d963SHong Zhang         }
19088a62d963SHong Zhang       }
19098a62d963SHong Zhang       bcol++;
19108a62d963SHong Zhang     }
19118a62d963SHong Zhang     *norm = 0.0;
1912d0f46423SBarry Smith     for (j=0; j<A->cmap->n; j++) {
19138a62d963SHong Zhang       if (tmp[j] > *norm) *norm = tmp[j];
19148a62d963SHong Zhang     }
19158a62d963SHong Zhang     ierr = PetscFree(tmp);CHKERRQ(ierr);
1916596552b5SBarry Smith   } else if (type == NORM_INFINITY) { /* maximum row sum */
1917596552b5SBarry Smith     *norm = 0.0;
1918596552b5SBarry Smith     for (k=0; k<bs; k++) {
191974f84c7bSSatish Balay       for (j=0; j<a->mbs; j++) {
1920596552b5SBarry Smith         v = a->a + bs2*a->i[j] + k;
1921596552b5SBarry Smith         sum = 0.0;
1922596552b5SBarry Smith         for (i=0; i<a->i[j+1]-a->i[j]; i++) {
19230e90e235SBarry Smith           for (k1=0; k1<bs; k1++){
1924596552b5SBarry Smith             sum += PetscAbsScalar(*v);
1925596552b5SBarry Smith             v   += bs;
19262d61bbb3SSatish Balay           }
19270e90e235SBarry Smith         }
1928596552b5SBarry Smith         if (sum > *norm) *norm = sum;
1929596552b5SBarry Smith       }
1930596552b5SBarry Smith     }
1931e7e72b3dSBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for this norm yet");
19322d61bbb3SSatish Balay   PetscFunctionReturn(0);
19332d61bbb3SSatish Balay }
19342d61bbb3SSatish Balay 
19352d61bbb3SSatish Balay 
19364a2ae208SSatish Balay #undef __FUNCT__
19374a2ae208SSatish Balay #define __FUNCT__ "MatEqual_SeqBAIJ"
1938dfbe8321SBarry Smith PetscErrorCode MatEqual_SeqBAIJ(Mat A,Mat B,PetscTruth* flg)
19392d61bbb3SSatish Balay {
19402d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ *)A->data,*b = (Mat_SeqBAIJ *)B->data;
1941dfbe8321SBarry Smith   PetscErrorCode ierr;
19422d61bbb3SSatish Balay 
19432d61bbb3SSatish Balay   PetscFunctionBegin;
19442d61bbb3SSatish Balay   /* If the  matrix/block dimensions are not equal, or no of nonzeros or shift */
1945d0f46423SBarry Smith   if ((A->rmap->N != B->rmap->N) || (A->cmap->n != B->cmap->n) || (A->rmap->bs != B->rmap->bs)|| (a->nz != b->nz)) {
1946273d9f13SBarry Smith     *flg = PETSC_FALSE;
1947273d9f13SBarry Smith     PetscFunctionReturn(0);
19482d61bbb3SSatish Balay   }
19492d61bbb3SSatish Balay 
19502d61bbb3SSatish Balay   /* if the a->i are the same */
1951690b6cddSBarry Smith   ierr = PetscMemcmp(a->i,b->i,(a->mbs+1)*sizeof(PetscInt),flg);CHKERRQ(ierr);
1952abc0a331SBarry Smith   if (!*flg) {
19530f5bd95cSBarry Smith     PetscFunctionReturn(0);
19542d61bbb3SSatish Balay   }
19552d61bbb3SSatish Balay 
19562d61bbb3SSatish Balay   /* if a->j are the same */
1957690b6cddSBarry Smith   ierr = PetscMemcmp(a->j,b->j,(a->nz)*sizeof(PetscInt),flg);CHKERRQ(ierr);
1958abc0a331SBarry Smith   if (!*flg) {
19590f5bd95cSBarry Smith     PetscFunctionReturn(0);
19602d61bbb3SSatish Balay   }
19612d61bbb3SSatish Balay   /* if a->a are the same */
1962d0f46423SBarry Smith   ierr = PetscMemcmp(a->a,b->a,(a->nz)*(A->rmap->bs)*(B->rmap->bs)*sizeof(PetscScalar),flg);CHKERRQ(ierr);
19632d61bbb3SSatish Balay   PetscFunctionReturn(0);
19642d61bbb3SSatish Balay 
19652d61bbb3SSatish Balay }
19662d61bbb3SSatish Balay 
19674a2ae208SSatish Balay #undef __FUNCT__
19684a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_SeqBAIJ"
1969dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_SeqBAIJ(Mat A,Vec v)
19702d61bbb3SSatish Balay {
19712d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
1972dfbe8321SBarry Smith   PetscErrorCode ierr;
1973690b6cddSBarry Smith   PetscInt       i,j,k,n,row,bs,*ai,*aj,ambs,bs2;
197487828ca2SBarry Smith   PetscScalar    *x,zero = 0.0;
19753f1db9ecSBarry Smith   MatScalar      *aa,*aa_j;
19762d61bbb3SSatish Balay 
19772d61bbb3SSatish Balay   PetscFunctionBegin;
1978e32f2f54SBarry Smith   if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1979d0f46423SBarry Smith   bs   = A->rmap->bs;
19802d61bbb3SSatish Balay   aa   = a->a;
19812d61bbb3SSatish Balay   ai   = a->i;
19822d61bbb3SSatish Balay   aj   = a->j;
19832d61bbb3SSatish Balay   ambs = a->mbs;
19842d61bbb3SSatish Balay   bs2  = a->bs2;
19852d61bbb3SSatish Balay 
19862dcb1b2aSMatthew Knepley   ierr = VecSet(v,zero);CHKERRQ(ierr);
19871ebc52fbSHong Zhang   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1988e1311b90SBarry Smith   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
1989e32f2f54SBarry Smith   if (n != A->rmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
19902d61bbb3SSatish Balay   for (i=0; i<ambs; i++) {
19912d61bbb3SSatish Balay     for (j=ai[i]; j<ai[i+1]; j++) {
19922d61bbb3SSatish Balay       if (aj[j] == i) {
19932d61bbb3SSatish Balay         row  = i*bs;
19942d61bbb3SSatish Balay         aa_j = aa+j*bs2;
19952d61bbb3SSatish Balay         for (k=0; k<bs2; k+=(bs+1),row++) x[row] = aa_j[k];
19962d61bbb3SSatish Balay         break;
19972d61bbb3SSatish Balay       }
19982d61bbb3SSatish Balay     }
19992d61bbb3SSatish Balay   }
20001ebc52fbSHong Zhang   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
20012d61bbb3SSatish Balay   PetscFunctionReturn(0);
20022d61bbb3SSatish Balay }
20032d61bbb3SSatish Balay 
20044a2ae208SSatish Balay #undef __FUNCT__
20054a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_SeqBAIJ"
2006dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_SeqBAIJ(Mat A,Vec ll,Vec rr)
20072d61bbb3SSatish Balay {
20082d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
200987828ca2SBarry Smith   PetscScalar    *l,*r,x,*li,*ri;
20103f1db9ecSBarry Smith   MatScalar      *aa,*v;
2011dfbe8321SBarry Smith   PetscErrorCode ierr;
2012690b6cddSBarry Smith   PetscInt       i,j,k,lm,rn,M,m,n,*ai,*aj,mbs,tmp,bs,bs2;
20132d61bbb3SSatish Balay 
20142d61bbb3SSatish Balay   PetscFunctionBegin;
20152d61bbb3SSatish Balay   ai  = a->i;
20162d61bbb3SSatish Balay   aj  = a->j;
20172d61bbb3SSatish Balay   aa  = a->a;
2018d0f46423SBarry Smith   m   = A->rmap->n;
2019d0f46423SBarry Smith   n   = A->cmap->n;
2020d0f46423SBarry Smith   bs  = A->rmap->bs;
20212d61bbb3SSatish Balay   mbs = a->mbs;
20222d61bbb3SSatish Balay   bs2 = a->bs2;
20232d61bbb3SSatish Balay   if (ll) {
20241ebc52fbSHong Zhang     ierr = VecGetArray(ll,&l);CHKERRQ(ierr);
2025e1311b90SBarry Smith     ierr = VecGetLocalSize(ll,&lm);CHKERRQ(ierr);
2026e32f2f54SBarry Smith     if (lm != m) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Left scaling vector wrong length");
20272d61bbb3SSatish Balay     for (i=0; i<mbs; i++) { /* for each block row */
20282d61bbb3SSatish Balay       M  = ai[i+1] - ai[i];
20292d61bbb3SSatish Balay       li = l + i*bs;
20302d61bbb3SSatish Balay       v  = aa + bs2*ai[i];
20312d61bbb3SSatish Balay       for (j=0; j<M; j++) { /* for each block */
20322d61bbb3SSatish Balay         for (k=0; k<bs2; k++) {
20332d61bbb3SSatish Balay           (*v++) *= li[k%bs];
20342d61bbb3SSatish Balay         }
20352d61bbb3SSatish Balay       }
20362d61bbb3SSatish Balay     }
20371ebc52fbSHong Zhang     ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr);
2038efee365bSSatish Balay     ierr = PetscLogFlops(a->nz);CHKERRQ(ierr);
20392d61bbb3SSatish Balay   }
20402d61bbb3SSatish Balay 
20412d61bbb3SSatish Balay   if (rr) {
20421ebc52fbSHong Zhang     ierr = VecGetArray(rr,&r);CHKERRQ(ierr);
2043e1311b90SBarry Smith     ierr = VecGetLocalSize(rr,&rn);CHKERRQ(ierr);
2044e32f2f54SBarry Smith     if (rn != n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Right scaling vector wrong length");
20452d61bbb3SSatish Balay     for (i=0; i<mbs; i++) { /* for each block row */
20462d61bbb3SSatish Balay       M  = ai[i+1] - ai[i];
20472d61bbb3SSatish Balay       v  = aa + bs2*ai[i];
20482d61bbb3SSatish Balay       for (j=0; j<M; j++) { /* for each block */
20492d61bbb3SSatish Balay         ri = r + bs*aj[ai[i]+j];
20502d61bbb3SSatish Balay         for (k=0; k<bs; k++) {
20512d61bbb3SSatish Balay           x = ri[k];
20522d61bbb3SSatish Balay           for (tmp=0; tmp<bs; tmp++) (*v++) *= x;
20532d61bbb3SSatish Balay         }
20542d61bbb3SSatish Balay       }
20552d61bbb3SSatish Balay     }
20561ebc52fbSHong Zhang     ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr);
2057efee365bSSatish Balay     ierr = PetscLogFlops(a->nz);CHKERRQ(ierr);
20582d61bbb3SSatish Balay   }
20592d61bbb3SSatish Balay   PetscFunctionReturn(0);
20602d61bbb3SSatish Balay }
20612d61bbb3SSatish Balay 
20622d61bbb3SSatish Balay 
20634a2ae208SSatish Balay #undef __FUNCT__
20644a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_SeqBAIJ"
2065dfbe8321SBarry Smith PetscErrorCode MatGetInfo_SeqBAIJ(Mat A,MatInfoType flag,MatInfo *info)
20662d61bbb3SSatish Balay {
20672d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
20682d61bbb3SSatish Balay 
20692d61bbb3SSatish Balay   PetscFunctionBegin;
20702d61bbb3SSatish Balay   info->block_size     = a->bs2;
20712d61bbb3SSatish Balay   info->nz_allocated   = a->maxnz;
20722d61bbb3SSatish Balay   info->nz_used        = a->bs2*a->nz;
20732d61bbb3SSatish Balay   info->nz_unneeded    = (double)(info->nz_allocated - info->nz_used);
20742d61bbb3SSatish Balay   info->assemblies   = A->num_ass;
20758e58a170SBarry Smith   info->mallocs      = A->info.mallocs;
20767adad957SLisandro Dalcin   info->memory       = ((PetscObject)A)->mem;
2077d5f3da31SBarry Smith   if (A->factortype) {
20782d61bbb3SSatish Balay     info->fill_ratio_given  = A->info.fill_ratio_given;
20792d61bbb3SSatish Balay     info->fill_ratio_needed = A->info.fill_ratio_needed;
20802d61bbb3SSatish Balay     info->factor_mallocs    = A->info.factor_mallocs;
20812d61bbb3SSatish Balay   } else {
20822d61bbb3SSatish Balay     info->fill_ratio_given  = 0;
20832d61bbb3SSatish Balay     info->fill_ratio_needed = 0;
20842d61bbb3SSatish Balay     info->factor_mallocs    = 0;
20852d61bbb3SSatish Balay   }
20862d61bbb3SSatish Balay   PetscFunctionReturn(0);
20872d61bbb3SSatish Balay }
20882d61bbb3SSatish Balay 
20892d61bbb3SSatish Balay 
20904a2ae208SSatish Balay #undef __FUNCT__
20914a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_SeqBAIJ"
2092dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_SeqBAIJ(Mat A)
20932d61bbb3SSatish Balay {
20942d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
2095dfbe8321SBarry Smith   PetscErrorCode ierr;
20962d61bbb3SSatish Balay 
20972d61bbb3SSatish Balay   PetscFunctionBegin;
2098549d3d68SSatish Balay   ierr = PetscMemzero(a->a,a->bs2*a->i[a->mbs]*sizeof(MatScalar));CHKERRQ(ierr);
20992d61bbb3SSatish Balay   PetscFunctionReturn(0);
21002d61bbb3SSatish Balay }
2101