xref: /petsc/src/mat/impls/baij/seq/baij2.c (revision 5ab264f3a2c5f7b4402bc8da2ae9fefa99c8eabd)
1 
2 #include <../src/mat/impls/baij/seq/baij.h>
3 #include <petsc/private/kernels/blockinvert.h>
4 #include <petscbt.h>
5 #include <petscblaslapack.h>
6 
7 #if defined(PETSC_HAVE_IMMINTRIN_H)
8 #include <immintrin.h>
9 #endif
10 
11 PetscErrorCode MatIncreaseOverlap_SeqBAIJ(Mat A,PetscInt is_max,IS is[],PetscInt ov)
12 {
13   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
14   PetscErrorCode ierr;
15   PetscInt       row,i,j,k,l,m,n,*nidx,isz,val,ival;
16   const PetscInt *idx;
17   PetscInt       start,end,*ai,*aj,bs,*nidx2;
18   PetscBT        table;
19 
20   PetscFunctionBegin;
21   m  = a->mbs;
22   ai = a->i;
23   aj = a->j;
24   bs = A->rmap->bs;
25 
26   if (ov < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap specified");
27 
28   ierr = PetscBTCreate(m,&table);CHKERRQ(ierr);
29   ierr = PetscMalloc1(m+1,&nidx);CHKERRQ(ierr);
30   ierr = PetscMalloc1(A->rmap->N+1,&nidx2);CHKERRQ(ierr);
31 
32   for (i=0; i<is_max; i++) {
33     /* Initialise the two local arrays */
34     isz  = 0;
35     ierr = PetscBTMemzero(m,table);CHKERRQ(ierr);
36 
37     /* Extract the indices, assume there can be duplicate entries */
38     ierr = ISGetIndices(is[i],&idx);CHKERRQ(ierr);
39     ierr = ISGetLocalSize(is[i],&n);CHKERRQ(ierr);
40 
41     /* Enter these into the temp arrays i.e mark table[row], enter row into new index */
42     for (j=0; j<n; ++j) {
43       ival = idx[j]/bs; /* convert the indices into block indices */
44       if (ival>=m) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"index greater than mat-dim");
45       if (!PetscBTLookupSet(table,ival)) nidx[isz++] = ival;
46     }
47     ierr = ISRestoreIndices(is[i],&idx);CHKERRQ(ierr);
48     ierr = ISDestroy(&is[i]);CHKERRQ(ierr);
49 
50     k = 0;
51     for (j=0; j<ov; j++) { /* for each overlap*/
52       n = isz;
53       for (; k<n; k++) {  /* do only those rows in nidx[k], which are not done yet */
54         row   = nidx[k];
55         start = ai[row];
56         end   = ai[row+1];
57         for (l = start; l<end; l++) {
58           val = aj[l];
59           if (!PetscBTLookupSet(table,val)) nidx[isz++] = val;
60         }
61       }
62     }
63     /* expand the Index Set */
64     for (j=0; j<isz; j++) {
65       for (k=0; k<bs; k++) nidx2[j*bs+k] = nidx[j]*bs+k;
66     }
67     ierr = ISCreateGeneral(PETSC_COMM_SELF,isz*bs,nidx2,PETSC_COPY_VALUES,is+i);CHKERRQ(ierr);
68   }
69   ierr = PetscBTDestroy(&table);CHKERRQ(ierr);
70   ierr = PetscFree(nidx);CHKERRQ(ierr);
71   ierr = PetscFree(nidx2);CHKERRQ(ierr);
72   PetscFunctionReturn(0);
73 }
74 
75 PetscErrorCode MatCreateSubMatrix_SeqBAIJ_Private(Mat A,IS isrow,IS iscol,MatReuse scall,Mat *B)
76 {
77   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*c;
78   PetscErrorCode ierr;
79   PetscInt       *smap,i,k,kstart,kend,oldcols = a->nbs,*lens;
80   PetscInt       row,mat_i,*mat_j,tcol,*mat_ilen;
81   const PetscInt *irow,*icol;
82   PetscInt       nrows,ncols,*ssmap,bs=A->rmap->bs,bs2=a->bs2;
83   PetscInt       *aj = a->j,*ai = a->i;
84   MatScalar      *mat_a;
85   Mat            C;
86   PetscBool      flag;
87 
88   PetscFunctionBegin;
89   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
90   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
91   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
92   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
93 
94   ierr  = PetscCalloc1(1+oldcols,&smap);CHKERRQ(ierr);
95   ssmap = smap;
96   ierr  = PetscMalloc1(1+nrows,&lens);CHKERRQ(ierr);
97   for (i=0; i<ncols; i++) smap[icol[i]] = i+1;
98   /* determine lens of each row */
99   for (i=0; i<nrows; i++) {
100     kstart  = ai[irow[i]];
101     kend    = kstart + a->ilen[irow[i]];
102     lens[i] = 0;
103     for (k=kstart; k<kend; k++) {
104       if (ssmap[aj[k]]) lens[i]++;
105     }
106   }
107   /* Create and fill new matrix */
108   if (scall == MAT_REUSE_MATRIX) {
109     c = (Mat_SeqBAIJ*)((*B)->data);
110 
111     if (c->mbs!=nrows || c->nbs!=ncols || (*B)->rmap->bs!=bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Submatrix wrong size");
112     ierr = PetscMemcmp(c->ilen,lens,c->mbs *sizeof(PetscInt),&flag);CHKERRQ(ierr);
113     if (!flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong no of nonzeros");
114     ierr = PetscMemzero(c->ilen,c->mbs*sizeof(PetscInt));CHKERRQ(ierr);
115     C    = *B;
116   } else {
117     ierr = MatCreate(PetscObjectComm((PetscObject)A),&C);CHKERRQ(ierr);
118     ierr = MatSetSizes(C,nrows*bs,ncols*bs,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
119     ierr = MatSetType(C,((PetscObject)A)->type_name);CHKERRQ(ierr);
120     ierr = MatSeqBAIJSetPreallocation(C,bs,0,lens);CHKERRQ(ierr);
121   }
122   c = (Mat_SeqBAIJ*)(C->data);
123   for (i=0; i<nrows; i++) {
124     row      = irow[i];
125     kstart   = ai[row];
126     kend     = kstart + a->ilen[row];
127     mat_i    = c->i[i];
128     mat_j    = c->j + mat_i;
129     mat_a    = c->a + mat_i*bs2;
130     mat_ilen = c->ilen + i;
131     for (k=kstart; k<kend; k++) {
132       if ((tcol=ssmap[a->j[k]])) {
133         *mat_j++ = tcol - 1;
134         ierr     = PetscMemcpy(mat_a,a->a+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
135         mat_a   += bs2;
136         (*mat_ilen)++;
137       }
138     }
139   }
140   /* sort */
141   {
142     MatScalar *work;
143     ierr = PetscMalloc1(bs2,&work);CHKERRQ(ierr);
144     for (i=0; i<nrows; i++) {
145       PetscInt ilen;
146       mat_i = c->i[i];
147       mat_j = c->j + mat_i;
148       mat_a = c->a + mat_i*bs2;
149       ilen  = c->ilen[i];
150       ierr  = PetscSortIntWithDataArray(ilen,mat_j,mat_a,bs2*sizeof(MatScalar),work);CHKERRQ(ierr);
151     }
152     ierr = PetscFree(work);CHKERRQ(ierr);
153   }
154 
155   /* Free work space */
156   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
157   ierr = PetscFree(smap);CHKERRQ(ierr);
158   ierr = PetscFree(lens);CHKERRQ(ierr);
159   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
160   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
161 
162   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
163   *B   = C;
164   PetscFunctionReturn(0);
165 }
166 
167 PetscErrorCode MatCreateSubMatrix_SeqBAIJ(Mat A,IS isrow,IS iscol,MatReuse scall,Mat *B)
168 {
169   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
170   IS             is1,is2;
171   PetscErrorCode ierr;
172   PetscInt       *vary,*iary,nrows,ncols,i,bs=A->rmap->bs,count,maxmnbs,j;
173   const PetscInt *irow,*icol;
174 
175   PetscFunctionBegin;
176   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
177   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
178   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
179   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
180 
181   /* Verify if the indices corespond to each element in a block
182    and form the IS with compressed IS */
183   maxmnbs = PetscMax(a->mbs,a->nbs);
184   ierr = PetscMalloc2(maxmnbs,&vary,maxmnbs,&iary);CHKERRQ(ierr);
185   ierr = PetscMemzero(vary,a->mbs*sizeof(PetscInt));CHKERRQ(ierr);
186   for (i=0; i<nrows; i++) vary[irow[i]/bs]++;
187   for (i=0; i<a->mbs; i++) {
188     if (vary[i]!=0 && vary[i]!=bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Index set does not match blocks");
189   }
190   count = 0;
191   for (i=0; i<nrows; i++) {
192     j = irow[i] / bs;
193     if ((vary[j]--)==bs) iary[count++] = j;
194   }
195   ierr = ISCreateGeneral(PETSC_COMM_SELF,count,iary,PETSC_COPY_VALUES,&is1);CHKERRQ(ierr);
196 
197   ierr = PetscMemzero(vary,(a->nbs)*sizeof(PetscInt));CHKERRQ(ierr);
198   for (i=0; i<ncols; i++) vary[icol[i]/bs]++;
199   for (i=0; i<a->nbs; i++) {
200     if (vary[i]!=0 && vary[i]!=bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Internal error in PETSc");
201   }
202   count = 0;
203   for (i=0; i<ncols; i++) {
204     j = icol[i] / bs;
205     if ((vary[j]--)==bs) iary[count++] = j;
206   }
207   ierr = ISCreateGeneral(PETSC_COMM_SELF,count,iary,PETSC_COPY_VALUES,&is2);CHKERRQ(ierr);
208   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
209   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
210   ierr = PetscFree2(vary,iary);CHKERRQ(ierr);
211 
212   ierr = MatCreateSubMatrix_SeqBAIJ_Private(A,is1,is2,scall,B);CHKERRQ(ierr);
213   ierr = ISDestroy(&is1);CHKERRQ(ierr);
214   ierr = ISDestroy(&is2);CHKERRQ(ierr);
215   PetscFunctionReturn(0);
216 }
217 
218 PetscErrorCode MatDestroySubMatrix_SeqBAIJ(Mat C)
219 {
220   PetscErrorCode ierr;
221   Mat_SeqBAIJ    *c = (Mat_SeqBAIJ*)C->data;
222   Mat_SubSppt    *submatj = c->submatis1;
223 
224   PetscFunctionBegin;
225   ierr = submatj->destroy(C);CHKERRQ(ierr);
226   ierr = MatDestroySubMatrix_Private(submatj);CHKERRQ(ierr);
227   PetscFunctionReturn(0);
228 }
229 
230 PetscErrorCode MatDestroySubMatrices_SeqBAIJ(PetscInt n,Mat *mat[])
231 {
232   PetscErrorCode ierr;
233   PetscInt       i;
234   Mat            C;
235   Mat_SeqBAIJ    *c;
236   Mat_SubSppt    *submatj;
237 
238   PetscFunctionBegin;
239   for (i=0; i<n; i++) {
240     C       = (*mat)[i];
241     c       = (Mat_SeqBAIJ*)C->data;
242     submatj = c->submatis1;
243     if (submatj) {
244       ierr = submatj->destroy(C);CHKERRQ(ierr);
245       ierr = MatDestroySubMatrix_Private(submatj);CHKERRQ(ierr);
246       ierr = PetscLayoutDestroy(&C->rmap);CHKERRQ(ierr);
247       ierr = PetscLayoutDestroy(&C->cmap);CHKERRQ(ierr);
248       ierr = PetscHeaderDestroy(&C);CHKERRQ(ierr);
249     } else {
250       ierr = MatDestroy(&C);CHKERRQ(ierr);
251     }
252   }
253   ierr = PetscFree(*mat);CHKERRQ(ierr);
254   PetscFunctionReturn(0);
255 }
256 
257 PetscErrorCode MatCreateSubMatrices_SeqBAIJ(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
258 {
259   PetscErrorCode ierr;
260   PetscInt       i;
261 
262   PetscFunctionBegin;
263   if (scall == MAT_INITIAL_MATRIX) {
264     ierr = PetscCalloc1(n+1,B);CHKERRQ(ierr);
265   }
266 
267   for (i=0; i<n; i++) {
268     ierr = MatCreateSubMatrix_SeqBAIJ(A,irow[i],icol[i],scall,&(*B)[i]);CHKERRQ(ierr);
269   }
270   PetscFunctionReturn(0);
271 }
272 
273 /* -------------------------------------------------------*/
274 /* Should check that shapes of vectors and matrices match */
275 /* -------------------------------------------------------*/
276 
277 PetscErrorCode MatMult_SeqBAIJ_1(Mat A,Vec xx,Vec zz)
278 {
279   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
280   PetscScalar       *z,sum;
281   const PetscScalar *x;
282   const MatScalar   *v;
283   PetscErrorCode    ierr;
284   PetscInt          mbs,i,n;
285   const PetscInt    *idx,*ii,*ridx=NULL;
286   PetscBool         usecprow=a->compressedrow.use;
287 
288   PetscFunctionBegin;
289   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
290   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
291 
292   if (usecprow) {
293     mbs  = a->compressedrow.nrows;
294     ii   = a->compressedrow.i;
295     ridx = a->compressedrow.rindex;
296     ierr = PetscMemzero(z,a->mbs*sizeof(PetscScalar));CHKERRQ(ierr);
297   } else {
298     mbs = a->mbs;
299     ii  = a->i;
300   }
301 
302   for (i=0; i<mbs; i++) {
303     n   = ii[1] - ii[0];
304     v   = a->a + ii[0];
305     idx = a->j + ii[0];
306     ii++;
307     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);   /* Indices for the next row (assumes same size as this one) */
308     PetscPrefetchBlock(v+1*n,1*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
309     sum = 0.0;
310     PetscSparseDensePlusDot(sum,x,v,idx,n);
311     if (usecprow) {
312       z[ridx[i]] = sum;
313     } else {
314       z[i]       = sum;
315     }
316   }
317   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
318   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
319   ierr = PetscLogFlops(2.0*a->nz - a->nonzerorowcnt);CHKERRQ(ierr);
320   PetscFunctionReturn(0);
321 }
322 
323 PetscErrorCode MatMult_SeqBAIJ_2(Mat A,Vec xx,Vec zz)
324 {
325   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
326   PetscScalar       *z = 0,sum1,sum2,*zarray;
327   const PetscScalar *x,*xb;
328   PetscScalar       x1,x2;
329   const MatScalar   *v;
330   PetscErrorCode    ierr;
331   PetscInt          mbs,i,*idx,*ii,j,n,*ridx=NULL;
332   PetscBool         usecprow=a->compressedrow.use;
333 
334   PetscFunctionBegin;
335   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
336   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
337 
338   idx = a->j;
339   v   = a->a;
340   if (usecprow) {
341     mbs  = a->compressedrow.nrows;
342     ii   = a->compressedrow.i;
343     ridx = a->compressedrow.rindex;
344     ierr = PetscMemzero(zarray,2*a->mbs*sizeof(PetscScalar));CHKERRQ(ierr);
345   } else {
346     mbs = a->mbs;
347     ii  = a->i;
348     z   = zarray;
349   }
350 
351   for (i=0; i<mbs; i++) {
352     n           = ii[1] - ii[0]; ii++;
353     sum1        = 0.0; sum2 = 0.0;
354     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);   /* Indices for the next row (assumes same size as this one) */
355     PetscPrefetchBlock(v+4*n,4*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
356     for (j=0; j<n; j++) {
357       xb    = x + 2*(*idx++); x1 = xb[0]; x2 = xb[1];
358       sum1 += v[0]*x1 + v[2]*x2;
359       sum2 += v[1]*x1 + v[3]*x2;
360       v    += 4;
361     }
362     if (usecprow) z = zarray + 2*ridx[i];
363     z[0] = sum1; z[1] = sum2;
364     if (!usecprow) z += 2;
365   }
366   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
367   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
368   ierr = PetscLogFlops(8.0*a->nz - 2.0*a->nonzerorowcnt);CHKERRQ(ierr);
369   PetscFunctionReturn(0);
370 }
371 
372 PetscErrorCode MatMult_SeqBAIJ_3(Mat A,Vec xx,Vec zz)
373 {
374   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
375   PetscScalar       *z = 0,sum1,sum2,sum3,x1,x2,x3,*zarray;
376   const PetscScalar *x,*xb;
377   const MatScalar   *v;
378   PetscErrorCode    ierr;
379   PetscInt          mbs,i,*idx,*ii,j,n,*ridx=NULL;
380   PetscBool         usecprow=a->compressedrow.use;
381 
382 #if defined(PETSC_HAVE_PRAGMA_DISJOINT)
383 #pragma disjoint(*v,*z,*xb)
384 #endif
385 
386   PetscFunctionBegin;
387   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
388   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
389 
390   idx = a->j;
391   v   = a->a;
392   if (usecprow) {
393     mbs  = a->compressedrow.nrows;
394     ii   = a->compressedrow.i;
395     ridx = a->compressedrow.rindex;
396     ierr = PetscMemzero(zarray,3*a->mbs*sizeof(PetscScalar));CHKERRQ(ierr);
397   } else {
398     mbs = a->mbs;
399     ii  = a->i;
400     z   = zarray;
401   }
402 
403   for (i=0; i<mbs; i++) {
404     n           = ii[1] - ii[0]; ii++;
405     sum1        = 0.0; sum2 = 0.0; sum3 = 0.0;
406     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);   /* Indices for the next row (assumes same size as this one) */
407     PetscPrefetchBlock(v+9*n,9*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
408     for (j=0; j<n; j++) {
409       xb = x + 3*(*idx++);
410       x1 = xb[0];
411       x2 = xb[1];
412       x3 = xb[2];
413 
414       sum1 += v[0]*x1 + v[3]*x2 + v[6]*x3;
415       sum2 += v[1]*x1 + v[4]*x2 + v[7]*x3;
416       sum3 += v[2]*x1 + v[5]*x2 + v[8]*x3;
417       v    += 9;
418     }
419     if (usecprow) z = zarray + 3*ridx[i];
420     z[0] = sum1; z[1] = sum2; z[2] = sum3;
421     if (!usecprow) z += 3;
422   }
423   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
424   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
425   ierr = PetscLogFlops(18.0*a->nz - 3.0*a->nonzerorowcnt);CHKERRQ(ierr);
426   PetscFunctionReturn(0);
427 }
428 
429 PetscErrorCode MatMult_SeqBAIJ_4(Mat A,Vec xx,Vec zz)
430 {
431   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
432   PetscScalar       *z = 0,sum1,sum2,sum3,sum4,x1,x2,x3,x4,*zarray;
433   const PetscScalar *x,*xb;
434   const MatScalar   *v;
435   PetscErrorCode    ierr;
436   PetscInt          mbs,i,*idx,*ii,j,n,*ridx=NULL;
437   PetscBool         usecprow=a->compressedrow.use;
438 
439   PetscFunctionBegin;
440   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
441   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
442 
443   idx = a->j;
444   v   = a->a;
445   if (usecprow) {
446     mbs  = a->compressedrow.nrows;
447     ii   = a->compressedrow.i;
448     ridx = a->compressedrow.rindex;
449     ierr = PetscMemzero(zarray,4*a->mbs*sizeof(PetscScalar));CHKERRQ(ierr);
450   } else {
451     mbs = a->mbs;
452     ii  = a->i;
453     z   = zarray;
454   }
455 
456   for (i=0; i<mbs; i++) {
457     n = ii[1] - ii[0];
458     ii++;
459     sum1 = 0.0;
460     sum2 = 0.0;
461     sum3 = 0.0;
462     sum4 = 0.0;
463 
464     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);     /* Indices for the next row (assumes same size as this one) */
465     PetscPrefetchBlock(v+16*n,16*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
466     for (j=0; j<n; j++) {
467       xb    = x + 4*(*idx++);
468       x1    = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
469       sum1 += v[0]*x1 + v[4]*x2 + v[8]*x3   + v[12]*x4;
470       sum2 += v[1]*x1 + v[5]*x2 + v[9]*x3   + v[13]*x4;
471       sum3 += v[2]*x1 + v[6]*x2 + v[10]*x3  + v[14]*x4;
472       sum4 += v[3]*x1 + v[7]*x2 + v[11]*x3  + v[15]*x4;
473       v    += 16;
474     }
475     if (usecprow) z = zarray + 4*ridx[i];
476     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4;
477     if (!usecprow) z += 4;
478   }
479   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
480   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
481   ierr = PetscLogFlops(32.0*a->nz - 4.0*a->nonzerorowcnt);CHKERRQ(ierr);
482   PetscFunctionReturn(0);
483 }
484 
485 PetscErrorCode MatMult_SeqBAIJ_5(Mat A,Vec xx,Vec zz)
486 {
487   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
488   PetscScalar       sum1,sum2,sum3,sum4,sum5,x1,x2,x3,x4,x5,*z = 0,*zarray;
489   const PetscScalar *xb,*x;
490   const MatScalar   *v;
491   PetscErrorCode    ierr;
492   const PetscInt    *idx,*ii,*ridx=NULL;
493   PetscInt          mbs,i,j,n;
494   PetscBool         usecprow=a->compressedrow.use;
495 
496   PetscFunctionBegin;
497   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
498   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
499 
500   idx = a->j;
501   v   = a->a;
502   if (usecprow) {
503     mbs  = a->compressedrow.nrows;
504     ii   = a->compressedrow.i;
505     ridx = a->compressedrow.rindex;
506     ierr = PetscMemzero(zarray,5*a->mbs*sizeof(PetscScalar));CHKERRQ(ierr);
507   } else {
508     mbs = a->mbs;
509     ii  = a->i;
510     z   = zarray;
511   }
512 
513   for (i=0; i<mbs; i++) {
514     n           = ii[1] - ii[0]; ii++;
515     sum1        = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0;
516     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);     /* Indices for the next row (assumes same size as this one) */
517     PetscPrefetchBlock(v+25*n,25*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
518     for (j=0; j<n; j++) {
519       xb    = x + 5*(*idx++);
520       x1    = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4];
521       sum1 += v[0]*x1 + v[5]*x2 + v[10]*x3  + v[15]*x4 + v[20]*x5;
522       sum2 += v[1]*x1 + v[6]*x2 + v[11]*x3  + v[16]*x4 + v[21]*x5;
523       sum3 += v[2]*x1 + v[7]*x2 + v[12]*x3  + v[17]*x4 + v[22]*x5;
524       sum4 += v[3]*x1 + v[8]*x2 + v[13]*x3  + v[18]*x4 + v[23]*x5;
525       sum5 += v[4]*x1 + v[9]*x2 + v[14]*x3  + v[19]*x4 + v[24]*x5;
526       v    += 25;
527     }
528     if (usecprow) z = zarray + 5*ridx[i];
529     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5;
530     if (!usecprow) z += 5;
531   }
532   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
533   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
534   ierr = PetscLogFlops(50.0*a->nz - 5.0*a->nonzerorowcnt);CHKERRQ(ierr);
535   PetscFunctionReturn(0);
536 }
537 
538 
539 
540 PetscErrorCode MatMult_SeqBAIJ_6(Mat A,Vec xx,Vec zz)
541 {
542   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
543   PetscScalar       *z = 0,sum1,sum2,sum3,sum4,sum5,sum6;
544   const PetscScalar *x,*xb;
545   PetscScalar       x1,x2,x3,x4,x5,x6,*zarray;
546   const MatScalar   *v;
547   PetscErrorCode    ierr;
548   PetscInt          mbs,i,*idx,*ii,j,n,*ridx=NULL;
549   PetscBool         usecprow=a->compressedrow.use;
550 
551   PetscFunctionBegin;
552   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
553   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
554 
555   idx = a->j;
556   v   = a->a;
557   if (usecprow) {
558     mbs  = a->compressedrow.nrows;
559     ii   = a->compressedrow.i;
560     ridx = a->compressedrow.rindex;
561     ierr = PetscMemzero(zarray,6*a->mbs*sizeof(PetscScalar));CHKERRQ(ierr);
562   } else {
563     mbs = a->mbs;
564     ii  = a->i;
565     z   = zarray;
566   }
567 
568   for (i=0; i<mbs; i++) {
569     n  = ii[1] - ii[0];
570     ii++;
571     sum1 = 0.0;
572     sum2 = 0.0;
573     sum3 = 0.0;
574     sum4 = 0.0;
575     sum5 = 0.0;
576     sum6 = 0.0;
577 
578     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);     /* Indices for the next row (assumes same size as this one) */
579     PetscPrefetchBlock(v+36*n,36*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
580     for (j=0; j<n; j++) {
581       xb    = x + 6*(*idx++);
582       x1    = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5];
583       sum1 += v[0]*x1 + v[6]*x2  + v[12]*x3  + v[18]*x4 + v[24]*x5 + v[30]*x6;
584       sum2 += v[1]*x1 + v[7]*x2  + v[13]*x3  + v[19]*x4 + v[25]*x5 + v[31]*x6;
585       sum3 += v[2]*x1 + v[8]*x2  + v[14]*x3  + v[20]*x4 + v[26]*x5 + v[32]*x6;
586       sum4 += v[3]*x1 + v[9]*x2  + v[15]*x3  + v[21]*x4 + v[27]*x5 + v[33]*x6;
587       sum5 += v[4]*x1 + v[10]*x2 + v[16]*x3  + v[22]*x4 + v[28]*x5 + v[34]*x6;
588       sum6 += v[5]*x1 + v[11]*x2 + v[17]*x3  + v[23]*x4 + v[29]*x5 + v[35]*x6;
589       v    += 36;
590     }
591     if (usecprow) z = zarray + 6*ridx[i];
592     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6;
593     if (!usecprow) z += 6;
594   }
595 
596   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
597   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
598   ierr = PetscLogFlops(72.0*a->nz - 6.0*a->nonzerorowcnt);CHKERRQ(ierr);
599   PetscFunctionReturn(0);
600 }
601 
602 PetscErrorCode MatMult_SeqBAIJ_7(Mat A,Vec xx,Vec zz)
603 {
604   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
605   PetscScalar       *z = 0,sum1,sum2,sum3,sum4,sum5,sum6,sum7;
606   const PetscScalar *x,*xb;
607   PetscScalar       x1,x2,x3,x4,x5,x6,x7,*zarray;
608   const MatScalar   *v;
609   PetscErrorCode    ierr;
610   PetscInt          mbs,i,*idx,*ii,j,n,*ridx=NULL;
611   PetscBool         usecprow=a->compressedrow.use;
612 
613   PetscFunctionBegin;
614   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
615   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
616 
617   idx = a->j;
618   v   = a->a;
619   if (usecprow) {
620     mbs  = a->compressedrow.nrows;
621     ii   = a->compressedrow.i;
622     ridx = a->compressedrow.rindex;
623     ierr = PetscMemzero(zarray,7*a->mbs*sizeof(PetscScalar));CHKERRQ(ierr);
624   } else {
625     mbs = a->mbs;
626     ii  = a->i;
627     z   = zarray;
628   }
629 
630   for (i=0; i<mbs; i++) {
631     n  = ii[1] - ii[0];
632     ii++;
633     sum1 = 0.0;
634     sum2 = 0.0;
635     sum3 = 0.0;
636     sum4 = 0.0;
637     sum5 = 0.0;
638     sum6 = 0.0;
639     sum7 = 0.0;
640 
641     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);     /* Indices for the next row (assumes same size as this one) */
642     PetscPrefetchBlock(v+49*n,49*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
643     for (j=0; j<n; j++) {
644       xb    = x + 7*(*idx++);
645       x1    = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6];
646       sum1 += v[0]*x1 + v[7]*x2  + v[14]*x3  + v[21]*x4 + v[28]*x5 + v[35]*x6 + v[42]*x7;
647       sum2 += v[1]*x1 + v[8]*x2  + v[15]*x3  + v[22]*x4 + v[29]*x5 + v[36]*x6 + v[43]*x7;
648       sum3 += v[2]*x1 + v[9]*x2  + v[16]*x3  + v[23]*x4 + v[30]*x5 + v[37]*x6 + v[44]*x7;
649       sum4 += v[3]*x1 + v[10]*x2 + v[17]*x3  + v[24]*x4 + v[31]*x5 + v[38]*x6 + v[45]*x7;
650       sum5 += v[4]*x1 + v[11]*x2 + v[18]*x3  + v[25]*x4 + v[32]*x5 + v[39]*x6 + v[46]*x7;
651       sum6 += v[5]*x1 + v[12]*x2 + v[19]*x3  + v[26]*x4 + v[33]*x5 + v[40]*x6 + v[47]*x7;
652       sum7 += v[6]*x1 + v[13]*x2 + v[20]*x3  + v[27]*x4 + v[34]*x5 + v[41]*x6 + v[48]*x7;
653       v    += 49;
654     }
655     if (usecprow) z = zarray + 7*ridx[i];
656     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
657     if (!usecprow) z += 7;
658   }
659 
660   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
661   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
662   ierr = PetscLogFlops(98.0*a->nz - 7.0*a->nonzerorowcnt);CHKERRQ(ierr);
663   PetscFunctionReturn(0);
664 }
665 
666 #if defined(PETSC_HAVE_IMMINTRIN_H) && defined(__AVX2__) && defined(PETSC_USE_REAL_DOUBLE) && !defined(PETSC_USE_COMPLEX)
667 PetscErrorCode MatMult_SeqBAIJ_9_AVX2(Mat A,Vec xx,Vec zz)
668 {
669   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
670   PetscScalar    *z = 0,*work,*workt,*zarray;
671   const PetscScalar *x,*xb;
672   const MatScalar   *v;
673   PetscErrorCode ierr;
674   PetscInt       mbs,i,bs=A->rmap->bs,j,n,bs2=a->bs2;
675   const PetscInt *idx,*ii,*ridx=NULL;
676   PetscInt       ncols,k;
677   PetscBool      usecprow=a->compressedrow.use;
678 
679   __m256d a0,a1,a2,a3,a4,a5;
680   __m256d w0,w1,w2,w3,w4,w5;
681   __m256d z0,z1,z2;
682   __m256i mask1 = _mm256_set_epi64x(0LL, 0LL, 0LL, 1LL<<63);
683 
684   PetscFunctionBegin;
685   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
686   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
687 
688   idx = a->j;
689   v   = a->a;
690   if (usecprow) {
691     mbs  = a->compressedrow.nrows;
692     ii   = a->compressedrow.i;
693     ridx = a->compressedrow.rindex;
694     ierr = PetscMemzero(zarray,bs*a->mbs*sizeof(PetscScalar));CHKERRQ(ierr);
695   } else {
696     mbs = a->mbs;
697     ii  = a->i;
698     z   = zarray;
699   }
700 
701   if (!a->mult_work) {
702     k    = PetscMax(A->rmap->n,A->cmap->n);
703     ierr = PetscMalloc1(k+1,&a->mult_work);CHKERRQ(ierr);
704   }
705 
706   work = a->mult_work;
707   for (i=0; i<mbs; i++) {
708     n           = ii[1] - ii[0]; ii++;
709     ncols       = n*bs;
710     workt       = work;
711     for (j=0; j<n; j++) {
712       xb = x + bs*(*idx++);
713       for (k=0; k<bs; k++) workt[k] = xb[k];
714       workt += bs;
715     }
716     if (usecprow) z = zarray + bs*ridx[i];
717 
718     z0 = _mm256_setzero_pd(); z1 = _mm256_setzero_pd(); z2 = _mm256_setzero_pd();
719 
720     for (j=0; j<n; j++) {
721       // first column of a
722       w0 = _mm256_set1_pd(work[j*9  ]);
723       a0 = _mm256_loadu_pd(&v[j*81  ]); z0 = _mm256_fmadd_pd(a0,w0,z0);
724       a1 = _mm256_loadu_pd(&v[j*81+4]); z1 = _mm256_fmadd_pd(a1,w0,z1);
725       a2 = _mm256_loadu_pd(&v[j*81+8]); z2 = _mm256_fmadd_pd(a2,w0,z2);
726 
727       // second column of a
728       w1 = _mm256_set1_pd(work[j*9+ 1]);
729       a0 = _mm256_loadu_pd(&v[j*81+ 9]); z0 = _mm256_fmadd_pd(a0,w1,z0);
730       a1 = _mm256_loadu_pd(&v[j*81+13]); z1 = _mm256_fmadd_pd(a1,w1,z1);
731       a2 = _mm256_loadu_pd(&v[j*81+17]); z2 = _mm256_fmadd_pd(a2,w1,z2);
732 
733       // third column of a
734       w2 = _mm256_set1_pd(work[j*9 +2]);
735       a3 = _mm256_loadu_pd(&v[j*81+18]); z0 = _mm256_fmadd_pd(a3,w2,z0);
736       a4 = _mm256_loadu_pd(&v[j*81+22]); z1 = _mm256_fmadd_pd(a4,w2,z1);
737       a5 = _mm256_loadu_pd(&v[j*81+26]); z2 = _mm256_fmadd_pd(a5,w2,z2);
738 
739       // fourth column of a
740       w3 = _mm256_set1_pd(work[j*9+ 3]);
741       a0 = _mm256_loadu_pd(&v[j*81+27]); z0 = _mm256_fmadd_pd(a0,w3,z0);
742       a1 = _mm256_loadu_pd(&v[j*81+31]); z1 = _mm256_fmadd_pd(a1,w3,z1);
743       a2 = _mm256_loadu_pd(&v[j*81+35]); z2 = _mm256_fmadd_pd(a2,w3,z2);
744 
745       // fifth column of a
746       w0 = _mm256_set1_pd(work[j*9+ 4]);
747       a3 = _mm256_loadu_pd(&v[j*81+36]); z0 = _mm256_fmadd_pd(a3,w0,z0);
748       a4 = _mm256_loadu_pd(&v[j*81+40]); z1 = _mm256_fmadd_pd(a4,w0,z1);
749       a5 = _mm256_loadu_pd(&v[j*81+44]); z2 = _mm256_fmadd_pd(a5,w0,z2);
750 
751       // sixth column of a
752       w1 = _mm256_set1_pd(work[j*9+ 5]);
753       a0 = _mm256_loadu_pd(&v[j*81+45]); z0 = _mm256_fmadd_pd(a0,w1,z0);
754       a1 = _mm256_loadu_pd(&v[j*81+49]); z1 = _mm256_fmadd_pd(a1,w1,z1);
755       a2 = _mm256_loadu_pd(&v[j*81+53]); z2 = _mm256_fmadd_pd(a2,w1,z2);
756 
757       // seventh column of a
758       w2 = _mm256_set1_pd(work[j*9+ 6]);
759       a0 = _mm256_loadu_pd(&v[j*81+54]); z0 = _mm256_fmadd_pd(a0,w2,z0);
760       a1 = _mm256_loadu_pd(&v[j*81+58]); z1 = _mm256_fmadd_pd(a1,w2,z1);
761       a2 = _mm256_loadu_pd(&v[j*81+62]); z2 = _mm256_fmadd_pd(a2,w2,z2);
762 
763       // eigth column of a
764       w3 = _mm256_set1_pd(work[j*9+ 7]);
765       a3 = _mm256_loadu_pd(&v[j*81+63]); z0 = _mm256_fmadd_pd(a3,w3,z0);
766       a4 = _mm256_loadu_pd(&v[j*81+67]); z1 = _mm256_fmadd_pd(a4,w3,z1);
767       a5 = _mm256_loadu_pd(&v[j*81+71]); z2 = _mm256_fmadd_pd(a5,w3,z2);
768 
769       // ninth column of a
770       w0 = _mm256_set1_pd(work[j*9+ 8]);
771       a0 = _mm256_loadu_pd(&v[j*81+72]); z0 = _mm256_fmadd_pd(a0,w0,z0);
772       a1 = _mm256_loadu_pd(&v[j*81+76]); z1 = _mm256_fmadd_pd(a1,w0,z1);
773       a2 = _mm256_maskload_pd(&v[j*81+80],mask1); z2 = _mm256_fmadd_pd(a2,w0,z2);
774     }
775 
776     _mm256_storeu_pd(&z[ 0], z0); _mm256_storeu_pd(&z[ 4], z1); _mm256_maskstore_pd(&z[8], mask1, z2);
777 
778     v += n*bs2;
779     if (!usecprow) z += bs;
780   }
781   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
782   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
783   ierr = PetscLogFlops(2.0*a->nz*bs2 - bs*a->nonzerorowcnt);CHKERRQ(ierr);
784   PetscFunctionReturn(0);
785 }
786 #endif
787 
788 PetscErrorCode MatMult_SeqBAIJ_11(Mat A,Vec xx,Vec zz)
789 {
790   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
791   PetscScalar       *z = 0,sum1,sum2,sum3,sum4,sum5,sum6,sum7,sum8,sum9,sum10,sum11;
792   const PetscScalar *x,*xb;
793   PetscScalar       *zarray,xv;
794   const MatScalar   *v;
795   PetscErrorCode    ierr;
796   const PetscInt    *ii,*ij=a->j,*idx;
797   PetscInt          mbs,i,j,k,n,*ridx=NULL;
798   PetscBool         usecprow=a->compressedrow.use;
799 
800   PetscFunctionBegin;
801   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
802   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
803 
804   v = a->a;
805   if (usecprow) {
806     mbs  = a->compressedrow.nrows;
807     ii   = a->compressedrow.i;
808     ridx = a->compressedrow.rindex;
809     ierr = PetscMemzero(zarray,11*a->mbs*sizeof(PetscScalar));CHKERRQ(ierr);
810   } else {
811     mbs = a->mbs;
812     ii  = a->i;
813     z   = zarray;
814   }
815 
816   for (i=0; i<mbs; i++) {
817     n    = ii[i+1] - ii[i];
818     idx  = ij + ii[i];
819     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0; sum7 = 0.0;
820     sum8 = 0.0; sum9 = 0.0; sum10 = 0.0; sum11 = 0.0;
821 
822     for (j=0; j<n; j++) {
823       xb = x + 11*(idx[j]);
824 
825       for (k=0; k<11; k++) {
826         xv     =  xb[k];
827         sum1  += v[0]*xv;
828         sum2  += v[1]*xv;
829         sum3  += v[2]*xv;
830         sum4  += v[3]*xv;
831         sum5  += v[4]*xv;
832         sum6  += v[5]*xv;
833         sum7  += v[6]*xv;
834         sum8  += v[7]*xv;
835         sum9  += v[8]*xv;
836         sum10 += v[9]*xv;
837         sum11 += v[10]*xv;
838         v     += 11;
839       }
840     }
841     if (usecprow) z = zarray + 11*ridx[i];
842     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
843     z[7] = sum8; z[8] = sum9; z[9] = sum10; z[10] = sum11;
844 
845     if (!usecprow) z += 11;
846   }
847 
848   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
849   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
850   ierr = PetscLogFlops(242.0*a->nz - 11.0*a->nonzerorowcnt);CHKERRQ(ierr);
851   PetscFunctionReturn(0);
852 }
853 
854 /* MatMult_SeqBAIJ_15 version 1: Columns in the block are accessed one at a time */
855 /* Default MatMult for block size 15 */
856 
857 PetscErrorCode MatMult_SeqBAIJ_15_ver1(Mat A,Vec xx,Vec zz)
858 {
859   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
860   PetscScalar       *z = 0,sum1,sum2,sum3,sum4,sum5,sum6,sum7,sum8,sum9,sum10,sum11,sum12,sum13,sum14,sum15;
861   const PetscScalar *x,*xb;
862   PetscScalar       *zarray,xv;
863   const MatScalar   *v;
864   PetscErrorCode    ierr;
865   const PetscInt    *ii,*ij=a->j,*idx;
866   PetscInt          mbs,i,j,k,n,*ridx=NULL;
867   PetscBool         usecprow=a->compressedrow.use;
868 
869   PetscFunctionBegin;
870   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
871   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
872 
873   v = a->a;
874   if (usecprow) {
875     mbs  = a->compressedrow.nrows;
876     ii   = a->compressedrow.i;
877     ridx = a->compressedrow.rindex;
878     ierr = PetscMemzero(zarray,15*a->mbs*sizeof(PetscScalar));CHKERRQ(ierr);
879   } else {
880     mbs = a->mbs;
881     ii  = a->i;
882     z   = zarray;
883   }
884 
885   for (i=0; i<mbs; i++) {
886     n    = ii[i+1] - ii[i];
887     idx  = ij + ii[i];
888     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0; sum7 = 0.0;
889     sum8 = 0.0; sum9 = 0.0; sum10 = 0.0; sum11 = 0.0; sum12 = 0.0; sum13 = 0.0; sum14 = 0.0;sum15 = 0.0;
890 
891     for (j=0; j<n; j++) {
892       xb = x + 15*(idx[j]);
893 
894       for (k=0; k<15; k++) {
895         xv     =  xb[k];
896         sum1  += v[0]*xv;
897         sum2  += v[1]*xv;
898         sum3  += v[2]*xv;
899         sum4  += v[3]*xv;
900         sum5  += v[4]*xv;
901         sum6  += v[5]*xv;
902         sum7  += v[6]*xv;
903         sum8  += v[7]*xv;
904         sum9  += v[8]*xv;
905         sum10 += v[9]*xv;
906         sum11 += v[10]*xv;
907         sum12 += v[11]*xv;
908         sum13 += v[12]*xv;
909         sum14 += v[13]*xv;
910         sum15 += v[14]*xv;
911         v     += 15;
912       }
913     }
914     if (usecprow) z = zarray + 15*ridx[i];
915     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
916     z[7] = sum8; z[8] = sum9; z[9] = sum10; z[10] = sum11; z[11] = sum12; z[12] = sum13; z[13] = sum14;z[14] = sum15;
917 
918     if (!usecprow) z += 15;
919   }
920 
921   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
922   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
923   ierr = PetscLogFlops(450.0*a->nz - 15.0*a->nonzerorowcnt);CHKERRQ(ierr);
924   PetscFunctionReturn(0);
925 }
926 
927 /* MatMult_SeqBAIJ_15_ver2 : Columns in the block are accessed in sets of 4,4,4,3 */
928 PetscErrorCode MatMult_SeqBAIJ_15_ver2(Mat A,Vec xx,Vec zz)
929 {
930   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
931   PetscScalar       *z = 0,sum1,sum2,sum3,sum4,sum5,sum6,sum7,sum8,sum9,sum10,sum11,sum12,sum13,sum14,sum15;
932   const PetscScalar *x,*xb;
933   PetscScalar       x1,x2,x3,x4,*zarray;
934   const MatScalar   *v;
935   PetscErrorCode    ierr;
936   const PetscInt    *ii,*ij=a->j,*idx;
937   PetscInt          mbs,i,j,n,*ridx=NULL;
938   PetscBool         usecprow=a->compressedrow.use;
939 
940   PetscFunctionBegin;
941   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
942   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
943 
944   v = a->a;
945   if (usecprow) {
946     mbs  = a->compressedrow.nrows;
947     ii   = a->compressedrow.i;
948     ridx = a->compressedrow.rindex;
949     ierr = PetscMemzero(zarray,15*a->mbs*sizeof(PetscScalar));CHKERRQ(ierr);
950   } else {
951     mbs = a->mbs;
952     ii  = a->i;
953     z   = zarray;
954   }
955 
956   for (i=0; i<mbs; i++) {
957     n    = ii[i+1] - ii[i];
958     idx  = ij + ii[i];
959     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0; sum7 = 0.0;
960     sum8 = 0.0; sum9 = 0.0; sum10 = 0.0; sum11 = 0.0; sum12 = 0.0; sum13 = 0.0; sum14 = 0.0;sum15 = 0.0;
961 
962     for (j=0; j<n; j++) {
963       xb = x + 15*(idx[j]);
964       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
965 
966       sum1  += v[0]*x1 + v[15]*x2 + v[30]*x3   + v[45]*x4;
967       sum2  += v[1]*x1 + v[16]*x2 + v[31]*x3   + v[46]*x4;
968       sum3  += v[2]*x1 + v[17]*x2 + v[32]*x3  + v[47]*x4;
969       sum4  += v[3]*x1 + v[18]*x2 + v[33]*x3  + v[48]*x4;
970       sum5  += v[4]*x1 + v[19]*x2 + v[34]*x3   + v[49]*x4;
971       sum6  += v[5]*x1 + v[20]*x2 + v[35]*x3   + v[50]*x4;
972       sum7  += v[6]*x1 + v[21]*x2 + v[36]*x3  + v[51]*x4;
973       sum8  += v[7]*x1 + v[22]*x2 + v[37]*x3  + v[52]*x4;
974       sum9  += v[8]*x1 + v[23]*x2 + v[38]*x3   + v[53]*x4;
975       sum10 += v[9]*x1 + v[24]*x2 + v[39]*x3   + v[54]*x4;
976       sum11 += v[10]*x1 + v[25]*x2 + v[40]*x3  + v[55]*x4;
977       sum12 += v[11]*x1 + v[26]*x2 + v[41]*x3  + v[56]*x4;
978       sum13 += v[12]*x1 + v[27]*x2 + v[42]*x3   + v[57]*x4;
979       sum14 += v[13]*x1 + v[28]*x2 + v[43]*x3   + v[58]*x4;
980       sum15 += v[14]*x1 + v[29]*x2 + v[44]*x3  + v[59]*x4;
981 
982       v += 60;
983 
984       x1 = xb[4]; x2 = xb[5]; x3 = xb[6]; x4 = xb[7];
985 
986       sum1  += v[0]*x1 + v[15]*x2 + v[30]*x3   + v[45]*x4;
987       sum2  += v[1]*x1 + v[16]*x2 + v[31]*x3   + v[46]*x4;
988       sum3  += v[2]*x1 + v[17]*x2 + v[32]*x3  + v[47]*x4;
989       sum4  += v[3]*x1 + v[18]*x2 + v[33]*x3  + v[48]*x4;
990       sum5  += v[4]*x1 + v[19]*x2 + v[34]*x3   + v[49]*x4;
991       sum6  += v[5]*x1 + v[20]*x2 + v[35]*x3   + v[50]*x4;
992       sum7  += v[6]*x1 + v[21]*x2 + v[36]*x3  + v[51]*x4;
993       sum8  += v[7]*x1 + v[22]*x2 + v[37]*x3  + v[52]*x4;
994       sum9  += v[8]*x1 + v[23]*x2 + v[38]*x3   + v[53]*x4;
995       sum10 += v[9]*x1 + v[24]*x2 + v[39]*x3   + v[54]*x4;
996       sum11 += v[10]*x1 + v[25]*x2 + v[40]*x3  + v[55]*x4;
997       sum12 += v[11]*x1 + v[26]*x2 + v[41]*x3  + v[56]*x4;
998       sum13 += v[12]*x1 + v[27]*x2 + v[42]*x3   + v[57]*x4;
999       sum14 += v[13]*x1 + v[28]*x2 + v[43]*x3   + v[58]*x4;
1000       sum15 += v[14]*x1 + v[29]*x2 + v[44]*x3  + v[59]*x4;
1001       v     += 60;
1002 
1003       x1     = xb[8]; x2 = xb[9]; x3 = xb[10]; x4 = xb[11];
1004       sum1  += v[0]*x1 + v[15]*x2 + v[30]*x3   + v[45]*x4;
1005       sum2  += v[1]*x1 + v[16]*x2 + v[31]*x3   + v[46]*x4;
1006       sum3  += v[2]*x1 + v[17]*x2 + v[32]*x3  + v[47]*x4;
1007       sum4  += v[3]*x1 + v[18]*x2 + v[33]*x3  + v[48]*x4;
1008       sum5  += v[4]*x1 + v[19]*x2 + v[34]*x3   + v[49]*x4;
1009       sum6  += v[5]*x1 + v[20]*x2 + v[35]*x3   + v[50]*x4;
1010       sum7  += v[6]*x1 + v[21]*x2 + v[36]*x3  + v[51]*x4;
1011       sum8  += v[7]*x1 + v[22]*x2 + v[37]*x3  + v[52]*x4;
1012       sum9  += v[8]*x1 + v[23]*x2 + v[38]*x3   + v[53]*x4;
1013       sum10 += v[9]*x1 + v[24]*x2 + v[39]*x3   + v[54]*x4;
1014       sum11 += v[10]*x1 + v[25]*x2 + v[40]*x3  + v[55]*x4;
1015       sum12 += v[11]*x1 + v[26]*x2 + v[41]*x3  + v[56]*x4;
1016       sum13 += v[12]*x1 + v[27]*x2 + v[42]*x3   + v[57]*x4;
1017       sum14 += v[13]*x1 + v[28]*x2 + v[43]*x3   + v[58]*x4;
1018       sum15 += v[14]*x1 + v[29]*x2 + v[44]*x3  + v[59]*x4;
1019       v     += 60;
1020 
1021       x1     = xb[12]; x2 = xb[13]; x3 = xb[14];
1022       sum1  += v[0]*x1 + v[15]*x2 + v[30]*x3;
1023       sum2  += v[1]*x1 + v[16]*x2 + v[31]*x3;
1024       sum3  += v[2]*x1 + v[17]*x2 + v[32]*x3;
1025       sum4  += v[3]*x1 + v[18]*x2 + v[33]*x3;
1026       sum5  += v[4]*x1 + v[19]*x2 + v[34]*x3;
1027       sum6  += v[5]*x1 + v[20]*x2 + v[35]*x3;
1028       sum7  += v[6]*x1 + v[21]*x2 + v[36]*x3;
1029       sum8  += v[7]*x1 + v[22]*x2 + v[37]*x3;
1030       sum9  += v[8]*x1 + v[23]*x2 + v[38]*x3;
1031       sum10 += v[9]*x1 + v[24]*x2 + v[39]*x3;
1032       sum11 += v[10]*x1 + v[25]*x2 + v[40]*x3;
1033       sum12 += v[11]*x1 + v[26]*x2 + v[41]*x3;
1034       sum13 += v[12]*x1 + v[27]*x2 + v[42]*x3;
1035       sum14 += v[13]*x1 + v[28]*x2 + v[43]*x3;
1036       sum15 += v[14]*x1 + v[29]*x2 + v[44]*x3;
1037       v     += 45;
1038     }
1039     if (usecprow) z = zarray + 15*ridx[i];
1040     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
1041     z[7] = sum8; z[8] = sum9; z[9] = sum10; z[10] = sum11; z[11] = sum12; z[12] = sum13; z[13] = sum14;z[14] = sum15;
1042 
1043     if (!usecprow) z += 15;
1044   }
1045 
1046   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
1047   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
1048   ierr = PetscLogFlops(450.0*a->nz - 15.0*a->nonzerorowcnt);CHKERRQ(ierr);
1049   PetscFunctionReturn(0);
1050 }
1051 
1052 /* MatMult_SeqBAIJ_15_ver3 : Columns in the block are accessed in sets of 8,7 */
1053 PetscErrorCode MatMult_SeqBAIJ_15_ver3(Mat A,Vec xx,Vec zz)
1054 {
1055   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
1056   PetscScalar       *z = 0,sum1,sum2,sum3,sum4,sum5,sum6,sum7,sum8,sum9,sum10,sum11,sum12,sum13,sum14,sum15;
1057   const PetscScalar *x,*xb;
1058   PetscScalar       x1,x2,x3,x4,x5,x6,x7,x8,*zarray;
1059   const MatScalar   *v;
1060   PetscErrorCode    ierr;
1061   const PetscInt    *ii,*ij=a->j,*idx;
1062   PetscInt          mbs,i,j,n,*ridx=NULL;
1063   PetscBool         usecprow=a->compressedrow.use;
1064 
1065   PetscFunctionBegin;
1066   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
1067   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
1068 
1069   v = a->a;
1070   if (usecprow) {
1071     mbs  = a->compressedrow.nrows;
1072     ii   = a->compressedrow.i;
1073     ridx = a->compressedrow.rindex;
1074     ierr = PetscMemzero(zarray,15*a->mbs*sizeof(PetscScalar));CHKERRQ(ierr);
1075   } else {
1076     mbs = a->mbs;
1077     ii  = a->i;
1078     z   = zarray;
1079   }
1080 
1081   for (i=0; i<mbs; i++) {
1082     n    = ii[i+1] - ii[i];
1083     idx  = ij + ii[i];
1084     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0; sum7 = 0.0;
1085     sum8 = 0.0; sum9 = 0.0; sum10 = 0.0; sum11 = 0.0; sum12 = 0.0; sum13 = 0.0; sum14 = 0.0;sum15 = 0.0;
1086 
1087     for (j=0; j<n; j++) {
1088       xb = x + 15*(idx[j]);
1089       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6];
1090       x8 = xb[7];
1091 
1092       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;
1093       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;
1094       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;
1095       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;
1096       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;
1097       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;
1098       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;
1099       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;
1100       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;
1101       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;
1102       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;
1103       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;
1104       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;
1105       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;
1106       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;
1107       v     += 120;
1108 
1109       x1 = xb[8]; x2 = xb[9]; x3 = xb[10]; x4 = xb[11]; x5 = xb[12]; x6 = xb[13]; x7 = xb[14];
1110 
1111       sum1  += v[0]*x1 + v[15]*x2  + v[30]*x3  + v[45]*x4 + v[60]*x5 + v[75]*x6 + v[90]*x7;
1112       sum2  += v[1]*x1 + v[16]*x2  + v[31]*x3  + v[46]*x4 + v[61]*x5 + v[76]*x6 + v[91]*x7;
1113       sum3  += v[2]*x1 + v[17]*x2  + v[32]*x3  + v[47]*x4 + v[62]*x5 + v[77]*x6 + v[92]*x7;
1114       sum4  += v[3]*x1 + v[18]*x2 + v[33]*x3  + v[48]*x4 + v[63]*x5 + v[78]*x6 + v[93]*x7;
1115       sum5  += v[4]*x1 + v[19]*x2 + v[34]*x3  + v[49]*x4 + v[64]*x5 + v[79]*x6 + v[94]*x7;
1116       sum6  += v[5]*x1 + v[20]*x2 + v[35]*x3  + v[50]*x4 + v[65]*x5 + v[80]*x6 + v[95]*x7;
1117       sum7  += v[6]*x1 + v[21]*x2 + v[36]*x3  + v[51]*x4 + v[66]*x5 + v[81]*x6 + v[96]*x7;
1118       sum8  += v[7]*x1 + v[22]*x2  + v[37]*x3  + v[52]*x4 + v[67]*x5 + v[82]*x6 + v[97]*x7;
1119       sum9  += v[8]*x1 + v[23]*x2  + v[38]*x3  + v[53]*x4 + v[68]*x5 + v[83]*x6 + v[98]*x7;
1120       sum10 += v[9]*x1 + v[24]*x2  + v[39]*x3  + v[54]*x4 + v[69]*x5 + v[84]*x6 + v[99]*x7;
1121       sum11 += v[10]*x1 + v[25]*x2 + v[40]*x3  + v[55]*x4 + v[70]*x5 + v[85]*x6 + v[100]*x7;
1122       sum12 += v[11]*x1 + v[26]*x2 + v[41]*x3  + v[56]*x4 + v[71]*x5 + v[86]*x6 + v[101]*x7;
1123       sum13 += v[12]*x1 + v[27]*x2 + v[42]*x3  + v[57]*x4 + v[72]*x5 + v[87]*x6 + v[102]*x7;
1124       sum14 += v[13]*x1 + v[28]*x2 + v[43]*x3  + v[58]*x4 + v[73]*x5 + v[88]*x6 + v[103]*x7;
1125       sum15 += v[14]*x1 + v[29]*x2 + v[44]*x3  + v[59]*x4 + v[74]*x5 + v[89]*x6 + v[104]*x7;
1126       v     += 105;
1127     }
1128     if (usecprow) z = zarray + 15*ridx[i];
1129     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
1130     z[7] = sum8; z[8] = sum9; z[9] = sum10; z[10] = sum11; z[11] = sum12; z[12] = sum13; z[13] = sum14;z[14] = sum15;
1131 
1132     if (!usecprow) z += 15;
1133   }
1134 
1135   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
1136   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
1137   ierr = PetscLogFlops(450.0*a->nz - 15.0*a->nonzerorowcnt);CHKERRQ(ierr);
1138   PetscFunctionReturn(0);
1139 }
1140 
1141 /* MatMult_SeqBAIJ_15_ver4 : All columns in the block are accessed at once */
1142 
1143 PetscErrorCode MatMult_SeqBAIJ_15_ver4(Mat A,Vec xx,Vec zz)
1144 {
1145   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
1146   PetscScalar       *z = 0,sum1,sum2,sum3,sum4,sum5,sum6,sum7,sum8,sum9,sum10,sum11,sum12,sum13,sum14,sum15;
1147   const PetscScalar *x,*xb;
1148   PetscScalar       x1,x2,x3,x4,x5,x6,x7,x8,x9,x10,x11,x12,x13,x14,x15,*zarray;
1149   const MatScalar   *v;
1150   PetscErrorCode    ierr;
1151   const PetscInt    *ii,*ij=a->j,*idx;
1152   PetscInt          mbs,i,j,n,*ridx=NULL;
1153   PetscBool         usecprow=a->compressedrow.use;
1154 
1155   PetscFunctionBegin;
1156   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
1157   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
1158 
1159   v = a->a;
1160   if (usecprow) {
1161     mbs  = a->compressedrow.nrows;
1162     ii   = a->compressedrow.i;
1163     ridx = a->compressedrow.rindex;
1164     ierr = PetscMemzero(zarray,15*a->mbs*sizeof(PetscScalar));CHKERRQ(ierr);
1165   } else {
1166     mbs = a->mbs;
1167     ii  = a->i;
1168     z   = zarray;
1169   }
1170 
1171   for (i=0; i<mbs; i++) {
1172     n    = ii[i+1] - ii[i];
1173     idx  = ij + ii[i];
1174     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0; sum7 = 0.0;
1175     sum8 = 0.0; sum9 = 0.0; sum10 = 0.0; sum11 = 0.0; sum12 = 0.0; sum13 = 0.0; sum14 = 0.0;sum15 = 0.0;
1176 
1177     for (j=0; j<n; j++) {
1178       xb = x + 15*(idx[j]);
1179       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6];
1180       x8 = xb[7]; x9 = xb[8]; x10 = xb[9]; x11 = xb[10]; x12 = xb[11]; x13 = xb[12]; x14 = xb[13];x15 = xb[14];
1181 
1182       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;
1183       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;
1184       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;
1185       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;
1186       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;
1187       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;
1188       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;
1189       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;
1190       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;
1191       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;
1192       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;
1193       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;
1194       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;
1195       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;
1196       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;
1197       v     += 225;
1198     }
1199     if (usecprow) z = zarray + 15*ridx[i];
1200     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
1201     z[7] = sum8; z[8] = sum9; z[9] = sum10; z[10] = sum11; z[11] = sum12; z[12] = sum13; z[13] = sum14;z[14] = sum15;
1202 
1203     if (!usecprow) z += 15;
1204   }
1205 
1206   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
1207   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
1208   ierr = PetscLogFlops(450.0*a->nz - 15.0*a->nonzerorowcnt);CHKERRQ(ierr);
1209   PetscFunctionReturn(0);
1210 }
1211 
1212 
1213 /*
1214     This will not work with MatScalar == float because it calls the BLAS
1215 */
1216 PetscErrorCode MatMult_SeqBAIJ_N(Mat A,Vec xx,Vec zz)
1217 {
1218   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
1219   PetscScalar       *z = 0,*work,*workt,*zarray;
1220   const PetscScalar *x,*xb;
1221   const MatScalar   *v;
1222   PetscErrorCode    ierr;
1223   PetscInt          mbs,i,bs=A->rmap->bs,j,n,bs2=a->bs2;
1224   const PetscInt    *idx,*ii,*ridx=NULL;
1225   PetscInt          ncols,k;
1226   PetscBool         usecprow=a->compressedrow.use;
1227 
1228   PetscFunctionBegin;
1229   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
1230   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
1231 
1232   idx = a->j;
1233   v   = a->a;
1234   if (usecprow) {
1235     mbs  = a->compressedrow.nrows;
1236     ii   = a->compressedrow.i;
1237     ridx = a->compressedrow.rindex;
1238     ierr = PetscMemzero(zarray,bs*a->mbs*sizeof(PetscScalar));CHKERRQ(ierr);
1239   } else {
1240     mbs = a->mbs;
1241     ii  = a->i;
1242     z   = zarray;
1243   }
1244 
1245   if (!a->mult_work) {
1246     k    = PetscMax(A->rmap->n,A->cmap->n);
1247     ierr = PetscMalloc1(k+1,&a->mult_work);CHKERRQ(ierr);
1248   }
1249   work = a->mult_work;
1250   for (i=0; i<mbs; i++) {
1251     n           = ii[1] - ii[0]; ii++;
1252     ncols       = n*bs;
1253     workt       = work;
1254     for (j=0; j<n; j++) {
1255       xb = x + bs*(*idx++);
1256       for (k=0; k<bs; k++) workt[k] = xb[k];
1257       workt += bs;
1258     }
1259     if (usecprow) z = zarray + bs*ridx[i];
1260     PetscKernel_w_gets_Ar_times_v(bs,ncols,work,v,z);
1261     /* BLASgemv_("N",&bs,&ncols,&_DOne,v,&bs,work,&_One,&_DZero,z,&_One); */
1262     v += n*bs2;
1263     if (!usecprow) z += bs;
1264   }
1265   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
1266   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
1267   ierr = PetscLogFlops(2.0*a->nz*bs2 - bs*a->nonzerorowcnt);CHKERRQ(ierr);
1268   PetscFunctionReturn(0);
1269 }
1270 
1271 PetscErrorCode MatMultAdd_SeqBAIJ_1(Mat A,Vec xx,Vec yy,Vec zz)
1272 {
1273   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
1274   const PetscScalar *x;
1275   PetscScalar       *y,*z,sum;
1276   const MatScalar   *v;
1277   PetscErrorCode    ierr;
1278   PetscInt          mbs=a->mbs,i,n,*ridx=NULL;
1279   const PetscInt    *idx,*ii;
1280   PetscBool         usecprow=a->compressedrow.use;
1281 
1282   PetscFunctionBegin;
1283   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
1284   ierr = VecGetArrayPair(yy,zz,&y,&z);CHKERRQ(ierr);
1285 
1286   idx = a->j;
1287   v   = a->a;
1288   if (usecprow) {
1289     if (zz != yy) {
1290       ierr = PetscMemcpy(z,y,mbs*sizeof(PetscScalar));CHKERRQ(ierr);
1291     }
1292     mbs  = a->compressedrow.nrows;
1293     ii   = a->compressedrow.i;
1294     ridx = a->compressedrow.rindex;
1295   } else {
1296     ii = a->i;
1297   }
1298 
1299   for (i=0; i<mbs; i++) {
1300     n = ii[1] - ii[0];
1301     ii++;
1302     if (!usecprow) {
1303       sum         = y[i];
1304     } else {
1305       sum = y[ridx[i]];
1306     }
1307     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA); /* Indices for the next row (assumes same size as this one) */
1308     PetscPrefetchBlock(v+n,n,0,PETSC_PREFETCH_HINT_NTA);   /* Entries for the next row */
1309     PetscSparseDensePlusDot(sum,x,v,idx,n);
1310     v   += n;
1311     idx += n;
1312     if (usecprow) {
1313       z[ridx[i]] = sum;
1314     } else {
1315       z[i] = sum;
1316     }
1317   }
1318   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
1319   ierr = VecRestoreArrayPair(yy,zz,&y,&z);CHKERRQ(ierr);
1320   ierr = PetscLogFlops(2.0*a->nz);CHKERRQ(ierr);
1321   PetscFunctionReturn(0);
1322 }
1323 
1324 PetscErrorCode MatMultAdd_SeqBAIJ_2(Mat A,Vec xx,Vec yy,Vec zz)
1325 {
1326   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
1327   PetscScalar       *y = 0,*z = 0,sum1,sum2;
1328   const PetscScalar *x,*xb;
1329   PetscScalar       x1,x2,*yarray,*zarray;
1330   const MatScalar   *v;
1331   PetscErrorCode    ierr;
1332   PetscInt          mbs = a->mbs,i,n,j;
1333   const PetscInt    *idx,*ii,*ridx = NULL;
1334   PetscBool         usecprow = a->compressedrow.use;
1335 
1336   PetscFunctionBegin;
1337   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
1338   ierr = VecGetArrayPair(yy,zz,&yarray,&zarray);CHKERRQ(ierr);
1339 
1340   idx = a->j;
1341   v   = a->a;
1342   if (usecprow) {
1343     if (zz != yy) {
1344       ierr = PetscMemcpy(zarray,yarray,2*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
1345     }
1346     mbs  = a->compressedrow.nrows;
1347     ii   = a->compressedrow.i;
1348     ridx = a->compressedrow.rindex;
1349   } else {
1350     ii = a->i;
1351     y  = yarray;
1352     z  = zarray;
1353   }
1354 
1355   for (i=0; i<mbs; i++) {
1356     n = ii[1] - ii[0]; ii++;
1357     if (usecprow) {
1358       z = zarray + 2*ridx[i];
1359       y = yarray + 2*ridx[i];
1360     }
1361     sum1 = y[0]; sum2 = y[1];
1362     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);   /* Indices for the next row (assumes same size as this one) */
1363     PetscPrefetchBlock(v+4*n,4*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
1364     for (j=0; j<n; j++) {
1365       xb = x + 2*(*idx++);
1366       x1 = xb[0];
1367       x2 = xb[1];
1368 
1369       sum1 += v[0]*x1 + v[2]*x2;
1370       sum2 += v[1]*x1 + v[3]*x2;
1371       v    += 4;
1372     }
1373     z[0] = sum1; z[1] = sum2;
1374     if (!usecprow) {
1375       z += 2; y += 2;
1376     }
1377   }
1378   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
1379   ierr = VecRestoreArrayPair(yy,zz,&yarray,&zarray);CHKERRQ(ierr);
1380   ierr = PetscLogFlops(4.0*a->nz);CHKERRQ(ierr);
1381   PetscFunctionReturn(0);
1382 }
1383 
1384 PetscErrorCode MatMultAdd_SeqBAIJ_3(Mat A,Vec xx,Vec yy,Vec zz)
1385 {
1386   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
1387   PetscScalar       *y = 0,*z = 0,sum1,sum2,sum3,x1,x2,x3,*yarray,*zarray;
1388   const PetscScalar *x,*xb;
1389   const MatScalar   *v;
1390   PetscErrorCode    ierr;
1391   PetscInt          mbs = a->mbs,i,j,n;
1392   const PetscInt    *idx,*ii,*ridx = NULL;
1393   PetscBool         usecprow = a->compressedrow.use;
1394 
1395   PetscFunctionBegin;
1396   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
1397   ierr = VecGetArrayPair(yy,zz,&yarray,&zarray);CHKERRQ(ierr);
1398 
1399   idx = a->j;
1400   v   = a->a;
1401   if (usecprow) {
1402     if (zz != yy) {
1403       ierr = PetscMemcpy(zarray,yarray,3*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
1404     }
1405     mbs  = a->compressedrow.nrows;
1406     ii   = a->compressedrow.i;
1407     ridx = a->compressedrow.rindex;
1408   } else {
1409     ii = a->i;
1410     y  = yarray;
1411     z  = zarray;
1412   }
1413 
1414   for (i=0; i<mbs; i++) {
1415     n = ii[1] - ii[0]; ii++;
1416     if (usecprow) {
1417       z = zarray + 3*ridx[i];
1418       y = yarray + 3*ridx[i];
1419     }
1420     sum1 = y[0]; sum2 = y[1]; sum3 = y[2];
1421     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);   /* Indices for the next row (assumes same size as this one) */
1422     PetscPrefetchBlock(v+9*n,9*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
1423     for (j=0; j<n; j++) {
1424       xb    = x + 3*(*idx++); x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
1425       sum1 += v[0]*x1 + v[3]*x2 + v[6]*x3;
1426       sum2 += v[1]*x1 + v[4]*x2 + v[7]*x3;
1427       sum3 += v[2]*x1 + v[5]*x2 + v[8]*x3;
1428       v    += 9;
1429     }
1430     z[0] = sum1; z[1] = sum2; z[2] = sum3;
1431     if (!usecprow) {
1432       z += 3; y += 3;
1433     }
1434   }
1435   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
1436   ierr = VecRestoreArrayPair(yy,zz,&yarray,&zarray);CHKERRQ(ierr);
1437   ierr = PetscLogFlops(18.0*a->nz);CHKERRQ(ierr);
1438   PetscFunctionReturn(0);
1439 }
1440 
1441 PetscErrorCode MatMultAdd_SeqBAIJ_4(Mat A,Vec xx,Vec yy,Vec zz)
1442 {
1443   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
1444   PetscScalar       *y = 0,*z = 0,sum1,sum2,sum3,sum4,x1,x2,x3,x4,*yarray,*zarray;
1445   const PetscScalar *x,*xb;
1446   const MatScalar   *v;
1447   PetscErrorCode    ierr;
1448   PetscInt          mbs = a->mbs,i,j,n;
1449   const PetscInt    *idx,*ii,*ridx=NULL;
1450   PetscBool         usecprow=a->compressedrow.use;
1451 
1452   PetscFunctionBegin;
1453   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
1454   ierr = VecGetArrayPair(yy,zz,&yarray,&zarray);CHKERRQ(ierr);
1455 
1456   idx = a->j;
1457   v   = a->a;
1458   if (usecprow) {
1459     if (zz != yy) {
1460       ierr = PetscMemcpy(zarray,yarray,4*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
1461     }
1462     mbs  = a->compressedrow.nrows;
1463     ii   = a->compressedrow.i;
1464     ridx = a->compressedrow.rindex;
1465   } else {
1466     ii = a->i;
1467     y  = yarray;
1468     z  = zarray;
1469   }
1470 
1471   for (i=0; i<mbs; i++) {
1472     n = ii[1] - ii[0]; ii++;
1473     if (usecprow) {
1474       z = zarray + 4*ridx[i];
1475       y = yarray + 4*ridx[i];
1476     }
1477     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3];
1478     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);     /* Indices for the next row (assumes same size as this one) */
1479     PetscPrefetchBlock(v+16*n,16*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
1480     for (j=0; j<n; j++) {
1481       xb    = x + 4*(*idx++);
1482       x1    = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
1483       sum1 += v[0]*x1 + v[4]*x2 + v[8]*x3   + v[12]*x4;
1484       sum2 += v[1]*x1 + v[5]*x2 + v[9]*x3   + v[13]*x4;
1485       sum3 += v[2]*x1 + v[6]*x2 + v[10]*x3  + v[14]*x4;
1486       sum4 += v[3]*x1 + v[7]*x2 + v[11]*x3  + v[15]*x4;
1487       v    += 16;
1488     }
1489     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4;
1490     if (!usecprow) {
1491       z += 4; y += 4;
1492     }
1493   }
1494   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
1495   ierr = VecRestoreArrayPair(yy,zz,&yarray,&zarray);CHKERRQ(ierr);
1496   ierr = PetscLogFlops(32.0*a->nz);CHKERRQ(ierr);
1497   PetscFunctionReturn(0);
1498 }
1499 
1500 PetscErrorCode MatMultAdd_SeqBAIJ_5(Mat A,Vec xx,Vec yy,Vec zz)
1501 {
1502   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
1503   PetscScalar       *y = 0,*z = 0,sum1,sum2,sum3,sum4,sum5,x1,x2,x3,x4,x5;
1504   const PetscScalar *x,*xb;
1505   PetscScalar       *yarray,*zarray;
1506   const MatScalar   *v;
1507   PetscErrorCode    ierr;
1508   PetscInt          mbs = a->mbs,i,j,n;
1509   const PetscInt    *idx,*ii,*ridx = NULL;
1510   PetscBool         usecprow=a->compressedrow.use;
1511 
1512   PetscFunctionBegin;
1513   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
1514   ierr = VecGetArrayPair(yy,zz,&yarray,&zarray);CHKERRQ(ierr);
1515 
1516   idx = a->j;
1517   v   = a->a;
1518   if (usecprow) {
1519     if (zz != yy) {
1520       ierr = PetscMemcpy(zarray,yarray,5*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
1521     }
1522     mbs  = a->compressedrow.nrows;
1523     ii   = a->compressedrow.i;
1524     ridx = a->compressedrow.rindex;
1525   } else {
1526     ii = a->i;
1527     y  = yarray;
1528     z  = zarray;
1529   }
1530 
1531   for (i=0; i<mbs; i++) {
1532     n = ii[1] - ii[0]; ii++;
1533     if (usecprow) {
1534       z = zarray + 5*ridx[i];
1535       y = yarray + 5*ridx[i];
1536     }
1537     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4];
1538     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);     /* Indices for the next row (assumes same size as this one) */
1539     PetscPrefetchBlock(v+25*n,25*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
1540     for (j=0; j<n; j++) {
1541       xb    = x + 5*(*idx++);
1542       x1    = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4];
1543       sum1 += v[0]*x1 + v[5]*x2 + v[10]*x3  + v[15]*x4 + v[20]*x5;
1544       sum2 += v[1]*x1 + v[6]*x2 + v[11]*x3  + v[16]*x4 + v[21]*x5;
1545       sum3 += v[2]*x1 + v[7]*x2 + v[12]*x3  + v[17]*x4 + v[22]*x5;
1546       sum4 += v[3]*x1 + v[8]*x2 + v[13]*x3  + v[18]*x4 + v[23]*x5;
1547       sum5 += v[4]*x1 + v[9]*x2 + v[14]*x3  + v[19]*x4 + v[24]*x5;
1548       v    += 25;
1549     }
1550     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5;
1551     if (!usecprow) {
1552       z += 5; y += 5;
1553     }
1554   }
1555   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
1556   ierr = VecRestoreArrayPair(yy,zz,&yarray,&zarray);CHKERRQ(ierr);
1557   ierr = PetscLogFlops(50.0*a->nz);CHKERRQ(ierr);
1558   PetscFunctionReturn(0);
1559 }
1560 PetscErrorCode MatMultAdd_SeqBAIJ_6(Mat A,Vec xx,Vec yy,Vec zz)
1561 {
1562   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
1563   PetscScalar       *y = 0,*z = 0,sum1,sum2,sum3,sum4,sum5,sum6;
1564   const PetscScalar *x,*xb;
1565   PetscScalar       x1,x2,x3,x4,x5,x6,*yarray,*zarray;
1566   const MatScalar   *v;
1567   PetscErrorCode    ierr;
1568   PetscInt          mbs = a->mbs,i,j,n;
1569   const PetscInt    *idx,*ii,*ridx=NULL;
1570   PetscBool         usecprow=a->compressedrow.use;
1571 
1572   PetscFunctionBegin;
1573   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
1574   ierr = VecGetArrayPair(yy,zz,&yarray,&zarray);CHKERRQ(ierr);
1575 
1576   idx = a->j;
1577   v   = a->a;
1578   if (usecprow) {
1579     if (zz != yy) {
1580       ierr = PetscMemcpy(zarray,yarray,6*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
1581     }
1582     mbs  = a->compressedrow.nrows;
1583     ii   = a->compressedrow.i;
1584     ridx = a->compressedrow.rindex;
1585   } else {
1586     ii = a->i;
1587     y  = yarray;
1588     z  = zarray;
1589   }
1590 
1591   for (i=0; i<mbs; i++) {
1592     n = ii[1] - ii[0]; ii++;
1593     if (usecprow) {
1594       z = zarray + 6*ridx[i];
1595       y = yarray + 6*ridx[i];
1596     }
1597     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4]; sum6 = y[5];
1598     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);     /* Indices for the next row (assumes same size as this one) */
1599     PetscPrefetchBlock(v+36*n,36*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
1600     for (j=0; j<n; j++) {
1601       xb    = x + 6*(*idx++);
1602       x1    = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5];
1603       sum1 += v[0]*x1 + v[6]*x2  + v[12]*x3  + v[18]*x4 + v[24]*x5 + v[30]*x6;
1604       sum2 += v[1]*x1 + v[7]*x2  + v[13]*x3  + v[19]*x4 + v[25]*x5 + v[31]*x6;
1605       sum3 += v[2]*x1 + v[8]*x2  + v[14]*x3  + v[20]*x4 + v[26]*x5 + v[32]*x6;
1606       sum4 += v[3]*x1 + v[9]*x2  + v[15]*x3  + v[21]*x4 + v[27]*x5 + v[33]*x6;
1607       sum5 += v[4]*x1 + v[10]*x2 + v[16]*x3  + v[22]*x4 + v[28]*x5 + v[34]*x6;
1608       sum6 += v[5]*x1 + v[11]*x2 + v[17]*x3  + v[23]*x4 + v[29]*x5 + v[35]*x6;
1609       v    += 36;
1610     }
1611     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6;
1612     if (!usecprow) {
1613       z += 6; y += 6;
1614     }
1615   }
1616   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
1617   ierr = VecRestoreArrayPair(yy,zz,&yarray,&zarray);CHKERRQ(ierr);
1618   ierr = PetscLogFlops(72.0*a->nz);CHKERRQ(ierr);
1619   PetscFunctionReturn(0);
1620 }
1621 
1622 PetscErrorCode MatMultAdd_SeqBAIJ_7(Mat A,Vec xx,Vec yy,Vec zz)
1623 {
1624   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
1625   PetscScalar       *y = 0,*z = 0,sum1,sum2,sum3,sum4,sum5,sum6,sum7;
1626   const PetscScalar *x,*xb;
1627   PetscScalar       x1,x2,x3,x4,x5,x6,x7,*yarray,*zarray;
1628   const MatScalar   *v;
1629   PetscErrorCode    ierr;
1630   PetscInt          mbs = a->mbs,i,j,n;
1631   const PetscInt    *idx,*ii,*ridx = NULL;
1632   PetscBool         usecprow=a->compressedrow.use;
1633 
1634   PetscFunctionBegin;
1635   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
1636   ierr = VecGetArrayPair(yy,zz,&yarray,&zarray);CHKERRQ(ierr);
1637 
1638   idx = a->j;
1639   v   = a->a;
1640   if (usecprow) {
1641     if (zz != yy) {
1642       ierr = PetscMemcpy(zarray,yarray,7*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
1643     }
1644     mbs  = a->compressedrow.nrows;
1645     ii   = a->compressedrow.i;
1646     ridx = a->compressedrow.rindex;
1647   } else {
1648     ii = a->i;
1649     y  = yarray;
1650     z  = zarray;
1651   }
1652 
1653   for (i=0; i<mbs; i++) {
1654     n = ii[1] - ii[0]; ii++;
1655     if (usecprow) {
1656       z = zarray + 7*ridx[i];
1657       y = yarray + 7*ridx[i];
1658     }
1659     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4]; sum6 = y[5]; sum7 = y[6];
1660     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);     /* Indices for the next row (assumes same size as this one) */
1661     PetscPrefetchBlock(v+49*n,49*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
1662     for (j=0; j<n; j++) {
1663       xb    = x + 7*(*idx++);
1664       x1    = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6];
1665       sum1 += v[0]*x1 + v[7]*x2  + v[14]*x3  + v[21]*x4 + v[28]*x5 + v[35]*x6 + v[42]*x7;
1666       sum2 += v[1]*x1 + v[8]*x2  + v[15]*x3  + v[22]*x4 + v[29]*x5 + v[36]*x6 + v[43]*x7;
1667       sum3 += v[2]*x1 + v[9]*x2  + v[16]*x3  + v[23]*x4 + v[30]*x5 + v[37]*x6 + v[44]*x7;
1668       sum4 += v[3]*x1 + v[10]*x2 + v[17]*x3  + v[24]*x4 + v[31]*x5 + v[38]*x6 + v[45]*x7;
1669       sum5 += v[4]*x1 + v[11]*x2 + v[18]*x3  + v[25]*x4 + v[32]*x5 + v[39]*x6 + v[46]*x7;
1670       sum6 += v[5]*x1 + v[12]*x2 + v[19]*x3  + v[26]*x4 + v[33]*x5 + v[40]*x6 + v[47]*x7;
1671       sum7 += v[6]*x1 + v[13]*x2 + v[20]*x3  + v[27]*x4 + v[34]*x5 + v[41]*x6 + v[48]*x7;
1672       v    += 49;
1673     }
1674     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
1675     if (!usecprow) {
1676       z += 7; y += 7;
1677     }
1678   }
1679   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
1680   ierr = VecRestoreArrayPair(yy,zz,&yarray,&zarray);CHKERRQ(ierr);
1681   ierr = PetscLogFlops(98.0*a->nz);CHKERRQ(ierr);
1682   PetscFunctionReturn(0);
1683 }
1684 
1685 #if defined(PETSC_HAVE_IMMINTRIN_H) && defined(__AVX2__) && defined(PETSC_USE_REAL_DOUBLE) && !defined(PETSC_USE_COMPLEX)
1686 PetscErrorCode MatMultAdd_SeqBAIJ_9_AVX2(Mat A,Vec xx,Vec yy,Vec zz)
1687 {
1688   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
1689   PetscScalar    *z = 0,*work,*workt,*zarray;
1690   const PetscScalar *x,*xb;
1691   const MatScalar   *v;
1692   PetscErrorCode ierr;
1693   PetscInt       mbs,i,bs=A->rmap->bs,j,n,bs2=a->bs2;
1694   PetscInt       ncols,k;
1695   PetscBool      usecprow=a->compressedrow.use;
1696   const PetscInt *idx,*ii,*ridx=NULL;
1697 
1698   __m256d a0,a1,a2,a3,a4,a5;
1699   __m256d w0,w1,w2,w3,w4,w5;
1700   __m256d z0,z1,z2;
1701   __m256i mask1 = _mm256_set_epi64x(0LL, 0LL, 0LL, 1LL<<63);
1702 
1703   PetscFunctionBegin;
1704   ierr = VecCopy(yy,zz);CHKERRQ(ierr);
1705   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
1706   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
1707 
1708   idx = a->j;
1709   v   = a->a;
1710   if (usecprow) {
1711     mbs  = a->compressedrow.nrows;
1712     ii   = a->compressedrow.i;
1713     ridx = a->compressedrow.rindex;
1714   } else {
1715     mbs = a->mbs;
1716     ii  = a->i;
1717     z   = zarray;
1718   }
1719 
1720   if (!a->mult_work) {
1721     k    = PetscMax(A->rmap->n,A->cmap->n);
1722     ierr = PetscMalloc1(k+1,&a->mult_work);CHKERRQ(ierr);
1723   }
1724 
1725   work = a->mult_work;
1726   for (i=0; i<mbs; i++) {
1727     n           = ii[1] - ii[0]; ii++;
1728     ncols       = n*bs;
1729     workt       = work;
1730     for (j=0; j<n; j++) {
1731       xb = x + bs*(*idx++);
1732       for (k=0; k<bs; k++) workt[k] = xb[k];
1733       workt += bs;
1734     }
1735     if (usecprow) z = zarray + bs*ridx[i];
1736 
1737     z0 = _mm256_loadu_pd(&z[ 0]); z1 = _mm256_loadu_pd(&z[ 4]); z2 = _mm256_set1_pd(z[ 8]);
1738 
1739     for (j=0; j<n; j++) {
1740       // first column of a
1741       w0 = _mm256_set1_pd(work[j*9  ]);
1742       a0 = _mm256_loadu_pd(&v[j*81  ]); z0 = _mm256_fmadd_pd(a0,w0,z0);
1743       a1 = _mm256_loadu_pd(&v[j*81+4]); z1 = _mm256_fmadd_pd(a1,w0,z1);
1744       a2 = _mm256_loadu_pd(&v[j*81+8]); z2 = _mm256_fmadd_pd(a2,w0,z2);
1745 
1746       // second column of a
1747       w1 = _mm256_set1_pd(work[j*9+ 1]);
1748       a0 = _mm256_loadu_pd(&v[j*81+ 9]); z0 = _mm256_fmadd_pd(a0,w1,z0);
1749       a1 = _mm256_loadu_pd(&v[j*81+13]); z1 = _mm256_fmadd_pd(a1,w1,z1);
1750       a2 = _mm256_loadu_pd(&v[j*81+17]); z2 = _mm256_fmadd_pd(a2,w1,z2);
1751 
1752       // third column of a
1753       w2 = _mm256_set1_pd(work[j*9+ 2]);
1754       a3 = _mm256_loadu_pd(&v[j*81+18]); z0 = _mm256_fmadd_pd(a3,w2,z0);
1755       a4 = _mm256_loadu_pd(&v[j*81+22]); z1 = _mm256_fmadd_pd(a4,w2,z1);
1756       a5 = _mm256_loadu_pd(&v[j*81+26]); z2 = _mm256_fmadd_pd(a5,w2,z2);
1757 
1758       // fourth column of a
1759       w3 = _mm256_set1_pd(work[j*9+ 3]);
1760       a0 = _mm256_loadu_pd(&v[j*81+27]); z0 = _mm256_fmadd_pd(a0,w3,z0);
1761       a1 = _mm256_loadu_pd(&v[j*81+31]); z1 = _mm256_fmadd_pd(a1,w3,z1);
1762       a2 = _mm256_loadu_pd(&v[j*81+35]); z2 = _mm256_fmadd_pd(a2,w3,z2);
1763 
1764       // fifth column of a
1765       w0 = _mm256_set1_pd(work[j*9+ 4]);
1766       a3 = _mm256_loadu_pd(&v[j*81+36]); z0 = _mm256_fmadd_pd(a3,w0,z0);
1767       a4 = _mm256_loadu_pd(&v[j*81+40]); z1 = _mm256_fmadd_pd(a4,w0,z1);
1768       a5 = _mm256_loadu_pd(&v[j*81+44]); z2 = _mm256_fmadd_pd(a5,w0,z2);
1769 
1770       // sixth column of a
1771       w1 = _mm256_set1_pd(work[j*9+ 5]);
1772       a0 = _mm256_loadu_pd(&v[j*81+45]); z0 = _mm256_fmadd_pd(a0,w1,z0);
1773       a1 = _mm256_loadu_pd(&v[j*81+49]); z1 = _mm256_fmadd_pd(a1,w1,z1);
1774       a2 = _mm256_loadu_pd(&v[j*81+53]); z2 = _mm256_fmadd_pd(a2,w1,z2);
1775 
1776       // seventh column of a
1777       w2 = _mm256_set1_pd(work[j*9+ 6]);
1778       a0 = _mm256_loadu_pd(&v[j*81+54]); z0 = _mm256_fmadd_pd(a0,w2,z0);
1779       a1 = _mm256_loadu_pd(&v[j*81+58]); z1 = _mm256_fmadd_pd(a1,w2,z1);
1780       a2 = _mm256_loadu_pd(&v[j*81+62]); z2 = _mm256_fmadd_pd(a2,w2,z2);
1781 
1782       // eigth column of a
1783       w3 = _mm256_set1_pd(work[j*9+ 7]);
1784       a3 = _mm256_loadu_pd(&v[j*81+63]); z0 = _mm256_fmadd_pd(a3,w3,z0);
1785       a4 = _mm256_loadu_pd(&v[j*81+67]); z1 = _mm256_fmadd_pd(a4,w3,z1);
1786       a5 = _mm256_loadu_pd(&v[j*81+71]); z2 = _mm256_fmadd_pd(a5,w3,z2);
1787 
1788       // ninth column of a
1789       w0 = _mm256_set1_pd(work[j*9+ 8]);
1790       a0 = _mm256_loadu_pd(&v[j*81+72]); z0 = _mm256_fmadd_pd(a0,w0,z0);
1791       a1 = _mm256_loadu_pd(&v[j*81+76]); z1 = _mm256_fmadd_pd(a1,w0,z1);
1792       a2 = _mm256_maskload_pd(&v[j*81+80],mask1); z2 = _mm256_fmadd_pd(a2,w0,z2);
1793     }
1794 
1795     _mm256_storeu_pd(&z[ 0], z0); _mm256_storeu_pd(&z[ 4], z1); _mm256_maskstore_pd(&z[8], mask1, z2);
1796 
1797     v += n*bs2;
1798     if (!usecprow) z += bs;
1799   }
1800   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
1801   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
1802   ierr = PetscLogFlops(162.0*a->nz);CHKERRQ(ierr);
1803   PetscFunctionReturn(0);
1804 }
1805 #endif
1806 
1807 PetscErrorCode MatMultAdd_SeqBAIJ_11(Mat A,Vec xx,Vec yy,Vec zz)
1808 {
1809   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
1810   PetscScalar       *y = 0,*z = 0,sum1,sum2,sum3,sum4,sum5,sum6,sum7,sum8,sum9,sum10,sum11;
1811   const PetscScalar *x,*xb;
1812   PetscScalar       x1,x2,x3,x4,x5,x6,x7,x8,x9,x10,x11,*yarray,*zarray;
1813   const MatScalar   *v;
1814   PetscErrorCode    ierr;
1815   PetscInt          mbs = a->mbs,i,j,n;
1816   const PetscInt    *idx,*ii,*ridx = NULL;
1817   PetscBool         usecprow=a->compressedrow.use;
1818 
1819   PetscFunctionBegin;
1820   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
1821   ierr = VecGetArrayPair(yy,zz,&yarray,&zarray);CHKERRQ(ierr);
1822 
1823   idx = a->j;
1824   v   = a->a;
1825   if (usecprow) {
1826     if (zz != yy) {
1827       ierr = PetscMemcpy(zarray,yarray,7*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
1828     }
1829     mbs  = a->compressedrow.nrows;
1830     ii   = a->compressedrow.i;
1831     ridx = a->compressedrow.rindex;
1832   } else {
1833     ii = a->i;
1834     y  = yarray;
1835     z  = zarray;
1836   }
1837 
1838   for (i=0; i<mbs; i++) {
1839     n = ii[1] - ii[0]; ii++;
1840     if (usecprow) {
1841       z = zarray + 11*ridx[i];
1842       y = yarray + 11*ridx[i];
1843     }
1844     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4]; sum6 = y[5]; sum7 = y[6];
1845     sum8 = y[7]; sum9 = y[8]; sum10 = y[9]; sum11 = y[10];
1846     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);     /* Indices for the next row (assumes same size as this one) */
1847     PetscPrefetchBlock(v+121*n,121*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
1848     for (j=0; j<n; j++) {
1849       xb    = x + 11*(*idx++);
1850       x1    = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6];x8 = xb[7]; x9 = xb[8]; x10 = xb[9]; x11 = xb[10];
1851       sum1 += v[  0]*x1 + v[  11]*x2  + v[2*11]*x3  + v[3*11]*x4 + v[4*11]*x5 + v[5*11]*x6 + v[6*11]*x7+ v[7*11]*x8 + v[8*11]*x9  + v[9*11]*x10  + v[10*11]*x11;
1852       sum2 += v[1+0]*x1 + v[1+11]*x2  + v[1+2*11]*x3  + v[1+3*11]*x4 + v[1+4*11]*x5 + v[1+5*11]*x6 + v[1+6*11]*x7+ v[1+7*11]*x8 + v[1+8*11]*x9  + v[1+9*11]*x10  + v[1+10*11]*x11;
1853       sum3 += v[2+0]*x1 + v[2+11]*x2  + v[2+2*11]*x3  + v[2+3*11]*x4 + v[2+4*11]*x5 + v[2+5*11]*x6 + v[2+6*11]*x7+ v[2+7*11]*x8 + v[2+8*11]*x9  + v[2+9*11]*x10  + v[2+10*11]*x11;
1854       sum4 += v[3+0]*x1 + v[3+11]*x2  + v[3+2*11]*x3  + v[3+3*11]*x4 + v[3+4*11]*x5 + v[3+5*11]*x6 + v[3+6*11]*x7+ v[3+7*11]*x8 + v[3+8*11]*x9  + v[3+9*11]*x10  + v[3+10*11]*x11;
1855       sum5 += v[4+0]*x1 + v[4+11]*x2  + v[4+2*11]*x3  + v[4+3*11]*x4 + v[4+4*11]*x5 + v[4+5*11]*x6 + v[4+6*11]*x7+ v[4+7*11]*x8 + v[4+8*11]*x9  + v[4+9*11]*x10  + v[4+10*11]*x11;
1856       sum6 += v[5+0]*x1 + v[5+11]*x2  + v[5+2*11]*x3  + v[5+3*11]*x4 + v[5+4*11]*x5 + v[5+5*11]*x6 + v[5+6*11]*x7+ v[5+7*11]*x8 + v[5+8*11]*x9  + v[5+9*11]*x10  + v[5+10*11]*x11;
1857       sum7 += v[6+0]*x1 + v[6+11]*x2  + v[6+2*11]*x3  + v[6+3*11]*x4 + v[6+4*11]*x5 + v[6+5*11]*x6 + v[6+6*11]*x7+ v[6+7*11]*x8 + v[6+8*11]*x9  + v[6+9*11]*x10  + v[6+10*11]*x11;
1858       sum8 += v[7+0]*x1 + v[7+11]*x2  + v[7+2*11]*x3  + v[7+3*11]*x4 + v[7+4*11]*x5 + v[7+5*11]*x6 + v[7+6*11]*x7+ v[7+7*11]*x8 + v[7+8*11]*x9  + v[7+9*11]*x10  + v[7+10*11]*x11;
1859       sum9 += v[8+0]*x1 + v[8+11]*x2  + v[8+2*11]*x3  + v[8+3*11]*x4 + v[8+4*11]*x5 + v[8+5*11]*x6 + v[8+6*11]*x7+ v[8+7*11]*x8 + v[8+8*11]*x9  + v[8+9*11]*x10  + v[8+10*11]*x11;
1860       sum10 += v[9+0]*x1 + v[9+11]*x2  + v[9+2*11]*x3  + v[9+3*11]*x4 + v[9+4*11]*x5 + v[9+5*11]*x6 + v[9+6*11]*x7+ v[9+7*11]*x8 + v[9+8*11]*x9  + v[9+9*11]*x10  + v[9+10*11]*x11;
1861       sum11 += v[10+0]*x1 + v[10+11]*x2  + v[10+2*11]*x3  + v[10+3*11]*x4 + v[10+4*11]*x5 + v[10+5*11]*x6 + v[10+6*11]*x7+ v[10+7*11]*x8 + v[10+8*11]*x9  + v[10+9*11]*x10  + v[10+10*11]*x11;
1862       v    += 121;
1863     }
1864     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
1865     z[7] = sum6; z[8] = sum7; z[9] = sum8; z[10] = sum9; z[11] = sum10;
1866     if (!usecprow) {
1867       z += 11; y += 11;
1868     }
1869   }
1870   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
1871   ierr = VecRestoreArrayPair(yy,zz,&yarray,&zarray);CHKERRQ(ierr);
1872   ierr = PetscLogFlops(242.0*a->nz);CHKERRQ(ierr);
1873   PetscFunctionReturn(0);
1874 }
1875 
1876 PetscErrorCode MatMultAdd_SeqBAIJ_N(Mat A,Vec xx,Vec yy,Vec zz)
1877 {
1878   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
1879   PetscScalar       *z = 0,*work,*workt,*zarray;
1880   const PetscScalar *x,*xb;
1881   const MatScalar   *v;
1882   PetscErrorCode    ierr;
1883   PetscInt          mbs,i,bs=A->rmap->bs,j,n,bs2=a->bs2;
1884   PetscInt          ncols,k;
1885   const PetscInt    *ridx = NULL,*idx,*ii;
1886   PetscBool         usecprow = a->compressedrow.use;
1887 
1888   PetscFunctionBegin;
1889   ierr = VecCopy(yy,zz);CHKERRQ(ierr);
1890   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
1891   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
1892 
1893   idx = a->j;
1894   v   = a->a;
1895   if (usecprow) {
1896     mbs  = a->compressedrow.nrows;
1897     ii   = a->compressedrow.i;
1898     ridx = a->compressedrow.rindex;
1899   } else {
1900     mbs = a->mbs;
1901     ii  = a->i;
1902     z   = zarray;
1903   }
1904 
1905   if (!a->mult_work) {
1906     k    = PetscMax(A->rmap->n,A->cmap->n);
1907     ierr = PetscMalloc1(k+1,&a->mult_work);CHKERRQ(ierr);
1908   }
1909   work = a->mult_work;
1910   for (i=0; i<mbs; i++) {
1911     n     = ii[1] - ii[0]; ii++;
1912     ncols = n*bs;
1913     workt = work;
1914     for (j=0; j<n; j++) {
1915       xb = x + bs*(*idx++);
1916       for (k=0; k<bs; k++) workt[k] = xb[k];
1917       workt += bs;
1918     }
1919     if (usecprow) z = zarray + bs*ridx[i];
1920     PetscKernel_w_gets_w_plus_Ar_times_v(bs,ncols,work,v,z);
1921     /* BLASgemv_("N",&bs,&ncols,&_DOne,v,&bs,work,&_One,&_DOne,z,&_One); */
1922     v += n*bs2;
1923     if (!usecprow) z += bs;
1924   }
1925   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
1926   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
1927   ierr = PetscLogFlops(2.0*a->nz*bs2);CHKERRQ(ierr);
1928   PetscFunctionReturn(0);
1929 }
1930 
1931 PetscErrorCode MatMultHermitianTranspose_SeqBAIJ(Mat A,Vec xx,Vec zz)
1932 {
1933   PetscScalar    zero = 0.0;
1934   PetscErrorCode ierr;
1935 
1936   PetscFunctionBegin;
1937   ierr = VecSet(zz,zero);CHKERRQ(ierr);
1938   ierr = MatMultHermitianTransposeAdd_SeqBAIJ(A,xx,zz,zz);CHKERRQ(ierr);
1939   PetscFunctionReturn(0);
1940 }
1941 
1942 PetscErrorCode MatMultTranspose_SeqBAIJ(Mat A,Vec xx,Vec zz)
1943 {
1944   PetscScalar    zero = 0.0;
1945   PetscErrorCode ierr;
1946 
1947   PetscFunctionBegin;
1948   ierr = VecSet(zz,zero);CHKERRQ(ierr);
1949   ierr = MatMultTransposeAdd_SeqBAIJ(A,xx,zz,zz);CHKERRQ(ierr);
1950   PetscFunctionReturn(0);
1951 }
1952 
1953 PetscErrorCode MatMultHermitianTransposeAdd_SeqBAIJ(Mat A,Vec xx,Vec yy,Vec zz)
1954 {
1955   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
1956   PetscScalar       *z,x1,x2,x3,x4,x5;
1957   const PetscScalar *x,*xb = NULL;
1958   const MatScalar   *v;
1959   PetscErrorCode    ierr;
1960   PetscInt          mbs,i,rval,bs=A->rmap->bs,j,n;
1961   const PetscInt    *idx,*ii,*ib,*ridx = NULL;
1962   Mat_CompressedRow cprow = a->compressedrow;
1963   PetscBool         usecprow = cprow.use;
1964 
1965   PetscFunctionBegin;
1966   if (yy != zz) { ierr = VecCopy(yy,zz);CHKERRQ(ierr); }
1967   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
1968   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
1969 
1970   idx = a->j;
1971   v   = a->a;
1972   if (usecprow) {
1973     mbs  = cprow.nrows;
1974     ii   = cprow.i;
1975     ridx = cprow.rindex;
1976   } else {
1977     mbs=a->mbs;
1978     ii = a->i;
1979     xb = x;
1980   }
1981 
1982   switch (bs) {
1983   case 1:
1984     for (i=0; i<mbs; i++) {
1985       if (usecprow) xb = x + ridx[i];
1986       x1 = xb[0];
1987       ib = idx + ii[0];
1988       n  = ii[1] - ii[0]; ii++;
1989       for (j=0; j<n; j++) {
1990         rval     = ib[j];
1991         z[rval] += PetscConj(*v) * x1;
1992         v++;
1993       }
1994       if (!usecprow) xb++;
1995     }
1996     break;
1997   case 2:
1998     for (i=0; i<mbs; i++) {
1999       if (usecprow) xb = x + 2*ridx[i];
2000       x1 = xb[0]; x2 = xb[1];
2001       ib = idx + ii[0];
2002       n  = ii[1] - ii[0]; ii++;
2003       for (j=0; j<n; j++) {
2004         rval       = ib[j]*2;
2005         z[rval++] += PetscConj(v[0])*x1 + PetscConj(v[1])*x2;
2006         z[rval++] += PetscConj(v[2])*x1 + PetscConj(v[3])*x2;
2007         v         += 4;
2008       }
2009       if (!usecprow) xb += 2;
2010     }
2011     break;
2012   case 3:
2013     for (i=0; i<mbs; i++) {
2014       if (usecprow) xb = x + 3*ridx[i];
2015       x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
2016       ib = idx + ii[0];
2017       n  = ii[1] - ii[0]; ii++;
2018       for (j=0; j<n; j++) {
2019         rval       = ib[j]*3;
2020         z[rval++] += PetscConj(v[0])*x1 + PetscConj(v[1])*x2 + PetscConj(v[2])*x3;
2021         z[rval++] += PetscConj(v[3])*x1 + PetscConj(v[4])*x2 + PetscConj(v[5])*x3;
2022         z[rval++] += PetscConj(v[6])*x1 + PetscConj(v[7])*x2 + PetscConj(v[8])*x3;
2023         v         += 9;
2024       }
2025       if (!usecprow) xb += 3;
2026     }
2027     break;
2028   case 4:
2029     for (i=0; i<mbs; i++) {
2030       if (usecprow) xb = x + 4*ridx[i];
2031       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
2032       ib = idx + ii[0];
2033       n  = ii[1] - ii[0]; ii++;
2034       for (j=0; j<n; j++) {
2035         rval       = ib[j]*4;
2036         z[rval++] +=  PetscConj(v[0])*x1 + PetscConj(v[1])*x2  + PetscConj(v[2])*x3  + PetscConj(v[3])*x4;
2037         z[rval++] +=  PetscConj(v[4])*x1 + PetscConj(v[5])*x2  + PetscConj(v[6])*x3  + PetscConj(v[7])*x4;
2038         z[rval++] +=  PetscConj(v[8])*x1 + PetscConj(v[9])*x2  + PetscConj(v[10])*x3 + PetscConj(v[11])*x4;
2039         z[rval++] += PetscConj(v[12])*x1 + PetscConj(v[13])*x2 + PetscConj(v[14])*x3 + PetscConj(v[15])*x4;
2040         v         += 16;
2041       }
2042       if (!usecprow) xb += 4;
2043     }
2044     break;
2045   case 5:
2046     for (i=0; i<mbs; i++) {
2047       if (usecprow) xb = x + 5*ridx[i];
2048       x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
2049       x4 = xb[3]; x5 = xb[4];
2050       ib = idx + ii[0];
2051       n  = ii[1] - ii[0]; ii++;
2052       for (j=0; j<n; j++) {
2053         rval       = ib[j]*5;
2054         z[rval++] +=  PetscConj(v[0])*x1 +  PetscConj(v[1])*x2 +  PetscConj(v[2])*x3 +  PetscConj(v[3])*x4 +  PetscConj(v[4])*x5;
2055         z[rval++] +=  PetscConj(v[5])*x1 +  PetscConj(v[6])*x2 +  PetscConj(v[7])*x3 +  PetscConj(v[8])*x4 +  PetscConj(v[9])*x5;
2056         z[rval++] += PetscConj(v[10])*x1 + PetscConj(v[11])*x2 + PetscConj(v[12])*x3 + PetscConj(v[13])*x4 + PetscConj(v[14])*x5;
2057         z[rval++] += PetscConj(v[15])*x1 + PetscConj(v[16])*x2 + PetscConj(v[17])*x3 + PetscConj(v[18])*x4 + PetscConj(v[19])*x5;
2058         z[rval++] += PetscConj(v[20])*x1 + PetscConj(v[21])*x2 + PetscConj(v[22])*x3 + PetscConj(v[23])*x4 + PetscConj(v[24])*x5;
2059         v         += 25;
2060       }
2061       if (!usecprow) xb += 5;
2062     }
2063     break;
2064   default: /* block sizes larger than 5 by 5 are handled by BLAS */
2065     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"block size larger than 5 is not supported yet");
2066 #if 0
2067     {
2068       PetscInt          ncols,k,bs2=a->bs2;
2069       PetscScalar       *work,*workt,zb;
2070       const PetscScalar *xtmp;
2071       if (!a->mult_work) {
2072         k    = PetscMax(A->rmap->n,A->cmap->n);
2073         ierr = PetscMalloc1(k+1,&a->mult_work);CHKERRQ(ierr);
2074       }
2075       work = a->mult_work;
2076       xtmp = x;
2077       for (i=0; i<mbs; i++) {
2078         n     = ii[1] - ii[0]; ii++;
2079         ncols = n*bs;
2080         ierr  = PetscMemzero(work,ncols*sizeof(PetscScalar));CHKERRQ(ierr);
2081         if (usecprow) xtmp = x + bs*ridx[i];
2082         PetscKernel_w_gets_w_plus_trans_Ar_times_v(bs,ncols,xtmp,v,work);
2083         /* BLASgemv_("T",&bs,&ncols,&_DOne,v,&bs,xtmp,&_One,&_DOne,work,&_One); */
2084         v += n*bs2;
2085         if (!usecprow) xtmp += bs;
2086         workt = work;
2087         for (j=0; j<n; j++) {
2088           zb = z + bs*(*idx++);
2089           for (k=0; k<bs; k++) zb[k] += workt[k] ;
2090           workt += bs;
2091         }
2092       }
2093     }
2094 #endif
2095   }
2096   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
2097   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
2098   ierr = PetscLogFlops(2.0*a->nz*a->bs2);CHKERRQ(ierr);
2099   PetscFunctionReturn(0);
2100 }
2101 
2102 PetscErrorCode MatMultTransposeAdd_SeqBAIJ(Mat A,Vec xx,Vec yy,Vec zz)
2103 {
2104   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
2105   PetscScalar       *zb,*z,x1,x2,x3,x4,x5;
2106   const PetscScalar *x,*xb = 0;
2107   const MatScalar   *v;
2108   PetscErrorCode    ierr;
2109   PetscInt          mbs,i,rval,bs=A->rmap->bs,j,n,bs2=a->bs2;
2110   const PetscInt    *idx,*ii,*ib,*ridx = NULL;
2111   Mat_CompressedRow cprow   = a->compressedrow;
2112   PetscBool         usecprow=cprow.use;
2113 
2114   PetscFunctionBegin;
2115   if (yy != zz) { ierr = VecCopy(yy,zz);CHKERRQ(ierr); }
2116   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
2117   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
2118 
2119   idx = a->j;
2120   v   = a->a;
2121   if (usecprow) {
2122     mbs  = cprow.nrows;
2123     ii   = cprow.i;
2124     ridx = cprow.rindex;
2125   } else {
2126     mbs=a->mbs;
2127     ii = a->i;
2128     xb = x;
2129   }
2130 
2131   switch (bs) {
2132   case 1:
2133     for (i=0; i<mbs; i++) {
2134       if (usecprow) xb = x + ridx[i];
2135       x1 = xb[0];
2136       ib = idx + ii[0];
2137       n  = ii[1] - ii[0]; ii++;
2138       for (j=0; j<n; j++) {
2139         rval     = ib[j];
2140         z[rval] += *v * x1;
2141         v++;
2142       }
2143       if (!usecprow) xb++;
2144     }
2145     break;
2146   case 2:
2147     for (i=0; i<mbs; i++) {
2148       if (usecprow) xb = x + 2*ridx[i];
2149       x1 = xb[0]; x2 = xb[1];
2150       ib = idx + ii[0];
2151       n  = ii[1] - ii[0]; ii++;
2152       for (j=0; j<n; j++) {
2153         rval       = ib[j]*2;
2154         z[rval++] += v[0]*x1 + v[1]*x2;
2155         z[rval++] += v[2]*x1 + v[3]*x2;
2156         v         += 4;
2157       }
2158       if (!usecprow) xb += 2;
2159     }
2160     break;
2161   case 3:
2162     for (i=0; i<mbs; i++) {
2163       if (usecprow) xb = x + 3*ridx[i];
2164       x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
2165       ib = idx + ii[0];
2166       n  = ii[1] - ii[0]; ii++;
2167       for (j=0; j<n; j++) {
2168         rval       = ib[j]*3;
2169         z[rval++] += v[0]*x1 + v[1]*x2 + v[2]*x3;
2170         z[rval++] += v[3]*x1 + v[4]*x2 + v[5]*x3;
2171         z[rval++] += v[6]*x1 + v[7]*x2 + v[8]*x3;
2172         v         += 9;
2173       }
2174       if (!usecprow) xb += 3;
2175     }
2176     break;
2177   case 4:
2178     for (i=0; i<mbs; i++) {
2179       if (usecprow) xb = x + 4*ridx[i];
2180       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
2181       ib = idx + ii[0];
2182       n  = ii[1] - ii[0]; ii++;
2183       for (j=0; j<n; j++) {
2184         rval       = ib[j]*4;
2185         z[rval++] +=  v[0]*x1 +  v[1]*x2 +  v[2]*x3 +  v[3]*x4;
2186         z[rval++] +=  v[4]*x1 +  v[5]*x2 +  v[6]*x3 +  v[7]*x4;
2187         z[rval++] +=  v[8]*x1 +  v[9]*x2 + v[10]*x3 + v[11]*x4;
2188         z[rval++] += v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4;
2189         v         += 16;
2190       }
2191       if (!usecprow) xb += 4;
2192     }
2193     break;
2194   case 5:
2195     for (i=0; i<mbs; i++) {
2196       if (usecprow) xb = x + 5*ridx[i];
2197       x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
2198       x4 = xb[3]; x5 = xb[4];
2199       ib = idx + ii[0];
2200       n  = ii[1] - ii[0]; ii++;
2201       for (j=0; j<n; j++) {
2202         rval       = ib[j]*5;
2203         z[rval++] +=  v[0]*x1 +  v[1]*x2 +  v[2]*x3 +  v[3]*x4 +  v[4]*x5;
2204         z[rval++] +=  v[5]*x1 +  v[6]*x2 +  v[7]*x3 +  v[8]*x4 +  v[9]*x5;
2205         z[rval++] += v[10]*x1 + v[11]*x2 + v[12]*x3 + v[13]*x4 + v[14]*x5;
2206         z[rval++] += v[15]*x1 + v[16]*x2 + v[17]*x3 + v[18]*x4 + v[19]*x5;
2207         z[rval++] += v[20]*x1 + v[21]*x2 + v[22]*x3 + v[23]*x4 + v[24]*x5;
2208         v         += 25;
2209       }
2210       if (!usecprow) xb += 5;
2211     }
2212     break;
2213   default: {      /* block sizes larger then 5 by 5 are handled by BLAS */
2214     PetscInt          ncols,k;
2215     PetscScalar       *work,*workt;
2216     const PetscScalar *xtmp;
2217     if (!a->mult_work) {
2218       k    = PetscMax(A->rmap->n,A->cmap->n);
2219       ierr = PetscMalloc1(k+1,&a->mult_work);CHKERRQ(ierr);
2220     }
2221     work = a->mult_work;
2222     xtmp = x;
2223     for (i=0; i<mbs; i++) {
2224       n     = ii[1] - ii[0]; ii++;
2225       ncols = n*bs;
2226       ierr  = PetscMemzero(work,ncols*sizeof(PetscScalar));CHKERRQ(ierr);
2227       if (usecprow) xtmp = x + bs*ridx[i];
2228       PetscKernel_w_gets_w_plus_trans_Ar_times_v(bs,ncols,xtmp,v,work);
2229       /* BLASgemv_("T",&bs,&ncols,&_DOne,v,&bs,xtmp,&_One,&_DOne,work,&_One); */
2230       v += n*bs2;
2231       if (!usecprow) xtmp += bs;
2232       workt = work;
2233       for (j=0; j<n; j++) {
2234         zb = z + bs*(*idx++);
2235         for (k=0; k<bs; k++) zb[k] += workt[k];
2236         workt += bs;
2237       }
2238     }
2239     }
2240   }
2241   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
2242   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
2243   ierr = PetscLogFlops(2.0*a->nz*a->bs2);CHKERRQ(ierr);
2244   PetscFunctionReturn(0);
2245 }
2246 
2247 PetscErrorCode MatScale_SeqBAIJ(Mat inA,PetscScalar alpha)
2248 {
2249   Mat_SeqBAIJ    *a      = (Mat_SeqBAIJ*)inA->data;
2250   PetscInt       totalnz = a->bs2*a->nz;
2251   PetscScalar    oalpha  = alpha;
2252   PetscErrorCode ierr;
2253   PetscBLASInt   one = 1,tnz;
2254 
2255   PetscFunctionBegin;
2256   ierr = PetscBLASIntCast(totalnz,&tnz);CHKERRQ(ierr);
2257   PetscStackCallBLAS("BLASscal",BLASscal_(&tnz,&oalpha,a->a,&one));
2258   ierr = PetscLogFlops(totalnz);CHKERRQ(ierr);
2259   PetscFunctionReturn(0);
2260 }
2261 
2262 PetscErrorCode MatNorm_SeqBAIJ(Mat A,NormType type,PetscReal *norm)
2263 {
2264   PetscErrorCode ierr;
2265   Mat_SeqBAIJ    *a  = (Mat_SeqBAIJ*)A->data;
2266   MatScalar      *v  = a->a;
2267   PetscReal      sum = 0.0;
2268   PetscInt       i,j,k,bs=A->rmap->bs,nz=a->nz,bs2=a->bs2,k1;
2269 
2270   PetscFunctionBegin;
2271   if (type == NORM_FROBENIUS) {
2272 #if defined(PETSC_USE_REAL___FP16)
2273     PetscBLASInt one = 1,cnt = bs2*nz;
2274     *norm = BLASnrm2_(&cnt,v,&one);
2275 #else
2276     for (i=0; i<bs2*nz; i++) {
2277       sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
2278     }
2279 #endif
2280     *norm = PetscSqrtReal(sum);
2281     ierr = PetscLogFlops(2*bs2*nz);CHKERRQ(ierr);
2282   } else if (type == NORM_1) { /* maximum column sum */
2283     PetscReal *tmp;
2284     PetscInt  *bcol = a->j;
2285     ierr = PetscCalloc1(A->cmap->n+1,&tmp);CHKERRQ(ierr);
2286     for (i=0; i<nz; i++) {
2287       for (j=0; j<bs; j++) {
2288         k1 = bs*(*bcol) + j; /* column index */
2289         for (k=0; k<bs; k++) {
2290           tmp[k1] += PetscAbsScalar(*v); v++;
2291         }
2292       }
2293       bcol++;
2294     }
2295     *norm = 0.0;
2296     for (j=0; j<A->cmap->n; j++) {
2297       if (tmp[j] > *norm) *norm = tmp[j];
2298     }
2299     ierr = PetscFree(tmp);CHKERRQ(ierr);
2300     ierr = PetscLogFlops(PetscMax(bs2*nz-1,0));CHKERRQ(ierr);
2301   } else if (type == NORM_INFINITY) { /* maximum row sum */
2302     *norm = 0.0;
2303     for (k=0; k<bs; k++) {
2304       for (j=0; j<a->mbs; j++) {
2305         v   = a->a + bs2*a->i[j] + k;
2306         sum = 0.0;
2307         for (i=0; i<a->i[j+1]-a->i[j]; i++) {
2308           for (k1=0; k1<bs; k1++) {
2309             sum += PetscAbsScalar(*v);
2310             v   += bs;
2311           }
2312         }
2313         if (sum > *norm) *norm = sum;
2314       }
2315     }
2316     ierr = PetscLogFlops(PetscMax(bs2*nz-1,0));CHKERRQ(ierr);
2317   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for this norm yet");
2318   PetscFunctionReturn(0);
2319 }
2320 
2321 
2322 PetscErrorCode MatEqual_SeqBAIJ(Mat A,Mat B,PetscBool * flg)
2323 {
2324   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)B->data;
2325   PetscErrorCode ierr;
2326 
2327   PetscFunctionBegin;
2328   /* If the  matrix/block dimensions are not equal, or no of nonzeros or shift */
2329   if ((A->rmap->N != B->rmap->N) || (A->cmap->n != B->cmap->n) || (A->rmap->bs != B->rmap->bs)|| (a->nz != b->nz)) {
2330     *flg = PETSC_FALSE;
2331     PetscFunctionReturn(0);
2332   }
2333 
2334   /* if the a->i are the same */
2335   ierr = PetscMemcmp(a->i,b->i,(a->mbs+1)*sizeof(PetscInt),flg);CHKERRQ(ierr);
2336   if (!*flg) PetscFunctionReturn(0);
2337 
2338   /* if a->j are the same */
2339   ierr = PetscMemcmp(a->j,b->j,(a->nz)*sizeof(PetscInt),flg);CHKERRQ(ierr);
2340   if (!*flg) PetscFunctionReturn(0);
2341 
2342   /* if a->a are the same */
2343   ierr = PetscMemcmp(a->a,b->a,(a->nz)*(A->rmap->bs)*(B->rmap->bs)*sizeof(PetscScalar),flg);CHKERRQ(ierr);
2344   PetscFunctionReturn(0);
2345 
2346 }
2347 
2348 PetscErrorCode MatGetDiagonal_SeqBAIJ(Mat A,Vec v)
2349 {
2350   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
2351   PetscErrorCode ierr;
2352   PetscInt       i,j,k,n,row,bs,*ai,*aj,ambs,bs2;
2353   PetscScalar    *x,zero = 0.0;
2354   MatScalar      *aa,*aa_j;
2355 
2356   PetscFunctionBegin;
2357   if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2358   bs   = A->rmap->bs;
2359   aa   = a->a;
2360   ai   = a->i;
2361   aj   = a->j;
2362   ambs = a->mbs;
2363   bs2  = a->bs2;
2364 
2365   ierr = VecSet(v,zero);CHKERRQ(ierr);
2366   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
2367   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
2368   if (n != A->rmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
2369   for (i=0; i<ambs; i++) {
2370     for (j=ai[i]; j<ai[i+1]; j++) {
2371       if (aj[j] == i) {
2372         row  = i*bs;
2373         aa_j = aa+j*bs2;
2374         for (k=0; k<bs2; k+=(bs+1),row++) x[row] = aa_j[k];
2375         break;
2376       }
2377     }
2378   }
2379   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
2380   PetscFunctionReturn(0);
2381 }
2382 
2383 PetscErrorCode MatDiagonalScale_SeqBAIJ(Mat A,Vec ll,Vec rr)
2384 {
2385   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
2386   const PetscScalar *l,*r,*li,*ri;
2387   PetscScalar       x;
2388   MatScalar         *aa, *v;
2389   PetscErrorCode    ierr;
2390   PetscInt          i,j,k,lm,rn,M,m,n,mbs,tmp,bs,bs2,iai;
2391   const PetscInt    *ai,*aj;
2392 
2393   PetscFunctionBegin;
2394   ai  = a->i;
2395   aj  = a->j;
2396   aa  = a->a;
2397   m   = A->rmap->n;
2398   n   = A->cmap->n;
2399   bs  = A->rmap->bs;
2400   mbs = a->mbs;
2401   bs2 = a->bs2;
2402   if (ll) {
2403     ierr = VecGetArrayRead(ll,&l);CHKERRQ(ierr);
2404     ierr = VecGetLocalSize(ll,&lm);CHKERRQ(ierr);
2405     if (lm != m) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Left scaling vector wrong length");
2406     for (i=0; i<mbs; i++) { /* for each block row */
2407       M  = ai[i+1] - ai[i];
2408       li = l + i*bs;
2409       v  = aa + bs2*ai[i];
2410       for (j=0; j<M; j++) { /* for each block */
2411         for (k=0; k<bs2; k++) {
2412           (*v++) *= li[k%bs];
2413         }
2414       }
2415     }
2416     ierr = VecRestoreArrayRead(ll,&l);CHKERRQ(ierr);
2417     ierr = PetscLogFlops(a->nz);CHKERRQ(ierr);
2418   }
2419 
2420   if (rr) {
2421     ierr = VecGetArrayRead(rr,&r);CHKERRQ(ierr);
2422     ierr = VecGetLocalSize(rr,&rn);CHKERRQ(ierr);
2423     if (rn != n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Right scaling vector wrong length");
2424     for (i=0; i<mbs; i++) { /* for each block row */
2425       iai = ai[i];
2426       M   = ai[i+1] - iai;
2427       v   = aa + bs2*iai;
2428       for (j=0; j<M; j++) { /* for each block */
2429         ri = r + bs*aj[iai+j];
2430         for (k=0; k<bs; k++) {
2431           x = ri[k];
2432           for (tmp=0; tmp<bs; tmp++) v[tmp] *= x;
2433           v += bs;
2434         }
2435       }
2436     }
2437     ierr = VecRestoreArrayRead(rr,&r);CHKERRQ(ierr);
2438     ierr = PetscLogFlops(a->nz);CHKERRQ(ierr);
2439   }
2440   PetscFunctionReturn(0);
2441 }
2442 
2443 
2444 PetscErrorCode MatGetInfo_SeqBAIJ(Mat A,MatInfoType flag,MatInfo *info)
2445 {
2446   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
2447 
2448   PetscFunctionBegin;
2449   info->block_size   = a->bs2;
2450   info->nz_allocated = a->bs2*a->maxnz;
2451   info->nz_used      = a->bs2*a->nz;
2452   info->nz_unneeded  = (double)(info->nz_allocated - info->nz_used);
2453   info->assemblies   = A->num_ass;
2454   info->mallocs      = A->info.mallocs;
2455   info->memory       = ((PetscObject)A)->mem;
2456   if (A->factortype) {
2457     info->fill_ratio_given  = A->info.fill_ratio_given;
2458     info->fill_ratio_needed = A->info.fill_ratio_needed;
2459     info->factor_mallocs    = A->info.factor_mallocs;
2460   } else {
2461     info->fill_ratio_given  = 0;
2462     info->fill_ratio_needed = 0;
2463     info->factor_mallocs    = 0;
2464   }
2465   PetscFunctionReturn(0);
2466 }
2467 
2468 PetscErrorCode MatZeroEntries_SeqBAIJ(Mat A)
2469 {
2470   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
2471   PetscErrorCode ierr;
2472 
2473   PetscFunctionBegin;
2474   ierr = PetscMemzero(a->a,a->bs2*a->i[a->mbs]*sizeof(MatScalar));CHKERRQ(ierr);
2475   PetscFunctionReturn(0);
2476 }
2477