xref: /petsc/src/mat/impls/baij/seq/baij2.c (revision 16b643558ddb2850265bc9ef423683e46e446b85)
1cac129eeSSatish Balay 
2c6db04a5SJed Brown #include <../src/mat/impls/baij/seq/baij.h>
3af0996ceSBarry Smith #include <petsc/private/kernels/blockinvert.h>
4c6db04a5SJed Brown #include <petscbt.h>
5c6db04a5SJed Brown #include <petscblaslapack.h>
6cac129eeSSatish Balay 
7690b6cddSBarry Smith PetscErrorCode MatIncreaseOverlap_SeqBAIJ(Mat A,PetscInt is_max,IS is[],PetscInt ov)
8a3192f15SSatish Balay {
9a3192f15SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
106849ba73SBarry Smith   PetscErrorCode ierr;
115d0c19d7SBarry Smith   PetscInt       row,i,j,k,l,m,n,*nidx,isz,val,ival;
125d0c19d7SBarry Smith   const PetscInt *idx;
13690b6cddSBarry Smith   PetscInt       start,end,*ai,*aj,bs,*nidx2;
14f1af5d2fSBarry Smith   PetscBT        table;
15a3192f15SSatish Balay 
163a40ed3dSBarry Smith   PetscFunctionBegin;
17a3192f15SSatish Balay   m  = a->mbs;
18a3192f15SSatish Balay   ai = a->i;
19a3192f15SSatish Balay   aj = a->j;
20d0f46423SBarry Smith   bs = A->rmap->bs;
21a3192f15SSatish Balay 
22e32f2f54SBarry Smith   if (ov < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap specified");
23a3192f15SSatish Balay 
2453b8de81SBarry Smith   ierr = PetscBTCreate(m,&table);CHKERRQ(ierr);
25854ce69bSBarry Smith   ierr = PetscMalloc1(m+1,&nidx);CHKERRQ(ierr);
26854ce69bSBarry Smith   ierr = PetscMalloc1(A->rmap->N+1,&nidx2);CHKERRQ(ierr);
27a3192f15SSatish Balay 
28a3192f15SSatish Balay   for (i=0; i<is_max; i++) {
29a3192f15SSatish Balay     /* Initialise the two local arrays */
30a3192f15SSatish Balay     isz  = 0;
316831982aSBarry Smith     ierr = PetscBTMemzero(m,table);CHKERRQ(ierr);
32a3192f15SSatish Balay 
33a3192f15SSatish Balay     /* Extract the indices, assume there can be duplicate entries */
34a3192f15SSatish Balay     ierr = ISGetIndices(is[i],&idx);CHKERRQ(ierr);
35b9b97703SBarry Smith     ierr = ISGetLocalSize(is[i],&n);CHKERRQ(ierr);
36a3192f15SSatish Balay 
37a3192f15SSatish Balay     /* Enter these into the temp arrays i.e mark table[row], enter row into new index */
38a3192f15SSatish Balay     for (j=0; j<n; ++j) {
39218c64b6SSatish Balay       ival = idx[j]/bs; /* convert the indices into block indices */
40e32f2f54SBarry Smith       if (ival>=m) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"index greater than mat-dim");
4126fbe8dcSKarl Rupp       if (!PetscBTLookupSet(table,ival)) nidx[isz++] = ival;
42a3192f15SSatish Balay     }
43a3192f15SSatish Balay     ierr = ISRestoreIndices(is[i],&idx);CHKERRQ(ierr);
446bf464f9SBarry Smith     ierr = ISDestroy(&is[i]);CHKERRQ(ierr);
45a3192f15SSatish Balay 
46a3192f15SSatish Balay     k = 0;
47a3192f15SSatish Balay     for (j=0; j<ov; j++) { /* for each overlap*/
48a3192f15SSatish Balay       n = isz;
49a3192f15SSatish Balay       for (; k<n; k++) {  /* do only those rows in nidx[k], which are not done yet */
50a3192f15SSatish Balay         row   = nidx[k];
51a3192f15SSatish Balay         start = ai[row];
52a3192f15SSatish Balay         end   = ai[row+1];
53a3192f15SSatish Balay         for (l = start; l<end; l++) {
54a3192f15SSatish Balay           val = aj[l];
5526fbe8dcSKarl Rupp           if (!PetscBTLookupSet(table,val)) nidx[isz++] = val;
56a3192f15SSatish Balay         }
57a3192f15SSatish Balay       }
58a3192f15SSatish Balay     }
59218c64b6SSatish Balay     /* expand the Index Set */
60218c64b6SSatish Balay     for (j=0; j<isz; j++) {
6126fbe8dcSKarl Rupp       for (k=0; k<bs; k++) nidx2[j*bs+k] = nidx[j]*bs+k;
62218c64b6SSatish Balay     }
6370b3c8c7SBarry Smith     ierr = ISCreateGeneral(PETSC_COMM_SELF,isz*bs,nidx2,PETSC_COPY_VALUES,is+i);CHKERRQ(ierr);
64a3192f15SSatish Balay   }
6594bacf5dSBarry Smith   ierr = PetscBTDestroy(&table);CHKERRQ(ierr);
66606d414cSSatish Balay   ierr = PetscFree(nidx);CHKERRQ(ierr);
67606d414cSSatish Balay   ierr = PetscFree(nidx2);CHKERRQ(ierr);
683a40ed3dSBarry Smith   PetscFunctionReturn(0);
69a3192f15SSatish Balay }
701c351548SSatish Balay 
714aa3045dSJed Brown PetscErrorCode MatGetSubMatrix_SeqBAIJ_Private(Mat A,IS isrow,IS iscol,MatReuse scall,Mat *B)
72736121d4SSatish Balay {
73736121d4SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*c;
746849ba73SBarry Smith   PetscErrorCode ierr;
75690b6cddSBarry Smith   PetscInt       *smap,i,k,kstart,kend,oldcols = a->nbs,*lens;
76690b6cddSBarry Smith   PetscInt       row,mat_i,*mat_j,tcol,*mat_ilen;
775d0c19d7SBarry Smith   const PetscInt *irow,*icol;
785d0c19d7SBarry Smith   PetscInt       nrows,ncols,*ssmap,bs=A->rmap->bs,bs2=a->bs2;
79690b6cddSBarry Smith   PetscInt       *aj = a->j,*ai = a->i;
803f1db9ecSBarry Smith   MatScalar      *mat_a;
81736121d4SSatish Balay   Mat            C;
826041f1b1SToby Isaac   PetscBool      flag;
83736121d4SSatish Balay 
843a40ed3dSBarry Smith   PetscFunctionBegin;
85736121d4SSatish Balay   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
86218c64b6SSatish Balay   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
87b9b97703SBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
88b9b97703SBarry Smith   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
89736121d4SSatish Balay 
90854ce69bSBarry Smith   ierr  = PetscCalloc1(1+oldcols,&smap);CHKERRQ(ierr);
91736121d4SSatish Balay   ssmap = smap;
92854ce69bSBarry Smith   ierr  = PetscMalloc1(1+nrows,&lens);CHKERRQ(ierr);
93736121d4SSatish Balay   for (i=0; i<ncols; i++) smap[icol[i]] = i+1;
94736121d4SSatish Balay   /* determine lens of each row */
95736121d4SSatish Balay   for (i=0; i<nrows; i++) {
96736121d4SSatish Balay     kstart  = ai[irow[i]];
97736121d4SSatish Balay     kend    = kstart + a->ilen[irow[i]];
98736121d4SSatish Balay     lens[i] = 0;
99736121d4SSatish Balay     for (k=kstart; k<kend; k++) {
10026fbe8dcSKarl Rupp       if (ssmap[aj[k]]) lens[i]++;
101736121d4SSatish Balay     }
102736121d4SSatish Balay   }
103736121d4SSatish Balay   /* Create and fill new matrix */
104736121d4SSatish Balay   if (scall == MAT_REUSE_MATRIX) {
105736121d4SSatish Balay     c = (Mat_SeqBAIJ*)((*B)->data);
106736121d4SSatish Balay 
107e32f2f54SBarry Smith     if (c->mbs!=nrows || c->nbs!=ncols || (*B)->rmap->bs!=bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Submatrix wrong size");
108690b6cddSBarry Smith     ierr = PetscMemcmp(c->ilen,lens,c->mbs *sizeof(PetscInt),&flag);CHKERRQ(ierr);
109e7e72b3dSBarry Smith     if (!flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong no of nonzeros");
110690b6cddSBarry Smith     ierr = PetscMemzero(c->ilen,c->mbs*sizeof(PetscInt));CHKERRQ(ierr);
111736121d4SSatish Balay     C    = *B;
1123a40ed3dSBarry Smith   } else {
113ce94432eSBarry Smith     ierr = MatCreate(PetscObjectComm((PetscObject)A),&C);CHKERRQ(ierr);
114f69a0ea3SMatthew Knepley     ierr = MatSetSizes(C,nrows*bs,ncols*bs,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1157adad957SLisandro Dalcin     ierr = MatSetType(C,((PetscObject)A)->type_name);CHKERRQ(ierr);
116367daffbSBarry Smith     ierr = MatSeqBAIJSetPreallocation(C,bs,0,lens);CHKERRQ(ierr);
117736121d4SSatish Balay   }
118736121d4SSatish Balay   c = (Mat_SeqBAIJ*)(C->data);
119736121d4SSatish Balay   for (i=0; i<nrows; i++) {
120736121d4SSatish Balay     row      = irow[i];
121736121d4SSatish Balay     kstart   = ai[row];
122736121d4SSatish Balay     kend     = kstart + a->ilen[row];
123736121d4SSatish Balay     mat_i    = c->i[i];
124736121d4SSatish Balay     mat_j    = c->j + mat_i;
125218c64b6SSatish Balay     mat_a    = c->a + mat_i*bs2;
126736121d4SSatish Balay     mat_ilen = c->ilen + i;
127736121d4SSatish Balay     for (k=kstart; k<kend; k++) {
128736121d4SSatish Balay       if ((tcol=ssmap[a->j[k]])) {
129736121d4SSatish Balay         *mat_j++ = tcol - 1;
130549d3d68SSatish Balay         ierr     = PetscMemcpy(mat_a,a->a+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
131549d3d68SSatish Balay         mat_a   += bs2;
132736121d4SSatish Balay         (*mat_ilen)++;
133736121d4SSatish Balay       }
134736121d4SSatish Balay     }
135736121d4SSatish Balay   }
136cdc6f3adSToby Isaac   /* sort */
137cdc6f3adSToby Isaac   {
138cdc6f3adSToby Isaac     MatScalar *work;
139cdc6f3adSToby Isaac     ierr = PetscMalloc1(bs2,&work);CHKERRQ(ierr);
140cdc6f3adSToby Isaac     for (i=0; i<nrows; i++) {
141cdc6f3adSToby Isaac       PetscInt ilen;
142cdc6f3adSToby Isaac       mat_i = c->i[i];
143cdc6f3adSToby Isaac       mat_j = c->j + mat_i;
144cdc6f3adSToby Isaac       mat_a = c->a + mat_i*bs2;
145cdc6f3adSToby Isaac       ilen  = c->ilen[i];
146cdc6f3adSToby Isaac       ierr  = PetscSortIntWithDataArray(ilen,mat_j,mat_a,bs2*sizeof(MatScalar),work);CHKERRQ(ierr);
147cdc6f3adSToby Isaac     }
148cdc6f3adSToby Isaac     ierr = PetscFree(work);CHKERRQ(ierr);
149cdc6f3adSToby Isaac   }
150218c64b6SSatish Balay 
151736121d4SSatish Balay   /* Free work space */
152736121d4SSatish Balay   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
153606d414cSSatish Balay   ierr = PetscFree(smap);CHKERRQ(ierr);
154606d414cSSatish Balay   ierr = PetscFree(lens);CHKERRQ(ierr);
1556d4a8577SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1566d4a8577SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
157736121d4SSatish Balay 
158736121d4SSatish Balay   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
159736121d4SSatish Balay   *B   = C;
1603a40ed3dSBarry Smith   PetscFunctionReturn(0);
161736121d4SSatish Balay }
162736121d4SSatish Balay 
1634aa3045dSJed Brown PetscErrorCode MatGetSubMatrix_SeqBAIJ(Mat A,IS isrow,IS iscol,MatReuse scall,Mat *B)
164218c64b6SSatish Balay {
165218c64b6SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
166218c64b6SSatish Balay   IS             is1,is2;
1676849ba73SBarry Smith   PetscErrorCode ierr;
168afebec48SHong Zhang   PetscInt       *vary,*iary,nrows,ncols,i,bs=A->rmap->bs,count,maxmnbs,j;
1695d0c19d7SBarry Smith   const PetscInt *irow,*icol;
170218c64b6SSatish Balay 
1713a40ed3dSBarry Smith   PetscFunctionBegin;
172218c64b6SSatish Balay   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
173218c64b6SSatish Balay   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
174b9b97703SBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
175b9b97703SBarry Smith   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
176218c64b6SSatish Balay 
177218c64b6SSatish Balay   /* Verify if the indices corespond to each element in a block
178218c64b6SSatish Balay    and form the IS with compressed IS */
179f8ecb639SStefano Zampini   maxmnbs = PetscMax(a->mbs,a->nbs);
180f8ecb639SStefano Zampini   ierr = PetscMalloc2(maxmnbs,&vary,maxmnbs,&iary);CHKERRQ(ierr);
181fca92195SBarry Smith   ierr = PetscMemzero(vary,a->mbs*sizeof(PetscInt));CHKERRQ(ierr);
182218c64b6SSatish Balay   for (i=0; i<nrows; i++) vary[irow[i]/bs]++;
183218c64b6SSatish Balay   for (i=0; i<a->mbs; i++) {
184e32f2f54SBarry Smith     if (vary[i]!=0 && vary[i]!=bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Index set does not match blocks");
1856041f1b1SToby Isaac   }
1866041f1b1SToby Isaac   count = 0;
1876041f1b1SToby Isaac   for (i=0; i<nrows; i++) {
188afebec48SHong Zhang     j = irow[i] / bs;
1896041f1b1SToby Isaac     if ((vary[j]--)==bs) iary[count++] = j;
190218c64b6SSatish Balay   }
19170b3c8c7SBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,count,iary,PETSC_COPY_VALUES,&is1);CHKERRQ(ierr);
192218c64b6SSatish Balay 
193f8ecb639SStefano Zampini   ierr = PetscMemzero(vary,(a->nbs)*sizeof(PetscInt));CHKERRQ(ierr);
194218c64b6SSatish Balay   for (i=0; i<ncols; i++) vary[icol[i]/bs]++;
195f8ecb639SStefano Zampini   for (i=0; i<a->nbs; i++) {
196e32f2f54SBarry Smith     if (vary[i]!=0 && vary[i]!=bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Internal error in PETSc");
197218c64b6SSatish Balay   }
1986041f1b1SToby Isaac   count = 0;
1996041f1b1SToby Isaac   for (i=0; i<ncols; i++) {
200afebec48SHong Zhang     j = icol[i] / bs;
2016041f1b1SToby Isaac     if ((vary[j]--)==bs) iary[count++] = j;
2026041f1b1SToby Isaac   }
20370b3c8c7SBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,count,iary,PETSC_COPY_VALUES,&is2);CHKERRQ(ierr);
204218c64b6SSatish Balay   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
205218c64b6SSatish Balay   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
206fca92195SBarry Smith   ierr = PetscFree2(vary,iary);CHKERRQ(ierr);
207218c64b6SSatish Balay 
2084aa3045dSJed Brown   ierr = MatGetSubMatrix_SeqBAIJ_Private(A,is1,is2,scall,B);CHKERRQ(ierr);
2096bf464f9SBarry Smith   ierr = ISDestroy(&is1);CHKERRQ(ierr);
2106bf464f9SBarry Smith   ierr = ISDestroy(&is2);CHKERRQ(ierr);
2113a40ed3dSBarry Smith   PetscFunctionReturn(0);
212218c64b6SSatish Balay }
213218c64b6SSatish Balay 
214*16b64355SHong Zhang PetscErrorCode MatDestroy_SeqBAIJ_Submatrices(Mat C)
215*16b64355SHong Zhang {
216*16b64355SHong Zhang   PetscErrorCode ierr;
217*16b64355SHong Zhang   Mat_SeqBAIJ    *c = (Mat_SeqBAIJ*)C->data;
218*16b64355SHong Zhang   Mat_SubMat     *submatj = c->submatis1;
219*16b64355SHong Zhang 
220*16b64355SHong Zhang   PetscFunctionBegin;
221*16b64355SHong Zhang   ierr = submatj->destroy(C);CHKERRQ(ierr);
222*16b64355SHong Zhang   ierr = MatDestroySubMatrices_Private(submatj);
223*16b64355SHong Zhang   PetscFunctionReturn(0);
224*16b64355SHong Zhang }
225*16b64355SHong Zhang 
226690b6cddSBarry Smith PetscErrorCode MatGetSubMatrices_SeqBAIJ(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
227736121d4SSatish Balay {
2286849ba73SBarry Smith   PetscErrorCode ierr;
229690b6cddSBarry Smith   PetscInt       i;
230736121d4SSatish Balay 
2313a40ed3dSBarry Smith   PetscFunctionBegin;
232736121d4SSatish Balay   if (scall == MAT_INITIAL_MATRIX) {
233df750dc8SHong Zhang     ierr = PetscCalloc1(n+1,B);CHKERRQ(ierr);
234736121d4SSatish Balay   }
235736121d4SSatish Balay 
236736121d4SSatish Balay   for (i=0; i<n; i++) {
2374aa3045dSJed Brown     ierr = MatGetSubMatrix_SeqBAIJ(A,irow[i],icol[i],scall,&(*B)[i]);CHKERRQ(ierr);
238736121d4SSatish Balay   }
2393a40ed3dSBarry Smith   PetscFunctionReturn(0);
240736121d4SSatish Balay }
241218c64b6SSatish Balay 
242218c64b6SSatish Balay 
2432d61bbb3SSatish Balay /* -------------------------------------------------------*/
2442d61bbb3SSatish Balay /* Should check that shapes of vectors and matrices match */
2452d61bbb3SSatish Balay /* -------------------------------------------------------*/
2462d61bbb3SSatish Balay 
247dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_1(Mat A,Vec xx,Vec zz)
2482d61bbb3SSatish Balay {
2492d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
250d9fead3dSBarry Smith   PetscScalar       *z,sum;
251d9fead3dSBarry Smith   const PetscScalar *x;
252d9fead3dSBarry Smith   const MatScalar   *v;
2536849ba73SBarry Smith   PetscErrorCode    ierr;
2547c565772SBarry Smith   PetscInt          mbs,i,n;
2550298fd71SBarry Smith   const PetscInt    *idx,*ii,*ridx=NULL;
256ace3abfcSBarry Smith   PetscBool         usecprow=a->compressedrow.use;
2572d61bbb3SSatish Balay 
2582d61bbb3SSatish Balay   PetscFunctionBegin;
2593649974fSBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
2601ebc52fbSHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
2612d61bbb3SSatish Balay 
26226e093fcSHong Zhang   if (usecprow) {
26326e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
26426e093fcSHong Zhang     ii   = a->compressedrow.i;
2657b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
266ef16e671SStefano Zampini     ierr = PetscMemzero(z,a->mbs*sizeof(PetscScalar));CHKERRQ(ierr);
26726e093fcSHong Zhang   } else {
26826e093fcSHong Zhang     mbs = a->mbs;
2692d61bbb3SSatish Balay     ii  = a->i;
27026e093fcSHong Zhang   }
2712d61bbb3SSatish Balay 
2722d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
273ee54c7eeSHong Zhang     n   = ii[1] - ii[0];
274ee54c7eeSHong Zhang     v   = a->a + ii[0];
275ee54c7eeSHong Zhang     idx = a->j + ii[0];
276ee54c7eeSHong Zhang     ii++;
277444d8c10SJed Brown     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);   /* Indices for the next row (assumes same size as this one) */
278444d8c10SJed Brown     PetscPrefetchBlock(v+1*n,1*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
2792d61bbb3SSatish Balay     sum = 0.0;
2802162cab8SBarry Smith     PetscSparseDensePlusDot(sum,x,v,idx,n);
28126e093fcSHong Zhang     if (usecprow) {
2827b2bb3b9SHong Zhang       z[ridx[i]] = sum;
28326e093fcSHong Zhang     } else {
2842d61bbb3SSatish Balay       z[i]       = sum;
2852d61bbb3SSatish Balay     }
28626e093fcSHong Zhang   }
2873649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
2881ebc52fbSHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
2897c565772SBarry Smith   ierr = PetscLogFlops(2.0*a->nz - a->nonzerorowcnt);CHKERRQ(ierr);
2902d61bbb3SSatish Balay   PetscFunctionReturn(0);
2912d61bbb3SSatish Balay }
2922d61bbb3SSatish Balay 
293dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_2(Mat A,Vec xx,Vec zz)
2942d61bbb3SSatish Balay {
2952d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
296d9fead3dSBarry Smith   PetscScalar       *z = 0,sum1,sum2,*zarray;
297d9fead3dSBarry Smith   const PetscScalar *x,*xb;
29887828ca2SBarry Smith   PetscScalar       x1,x2;
299d9fead3dSBarry Smith   const MatScalar   *v;
300dfbe8321SBarry Smith   PetscErrorCode    ierr;
3017c565772SBarry Smith   PetscInt          mbs,i,*idx,*ii,j,n,*ridx=NULL;
302ace3abfcSBarry Smith   PetscBool         usecprow=a->compressedrow.use;
3032d61bbb3SSatish Balay 
3042d61bbb3SSatish Balay   PetscFunctionBegin;
3053649974fSBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
30626e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
3072d61bbb3SSatish Balay 
3082d61bbb3SSatish Balay   idx = a->j;
3092d61bbb3SSatish Balay   v   = a->a;
31026e093fcSHong Zhang   if (usecprow) {
31126e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
31226e093fcSHong Zhang     ii   = a->compressedrow.i;
3137b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
314ef16e671SStefano Zampini     ierr = PetscMemzero(zarray,2*a->mbs*sizeof(PetscScalar));CHKERRQ(ierr);
31526e093fcSHong Zhang   } else {
31626e093fcSHong Zhang     mbs = a->mbs;
3172d61bbb3SSatish Balay     ii  = a->i;
31826e093fcSHong Zhang     z   = zarray;
31926e093fcSHong Zhang   }
3202d61bbb3SSatish Balay 
3212d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
3222d61bbb3SSatish Balay     n           = ii[1] - ii[0]; ii++;
3232d61bbb3SSatish Balay     sum1        = 0.0; sum2 = 0.0;
324444d8c10SJed Brown     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);   /* Indices for the next row (assumes same size as this one) */
325444d8c10SJed Brown     PetscPrefetchBlock(v+4*n,4*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
3262d61bbb3SSatish Balay     for (j=0; j<n; j++) {
3272d61bbb3SSatish Balay       xb    = x + 2*(*idx++); x1 = xb[0]; x2 = xb[1];
3282d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[2]*x2;
3292d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[3]*x2;
3302d61bbb3SSatish Balay       v    += 4;
3312d61bbb3SSatish Balay     }
3327b2bb3b9SHong Zhang     if (usecprow) z = zarray + 2*ridx[i];
3332d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2;
33426e093fcSHong Zhang     if (!usecprow) z += 2;
3352d61bbb3SSatish Balay   }
3363649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
33726e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
3387c565772SBarry Smith   ierr = PetscLogFlops(8.0*a->nz - 2.0*a->nonzerorowcnt);CHKERRQ(ierr);
3392d61bbb3SSatish Balay   PetscFunctionReturn(0);
3402d61bbb3SSatish Balay }
3412d61bbb3SSatish Balay 
342dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_3(Mat A,Vec xx,Vec zz)
3432d61bbb3SSatish Balay {
3442d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
345d9fead3dSBarry Smith   PetscScalar       *z = 0,sum1,sum2,sum3,x1,x2,x3,*zarray;
346d9fead3dSBarry Smith   const PetscScalar *x,*xb;
347d9fead3dSBarry Smith   const MatScalar   *v;
348dfbe8321SBarry Smith   PetscErrorCode    ierr;
3497c565772SBarry Smith   PetscInt          mbs,i,*idx,*ii,j,n,*ridx=NULL;
350ace3abfcSBarry Smith   PetscBool         usecprow=a->compressedrow.use;
35126e093fcSHong Zhang 
352b6410449SSatish Balay #if defined(PETSC_HAVE_PRAGMA_DISJOINT)
353fee21e36SBarry Smith #pragma disjoint(*v,*z,*xb)
354fee21e36SBarry Smith #endif
355fee21e36SBarry Smith 
3562d61bbb3SSatish Balay   PetscFunctionBegin;
3573649974fSBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
35826e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
3592d61bbb3SSatish Balay 
3602d61bbb3SSatish Balay   idx = a->j;
3612d61bbb3SSatish Balay   v   = a->a;
36226e093fcSHong Zhang   if (usecprow) {
36326e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
36426e093fcSHong Zhang     ii   = a->compressedrow.i;
3657b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
366ef16e671SStefano Zampini     ierr = PetscMemzero(zarray,3*a->mbs*sizeof(PetscScalar));CHKERRQ(ierr);
36726e093fcSHong Zhang   } else {
36826e093fcSHong Zhang     mbs = a->mbs;
3692d61bbb3SSatish Balay     ii  = a->i;
37026e093fcSHong Zhang     z   = zarray;
37126e093fcSHong Zhang   }
3722d61bbb3SSatish Balay 
3732d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
3742d61bbb3SSatish Balay     n           = ii[1] - ii[0]; ii++;
3752d61bbb3SSatish Balay     sum1        = 0.0; sum2 = 0.0; sum3 = 0.0;
376444d8c10SJed Brown     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);   /* Indices for the next row (assumes same size as this one) */
377444d8c10SJed Brown     PetscPrefetchBlock(v+9*n,9*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
3782d61bbb3SSatish Balay     for (j=0; j<n; j++) {
37926fbe8dcSKarl Rupp       xb = x + 3*(*idx++);
38026fbe8dcSKarl Rupp       x1 = xb[0];
38126fbe8dcSKarl Rupp       x2 = xb[1];
38226fbe8dcSKarl Rupp       x3 = xb[2];
38326fbe8dcSKarl Rupp 
3842d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[3]*x2 + v[6]*x3;
3852d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[4]*x2 + v[7]*x3;
3862d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[5]*x2 + v[8]*x3;
3872d61bbb3SSatish Balay       v    += 9;
3882d61bbb3SSatish Balay     }
3897b2bb3b9SHong Zhang     if (usecprow) z = zarray + 3*ridx[i];
3902d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3;
39126e093fcSHong Zhang     if (!usecprow) z += 3;
3922d61bbb3SSatish Balay   }
3933649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
39426e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
3957c565772SBarry Smith   ierr = PetscLogFlops(18.0*a->nz - 3.0*a->nonzerorowcnt);CHKERRQ(ierr);
3962d61bbb3SSatish Balay   PetscFunctionReturn(0);
3972d61bbb3SSatish Balay }
3982d61bbb3SSatish Balay 
399dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_4(Mat A,Vec xx,Vec zz)
4002d61bbb3SSatish Balay {
4012d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
402d9fead3dSBarry Smith   PetscScalar       *z = 0,sum1,sum2,sum3,sum4,x1,x2,x3,x4,*zarray;
403d9fead3dSBarry Smith   const PetscScalar *x,*xb;
404d9fead3dSBarry Smith   const MatScalar   *v;
405dfbe8321SBarry Smith   PetscErrorCode    ierr;
4067c565772SBarry Smith   PetscInt          mbs,i,*idx,*ii,j,n,*ridx=NULL;
407ace3abfcSBarry Smith   PetscBool         usecprow=a->compressedrow.use;
4082d61bbb3SSatish Balay 
4092d61bbb3SSatish Balay   PetscFunctionBegin;
4103649974fSBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
41126e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
4122d61bbb3SSatish Balay 
4132d61bbb3SSatish Balay   idx = a->j;
4142d61bbb3SSatish Balay   v   = a->a;
41526e093fcSHong Zhang   if (usecprow) {
41626e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
41726e093fcSHong Zhang     ii   = a->compressedrow.i;
4187b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
419ef16e671SStefano Zampini     ierr = PetscMemzero(zarray,4*a->mbs*sizeof(PetscScalar));CHKERRQ(ierr);
42026e093fcSHong Zhang   } else {
42126e093fcSHong Zhang     mbs = a->mbs;
4222d61bbb3SSatish Balay     ii  = a->i;
42326e093fcSHong Zhang     z   = zarray;
42426e093fcSHong Zhang   }
4252d61bbb3SSatish Balay 
4262d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
42726fbe8dcSKarl Rupp     n = ii[1] - ii[0];
42826fbe8dcSKarl Rupp     ii++;
42926fbe8dcSKarl Rupp     sum1 = 0.0;
43026fbe8dcSKarl Rupp     sum2 = 0.0;
43126fbe8dcSKarl Rupp     sum3 = 0.0;
43226fbe8dcSKarl Rupp     sum4 = 0.0;
43326fbe8dcSKarl Rupp 
434444d8c10SJed Brown     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);     /* Indices for the next row (assumes same size as this one) */
435444d8c10SJed Brown     PetscPrefetchBlock(v+16*n,16*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
4362d61bbb3SSatish Balay     for (j=0; j<n; j++) {
4372d61bbb3SSatish Balay       xb    = x + 4*(*idx++);
4382d61bbb3SSatish Balay       x1    = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
4392d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[4]*x2 + v[8]*x3   + v[12]*x4;
4402d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[5]*x2 + v[9]*x3   + v[13]*x4;
4412d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[6]*x2 + v[10]*x3  + v[14]*x4;
4422d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[7]*x2 + v[11]*x3  + v[15]*x4;
4432d61bbb3SSatish Balay       v    += 16;
4442d61bbb3SSatish Balay     }
4457b2bb3b9SHong Zhang     if (usecprow) z = zarray + 4*ridx[i];
4462d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4;
44726e093fcSHong Zhang     if (!usecprow) z += 4;
4482d61bbb3SSatish Balay   }
4493649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
45026e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
4517c565772SBarry Smith   ierr = PetscLogFlops(32.0*a->nz - 4.0*a->nonzerorowcnt);CHKERRQ(ierr);
4522d61bbb3SSatish Balay   PetscFunctionReturn(0);
4532d61bbb3SSatish Balay }
4542d61bbb3SSatish Balay 
455dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_5(Mat A,Vec xx,Vec zz)
4562d61bbb3SSatish Balay {
4572d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
458d9fead3dSBarry Smith   PetscScalar       sum1,sum2,sum3,sum4,sum5,x1,x2,x3,x4,x5,*z = 0,*zarray;
459d9fead3dSBarry Smith   const PetscScalar *xb,*x;
460d9fead3dSBarry Smith   const MatScalar   *v;
461dfbe8321SBarry Smith   PetscErrorCode    ierr;
4620298fd71SBarry Smith   const PetscInt    *idx,*ii,*ridx=NULL;
4637c565772SBarry Smith   PetscInt          mbs,i,j,n;
464ace3abfcSBarry Smith   PetscBool         usecprow=a->compressedrow.use;
4652d61bbb3SSatish Balay 
466433994e6SBarry Smith   PetscFunctionBegin;
4673649974fSBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
46826e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
4692d61bbb3SSatish Balay 
4702d61bbb3SSatish Balay   idx = a->j;
4712d61bbb3SSatish Balay   v   = a->a;
47226e093fcSHong Zhang   if (usecprow) {
47326e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
47426e093fcSHong Zhang     ii   = a->compressedrow.i;
4757b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
476ef16e671SStefano Zampini     ierr = PetscMemzero(zarray,5*a->mbs*sizeof(PetscScalar));CHKERRQ(ierr);
47726e093fcSHong Zhang   } else {
47826e093fcSHong Zhang     mbs = a->mbs;
4792d61bbb3SSatish Balay     ii  = a->i;
48026e093fcSHong Zhang     z   = zarray;
48126e093fcSHong Zhang   }
4822d61bbb3SSatish Balay 
4832d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
4842d61bbb3SSatish Balay     n           = ii[1] - ii[0]; ii++;
4852d61bbb3SSatish Balay     sum1        = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0;
486444d8c10SJed Brown     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);     /* Indices for the next row (assumes same size as this one) */
487444d8c10SJed Brown     PetscPrefetchBlock(v+25*n,25*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
4882d61bbb3SSatish Balay     for (j=0; j<n; j++) {
4892d61bbb3SSatish Balay       xb    = x + 5*(*idx++);
4902d61bbb3SSatish Balay       x1    = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4];
4912d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[5]*x2 + v[10]*x3  + v[15]*x4 + v[20]*x5;
4922d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[6]*x2 + v[11]*x3  + v[16]*x4 + v[21]*x5;
4932d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[7]*x2 + v[12]*x3  + v[17]*x4 + v[22]*x5;
4942d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[8]*x2 + v[13]*x3  + v[18]*x4 + v[23]*x5;
4952d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[9]*x2 + v[14]*x3  + v[19]*x4 + v[24]*x5;
4962d61bbb3SSatish Balay       v    += 25;
4972d61bbb3SSatish Balay     }
4987b2bb3b9SHong Zhang     if (usecprow) z = zarray + 5*ridx[i];
4992d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5;
50026e093fcSHong Zhang     if (!usecprow) z += 5;
5012d61bbb3SSatish Balay   }
5023649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
50326e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
5047c565772SBarry Smith   ierr = PetscLogFlops(50.0*a->nz - 5.0*a->nonzerorowcnt);CHKERRQ(ierr);
5052d61bbb3SSatish Balay   PetscFunctionReturn(0);
5062d61bbb3SSatish Balay }
5072d61bbb3SSatish Balay 
50815091d37SBarry Smith 
509d6232bceSMichael Lange 
510dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_6(Mat A,Vec xx,Vec zz)
51115091d37SBarry Smith {
51215091d37SBarry Smith   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
513d9fead3dSBarry Smith   PetscScalar       *z = 0,sum1,sum2,sum3,sum4,sum5,sum6;
514d9fead3dSBarry Smith   const PetscScalar *x,*xb;
51526e093fcSHong Zhang   PetscScalar       x1,x2,x3,x4,x5,x6,*zarray;
516d9fead3dSBarry Smith   const MatScalar   *v;
517dfbe8321SBarry Smith   PetscErrorCode    ierr;
5187c565772SBarry Smith   PetscInt          mbs,i,*idx,*ii,j,n,*ridx=NULL;
519ace3abfcSBarry Smith   PetscBool         usecprow=a->compressedrow.use;
52015091d37SBarry Smith 
521433994e6SBarry Smith   PetscFunctionBegin;
5223649974fSBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
52326e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
52415091d37SBarry Smith 
52515091d37SBarry Smith   idx = a->j;
52615091d37SBarry Smith   v   = a->a;
52726e093fcSHong Zhang   if (usecprow) {
52826e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
52926e093fcSHong Zhang     ii   = a->compressedrow.i;
5307b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
531ef16e671SStefano Zampini     ierr = PetscMemzero(zarray,6*a->mbs*sizeof(PetscScalar));CHKERRQ(ierr);
53226e093fcSHong Zhang   } else {
53326e093fcSHong Zhang     mbs = a->mbs;
53415091d37SBarry Smith     ii  = a->i;
53526e093fcSHong Zhang     z   = zarray;
53626e093fcSHong Zhang   }
53715091d37SBarry Smith 
53815091d37SBarry Smith   for (i=0; i<mbs; i++) {
53926fbe8dcSKarl Rupp     n  = ii[1] - ii[0];
54026fbe8dcSKarl Rupp     ii++;
54126fbe8dcSKarl Rupp     sum1 = 0.0;
54226fbe8dcSKarl Rupp     sum2 = 0.0;
54326fbe8dcSKarl Rupp     sum3 = 0.0;
54426fbe8dcSKarl Rupp     sum4 = 0.0;
54526fbe8dcSKarl Rupp     sum5 = 0.0;
54626fbe8dcSKarl Rupp     sum6 = 0.0;
54726fbe8dcSKarl Rupp 
548444d8c10SJed Brown     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);     /* Indices for the next row (assumes same size as this one) */
549444d8c10SJed Brown     PetscPrefetchBlock(v+36*n,36*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
55015091d37SBarry Smith     for (j=0; j<n; j++) {
55115091d37SBarry Smith       xb    = x + 6*(*idx++);
55215091d37SBarry Smith       x1    = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5];
55315091d37SBarry Smith       sum1 += v[0]*x1 + v[6]*x2  + v[12]*x3  + v[18]*x4 + v[24]*x5 + v[30]*x6;
55415091d37SBarry Smith       sum2 += v[1]*x1 + v[7]*x2  + v[13]*x3  + v[19]*x4 + v[25]*x5 + v[31]*x6;
55515091d37SBarry Smith       sum3 += v[2]*x1 + v[8]*x2  + v[14]*x3  + v[20]*x4 + v[26]*x5 + v[32]*x6;
55615091d37SBarry Smith       sum4 += v[3]*x1 + v[9]*x2  + v[15]*x3  + v[21]*x4 + v[27]*x5 + v[33]*x6;
55715091d37SBarry Smith       sum5 += v[4]*x1 + v[10]*x2 + v[16]*x3  + v[22]*x4 + v[28]*x5 + v[34]*x6;
55815091d37SBarry Smith       sum6 += v[5]*x1 + v[11]*x2 + v[17]*x3  + v[23]*x4 + v[29]*x5 + v[35]*x6;
55915091d37SBarry Smith       v    += 36;
56015091d37SBarry Smith     }
5617b2bb3b9SHong Zhang     if (usecprow) z = zarray + 6*ridx[i];
56215091d37SBarry Smith     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6;
56326e093fcSHong Zhang     if (!usecprow) z += 6;
56415091d37SBarry Smith   }
56515091d37SBarry Smith 
5663649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
56726e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
5687c565772SBarry Smith   ierr = PetscLogFlops(72.0*a->nz - 6.0*a->nonzerorowcnt);CHKERRQ(ierr);
56915091d37SBarry Smith   PetscFunctionReturn(0);
57015091d37SBarry Smith }
5718ab949d8SShri Abhyankar 
572dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_7(Mat A,Vec xx,Vec zz)
5732d61bbb3SSatish Balay {
5742d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
575d9fead3dSBarry Smith   PetscScalar       *z = 0,sum1,sum2,sum3,sum4,sum5,sum6,sum7;
576d9fead3dSBarry Smith   const PetscScalar *x,*xb;
57726e093fcSHong Zhang   PetscScalar       x1,x2,x3,x4,x5,x6,x7,*zarray;
578d9fead3dSBarry Smith   const MatScalar   *v;
579dfbe8321SBarry Smith   PetscErrorCode    ierr;
5807c565772SBarry Smith   PetscInt          mbs,i,*idx,*ii,j,n,*ridx=NULL;
581ace3abfcSBarry Smith   PetscBool         usecprow=a->compressedrow.use;
5822d61bbb3SSatish Balay 
583433994e6SBarry Smith   PetscFunctionBegin;
5843649974fSBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
58526e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
5862d61bbb3SSatish Balay 
5872d61bbb3SSatish Balay   idx = a->j;
5882d61bbb3SSatish Balay   v   = a->a;
58926e093fcSHong Zhang   if (usecprow) {
59026e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
59126e093fcSHong Zhang     ii   = a->compressedrow.i;
5927b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
593ef16e671SStefano Zampini     ierr = PetscMemzero(zarray,7*a->mbs*sizeof(PetscScalar));CHKERRQ(ierr);
59426e093fcSHong Zhang   } else {
59526e093fcSHong Zhang     mbs = a->mbs;
5962d61bbb3SSatish Balay     ii  = a->i;
59726e093fcSHong Zhang     z   = zarray;
59826e093fcSHong Zhang   }
5992d61bbb3SSatish Balay 
6002d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
60126fbe8dcSKarl Rupp     n  = ii[1] - ii[0];
60226fbe8dcSKarl Rupp     ii++;
60326fbe8dcSKarl Rupp     sum1 = 0.0;
60426fbe8dcSKarl Rupp     sum2 = 0.0;
60526fbe8dcSKarl Rupp     sum3 = 0.0;
60626fbe8dcSKarl Rupp     sum4 = 0.0;
60726fbe8dcSKarl Rupp     sum5 = 0.0;
60826fbe8dcSKarl Rupp     sum6 = 0.0;
60926fbe8dcSKarl Rupp     sum7 = 0.0;
61026fbe8dcSKarl Rupp 
611444d8c10SJed Brown     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);     /* Indices for the next row (assumes same size as this one) */
612444d8c10SJed Brown     PetscPrefetchBlock(v+49*n,49*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
6132d61bbb3SSatish Balay     for (j=0; j<n; j++) {
6142d61bbb3SSatish Balay       xb    = x + 7*(*idx++);
6152d61bbb3SSatish Balay       x1    = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6];
6162d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[7]*x2  + v[14]*x3  + v[21]*x4 + v[28]*x5 + v[35]*x6 + v[42]*x7;
6172d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[8]*x2  + v[15]*x3  + v[22]*x4 + v[29]*x5 + v[36]*x6 + v[43]*x7;
6182d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[9]*x2  + v[16]*x3  + v[23]*x4 + v[30]*x5 + v[37]*x6 + v[44]*x7;
6192d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[10]*x2 + v[17]*x3  + v[24]*x4 + v[31]*x5 + v[38]*x6 + v[45]*x7;
6202d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[11]*x2 + v[18]*x3  + v[25]*x4 + v[32]*x5 + v[39]*x6 + v[46]*x7;
6212d61bbb3SSatish Balay       sum6 += v[5]*x1 + v[12]*x2 + v[19]*x3  + v[26]*x4 + v[33]*x5 + v[40]*x6 + v[47]*x7;
6222d61bbb3SSatish Balay       sum7 += v[6]*x1 + v[13]*x2 + v[20]*x3  + v[27]*x4 + v[34]*x5 + v[41]*x6 + v[48]*x7;
6232d61bbb3SSatish Balay       v    += 49;
6242d61bbb3SSatish Balay     }
6257b2bb3b9SHong Zhang     if (usecprow) z = zarray + 7*ridx[i];
6262d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
62726e093fcSHong Zhang     if (!usecprow) z += 7;
6282d61bbb3SSatish Balay   }
6292d61bbb3SSatish Balay 
6303649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
63126e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
6327c565772SBarry Smith   ierr = PetscLogFlops(98.0*a->nz - 7.0*a->nonzerorowcnt);CHKERRQ(ierr);
6332d61bbb3SSatish Balay   PetscFunctionReturn(0);
6342d61bbb3SSatish Balay }
6352d61bbb3SSatish Balay 
6368ab949d8SShri Abhyankar /* MatMult_SeqBAIJ_15 version 1: Columns in the block are accessed one at a time */
637832cc040SShri Abhyankar /* Default MatMult for block size 15 */
6388ab949d8SShri Abhyankar 
6398ab949d8SShri Abhyankar PetscErrorCode MatMult_SeqBAIJ_15_ver1(Mat A,Vec xx,Vec zz)
6408ab949d8SShri Abhyankar {
6418ab949d8SShri Abhyankar   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
6428ab949d8SShri Abhyankar   PetscScalar       *z = 0,sum1,sum2,sum3,sum4,sum5,sum6,sum7,sum8,sum9,sum10,sum11,sum12,sum13,sum14,sum15;
6438ab949d8SShri Abhyankar   const PetscScalar *x,*xb;
64453ef36baSBarry Smith   PetscScalar       *zarray,xv;
6458ab949d8SShri Abhyankar   const MatScalar   *v;
6468ab949d8SShri Abhyankar   PetscErrorCode    ierr;
6478ab949d8SShri Abhyankar   const PetscInt    *ii,*ij=a->j,*idx;
6487c565772SBarry Smith   PetscInt          mbs,i,j,k,n,*ridx=NULL;
649ace3abfcSBarry Smith   PetscBool         usecprow=a->compressedrow.use;
6508ab949d8SShri Abhyankar 
6518ab949d8SShri Abhyankar   PetscFunctionBegin;
6523649974fSBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
6538ab949d8SShri Abhyankar   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
6548ab949d8SShri Abhyankar 
6558ab949d8SShri Abhyankar   v = a->a;
6568ab949d8SShri Abhyankar   if (usecprow) {
6578ab949d8SShri Abhyankar     mbs  = a->compressedrow.nrows;
6588ab949d8SShri Abhyankar     ii   = a->compressedrow.i;
6598ab949d8SShri Abhyankar     ridx = a->compressedrow.rindex;
660ef16e671SStefano Zampini     ierr = PetscMemzero(zarray,15*a->mbs*sizeof(PetscScalar));CHKERRQ(ierr);
6618ab949d8SShri Abhyankar   } else {
6628ab949d8SShri Abhyankar     mbs = a->mbs;
6638ab949d8SShri Abhyankar     ii  = a->i;
6648ab949d8SShri Abhyankar     z   = zarray;
6658ab949d8SShri Abhyankar   }
6668ab949d8SShri Abhyankar 
6678ab949d8SShri Abhyankar   for (i=0; i<mbs; i++) {
6688ab949d8SShri Abhyankar     n    = ii[i+1] - ii[i];
6698ab949d8SShri Abhyankar     idx  = ij + ii[i];
6708ab949d8SShri Abhyankar     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0; sum7 = 0.0;
6718ab949d8SShri 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;
6728ab949d8SShri Abhyankar 
6738ab949d8SShri Abhyankar     for (j=0; j<n; j++) {
6748ab949d8SShri Abhyankar       xb = x + 15*(idx[j]);
6758ab949d8SShri Abhyankar 
6768ab949d8SShri Abhyankar       for (k=0; k<15; k++) {
67753ef36baSBarry Smith         xv     =  xb[k];
67853ef36baSBarry Smith         sum1  += v[0]*xv;
67953ef36baSBarry Smith         sum2  += v[1]*xv;
68053ef36baSBarry Smith         sum3  += v[2]*xv;
68153ef36baSBarry Smith         sum4  += v[3]*xv;
68253ef36baSBarry Smith         sum5  += v[4]*xv;
68353ef36baSBarry Smith         sum6  += v[5]*xv;
68453ef36baSBarry Smith         sum7  += v[6]*xv;
68553ef36baSBarry Smith         sum8  += v[7]*xv;
68653ef36baSBarry Smith         sum9  += v[8]*xv;
68753ef36baSBarry Smith         sum10 += v[9]*xv;
68853ef36baSBarry Smith         sum11 += v[10]*xv;
68953ef36baSBarry Smith         sum12 += v[11]*xv;
69053ef36baSBarry Smith         sum13 += v[12]*xv;
69153ef36baSBarry Smith         sum14 += v[13]*xv;
69253ef36baSBarry Smith         sum15 += v[14]*xv;
6938ab949d8SShri Abhyankar         v     += 15;
6948ab949d8SShri Abhyankar       }
6958ab949d8SShri Abhyankar     }
6968ab949d8SShri Abhyankar     if (usecprow) z = zarray + 15*ridx[i];
6978ab949d8SShri Abhyankar     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
6988ab949d8SShri 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;
6998ab949d8SShri Abhyankar 
7008ab949d8SShri Abhyankar     if (!usecprow) z += 15;
7018ab949d8SShri Abhyankar   }
7028ab949d8SShri Abhyankar 
7033649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
7048ab949d8SShri Abhyankar   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
7057c565772SBarry Smith   ierr = PetscLogFlops(450.0*a->nz - 15.0*a->nonzerorowcnt);CHKERRQ(ierr);
7068ab949d8SShri Abhyankar   PetscFunctionReturn(0);
7078ab949d8SShri Abhyankar }
7088ab949d8SShri Abhyankar 
7098ab949d8SShri Abhyankar /* MatMult_SeqBAIJ_15_ver2 : Columns in the block are accessed in sets of 4,4,4,3 */
7108ab949d8SShri Abhyankar PetscErrorCode MatMult_SeqBAIJ_15_ver2(Mat A,Vec xx,Vec zz)
7118ab949d8SShri Abhyankar {
7128ab949d8SShri Abhyankar   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
7138ab949d8SShri Abhyankar   PetscScalar       *z = 0,sum1,sum2,sum3,sum4,sum5,sum6,sum7,sum8,sum9,sum10,sum11,sum12,sum13,sum14,sum15;
7148ab949d8SShri Abhyankar   const PetscScalar *x,*xb;
7150b8f6341SShri Abhyankar   PetscScalar       x1,x2,x3,x4,*zarray;
7168ab949d8SShri Abhyankar   const MatScalar   *v;
7178ab949d8SShri Abhyankar   PetscErrorCode    ierr;
7188ab949d8SShri Abhyankar   const PetscInt    *ii,*ij=a->j,*idx;
7197c565772SBarry Smith   PetscInt          mbs,i,j,n,*ridx=NULL;
720ace3abfcSBarry Smith   PetscBool         usecprow=a->compressedrow.use;
7218ab949d8SShri Abhyankar 
7228ab949d8SShri Abhyankar   PetscFunctionBegin;
7233649974fSBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
7248ab949d8SShri Abhyankar   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
7258ab949d8SShri Abhyankar 
7268ab949d8SShri Abhyankar   v = a->a;
7278ab949d8SShri Abhyankar   if (usecprow) {
7288ab949d8SShri Abhyankar     mbs  = a->compressedrow.nrows;
7298ab949d8SShri Abhyankar     ii   = a->compressedrow.i;
7308ab949d8SShri Abhyankar     ridx = a->compressedrow.rindex;
731ef16e671SStefano Zampini     ierr = PetscMemzero(zarray,15*a->mbs*sizeof(PetscScalar));CHKERRQ(ierr);
7328ab949d8SShri Abhyankar   } else {
7338ab949d8SShri Abhyankar     mbs = a->mbs;
7348ab949d8SShri Abhyankar     ii  = a->i;
7358ab949d8SShri Abhyankar     z   = zarray;
7368ab949d8SShri Abhyankar   }
7378ab949d8SShri Abhyankar 
7388ab949d8SShri Abhyankar   for (i=0; i<mbs; i++) {
7398ab949d8SShri Abhyankar     n    = ii[i+1] - ii[i];
7408ab949d8SShri Abhyankar     idx  = ij + ii[i];
7418ab949d8SShri Abhyankar     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0; sum7 = 0.0;
7428ab949d8SShri 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;
7438ab949d8SShri Abhyankar 
7448ab949d8SShri Abhyankar     for (j=0; j<n; j++) {
7458ab949d8SShri Abhyankar       xb = x + 15*(idx[j]);
7460b8f6341SShri Abhyankar       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
7478ab949d8SShri Abhyankar 
7488ab949d8SShri Abhyankar       sum1  += v[0]*x1 + v[15]*x2 + v[30]*x3   + v[45]*x4;
7498ab949d8SShri Abhyankar       sum2  += v[1]*x1 + v[16]*x2 + v[31]*x3   + v[46]*x4;
7508ab949d8SShri Abhyankar       sum3  += v[2]*x1 + v[17]*x2 + v[32]*x3  + v[47]*x4;
7518ab949d8SShri Abhyankar       sum4  += v[3]*x1 + v[18]*x2 + v[33]*x3  + v[48]*x4;
7528ab949d8SShri Abhyankar       sum5  += v[4]*x1 + v[19]*x2 + v[34]*x3   + v[49]*x4;
7538ab949d8SShri Abhyankar       sum6  += v[5]*x1 + v[20]*x2 + v[35]*x3   + v[50]*x4;
7548ab949d8SShri Abhyankar       sum7  += v[6]*x1 + v[21]*x2 + v[36]*x3  + v[51]*x4;
7558ab949d8SShri Abhyankar       sum8  += v[7]*x1 + v[22]*x2 + v[37]*x3  + v[52]*x4;
7568ab949d8SShri Abhyankar       sum9  += v[8]*x1 + v[23]*x2 + v[38]*x3   + v[53]*x4;
7578ab949d8SShri Abhyankar       sum10 += v[9]*x1 + v[24]*x2 + v[39]*x3   + v[54]*x4;
7588ab949d8SShri Abhyankar       sum11 += v[10]*x1 + v[25]*x2 + v[40]*x3  + v[55]*x4;
7598ab949d8SShri Abhyankar       sum12 += v[11]*x1 + v[26]*x2 + v[41]*x3  + v[56]*x4;
7608ab949d8SShri Abhyankar       sum13 += v[12]*x1 + v[27]*x2 + v[42]*x3   + v[57]*x4;
7618ab949d8SShri Abhyankar       sum14 += v[13]*x1 + v[28]*x2 + v[43]*x3   + v[58]*x4;
7628ab949d8SShri Abhyankar       sum15 += v[14]*x1 + v[29]*x2 + v[44]*x3  + v[59]*x4;
7638ab949d8SShri Abhyankar 
7648ab949d8SShri Abhyankar       v += 60;
7658ab949d8SShri Abhyankar 
7660b8f6341SShri Abhyankar       x1 = xb[4]; x2 = xb[5]; x3 = xb[6]; x4 = xb[7];
7678ab949d8SShri Abhyankar 
7688ab949d8SShri Abhyankar       sum1  += v[0]*x1 + v[15]*x2 + v[30]*x3   + v[45]*x4;
7698ab949d8SShri Abhyankar       sum2  += v[1]*x1 + v[16]*x2 + v[31]*x3   + v[46]*x4;
7708ab949d8SShri Abhyankar       sum3  += v[2]*x1 + v[17]*x2 + v[32]*x3  + v[47]*x4;
7718ab949d8SShri Abhyankar       sum4  += v[3]*x1 + v[18]*x2 + v[33]*x3  + v[48]*x4;
7728ab949d8SShri Abhyankar       sum5  += v[4]*x1 + v[19]*x2 + v[34]*x3   + v[49]*x4;
7738ab949d8SShri Abhyankar       sum6  += v[5]*x1 + v[20]*x2 + v[35]*x3   + v[50]*x4;
7748ab949d8SShri Abhyankar       sum7  += v[6]*x1 + v[21]*x2 + v[36]*x3  + v[51]*x4;
7758ab949d8SShri Abhyankar       sum8  += v[7]*x1 + v[22]*x2 + v[37]*x3  + v[52]*x4;
7768ab949d8SShri Abhyankar       sum9  += v[8]*x1 + v[23]*x2 + v[38]*x3   + v[53]*x4;
7778ab949d8SShri Abhyankar       sum10 += v[9]*x1 + v[24]*x2 + v[39]*x3   + v[54]*x4;
7788ab949d8SShri Abhyankar       sum11 += v[10]*x1 + v[25]*x2 + v[40]*x3  + v[55]*x4;
7798ab949d8SShri Abhyankar       sum12 += v[11]*x1 + v[26]*x2 + v[41]*x3  + v[56]*x4;
7808ab949d8SShri Abhyankar       sum13 += v[12]*x1 + v[27]*x2 + v[42]*x3   + v[57]*x4;
7818ab949d8SShri Abhyankar       sum14 += v[13]*x1 + v[28]*x2 + v[43]*x3   + v[58]*x4;
7828ab949d8SShri Abhyankar       sum15 += v[14]*x1 + v[29]*x2 + v[44]*x3  + v[59]*x4;
7838ab949d8SShri Abhyankar       v     += 60;
7848ab949d8SShri Abhyankar 
7850b8f6341SShri Abhyankar       x1     = xb[8]; x2 = xb[9]; x3 = xb[10]; x4 = xb[11];
7860b8f6341SShri Abhyankar       sum1  += v[0]*x1 + v[15]*x2 + v[30]*x3   + v[45]*x4;
7870b8f6341SShri Abhyankar       sum2  += v[1]*x1 + v[16]*x2 + v[31]*x3   + v[46]*x4;
7880b8f6341SShri Abhyankar       sum3  += v[2]*x1 + v[17]*x2 + v[32]*x3  + v[47]*x4;
7890b8f6341SShri Abhyankar       sum4  += v[3]*x1 + v[18]*x2 + v[33]*x3  + v[48]*x4;
7900b8f6341SShri Abhyankar       sum5  += v[4]*x1 + v[19]*x2 + v[34]*x3   + v[49]*x4;
7910b8f6341SShri Abhyankar       sum6  += v[5]*x1 + v[20]*x2 + v[35]*x3   + v[50]*x4;
7920b8f6341SShri Abhyankar       sum7  += v[6]*x1 + v[21]*x2 + v[36]*x3  + v[51]*x4;
7930b8f6341SShri Abhyankar       sum8  += v[7]*x1 + v[22]*x2 + v[37]*x3  + v[52]*x4;
7940b8f6341SShri Abhyankar       sum9  += v[8]*x1 + v[23]*x2 + v[38]*x3   + v[53]*x4;
7950b8f6341SShri Abhyankar       sum10 += v[9]*x1 + v[24]*x2 + v[39]*x3   + v[54]*x4;
7960b8f6341SShri Abhyankar       sum11 += v[10]*x1 + v[25]*x2 + v[40]*x3  + v[55]*x4;
7970b8f6341SShri Abhyankar       sum12 += v[11]*x1 + v[26]*x2 + v[41]*x3  + v[56]*x4;
7980b8f6341SShri Abhyankar       sum13 += v[12]*x1 + v[27]*x2 + v[42]*x3   + v[57]*x4;
7990b8f6341SShri Abhyankar       sum14 += v[13]*x1 + v[28]*x2 + v[43]*x3   + v[58]*x4;
8000b8f6341SShri Abhyankar       sum15 += v[14]*x1 + v[29]*x2 + v[44]*x3  + v[59]*x4;
8010b8f6341SShri Abhyankar       v     += 60;
8020b8f6341SShri Abhyankar 
8030b8f6341SShri Abhyankar       x1     = xb[12]; x2 = xb[13]; x3 = xb[14];
8048ab949d8SShri Abhyankar       sum1  += v[0]*x1 + v[15]*x2 + v[30]*x3;
8058ab949d8SShri Abhyankar       sum2  += v[1]*x1 + v[16]*x2 + v[31]*x3;
8068ab949d8SShri Abhyankar       sum3  += v[2]*x1 + v[17]*x2 + v[32]*x3;
8078ab949d8SShri Abhyankar       sum4  += v[3]*x1 + v[18]*x2 + v[33]*x3;
8088ab949d8SShri Abhyankar       sum5  += v[4]*x1 + v[19]*x2 + v[34]*x3;
8098ab949d8SShri Abhyankar       sum6  += v[5]*x1 + v[20]*x2 + v[35]*x3;
8108ab949d8SShri Abhyankar       sum7  += v[6]*x1 + v[21]*x2 + v[36]*x3;
8118ab949d8SShri Abhyankar       sum8  += v[7]*x1 + v[22]*x2 + v[37]*x3;
8128ab949d8SShri Abhyankar       sum9  += v[8]*x1 + v[23]*x2 + v[38]*x3;
8138ab949d8SShri Abhyankar       sum10 += v[9]*x1 + v[24]*x2 + v[39]*x3;
8148ab949d8SShri Abhyankar       sum11 += v[10]*x1 + v[25]*x2 + v[40]*x3;
8158ab949d8SShri Abhyankar       sum12 += v[11]*x1 + v[26]*x2 + v[41]*x3;
8168ab949d8SShri Abhyankar       sum13 += v[12]*x1 + v[27]*x2 + v[42]*x3;
8178ab949d8SShri Abhyankar       sum14 += v[13]*x1 + v[28]*x2 + v[43]*x3;
8188ab949d8SShri Abhyankar       sum15 += v[14]*x1 + v[29]*x2 + v[44]*x3;
8198ab949d8SShri Abhyankar       v     += 45;
8208ab949d8SShri Abhyankar     }
8218ab949d8SShri Abhyankar     if (usecprow) z = zarray + 15*ridx[i];
8228ab949d8SShri Abhyankar     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
8238ab949d8SShri 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;
8248ab949d8SShri Abhyankar 
8258ab949d8SShri Abhyankar     if (!usecprow) z += 15;
8268ab949d8SShri Abhyankar   }
8278ab949d8SShri Abhyankar 
8283649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
8298ab949d8SShri Abhyankar   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
8307c565772SBarry Smith   ierr = PetscLogFlops(450.0*a->nz - 15.0*a->nonzerorowcnt);CHKERRQ(ierr);
8318ab949d8SShri Abhyankar   PetscFunctionReturn(0);
8328ab949d8SShri Abhyankar }
8338ab949d8SShri Abhyankar 
8348ab949d8SShri Abhyankar /* MatMult_SeqBAIJ_15_ver3 : Columns in the block are accessed in sets of 8,7 */
8358ab949d8SShri Abhyankar PetscErrorCode MatMult_SeqBAIJ_15_ver3(Mat A,Vec xx,Vec zz)
8368ab949d8SShri Abhyankar {
8378ab949d8SShri Abhyankar   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
8388ab949d8SShri Abhyankar   PetscScalar       *z = 0,sum1,sum2,sum3,sum4,sum5,sum6,sum7,sum8,sum9,sum10,sum11,sum12,sum13,sum14,sum15;
8398ab949d8SShri Abhyankar   const PetscScalar *x,*xb;
8400b8f6341SShri Abhyankar   PetscScalar       x1,x2,x3,x4,x5,x6,x7,x8,*zarray;
8418ab949d8SShri Abhyankar   const MatScalar   *v;
8428ab949d8SShri Abhyankar   PetscErrorCode    ierr;
8438ab949d8SShri Abhyankar   const PetscInt    *ii,*ij=a->j,*idx;
8447c565772SBarry Smith   PetscInt          mbs,i,j,n,*ridx=NULL;
845ace3abfcSBarry Smith   PetscBool         usecprow=a->compressedrow.use;
8468ab949d8SShri Abhyankar 
8478ab949d8SShri Abhyankar   PetscFunctionBegin;
8483649974fSBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
8498ab949d8SShri Abhyankar   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
8508ab949d8SShri Abhyankar 
8518ab949d8SShri Abhyankar   v = a->a;
8528ab949d8SShri Abhyankar   if (usecprow) {
8538ab949d8SShri Abhyankar     mbs  = a->compressedrow.nrows;
8548ab949d8SShri Abhyankar     ii   = a->compressedrow.i;
8558ab949d8SShri Abhyankar     ridx = a->compressedrow.rindex;
856ef16e671SStefano Zampini     ierr = PetscMemzero(zarray,15*a->mbs*sizeof(PetscScalar));CHKERRQ(ierr);
8578ab949d8SShri Abhyankar   } else {
8588ab949d8SShri Abhyankar     mbs = a->mbs;
8598ab949d8SShri Abhyankar     ii  = a->i;
8608ab949d8SShri Abhyankar     z   = zarray;
8618ab949d8SShri Abhyankar   }
8628ab949d8SShri Abhyankar 
8638ab949d8SShri Abhyankar   for (i=0; i<mbs; i++) {
8648ab949d8SShri Abhyankar     n    = ii[i+1] - ii[i];
8658ab949d8SShri Abhyankar     idx  = ij + ii[i];
8668ab949d8SShri Abhyankar     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0; sum7 = 0.0;
8678ab949d8SShri 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;
8688ab949d8SShri Abhyankar 
8698ab949d8SShri Abhyankar     for (j=0; j<n; j++) {
8708ab949d8SShri Abhyankar       xb = x + 15*(idx[j]);
8718ab949d8SShri Abhyankar       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6];
8720b8f6341SShri Abhyankar       x8 = xb[7];
8738ab949d8SShri Abhyankar 
8748ab949d8SShri 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;
8758ab949d8SShri 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;
8768ab949d8SShri 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;
8778ab949d8SShri 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;
8788ab949d8SShri 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;
8798ab949d8SShri 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;
8808ab949d8SShri 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;
8818ab949d8SShri 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;
8828ab949d8SShri 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;
8838ab949d8SShri 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;
8848ab949d8SShri 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;
8858ab949d8SShri 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;
8868ab949d8SShri 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;
8878ab949d8SShri 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;
8888ab949d8SShri 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;
8898ab949d8SShri Abhyankar       v     += 120;
8908ab949d8SShri Abhyankar 
8910b8f6341SShri Abhyankar       x1 = xb[8]; x2 = xb[9]; x3 = xb[10]; x4 = xb[11]; x5 = xb[12]; x6 = xb[13]; x7 = xb[14];
8920b8f6341SShri Abhyankar 
8938ab949d8SShri Abhyankar       sum1  += v[0]*x1 + v[15]*x2  + v[30]*x3  + v[45]*x4 + v[60]*x5 + v[75]*x6 + v[90]*x7;
8948ab949d8SShri Abhyankar       sum2  += v[1]*x1 + v[16]*x2  + v[31]*x3  + v[46]*x4 + v[61]*x5 + v[76]*x6 + v[91]*x7;
8958ab949d8SShri Abhyankar       sum3  += v[2]*x1 + v[17]*x2  + v[32]*x3  + v[47]*x4 + v[62]*x5 + v[77]*x6 + v[92]*x7;
8968ab949d8SShri Abhyankar       sum4  += v[3]*x1 + v[18]*x2 + v[33]*x3  + v[48]*x4 + v[63]*x5 + v[78]*x6 + v[93]*x7;
8978ab949d8SShri Abhyankar       sum5  += v[4]*x1 + v[19]*x2 + v[34]*x3  + v[49]*x4 + v[64]*x5 + v[79]*x6 + v[94]*x7;
8988ab949d8SShri Abhyankar       sum6  += v[5]*x1 + v[20]*x2 + v[35]*x3  + v[50]*x4 + v[65]*x5 + v[80]*x6 + v[95]*x7;
8998ab949d8SShri Abhyankar       sum7  += v[6]*x1 + v[21]*x2 + v[36]*x3  + v[51]*x4 + v[66]*x5 + v[81]*x6 + v[96]*x7;
9008ab949d8SShri Abhyankar       sum8  += v[7]*x1 + v[22]*x2  + v[37]*x3  + v[52]*x4 + v[67]*x5 + v[82]*x6 + v[97]*x7;
9018ab949d8SShri Abhyankar       sum9  += v[8]*x1 + v[23]*x2  + v[38]*x3  + v[53]*x4 + v[68]*x5 + v[83]*x6 + v[98]*x7;
9028ab949d8SShri Abhyankar       sum10 += v[9]*x1 + v[24]*x2  + v[39]*x3  + v[54]*x4 + v[69]*x5 + v[84]*x6 + v[99]*x7;
9038ab949d8SShri Abhyankar       sum11 += v[10]*x1 + v[25]*x2 + v[40]*x3  + v[55]*x4 + v[70]*x5 + v[85]*x6 + v[100]*x7;
9048ab949d8SShri Abhyankar       sum12 += v[11]*x1 + v[26]*x2 + v[41]*x3  + v[56]*x4 + v[71]*x5 + v[86]*x6 + v[101]*x7;
9058ab949d8SShri Abhyankar       sum13 += v[12]*x1 + v[27]*x2 + v[42]*x3  + v[57]*x4 + v[72]*x5 + v[87]*x6 + v[102]*x7;
9068ab949d8SShri Abhyankar       sum14 += v[13]*x1 + v[28]*x2 + v[43]*x3  + v[58]*x4 + v[73]*x5 + v[88]*x6 + v[103]*x7;
9078ab949d8SShri Abhyankar       sum15 += v[14]*x1 + v[29]*x2 + v[44]*x3  + v[59]*x4 + v[74]*x5 + v[89]*x6 + v[104]*x7;
9088ab949d8SShri Abhyankar       v     += 105;
9098ab949d8SShri Abhyankar     }
9108ab949d8SShri Abhyankar     if (usecprow) z = zarray + 15*ridx[i];
9118ab949d8SShri Abhyankar     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
9128ab949d8SShri 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;
9138ab949d8SShri Abhyankar 
9148ab949d8SShri Abhyankar     if (!usecprow) z += 15;
9158ab949d8SShri Abhyankar   }
9168ab949d8SShri Abhyankar 
9173649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
9188ab949d8SShri Abhyankar   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
9197c565772SBarry Smith   ierr = PetscLogFlops(450.0*a->nz - 15.0*a->nonzerorowcnt);CHKERRQ(ierr);
9208ab949d8SShri Abhyankar   PetscFunctionReturn(0);
9218ab949d8SShri Abhyankar }
9228ab949d8SShri Abhyankar 
9238ab949d8SShri Abhyankar /* MatMult_SeqBAIJ_15_ver4 : All columns in the block are accessed at once */
9248ab949d8SShri Abhyankar 
9258ab949d8SShri Abhyankar PetscErrorCode MatMult_SeqBAIJ_15_ver4(Mat A,Vec xx,Vec zz)
9268ab949d8SShri Abhyankar {
9278ab949d8SShri Abhyankar   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
9288ab949d8SShri Abhyankar   PetscScalar       *z = 0,sum1,sum2,sum3,sum4,sum5,sum6,sum7,sum8,sum9,sum10,sum11,sum12,sum13,sum14,sum15;
9298ab949d8SShri Abhyankar   const PetscScalar *x,*xb;
9308ab949d8SShri Abhyankar   PetscScalar       x1,x2,x3,x4,x5,x6,x7,x8,x9,x10,x11,x12,x13,x14,x15,*zarray;
9318ab949d8SShri Abhyankar   const MatScalar   *v;
9328ab949d8SShri Abhyankar   PetscErrorCode    ierr;
9338ab949d8SShri Abhyankar   const PetscInt    *ii,*ij=a->j,*idx;
9347c565772SBarry Smith   PetscInt          mbs,i,j,n,*ridx=NULL;
935ace3abfcSBarry Smith   PetscBool         usecprow=a->compressedrow.use;
9368ab949d8SShri Abhyankar 
9378ab949d8SShri Abhyankar   PetscFunctionBegin;
9383649974fSBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
9398ab949d8SShri Abhyankar   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
9408ab949d8SShri Abhyankar 
9418ab949d8SShri Abhyankar   v = a->a;
9428ab949d8SShri Abhyankar   if (usecprow) {
9438ab949d8SShri Abhyankar     mbs  = a->compressedrow.nrows;
9448ab949d8SShri Abhyankar     ii   = a->compressedrow.i;
9458ab949d8SShri Abhyankar     ridx = a->compressedrow.rindex;
946ef16e671SStefano Zampini     ierr = PetscMemzero(zarray,15*a->mbs*sizeof(PetscScalar));CHKERRQ(ierr);
9478ab949d8SShri Abhyankar   } else {
9488ab949d8SShri Abhyankar     mbs = a->mbs;
9498ab949d8SShri Abhyankar     ii  = a->i;
9508ab949d8SShri Abhyankar     z   = zarray;
9518ab949d8SShri Abhyankar   }
9528ab949d8SShri Abhyankar 
9538ab949d8SShri Abhyankar   for (i=0; i<mbs; i++) {
9548ab949d8SShri Abhyankar     n    = ii[i+1] - ii[i];
9558ab949d8SShri Abhyankar     idx  = ij + ii[i];
9568ab949d8SShri Abhyankar     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0; sum7 = 0.0;
9578ab949d8SShri 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;
9588ab949d8SShri Abhyankar 
9598ab949d8SShri Abhyankar     for (j=0; j<n; j++) {
9608ab949d8SShri Abhyankar       xb = x + 15*(idx[j]);
9618ab949d8SShri Abhyankar       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6];
9628ab949d8SShri 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];
9638ab949d8SShri Abhyankar 
9648ab949d8SShri 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;
9658ab949d8SShri 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;
9668ab949d8SShri 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;
9678ab949d8SShri 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;
9688ab949d8SShri 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;
9698ab949d8SShri 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;
9708ab949d8SShri 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;
9718ab949d8SShri 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;
9728ab949d8SShri 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;
9738ab949d8SShri 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;
9748ab949d8SShri 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;
9758ab949d8SShri 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;
9768ab949d8SShri 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;
9778ab949d8SShri 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;
9788ab949d8SShri 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;
9798ab949d8SShri Abhyankar       v     += 225;
9808ab949d8SShri Abhyankar     }
9818ab949d8SShri Abhyankar     if (usecprow) z = zarray + 15*ridx[i];
9828ab949d8SShri Abhyankar     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
9838ab949d8SShri 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;
9848ab949d8SShri Abhyankar 
9858ab949d8SShri Abhyankar     if (!usecprow) z += 15;
9868ab949d8SShri Abhyankar   }
9878ab949d8SShri Abhyankar 
9883649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
9898ab949d8SShri Abhyankar   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
9907c565772SBarry Smith   ierr = PetscLogFlops(450.0*a->nz - 15.0*a->nonzerorowcnt);CHKERRQ(ierr);
9918ab949d8SShri Abhyankar   PetscFunctionReturn(0);
9928ab949d8SShri Abhyankar }
9938ab949d8SShri Abhyankar 
9948ab949d8SShri Abhyankar 
9953f1db9ecSBarry Smith /*
9963f1db9ecSBarry Smith     This will not work with MatScalar == float because it calls the BLAS
9973f1db9ecSBarry Smith */
998dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_N(Mat A,Vec xx,Vec zz)
9992d61bbb3SSatish Balay {
10002d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
1001d9ca1df4SBarry Smith   PetscScalar       *z = 0,*work,*workt,*zarray;
1002d9ca1df4SBarry Smith   const PetscScalar *x,*xb;
1003d9ca1df4SBarry Smith   const MatScalar   *v;
1004dfbe8321SBarry Smith   PetscErrorCode    ierr;
1005d9ca1df4SBarry Smith   PetscInt          mbs,i,bs=A->rmap->bs,j,n,bs2=a->bs2;
1006d9ca1df4SBarry Smith   const PetscInt    *idx,*ii,*ridx=NULL;
1007d9ca1df4SBarry Smith   PetscInt          ncols,k;
1008ace3abfcSBarry Smith   PetscBool         usecprow=a->compressedrow.use;
10092d61bbb3SSatish Balay 
10102d61bbb3SSatish Balay   PetscFunctionBegin;
1011d9ca1df4SBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
101226e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
10132d61bbb3SSatish Balay 
10142d61bbb3SSatish Balay   idx = a->j;
10152d61bbb3SSatish Balay   v   = a->a;
101626e093fcSHong Zhang   if (usecprow) {
101726e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
101826e093fcSHong Zhang     ii   = a->compressedrow.i;
10197b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
1020ef16e671SStefano Zampini     ierr = PetscMemzero(zarray,bs*a->mbs*sizeof(PetscScalar));CHKERRQ(ierr);
102126e093fcSHong Zhang   } else {
102226e093fcSHong Zhang     mbs = a->mbs;
10232d61bbb3SSatish Balay     ii  = a->i;
102426e093fcSHong Zhang     z   = zarray;
102526e093fcSHong Zhang   }
1026218c64b6SSatish Balay 
10272d61bbb3SSatish Balay   if (!a->mult_work) {
1028d0f46423SBarry Smith     k    = PetscMax(A->rmap->n,A->cmap->n);
1029854ce69bSBarry Smith     ierr = PetscMalloc1(k+1,&a->mult_work);CHKERRQ(ierr);
10302d61bbb3SSatish Balay   }
10312d61bbb3SSatish Balay   work = a->mult_work;
10322d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
10332d61bbb3SSatish Balay     n           = ii[1] - ii[0]; ii++;
10342d61bbb3SSatish Balay     ncols       = n*bs;
10352d61bbb3SSatish Balay     workt       = work;
10362d61bbb3SSatish Balay     for (j=0; j<n; j++) {
10372d61bbb3SSatish Balay       xb = x + bs*(*idx++);
10382d61bbb3SSatish Balay       for (k=0; k<bs; k++) workt[k] = xb[k];
10392d61bbb3SSatish Balay       workt += bs;
10402d61bbb3SSatish Balay     }
10417b2bb3b9SHong Zhang     if (usecprow) z = zarray + bs*ridx[i];
104296b95a6bSBarry Smith     PetscKernel_w_gets_Ar_times_v(bs,ncols,work,v,z);
104371044d3cSBarry Smith     /* BLASgemv_("N",&bs,&ncols,&_DOne,v,&bs,work,&_One,&_DZero,z,&_One); */
10442d61bbb3SSatish Balay     v += n*bs2;
104526e093fcSHong Zhang     if (!usecprow) z += bs;
10462d61bbb3SSatish Balay   }
1047d9ca1df4SBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
104826e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
10497c565772SBarry Smith   ierr = PetscLogFlops(2.0*a->nz*bs2 - bs*a->nonzerorowcnt);CHKERRQ(ierr);
10502d61bbb3SSatish Balay   PetscFunctionReturn(0);
10512d61bbb3SSatish Balay }
10522d61bbb3SSatish Balay 
1053dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_1(Mat A,Vec xx,Vec yy,Vec zz)
10542d61bbb3SSatish Balay {
10552d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
1056122f12eaSBarry Smith   const PetscScalar *x;
1057122f12eaSBarry Smith   PetscScalar       *y,*z,sum;
1058122f12eaSBarry Smith   const MatScalar   *v;
1059dfbe8321SBarry Smith   PetscErrorCode    ierr;
10607c565772SBarry Smith   PetscInt          mbs=a->mbs,i,n,*ridx=NULL;
1061122f12eaSBarry Smith   const PetscInt    *idx,*ii;
1062ace3abfcSBarry Smith   PetscBool         usecprow=a->compressedrow.use;
10632d61bbb3SSatish Balay 
10642d61bbb3SSatish Balay   PetscFunctionBegin;
10653649974fSBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
1066d9ca1df4SBarry Smith   ierr = VecGetArrayPair(yy,zz,&y,&z);CHKERRQ(ierr);
10672d61bbb3SSatish Balay 
10682d61bbb3SSatish Balay   idx = a->j;
10692d61bbb3SSatish Balay   v   = a->a;
107026e093fcSHong Zhang   if (usecprow) {
10714eb6d288SHong Zhang     if (zz != yy) {
1072122f12eaSBarry Smith       ierr = PetscMemcpy(z,y,mbs*sizeof(PetscScalar));CHKERRQ(ierr);
10734eb6d288SHong Zhang     }
107426e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
107526e093fcSHong Zhang     ii   = a->compressedrow.i;
10767b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
107726e093fcSHong Zhang   } else {
10782d61bbb3SSatish Balay     ii = a->i;
107926e093fcSHong Zhang   }
10802d61bbb3SSatish Balay 
10812d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
1082122f12eaSBarry Smith     n = ii[1] - ii[0];
1083122f12eaSBarry Smith     ii++;
108426e093fcSHong Zhang     if (!usecprow) {
1085122f12eaSBarry Smith       sum         = y[i];
1086122f12eaSBarry Smith     } else {
1087122f12eaSBarry Smith       sum = y[ridx[i]];
1088122f12eaSBarry Smith     }
1089444d8c10SJed Brown     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA); /* Indices for the next row (assumes same size as this one) */
1090444d8c10SJed Brown     PetscPrefetchBlock(v+n,n,0,PETSC_PREFETCH_HINT_NTA);   /* Entries for the next row */
1091122f12eaSBarry Smith     PetscSparseDensePlusDot(sum,x,v,idx,n);
1092122f12eaSBarry Smith     v   += n;
1093122f12eaSBarry Smith     idx += n;
1094122f12eaSBarry Smith     if (usecprow) {
1095122f12eaSBarry Smith       z[ridx[i]] = sum;
1096122f12eaSBarry Smith     } else {
1097122f12eaSBarry Smith       z[i] = sum;
109826e093fcSHong Zhang     }
10992d61bbb3SSatish Balay   }
11003649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
1101d9ca1df4SBarry Smith   ierr = VecRestoreArrayPair(yy,zz,&y,&z);CHKERRQ(ierr);
11027c565772SBarry Smith   ierr = PetscLogFlops(2.0*a->nz);CHKERRQ(ierr);
11032d61bbb3SSatish Balay   PetscFunctionReturn(0);
11042d61bbb3SSatish Balay }
11052d61bbb3SSatish Balay 
1106dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_2(Mat A,Vec xx,Vec yy,Vec zz)
11072d61bbb3SSatish Balay {
11082d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
1109d9ca1df4SBarry Smith   PetscScalar       *y = 0,*z = 0,sum1,sum2;
1110d9ca1df4SBarry Smith   const PetscScalar *x,*xb;
111126e093fcSHong Zhang   PetscScalar       x1,x2,*yarray,*zarray;
1112d9ca1df4SBarry Smith   const MatScalar   *v;
1113dfbe8321SBarry Smith   PetscErrorCode    ierr;
1114d9ca1df4SBarry Smith   PetscInt          mbs = a->mbs,i,n,j;
1115d9ca1df4SBarry Smith   const PetscInt    *idx,*ii,*ridx = NULL;
1116ace3abfcSBarry Smith   PetscBool         usecprow = a->compressedrow.use;
11172d61bbb3SSatish Balay 
11182d61bbb3SSatish Balay   PetscFunctionBegin;
1119d9ca1df4SBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
1120d9ca1df4SBarry Smith   ierr = VecGetArrayPair(yy,zz,&yarray,&zarray);CHKERRQ(ierr);
11212d61bbb3SSatish Balay 
11222d61bbb3SSatish Balay   idx = a->j;
11232d61bbb3SSatish Balay   v   = a->a;
112426e093fcSHong Zhang   if (usecprow) {
11254eb6d288SHong Zhang     if (zz != yy) {
11264eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,2*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
11274eb6d288SHong Zhang     }
112826e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
112926e093fcSHong Zhang     ii   = a->compressedrow.i;
11307b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
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++) {
114726fbe8dcSKarl Rupp       xb = x + 2*(*idx++);
114826fbe8dcSKarl Rupp       x1 = xb[0];
114926fbe8dcSKarl Rupp       x2 = xb[1];
115026fbe8dcSKarl Rupp 
11512d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[2]*x2;
11522d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[3]*x2;
11532d61bbb3SSatish Balay       v    += 4;
11542d61bbb3SSatish Balay     }
11552d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2;
115626e093fcSHong Zhang     if (!usecprow) {
11572d61bbb3SSatish Balay       z += 2; y += 2;
11582d61bbb3SSatish Balay     }
115926e093fcSHong Zhang   }
1160d9ca1df4SBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
1161d9ca1df4SBarry Smith   ierr = VecRestoreArrayPair(yy,zz,&yarray,&zarray);CHKERRQ(ierr);
1162dc0b31edSSatish Balay   ierr = PetscLogFlops(4.0*a->nz);CHKERRQ(ierr);
11632d61bbb3SSatish Balay   PetscFunctionReturn(0);
11642d61bbb3SSatish Balay }
11652d61bbb3SSatish Balay 
1166dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_3(Mat A,Vec xx,Vec yy,Vec zz)
11672d61bbb3SSatish Balay {
11682d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
1169d9ca1df4SBarry Smith   PetscScalar       *y = 0,*z = 0,sum1,sum2,sum3,x1,x2,x3,*yarray,*zarray;
1170d9ca1df4SBarry Smith   const PetscScalar *x,*xb;
1171d9ca1df4SBarry Smith   const MatScalar   *v;
1172dfbe8321SBarry Smith   PetscErrorCode    ierr;
1173d9ca1df4SBarry Smith   PetscInt          mbs = a->mbs,i,j,n;
1174d9ca1df4SBarry Smith   const PetscInt    *idx,*ii,*ridx = NULL;
1175ace3abfcSBarry Smith   PetscBool         usecprow = a->compressedrow.use;
11762d61bbb3SSatish Balay 
11772d61bbb3SSatish Balay   PetscFunctionBegin;
1178d9ca1df4SBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
1179d9ca1df4SBarry Smith   ierr = VecGetArrayPair(yy,zz,&yarray,&zarray);CHKERRQ(ierr);
11802d61bbb3SSatish Balay 
11812d61bbb3SSatish Balay   idx = a->j;
11822d61bbb3SSatish Balay   v   = a->a;
118326e093fcSHong Zhang   if (usecprow) {
11844eb6d288SHong Zhang     if (zz != yy) {
11854eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,3*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
11864eb6d288SHong Zhang     }
118726e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
118826e093fcSHong Zhang     ii   = a->compressedrow.i;
11897b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
119026e093fcSHong Zhang   } else {
11912d61bbb3SSatish Balay     ii = a->i;
119226e093fcSHong Zhang     y  = yarray;
119326e093fcSHong Zhang     z  = zarray;
119426e093fcSHong Zhang   }
11952d61bbb3SSatish Balay 
11962d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
11972d61bbb3SSatish Balay     n = ii[1] - ii[0]; ii++;
119826e093fcSHong Zhang     if (usecprow) {
11997b2bb3b9SHong Zhang       z = zarray + 3*ridx[i];
12007b2bb3b9SHong Zhang       y = yarray + 3*ridx[i];
120126e093fcSHong Zhang     }
12022d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2];
1203444d8c10SJed Brown     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);   /* Indices for the next row (assumes same size as this one) */
1204444d8c10SJed Brown     PetscPrefetchBlock(v+9*n,9*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
12052d61bbb3SSatish Balay     for (j=0; j<n; j++) {
12062d61bbb3SSatish Balay       xb    = x + 3*(*idx++); x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
12072d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[3]*x2 + v[6]*x3;
12082d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[4]*x2 + v[7]*x3;
12092d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[5]*x2 + v[8]*x3;
12102d61bbb3SSatish Balay       v    += 9;
12112d61bbb3SSatish Balay     }
12122d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3;
121326e093fcSHong Zhang     if (!usecprow) {
12142d61bbb3SSatish Balay       z += 3; y += 3;
12152d61bbb3SSatish Balay     }
121626e093fcSHong Zhang   }
1217d9ca1df4SBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
1218d9ca1df4SBarry Smith   ierr = VecRestoreArrayPair(yy,zz,&yarray,&zarray);CHKERRQ(ierr);
1219dc0b31edSSatish Balay   ierr = PetscLogFlops(18.0*a->nz);CHKERRQ(ierr);
12202d61bbb3SSatish Balay   PetscFunctionReturn(0);
12212d61bbb3SSatish Balay }
12222d61bbb3SSatish Balay 
1223dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_4(Mat A,Vec xx,Vec yy,Vec zz)
12242d61bbb3SSatish Balay {
12252d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
1226d9ca1df4SBarry Smith   PetscScalar       *y = 0,*z = 0,sum1,sum2,sum3,sum4,x1,x2,x3,x4,*yarray,*zarray;
1227d9ca1df4SBarry Smith   const PetscScalar *x,*xb;
1228d9ca1df4SBarry Smith   const MatScalar   *v;
1229dfbe8321SBarry Smith   PetscErrorCode    ierr;
1230d9ca1df4SBarry Smith   PetscInt          mbs = a->mbs,i,j,n;
1231d9ca1df4SBarry Smith   const PetscInt    *idx,*ii,*ridx=NULL;
1232ace3abfcSBarry Smith   PetscBool         usecprow=a->compressedrow.use;
12332d61bbb3SSatish Balay 
12342d61bbb3SSatish Balay   PetscFunctionBegin;
1235d9ca1df4SBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
1236d9ca1df4SBarry Smith   ierr = VecGetArrayPair(yy,zz,&yarray,&zarray);CHKERRQ(ierr);
12372d61bbb3SSatish Balay 
12382d61bbb3SSatish Balay   idx = a->j;
12392d61bbb3SSatish Balay   v   = a->a;
124026e093fcSHong Zhang   if (usecprow) {
12414eb6d288SHong Zhang     if (zz != yy) {
12424eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,4*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
12434eb6d288SHong Zhang     }
124426e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
124526e093fcSHong Zhang     ii   = a->compressedrow.i;
12467b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
124726e093fcSHong Zhang   } else {
12482d61bbb3SSatish Balay     ii = a->i;
124926e093fcSHong Zhang     y  = yarray;
125026e093fcSHong Zhang     z  = zarray;
125126e093fcSHong Zhang   }
12522d61bbb3SSatish Balay 
12532d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
12542d61bbb3SSatish Balay     n = ii[1] - ii[0]; ii++;
125526e093fcSHong Zhang     if (usecprow) {
12567b2bb3b9SHong Zhang       z = zarray + 4*ridx[i];
12577b2bb3b9SHong Zhang       y = yarray + 4*ridx[i];
125826e093fcSHong Zhang     }
12592d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3];
1260444d8c10SJed Brown     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);     /* Indices for the next row (assumes same size as this one) */
1261444d8c10SJed Brown     PetscPrefetchBlock(v+16*n,16*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
12622d61bbb3SSatish Balay     for (j=0; j<n; j++) {
12632d61bbb3SSatish Balay       xb    = x + 4*(*idx++);
12642d61bbb3SSatish Balay       x1    = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
12652d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[4]*x2 + v[8]*x3   + v[12]*x4;
12662d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[5]*x2 + v[9]*x3   + v[13]*x4;
12672d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[6]*x2 + v[10]*x3  + v[14]*x4;
12682d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[7]*x2 + v[11]*x3  + v[15]*x4;
12692d61bbb3SSatish Balay       v    += 16;
12702d61bbb3SSatish Balay     }
12712d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4;
127226e093fcSHong Zhang     if (!usecprow) {
12732d61bbb3SSatish Balay       z += 4; y += 4;
12742d61bbb3SSatish Balay     }
127526e093fcSHong Zhang   }
1276d9ca1df4SBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
1277d9ca1df4SBarry Smith   ierr = VecRestoreArrayPair(yy,zz,&yarray,&zarray);CHKERRQ(ierr);
1278dc0b31edSSatish Balay   ierr = PetscLogFlops(32.0*a->nz);CHKERRQ(ierr);
12792d61bbb3SSatish Balay   PetscFunctionReturn(0);
12802d61bbb3SSatish Balay }
12812d61bbb3SSatish Balay 
1282dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_5(Mat A,Vec xx,Vec yy,Vec zz)
12832d61bbb3SSatish Balay {
12842d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
1285d9ca1df4SBarry Smith   PetscScalar       *y = 0,*z = 0,sum1,sum2,sum3,sum4,sum5,x1,x2,x3,x4,x5;
1286d9ca1df4SBarry Smith   const PetscScalar *x,*xb;
128726e093fcSHong Zhang   PetscScalar       *yarray,*zarray;
1288d9ca1df4SBarry Smith   const MatScalar   *v;
1289dfbe8321SBarry Smith   PetscErrorCode    ierr;
1290d9ca1df4SBarry Smith   PetscInt          mbs = a->mbs,i,j,n;
1291d9ca1df4SBarry Smith   const PetscInt    *idx,*ii,*ridx = NULL;
1292ace3abfcSBarry Smith   PetscBool         usecprow=a->compressedrow.use;
12932d61bbb3SSatish Balay 
12942d61bbb3SSatish Balay   PetscFunctionBegin;
1295d9ca1df4SBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
1296d9ca1df4SBarry Smith   ierr = VecGetArrayPair(yy,zz,&yarray,&zarray);CHKERRQ(ierr);
12972d61bbb3SSatish Balay 
12982d61bbb3SSatish Balay   idx = a->j;
12992d61bbb3SSatish Balay   v   = a->a;
130026e093fcSHong Zhang   if (usecprow) {
13014eb6d288SHong Zhang     if (zz != yy) {
13024eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,5*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
13034eb6d288SHong Zhang     }
130426e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
130526e093fcSHong Zhang     ii   = a->compressedrow.i;
13067b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
130726e093fcSHong Zhang   } else {
13082d61bbb3SSatish Balay     ii = a->i;
130926e093fcSHong Zhang     y  = yarray;
131026e093fcSHong Zhang     z  = zarray;
131126e093fcSHong Zhang   }
13122d61bbb3SSatish Balay 
13132d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
13142d61bbb3SSatish Balay     n = ii[1] - ii[0]; ii++;
131526e093fcSHong Zhang     if (usecprow) {
13167b2bb3b9SHong Zhang       z = zarray + 5*ridx[i];
13177b2bb3b9SHong Zhang       y = yarray + 5*ridx[i];
131826e093fcSHong Zhang     }
13192d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4];
1320444d8c10SJed Brown     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);     /* Indices for the next row (assumes same size as this one) */
1321444d8c10SJed Brown     PetscPrefetchBlock(v+25*n,25*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
13222d61bbb3SSatish Balay     for (j=0; j<n; j++) {
13232d61bbb3SSatish Balay       xb    = x + 5*(*idx++);
13242d61bbb3SSatish Balay       x1    = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4];
13252d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[5]*x2 + v[10]*x3  + v[15]*x4 + v[20]*x5;
13262d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[6]*x2 + v[11]*x3  + v[16]*x4 + v[21]*x5;
13272d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[7]*x2 + v[12]*x3  + v[17]*x4 + v[22]*x5;
13282d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[8]*x2 + v[13]*x3  + v[18]*x4 + v[23]*x5;
13292d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[9]*x2 + v[14]*x3  + v[19]*x4 + v[24]*x5;
13302d61bbb3SSatish Balay       v    += 25;
13312d61bbb3SSatish Balay     }
13322d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5;
133326e093fcSHong Zhang     if (!usecprow) {
13342d61bbb3SSatish Balay       z += 5; y += 5;
13352d61bbb3SSatish Balay     }
133626e093fcSHong Zhang   }
1337d9ca1df4SBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
1338d9ca1df4SBarry Smith   ierr = VecRestoreArrayPair(yy,zz,&yarray,&zarray);CHKERRQ(ierr);
1339dc0b31edSSatish Balay   ierr = PetscLogFlops(50.0*a->nz);CHKERRQ(ierr);
13402d61bbb3SSatish Balay   PetscFunctionReturn(0);
13412d61bbb3SSatish Balay }
1342dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_6(Mat A,Vec xx,Vec yy,Vec zz)
134315091d37SBarry Smith {
134415091d37SBarry Smith   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
1345d9ca1df4SBarry Smith   PetscScalar       *y = 0,*z = 0,sum1,sum2,sum3,sum4,sum5,sum6;
1346d9ca1df4SBarry Smith   const PetscScalar *x,*xb;
134726e093fcSHong Zhang   PetscScalar       x1,x2,x3,x4,x5,x6,*yarray,*zarray;
1348d9ca1df4SBarry Smith   const MatScalar   *v;
1349dfbe8321SBarry Smith   PetscErrorCode    ierr;
1350d9ca1df4SBarry Smith   PetscInt          mbs = a->mbs,i,j,n;
1351d9ca1df4SBarry Smith   const PetscInt    *idx,*ii,*ridx=NULL;
1352ace3abfcSBarry Smith   PetscBool         usecprow=a->compressedrow.use;
135315091d37SBarry Smith 
135415091d37SBarry Smith   PetscFunctionBegin;
1355d9ca1df4SBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
1356d9ca1df4SBarry Smith   ierr = VecGetArrayPair(yy,zz,&yarray,&zarray);CHKERRQ(ierr);
135715091d37SBarry Smith 
135815091d37SBarry Smith   idx = a->j;
135915091d37SBarry Smith   v   = a->a;
136026e093fcSHong Zhang   if (usecprow) {
13614eb6d288SHong Zhang     if (zz != yy) {
13624eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,6*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
13634eb6d288SHong Zhang     }
136426e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
136526e093fcSHong Zhang     ii   = a->compressedrow.i;
13667b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
136726e093fcSHong Zhang   } else {
136815091d37SBarry Smith     ii = a->i;
136926e093fcSHong Zhang     y  = yarray;
137026e093fcSHong Zhang     z  = zarray;
137126e093fcSHong Zhang   }
137215091d37SBarry Smith 
137315091d37SBarry Smith   for (i=0; i<mbs; i++) {
137415091d37SBarry Smith     n = ii[1] - ii[0]; ii++;
137526e093fcSHong Zhang     if (usecprow) {
13767b2bb3b9SHong Zhang       z = zarray + 6*ridx[i];
13777b2bb3b9SHong Zhang       y = yarray + 6*ridx[i];
137826e093fcSHong Zhang     }
137915091d37SBarry Smith     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4]; sum6 = y[5];
1380444d8c10SJed Brown     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);     /* Indices for the next row (assumes same size as this one) */
1381444d8c10SJed Brown     PetscPrefetchBlock(v+36*n,36*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
138215091d37SBarry Smith     for (j=0; j<n; j++) {
13833b95cb0eSSatish Balay       xb    = x + 6*(*idx++);
138415091d37SBarry Smith       x1    = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5];
138515091d37SBarry Smith       sum1 += v[0]*x1 + v[6]*x2  + v[12]*x3  + v[18]*x4 + v[24]*x5 + v[30]*x6;
138615091d37SBarry Smith       sum2 += v[1]*x1 + v[7]*x2  + v[13]*x3  + v[19]*x4 + v[25]*x5 + v[31]*x6;
138715091d37SBarry Smith       sum3 += v[2]*x1 + v[8]*x2  + v[14]*x3  + v[20]*x4 + v[26]*x5 + v[32]*x6;
138815091d37SBarry Smith       sum4 += v[3]*x1 + v[9]*x2  + v[15]*x3  + v[21]*x4 + v[27]*x5 + v[33]*x6;
138915091d37SBarry Smith       sum5 += v[4]*x1 + v[10]*x2 + v[16]*x3  + v[22]*x4 + v[28]*x5 + v[34]*x6;
139015091d37SBarry Smith       sum6 += v[5]*x1 + v[11]*x2 + v[17]*x3  + v[23]*x4 + v[29]*x5 + v[35]*x6;
139115091d37SBarry Smith       v    += 36;
139215091d37SBarry Smith     }
139315091d37SBarry Smith     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6;
139426e093fcSHong Zhang     if (!usecprow) {
139515091d37SBarry Smith       z += 6; y += 6;
139615091d37SBarry Smith     }
139726e093fcSHong Zhang   }
1398d9ca1df4SBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
1399d9ca1df4SBarry Smith   ierr = VecRestoreArrayPair(yy,zz,&yarray,&zarray);CHKERRQ(ierr);
1400dc0b31edSSatish Balay   ierr = PetscLogFlops(72.0*a->nz);CHKERRQ(ierr);
140115091d37SBarry Smith   PetscFunctionReturn(0);
140215091d37SBarry Smith }
14032d61bbb3SSatish Balay 
1404dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_7(Mat A,Vec xx,Vec yy,Vec zz)
14052d61bbb3SSatish Balay {
14062d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
1407d9ca1df4SBarry Smith   PetscScalar       *y = 0,*z = 0,sum1,sum2,sum3,sum4,sum5,sum6,sum7;
1408d9ca1df4SBarry Smith   const PetscScalar *x,*xb;
140926e093fcSHong Zhang   PetscScalar       x1,x2,x3,x4,x5,x6,x7,*yarray,*zarray;
1410d9ca1df4SBarry Smith   const MatScalar   *v;
1411dfbe8321SBarry Smith   PetscErrorCode    ierr;
1412d9ca1df4SBarry Smith   PetscInt          mbs = a->mbs,i,j,n;
1413d9ca1df4SBarry Smith   const PetscInt    *idx,*ii,*ridx = NULL;
1414ace3abfcSBarry Smith   PetscBool         usecprow=a->compressedrow.use;
14152d61bbb3SSatish Balay 
14162d61bbb3SSatish Balay   PetscFunctionBegin;
1417d9ca1df4SBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
1418d9ca1df4SBarry Smith   ierr = VecGetArrayPair(yy,zz,&yarray,&zarray);CHKERRQ(ierr);
14192d61bbb3SSatish Balay 
14202d61bbb3SSatish Balay   idx = a->j;
14212d61bbb3SSatish Balay   v   = a->a;
142226e093fcSHong Zhang   if (usecprow) {
14234eb6d288SHong Zhang     if (zz != yy) {
14244eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,7*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
14254eb6d288SHong Zhang     }
142626e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
142726e093fcSHong Zhang     ii   = a->compressedrow.i;
14287b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
142926e093fcSHong Zhang   } else {
14302d61bbb3SSatish Balay     ii = a->i;
143126e093fcSHong Zhang     y  = yarray;
143226e093fcSHong Zhang     z  = zarray;
143326e093fcSHong Zhang   }
14342d61bbb3SSatish Balay 
14352d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
14362d61bbb3SSatish Balay     n = ii[1] - ii[0]; ii++;
143726e093fcSHong Zhang     if (usecprow) {
14387b2bb3b9SHong Zhang       z = zarray + 7*ridx[i];
14397b2bb3b9SHong Zhang       y = yarray + 7*ridx[i];
144026e093fcSHong Zhang     }
14412d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4]; sum6 = y[5]; sum7 = y[6];
1442444d8c10SJed Brown     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);     /* Indices for the next row (assumes same size as this one) */
1443444d8c10SJed Brown     PetscPrefetchBlock(v+49*n,49*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
14442d61bbb3SSatish Balay     for (j=0; j<n; j++) {
14452d61bbb3SSatish Balay       xb    = x + 7*(*idx++);
14462d61bbb3SSatish Balay       x1    = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6];
14472d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[7]*x2  + v[14]*x3  + v[21]*x4 + v[28]*x5 + v[35]*x6 + v[42]*x7;
14482d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[8]*x2  + v[15]*x3  + v[22]*x4 + v[29]*x5 + v[36]*x6 + v[43]*x7;
14492d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[9]*x2  + v[16]*x3  + v[23]*x4 + v[30]*x5 + v[37]*x6 + v[44]*x7;
14502d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[10]*x2 + v[17]*x3  + v[24]*x4 + v[31]*x5 + v[38]*x6 + v[45]*x7;
14512d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[11]*x2 + v[18]*x3  + v[25]*x4 + v[32]*x5 + v[39]*x6 + v[46]*x7;
14522d61bbb3SSatish Balay       sum6 += v[5]*x1 + v[12]*x2 + v[19]*x3  + v[26]*x4 + v[33]*x5 + v[40]*x6 + v[47]*x7;
14532d61bbb3SSatish Balay       sum7 += v[6]*x1 + v[13]*x2 + v[20]*x3  + v[27]*x4 + v[34]*x5 + v[41]*x6 + v[48]*x7;
14542d61bbb3SSatish Balay       v    += 49;
14552d61bbb3SSatish Balay     }
14562d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
145726e093fcSHong Zhang     if (!usecprow) {
14582d61bbb3SSatish Balay       z += 7; y += 7;
14592d61bbb3SSatish Balay     }
146026e093fcSHong Zhang   }
1461d9ca1df4SBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
1462d9ca1df4SBarry Smith   ierr = VecRestoreArrayPair(yy,zz,&yarray,&zarray);CHKERRQ(ierr);
1463dc0b31edSSatish Balay   ierr = PetscLogFlops(98.0*a->nz);CHKERRQ(ierr);
14642d61bbb3SSatish Balay   PetscFunctionReturn(0);
14652d61bbb3SSatish Balay }
1466218c64b6SSatish Balay 
1467dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_N(Mat A,Vec xx,Vec yy,Vec zz)
14682d61bbb3SSatish Balay {
14692d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
1470d9ca1df4SBarry Smith   PetscScalar       *z = 0,*work,*workt,*zarray;
1471d9ca1df4SBarry Smith   const PetscScalar *x,*xb;
1472d9ca1df4SBarry Smith   const MatScalar   *v;
14736849ba73SBarry Smith   PetscErrorCode    ierr;
1474d9ca1df4SBarry Smith   PetscInt          mbs,i,bs=A->rmap->bs,j,n,bs2=a->bs2;
1475d9ca1df4SBarry Smith   PetscInt          ncols,k;
1476d9ca1df4SBarry Smith   const PetscInt    *ridx = NULL,*idx,*ii;
1477ace3abfcSBarry Smith   PetscBool         usecprow = a->compressedrow.use;
1478218c64b6SSatish Balay 
14792d61bbb3SSatish Balay   PetscFunctionBegin;
1480343551c4SBarry Smith   ierr = VecCopy(yy,zz);CHKERRQ(ierr);
1481d9ca1df4SBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
148226e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
14832d61bbb3SSatish Balay 
14842d61bbb3SSatish Balay   idx = a->j;
14852d61bbb3SSatish Balay   v   = a->a;
148626e093fcSHong Zhang   if (usecprow) {
148726e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
148826e093fcSHong Zhang     ii   = a->compressedrow.i;
14897b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
149026e093fcSHong Zhang   } else {
149126e093fcSHong Zhang     mbs = a->mbs;
14922d61bbb3SSatish Balay     ii  = a->i;
149326e093fcSHong Zhang     z   = zarray;
149426e093fcSHong Zhang   }
14952d61bbb3SSatish Balay 
14962d61bbb3SSatish Balay   if (!a->mult_work) {
1497d0f46423SBarry Smith     k    = PetscMax(A->rmap->n,A->cmap->n);
1498854ce69bSBarry Smith     ierr = PetscMalloc1(k+1,&a->mult_work);CHKERRQ(ierr);
14992d61bbb3SSatish Balay   }
15002d61bbb3SSatish Balay   work = a->mult_work;
15012d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
15022d61bbb3SSatish Balay     n     = ii[1] - ii[0]; ii++;
15032d61bbb3SSatish Balay     ncols = n*bs;
15042d61bbb3SSatish Balay     workt = work;
15052d61bbb3SSatish Balay     for (j=0; j<n; j++) {
15062d61bbb3SSatish Balay       xb = x + bs*(*idx++);
15072d61bbb3SSatish Balay       for (k=0; k<bs; k++) workt[k] = xb[k];
15082d61bbb3SSatish Balay       workt += bs;
15092d61bbb3SSatish Balay     }
15107b2bb3b9SHong Zhang     if (usecprow) z = zarray + bs*ridx[i];
151196b95a6bSBarry Smith     PetscKernel_w_gets_w_plus_Ar_times_v(bs,ncols,work,v,z);
151271044d3cSBarry Smith     /* BLASgemv_("N",&bs,&ncols,&_DOne,v,&bs,work,&_One,&_DOne,z,&_One); */
15132d61bbb3SSatish Balay     v += n*bs2;
151426fbe8dcSKarl Rupp     if (!usecprow) z += bs;
151526e093fcSHong Zhang   }
1516d9ca1df4SBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
151726e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
1518dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*a->nz*bs2);CHKERRQ(ierr);
15192d61bbb3SSatish Balay   PetscFunctionReturn(0);
15202d61bbb3SSatish Balay }
15212d61bbb3SSatish Balay 
1522547795f9SHong Zhang PetscErrorCode MatMultHermitianTranspose_SeqBAIJ(Mat A,Vec xx,Vec zz)
1523547795f9SHong Zhang {
1524547795f9SHong Zhang   PetscScalar    zero = 0.0;
1525547795f9SHong Zhang   PetscErrorCode ierr;
1526547795f9SHong Zhang 
1527547795f9SHong Zhang   PetscFunctionBegin;
1528547795f9SHong Zhang   ierr = VecSet(zz,zero);CHKERRQ(ierr);
1529547795f9SHong Zhang   ierr = MatMultHermitianTransposeAdd_SeqBAIJ(A,xx,zz,zz);CHKERRQ(ierr);
1530547795f9SHong Zhang   PetscFunctionReturn(0);
1531547795f9SHong Zhang }
1532547795f9SHong Zhang 
1533dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_SeqBAIJ(Mat A,Vec xx,Vec zz)
15342d61bbb3SSatish Balay {
15353447b6efSHong Zhang   PetscScalar    zero = 0.0;
15366849ba73SBarry Smith   PetscErrorCode ierr;
15372d61bbb3SSatish Balay 
15382d61bbb3SSatish Balay   PetscFunctionBegin;
15392dcb1b2aSMatthew Knepley   ierr = VecSet(zz,zero);CHKERRQ(ierr);
15403447b6efSHong Zhang   ierr = MatMultTransposeAdd_SeqBAIJ(A,xx,zz,zz);CHKERRQ(ierr);
15412d61bbb3SSatish Balay   PetscFunctionReturn(0);
15422d61bbb3SSatish Balay }
15432d61bbb3SSatish Balay 
1544547795f9SHong Zhang PetscErrorCode MatMultHermitianTransposeAdd_SeqBAIJ(Mat A,Vec xx,Vec yy,Vec zz)
1545547795f9SHong Zhang {
1546547795f9SHong Zhang   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
1547b8c08b77SHong Zhang   PetscScalar       *z,x1,x2,x3,x4,x5;
1548d9ca1df4SBarry Smith   const PetscScalar *x,*xb = NULL;
1549d9ca1df4SBarry Smith   const MatScalar   *v;
1550547795f9SHong Zhang   PetscErrorCode    ierr;
1551b8c08b77SHong Zhang   PetscInt          mbs,i,rval,bs=A->rmap->bs,j,n;
1552d9ca1df4SBarry Smith   const PetscInt    *idx,*ii,*ib,*ridx = NULL;
1553547795f9SHong Zhang   Mat_CompressedRow cprow = a->compressedrow;
1554ace3abfcSBarry Smith   PetscBool         usecprow = cprow.use;
1555547795f9SHong Zhang 
1556547795f9SHong Zhang   PetscFunctionBegin;
1557343551c4SBarry Smith   if (yy != zz) { ierr = VecCopy(yy,zz);CHKERRQ(ierr); }
1558d9ca1df4SBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
1559547795f9SHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
1560547795f9SHong Zhang 
1561547795f9SHong Zhang   idx = a->j;
1562547795f9SHong Zhang   v   = a->a;
1563547795f9SHong Zhang   if (usecprow) {
1564547795f9SHong Zhang     mbs  = cprow.nrows;
1565547795f9SHong Zhang     ii   = cprow.i;
1566547795f9SHong Zhang     ridx = cprow.rindex;
1567547795f9SHong Zhang   } else {
1568547795f9SHong Zhang     mbs=a->mbs;
1569547795f9SHong Zhang     ii = a->i;
1570547795f9SHong Zhang     xb = x;
1571547795f9SHong Zhang   }
1572547795f9SHong Zhang 
1573547795f9SHong Zhang   switch (bs) {
1574547795f9SHong Zhang   case 1:
1575547795f9SHong Zhang     for (i=0; i<mbs; i++) {
1576547795f9SHong Zhang       if (usecprow) xb = x + ridx[i];
1577547795f9SHong Zhang       x1 = xb[0];
1578547795f9SHong Zhang       ib = idx + ii[0];
1579547795f9SHong Zhang       n  = ii[1] - ii[0]; ii++;
1580547795f9SHong Zhang       for (j=0; j<n; j++) {
1581547795f9SHong Zhang         rval     = ib[j];
1582547795f9SHong Zhang         z[rval] += PetscConj(*v) * x1;
1583547795f9SHong Zhang         v++;
1584547795f9SHong Zhang       }
1585547795f9SHong Zhang       if (!usecprow) xb++;
1586547795f9SHong Zhang     }
1587547795f9SHong Zhang     break;
1588547795f9SHong Zhang   case 2:
1589547795f9SHong Zhang     for (i=0; i<mbs; i++) {
1590547795f9SHong Zhang       if (usecprow) xb = x + 2*ridx[i];
1591547795f9SHong Zhang       x1 = xb[0]; x2 = xb[1];
1592547795f9SHong Zhang       ib = idx + ii[0];
1593547795f9SHong Zhang       n  = ii[1] - ii[0]; ii++;
1594547795f9SHong Zhang       for (j=0; j<n; j++) {
1595547795f9SHong Zhang         rval       = ib[j]*2;
1596547795f9SHong Zhang         z[rval++] += PetscConj(v[0])*x1 + PetscConj(v[1])*x2;
1597547795f9SHong Zhang         z[rval++] += PetscConj(v[2])*x1 + PetscConj(v[3])*x2;
1598547795f9SHong Zhang         v         += 4;
1599547795f9SHong Zhang       }
1600547795f9SHong Zhang       if (!usecprow) xb += 2;
1601547795f9SHong Zhang     }
1602547795f9SHong Zhang     break;
1603547795f9SHong Zhang   case 3:
1604547795f9SHong Zhang     for (i=0; i<mbs; i++) {
1605547795f9SHong Zhang       if (usecprow) xb = x + 3*ridx[i];
1606547795f9SHong Zhang       x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
1607547795f9SHong Zhang       ib = idx + ii[0];
1608547795f9SHong Zhang       n  = ii[1] - ii[0]; ii++;
1609547795f9SHong Zhang       for (j=0; j<n; j++) {
1610547795f9SHong Zhang         rval       = ib[j]*3;
1611547795f9SHong Zhang         z[rval++] += PetscConj(v[0])*x1 + PetscConj(v[1])*x2 + PetscConj(v[2])*x3;
1612547795f9SHong Zhang         z[rval++] += PetscConj(v[3])*x1 + PetscConj(v[4])*x2 + PetscConj(v[5])*x3;
1613547795f9SHong Zhang         z[rval++] += PetscConj(v[6])*x1 + PetscConj(v[7])*x2 + PetscConj(v[8])*x3;
1614547795f9SHong Zhang         v         += 9;
1615547795f9SHong Zhang       }
1616547795f9SHong Zhang       if (!usecprow) xb += 3;
1617547795f9SHong Zhang     }
1618547795f9SHong Zhang     break;
1619547795f9SHong Zhang   case 4:
1620547795f9SHong Zhang     for (i=0; i<mbs; i++) {
1621547795f9SHong Zhang       if (usecprow) xb = x + 4*ridx[i];
1622547795f9SHong Zhang       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
1623547795f9SHong Zhang       ib = idx + ii[0];
1624547795f9SHong Zhang       n  = ii[1] - ii[0]; ii++;
1625547795f9SHong Zhang       for (j=0; j<n; j++) {
1626547795f9SHong Zhang         rval       = ib[j]*4;
1627547795f9SHong Zhang         z[rval++] +=  PetscConj(v[0])*x1 + PetscConj(v[1])*x2  + PetscConj(v[2])*x3  + PetscConj(v[3])*x4;
1628547795f9SHong Zhang         z[rval++] +=  PetscConj(v[4])*x1 + PetscConj(v[5])*x2  + PetscConj(v[6])*x3  + PetscConj(v[7])*x4;
1629547795f9SHong Zhang         z[rval++] +=  PetscConj(v[8])*x1 + PetscConj(v[9])*x2  + PetscConj(v[10])*x3 + PetscConj(v[11])*x4;
1630547795f9SHong Zhang         z[rval++] += PetscConj(v[12])*x1 + PetscConj(v[13])*x2 + PetscConj(v[14])*x3 + PetscConj(v[15])*x4;
1631547795f9SHong Zhang         v         += 16;
1632547795f9SHong Zhang       }
1633547795f9SHong Zhang       if (!usecprow) xb += 4;
1634547795f9SHong Zhang     }
1635547795f9SHong Zhang     break;
1636547795f9SHong Zhang   case 5:
1637547795f9SHong Zhang     for (i=0; i<mbs; i++) {
1638547795f9SHong Zhang       if (usecprow) xb = x + 5*ridx[i];
1639547795f9SHong Zhang       x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
1640547795f9SHong Zhang       x4 = xb[3]; x5 = xb[4];
1641547795f9SHong Zhang       ib = idx + ii[0];
1642547795f9SHong Zhang       n  = ii[1] - ii[0]; ii++;
1643547795f9SHong Zhang       for (j=0; j<n; j++) {
1644547795f9SHong Zhang         rval       = ib[j]*5;
1645547795f9SHong Zhang         z[rval++] +=  PetscConj(v[0])*x1 +  PetscConj(v[1])*x2 +  PetscConj(v[2])*x3 +  PetscConj(v[3])*x4 +  PetscConj(v[4])*x5;
1646547795f9SHong Zhang         z[rval++] +=  PetscConj(v[5])*x1 +  PetscConj(v[6])*x2 +  PetscConj(v[7])*x3 +  PetscConj(v[8])*x4 +  PetscConj(v[9])*x5;
1647547795f9SHong Zhang         z[rval++] += PetscConj(v[10])*x1 + PetscConj(v[11])*x2 + PetscConj(v[12])*x3 + PetscConj(v[13])*x4 + PetscConj(v[14])*x5;
1648547795f9SHong Zhang         z[rval++] += PetscConj(v[15])*x1 + PetscConj(v[16])*x2 + PetscConj(v[17])*x3 + PetscConj(v[18])*x4 + PetscConj(v[19])*x5;
1649547795f9SHong Zhang         z[rval++] += PetscConj(v[20])*x1 + PetscConj(v[21])*x2 + PetscConj(v[22])*x3 + PetscConj(v[23])*x4 + PetscConj(v[24])*x5;
1650547795f9SHong Zhang         v         += 25;
1651547795f9SHong Zhang       }
1652547795f9SHong Zhang       if (!usecprow) xb += 5;
1653547795f9SHong Zhang     }
1654547795f9SHong Zhang     break;
1655968ae2c8SSatish Balay   default: /* block sizes larger than 5 by 5 are handled by BLAS */
1656968ae2c8SSatish Balay     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"block size larger than 5 is not supported yet");
1657968ae2c8SSatish Balay #if 0
1658968ae2c8SSatish Balay     {
1659b8c08b77SHong Zhang       PetscInt          ncols,k,bs2=a->bs2;
1660b8c08b77SHong Zhang       PetscScalar       *work,*workt,zb;
1661d9ca1df4SBarry Smith       const PetscScalar *xtmp;
1662547795f9SHong Zhang       if (!a->mult_work) {
1663547795f9SHong Zhang         k    = PetscMax(A->rmap->n,A->cmap->n);
1664854ce69bSBarry Smith         ierr = PetscMalloc1(k+1,&a->mult_work);CHKERRQ(ierr);
1665547795f9SHong Zhang       }
1666547795f9SHong Zhang       work = a->mult_work;
1667547795f9SHong Zhang       xtmp = x;
1668547795f9SHong Zhang       for (i=0; i<mbs; i++) {
1669547795f9SHong Zhang         n     = ii[1] - ii[0]; ii++;
1670547795f9SHong Zhang         ncols = n*bs;
1671547795f9SHong Zhang         ierr  = PetscMemzero(work,ncols*sizeof(PetscScalar));CHKERRQ(ierr);
167226fbe8dcSKarl Rupp         if (usecprow) xtmp = x + bs*ridx[i];
167396b95a6bSBarry Smith         PetscKernel_w_gets_w_plus_trans_Ar_times_v(bs,ncols,xtmp,v,work);
1674547795f9SHong Zhang         /* BLASgemv_("T",&bs,&ncols,&_DOne,v,&bs,xtmp,&_One,&_DOne,work,&_One); */
1675547795f9SHong Zhang         v += n*bs2;
1676547795f9SHong Zhang         if (!usecprow) xtmp += bs;
1677547795f9SHong Zhang         workt = work;
1678547795f9SHong Zhang         for (j=0; j<n; j++) {
1679547795f9SHong Zhang           zb = z + bs*(*idx++);
1680547795f9SHong Zhang           for (k=0; k<bs; k++) zb[k] += workt[k] ;
1681547795f9SHong Zhang           workt += bs;
1682547795f9SHong Zhang         }
1683547795f9SHong Zhang       }
1684547795f9SHong Zhang     }
1685968ae2c8SSatish Balay #endif
1686547795f9SHong Zhang   }
1687d9ca1df4SBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
1688547795f9SHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
1689547795f9SHong Zhang   ierr = PetscLogFlops(2.0*a->nz*a->bs2);CHKERRQ(ierr);
1690547795f9SHong Zhang   PetscFunctionReturn(0);
1691547795f9SHong Zhang }
1692547795f9SHong Zhang 
1693dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_SeqBAIJ(Mat A,Vec xx,Vec yy,Vec zz)
16942d61bbb3SSatish Balay {
16952d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
1696d9ca1df4SBarry Smith   PetscScalar       *zb,*z,x1,x2,x3,x4,x5;
1697d9ca1df4SBarry Smith   const PetscScalar *x,*xb = 0;
1698d9ca1df4SBarry Smith   const MatScalar   *v;
16996849ba73SBarry Smith   PetscErrorCode    ierr;
1700d9ca1df4SBarry Smith   PetscInt          mbs,i,rval,bs=A->rmap->bs,j,n,bs2=a->bs2;
1701d9ca1df4SBarry Smith   const PetscInt    *idx,*ii,*ib,*ridx = NULL;
17023447b6efSHong Zhang   Mat_CompressedRow cprow   = a->compressedrow;
1703ace3abfcSBarry Smith   PetscBool         usecprow=cprow.use;
17042d61bbb3SSatish Balay 
17052d61bbb3SSatish Balay   PetscFunctionBegin;
1706343551c4SBarry Smith   if (yy != zz) { ierr = VecCopy(yy,zz);CHKERRQ(ierr); }
1707d9ca1df4SBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
17083447b6efSHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
17092d61bbb3SSatish Balay 
17102d61bbb3SSatish Balay   idx = a->j;
17112d61bbb3SSatish Balay   v   = a->a;
17123447b6efSHong Zhang   if (usecprow) {
17133447b6efSHong Zhang     mbs  = cprow.nrows;
17143447b6efSHong Zhang     ii   = cprow.i;
17157b2bb3b9SHong Zhang     ridx = cprow.rindex;
17163447b6efSHong Zhang   } else {
17173447b6efSHong Zhang     mbs=a->mbs;
17182d61bbb3SSatish Balay     ii = a->i;
1719f1af5d2fSBarry Smith     xb = x;
17203447b6efSHong Zhang   }
17212d61bbb3SSatish Balay 
17222d61bbb3SSatish Balay   switch (bs) {
17232d61bbb3SSatish Balay   case 1:
17242d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
17257b2bb3b9SHong Zhang       if (usecprow) xb = x + ridx[i];
1726f1af5d2fSBarry Smith       x1 = xb[0];
17273447b6efSHong Zhang       ib = idx + ii[0];
17283447b6efSHong Zhang       n  = ii[1] - ii[0]; ii++;
17292d61bbb3SSatish Balay       for (j=0; j<n; j++) {
17302d61bbb3SSatish Balay         rval     = ib[j];
1731f1af5d2fSBarry Smith         z[rval] += *v * x1;
1732f1af5d2fSBarry Smith         v++;
17332d61bbb3SSatish Balay       }
17343447b6efSHong Zhang       if (!usecprow) xb++;
17352d61bbb3SSatish Balay     }
17362d61bbb3SSatish Balay     break;
17372d61bbb3SSatish Balay   case 2:
17382d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
17397b2bb3b9SHong Zhang       if (usecprow) xb = x + 2*ridx[i];
1740f1af5d2fSBarry Smith       x1 = xb[0]; x2 = xb[1];
17413447b6efSHong Zhang       ib = idx + ii[0];
17423447b6efSHong Zhang       n  = ii[1] - ii[0]; ii++;
17432d61bbb3SSatish Balay       for (j=0; j<n; j++) {
17442d61bbb3SSatish Balay         rval       = ib[j]*2;
17452d61bbb3SSatish Balay         z[rval++] += v[0]*x1 + v[1]*x2;
17462d61bbb3SSatish Balay         z[rval++] += v[2]*x1 + v[3]*x2;
17472d61bbb3SSatish Balay         v         += 4;
17482d61bbb3SSatish Balay       }
17493447b6efSHong Zhang       if (!usecprow) xb += 2;
17502d61bbb3SSatish Balay     }
17512d61bbb3SSatish Balay     break;
17522d61bbb3SSatish Balay   case 3:
17532d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
17547b2bb3b9SHong Zhang       if (usecprow) xb = x + 3*ridx[i];
1755f1af5d2fSBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
17563447b6efSHong Zhang       ib = idx + ii[0];
17573447b6efSHong Zhang       n  = ii[1] - ii[0]; ii++;
17582d61bbb3SSatish Balay       for (j=0; j<n; j++) {
17592d61bbb3SSatish Balay         rval       = ib[j]*3;
17602d61bbb3SSatish Balay         z[rval++] += v[0]*x1 + v[1]*x2 + v[2]*x3;
17612d61bbb3SSatish Balay         z[rval++] += v[3]*x1 + v[4]*x2 + v[5]*x3;
17622d61bbb3SSatish Balay         z[rval++] += v[6]*x1 + v[7]*x2 + v[8]*x3;
17632d61bbb3SSatish Balay         v         += 9;
17642d61bbb3SSatish Balay       }
17653447b6efSHong Zhang       if (!usecprow) xb += 3;
17662d61bbb3SSatish Balay     }
17672d61bbb3SSatish Balay     break;
17682d61bbb3SSatish Balay   case 4:
17692d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
17707b2bb3b9SHong Zhang       if (usecprow) xb = x + 4*ridx[i];
1771f1af5d2fSBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
17723447b6efSHong Zhang       ib = idx + ii[0];
17733447b6efSHong Zhang       n  = ii[1] - ii[0]; ii++;
17742d61bbb3SSatish Balay       for (j=0; j<n; j++) {
17752d61bbb3SSatish Balay         rval       = ib[j]*4;
17762d61bbb3SSatish Balay         z[rval++] +=  v[0]*x1 +  v[1]*x2 +  v[2]*x3 +  v[3]*x4;
17772d61bbb3SSatish Balay         z[rval++] +=  v[4]*x1 +  v[5]*x2 +  v[6]*x3 +  v[7]*x4;
17782d61bbb3SSatish Balay         z[rval++] +=  v[8]*x1 +  v[9]*x2 + v[10]*x3 + v[11]*x4;
17792d61bbb3SSatish Balay         z[rval++] += v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4;
17802d61bbb3SSatish Balay         v         += 16;
17812d61bbb3SSatish Balay       }
17823447b6efSHong Zhang       if (!usecprow) xb += 4;
17832d61bbb3SSatish Balay     }
17842d61bbb3SSatish Balay     break;
17852d61bbb3SSatish Balay   case 5:
17862d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
17877b2bb3b9SHong Zhang       if (usecprow) xb = x + 5*ridx[i];
1788f1af5d2fSBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
17892d61bbb3SSatish Balay       x4 = xb[3]; x5 = xb[4];
17903447b6efSHong Zhang       ib = idx + ii[0];
17913447b6efSHong Zhang       n  = ii[1] - ii[0]; ii++;
17922d61bbb3SSatish Balay       for (j=0; j<n; j++) {
17932d61bbb3SSatish Balay         rval       = ib[j]*5;
17942d61bbb3SSatish Balay         z[rval++] +=  v[0]*x1 +  v[1]*x2 +  v[2]*x3 +  v[3]*x4 +  v[4]*x5;
17952d61bbb3SSatish Balay         z[rval++] +=  v[5]*x1 +  v[6]*x2 +  v[7]*x3 +  v[8]*x4 +  v[9]*x5;
17962d61bbb3SSatish Balay         z[rval++] += v[10]*x1 + v[11]*x2 + v[12]*x3 + v[13]*x4 + v[14]*x5;
17972d61bbb3SSatish Balay         z[rval++] += v[15]*x1 + v[16]*x2 + v[17]*x3 + v[18]*x4 + v[19]*x5;
17982d61bbb3SSatish Balay         z[rval++] += v[20]*x1 + v[21]*x2 + v[22]*x3 + v[23]*x4 + v[24]*x5;
17992d61bbb3SSatish Balay         v         += 25;
18002d61bbb3SSatish Balay       }
18013447b6efSHong Zhang       if (!usecprow) xb += 5;
18022d61bbb3SSatish Balay     }
18032d61bbb3SSatish Balay     break;
1804f1af5d2fSBarry Smith   default: {      /* block sizes larger then 5 by 5 are handled by BLAS */
1805690b6cddSBarry Smith     PetscInt          ncols,k;
1806d9ca1df4SBarry Smith     PetscScalar       *work,*workt;
1807d9ca1df4SBarry Smith     const PetscScalar *xtmp;
18082d61bbb3SSatish Balay     if (!a->mult_work) {
1809d0f46423SBarry Smith       k    = PetscMax(A->rmap->n,A->cmap->n);
1810854ce69bSBarry Smith       ierr = PetscMalloc1(k+1,&a->mult_work);CHKERRQ(ierr);
18112d61bbb3SSatish Balay     }
18122d61bbb3SSatish Balay     work = a->mult_work;
18133447b6efSHong Zhang     xtmp = x;
18142d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
18152d61bbb3SSatish Balay       n     = ii[1] - ii[0]; ii++;
18162d61bbb3SSatish Balay       ncols = n*bs;
181787828ca2SBarry Smith       ierr  = PetscMemzero(work,ncols*sizeof(PetscScalar));CHKERRQ(ierr);
181826fbe8dcSKarl Rupp       if (usecprow) xtmp = x + bs*ridx[i];
181996b95a6bSBarry Smith       PetscKernel_w_gets_w_plus_trans_Ar_times_v(bs,ncols,xtmp,v,work);
182071044d3cSBarry Smith       /* BLASgemv_("T",&bs,&ncols,&_DOne,v,&bs,xtmp,&_One,&_DOne,work,&_One); */
18212d61bbb3SSatish Balay       v += n*bs2;
18223447b6efSHong Zhang       if (!usecprow) xtmp += bs;
18232d61bbb3SSatish Balay       workt = work;
18242d61bbb3SSatish Balay       for (j=0; j<n; j++) {
18252d61bbb3SSatish Balay         zb = z + bs*(*idx++);
18262d61bbb3SSatish Balay         for (k=0; k<bs; k++) zb[k] += workt[k];
18272d61bbb3SSatish Balay         workt += bs;
18282d61bbb3SSatish Balay       }
18292d61bbb3SSatish Balay     }
18302d61bbb3SSatish Balay     }
18312d61bbb3SSatish Balay   }
1832d9ca1df4SBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
18333447b6efSHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
1834dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*a->nz*a->bs2);CHKERRQ(ierr);
18352d61bbb3SSatish Balay   PetscFunctionReturn(0);
18362d61bbb3SSatish Balay }
18372d61bbb3SSatish Balay 
1838f4df32b1SMatthew Knepley PetscErrorCode MatScale_SeqBAIJ(Mat inA,PetscScalar alpha)
18392d61bbb3SSatish Balay {
18402d61bbb3SSatish Balay   Mat_SeqBAIJ    *a      = (Mat_SeqBAIJ*)inA->data;
1841690b6cddSBarry Smith   PetscInt       totalnz = a->bs2*a->nz;
1842f4df32b1SMatthew Knepley   PetscScalar    oalpha  = alpha;
1843efee365bSSatish Balay   PetscErrorCode ierr;
1844c5df96a5SBarry Smith   PetscBLASInt   one = 1,tnz;
18452d61bbb3SSatish Balay 
18462d61bbb3SSatish Balay   PetscFunctionBegin;
1847c5df96a5SBarry Smith   ierr = PetscBLASIntCast(totalnz,&tnz);CHKERRQ(ierr);
18488b83055fSJed Brown   PetscStackCallBLAS("BLASscal",BLASscal_(&tnz,&oalpha,a->a,&one));
1849efee365bSSatish Balay   ierr = PetscLogFlops(totalnz);CHKERRQ(ierr);
18502d61bbb3SSatish Balay   PetscFunctionReturn(0);
18512d61bbb3SSatish Balay }
18522d61bbb3SSatish Balay 
1853dfbe8321SBarry Smith PetscErrorCode MatNorm_SeqBAIJ(Mat A,NormType type,PetscReal *norm)
18542d61bbb3SSatish Balay {
18558a62d963SHong Zhang   PetscErrorCode ierr;
18562d61bbb3SSatish Balay   Mat_SeqBAIJ    *a  = (Mat_SeqBAIJ*)A->data;
18573f1db9ecSBarry Smith   MatScalar      *v  = a->a;
1858329f5518SBarry Smith   PetscReal      sum = 0.0;
1859d0f46423SBarry Smith   PetscInt       i,j,k,bs=A->rmap->bs,nz=a->nz,bs2=a->bs2,k1;
18602d61bbb3SSatish Balay 
18612d61bbb3SSatish Balay   PetscFunctionBegin;
18622d61bbb3SSatish Balay   if (type == NORM_FROBENIUS) {
1863570b7f6dSBarry Smith #if defined(PETSC_USE_REAL___FP16)
1864570b7f6dSBarry Smith     PetscBLASInt one = 1,cnt = bs2*nz;
1865570b7f6dSBarry Smith     *norm = BLASnrm2_(&cnt,v,&one);
1866570b7f6dSBarry Smith #else
18672d61bbb3SSatish Balay     for (i=0; i<bs2*nz; i++) {
1868329f5518SBarry Smith       sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
18692d61bbb3SSatish Balay     }
1870570b7f6dSBarry Smith #endif
18718f1a2a5eSBarry Smith     *norm = PetscSqrtReal(sum);
187251f70360SJed Brown     ierr = PetscLogFlops(2*bs2*nz);CHKERRQ(ierr);
18738a62d963SHong Zhang   } else if (type == NORM_1) { /* maximum column sum */
18748a62d963SHong Zhang     PetscReal *tmp;
18758a62d963SHong Zhang     PetscInt  *bcol = a->j;
1876854ce69bSBarry Smith     ierr = PetscCalloc1(A->cmap->n+1,&tmp);CHKERRQ(ierr);
18778a62d963SHong Zhang     for (i=0; i<nz; i++) {
18788a62d963SHong Zhang       for (j=0; j<bs; j++) {
18798a62d963SHong Zhang         k1 = bs*(*bcol) + j; /* column index */
18808a62d963SHong Zhang         for (k=0; k<bs; k++) {
18818a62d963SHong Zhang           tmp[k1] += PetscAbsScalar(*v); v++;
18828a62d963SHong Zhang         }
18838a62d963SHong Zhang       }
18848a62d963SHong Zhang       bcol++;
18858a62d963SHong Zhang     }
18868a62d963SHong Zhang     *norm = 0.0;
1887d0f46423SBarry Smith     for (j=0; j<A->cmap->n; j++) {
18888a62d963SHong Zhang       if (tmp[j] > *norm) *norm = tmp[j];
18898a62d963SHong Zhang     }
18908a62d963SHong Zhang     ierr = PetscFree(tmp);CHKERRQ(ierr);
189151f70360SJed Brown     ierr = PetscLogFlops(PetscMax(bs2*nz-1,0));CHKERRQ(ierr);
1892596552b5SBarry Smith   } else if (type == NORM_INFINITY) { /* maximum row sum */
1893596552b5SBarry Smith     *norm = 0.0;
1894596552b5SBarry Smith     for (k=0; k<bs; k++) {
189574f84c7bSSatish Balay       for (j=0; j<a->mbs; j++) {
1896596552b5SBarry Smith         v   = a->a + bs2*a->i[j] + k;
1897596552b5SBarry Smith         sum = 0.0;
1898596552b5SBarry Smith         for (i=0; i<a->i[j+1]-a->i[j]; i++) {
18990e90e235SBarry Smith           for (k1=0; k1<bs; k1++) {
1900596552b5SBarry Smith             sum += PetscAbsScalar(*v);
1901596552b5SBarry Smith             v   += bs;
19022d61bbb3SSatish Balay           }
19030e90e235SBarry Smith         }
1904596552b5SBarry Smith         if (sum > *norm) *norm = sum;
1905596552b5SBarry Smith       }
1906596552b5SBarry Smith     }
190751f70360SJed Brown     ierr = PetscLogFlops(PetscMax(bs2*nz-1,0));CHKERRQ(ierr);
1908e7e72b3dSBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for this norm yet");
19092d61bbb3SSatish Balay   PetscFunctionReturn(0);
19102d61bbb3SSatish Balay }
19112d61bbb3SSatish Balay 
19122d61bbb3SSatish Balay 
1913ace3abfcSBarry Smith PetscErrorCode MatEqual_SeqBAIJ(Mat A,Mat B,PetscBool * flg)
19142d61bbb3SSatish Balay {
19152d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)B->data;
1916dfbe8321SBarry Smith   PetscErrorCode ierr;
19172d61bbb3SSatish Balay 
19182d61bbb3SSatish Balay   PetscFunctionBegin;
19192d61bbb3SSatish Balay   /* If the  matrix/block dimensions are not equal, or no of nonzeros or shift */
1920d0f46423SBarry Smith   if ((A->rmap->N != B->rmap->N) || (A->cmap->n != B->cmap->n) || (A->rmap->bs != B->rmap->bs)|| (a->nz != b->nz)) {
1921273d9f13SBarry Smith     *flg = PETSC_FALSE;
1922273d9f13SBarry Smith     PetscFunctionReturn(0);
19232d61bbb3SSatish Balay   }
19242d61bbb3SSatish Balay 
19252d61bbb3SSatish Balay   /* if the a->i are the same */
1926690b6cddSBarry Smith   ierr = PetscMemcmp(a->i,b->i,(a->mbs+1)*sizeof(PetscInt),flg);CHKERRQ(ierr);
192726fbe8dcSKarl Rupp   if (!*flg) PetscFunctionReturn(0);
19282d61bbb3SSatish Balay 
19292d61bbb3SSatish Balay   /* if a->j are the same */
1930690b6cddSBarry Smith   ierr = PetscMemcmp(a->j,b->j,(a->nz)*sizeof(PetscInt),flg);CHKERRQ(ierr);
193126fbe8dcSKarl Rupp   if (!*flg) PetscFunctionReturn(0);
193226fbe8dcSKarl Rupp 
19332d61bbb3SSatish Balay   /* if a->a are the same */
1934d0f46423SBarry Smith   ierr = PetscMemcmp(a->a,b->a,(a->nz)*(A->rmap->bs)*(B->rmap->bs)*sizeof(PetscScalar),flg);CHKERRQ(ierr);
19352d61bbb3SSatish Balay   PetscFunctionReturn(0);
19362d61bbb3SSatish Balay 
19372d61bbb3SSatish Balay }
19382d61bbb3SSatish Balay 
1939dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_SeqBAIJ(Mat A,Vec v)
19402d61bbb3SSatish Balay {
19412d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
1942dfbe8321SBarry Smith   PetscErrorCode ierr;
1943690b6cddSBarry Smith   PetscInt       i,j,k,n,row,bs,*ai,*aj,ambs,bs2;
194487828ca2SBarry Smith   PetscScalar    *x,zero = 0.0;
19453f1db9ecSBarry Smith   MatScalar      *aa,*aa_j;
19462d61bbb3SSatish Balay 
19472d61bbb3SSatish Balay   PetscFunctionBegin;
1948e32f2f54SBarry Smith   if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1949d0f46423SBarry Smith   bs   = A->rmap->bs;
19502d61bbb3SSatish Balay   aa   = a->a;
19512d61bbb3SSatish Balay   ai   = a->i;
19522d61bbb3SSatish Balay   aj   = a->j;
19532d61bbb3SSatish Balay   ambs = a->mbs;
19542d61bbb3SSatish Balay   bs2  = a->bs2;
19552d61bbb3SSatish Balay 
19562dcb1b2aSMatthew Knepley   ierr = VecSet(v,zero);CHKERRQ(ierr);
19571ebc52fbSHong Zhang   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1958e1311b90SBarry Smith   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
1959e32f2f54SBarry Smith   if (n != A->rmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
19602d61bbb3SSatish Balay   for (i=0; i<ambs; i++) {
19612d61bbb3SSatish Balay     for (j=ai[i]; j<ai[i+1]; j++) {
19622d61bbb3SSatish Balay       if (aj[j] == i) {
19632d61bbb3SSatish Balay         row  = i*bs;
19642d61bbb3SSatish Balay         aa_j = aa+j*bs2;
19652d61bbb3SSatish Balay         for (k=0; k<bs2; k+=(bs+1),row++) x[row] = aa_j[k];
19662d61bbb3SSatish Balay         break;
19672d61bbb3SSatish Balay       }
19682d61bbb3SSatish Balay     }
19692d61bbb3SSatish Balay   }
19701ebc52fbSHong Zhang   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
19712d61bbb3SSatish Balay   PetscFunctionReturn(0);
19722d61bbb3SSatish Balay }
19732d61bbb3SSatish Balay 
1974dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_SeqBAIJ(Mat A,Vec ll,Vec rr)
19752d61bbb3SSatish Balay {
19762d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
197753ef36baSBarry Smith   const PetscScalar *l,*r,*li,*ri;
197853ef36baSBarry Smith   PetscScalar       x;
19793f1db9ecSBarry Smith   MatScalar         *aa, *v;
1980dfbe8321SBarry Smith   PetscErrorCode    ierr;
198153ef36baSBarry Smith   PetscInt          i,j,k,lm,rn,M,m,n,mbs,tmp,bs,bs2,iai;
198253ef36baSBarry Smith   const PetscInt    *ai,*aj;
19832d61bbb3SSatish Balay 
19842d61bbb3SSatish Balay   PetscFunctionBegin;
19852d61bbb3SSatish Balay   ai  = a->i;
19862d61bbb3SSatish Balay   aj  = a->j;
19872d61bbb3SSatish Balay   aa  = a->a;
1988d0f46423SBarry Smith   m   = A->rmap->n;
1989d0f46423SBarry Smith   n   = A->cmap->n;
1990d0f46423SBarry Smith   bs  = A->rmap->bs;
19912d61bbb3SSatish Balay   mbs = a->mbs;
19922d61bbb3SSatish Balay   bs2 = a->bs2;
19932d61bbb3SSatish Balay   if (ll) {
199453ef36baSBarry Smith     ierr = VecGetArrayRead(ll,&l);CHKERRQ(ierr);
1995e1311b90SBarry Smith     ierr = VecGetLocalSize(ll,&lm);CHKERRQ(ierr);
1996e32f2f54SBarry Smith     if (lm != m) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Left scaling vector wrong length");
19972d61bbb3SSatish Balay     for (i=0; i<mbs; i++) { /* for each block row */
19982d61bbb3SSatish Balay       M  = ai[i+1] - ai[i];
19992d61bbb3SSatish Balay       li = l + i*bs;
20002d61bbb3SSatish Balay       v  = aa + bs2*ai[i];
20012d61bbb3SSatish Balay       for (j=0; j<M; j++) { /* for each block */
20022d61bbb3SSatish Balay         for (k=0; k<bs2; k++) {
20032d61bbb3SSatish Balay           (*v++) *= li[k%bs];
20042d61bbb3SSatish Balay         }
20052d61bbb3SSatish Balay       }
20062d61bbb3SSatish Balay     }
200753ef36baSBarry Smith     ierr = VecRestoreArrayRead(ll,&l);CHKERRQ(ierr);
2008efee365bSSatish Balay     ierr = PetscLogFlops(a->nz);CHKERRQ(ierr);
20092d61bbb3SSatish Balay   }
20102d61bbb3SSatish Balay 
20112d61bbb3SSatish Balay   if (rr) {
201253ef36baSBarry Smith     ierr = VecGetArrayRead(rr,&r);CHKERRQ(ierr);
2013e1311b90SBarry Smith     ierr = VecGetLocalSize(rr,&rn);CHKERRQ(ierr);
2014e32f2f54SBarry Smith     if (rn != n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Right scaling vector wrong length");
20152d61bbb3SSatish Balay     for (i=0; i<mbs; i++) { /* for each block row */
201653ef36baSBarry Smith       iai = ai[i];
201753ef36baSBarry Smith       M   = ai[i+1] - iai;
201853ef36baSBarry Smith       v   = aa + bs2*iai;
20192d61bbb3SSatish Balay       for (j=0; j<M; j++) { /* for each block */
202053ef36baSBarry Smith         ri = r + bs*aj[iai+j];
20212d61bbb3SSatish Balay         for (k=0; k<bs; k++) {
20222d61bbb3SSatish Balay           x = ri[k];
202353ef36baSBarry Smith           for (tmp=0; tmp<bs; tmp++) v[tmp] *= x;
202453ef36baSBarry Smith           v += bs;
20252d61bbb3SSatish Balay         }
20262d61bbb3SSatish Balay       }
20272d61bbb3SSatish Balay     }
202853ef36baSBarry Smith     ierr = VecRestoreArrayRead(rr,&r);CHKERRQ(ierr);
2029efee365bSSatish Balay     ierr = PetscLogFlops(a->nz);CHKERRQ(ierr);
20302d61bbb3SSatish Balay   }
20312d61bbb3SSatish Balay   PetscFunctionReturn(0);
20322d61bbb3SSatish Balay }
20332d61bbb3SSatish Balay 
20342d61bbb3SSatish Balay 
2035dfbe8321SBarry Smith PetscErrorCode MatGetInfo_SeqBAIJ(Mat A,MatInfoType flag,MatInfo *info)
20362d61bbb3SSatish Balay {
20372d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
20382d61bbb3SSatish Balay 
20392d61bbb3SSatish Balay   PetscFunctionBegin;
20402d61bbb3SSatish Balay   info->block_size   = a->bs2;
2041ceed8ce5SJed Brown   info->nz_allocated = a->bs2*a->maxnz;
20422d61bbb3SSatish Balay   info->nz_used      = a->bs2*a->nz;
20432d61bbb3SSatish Balay   info->nz_unneeded  = (double)(info->nz_allocated - info->nz_used);
20442d61bbb3SSatish Balay   info->assemblies   = A->num_ass;
20458e58a170SBarry Smith   info->mallocs      = A->info.mallocs;
20467adad957SLisandro Dalcin   info->memory       = ((PetscObject)A)->mem;
2047d5f3da31SBarry Smith   if (A->factortype) {
20482d61bbb3SSatish Balay     info->fill_ratio_given  = A->info.fill_ratio_given;
20492d61bbb3SSatish Balay     info->fill_ratio_needed = A->info.fill_ratio_needed;
20502d61bbb3SSatish Balay     info->factor_mallocs    = A->info.factor_mallocs;
20512d61bbb3SSatish Balay   } else {
20522d61bbb3SSatish Balay     info->fill_ratio_given  = 0;
20532d61bbb3SSatish Balay     info->fill_ratio_needed = 0;
20542d61bbb3SSatish Balay     info->factor_mallocs    = 0;
20552d61bbb3SSatish Balay   }
20562d61bbb3SSatish Balay   PetscFunctionReturn(0);
20572d61bbb3SSatish Balay }
20582d61bbb3SSatish Balay 
2059dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_SeqBAIJ(Mat A)
20602d61bbb3SSatish Balay {
20612d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
2062dfbe8321SBarry Smith   PetscErrorCode ierr;
20632d61bbb3SSatish Balay 
20642d61bbb3SSatish Balay   PetscFunctionBegin;
2065549d3d68SSatish Balay   ierr = PetscMemzero(a->a,a->bs2*a->i[a->mbs]*sizeof(MatScalar));CHKERRQ(ierr);
20662d61bbb3SSatish Balay   PetscFunctionReturn(0);
20672d61bbb3SSatish Balay }
2068