xref: /petsc/src/mat/impls/baij/seq/baij2.c (revision 8f1a2a5e609bf3d0ec1c7e7bf09507165b68d274)
1cac129eeSSatish Balay 
2c6db04a5SJed Brown #include <../src/mat/impls/baij/seq/baij.h>
3c6db04a5SJed Brown #include <../src/mat/blockinvert.h>
4c6db04a5SJed Brown #include <petscbt.h>
5c6db04a5SJed Brown #include <petscblaslapack.h>
6cac129eeSSatish Balay 
74a2ae208SSatish Balay #undef __FUNCT__
84a2ae208SSatish Balay #define __FUNCT__ "MatIncreaseOverlap_SeqBAIJ"
9690b6cddSBarry Smith PetscErrorCode MatIncreaseOverlap_SeqBAIJ(Mat A,PetscInt is_max,IS is[],PetscInt ov)
10a3192f15SSatish Balay {
11a3192f15SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
126849ba73SBarry Smith   PetscErrorCode ierr;
135d0c19d7SBarry Smith   PetscInt       row,i,j,k,l,m,n,*nidx,isz,val,ival;
145d0c19d7SBarry Smith   const PetscInt *idx;
15690b6cddSBarry Smith   PetscInt       start,end,*ai,*aj,bs,*nidx2;
16f1af5d2fSBarry Smith   PetscBT        table;
17a3192f15SSatish Balay 
183a40ed3dSBarry Smith   PetscFunctionBegin;
19a3192f15SSatish Balay   m     = a->mbs;
20a3192f15SSatish Balay   ai    = a->i;
21a3192f15SSatish Balay   aj    = a->j;
22d0f46423SBarry Smith   bs    = A->rmap->bs;
23a3192f15SSatish Balay 
24e32f2f54SBarry Smith   if (ov < 0)  SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap specified");
25a3192f15SSatish Balay 
266831982aSBarry Smith   ierr = PetscBTCreate(m,table);CHKERRQ(ierr);
27690b6cddSBarry Smith   ierr = PetscMalloc((m+1)*sizeof(PetscInt),&nidx);CHKERRQ(ierr);
28d0f46423SBarry Smith   ierr = PetscMalloc((A->rmap->N+1)*sizeof(PetscInt),&nidx2);CHKERRQ(ierr);
29a3192f15SSatish Balay 
30a3192f15SSatish Balay   for (i=0; i<is_max; i++) {
31a3192f15SSatish Balay     /* Initialise the two local arrays */
32a3192f15SSatish Balay     isz  = 0;
336831982aSBarry Smith     ierr = PetscBTMemzero(m,table);CHKERRQ(ierr);
34a3192f15SSatish Balay 
35a3192f15SSatish Balay     /* Extract the indices, assume there can be duplicate entries */
36a3192f15SSatish Balay     ierr = ISGetIndices(is[i],&idx);CHKERRQ(ierr);
37b9b97703SBarry Smith     ierr = ISGetLocalSize(is[i],&n);CHKERRQ(ierr);
38a3192f15SSatish Balay 
39a3192f15SSatish Balay     /* Enter these into the temp arrays i.e mark table[row], enter row into new index */
40a3192f15SSatish Balay     for (j=0; j<n ; ++j){
41218c64b6SSatish Balay       ival = idx[j]/bs; /* convert the indices into block indices */
42e32f2f54SBarry Smith       if (ival>=m) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"index greater than mat-dim");
43f1af5d2fSBarry Smith       if(!PetscBTLookupSet(table,ival)) { nidx[isz++] = ival;}
44a3192f15SSatish Balay     }
45a3192f15SSatish Balay     ierr = ISRestoreIndices(is[i],&idx);CHKERRQ(ierr);
466bf464f9SBarry Smith     ierr = ISDestroy(&is[i]);CHKERRQ(ierr);
47a3192f15SSatish Balay 
48a3192f15SSatish Balay     k = 0;
49a3192f15SSatish Balay     for (j=0; j<ov; j++){ /* for each overlap*/
50a3192f15SSatish Balay       n = isz;
51a3192f15SSatish Balay       for (; k<n ; k++){ /* do only those rows in nidx[k], which are not done yet */
52a3192f15SSatish Balay         row   = nidx[k];
53a3192f15SSatish Balay         start = ai[row];
54a3192f15SSatish Balay         end   = ai[row+1];
55a3192f15SSatish Balay         for (l = start; l<end ; l++){
56a3192f15SSatish Balay           val = aj[l];
57f1af5d2fSBarry Smith           if (!PetscBTLookupSet(table,val)) {nidx[isz++] = val;}
58a3192f15SSatish Balay         }
59a3192f15SSatish Balay       }
60a3192f15SSatish Balay     }
61218c64b6SSatish Balay     /* expand the Index Set */
62218c64b6SSatish Balay     for (j=0; j<isz; j++) {
63218c64b6SSatish Balay       for (k=0; k<bs; k++)
64218c64b6SSatish Balay         nidx2[j*bs+k] = nidx[j]*bs+k;
65218c64b6SSatish Balay     }
6670b3c8c7SBarry Smith     ierr = ISCreateGeneral(PETSC_COMM_SELF,isz*bs,nidx2,PETSC_COPY_VALUES,is+i);CHKERRQ(ierr);
67a3192f15SSatish Balay   }
686831982aSBarry Smith   ierr = PetscBTDestroy(table);CHKERRQ(ierr);
69606d414cSSatish Balay   ierr = PetscFree(nidx);CHKERRQ(ierr);
70606d414cSSatish Balay   ierr = PetscFree(nidx2);CHKERRQ(ierr);
713a40ed3dSBarry Smith   PetscFunctionReturn(0);
72a3192f15SSatish Balay }
731c351548SSatish Balay 
744a2ae208SSatish Balay #undef __FUNCT__
754a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqBAIJ_Private"
764aa3045dSJed Brown PetscErrorCode MatGetSubMatrix_SeqBAIJ_Private(Mat A,IS isrow,IS iscol,MatReuse scall,Mat *B)
77736121d4SSatish Balay {
78736121d4SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*c;
796849ba73SBarry Smith   PetscErrorCode ierr;
80690b6cddSBarry Smith   PetscInt       *smap,i,k,kstart,kend,oldcols = a->nbs,*lens;
81690b6cddSBarry Smith   PetscInt       row,mat_i,*mat_j,tcol,*mat_ilen;
825d0c19d7SBarry Smith   const PetscInt *irow,*icol;
835d0c19d7SBarry Smith   PetscInt       nrows,ncols,*ssmap,bs=A->rmap->bs,bs2=a->bs2;
84690b6cddSBarry Smith   PetscInt       *aj = a->j,*ai = a->i;
853f1db9ecSBarry Smith   MatScalar      *mat_a;
86736121d4SSatish Balay   Mat            C;
87ace3abfcSBarry Smith   PetscBool      flag,sorted;
88736121d4SSatish Balay 
893a40ed3dSBarry Smith   PetscFunctionBegin;
9014ca34e6SBarry Smith   ierr = ISSorted(iscol,&sorted);CHKERRQ(ierr);
91e32f2f54SBarry Smith   if (!sorted) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"IS is not sorted");
92736121d4SSatish Balay 
93736121d4SSatish Balay   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
94218c64b6SSatish Balay   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
95b9b97703SBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
96b9b97703SBarry Smith   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
97736121d4SSatish Balay 
98690b6cddSBarry Smith   ierr = PetscMalloc((1+oldcols)*sizeof(PetscInt),&smap);CHKERRQ(ierr);
99736121d4SSatish Balay   ssmap = smap;
100690b6cddSBarry Smith   ierr = PetscMalloc((1+nrows)*sizeof(PetscInt),&lens);CHKERRQ(ierr);
101690b6cddSBarry Smith   ierr  = PetscMemzero(smap,oldcols*sizeof(PetscInt));CHKERRQ(ierr);
102736121d4SSatish Balay   for (i=0; i<ncols; i++) smap[icol[i]] = i+1;
103736121d4SSatish Balay   /* determine lens of each row */
104736121d4SSatish Balay   for (i=0; i<nrows; i++) {
105736121d4SSatish Balay     kstart  = ai[irow[i]];
106736121d4SSatish Balay     kend    = kstart + a->ilen[irow[i]];
107736121d4SSatish Balay     lens[i] = 0;
108736121d4SSatish Balay       for (k=kstart; k<kend; k++) {
109736121d4SSatish Balay         if (ssmap[aj[k]]) {
110736121d4SSatish Balay           lens[i]++;
111736121d4SSatish Balay         }
112736121d4SSatish Balay       }
113736121d4SSatish Balay     }
114736121d4SSatish Balay   /* Create and fill new matrix */
115736121d4SSatish Balay   if (scall == MAT_REUSE_MATRIX) {
116736121d4SSatish Balay     c = (Mat_SeqBAIJ *)((*B)->data);
117736121d4SSatish Balay 
118e32f2f54SBarry Smith     if (c->mbs!=nrows || c->nbs!=ncols || (*B)->rmap->bs!=bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Submatrix wrong size");
119690b6cddSBarry Smith     ierr = PetscMemcmp(c->ilen,lens,c->mbs *sizeof(PetscInt),&flag);CHKERRQ(ierr);
120e7e72b3dSBarry Smith     if (!flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong no of nonzeros");
121690b6cddSBarry Smith     ierr = PetscMemzero(c->ilen,c->mbs*sizeof(PetscInt));CHKERRQ(ierr);
122736121d4SSatish Balay     C = *B;
1233a40ed3dSBarry Smith   } else {
1247adad957SLisandro Dalcin     ierr = MatCreate(((PetscObject)A)->comm,&C);CHKERRQ(ierr);
125f69a0ea3SMatthew Knepley     ierr = MatSetSizes(C,nrows*bs,ncols*bs,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1267adad957SLisandro Dalcin     ierr = MatSetType(C,((PetscObject)A)->type_name);CHKERRQ(ierr);
127ab93d7beSBarry Smith     ierr = MatSeqBAIJSetPreallocation_SeqBAIJ(C,bs,0,lens);CHKERRQ(ierr);
128736121d4SSatish Balay   }
129736121d4SSatish Balay   c = (Mat_SeqBAIJ *)(C->data);
130736121d4SSatish Balay   for (i=0; i<nrows; i++) {
131736121d4SSatish Balay     row    = irow[i];
132736121d4SSatish Balay     kstart = ai[row];
133736121d4SSatish Balay     kend   = kstart + a->ilen[row];
134736121d4SSatish Balay     mat_i  = c->i[i];
135736121d4SSatish Balay     mat_j  = c->j + mat_i;
136218c64b6SSatish Balay     mat_a  = c->a + mat_i*bs2;
137736121d4SSatish Balay     mat_ilen = c->ilen + i;
138736121d4SSatish Balay     for (k=kstart; k<kend; k++) {
139736121d4SSatish Balay       if ((tcol=ssmap[a->j[k]])) {
140736121d4SSatish Balay         *mat_j++ = tcol - 1;
141549d3d68SSatish Balay         ierr     = PetscMemcpy(mat_a,a->a+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
142549d3d68SSatish Balay         mat_a   += bs2;
143736121d4SSatish Balay         (*mat_ilen)++;
144736121d4SSatish Balay       }
145736121d4SSatish Balay     }
146736121d4SSatish Balay   }
147218c64b6SSatish Balay 
148736121d4SSatish Balay   /* Free work space */
149736121d4SSatish Balay   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
150606d414cSSatish Balay   ierr = PetscFree(smap);CHKERRQ(ierr);
151606d414cSSatish Balay   ierr = PetscFree(lens);CHKERRQ(ierr);
1526d4a8577SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1536d4a8577SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
154736121d4SSatish Balay 
155736121d4SSatish Balay   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
156736121d4SSatish Balay   *B = C;
1573a40ed3dSBarry Smith   PetscFunctionReturn(0);
158736121d4SSatish Balay }
159736121d4SSatish Balay 
1604a2ae208SSatish Balay #undef __FUNCT__
1614a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqBAIJ"
1624aa3045dSJed Brown PetscErrorCode MatGetSubMatrix_SeqBAIJ(Mat A,IS isrow,IS iscol,MatReuse scall,Mat *B)
163218c64b6SSatish Balay {
164218c64b6SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
165218c64b6SSatish Balay   IS             is1,is2;
1666849ba73SBarry Smith   PetscErrorCode ierr;
1675d0c19d7SBarry Smith   PetscInt       *vary,*iary,nrows,ncols,i,bs=A->rmap->bs,count;
1685d0c19d7SBarry Smith   const PetscInt *irow,*icol;
169218c64b6SSatish Balay 
1703a40ed3dSBarry Smith   PetscFunctionBegin;
171218c64b6SSatish Balay   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
172218c64b6SSatish Balay   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
173b9b97703SBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
174b9b97703SBarry Smith   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
175218c64b6SSatish Balay 
176218c64b6SSatish Balay   /* Verify if the indices corespond to each element in a block
177218c64b6SSatish Balay    and form the IS with compressed IS */
178fca92195SBarry Smith   ierr = PetscMalloc2(a->mbs,PetscInt,&vary,a->mbs,PetscInt,&iary);CHKERRQ(ierr);
179fca92195SBarry Smith   ierr = PetscMemzero(vary,a->mbs*sizeof(PetscInt));CHKERRQ(ierr);
180218c64b6SSatish Balay   for (i=0; i<nrows; i++) vary[irow[i]/bs]++;
181218c64b6SSatish Balay   count = 0;
182218c64b6SSatish Balay   for (i=0; i<a->mbs; i++) {
183e32f2f54SBarry Smith     if (vary[i]!=0 && vary[i]!=bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Index set does not match blocks");
184218c64b6SSatish Balay     if (vary[i]==bs) iary[count++] = i;
185218c64b6SSatish Balay   }
18670b3c8c7SBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,count,iary,PETSC_COPY_VALUES,&is1);CHKERRQ(ierr);
187218c64b6SSatish Balay 
188690b6cddSBarry Smith   ierr = PetscMemzero(vary,(a->mbs)*sizeof(PetscInt));CHKERRQ(ierr);
189218c64b6SSatish Balay   for (i=0; i<ncols; i++) vary[icol[i]/bs]++;
190218c64b6SSatish Balay   count = 0;
191218c64b6SSatish Balay   for (i=0; i<a->mbs; i++) {
192e32f2f54SBarry Smith     if (vary[i]!=0 && vary[i]!=bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Internal error in PETSc");
193218c64b6SSatish Balay     if (vary[i]==bs) iary[count++] = i;
194218c64b6SSatish Balay   }
19570b3c8c7SBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,count,iary,PETSC_COPY_VALUES,&is2);CHKERRQ(ierr);
196218c64b6SSatish Balay   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
197218c64b6SSatish Balay   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
198fca92195SBarry Smith   ierr = PetscFree2(vary,iary);CHKERRQ(ierr);
199218c64b6SSatish Balay 
2004aa3045dSJed Brown   ierr = MatGetSubMatrix_SeqBAIJ_Private(A,is1,is2,scall,B);CHKERRQ(ierr);
2016bf464f9SBarry Smith   ierr = ISDestroy(&is1);CHKERRQ(ierr);
2026bf464f9SBarry Smith   ierr = ISDestroy(&is2);CHKERRQ(ierr);
2033a40ed3dSBarry Smith   PetscFunctionReturn(0);
204218c64b6SSatish Balay }
205218c64b6SSatish Balay 
2064a2ae208SSatish Balay #undef __FUNCT__
2074a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_SeqBAIJ"
208690b6cddSBarry Smith PetscErrorCode MatGetSubMatrices_SeqBAIJ(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
209736121d4SSatish Balay {
2106849ba73SBarry Smith   PetscErrorCode ierr;
211690b6cddSBarry Smith   PetscInt       i;
212736121d4SSatish Balay 
2133a40ed3dSBarry Smith   PetscFunctionBegin;
214736121d4SSatish Balay   if (scall == MAT_INITIAL_MATRIX) {
21582502324SSatish Balay     ierr = PetscMalloc((n+1)*sizeof(Mat),B);CHKERRQ(ierr);
216736121d4SSatish Balay   }
217736121d4SSatish Balay 
218736121d4SSatish Balay   for (i=0; i<n; i++) {
2194aa3045dSJed Brown     ierr = MatGetSubMatrix_SeqBAIJ(A,irow[i],icol[i],scall,&(*B)[i]);CHKERRQ(ierr);
220736121d4SSatish Balay   }
2213a40ed3dSBarry Smith   PetscFunctionReturn(0);
222736121d4SSatish Balay }
223218c64b6SSatish Balay 
224218c64b6SSatish Balay 
2252d61bbb3SSatish Balay /* -------------------------------------------------------*/
2262d61bbb3SSatish Balay /* Should check that shapes of vectors and matrices match */
2272d61bbb3SSatish Balay /* -------------------------------------------------------*/
2282d61bbb3SSatish Balay 
2294a2ae208SSatish Balay #undef __FUNCT__
2304a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_1"
231dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_1(Mat A,Vec xx,Vec zz)
2322d61bbb3SSatish Balay {
2332d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
234d9fead3dSBarry Smith   PetscScalar       *z,sum;
235d9fead3dSBarry Smith   const PetscScalar *x;
236d9fead3dSBarry Smith   const MatScalar   *v;
2376849ba73SBarry Smith   PetscErrorCode    ierr;
2382162cab8SBarry Smith   PetscInt          mbs,i,n,nonzerorow=0;
2392162cab8SBarry Smith   const PetscInt    *idx,*ii,*ridx=PETSC_NULL;
240ace3abfcSBarry Smith   PetscBool         usecprow=a->compressedrow.use;
2412d61bbb3SSatish Balay 
2422d61bbb3SSatish Balay   PetscFunctionBegin;
2433649974fSBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
2441ebc52fbSHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
2452d61bbb3SSatish Balay 
24626e093fcSHong Zhang   if (usecprow){
24726e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
24826e093fcSHong Zhang     ii   = a->compressedrow.i;
2497b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
250a1c3900fSBarry Smith     ierr = PetscMemzero(z,mbs*sizeof(PetscScalar));CHKERRQ(ierr);
25126e093fcSHong Zhang   } else {
25226e093fcSHong Zhang     mbs = a->mbs;
2532d61bbb3SSatish Balay     ii  = a->i;
25426e093fcSHong Zhang   }
2552d61bbb3SSatish Balay 
2562d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
257ee54c7eeSHong Zhang     n    = ii[1] - ii[0];
258ee54c7eeSHong Zhang     v    = a->a + ii[0];
259ee54c7eeSHong Zhang     idx  = a->j + ii[0];
260ee54c7eeSHong Zhang     ii++;
261444d8c10SJed Brown     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);   /* Indices for the next row (assumes same size as this one) */
262444d8c10SJed Brown     PetscPrefetchBlock(v+1*n,1*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
2632d61bbb3SSatish Balay     sum  = 0.0;
2642162cab8SBarry Smith     PetscSparseDensePlusDot(sum,x,v,idx,n);
26526e093fcSHong Zhang     if (usecprow){
2667b2bb3b9SHong Zhang       z[ridx[i]] = sum;
26726e093fcSHong Zhang     } else {
268122f12eaSBarry Smith       nonzerorow += (n>0);
2692d61bbb3SSatish Balay       z[i] = sum;
2702d61bbb3SSatish Balay     }
27126e093fcSHong Zhang   }
2723649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
2731ebc52fbSHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
274dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*a->nz - nonzerorow);CHKERRQ(ierr);
2752d61bbb3SSatish Balay   PetscFunctionReturn(0);
2762d61bbb3SSatish Balay }
2772d61bbb3SSatish Balay 
2784a2ae208SSatish Balay #undef __FUNCT__
2794a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_2"
280dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_2(Mat A,Vec xx,Vec zz)
2812d61bbb3SSatish Balay {
2822d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
283d9fead3dSBarry Smith   PetscScalar       *z = 0,sum1,sum2,*zarray;
284d9fead3dSBarry Smith   const PetscScalar *x,*xb;
28587828ca2SBarry Smith   PetscScalar       x1,x2;
286d9fead3dSBarry Smith   const MatScalar   *v;
287dfbe8321SBarry Smith   PetscErrorCode    ierr;
28898c9bda7SSatish Balay   PetscInt          mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL,nonzerorow=0;
289ace3abfcSBarry Smith   PetscBool         usecprow=a->compressedrow.use;
2902d61bbb3SSatish Balay 
2912d61bbb3SSatish Balay   PetscFunctionBegin;
2923649974fSBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
29326e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
2942d61bbb3SSatish Balay 
2952d61bbb3SSatish Balay   idx = a->j;
2962d61bbb3SSatish Balay   v   = a->a;
29726e093fcSHong Zhang   if (usecprow){
29826e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
29926e093fcSHong Zhang     ii   = a->compressedrow.i;
3007b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
30126e093fcSHong Zhang   } else {
30226e093fcSHong Zhang     mbs = a->mbs;
3032d61bbb3SSatish Balay     ii  = a->i;
30426e093fcSHong Zhang     z   = zarray;
30526e093fcSHong Zhang   }
3062d61bbb3SSatish Balay 
3072d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
3082d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
3092d61bbb3SSatish Balay     sum1 = 0.0; sum2 = 0.0;
31098c9bda7SSatish Balay     nonzerorow += (n>0);
311444d8c10SJed Brown     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);   /* Indices for the next row (assumes same size as this one) */
312444d8c10SJed Brown     PetscPrefetchBlock(v+4*n,4*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
3132d61bbb3SSatish Balay     for (j=0; j<n; j++) {
3142d61bbb3SSatish Balay       xb = x + 2*(*idx++); x1 = xb[0]; x2 = xb[1];
3152d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[2]*x2;
3162d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[3]*x2;
3172d61bbb3SSatish Balay       v += 4;
3182d61bbb3SSatish Balay     }
3197b2bb3b9SHong Zhang     if (usecprow) z = zarray + 2*ridx[i];
3202d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2;
32126e093fcSHong Zhang     if (!usecprow) z += 2;
3222d61bbb3SSatish Balay   }
3233649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
32426e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
325dc0b31edSSatish Balay   ierr = PetscLogFlops(8.0*a->nz - 2.0*nonzerorow);CHKERRQ(ierr);
3262d61bbb3SSatish Balay   PetscFunctionReturn(0);
3272d61bbb3SSatish Balay }
3282d61bbb3SSatish Balay 
3294a2ae208SSatish Balay #undef __FUNCT__
3304a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_3"
331dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_3(Mat A,Vec xx,Vec zz)
3322d61bbb3SSatish Balay {
3332d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
334d9fead3dSBarry Smith   PetscScalar       *z = 0,sum1,sum2,sum3,x1,x2,x3,*zarray;
335d9fead3dSBarry Smith   const PetscScalar *x,*xb;
336d9fead3dSBarry Smith   const MatScalar   *v;
337dfbe8321SBarry Smith   PetscErrorCode    ierr;
338028cd4eaSSatish Balay   PetscInt          mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL,nonzerorow=0;
339ace3abfcSBarry Smith   PetscBool         usecprow=a->compressedrow.use;
34026e093fcSHong Zhang 
3412d61bbb3SSatish Balay 
342b6410449SSatish Balay #if defined(PETSC_HAVE_PRAGMA_DISJOINT)
343fee21e36SBarry Smith #pragma disjoint(*v,*z,*xb)
344fee21e36SBarry Smith #endif
345fee21e36SBarry Smith 
3462d61bbb3SSatish Balay   PetscFunctionBegin;
3473649974fSBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
34826e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
3492d61bbb3SSatish Balay 
3502d61bbb3SSatish Balay   idx = a->j;
3512d61bbb3SSatish Balay   v   = a->a;
35226e093fcSHong Zhang   if (usecprow){
35326e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
35426e093fcSHong Zhang     ii   = a->compressedrow.i;
3557b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
35626e093fcSHong Zhang   } else {
35726e093fcSHong Zhang     mbs = a->mbs;
3582d61bbb3SSatish Balay     ii  = a->i;
35926e093fcSHong Zhang     z   = zarray;
36026e093fcSHong Zhang   }
3612d61bbb3SSatish Balay 
3622d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
3632d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
3642d61bbb3SSatish Balay     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0;
36598c9bda7SSatish Balay     nonzerorow += (n>0);
366444d8c10SJed Brown     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);   /* Indices for the next row (assumes same size as this one) */
367444d8c10SJed Brown     PetscPrefetchBlock(v+9*n,9*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
3682d61bbb3SSatish Balay     for (j=0; j<n; j++) {
3692d61bbb3SSatish Balay       xb = x + 3*(*idx++); x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
3702d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[3]*x2 + v[6]*x3;
3712d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[4]*x2 + v[7]*x3;
3722d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[5]*x2 + v[8]*x3;
3732d61bbb3SSatish Balay       v += 9;
3742d61bbb3SSatish Balay     }
3757b2bb3b9SHong Zhang     if (usecprow) z = zarray + 3*ridx[i];
3762d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3;
37726e093fcSHong Zhang     if (!usecprow) z += 3;
3782d61bbb3SSatish Balay   }
3793649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
38026e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
381dc0b31edSSatish Balay   ierr = PetscLogFlops(18.0*a->nz - 3.0*nonzerorow);CHKERRQ(ierr);
3822d61bbb3SSatish Balay   PetscFunctionReturn(0);
3832d61bbb3SSatish Balay }
3842d61bbb3SSatish Balay 
3854a2ae208SSatish Balay #undef __FUNCT__
3864a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_4"
387dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_4(Mat A,Vec xx,Vec zz)
3882d61bbb3SSatish Balay {
3892d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
390d9fead3dSBarry Smith   PetscScalar       *z = 0,sum1,sum2,sum3,sum4,x1,x2,x3,x4,*zarray;
391d9fead3dSBarry Smith   const PetscScalar *x,*xb;
392d9fead3dSBarry Smith   const MatScalar   *v;
393dfbe8321SBarry Smith   PetscErrorCode    ierr;
39498c9bda7SSatish Balay   PetscInt          mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL,nonzerorow=0;
395ace3abfcSBarry Smith   PetscBool         usecprow=a->compressedrow.use;
3962d61bbb3SSatish Balay 
3972d61bbb3SSatish Balay   PetscFunctionBegin;
3983649974fSBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
39926e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
4002d61bbb3SSatish Balay 
4012d61bbb3SSatish Balay   idx = a->j;
4022d61bbb3SSatish Balay   v   = a->a;
40326e093fcSHong Zhang   if (usecprow){
40426e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
40526e093fcSHong Zhang     ii   = a->compressedrow.i;
4067b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
40726e093fcSHong Zhang   } else {
40826e093fcSHong Zhang     mbs = a->mbs;
4092d61bbb3SSatish Balay     ii  = a->i;
41026e093fcSHong Zhang     z   = zarray;
41126e093fcSHong Zhang   }
4122d61bbb3SSatish Balay 
4132d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
4142d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
4152d61bbb3SSatish Balay     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0;
41698c9bda7SSatish Balay     nonzerorow += (n>0);
417444d8c10SJed Brown     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);     /* Indices for the next row (assumes same size as this one) */
418444d8c10SJed Brown     PetscPrefetchBlock(v+16*n,16*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
4192d61bbb3SSatish Balay     for (j=0; j<n; j++) {
4202d61bbb3SSatish Balay       xb = x + 4*(*idx++);
4212d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
4222d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[4]*x2 + v[8]*x3   + v[12]*x4;
4232d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[5]*x2 + v[9]*x3   + v[13]*x4;
4242d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[6]*x2 + v[10]*x3  + v[14]*x4;
4252d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[7]*x2 + v[11]*x3  + v[15]*x4;
4262d61bbb3SSatish Balay       v += 16;
4272d61bbb3SSatish Balay     }
4287b2bb3b9SHong Zhang     if (usecprow) z = zarray + 4*ridx[i];
4292d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4;
43026e093fcSHong Zhang     if (!usecprow) z += 4;
4312d61bbb3SSatish Balay   }
4323649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
43326e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
434dc0b31edSSatish Balay   ierr = PetscLogFlops(32.0*a->nz - 4.0*nonzerorow);CHKERRQ(ierr);
4352d61bbb3SSatish Balay   PetscFunctionReturn(0);
4362d61bbb3SSatish Balay }
4372d61bbb3SSatish Balay 
4384a2ae208SSatish Balay #undef __FUNCT__
4394a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_5"
440dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_5(Mat A,Vec xx,Vec zz)
4412d61bbb3SSatish Balay {
4422d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
443d9fead3dSBarry Smith   PetscScalar       sum1,sum2,sum3,sum4,sum5,x1,x2,x3,x4,x5,*z = 0,*zarray;
444d9fead3dSBarry Smith   const PetscScalar *xb,*x;
445d9fead3dSBarry Smith   const MatScalar   *v;
446dfbe8321SBarry Smith   PetscErrorCode    ierr;
447766f9fbaSBarry Smith   const PetscInt    *idx,*ii,*ridx=PETSC_NULL;
448766f9fbaSBarry Smith   PetscInt          mbs,i,j,n,nonzerorow=0;
449ace3abfcSBarry Smith   PetscBool         usecprow=a->compressedrow.use;
4502d61bbb3SSatish Balay 
451433994e6SBarry Smith   PetscFunctionBegin;
4523649974fSBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
45326e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
4542d61bbb3SSatish Balay 
4552d61bbb3SSatish Balay   idx = a->j;
4562d61bbb3SSatish Balay   v   = a->a;
45726e093fcSHong Zhang   if (usecprow){
45826e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
45926e093fcSHong Zhang     ii   = a->compressedrow.i;
4607b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
46126e093fcSHong Zhang   } else {
46226e093fcSHong Zhang     mbs = a->mbs;
4632d61bbb3SSatish Balay     ii  = a->i;
46426e093fcSHong Zhang     z   = zarray;
46526e093fcSHong Zhang   }
4662d61bbb3SSatish Balay 
4672d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
4682d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
4692d61bbb3SSatish Balay     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0;
47098c9bda7SSatish Balay     nonzerorow += (n>0);
471444d8c10SJed Brown     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);     /* Indices for the next row (assumes same size as this one) */
472444d8c10SJed Brown     PetscPrefetchBlock(v+25*n,25*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
4732d61bbb3SSatish Balay     for (j=0; j<n; j++) {
4742d61bbb3SSatish Balay       xb = x + 5*(*idx++);
4752d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4];
4762d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[5]*x2 + v[10]*x3  + v[15]*x4 + v[20]*x5;
4772d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[6]*x2 + v[11]*x3  + v[16]*x4 + v[21]*x5;
4782d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[7]*x2 + v[12]*x3  + v[17]*x4 + v[22]*x5;
4792d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[8]*x2 + v[13]*x3  + v[18]*x4 + v[23]*x5;
4802d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[9]*x2 + v[14]*x3  + v[19]*x4 + v[24]*x5;
4812d61bbb3SSatish Balay       v += 25;
4822d61bbb3SSatish Balay     }
4837b2bb3b9SHong Zhang     if (usecprow) z = zarray + 5*ridx[i];
4842d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5;
48526e093fcSHong Zhang     if (!usecprow) z += 5;
4862d61bbb3SSatish Balay   }
4873649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
48826e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
489dc0b31edSSatish Balay   ierr = PetscLogFlops(50.0*a->nz - 5.0*nonzerorow);CHKERRQ(ierr);
4902d61bbb3SSatish Balay   PetscFunctionReturn(0);
4912d61bbb3SSatish Balay }
4922d61bbb3SSatish Balay 
49315091d37SBarry Smith 
4944a2ae208SSatish Balay #undef __FUNCT__
4954a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_6"
496dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_6(Mat A,Vec xx,Vec zz)
49715091d37SBarry Smith {
49815091d37SBarry Smith   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
499d9fead3dSBarry Smith   PetscScalar       *z = 0,sum1,sum2,sum3,sum4,sum5,sum6;
500d9fead3dSBarry Smith   const PetscScalar *x,*xb;
50126e093fcSHong Zhang   PetscScalar       x1,x2,x3,x4,x5,x6,*zarray;
502d9fead3dSBarry Smith   const MatScalar   *v;
503dfbe8321SBarry Smith   PetscErrorCode    ierr;
50498c9bda7SSatish Balay   PetscInt          mbs=a->mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL,nonzerorow=0;
505ace3abfcSBarry Smith   PetscBool         usecprow=a->compressedrow.use;
50615091d37SBarry Smith 
507433994e6SBarry Smith   PetscFunctionBegin;
5083649974fSBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
50926e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
51015091d37SBarry Smith 
51115091d37SBarry Smith   idx = a->j;
51215091d37SBarry Smith   v   = a->a;
51326e093fcSHong Zhang   if (usecprow){
51426e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
51526e093fcSHong Zhang     ii   = a->compressedrow.i;
5167b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
51726e093fcSHong Zhang   } else {
51826e093fcSHong Zhang     mbs = a->mbs;
51915091d37SBarry Smith     ii  = a->i;
52026e093fcSHong Zhang     z   = zarray;
52126e093fcSHong Zhang   }
52215091d37SBarry Smith 
52315091d37SBarry Smith   for (i=0; i<mbs; i++) {
52415091d37SBarry Smith     n  = ii[1] - ii[0]; ii++;
52515091d37SBarry Smith     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0;
52698c9bda7SSatish Balay     nonzerorow += (n>0);
527444d8c10SJed Brown     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);     /* Indices for the next row (assumes same size as this one) */
528444d8c10SJed Brown     PetscPrefetchBlock(v+36*n,36*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
52915091d37SBarry Smith     for (j=0; j<n; j++) {
53015091d37SBarry Smith       xb = x + 6*(*idx++);
53115091d37SBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5];
53215091d37SBarry Smith       sum1 += v[0]*x1 + v[6]*x2  + v[12]*x3  + v[18]*x4 + v[24]*x5 + v[30]*x6;
53315091d37SBarry Smith       sum2 += v[1]*x1 + v[7]*x2  + v[13]*x3  + v[19]*x4 + v[25]*x5 + v[31]*x6;
53415091d37SBarry Smith       sum3 += v[2]*x1 + v[8]*x2  + v[14]*x3  + v[20]*x4 + v[26]*x5 + v[32]*x6;
53515091d37SBarry Smith       sum4 += v[3]*x1 + v[9]*x2  + v[15]*x3  + v[21]*x4 + v[27]*x5 + v[33]*x6;
53615091d37SBarry Smith       sum5 += v[4]*x1 + v[10]*x2 + v[16]*x3  + v[22]*x4 + v[28]*x5 + v[34]*x6;
53715091d37SBarry Smith       sum6 += v[5]*x1 + v[11]*x2 + v[17]*x3  + v[23]*x4 + v[29]*x5 + v[35]*x6;
53815091d37SBarry Smith       v += 36;
53915091d37SBarry Smith     }
5407b2bb3b9SHong Zhang     if (usecprow) z = zarray + 6*ridx[i];
54115091d37SBarry Smith     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6;
54226e093fcSHong Zhang     if (!usecprow) z += 6;
54315091d37SBarry Smith   }
54415091d37SBarry Smith 
5453649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
54626e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
547dc0b31edSSatish Balay   ierr = PetscLogFlops(72.0*a->nz - 6.0*nonzerorow);CHKERRQ(ierr);
54815091d37SBarry Smith   PetscFunctionReturn(0);
54915091d37SBarry Smith }
5508ab949d8SShri Abhyankar 
5514a2ae208SSatish Balay #undef __FUNCT__
5524a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_7"
553dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_7(Mat A,Vec xx,Vec zz)
5542d61bbb3SSatish Balay {
5552d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
556d9fead3dSBarry Smith   PetscScalar       *z = 0,sum1,sum2,sum3,sum4,sum5,sum6,sum7;
557d9fead3dSBarry Smith   const PetscScalar *x,*xb;
55826e093fcSHong Zhang   PetscScalar       x1,x2,x3,x4,x5,x6,x7,*zarray;
559d9fead3dSBarry Smith   const MatScalar   *v;
560dfbe8321SBarry Smith   PetscErrorCode    ierr;
56198c9bda7SSatish Balay   PetscInt          mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL,nonzerorow=0;
562ace3abfcSBarry Smith   PetscBool         usecprow=a->compressedrow.use;
5632d61bbb3SSatish Balay 
564433994e6SBarry Smith   PetscFunctionBegin;
5653649974fSBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
56626e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
5672d61bbb3SSatish Balay 
5682d61bbb3SSatish Balay   idx = a->j;
5692d61bbb3SSatish Balay   v   = a->a;
57026e093fcSHong Zhang   if (usecprow){
57126e093fcSHong Zhang     mbs    = a->compressedrow.nrows;
57226e093fcSHong Zhang     ii     = a->compressedrow.i;
5737b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
57426e093fcSHong Zhang   } else {
57526e093fcSHong Zhang     mbs = a->mbs;
5762d61bbb3SSatish Balay     ii  = a->i;
57726e093fcSHong Zhang     z   = zarray;
57826e093fcSHong Zhang   }
5792d61bbb3SSatish Balay 
5802d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
5812d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
5822d61bbb3SSatish Balay     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0; sum7 = 0.0;
58398c9bda7SSatish Balay     nonzerorow += (n>0);
584444d8c10SJed Brown     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);     /* Indices for the next row (assumes same size as this one) */
585444d8c10SJed Brown     PetscPrefetchBlock(v+49*n,49*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
5862d61bbb3SSatish Balay     for (j=0; j<n; j++) {
5872d61bbb3SSatish Balay       xb = x + 7*(*idx++);
5882d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6];
5892d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[7]*x2  + v[14]*x3  + v[21]*x4 + v[28]*x5 + v[35]*x6 + v[42]*x7;
5902d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[8]*x2  + v[15]*x3  + v[22]*x4 + v[29]*x5 + v[36]*x6 + v[43]*x7;
5912d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[9]*x2  + v[16]*x3  + v[23]*x4 + v[30]*x5 + v[37]*x6 + v[44]*x7;
5922d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[10]*x2 + v[17]*x3  + v[24]*x4 + v[31]*x5 + v[38]*x6 + v[45]*x7;
5932d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[11]*x2 + v[18]*x3  + v[25]*x4 + v[32]*x5 + v[39]*x6 + v[46]*x7;
5942d61bbb3SSatish Balay       sum6 += v[5]*x1 + v[12]*x2 + v[19]*x3  + v[26]*x4 + v[33]*x5 + v[40]*x6 + v[47]*x7;
5952d61bbb3SSatish Balay       sum7 += v[6]*x1 + v[13]*x2 + v[20]*x3  + v[27]*x4 + v[34]*x5 + v[41]*x6 + v[48]*x7;
5962d61bbb3SSatish Balay       v += 49;
5972d61bbb3SSatish Balay     }
5987b2bb3b9SHong Zhang     if (usecprow) z = zarray + 7*ridx[i];
5992d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
60026e093fcSHong Zhang     if (!usecprow) z += 7;
6012d61bbb3SSatish Balay   }
6022d61bbb3SSatish Balay 
6033649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
60426e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
605dc0b31edSSatish Balay   ierr = PetscLogFlops(98.0*a->nz - 7.0*nonzerorow);CHKERRQ(ierr);
6062d61bbb3SSatish Balay   PetscFunctionReturn(0);
6072d61bbb3SSatish Balay }
6082d61bbb3SSatish Balay 
6098ab949d8SShri Abhyankar /* MatMult_SeqBAIJ_15 version 1: Columns in the block are accessed one at a time */
610832cc040SShri Abhyankar /* Default MatMult for block size 15 */
6118ab949d8SShri Abhyankar 
6128ab949d8SShri Abhyankar #undef __FUNCT__
6138ab949d8SShri Abhyankar #define __FUNCT__ "MatMult_SeqBAIJ_15_ver1"
6148ab949d8SShri Abhyankar PetscErrorCode MatMult_SeqBAIJ_15_ver1(Mat A,Vec xx,Vec zz)
6158ab949d8SShri Abhyankar {
6168ab949d8SShri Abhyankar   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
6178ab949d8SShri Abhyankar   PetscScalar       *z = 0,sum1,sum2,sum3,sum4,sum5,sum6,sum7,sum8,sum9,sum10,sum11,sum12,sum13,sum14,sum15;
6188ab949d8SShri Abhyankar   const PetscScalar *x,*xb;
61953ef36baSBarry Smith   PetscScalar       *zarray,xv;
6208ab949d8SShri Abhyankar   const MatScalar   *v;
6218ab949d8SShri Abhyankar   PetscErrorCode    ierr;
6228ab949d8SShri Abhyankar   const PetscInt    *ii,*ij=a->j,*idx;
6238ab949d8SShri Abhyankar   PetscInt          mbs,i,j,k,n,*ridx=PETSC_NULL,nonzerorow=0;
624ace3abfcSBarry Smith   PetscBool         usecprow=a->compressedrow.use;
6258ab949d8SShri Abhyankar 
6268ab949d8SShri Abhyankar   PetscFunctionBegin;
6273649974fSBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
6288ab949d8SShri Abhyankar   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
6298ab949d8SShri Abhyankar 
6308ab949d8SShri Abhyankar   v   = a->a;
6318ab949d8SShri Abhyankar   if (usecprow){
6328ab949d8SShri Abhyankar     mbs    = a->compressedrow.nrows;
6338ab949d8SShri Abhyankar     ii     = a->compressedrow.i;
6348ab949d8SShri Abhyankar     ridx = a->compressedrow.rindex;
6358ab949d8SShri Abhyankar   } else {
6368ab949d8SShri Abhyankar     mbs = a->mbs;
6378ab949d8SShri Abhyankar     ii  = a->i;
6388ab949d8SShri Abhyankar     z   = zarray;
6398ab949d8SShri Abhyankar   }
6408ab949d8SShri Abhyankar 
6418ab949d8SShri Abhyankar   for (i=0; i<mbs; i++) {
6428ab949d8SShri Abhyankar     n  = ii[i+1] - ii[i];
6438ab949d8SShri Abhyankar     idx = ij + ii[i];
6448ab949d8SShri Abhyankar     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0; sum7 = 0.0;
6458ab949d8SShri 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;
6468ab949d8SShri Abhyankar 
6478ab949d8SShri Abhyankar     nonzerorow += (n>0);
6488ab949d8SShri Abhyankar     for (j=0; j<n; j++) {
6498ab949d8SShri Abhyankar       xb = x + 15*(idx[j]);
6508ab949d8SShri Abhyankar 
6518ab949d8SShri Abhyankar       for(k=0;k<15;k++){
65253ef36baSBarry Smith         xv    =  xb[k];
65353ef36baSBarry Smith 	sum1  += v[0]*xv;
65453ef36baSBarry Smith         sum2  += v[1]*xv;
65553ef36baSBarry Smith 	sum3  += v[2]*xv;
65653ef36baSBarry Smith 	sum4  += v[3]*xv;
65753ef36baSBarry Smith 	sum5  += v[4]*xv;
65853ef36baSBarry Smith         sum6  += v[5]*xv;
65953ef36baSBarry Smith 	sum7  += v[6]*xv;
66053ef36baSBarry Smith 	sum8  += v[7]*xv;
66153ef36baSBarry Smith         sum9  += v[8]*xv;
66253ef36baSBarry Smith         sum10 += v[9]*xv;
66353ef36baSBarry Smith 	sum11 += v[10]*xv;
66453ef36baSBarry Smith 	sum12 += v[11]*xv;
66553ef36baSBarry Smith 	sum13 += v[12]*xv;
66653ef36baSBarry Smith         sum14 += v[13]*xv;
66753ef36baSBarry Smith 	sum15 += v[14]*xv;
6688ab949d8SShri Abhyankar 	v += 15;
6698ab949d8SShri Abhyankar       }
6708ab949d8SShri Abhyankar     }
6718ab949d8SShri Abhyankar     if (usecprow) z = zarray + 15*ridx[i];
6728ab949d8SShri Abhyankar     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
6738ab949d8SShri 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;
6748ab949d8SShri Abhyankar 
6758ab949d8SShri Abhyankar     if (!usecprow) z += 15;
6768ab949d8SShri Abhyankar   }
6778ab949d8SShri Abhyankar 
6783649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
6798ab949d8SShri Abhyankar   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
6808ab949d8SShri Abhyankar   ierr = PetscLogFlops(450.0*a->nz - 15.0*nonzerorow);CHKERRQ(ierr);
6818ab949d8SShri Abhyankar   PetscFunctionReturn(0);
6828ab949d8SShri Abhyankar }
6838ab949d8SShri Abhyankar 
6848ab949d8SShri Abhyankar /* MatMult_SeqBAIJ_15_ver2 : Columns in the block are accessed in sets of 4,4,4,3 */
6858ab949d8SShri Abhyankar #undef __FUNCT__
6868ab949d8SShri Abhyankar #define __FUNCT__ "MatMult_SeqBAIJ_15_ver2"
6878ab949d8SShri Abhyankar PetscErrorCode MatMult_SeqBAIJ_15_ver2(Mat A,Vec xx,Vec zz)
6888ab949d8SShri Abhyankar {
6898ab949d8SShri Abhyankar   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
6908ab949d8SShri Abhyankar   PetscScalar       *z = 0,sum1,sum2,sum3,sum4,sum5,sum6,sum7,sum8,sum9,sum10,sum11,sum12,sum13,sum14,sum15;
6918ab949d8SShri Abhyankar   const PetscScalar *x,*xb;
6920b8f6341SShri Abhyankar   PetscScalar        x1,x2,x3,x4,*zarray;
6938ab949d8SShri Abhyankar   const MatScalar   *v;
6948ab949d8SShri Abhyankar   PetscErrorCode    ierr;
6958ab949d8SShri Abhyankar   const PetscInt    *ii,*ij=a->j,*idx;
6968ab949d8SShri Abhyankar   PetscInt          mbs,i,j,n,*ridx=PETSC_NULL,nonzerorow=0;
697ace3abfcSBarry Smith   PetscBool         usecprow=a->compressedrow.use;
6988ab949d8SShri Abhyankar 
6998ab949d8SShri Abhyankar   PetscFunctionBegin;
7003649974fSBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
7018ab949d8SShri Abhyankar   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
7028ab949d8SShri Abhyankar 
7038ab949d8SShri Abhyankar   v   = a->a;
7048ab949d8SShri Abhyankar   if (usecprow){
7058ab949d8SShri Abhyankar     mbs    = a->compressedrow.nrows;
7068ab949d8SShri Abhyankar     ii     = a->compressedrow.i;
7078ab949d8SShri Abhyankar     ridx = a->compressedrow.rindex;
7088ab949d8SShri Abhyankar   } else {
7098ab949d8SShri Abhyankar     mbs = a->mbs;
7108ab949d8SShri Abhyankar     ii  = a->i;
7118ab949d8SShri Abhyankar     z   = zarray;
7128ab949d8SShri Abhyankar   }
7138ab949d8SShri Abhyankar 
7148ab949d8SShri Abhyankar   for (i=0; i<mbs; i++) {
7158ab949d8SShri Abhyankar     n  = ii[i+1] - ii[i];
7168ab949d8SShri Abhyankar     idx = ij + ii[i];
7178ab949d8SShri Abhyankar     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0; sum7 = 0.0;
7188ab949d8SShri 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;
7198ab949d8SShri Abhyankar 
7208ab949d8SShri Abhyankar     nonzerorow += (n>0);
7218ab949d8SShri Abhyankar     for (j=0; j<n; j++) {
7228ab949d8SShri Abhyankar       xb = x + 15*(idx[j]);
7230b8f6341SShri Abhyankar       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
7248ab949d8SShri Abhyankar 
7258ab949d8SShri Abhyankar       sum1  += v[0]*x1 + v[15]*x2 + v[30]*x3   + v[45]*x4;
7268ab949d8SShri Abhyankar       sum2  += v[1]*x1 + v[16]*x2 + v[31]*x3   + v[46]*x4;
7278ab949d8SShri Abhyankar       sum3  += v[2]*x1 + v[17]*x2 + v[32]*x3  + v[47]*x4;
7288ab949d8SShri Abhyankar       sum4  += v[3]*x1 + v[18]*x2 + v[33]*x3  + v[48]*x4;
7298ab949d8SShri Abhyankar       sum5  += v[4]*x1 + v[19]*x2 + v[34]*x3   + v[49]*x4;
7308ab949d8SShri Abhyankar       sum6  += v[5]*x1 + v[20]*x2 + v[35]*x3   + v[50]*x4;
7318ab949d8SShri Abhyankar       sum7  += v[6]*x1 + v[21]*x2 + v[36]*x3  + v[51]*x4;
7328ab949d8SShri Abhyankar       sum8  += v[7]*x1 + v[22]*x2 + v[37]*x3  + v[52]*x4;
7338ab949d8SShri Abhyankar       sum9  += v[8]*x1 + v[23]*x2 + v[38]*x3   + v[53]*x4;
7348ab949d8SShri Abhyankar       sum10 += v[9]*x1 + v[24]*x2 + v[39]*x3   + v[54]*x4;
7358ab949d8SShri Abhyankar       sum11 += v[10]*x1 + v[25]*x2 + v[40]*x3  + v[55]*x4;
7368ab949d8SShri Abhyankar       sum12 += v[11]*x1 + v[26]*x2 + v[41]*x3  + v[56]*x4;
7378ab949d8SShri Abhyankar       sum13 += v[12]*x1 + v[27]*x2 + v[42]*x3   + v[57]*x4;
7388ab949d8SShri Abhyankar       sum14 += v[13]*x1 + v[28]*x2 + v[43]*x3   + v[58]*x4;
7398ab949d8SShri Abhyankar       sum15 += v[14]*x1 + v[29]*x2 + v[44]*x3  + v[59]*x4;
7408ab949d8SShri Abhyankar 
7418ab949d8SShri Abhyankar       v += 60;
7428ab949d8SShri Abhyankar 
7430b8f6341SShri Abhyankar       x1 = xb[4]; x2 = xb[5]; x3 = xb[6]; x4 = xb[7];
7448ab949d8SShri Abhyankar 
7458ab949d8SShri Abhyankar       sum1  += v[0]*x1 + v[15]*x2 + v[30]*x3   + v[45]*x4;
7468ab949d8SShri Abhyankar       sum2  += v[1]*x1 + v[16]*x2 + v[31]*x3   + v[46]*x4;
7478ab949d8SShri Abhyankar       sum3  += v[2]*x1 + v[17]*x2 + v[32]*x3  + v[47]*x4;
7488ab949d8SShri Abhyankar       sum4  += v[3]*x1 + v[18]*x2 + v[33]*x3  + v[48]*x4;
7498ab949d8SShri Abhyankar       sum5  += v[4]*x1 + v[19]*x2 + v[34]*x3   + v[49]*x4;
7508ab949d8SShri Abhyankar       sum6  += v[5]*x1 + v[20]*x2 + v[35]*x3   + v[50]*x4;
7518ab949d8SShri Abhyankar       sum7  += v[6]*x1 + v[21]*x2 + v[36]*x3  + v[51]*x4;
7528ab949d8SShri Abhyankar       sum8  += v[7]*x1 + v[22]*x2 + v[37]*x3  + v[52]*x4;
7538ab949d8SShri Abhyankar       sum9  += v[8]*x1 + v[23]*x2 + v[38]*x3   + v[53]*x4;
7548ab949d8SShri Abhyankar       sum10 += v[9]*x1 + v[24]*x2 + v[39]*x3   + v[54]*x4;
7558ab949d8SShri Abhyankar       sum11 += v[10]*x1 + v[25]*x2 + v[40]*x3  + v[55]*x4;
7568ab949d8SShri Abhyankar       sum12 += v[11]*x1 + v[26]*x2 + v[41]*x3  + v[56]*x4;
7578ab949d8SShri Abhyankar       sum13 += v[12]*x1 + v[27]*x2 + v[42]*x3   + v[57]*x4;
7588ab949d8SShri Abhyankar       sum14 += v[13]*x1 + v[28]*x2 + v[43]*x3   + v[58]*x4;
7598ab949d8SShri Abhyankar       sum15 += v[14]*x1 + v[29]*x2 + v[44]*x3  + v[59]*x4;
7608ab949d8SShri Abhyankar       v += 60;
7618ab949d8SShri Abhyankar 
7620b8f6341SShri Abhyankar       x1 = xb[8]; x2 = xb[9]; x3 = xb[10]; x4 = xb[11];
7630b8f6341SShri Abhyankar       sum1  += v[0]*x1 + v[15]*x2 + v[30]*x3   + v[45]*x4;
7640b8f6341SShri Abhyankar       sum2  += v[1]*x1 + v[16]*x2 + v[31]*x3   + v[46]*x4;
7650b8f6341SShri Abhyankar       sum3  += v[2]*x1 + v[17]*x2 + v[32]*x3  + v[47]*x4;
7660b8f6341SShri Abhyankar       sum4  += v[3]*x1 + v[18]*x2 + v[33]*x3  + v[48]*x4;
7670b8f6341SShri Abhyankar       sum5  += v[4]*x1 + v[19]*x2 + v[34]*x3   + v[49]*x4;
7680b8f6341SShri Abhyankar       sum6  += v[5]*x1 + v[20]*x2 + v[35]*x3   + v[50]*x4;
7690b8f6341SShri Abhyankar       sum7  += v[6]*x1 + v[21]*x2 + v[36]*x3  + v[51]*x4;
7700b8f6341SShri Abhyankar       sum8  += v[7]*x1 + v[22]*x2 + v[37]*x3  + v[52]*x4;
7710b8f6341SShri Abhyankar       sum9  += v[8]*x1 + v[23]*x2 + v[38]*x3   + v[53]*x4;
7720b8f6341SShri Abhyankar       sum10 += v[9]*x1 + v[24]*x2 + v[39]*x3   + v[54]*x4;
7730b8f6341SShri Abhyankar       sum11 += v[10]*x1 + v[25]*x2 + v[40]*x3  + v[55]*x4;
7740b8f6341SShri Abhyankar       sum12 += v[11]*x1 + v[26]*x2 + v[41]*x3  + v[56]*x4;
7750b8f6341SShri Abhyankar       sum13 += v[12]*x1 + v[27]*x2 + v[42]*x3   + v[57]*x4;
7760b8f6341SShri Abhyankar       sum14 += v[13]*x1 + v[28]*x2 + v[43]*x3   + v[58]*x4;
7770b8f6341SShri Abhyankar       sum15 += v[14]*x1 + v[29]*x2 + v[44]*x3  + v[59]*x4;
7780b8f6341SShri Abhyankar       v  += 60;
7790b8f6341SShri Abhyankar 
7800b8f6341SShri Abhyankar       x1 = xb[12]; x2 = xb[13]; x3 = xb[14];
7818ab949d8SShri Abhyankar       sum1  += v[0]*x1 + v[15]*x2 + v[30]*x3;
7828ab949d8SShri Abhyankar       sum2  += v[1]*x1 + v[16]*x2 + v[31]*x3;
7838ab949d8SShri Abhyankar       sum3  += v[2]*x1 + v[17]*x2 + v[32]*x3;
7848ab949d8SShri Abhyankar       sum4  += v[3]*x1 + v[18]*x2 + v[33]*x3;
7858ab949d8SShri Abhyankar       sum5  += v[4]*x1 + v[19]*x2 + v[34]*x3;
7868ab949d8SShri Abhyankar       sum6  += v[5]*x1 + v[20]*x2 + v[35]*x3;
7878ab949d8SShri Abhyankar       sum7  += v[6]*x1 + v[21]*x2 + v[36]*x3;
7888ab949d8SShri Abhyankar       sum8  += v[7]*x1 + v[22]*x2 + v[37]*x3;
7898ab949d8SShri Abhyankar       sum9  += v[8]*x1 + v[23]*x2 + v[38]*x3;
7908ab949d8SShri Abhyankar       sum10 += v[9]*x1 + v[24]*x2 + v[39]*x3;
7918ab949d8SShri Abhyankar       sum11 += v[10]*x1 + v[25]*x2 + v[40]*x3;
7928ab949d8SShri Abhyankar       sum12 += v[11]*x1 + v[26]*x2 + v[41]*x3;
7938ab949d8SShri Abhyankar       sum13 += v[12]*x1 + v[27]*x2 + v[42]*x3;
7948ab949d8SShri Abhyankar       sum14 += v[13]*x1 + v[28]*x2 + v[43]*x3;
7958ab949d8SShri Abhyankar       sum15 += v[14]*x1 + v[29]*x2 + v[44]*x3;
7968ab949d8SShri Abhyankar       v += 45;
7978ab949d8SShri Abhyankar     }
7988ab949d8SShri Abhyankar     if (usecprow) z = zarray + 15*ridx[i];
7998ab949d8SShri Abhyankar     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
8008ab949d8SShri 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;
8018ab949d8SShri Abhyankar 
8028ab949d8SShri Abhyankar     if (!usecprow) z += 15;
8038ab949d8SShri Abhyankar   }
8048ab949d8SShri Abhyankar 
8053649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
8068ab949d8SShri Abhyankar   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
8078ab949d8SShri Abhyankar   ierr = PetscLogFlops(450.0*a->nz - 15.0*nonzerorow);CHKERRQ(ierr);
8088ab949d8SShri Abhyankar   PetscFunctionReturn(0);
8098ab949d8SShri Abhyankar }
8108ab949d8SShri Abhyankar 
8118ab949d8SShri Abhyankar /* MatMult_SeqBAIJ_15_ver3 : Columns in the block are accessed in sets of 8,7 */
8128ab949d8SShri Abhyankar #undef __FUNCT__
8138ab949d8SShri Abhyankar #define __FUNCT__ "MatMult_SeqBAIJ_15_ver3"
8148ab949d8SShri Abhyankar PetscErrorCode MatMult_SeqBAIJ_15_ver3(Mat A,Vec xx,Vec zz)
8158ab949d8SShri Abhyankar {
8168ab949d8SShri Abhyankar   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
8178ab949d8SShri Abhyankar   PetscScalar       *z = 0,sum1,sum2,sum3,sum4,sum5,sum6,sum7,sum8,sum9,sum10,sum11,sum12,sum13,sum14,sum15;
8188ab949d8SShri Abhyankar   const PetscScalar *x,*xb;
8190b8f6341SShri Abhyankar   PetscScalar        x1,x2,x3,x4,x5,x6,x7,x8,*zarray;
8208ab949d8SShri Abhyankar   const MatScalar   *v;
8218ab949d8SShri Abhyankar   PetscErrorCode    ierr;
8228ab949d8SShri Abhyankar   const PetscInt    *ii,*ij=a->j,*idx;
8238ab949d8SShri Abhyankar   PetscInt          mbs,i,j,n,*ridx=PETSC_NULL,nonzerorow=0;
824ace3abfcSBarry Smith   PetscBool         usecprow=a->compressedrow.use;
8258ab949d8SShri Abhyankar 
8268ab949d8SShri Abhyankar   PetscFunctionBegin;
8273649974fSBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
8288ab949d8SShri Abhyankar   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
8298ab949d8SShri Abhyankar 
8308ab949d8SShri Abhyankar   v   = a->a;
8318ab949d8SShri Abhyankar   if (usecprow){
8328ab949d8SShri Abhyankar     mbs    = a->compressedrow.nrows;
8338ab949d8SShri Abhyankar     ii     = a->compressedrow.i;
8348ab949d8SShri Abhyankar     ridx = a->compressedrow.rindex;
8358ab949d8SShri Abhyankar   } else {
8368ab949d8SShri Abhyankar     mbs = a->mbs;
8378ab949d8SShri Abhyankar     ii  = a->i;
8388ab949d8SShri Abhyankar     z   = zarray;
8398ab949d8SShri Abhyankar   }
8408ab949d8SShri Abhyankar 
8418ab949d8SShri Abhyankar   for (i=0; i<mbs; i++) {
8428ab949d8SShri Abhyankar     n  = ii[i+1] - ii[i];
8438ab949d8SShri Abhyankar     idx = ij + ii[i];
8448ab949d8SShri Abhyankar     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0; sum7 = 0.0;
8458ab949d8SShri 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;
8468ab949d8SShri Abhyankar 
8478ab949d8SShri Abhyankar     nonzerorow += (n>0);
8488ab949d8SShri Abhyankar     for (j=0; j<n; j++) {
8498ab949d8SShri Abhyankar       xb = x + 15*(idx[j]);
8508ab949d8SShri Abhyankar       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6];
8510b8f6341SShri Abhyankar       x8 = xb[7];
8528ab949d8SShri Abhyankar 
8538ab949d8SShri 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;
8548ab949d8SShri 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;
8558ab949d8SShri 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;
8568ab949d8SShri 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;
8578ab949d8SShri 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;
8588ab949d8SShri 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;
8598ab949d8SShri 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;
8608ab949d8SShri 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;
8618ab949d8SShri 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;
8628ab949d8SShri 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;
8638ab949d8SShri 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;
8648ab949d8SShri 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;
8658ab949d8SShri 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;
8668ab949d8SShri 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;
8678ab949d8SShri 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;
8688ab949d8SShri Abhyankar       v += 120;
8698ab949d8SShri Abhyankar 
8700b8f6341SShri Abhyankar       x1 = xb[8]; x2 = xb[9]; x3 = xb[10]; x4 = xb[11]; x5 = xb[12]; x6 = xb[13]; x7 = xb[14];
8710b8f6341SShri Abhyankar 
8728ab949d8SShri Abhyankar       sum1  += v[0]*x1 + v[15]*x2  + v[30]*x3  + v[45]*x4 + v[60]*x5 + v[75]*x6 + v[90]*x7;
8738ab949d8SShri Abhyankar       sum2  += v[1]*x1 + v[16]*x2  + v[31]*x3  + v[46]*x4 + v[61]*x5 + v[76]*x6 + v[91]*x7;
8748ab949d8SShri Abhyankar       sum3  += v[2]*x1 + v[17]*x2  + v[32]*x3  + v[47]*x4 + v[62]*x5 + v[77]*x6 + v[92]*x7;
8758ab949d8SShri Abhyankar       sum4  += v[3]*x1 + v[18]*x2 + v[33]*x3  + v[48]*x4 + v[63]*x5 + v[78]*x6 + v[93]*x7;
8768ab949d8SShri Abhyankar       sum5  += v[4]*x1 + v[19]*x2 + v[34]*x3  + v[49]*x4 + v[64]*x5 + v[79]*x6 + v[94]*x7;
8778ab949d8SShri Abhyankar       sum6  += v[5]*x1 + v[20]*x2 + v[35]*x3  + v[50]*x4 + v[65]*x5 + v[80]*x6 + v[95]*x7;
8788ab949d8SShri Abhyankar       sum7  += v[6]*x1 + v[21]*x2 + v[36]*x3  + v[51]*x4 + v[66]*x5 + v[81]*x6 + v[96]*x7;
8798ab949d8SShri Abhyankar       sum8  += v[7]*x1 + v[22]*x2  + v[37]*x3  + v[52]*x4 + v[67]*x5 + v[82]*x6 + v[97]*x7;
8808ab949d8SShri Abhyankar       sum9  += v[8]*x1 + v[23]*x2  + v[38]*x3  + v[53]*x4 + v[68]*x5 + v[83]*x6 + v[98]*x7;
8818ab949d8SShri Abhyankar       sum10 += v[9]*x1 + v[24]*x2  + v[39]*x3  + v[54]*x4 + v[69]*x5 + v[84]*x6 + v[99]*x7;
8828ab949d8SShri Abhyankar       sum11 += v[10]*x1 + v[25]*x2 + v[40]*x3  + v[55]*x4 + v[70]*x5 + v[85]*x6 + v[100]*x7;
8838ab949d8SShri Abhyankar       sum12 += v[11]*x1 + v[26]*x2 + v[41]*x3  + v[56]*x4 + v[71]*x5 + v[86]*x6 + v[101]*x7;
8848ab949d8SShri Abhyankar       sum13 += v[12]*x1 + v[27]*x2 + v[42]*x3  + v[57]*x4 + v[72]*x5 + v[87]*x6 + v[102]*x7;
8858ab949d8SShri Abhyankar       sum14 += v[13]*x1 + v[28]*x2 + v[43]*x3  + v[58]*x4 + v[73]*x5 + v[88]*x6 + v[103]*x7;
8868ab949d8SShri Abhyankar       sum15 += v[14]*x1 + v[29]*x2 + v[44]*x3  + v[59]*x4 + v[74]*x5 + v[89]*x6 + v[104]*x7;
8878ab949d8SShri Abhyankar       v += 105;
8888ab949d8SShri Abhyankar     }
8898ab949d8SShri Abhyankar     if (usecprow) z = zarray + 15*ridx[i];
8908ab949d8SShri Abhyankar     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
8918ab949d8SShri 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;
8928ab949d8SShri Abhyankar 
8938ab949d8SShri Abhyankar     if (!usecprow) z += 15;
8948ab949d8SShri Abhyankar   }
8958ab949d8SShri Abhyankar 
8963649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
8978ab949d8SShri Abhyankar   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
8988ab949d8SShri Abhyankar   ierr = PetscLogFlops(450.0*a->nz - 15.0*nonzerorow);CHKERRQ(ierr);
8998ab949d8SShri Abhyankar   PetscFunctionReturn(0);
9008ab949d8SShri Abhyankar }
9018ab949d8SShri Abhyankar 
9028ab949d8SShri Abhyankar /* MatMult_SeqBAIJ_15_ver4 : All columns in the block are accessed at once */
9038ab949d8SShri Abhyankar 
9048ab949d8SShri Abhyankar #undef __FUNCT__
9058ab949d8SShri Abhyankar #define __FUNCT__ "MatMult_SeqBAIJ_15_ver4"
9068ab949d8SShri Abhyankar PetscErrorCode MatMult_SeqBAIJ_15_ver4(Mat A,Vec xx,Vec zz)
9078ab949d8SShri Abhyankar {
9088ab949d8SShri Abhyankar   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
9098ab949d8SShri Abhyankar   PetscScalar       *z = 0,sum1,sum2,sum3,sum4,sum5,sum6,sum7,sum8,sum9,sum10,sum11,sum12,sum13,sum14,sum15;
9108ab949d8SShri Abhyankar   const PetscScalar *x,*xb;
9118ab949d8SShri Abhyankar   PetscScalar        x1,x2,x3,x4,x5,x6,x7,x8,x9,x10,x11,x12,x13,x14,x15,*zarray;
9128ab949d8SShri Abhyankar   const MatScalar   *v;
9138ab949d8SShri Abhyankar   PetscErrorCode    ierr;
9148ab949d8SShri Abhyankar   const PetscInt    *ii,*ij=a->j,*idx;
9158ab949d8SShri Abhyankar   PetscInt          mbs,i,j,n,*ridx=PETSC_NULL,nonzerorow=0;
916ace3abfcSBarry Smith   PetscBool         usecprow=a->compressedrow.use;
9178ab949d8SShri Abhyankar 
9188ab949d8SShri Abhyankar   PetscFunctionBegin;
9193649974fSBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
9208ab949d8SShri Abhyankar   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
9218ab949d8SShri Abhyankar 
9228ab949d8SShri Abhyankar   v   = a->a;
9238ab949d8SShri Abhyankar   if (usecprow){
9248ab949d8SShri Abhyankar     mbs    = a->compressedrow.nrows;
9258ab949d8SShri Abhyankar     ii     = a->compressedrow.i;
9268ab949d8SShri Abhyankar     ridx = a->compressedrow.rindex;
9278ab949d8SShri Abhyankar   } else {
9288ab949d8SShri Abhyankar     mbs = a->mbs;
9298ab949d8SShri Abhyankar     ii  = a->i;
9308ab949d8SShri Abhyankar     z   = zarray;
9318ab949d8SShri Abhyankar   }
9328ab949d8SShri Abhyankar 
9338ab949d8SShri Abhyankar   for (i=0; i<mbs; i++) {
9348ab949d8SShri Abhyankar     n  = ii[i+1] - ii[i];
9358ab949d8SShri Abhyankar     idx = ij + ii[i];
9368ab949d8SShri Abhyankar     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0; sum7 = 0.0;
9378ab949d8SShri 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;
9388ab949d8SShri Abhyankar 
9398ab949d8SShri Abhyankar     nonzerorow += (n>0);
9408ab949d8SShri Abhyankar     for (j=0; j<n; j++) {
9418ab949d8SShri Abhyankar       xb = x + 15*(idx[j]);
9428ab949d8SShri Abhyankar       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6];
9438ab949d8SShri 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];
9448ab949d8SShri Abhyankar 
9458ab949d8SShri 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;
9468ab949d8SShri 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;
9478ab949d8SShri 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;
9488ab949d8SShri 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;
9498ab949d8SShri 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;
9508ab949d8SShri 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;
9518ab949d8SShri 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;
9528ab949d8SShri 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;
9538ab949d8SShri 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;
9548ab949d8SShri 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;
9558ab949d8SShri 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;
9568ab949d8SShri 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;
9578ab949d8SShri 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;
9588ab949d8SShri 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;
9598ab949d8SShri 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;
9608ab949d8SShri Abhyankar       v += 225;
9618ab949d8SShri Abhyankar     }
9628ab949d8SShri Abhyankar     if (usecprow) z = zarray + 15*ridx[i];
9638ab949d8SShri Abhyankar     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
9648ab949d8SShri 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;
9658ab949d8SShri Abhyankar 
9668ab949d8SShri Abhyankar     if (!usecprow) z += 15;
9678ab949d8SShri Abhyankar   }
9688ab949d8SShri Abhyankar 
9693649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
9708ab949d8SShri Abhyankar   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
9718ab949d8SShri Abhyankar   ierr = PetscLogFlops(450.0*a->nz - 15.0*nonzerorow);CHKERRQ(ierr);
9728ab949d8SShri Abhyankar   PetscFunctionReturn(0);
9738ab949d8SShri Abhyankar }
9748ab949d8SShri Abhyankar 
9758ab949d8SShri Abhyankar 
9763f1db9ecSBarry Smith /*
9773f1db9ecSBarry Smith     This will not work with MatScalar == float because it calls the BLAS
9783f1db9ecSBarry Smith */
9794a2ae208SSatish Balay #undef __FUNCT__
9804a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_N"
981dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_N(Mat A,Vec xx,Vec zz)
9822d61bbb3SSatish Balay {
9832d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
984dcf5cc72SBarry Smith   PetscScalar    *x,*z = 0,*xb,*work,*workt,*zarray;
9853f1db9ecSBarry Smith   MatScalar      *v;
986dfbe8321SBarry Smith   PetscErrorCode ierr;
987d0f46423SBarry Smith   PetscInt       mbs=a->mbs,i,*idx,*ii,bs=A->rmap->bs,j,n,bs2=a->bs2;
98898c9bda7SSatish Balay   PetscInt       ncols,k,*ridx=PETSC_NULL,nonzerorow=0;
989ace3abfcSBarry Smith   PetscBool      usecprow=a->compressedrow.use;
9902d61bbb3SSatish Balay 
9912d61bbb3SSatish Balay   PetscFunctionBegin;
9921ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
99326e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
9942d61bbb3SSatish Balay 
9952d61bbb3SSatish Balay   idx = a->j;
9962d61bbb3SSatish Balay   v   = a->a;
99726e093fcSHong Zhang   if (usecprow){
99826e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
99926e093fcSHong Zhang     ii   = a->compressedrow.i;
10007b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
100126e093fcSHong Zhang   } else {
100226e093fcSHong Zhang     mbs = a->mbs;
10032d61bbb3SSatish Balay     ii  = a->i;
100426e093fcSHong Zhang     z   = zarray;
100526e093fcSHong Zhang   }
1006218c64b6SSatish Balay 
10072d61bbb3SSatish Balay   if (!a->mult_work) {
1008d0f46423SBarry Smith     k    = PetscMax(A->rmap->n,A->cmap->n);
100987828ca2SBarry Smith     ierr = PetscMalloc((k+1)*sizeof(PetscScalar),&a->mult_work);CHKERRQ(ierr);
10102d61bbb3SSatish Balay   }
10112d61bbb3SSatish Balay   work = a->mult_work;
10122d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
10132d61bbb3SSatish Balay     n     = ii[1] - ii[0]; ii++;
10142d61bbb3SSatish Balay     ncols = n*bs;
10152d61bbb3SSatish Balay     workt = work;
101698c9bda7SSatish Balay     nonzerorow += (n>0);
10172d61bbb3SSatish Balay     for (j=0; j<n; j++) {
10182d61bbb3SSatish Balay       xb = x + bs*(*idx++);
10192d61bbb3SSatish Balay       for (k=0; k<bs; k++) workt[k] = xb[k];
10202d61bbb3SSatish Balay       workt += bs;
10212d61bbb3SSatish Balay     }
10227b2bb3b9SHong Zhang     if (usecprow) z = zarray + bs*ridx[i];
1023737d121aSSatish Balay     Kernel_w_gets_Ar_times_v(bs,ncols,work,v,z);
102471044d3cSBarry Smith     /* BLASgemv_("N",&bs,&ncols,&_DOne,v,&bs,work,&_One,&_DZero,z,&_One); */
10252d61bbb3SSatish Balay     v += n*bs2;
102626e093fcSHong Zhang     if (!usecprow) z += bs;
10272d61bbb3SSatish Balay   }
10281ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
102926e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
1030dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*a->nz*bs2 - bs*nonzerorow);CHKERRQ(ierr);
10312d61bbb3SSatish Balay   PetscFunctionReturn(0);
10322d61bbb3SSatish Balay }
10332d61bbb3SSatish Balay 
10344a2ae208SSatish Balay #undef __FUNCT__
10354a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_1"
1036dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_1(Mat A,Vec xx,Vec yy,Vec zz)
10372d61bbb3SSatish Balay {
10382d61bbb3SSatish Balay   Mat_SeqBAIJ        *a = (Mat_SeqBAIJ*)A->data;
1039122f12eaSBarry Smith   const PetscScalar  *x;
1040122f12eaSBarry Smith   PetscScalar        *y,*z,sum;
1041122f12eaSBarry Smith   const MatScalar    *v;
1042dfbe8321SBarry Smith   PetscErrorCode     ierr;
1043122f12eaSBarry Smith   PetscInt           mbs=a->mbs,i,n,*ridx=PETSC_NULL,nonzerorow=0;
1044122f12eaSBarry Smith   const PetscInt     *idx,*ii;
1045ace3abfcSBarry Smith   PetscBool          usecprow=a->compressedrow.use;
10462d61bbb3SSatish Balay 
10472d61bbb3SSatish Balay   PetscFunctionBegin;
10483649974fSBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
1049122f12eaSBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
10502e8a6d31SBarry Smith   if (zz != yy) {
1051122f12eaSBarry Smith     ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
10522e8a6d31SBarry Smith   } else {
1053122f12eaSBarry Smith     z = y;
10542e8a6d31SBarry Smith   }
10552d61bbb3SSatish Balay 
10562d61bbb3SSatish Balay   idx = a->j;
10572d61bbb3SSatish Balay   v   = a->a;
105826e093fcSHong Zhang   if (usecprow){
10594eb6d288SHong Zhang     if (zz != yy){
1060122f12eaSBarry Smith       ierr = PetscMemcpy(z,y,mbs*sizeof(PetscScalar));CHKERRQ(ierr);
10614eb6d288SHong Zhang     }
106226e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
106326e093fcSHong Zhang     ii   = a->compressedrow.i;
10647b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
106526e093fcSHong Zhang   } else {
10662d61bbb3SSatish Balay     ii  = a->i;
106726e093fcSHong Zhang   }
10682d61bbb3SSatish Balay 
10692d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
1070122f12eaSBarry Smith     n    = ii[1] - ii[0];
1071122f12eaSBarry Smith     ii++;
107226e093fcSHong Zhang     if (!usecprow){
1073122f12eaSBarry Smith       nonzerorow += (n>0);
1074122f12eaSBarry Smith       sum = y[i];
1075122f12eaSBarry Smith     } else {
1076122f12eaSBarry Smith       sum = y[ridx[i]];
1077122f12eaSBarry Smith     }
1078444d8c10SJed Brown     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA); /* Indices for the next row (assumes same size as this one) */
1079444d8c10SJed Brown     PetscPrefetchBlock(v+n,n,0,PETSC_PREFETCH_HINT_NTA);   /* Entries for the next row */
1080122f12eaSBarry Smith     PetscSparseDensePlusDot(sum,x,v,idx,n);
1081122f12eaSBarry Smith     v += n;
1082122f12eaSBarry Smith     idx += n;
1083122f12eaSBarry Smith     if (usecprow){
1084122f12eaSBarry Smith       z[ridx[i]] = sum;
1085122f12eaSBarry Smith     } else {
1086122f12eaSBarry Smith       z[i] = sum;
108726e093fcSHong Zhang     }
10882d61bbb3SSatish Balay   }
10893649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
1090122f12eaSBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
10912e8a6d31SBarry Smith   if (zz != yy) {
1092122f12eaSBarry Smith     ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
10932e8a6d31SBarry Smith   }
1094122f12eaSBarry Smith   ierr = PetscLogFlops(2.0*a->nz - nonzerorow);CHKERRQ(ierr);
10952d61bbb3SSatish Balay   PetscFunctionReturn(0);
10962d61bbb3SSatish Balay }
10972d61bbb3SSatish Balay 
10984a2ae208SSatish Balay #undef __FUNCT__
10994a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_2"
1100dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_2(Mat A,Vec xx,Vec yy,Vec zz)
11012d61bbb3SSatish Balay {
11022d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
1103dcf5cc72SBarry Smith   PetscScalar    *x,*y = 0,*z = 0,*xb,sum1,sum2;
110426e093fcSHong Zhang   PetscScalar    x1,x2,*yarray,*zarray;
11053f1db9ecSBarry Smith   MatScalar      *v;
1106dfbe8321SBarry Smith   PetscErrorCode ierr;
11074eb6d288SHong Zhang   PetscInt       mbs=a->mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL;
1108ace3abfcSBarry Smith   PetscBool      usecprow=a->compressedrow.use;
11092d61bbb3SSatish Balay 
11102d61bbb3SSatish Balay   PetscFunctionBegin;
11111ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
111226e093fcSHong Zhang   ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr);
11132e8a6d31SBarry Smith   if (zz != yy) {
111426e093fcSHong Zhang     ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
11152e8a6d31SBarry Smith   } else {
111626e093fcSHong Zhang     zarray = yarray;
11172e8a6d31SBarry Smith   }
11182d61bbb3SSatish Balay 
11192d61bbb3SSatish Balay   idx = a->j;
11202d61bbb3SSatish Balay   v   = a->a;
112126e093fcSHong Zhang   if (usecprow){
11224eb6d288SHong Zhang     if (zz != yy){
11234eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,2*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
11244eb6d288SHong Zhang     }
112526e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
112626e093fcSHong Zhang     ii   = a->compressedrow.i;
11277b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
11284eb6d288SHong Zhang     if (zz != yy){
11294eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,a->mbs*sizeof(PetscScalar));CHKERRQ(ierr);
11304eb6d288SHong Zhang     }
113126e093fcSHong Zhang   } else {
11322d61bbb3SSatish Balay     ii  = a->i;
113326e093fcSHong Zhang     y   = yarray;
113426e093fcSHong Zhang     z   = zarray;
113526e093fcSHong Zhang   }
11362d61bbb3SSatish Balay 
11372d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
11382d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
113926e093fcSHong Zhang     if (usecprow){
11407b2bb3b9SHong Zhang       z = zarray + 2*ridx[i];
11417b2bb3b9SHong Zhang       y = yarray + 2*ridx[i];
114226e093fcSHong Zhang     }
11432d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1];
1144444d8c10SJed Brown     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);   /* Indices for the next row (assumes same size as this one) */
1145444d8c10SJed Brown     PetscPrefetchBlock(v+4*n,4*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
11462d61bbb3SSatish Balay     for (j=0; j<n; j++) {
11472d61bbb3SSatish Balay       xb = x + 2*(*idx++); x1 = xb[0]; x2 = xb[1];
11482d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[2]*x2;
11492d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[3]*x2;
11502d61bbb3SSatish Balay       v += 4;
11512d61bbb3SSatish Balay     }
11522d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2;
115326e093fcSHong Zhang     if (!usecprow){
11542d61bbb3SSatish Balay       z += 2; y += 2;
11552d61bbb3SSatish Balay     }
115626e093fcSHong Zhang   }
11571ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
115826e093fcSHong Zhang   ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr);
11592e8a6d31SBarry Smith   if (zz != yy) {
116026e093fcSHong Zhang     ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
11612e8a6d31SBarry Smith   }
1162dc0b31edSSatish Balay   ierr = PetscLogFlops(4.0*a->nz);CHKERRQ(ierr);
11632d61bbb3SSatish Balay   PetscFunctionReturn(0);
11642d61bbb3SSatish Balay }
11652d61bbb3SSatish Balay 
11664a2ae208SSatish Balay #undef __FUNCT__
11674a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_3"
1168dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_3(Mat A,Vec xx,Vec yy,Vec zz)
11692d61bbb3SSatish Balay {
11702d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
1171dcf5cc72SBarry Smith   PetscScalar    *x,*y = 0,*z = 0,*xb,sum1,sum2,sum3,x1,x2,x3,*yarray,*zarray;
11723f1db9ecSBarry Smith   MatScalar      *v;
1173dfbe8321SBarry Smith   PetscErrorCode ierr;
11744eb6d288SHong Zhang   PetscInt       mbs=a->mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL;
1175ace3abfcSBarry Smith   PetscBool      usecprow=a->compressedrow.use;
11762d61bbb3SSatish Balay 
11772d61bbb3SSatish Balay   PetscFunctionBegin;
11781ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
117926e093fcSHong Zhang   ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr);
11802e8a6d31SBarry Smith   if (zz != yy) {
118126e093fcSHong Zhang     ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
11822e8a6d31SBarry Smith   } else {
118326e093fcSHong Zhang     zarray = yarray;
11842e8a6d31SBarry Smith   }
11852d61bbb3SSatish Balay 
11862d61bbb3SSatish Balay   idx = a->j;
11872d61bbb3SSatish Balay   v   = a->a;
118826e093fcSHong Zhang   if (usecprow){
11894eb6d288SHong Zhang     if (zz != yy){
11904eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,3*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
11914eb6d288SHong Zhang     }
119226e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
119326e093fcSHong Zhang     ii   = a->compressedrow.i;
11947b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
119526e093fcSHong Zhang   } else {
11962d61bbb3SSatish Balay     ii  = a->i;
119726e093fcSHong Zhang     y   = yarray;
119826e093fcSHong Zhang     z   = zarray;
119926e093fcSHong Zhang   }
12002d61bbb3SSatish Balay 
12012d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
12022d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
120326e093fcSHong Zhang     if (usecprow){
12047b2bb3b9SHong Zhang       z = zarray + 3*ridx[i];
12057b2bb3b9SHong Zhang       y = yarray + 3*ridx[i];
120626e093fcSHong Zhang     }
12072d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2];
1208444d8c10SJed Brown     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);   /* Indices for the next row (assumes same size as this one) */
1209444d8c10SJed Brown     PetscPrefetchBlock(v+9*n,9*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
12102d61bbb3SSatish Balay     for (j=0; j<n; j++) {
12112d61bbb3SSatish Balay       xb = x + 3*(*idx++); x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
12122d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[3]*x2 + v[6]*x3;
12132d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[4]*x2 + v[7]*x3;
12142d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[5]*x2 + v[8]*x3;
12152d61bbb3SSatish Balay       v += 9;
12162d61bbb3SSatish Balay     }
12172d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3;
121826e093fcSHong Zhang     if (!usecprow){
12192d61bbb3SSatish Balay       z += 3; y += 3;
12202d61bbb3SSatish Balay     }
122126e093fcSHong Zhang   }
12221ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
122326e093fcSHong Zhang   ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr);
12242e8a6d31SBarry Smith   if (zz != yy) {
122526e093fcSHong Zhang     ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
12262e8a6d31SBarry Smith   }
1227dc0b31edSSatish Balay   ierr = PetscLogFlops(18.0*a->nz);CHKERRQ(ierr);
12282d61bbb3SSatish Balay   PetscFunctionReturn(0);
12292d61bbb3SSatish Balay }
12302d61bbb3SSatish Balay 
12314a2ae208SSatish Balay #undef __FUNCT__
12324a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_4"
1233dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_4(Mat A,Vec xx,Vec yy,Vec zz)
12342d61bbb3SSatish Balay {
12352d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
1236dcf5cc72SBarry Smith   PetscScalar    *x,*y = 0,*z = 0,*xb,sum1,sum2,sum3,sum4,x1,x2,x3,x4,*yarray,*zarray;
12373f1db9ecSBarry Smith   MatScalar      *v;
1238dfbe8321SBarry Smith   PetscErrorCode ierr;
12394eb6d288SHong Zhang   PetscInt       mbs=a->mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL;
1240ace3abfcSBarry Smith   PetscBool      usecprow=a->compressedrow.use;
12412d61bbb3SSatish Balay 
12422d61bbb3SSatish Balay   PetscFunctionBegin;
12431ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
124426e093fcSHong Zhang   ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr);
12452e8a6d31SBarry Smith   if (zz != yy) {
124626e093fcSHong Zhang     ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
12472e8a6d31SBarry Smith   } else {
124826e093fcSHong Zhang     zarray = yarray;
12492e8a6d31SBarry Smith   }
12502d61bbb3SSatish Balay 
12512d61bbb3SSatish Balay   idx   = a->j;
12522d61bbb3SSatish Balay   v     = a->a;
125326e093fcSHong Zhang   if (usecprow){
12544eb6d288SHong Zhang     if (zz != yy){
12554eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,4*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
12564eb6d288SHong Zhang     }
125726e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
125826e093fcSHong Zhang     ii   = a->compressedrow.i;
12597b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
126026e093fcSHong Zhang   } else {
12612d61bbb3SSatish Balay     ii  = a->i;
126226e093fcSHong Zhang     y   = yarray;
126326e093fcSHong Zhang     z   = zarray;
126426e093fcSHong Zhang   }
12652d61bbb3SSatish Balay 
12662d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
12672d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
126826e093fcSHong Zhang     if (usecprow){
12697b2bb3b9SHong Zhang       z = zarray + 4*ridx[i];
12707b2bb3b9SHong Zhang       y = yarray + 4*ridx[i];
127126e093fcSHong Zhang     }
12722d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3];
1273444d8c10SJed Brown     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);     /* Indices for the next row (assumes same size as this one) */
1274444d8c10SJed Brown     PetscPrefetchBlock(v+16*n,16*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
12752d61bbb3SSatish Balay     for (j=0; j<n; j++) {
12762d61bbb3SSatish Balay       xb = x + 4*(*idx++);
12772d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
12782d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[4]*x2 + v[8]*x3   + v[12]*x4;
12792d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[5]*x2 + v[9]*x3   + v[13]*x4;
12802d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[6]*x2 + v[10]*x3  + v[14]*x4;
12812d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[7]*x2 + v[11]*x3  + v[15]*x4;
12822d61bbb3SSatish Balay       v += 16;
12832d61bbb3SSatish Balay     }
12842d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4;
128526e093fcSHong Zhang     if (!usecprow){
12862d61bbb3SSatish Balay       z += 4; y += 4;
12872d61bbb3SSatish Balay     }
128826e093fcSHong Zhang   }
12891ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
129026e093fcSHong Zhang   ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr);
12912e8a6d31SBarry Smith   if (zz != yy) {
129226e093fcSHong Zhang     ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
12932e8a6d31SBarry Smith   }
1294dc0b31edSSatish Balay   ierr = PetscLogFlops(32.0*a->nz);CHKERRQ(ierr);
12952d61bbb3SSatish Balay   PetscFunctionReturn(0);
12962d61bbb3SSatish Balay }
12972d61bbb3SSatish Balay 
12984a2ae208SSatish Balay #undef __FUNCT__
12994a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_5"
1300dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_5(Mat A,Vec xx,Vec yy,Vec zz)
13012d61bbb3SSatish Balay {
13022d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
1303dcf5cc72SBarry Smith   PetscScalar    *x,*y = 0,*z = 0,*xb,sum1,sum2,sum3,sum4,sum5,x1,x2,x3,x4,x5;
130426e093fcSHong Zhang   PetscScalar    *yarray,*zarray;
13053f1db9ecSBarry Smith   MatScalar      *v;
1306dfbe8321SBarry Smith   PetscErrorCode ierr;
13074eb6d288SHong Zhang   PetscInt       mbs=a->mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL;
1308ace3abfcSBarry Smith   PetscBool      usecprow=a->compressedrow.use;
13092d61bbb3SSatish Balay 
13102d61bbb3SSatish Balay   PetscFunctionBegin;
13111ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
131226e093fcSHong Zhang   ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr);
13132e8a6d31SBarry Smith   if (zz != yy) {
131426e093fcSHong Zhang     ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
13152e8a6d31SBarry Smith   } else {
131626e093fcSHong Zhang     zarray = yarray;
13172e8a6d31SBarry Smith   }
13182d61bbb3SSatish Balay 
13192d61bbb3SSatish Balay   idx = a->j;
13202d61bbb3SSatish Balay   v   = a->a;
132126e093fcSHong Zhang   if (usecprow){
13224eb6d288SHong Zhang     if (zz != yy){
13234eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,5*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
13244eb6d288SHong Zhang     }
132526e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
132626e093fcSHong Zhang     ii   = a->compressedrow.i;
13277b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
132826e093fcSHong Zhang   } else {
13292d61bbb3SSatish Balay     ii  = a->i;
133026e093fcSHong Zhang     y   = yarray;
133126e093fcSHong Zhang     z   = zarray;
133226e093fcSHong Zhang   }
13332d61bbb3SSatish Balay 
13342d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
13352d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
133626e093fcSHong Zhang     if (usecprow){
13377b2bb3b9SHong Zhang       z = zarray + 5*ridx[i];
13387b2bb3b9SHong Zhang       y = yarray + 5*ridx[i];
133926e093fcSHong Zhang     }
13402d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4];
1341444d8c10SJed Brown     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);     /* Indices for the next row (assumes same size as this one) */
1342444d8c10SJed Brown     PetscPrefetchBlock(v+25*n,25*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
13432d61bbb3SSatish Balay     for (j=0; j<n; j++) {
13442d61bbb3SSatish Balay       xb = x + 5*(*idx++);
13452d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4];
13462d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[5]*x2 + v[10]*x3  + v[15]*x4 + v[20]*x5;
13472d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[6]*x2 + v[11]*x3  + v[16]*x4 + v[21]*x5;
13482d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[7]*x2 + v[12]*x3  + v[17]*x4 + v[22]*x5;
13492d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[8]*x2 + v[13]*x3  + v[18]*x4 + v[23]*x5;
13502d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[9]*x2 + v[14]*x3  + v[19]*x4 + v[24]*x5;
13512d61bbb3SSatish Balay       v += 25;
13522d61bbb3SSatish Balay     }
13532d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5;
135426e093fcSHong Zhang     if (!usecprow){
13552d61bbb3SSatish Balay       z += 5; y += 5;
13562d61bbb3SSatish Balay     }
135726e093fcSHong Zhang   }
13581ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
135926e093fcSHong Zhang   ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr);
13602e8a6d31SBarry Smith   if (zz != yy) {
136126e093fcSHong Zhang     ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
13622e8a6d31SBarry Smith   }
1363dc0b31edSSatish Balay   ierr = PetscLogFlops(50.0*a->nz);CHKERRQ(ierr);
13642d61bbb3SSatish Balay   PetscFunctionReturn(0);
13652d61bbb3SSatish Balay }
13664a2ae208SSatish Balay #undef __FUNCT__
13674a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_6"
1368dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_6(Mat A,Vec xx,Vec yy,Vec zz)
136915091d37SBarry Smith {
137015091d37SBarry Smith   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
1371dcf5cc72SBarry Smith   PetscScalar    *x,*y = 0,*z = 0,*xb,sum1,sum2,sum3,sum4,sum5,sum6;
137226e093fcSHong Zhang   PetscScalar    x1,x2,x3,x4,x5,x6,*yarray,*zarray;
137315091d37SBarry Smith   MatScalar      *v;
1374dfbe8321SBarry Smith   PetscErrorCode ierr;
13754eb6d288SHong Zhang   PetscInt       mbs=a->mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL;
1376ace3abfcSBarry Smith   PetscBool      usecprow=a->compressedrow.use;
137715091d37SBarry Smith 
137815091d37SBarry Smith   PetscFunctionBegin;
13791ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
138026e093fcSHong Zhang   ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr);
138115091d37SBarry Smith   if (zz != yy) {
138226e093fcSHong Zhang     ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
138315091d37SBarry Smith   } else {
138426e093fcSHong Zhang     zarray = yarray;
138515091d37SBarry Smith   }
138615091d37SBarry Smith 
138715091d37SBarry Smith   idx = a->j;
138815091d37SBarry Smith   v   = a->a;
138926e093fcSHong Zhang   if (usecprow){
13904eb6d288SHong Zhang     if (zz != yy){
13914eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,6*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
13924eb6d288SHong Zhang     }
139326e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
139426e093fcSHong Zhang     ii   = a->compressedrow.i;
13957b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
139626e093fcSHong Zhang   } else {
139715091d37SBarry Smith     ii  = a->i;
139826e093fcSHong Zhang     y   = yarray;
139926e093fcSHong Zhang     z   = zarray;
140026e093fcSHong Zhang   }
140115091d37SBarry Smith 
140215091d37SBarry Smith   for (i=0; i<mbs; i++) {
140315091d37SBarry Smith     n  = ii[1] - ii[0]; ii++;
140426e093fcSHong Zhang     if (usecprow){
14057b2bb3b9SHong Zhang       z = zarray + 6*ridx[i];
14067b2bb3b9SHong Zhang       y = yarray + 6*ridx[i];
140726e093fcSHong Zhang     }
140815091d37SBarry Smith     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4]; sum6 = y[5];
1409444d8c10SJed Brown     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);     /* Indices for the next row (assumes same size as this one) */
1410444d8c10SJed Brown     PetscPrefetchBlock(v+36*n,36*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
141115091d37SBarry Smith     for (j=0; j<n; j++) {
14123b95cb0eSSatish Balay       xb = x + 6*(*idx++);
141315091d37SBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5];
141415091d37SBarry Smith       sum1 += v[0]*x1 + v[6]*x2  + v[12]*x3  + v[18]*x4 + v[24]*x5 + v[30]*x6;
141515091d37SBarry Smith       sum2 += v[1]*x1 + v[7]*x2  + v[13]*x3  + v[19]*x4 + v[25]*x5 + v[31]*x6;
141615091d37SBarry Smith       sum3 += v[2]*x1 + v[8]*x2  + v[14]*x3  + v[20]*x4 + v[26]*x5 + v[32]*x6;
141715091d37SBarry Smith       sum4 += v[3]*x1 + v[9]*x2  + v[15]*x3  + v[21]*x4 + v[27]*x5 + v[33]*x6;
141815091d37SBarry Smith       sum5 += v[4]*x1 + v[10]*x2 + v[16]*x3  + v[22]*x4 + v[28]*x5 + v[34]*x6;
141915091d37SBarry Smith       sum6 += v[5]*x1 + v[11]*x2 + v[17]*x3  + v[23]*x4 + v[29]*x5 + v[35]*x6;
142015091d37SBarry Smith       v += 36;
142115091d37SBarry Smith     }
142215091d37SBarry Smith     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6;
142326e093fcSHong Zhang     if (!usecprow){
142415091d37SBarry Smith       z += 6; y += 6;
142515091d37SBarry Smith     }
142626e093fcSHong Zhang   }
14271ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
142826e093fcSHong Zhang   ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr);
142915091d37SBarry Smith   if (zz != yy) {
143026e093fcSHong Zhang     ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
143115091d37SBarry Smith   }
1432dc0b31edSSatish Balay   ierr = PetscLogFlops(72.0*a->nz);CHKERRQ(ierr);
143315091d37SBarry Smith   PetscFunctionReturn(0);
143415091d37SBarry Smith }
14352d61bbb3SSatish Balay 
14364a2ae208SSatish Balay #undef __FUNCT__
14374a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_7"
1438dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_7(Mat A,Vec xx,Vec yy,Vec zz)
14392d61bbb3SSatish Balay {
14402d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
1441dcf5cc72SBarry Smith   PetscScalar    *x,*y = 0,*z = 0,*xb,sum1,sum2,sum3,sum4,sum5,sum6,sum7;
144226e093fcSHong Zhang   PetscScalar    x1,x2,x3,x4,x5,x6,x7,*yarray,*zarray;
14433f1db9ecSBarry Smith   MatScalar      *v;
1444dfbe8321SBarry Smith   PetscErrorCode ierr;
14457b2bb3b9SHong Zhang   PetscInt       mbs=a->mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL;
1446ace3abfcSBarry Smith   PetscBool      usecprow=a->compressedrow.use;
14472d61bbb3SSatish Balay 
14482d61bbb3SSatish Balay   PetscFunctionBegin;
14491ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
145026e093fcSHong Zhang   ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr);
14512e8a6d31SBarry Smith   if (zz != yy) {
145226e093fcSHong Zhang     ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
14532e8a6d31SBarry Smith   } else {
145426e093fcSHong Zhang     zarray = yarray;
14552e8a6d31SBarry Smith   }
14562d61bbb3SSatish Balay 
14572d61bbb3SSatish Balay   idx = a->j;
14582d61bbb3SSatish Balay   v   = a->a;
145926e093fcSHong Zhang   if (usecprow){
14604eb6d288SHong Zhang     if (zz != yy){
14614eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,7*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
14624eb6d288SHong Zhang     }
146326e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
146426e093fcSHong Zhang     ii   = a->compressedrow.i;
14657b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
146626e093fcSHong Zhang   } else {
14672d61bbb3SSatish Balay     ii  = a->i;
146826e093fcSHong Zhang     y   = yarray;
146926e093fcSHong Zhang     z   = zarray;
147026e093fcSHong Zhang   }
14712d61bbb3SSatish Balay 
14722d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
14732d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
147426e093fcSHong Zhang     if (usecprow){
14757b2bb3b9SHong Zhang       z = zarray + 7*ridx[i];
14767b2bb3b9SHong Zhang       y = yarray + 7*ridx[i];
147726e093fcSHong Zhang     }
14782d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4]; sum6 = y[5]; sum7 = y[6];
1479444d8c10SJed Brown     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);     /* Indices for the next row (assumes same size as this one) */
1480444d8c10SJed Brown     PetscPrefetchBlock(v+49*n,49*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
14812d61bbb3SSatish Balay     for (j=0; j<n; j++) {
14822d61bbb3SSatish Balay       xb = x + 7*(*idx++);
14832d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6];
14842d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[7]*x2  + v[14]*x3  + v[21]*x4 + v[28]*x5 + v[35]*x6 + v[42]*x7;
14852d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[8]*x2  + v[15]*x3  + v[22]*x4 + v[29]*x5 + v[36]*x6 + v[43]*x7;
14862d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[9]*x2  + v[16]*x3  + v[23]*x4 + v[30]*x5 + v[37]*x6 + v[44]*x7;
14872d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[10]*x2 + v[17]*x3  + v[24]*x4 + v[31]*x5 + v[38]*x6 + v[45]*x7;
14882d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[11]*x2 + v[18]*x3  + v[25]*x4 + v[32]*x5 + v[39]*x6 + v[46]*x7;
14892d61bbb3SSatish Balay       sum6 += v[5]*x1 + v[12]*x2 + v[19]*x3  + v[26]*x4 + v[33]*x5 + v[40]*x6 + v[47]*x7;
14902d61bbb3SSatish Balay       sum7 += v[6]*x1 + v[13]*x2 + v[20]*x3  + v[27]*x4 + v[34]*x5 + v[41]*x6 + v[48]*x7;
14912d61bbb3SSatish Balay       v += 49;
14922d61bbb3SSatish Balay     }
14932d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
149426e093fcSHong Zhang     if (!usecprow){
14952d61bbb3SSatish Balay       z += 7; y += 7;
14962d61bbb3SSatish Balay     }
149726e093fcSHong Zhang   }
14981ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
149926e093fcSHong Zhang   ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr);
15002e8a6d31SBarry Smith   if (zz != yy) {
150126e093fcSHong Zhang     ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
15022e8a6d31SBarry Smith   }
1503dc0b31edSSatish Balay   ierr = PetscLogFlops(98.0*a->nz);CHKERRQ(ierr);
15042d61bbb3SSatish Balay   PetscFunctionReturn(0);
15052d61bbb3SSatish Balay }
1506218c64b6SSatish Balay 
15074a2ae208SSatish Balay #undef __FUNCT__
15084a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_N"
1509dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_N(Mat A,Vec xx,Vec yy,Vec zz)
15102d61bbb3SSatish Balay {
15112d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
1512bba805e6SBarry Smith   PetscScalar    *x,*z = 0,*xb,*work,*workt,*zarray;
15133f1db9ecSBarry Smith   MatScalar      *v;
15146849ba73SBarry Smith   PetscErrorCode ierr;
1515d0f46423SBarry Smith   PetscInt       mbs,i,*idx,*ii,bs=A->rmap->bs,j,n,bs2=a->bs2;
15167b2bb3b9SHong Zhang   PetscInt       ncols,k,*ridx=PETSC_NULL;
1517ace3abfcSBarry Smith   PetscBool      usecprow=a->compressedrow.use;
1518218c64b6SSatish Balay 
15192d61bbb3SSatish Balay   PetscFunctionBegin;
1520343551c4SBarry Smith   ierr = VecCopy(yy,zz);CHKERRQ(ierr);
15211ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
152226e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
15232d61bbb3SSatish Balay 
15242d61bbb3SSatish Balay   idx = a->j;
15252d61bbb3SSatish Balay   v   = a->a;
152626e093fcSHong Zhang   if (usecprow){
152726e093fcSHong Zhang     mbs    = a->compressedrow.nrows;
152826e093fcSHong Zhang     ii     = a->compressedrow.i;
15297b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
153026e093fcSHong Zhang   } else {
153126e093fcSHong Zhang     mbs = a->mbs;
15322d61bbb3SSatish Balay     ii  = a->i;
153326e093fcSHong Zhang     z   = zarray;
153426e093fcSHong Zhang   }
15352d61bbb3SSatish Balay 
15362d61bbb3SSatish Balay   if (!a->mult_work) {
1537d0f46423SBarry Smith     k    = PetscMax(A->rmap->n,A->cmap->n);
153887828ca2SBarry Smith     ierr = PetscMalloc((k+1)*sizeof(PetscScalar),&a->mult_work);CHKERRQ(ierr);
15392d61bbb3SSatish Balay   }
15402d61bbb3SSatish Balay   work = a->mult_work;
15412d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
15422d61bbb3SSatish Balay     n     = ii[1] - ii[0]; ii++;
15432d61bbb3SSatish Balay     ncols = n*bs;
15442d61bbb3SSatish Balay     workt = work;
15452d61bbb3SSatish Balay     for (j=0; j<n; j++) {
15462d61bbb3SSatish Balay       xb = x + bs*(*idx++);
15472d61bbb3SSatish Balay       for (k=0; k<bs; k++) workt[k] = xb[k];
15482d61bbb3SSatish Balay       workt += bs;
15492d61bbb3SSatish Balay     }
15507b2bb3b9SHong Zhang     if (usecprow) z = zarray + bs*ridx[i];
15513f1db9ecSBarry Smith     Kernel_w_gets_w_plus_Ar_times_v(bs,ncols,work,v,z);
155271044d3cSBarry Smith     /* BLASgemv_("N",&bs,&ncols,&_DOne,v,&bs,work,&_One,&_DOne,z,&_One); */
15532d61bbb3SSatish Balay     v += n*bs2;
155426e093fcSHong Zhang     if (!usecprow){
15552d61bbb3SSatish Balay       z += bs;
15562d61bbb3SSatish Balay     }
155726e093fcSHong Zhang   }
15581ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
155926e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
1560dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*a->nz*bs2);CHKERRQ(ierr);
15612d61bbb3SSatish Balay   PetscFunctionReturn(0);
15622d61bbb3SSatish Balay }
15632d61bbb3SSatish Balay 
15644a2ae208SSatish Balay #undef __FUNCT__
1565547795f9SHong Zhang #define __FUNCT__ "MatMultHermitianTranspose_SeqBAIJ"
1566547795f9SHong Zhang PetscErrorCode MatMultHermitianTranspose_SeqBAIJ(Mat A,Vec xx,Vec zz)
1567547795f9SHong Zhang {
1568547795f9SHong Zhang   PetscScalar    zero = 0.0;
1569547795f9SHong Zhang   PetscErrorCode ierr;
1570547795f9SHong Zhang 
1571547795f9SHong Zhang   PetscFunctionBegin;
1572547795f9SHong Zhang   ierr = VecSet(zz,zero);CHKERRQ(ierr);
1573547795f9SHong Zhang   ierr = MatMultHermitianTransposeAdd_SeqBAIJ(A,xx,zz,zz);CHKERRQ(ierr);
1574547795f9SHong Zhang   PetscFunctionReturn(0);
1575547795f9SHong Zhang }
1576547795f9SHong Zhang 
1577547795f9SHong Zhang #undef __FUNCT__
15784a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_SeqBAIJ"
1579dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_SeqBAIJ(Mat A,Vec xx,Vec zz)
15802d61bbb3SSatish Balay {
15813447b6efSHong Zhang   PetscScalar    zero = 0.0;
15826849ba73SBarry Smith   PetscErrorCode ierr;
15832d61bbb3SSatish Balay 
15842d61bbb3SSatish Balay   PetscFunctionBegin;
15852dcb1b2aSMatthew Knepley   ierr = VecSet(zz,zero);CHKERRQ(ierr);
15863447b6efSHong Zhang   ierr = MatMultTransposeAdd_SeqBAIJ(A,xx,zz,zz);CHKERRQ(ierr);
15872d61bbb3SSatish Balay   PetscFunctionReturn(0);
15882d61bbb3SSatish Balay }
15892d61bbb3SSatish Balay 
15904a2ae208SSatish Balay #undef __FUNCT__
1591547795f9SHong Zhang #define __FUNCT__ "MatMultHermitianTransposeAdd_SeqBAIJ"
1592547795f9SHong Zhang PetscErrorCode MatMultHermitianTransposeAdd_SeqBAIJ(Mat A,Vec xx,Vec yy,Vec zz)
1593547795f9SHong Zhang 
1594547795f9SHong Zhang {
1595547795f9SHong Zhang   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
1596547795f9SHong Zhang   PetscScalar       *zb,*x,*z,*xb = 0,x1,x2,x3,x4,x5;
1597547795f9SHong Zhang   MatScalar         *v;
1598547795f9SHong Zhang   PetscErrorCode    ierr;
1599547795f9SHong Zhang   PetscInt          mbs,i,*idx,*ii,rval,bs=A->rmap->bs,j,n,bs2=a->bs2,*ib,*ridx=PETSC_NULL;
1600547795f9SHong Zhang   Mat_CompressedRow cprow = a->compressedrow;
1601ace3abfcSBarry Smith   PetscBool         usecprow=cprow.use;
1602547795f9SHong Zhang 
1603547795f9SHong Zhang   PetscFunctionBegin;
1604343551c4SBarry Smith   if (yy != zz) { ierr = VecCopy(yy,zz);CHKERRQ(ierr); }
1605547795f9SHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1606547795f9SHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
1607547795f9SHong Zhang 
1608547795f9SHong Zhang   idx = a->j;
1609547795f9SHong Zhang   v   = a->a;
1610547795f9SHong Zhang   if (usecprow){
1611547795f9SHong Zhang     mbs  = cprow.nrows;
1612547795f9SHong Zhang     ii   = cprow.i;
1613547795f9SHong Zhang     ridx = cprow.rindex;
1614547795f9SHong Zhang   } else {
1615547795f9SHong Zhang     mbs=a->mbs;
1616547795f9SHong Zhang     ii = a->i;
1617547795f9SHong Zhang     xb = x;
1618547795f9SHong Zhang   }
1619547795f9SHong Zhang 
1620547795f9SHong Zhang   switch (bs) {
1621547795f9SHong Zhang   case 1:
1622547795f9SHong Zhang     for (i=0; i<mbs; i++) {
1623547795f9SHong Zhang       if (usecprow) xb = x + ridx[i];
1624547795f9SHong Zhang       x1 = xb[0];
1625547795f9SHong Zhang       ib = idx + ii[0];
1626547795f9SHong Zhang       n  = ii[1] - ii[0]; ii++;
1627547795f9SHong Zhang       for (j=0; j<n; j++) {
1628547795f9SHong Zhang         rval    = ib[j];
1629547795f9SHong Zhang         z[rval] += PetscConj(*v) * x1;
1630547795f9SHong Zhang         v++;
1631547795f9SHong Zhang       }
1632547795f9SHong Zhang       if (!usecprow) xb++;
1633547795f9SHong Zhang     }
1634547795f9SHong Zhang     break;
1635547795f9SHong Zhang   case 2:
1636547795f9SHong Zhang     for (i=0; i<mbs; i++) {
1637547795f9SHong Zhang       if (usecprow) xb = x + 2*ridx[i];
1638547795f9SHong Zhang       x1 = xb[0]; x2 = xb[1];
1639547795f9SHong Zhang       ib = idx + ii[0];
1640547795f9SHong Zhang       n  = ii[1] - ii[0]; ii++;
1641547795f9SHong Zhang       for (j=0; j<n; j++) {
1642547795f9SHong Zhang         rval      = ib[j]*2;
1643547795f9SHong Zhang         z[rval++] += PetscConj(v[0])*x1 + PetscConj(v[1])*x2;
1644547795f9SHong Zhang         z[rval++] += PetscConj(v[2])*x1 + PetscConj(v[3])*x2;
1645547795f9SHong Zhang         v  += 4;
1646547795f9SHong Zhang       }
1647547795f9SHong Zhang       if (!usecprow) xb += 2;
1648547795f9SHong Zhang     }
1649547795f9SHong Zhang     break;
1650547795f9SHong Zhang   case 3:
1651547795f9SHong Zhang     for (i=0; i<mbs; i++) {
1652547795f9SHong Zhang       if (usecprow) xb = x + 3*ridx[i];
1653547795f9SHong Zhang       x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
1654547795f9SHong Zhang       ib = idx + ii[0];
1655547795f9SHong Zhang       n  = ii[1] - ii[0]; ii++;
1656547795f9SHong Zhang       for (j=0; j<n; j++) {
1657547795f9SHong Zhang         rval      = ib[j]*3;
1658547795f9SHong Zhang         z[rval++] += PetscConj(v[0])*x1 + PetscConj(v[1])*x2 + PetscConj(v[2])*x3;
1659547795f9SHong Zhang         z[rval++] += PetscConj(v[3])*x1 + PetscConj(v[4])*x2 + PetscConj(v[5])*x3;
1660547795f9SHong Zhang         z[rval++] += PetscConj(v[6])*x1 + PetscConj(v[7])*x2 + PetscConj(v[8])*x3;
1661547795f9SHong Zhang         v  += 9;
1662547795f9SHong Zhang       }
1663547795f9SHong Zhang       if (!usecprow) xb += 3;
1664547795f9SHong Zhang     }
1665547795f9SHong Zhang     break;
1666547795f9SHong Zhang   case 4:
1667547795f9SHong Zhang     for (i=0; i<mbs; i++) {
1668547795f9SHong Zhang       if (usecprow) xb = x + 4*ridx[i];
1669547795f9SHong Zhang       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
1670547795f9SHong Zhang       ib = idx + ii[0];
1671547795f9SHong Zhang       n  = ii[1] - ii[0]; ii++;
1672547795f9SHong Zhang       for (j=0; j<n; j++) {
1673547795f9SHong Zhang         rval      = ib[j]*4;
1674547795f9SHong Zhang         z[rval++] +=  PetscConj(v[0])*x1 + PetscConj(v[1])*x2  + PetscConj(v[2])*x3  + PetscConj(v[3])*x4;
1675547795f9SHong Zhang         z[rval++] +=  PetscConj(v[4])*x1 + PetscConj(v[5])*x2  + PetscConj(v[6])*x3  + PetscConj(v[7])*x4;
1676547795f9SHong Zhang         z[rval++] +=  PetscConj(v[8])*x1 + PetscConj(v[9])*x2  + PetscConj(v[10])*x3 + PetscConj(v[11])*x4;
1677547795f9SHong Zhang         z[rval++] += PetscConj(v[12])*x1 + PetscConj(v[13])*x2 + PetscConj(v[14])*x3 + PetscConj(v[15])*x4;
1678547795f9SHong Zhang         v  += 16;
1679547795f9SHong Zhang       }
1680547795f9SHong Zhang       if (!usecprow) xb += 4;
1681547795f9SHong Zhang     }
1682547795f9SHong Zhang     break;
1683547795f9SHong Zhang   case 5:
1684547795f9SHong Zhang     for (i=0; i<mbs; i++) {
1685547795f9SHong Zhang       if (usecprow) xb = x + 5*ridx[i];
1686547795f9SHong Zhang       x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
1687547795f9SHong Zhang       x4 = xb[3]; x5 = xb[4];
1688547795f9SHong Zhang       ib = idx + ii[0];
1689547795f9SHong Zhang       n  = ii[1] - ii[0]; ii++;
1690547795f9SHong Zhang       for (j=0; j<n; j++) {
1691547795f9SHong Zhang         rval      = ib[j]*5;
1692547795f9SHong Zhang         z[rval++] +=  PetscConj(v[0])*x1 +  PetscConj(v[1])*x2 +  PetscConj(v[2])*x3 +  PetscConj(v[3])*x4 +  PetscConj(v[4])*x5;
1693547795f9SHong Zhang         z[rval++] +=  PetscConj(v[5])*x1 +  PetscConj(v[6])*x2 +  PetscConj(v[7])*x3 +  PetscConj(v[8])*x4 +  PetscConj(v[9])*x5;
1694547795f9SHong Zhang         z[rval++] += PetscConj(v[10])*x1 + PetscConj(v[11])*x2 + PetscConj(v[12])*x3 + PetscConj(v[13])*x4 + PetscConj(v[14])*x5;
1695547795f9SHong Zhang         z[rval++] += PetscConj(v[15])*x1 + PetscConj(v[16])*x2 + PetscConj(v[17])*x3 + PetscConj(v[18])*x4 + PetscConj(v[19])*x5;
1696547795f9SHong Zhang         z[rval++] += PetscConj(v[20])*x1 + PetscConj(v[21])*x2 + PetscConj(v[22])*x3 + PetscConj(v[23])*x4 + PetscConj(v[24])*x5;
1697547795f9SHong Zhang         v  += 25;
1698547795f9SHong Zhang       }
1699547795f9SHong Zhang       if (!usecprow) xb += 5;
1700547795f9SHong Zhang     }
1701547795f9SHong Zhang     break;
1702547795f9SHong Zhang   default: {      /* block sizes larger than 5 by 5 are handled by BLAS */
1703547795f9SHong Zhang       PetscInt     ncols,k;
1704547795f9SHong Zhang       PetscScalar  *work,*workt,*xtmp;
1705547795f9SHong Zhang 
1706e32f2f54SBarry Smith       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"block size larger than 5 is not supported yet");
1707547795f9SHong Zhang       if (!a->mult_work) {
1708547795f9SHong Zhang         k = PetscMax(A->rmap->n,A->cmap->n);
1709547795f9SHong Zhang         ierr = PetscMalloc((k+1)*sizeof(PetscScalar),&a->mult_work);CHKERRQ(ierr);
1710547795f9SHong Zhang       }
1711547795f9SHong Zhang       work = a->mult_work;
1712547795f9SHong Zhang       xtmp = x;
1713547795f9SHong Zhang       for (i=0; i<mbs; i++) {
1714547795f9SHong Zhang         n     = ii[1] - ii[0]; ii++;
1715547795f9SHong Zhang         ncols = n*bs;
1716547795f9SHong Zhang         ierr  = PetscMemzero(work,ncols*sizeof(PetscScalar));CHKERRQ(ierr);
1717547795f9SHong Zhang         if (usecprow) {
1718547795f9SHong Zhang           xtmp = x + bs*ridx[i];
1719547795f9SHong Zhang         }
1720547795f9SHong Zhang         Kernel_w_gets_w_plus_trans_Ar_times_v(bs,ncols,xtmp,v,work);
1721547795f9SHong Zhang         /* BLASgemv_("T",&bs,&ncols,&_DOne,v,&bs,xtmp,&_One,&_DOne,work,&_One); */
1722547795f9SHong Zhang         v += n*bs2;
1723547795f9SHong Zhang         if (!usecprow) xtmp += bs;
1724547795f9SHong Zhang         workt = work;
1725547795f9SHong Zhang         for (j=0; j<n; j++) {
1726547795f9SHong Zhang           zb = z + bs*(*idx++);
1727547795f9SHong Zhang           for (k=0; k<bs; k++) zb[k] += workt[k] ;
1728547795f9SHong Zhang           workt += bs;
1729547795f9SHong Zhang         }
1730547795f9SHong Zhang       }
1731547795f9SHong Zhang     }
1732547795f9SHong Zhang   }
1733547795f9SHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
1734547795f9SHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
1735547795f9SHong Zhang   ierr = PetscLogFlops(2.0*a->nz*a->bs2);CHKERRQ(ierr);
1736547795f9SHong Zhang   PetscFunctionReturn(0);
1737547795f9SHong Zhang }
1738547795f9SHong Zhang 
1739547795f9SHong Zhang #undef __FUNCT__
17404a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_SeqBAIJ"
1741dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_SeqBAIJ(Mat A,Vec xx,Vec yy,Vec zz)
17422d61bbb3SSatish Balay {
17432d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
1744dcf5cc72SBarry Smith   PetscScalar       *zb,*x,*z,*xb = 0,x1,x2,x3,x4,x5;
17453f1db9ecSBarry Smith   MatScalar         *v;
17466849ba73SBarry Smith   PetscErrorCode    ierr;
1747d0f46423SBarry Smith   PetscInt          mbs,i,*idx,*ii,rval,bs=A->rmap->bs,j,n,bs2=a->bs2,*ib,*ridx=PETSC_NULL;
17483447b6efSHong Zhang   Mat_CompressedRow cprow = a->compressedrow;
1749ace3abfcSBarry Smith   PetscBool         usecprow=cprow.use;
17502d61bbb3SSatish Balay 
17512d61bbb3SSatish Balay   PetscFunctionBegin;
1752343551c4SBarry Smith   if (yy != zz) { ierr = VecCopy(yy,zz);CHKERRQ(ierr); }
17533447b6efSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
17543447b6efSHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
17552d61bbb3SSatish Balay 
17562d61bbb3SSatish Balay   idx = a->j;
17572d61bbb3SSatish Balay   v   = a->a;
17583447b6efSHong Zhang   if (usecprow){
17593447b6efSHong Zhang     mbs  = cprow.nrows;
17603447b6efSHong Zhang     ii   = cprow.i;
17617b2bb3b9SHong Zhang     ridx = cprow.rindex;
17623447b6efSHong Zhang   } else {
17633447b6efSHong Zhang     mbs=a->mbs;
17642d61bbb3SSatish Balay     ii = a->i;
1765f1af5d2fSBarry Smith     xb = x;
17663447b6efSHong Zhang   }
17672d61bbb3SSatish Balay 
17682d61bbb3SSatish Balay   switch (bs) {
17692d61bbb3SSatish Balay   case 1:
17702d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
17717b2bb3b9SHong Zhang       if (usecprow) xb = x + ridx[i];
1772f1af5d2fSBarry Smith       x1 = xb[0];
17733447b6efSHong Zhang       ib = idx + ii[0];
17743447b6efSHong Zhang       n  = ii[1] - ii[0]; ii++;
17752d61bbb3SSatish Balay       for (j=0; j<n; j++) {
17762d61bbb3SSatish Balay         rval    = ib[j];
1777f1af5d2fSBarry Smith         z[rval] += *v * x1;
1778f1af5d2fSBarry Smith         v++;
17792d61bbb3SSatish Balay       }
17803447b6efSHong Zhang       if (!usecprow) xb++;
17812d61bbb3SSatish Balay     }
17822d61bbb3SSatish Balay     break;
17832d61bbb3SSatish Balay   case 2:
17842d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
17857b2bb3b9SHong Zhang       if (usecprow) xb = x + 2*ridx[i];
1786f1af5d2fSBarry Smith       x1 = xb[0]; x2 = xb[1];
17873447b6efSHong Zhang       ib = idx + ii[0];
17883447b6efSHong Zhang       n  = ii[1] - ii[0]; ii++;
17892d61bbb3SSatish Balay       for (j=0; j<n; j++) {
17902d61bbb3SSatish Balay         rval      = ib[j]*2;
17912d61bbb3SSatish Balay         z[rval++] += v[0]*x1 + v[1]*x2;
17922d61bbb3SSatish Balay         z[rval++] += v[2]*x1 + v[3]*x2;
17932d61bbb3SSatish Balay         v  += 4;
17942d61bbb3SSatish Balay       }
17953447b6efSHong Zhang       if (!usecprow) xb += 2;
17962d61bbb3SSatish Balay     }
17972d61bbb3SSatish Balay     break;
17982d61bbb3SSatish Balay   case 3:
17992d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
18007b2bb3b9SHong Zhang       if (usecprow) xb = x + 3*ridx[i];
1801f1af5d2fSBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
18023447b6efSHong Zhang       ib = idx + ii[0];
18033447b6efSHong Zhang       n  = ii[1] - ii[0]; ii++;
18042d61bbb3SSatish Balay       for (j=0; j<n; j++) {
18052d61bbb3SSatish Balay         rval      = ib[j]*3;
18062d61bbb3SSatish Balay         z[rval++] += v[0]*x1 + v[1]*x2 + v[2]*x3;
18072d61bbb3SSatish Balay         z[rval++] += v[3]*x1 + v[4]*x2 + v[5]*x3;
18082d61bbb3SSatish Balay         z[rval++] += v[6]*x1 + v[7]*x2 + v[8]*x3;
18092d61bbb3SSatish Balay         v  += 9;
18102d61bbb3SSatish Balay       }
18113447b6efSHong Zhang       if (!usecprow) xb += 3;
18122d61bbb3SSatish Balay     }
18132d61bbb3SSatish Balay     break;
18142d61bbb3SSatish Balay   case 4:
18152d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
18167b2bb3b9SHong Zhang       if (usecprow) xb = x + 4*ridx[i];
1817f1af5d2fSBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
18183447b6efSHong Zhang       ib = idx + ii[0];
18193447b6efSHong Zhang       n  = ii[1] - ii[0]; ii++;
18202d61bbb3SSatish Balay       for (j=0; j<n; j++) {
18212d61bbb3SSatish Balay         rval      = ib[j]*4;
18222d61bbb3SSatish Balay         z[rval++] +=  v[0]*x1 +  v[1]*x2 +  v[2]*x3 +  v[3]*x4;
18232d61bbb3SSatish Balay         z[rval++] +=  v[4]*x1 +  v[5]*x2 +  v[6]*x3 +  v[7]*x4;
18242d61bbb3SSatish Balay         z[rval++] +=  v[8]*x1 +  v[9]*x2 + v[10]*x3 + v[11]*x4;
18252d61bbb3SSatish Balay         z[rval++] += v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4;
18262d61bbb3SSatish Balay         v  += 16;
18272d61bbb3SSatish Balay       }
18283447b6efSHong Zhang       if (!usecprow) xb += 4;
18292d61bbb3SSatish Balay     }
18302d61bbb3SSatish Balay     break;
18312d61bbb3SSatish Balay   case 5:
18322d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
18337b2bb3b9SHong Zhang       if (usecprow) xb = x + 5*ridx[i];
1834f1af5d2fSBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
18352d61bbb3SSatish Balay       x4 = xb[3]; x5 = xb[4];
18363447b6efSHong Zhang       ib = idx + ii[0];
18373447b6efSHong Zhang       n  = ii[1] - ii[0]; ii++;
18382d61bbb3SSatish Balay       for (j=0; j<n; j++) {
18392d61bbb3SSatish Balay         rval      = ib[j]*5;
18402d61bbb3SSatish Balay         z[rval++] +=  v[0]*x1 +  v[1]*x2 +  v[2]*x3 +  v[3]*x4 +  v[4]*x5;
18412d61bbb3SSatish Balay         z[rval++] +=  v[5]*x1 +  v[6]*x2 +  v[7]*x3 +  v[8]*x4 +  v[9]*x5;
18422d61bbb3SSatish Balay         z[rval++] += v[10]*x1 + v[11]*x2 + v[12]*x3 + v[13]*x4 + v[14]*x5;
18432d61bbb3SSatish Balay         z[rval++] += v[15]*x1 + v[16]*x2 + v[17]*x3 + v[18]*x4 + v[19]*x5;
18442d61bbb3SSatish Balay         z[rval++] += v[20]*x1 + v[21]*x2 + v[22]*x3 + v[23]*x4 + v[24]*x5;
18452d61bbb3SSatish Balay         v  += 25;
18462d61bbb3SSatish Balay       }
18473447b6efSHong Zhang       if (!usecprow) xb += 5;
18482d61bbb3SSatish Balay     }
18492d61bbb3SSatish Balay     break;
1850f1af5d2fSBarry Smith   default: {      /* block sizes larger then 5 by 5 are handled by BLAS */
1851690b6cddSBarry Smith       PetscInt     ncols,k;
18523447b6efSHong Zhang       PetscScalar  *work,*workt,*xtmp;
18533f1db9ecSBarry Smith 
18542d61bbb3SSatish Balay       if (!a->mult_work) {
1855d0f46423SBarry Smith         k = PetscMax(A->rmap->n,A->cmap->n);
185687828ca2SBarry Smith         ierr = PetscMalloc((k+1)*sizeof(PetscScalar),&a->mult_work);CHKERRQ(ierr);
18572d61bbb3SSatish Balay       }
18582d61bbb3SSatish Balay       work = a->mult_work;
18593447b6efSHong Zhang       xtmp = x;
18602d61bbb3SSatish Balay       for (i=0; i<mbs; i++) {
18612d61bbb3SSatish Balay         n     = ii[1] - ii[0]; ii++;
18622d61bbb3SSatish Balay         ncols = n*bs;
186387828ca2SBarry Smith         ierr  = PetscMemzero(work,ncols*sizeof(PetscScalar));CHKERRQ(ierr);
18643447b6efSHong Zhang         if (usecprow) {
18657b2bb3b9SHong Zhang           xtmp = x + bs*ridx[i];
18663447b6efSHong Zhang         }
18673447b6efSHong Zhang         Kernel_w_gets_w_plus_trans_Ar_times_v(bs,ncols,xtmp,v,work);
186871044d3cSBarry Smith         /* BLASgemv_("T",&bs,&ncols,&_DOne,v,&bs,xtmp,&_One,&_DOne,work,&_One); */
18692d61bbb3SSatish Balay         v += n*bs2;
18703447b6efSHong Zhang         if (!usecprow) xtmp += bs;
18712d61bbb3SSatish Balay         workt = work;
18722d61bbb3SSatish Balay         for (j=0; j<n; j++) {
18732d61bbb3SSatish Balay           zb = z + bs*(*idx++);
18742d61bbb3SSatish Balay           for (k=0; k<bs; k++) zb[k] += workt[k] ;
18752d61bbb3SSatish Balay           workt += bs;
18762d61bbb3SSatish Balay         }
18772d61bbb3SSatish Balay       }
18782d61bbb3SSatish Balay     }
18792d61bbb3SSatish Balay   }
18803447b6efSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
18813447b6efSHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
1882dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*a->nz*a->bs2);CHKERRQ(ierr);
18832d61bbb3SSatish Balay   PetscFunctionReturn(0);
18842d61bbb3SSatish Balay }
18852d61bbb3SSatish Balay 
18864a2ae208SSatish Balay #undef __FUNCT__
18874a2ae208SSatish Balay #define __FUNCT__ "MatScale_SeqBAIJ"
1888f4df32b1SMatthew Knepley PetscErrorCode MatScale_SeqBAIJ(Mat inA,PetscScalar alpha)
18892d61bbb3SSatish Balay {
18902d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)inA->data;
1891690b6cddSBarry Smith   PetscInt       totalnz = a->bs2*a->nz;
1892f4df32b1SMatthew Knepley   PetscScalar    oalpha = alpha;
1893efee365bSSatish Balay   PetscErrorCode ierr;
18940805154bSBarry Smith   PetscBLASInt   one = 1,tnz = PetscBLASIntCast(totalnz);
18952d61bbb3SSatish Balay 
18962d61bbb3SSatish Balay   PetscFunctionBegin;
1897f4df32b1SMatthew Knepley   BLASscal_(&tnz,&oalpha,a->a,&one);
1898efee365bSSatish Balay   ierr = PetscLogFlops(totalnz);CHKERRQ(ierr);
18992d61bbb3SSatish Balay   PetscFunctionReturn(0);
19002d61bbb3SSatish Balay }
19012d61bbb3SSatish Balay 
19024a2ae208SSatish Balay #undef __FUNCT__
19034a2ae208SSatish Balay #define __FUNCT__ "MatNorm_SeqBAIJ"
1904dfbe8321SBarry Smith PetscErrorCode MatNorm_SeqBAIJ(Mat A,NormType type,PetscReal *norm)
19052d61bbb3SSatish Balay {
19068a62d963SHong Zhang   PetscErrorCode ierr;
19072d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
19083f1db9ecSBarry Smith   MatScalar      *v = a->a;
1909329f5518SBarry Smith   PetscReal      sum = 0.0;
1910d0f46423SBarry Smith   PetscInt       i,j,k,bs=A->rmap->bs,nz=a->nz,bs2=a->bs2,k1;
19112d61bbb3SSatish Balay 
19122d61bbb3SSatish Balay   PetscFunctionBegin;
19132d61bbb3SSatish Balay   if (type == NORM_FROBENIUS) {
19142d61bbb3SSatish Balay     for (i=0; i< bs2*nz; i++) {
1915aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
1916329f5518SBarry Smith       sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
19172d61bbb3SSatish Balay #else
19182d61bbb3SSatish Balay       sum += (*v)*(*v); v++;
19192d61bbb3SSatish Balay #endif
19202d61bbb3SSatish Balay     }
1921*8f1a2a5eSBarry Smith     *norm = PetscSqrtReal(sum);
19228a62d963SHong Zhang   } else if (type == NORM_1) { /* maximum column sum */
19238a62d963SHong Zhang     PetscReal *tmp;
19248a62d963SHong Zhang     PetscInt  *bcol = a->j;
1925d0f46423SBarry Smith     ierr = PetscMalloc((A->cmap->n+1)*sizeof(PetscReal),&tmp);CHKERRQ(ierr);
1926d0f46423SBarry Smith     ierr = PetscMemzero(tmp,A->cmap->n*sizeof(PetscReal));CHKERRQ(ierr);
19278a62d963SHong Zhang     for (i=0; i<nz; i++){
19288a62d963SHong Zhang       for (j=0; j<bs; j++){
19298a62d963SHong Zhang         k1 = bs*(*bcol) + j; /* column index */
19308a62d963SHong Zhang         for (k=0; k<bs; k++){
19318a62d963SHong Zhang           tmp[k1] += PetscAbsScalar(*v); v++;
19328a62d963SHong Zhang         }
19338a62d963SHong Zhang       }
19348a62d963SHong Zhang       bcol++;
19358a62d963SHong Zhang     }
19368a62d963SHong Zhang     *norm = 0.0;
1937d0f46423SBarry Smith     for (j=0; j<A->cmap->n; j++) {
19388a62d963SHong Zhang       if (tmp[j] > *norm) *norm = tmp[j];
19398a62d963SHong Zhang     }
19408a62d963SHong Zhang     ierr = PetscFree(tmp);CHKERRQ(ierr);
1941596552b5SBarry Smith   } else if (type == NORM_INFINITY) { /* maximum row sum */
1942596552b5SBarry Smith     *norm = 0.0;
1943596552b5SBarry Smith     for (k=0; k<bs; k++) {
194474f84c7bSSatish Balay       for (j=0; j<a->mbs; j++) {
1945596552b5SBarry Smith         v = a->a + bs2*a->i[j] + k;
1946596552b5SBarry Smith         sum = 0.0;
1947596552b5SBarry Smith         for (i=0; i<a->i[j+1]-a->i[j]; i++) {
19480e90e235SBarry Smith           for (k1=0; k1<bs; k1++){
1949596552b5SBarry Smith             sum += PetscAbsScalar(*v);
1950596552b5SBarry Smith             v   += bs;
19512d61bbb3SSatish Balay           }
19520e90e235SBarry Smith         }
1953596552b5SBarry Smith         if (sum > *norm) *norm = sum;
1954596552b5SBarry Smith       }
1955596552b5SBarry Smith     }
1956e7e72b3dSBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for this norm yet");
19572d61bbb3SSatish Balay   PetscFunctionReturn(0);
19582d61bbb3SSatish Balay }
19592d61bbb3SSatish Balay 
19602d61bbb3SSatish Balay 
19614a2ae208SSatish Balay #undef __FUNCT__
19624a2ae208SSatish Balay #define __FUNCT__ "MatEqual_SeqBAIJ"
1963ace3abfcSBarry Smith PetscErrorCode MatEqual_SeqBAIJ(Mat A,Mat B,PetscBool * flg)
19642d61bbb3SSatish Balay {
19652d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ *)A->data,*b = (Mat_SeqBAIJ *)B->data;
1966dfbe8321SBarry Smith   PetscErrorCode ierr;
19672d61bbb3SSatish Balay 
19682d61bbb3SSatish Balay   PetscFunctionBegin;
19692d61bbb3SSatish Balay   /* If the  matrix/block dimensions are not equal, or no of nonzeros or shift */
1970d0f46423SBarry Smith   if ((A->rmap->N != B->rmap->N) || (A->cmap->n != B->cmap->n) || (A->rmap->bs != B->rmap->bs)|| (a->nz != b->nz)) {
1971273d9f13SBarry Smith     *flg = PETSC_FALSE;
1972273d9f13SBarry Smith     PetscFunctionReturn(0);
19732d61bbb3SSatish Balay   }
19742d61bbb3SSatish Balay 
19752d61bbb3SSatish Balay   /* if the a->i are the same */
1976690b6cddSBarry Smith   ierr = PetscMemcmp(a->i,b->i,(a->mbs+1)*sizeof(PetscInt),flg);CHKERRQ(ierr);
1977abc0a331SBarry Smith   if (!*flg) {
19780f5bd95cSBarry Smith     PetscFunctionReturn(0);
19792d61bbb3SSatish Balay   }
19802d61bbb3SSatish Balay 
19812d61bbb3SSatish Balay   /* if a->j are the same */
1982690b6cddSBarry Smith   ierr = PetscMemcmp(a->j,b->j,(a->nz)*sizeof(PetscInt),flg);CHKERRQ(ierr);
1983abc0a331SBarry Smith   if (!*flg) {
19840f5bd95cSBarry Smith     PetscFunctionReturn(0);
19852d61bbb3SSatish Balay   }
19862d61bbb3SSatish Balay   /* if a->a are the same */
1987d0f46423SBarry Smith   ierr = PetscMemcmp(a->a,b->a,(a->nz)*(A->rmap->bs)*(B->rmap->bs)*sizeof(PetscScalar),flg);CHKERRQ(ierr);
19882d61bbb3SSatish Balay   PetscFunctionReturn(0);
19892d61bbb3SSatish Balay 
19902d61bbb3SSatish Balay }
19912d61bbb3SSatish Balay 
19924a2ae208SSatish Balay #undef __FUNCT__
19934a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_SeqBAIJ"
1994dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_SeqBAIJ(Mat A,Vec v)
19952d61bbb3SSatish Balay {
19962d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
1997dfbe8321SBarry Smith   PetscErrorCode ierr;
1998690b6cddSBarry Smith   PetscInt       i,j,k,n,row,bs,*ai,*aj,ambs,bs2;
199987828ca2SBarry Smith   PetscScalar    *x,zero = 0.0;
20003f1db9ecSBarry Smith   MatScalar      *aa,*aa_j;
20012d61bbb3SSatish Balay 
20022d61bbb3SSatish Balay   PetscFunctionBegin;
2003e32f2f54SBarry Smith   if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2004d0f46423SBarry Smith   bs   = A->rmap->bs;
20052d61bbb3SSatish Balay   aa   = a->a;
20062d61bbb3SSatish Balay   ai   = a->i;
20072d61bbb3SSatish Balay   aj   = a->j;
20082d61bbb3SSatish Balay   ambs = a->mbs;
20092d61bbb3SSatish Balay   bs2  = a->bs2;
20102d61bbb3SSatish Balay 
20112dcb1b2aSMatthew Knepley   ierr = VecSet(v,zero);CHKERRQ(ierr);
20121ebc52fbSHong Zhang   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
2013e1311b90SBarry Smith   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
2014e32f2f54SBarry Smith   if (n != A->rmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
20152d61bbb3SSatish Balay   for (i=0; i<ambs; i++) {
20162d61bbb3SSatish Balay     for (j=ai[i]; j<ai[i+1]; j++) {
20172d61bbb3SSatish Balay       if (aj[j] == i) {
20182d61bbb3SSatish Balay         row  = i*bs;
20192d61bbb3SSatish Balay         aa_j = aa+j*bs2;
20202d61bbb3SSatish Balay         for (k=0; k<bs2; k+=(bs+1),row++) x[row] = aa_j[k];
20212d61bbb3SSatish Balay         break;
20222d61bbb3SSatish Balay       }
20232d61bbb3SSatish Balay     }
20242d61bbb3SSatish Balay   }
20251ebc52fbSHong Zhang   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
20262d61bbb3SSatish Balay   PetscFunctionReturn(0);
20272d61bbb3SSatish Balay }
20282d61bbb3SSatish Balay 
20294a2ae208SSatish Balay #undef __FUNCT__
20304a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_SeqBAIJ"
2031dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_SeqBAIJ(Mat A,Vec ll,Vec rr)
20322d61bbb3SSatish Balay {
20332d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
203453ef36baSBarry Smith   const PetscScalar *l,*r,*li,*ri;
203553ef36baSBarry Smith   PetscScalar       x;
20363f1db9ecSBarry Smith   MatScalar         *aa, *v;
2037dfbe8321SBarry Smith   PetscErrorCode    ierr;
203853ef36baSBarry Smith   PetscInt          i,j,k,lm,rn,M,m,n,mbs,tmp,bs,bs2,iai;
203953ef36baSBarry Smith   const PetscInt    *ai,*aj;
20402d61bbb3SSatish Balay 
20412d61bbb3SSatish Balay   PetscFunctionBegin;
20422d61bbb3SSatish Balay   ai  = a->i;
20432d61bbb3SSatish Balay   aj  = a->j;
20442d61bbb3SSatish Balay   aa  = a->a;
2045d0f46423SBarry Smith   m   = A->rmap->n;
2046d0f46423SBarry Smith   n   = A->cmap->n;
2047d0f46423SBarry Smith   bs  = A->rmap->bs;
20482d61bbb3SSatish Balay   mbs = a->mbs;
20492d61bbb3SSatish Balay   bs2 = a->bs2;
20502d61bbb3SSatish Balay   if (ll) {
205153ef36baSBarry Smith     ierr = VecGetArrayRead(ll,&l);CHKERRQ(ierr);
2052e1311b90SBarry Smith     ierr = VecGetLocalSize(ll,&lm);CHKERRQ(ierr);
2053e32f2f54SBarry Smith     if (lm != m) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Left scaling vector wrong length");
20542d61bbb3SSatish Balay     for (i=0; i<mbs; i++) { /* for each block row */
20552d61bbb3SSatish Balay       M  = ai[i+1] - ai[i];
20562d61bbb3SSatish Balay       li = l + i*bs;
20572d61bbb3SSatish Balay       v  = aa + bs2*ai[i];
20582d61bbb3SSatish Balay       for (j=0; j<M; j++) { /* for each block */
20592d61bbb3SSatish Balay         for (k=0; k<bs2; k++) {
20602d61bbb3SSatish Balay           (*v++) *= li[k%bs];
20612d61bbb3SSatish Balay         }
20622d61bbb3SSatish Balay       }
20632d61bbb3SSatish Balay     }
206453ef36baSBarry Smith     ierr = VecRestoreArrayRead(ll,&l);CHKERRQ(ierr);
2065efee365bSSatish Balay     ierr = PetscLogFlops(a->nz);CHKERRQ(ierr);
20662d61bbb3SSatish Balay   }
20672d61bbb3SSatish Balay 
20682d61bbb3SSatish Balay   if (rr) {
206953ef36baSBarry Smith     ierr = VecGetArrayRead(rr,&r);CHKERRQ(ierr);
2070e1311b90SBarry Smith     ierr = VecGetLocalSize(rr,&rn);CHKERRQ(ierr);
2071e32f2f54SBarry Smith     if (rn != n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Right scaling vector wrong length");
20722d61bbb3SSatish Balay     for (i=0; i<mbs; i++) { /* for each block row */
207353ef36baSBarry Smith       iai = ai[i];
207453ef36baSBarry Smith       M   = ai[i+1] - iai;
207553ef36baSBarry Smith       v   = aa + bs2*iai;
20762d61bbb3SSatish Balay       for (j=0; j<M; j++) { /* for each block */
207753ef36baSBarry Smith         ri = r + bs*aj[iai+j];
20782d61bbb3SSatish Balay         for (k=0; k<bs; k++) {
20792d61bbb3SSatish Balay           x = ri[k];
208053ef36baSBarry Smith           for (tmp=0; tmp<bs; tmp++) v[tmp] *= x;
208153ef36baSBarry Smith           v += bs;
20822d61bbb3SSatish Balay         }
20832d61bbb3SSatish Balay       }
20842d61bbb3SSatish Balay     }
208553ef36baSBarry Smith     ierr = VecRestoreArrayRead(rr,&r);CHKERRQ(ierr);
2086efee365bSSatish Balay     ierr = PetscLogFlops(a->nz);CHKERRQ(ierr);
20872d61bbb3SSatish Balay   }
20882d61bbb3SSatish Balay   PetscFunctionReturn(0);
20892d61bbb3SSatish Balay }
20902d61bbb3SSatish Balay 
20912d61bbb3SSatish Balay 
20924a2ae208SSatish Balay #undef __FUNCT__
20934a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_SeqBAIJ"
2094dfbe8321SBarry Smith PetscErrorCode MatGetInfo_SeqBAIJ(Mat A,MatInfoType flag,MatInfo *info)
20952d61bbb3SSatish Balay {
20962d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
20972d61bbb3SSatish Balay 
20982d61bbb3SSatish Balay   PetscFunctionBegin;
20992d61bbb3SSatish Balay   info->block_size     = a->bs2;
2100ceed8ce5SJed Brown   info->nz_allocated   = a->bs2*a->maxnz;
21012d61bbb3SSatish Balay   info->nz_used        = a->bs2*a->nz;
21022d61bbb3SSatish Balay   info->nz_unneeded    = (double)(info->nz_allocated - info->nz_used);
21032d61bbb3SSatish Balay   info->assemblies   = A->num_ass;
21048e58a170SBarry Smith   info->mallocs      = A->info.mallocs;
21057adad957SLisandro Dalcin   info->memory       = ((PetscObject)A)->mem;
2106d5f3da31SBarry Smith   if (A->factortype) {
21072d61bbb3SSatish Balay     info->fill_ratio_given  = A->info.fill_ratio_given;
21082d61bbb3SSatish Balay     info->fill_ratio_needed = A->info.fill_ratio_needed;
21092d61bbb3SSatish Balay     info->factor_mallocs    = A->info.factor_mallocs;
21102d61bbb3SSatish Balay   } else {
21112d61bbb3SSatish Balay     info->fill_ratio_given  = 0;
21122d61bbb3SSatish Balay     info->fill_ratio_needed = 0;
21132d61bbb3SSatish Balay     info->factor_mallocs    = 0;
21142d61bbb3SSatish Balay   }
21152d61bbb3SSatish Balay   PetscFunctionReturn(0);
21162d61bbb3SSatish Balay }
21172d61bbb3SSatish Balay 
21182d61bbb3SSatish Balay 
21194a2ae208SSatish Balay #undef __FUNCT__
21204a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_SeqBAIJ"
2121dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_SeqBAIJ(Mat A)
21222d61bbb3SSatish Balay {
21232d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
2124dfbe8321SBarry Smith   PetscErrorCode ierr;
21252d61bbb3SSatish Balay 
21262d61bbb3SSatish Balay   PetscFunctionBegin;
2127549d3d68SSatish Balay   ierr = PetscMemzero(a->a,a->bs2*a->i[a->mbs]*sizeof(MatScalar));CHKERRQ(ierr);
21282d61bbb3SSatish Balay   PetscFunctionReturn(0);
21292d61bbb3SSatish Balay }
2130