xref: /petsc/src/mat/impls/sbaij/seq/sbaij2.c (revision 2b8d69ca7ea5fe9190df62c1dce3bbd66fce84dd)
1 
2 #include <../src/mat/impls/baij/seq/baij.h>
3 #include <petsc/private/kernels/blockinvert.h>
4 #include <petscbt.h>
5 #include <../src/mat/impls/sbaij/seq/sbaij.h>
6 #include <petscblaslapack.h>
7 
8 #undef __FUNCT__
9 #define __FUNCT__ "MatIncreaseOverlap_SeqSBAIJ"
10 PetscErrorCode MatIncreaseOverlap_SeqSBAIJ(Mat A,PetscInt is_max,IS is[],PetscInt ov)
11 {
12   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
13   PetscErrorCode ierr;
14   PetscInt       brow,i,j,k,l,mbs,n,*nidx,isz,bcol,bcol_max,start,end,*ai,*aj,bs,*nidx2;
15   const PetscInt *idx;
16   PetscBT        table_out,table_in;
17 
18   PetscFunctionBegin;
19   if (ov < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap specified");
20   mbs  = a->mbs;
21   ai   = a->i;
22   aj   = a->j;
23   bs   = A->rmap->bs;
24   ierr = PetscBTCreate(mbs,&table_out);CHKERRQ(ierr);
25   ierr = PetscMalloc1(mbs+1,&nidx);CHKERRQ(ierr);
26   ierr = PetscMalloc1(A->rmap->N+1,&nidx2);CHKERRQ(ierr);
27   ierr = PetscBTCreate(mbs,&table_in);CHKERRQ(ierr);
28 
29   for (i=0; i<is_max; i++) { /* for each is */
30     isz  = 0;
31     ierr = PetscBTMemzero(mbs,table_out);CHKERRQ(ierr);
32 
33     /* Extract the indices, assume there can be duplicate entries */
34     ierr = ISGetIndices(is[i],&idx);CHKERRQ(ierr);
35     ierr = ISGetLocalSize(is[i],&n);CHKERRQ(ierr);
36 
37     /* Enter these into the temp arrays i.e mark table_out[brow], enter brow into new index */
38     bcol_max = 0;
39     for (j=0; j<n; ++j) {
40       brow = idx[j]/bs; /* convert the indices into block indices */
41       if (brow >= mbs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"index greater than mat-dim");
42       if (!PetscBTLookupSet(table_out,brow)) {
43         nidx[isz++] = brow;
44         if (bcol_max < brow) bcol_max = brow;
45       }
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       /* set table_in for lookup - only mark entries that are added onto nidx in (j-1)-th overlap */
53       ierr = PetscBTMemzero(mbs,table_in);CHKERRQ(ierr);
54       for (l=k; l<isz; l++) { ierr = PetscBTSet(table_in,nidx[l]);CHKERRQ(ierr); }
55 
56       n = isz;  /* length of the updated is[i] */
57       for (brow=0; brow<mbs; brow++) {
58         start = ai[brow]; end   = ai[brow+1];
59         if (PetscBTLookup(table_in,brow)) { /* brow is on nidx - row search: collect all bcol in this brow */
60           for (l = start; l<end; l++) {
61             bcol = aj[l];
62             if (!PetscBTLookupSet(table_out,bcol)) {
63               nidx[isz++] = bcol;
64               if (bcol_max < bcol) bcol_max = bcol;
65             }
66           }
67           k++;
68           if (k >= n) break; /* for (brow=0; brow<mbs; brow++) */
69         } else { /* brow is not on nidx - col serach: add brow onto nidx if there is a bcol in nidx */
70           for (l = start; l<end; l++) {
71             bcol = aj[l];
72             if (bcol > bcol_max) break;
73             if (PetscBTLookup(table_in,bcol)) {
74               if (!PetscBTLookupSet(table_out,brow)) nidx[isz++] = brow;
75               break; /* for l = start; l<end ; l++) */
76             }
77           }
78         }
79       }
80     } /* for each overlap */
81 
82     /* expand the Index Set */
83     for (j=0; j<isz; j++) {
84       for (k=0; k<bs; k++) nidx2[j*bs+k] = nidx[j]*bs+k;
85     }
86     ierr = ISCreateGeneral(PETSC_COMM_SELF,isz*bs,nidx2,PETSC_COPY_VALUES,is+i);CHKERRQ(ierr);
87   }
88   ierr = PetscBTDestroy(&table_out);CHKERRQ(ierr);
89   ierr = PetscFree(nidx);CHKERRQ(ierr);
90   ierr = PetscFree(nidx2);CHKERRQ(ierr);
91   ierr = PetscBTDestroy(&table_in);CHKERRQ(ierr);
92   PetscFunctionReturn(0);
93 }
94 
95 #undef __FUNCT__
96 #define __FUNCT__ "MatGetSubMatrix_SeqSBAIJ_Private"
97 PetscErrorCode MatGetSubMatrix_SeqSBAIJ_Private(Mat A,IS isrow,MatReuse scall,Mat *B)
98 {
99   Mat_SeqSBAIJ    *a = (Mat_SeqSBAIJ*)A->data,*c;
100   PetscErrorCode  ierr;
101   PetscInt        *smap,i,k,kstart,kend,oldcols = a->mbs,*lens;
102   PetscInt        row,mat_i,*mat_j,tcol,*mat_ilen;
103   PetscInt        nrows,*ssmap,bs=A->rmap->bs,bs2=a->bs2;
104   const PetscInt  *irow,*aj = a->j,*ai = a->i;
105   MatScalar       *mat_a;
106   Mat             C;
107   PetscBool       flag,sorted;
108 
109   PetscFunctionBegin;
110   ierr = ISSorted(isrow,&sorted);CHKERRQ(ierr);
111   if (!sorted) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"IS is not sorted");
112 
113   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
114   ierr = ISGetSize(isrow,&nrows);CHKERRQ(ierr);
115 
116   ierr  = PetscMalloc1(oldcols,&smap);CHKERRQ(ierr);
117   ierr  = PetscMemzero(smap,oldcols*sizeof(PetscInt));CHKERRQ(ierr);
118   ssmap = smap;
119   ierr  = PetscMalloc1(1+nrows,&lens);CHKERRQ(ierr);
120   for (i=0; i<nrows; i++) smap[irow[i]] = i+1; /* nrows = ncols */
121   /* determine lens of each row */
122   for (i=0; i<nrows; i++) {
123     kstart  = ai[irow[i]];
124     kend    = kstart + a->ilen[irow[i]];
125     lens[i] = 0;
126     for (k=kstart; k<kend; k++) {
127       if (ssmap[aj[k]]) lens[i]++;
128     }
129   }
130   /* Create and fill new matrix */
131   if (scall == MAT_REUSE_MATRIX) {
132     c = (Mat_SeqSBAIJ*)((*B)->data);
133 
134     if (c->mbs!=nrows || (*B)->rmap->bs!=bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Submatrix wrong size");
135     ierr = PetscMemcmp(c->ilen,lens,c->mbs *sizeof(PetscInt),&flag);CHKERRQ(ierr);
136     if (!flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong no of nonzeros");
137     ierr = PetscMemzero(c->ilen,c->mbs*sizeof(PetscInt));CHKERRQ(ierr);
138     C    = *B;
139   } else {
140     ierr = MatCreate(PetscObjectComm((PetscObject)A),&C);CHKERRQ(ierr);
141     ierr = MatSetSizes(C,nrows*bs,nrows*bs,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
142     ierr = MatSetType(C,((PetscObject)A)->type_name);CHKERRQ(ierr);
143     ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(C,bs,0,lens);CHKERRQ(ierr);
144   }
145   c = (Mat_SeqSBAIJ*)(C->data);
146   for (i=0; i<nrows; i++) {
147     row      = irow[i];
148     kstart   = ai[row];
149     kend     = kstart + a->ilen[row];
150     mat_i    = c->i[i];
151     mat_j    = c->j + mat_i;
152     mat_a    = c->a + mat_i*bs2;
153     mat_ilen = c->ilen + i;
154     for (k=kstart; k<kend; k++) {
155       if ((tcol=ssmap[a->j[k]])) {
156         *mat_j++ = tcol - 1;
157         ierr     = PetscMemcpy(mat_a,a->a+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
158         mat_a   += bs2;
159         (*mat_ilen)++;
160       }
161     }
162   }
163 
164   /* Free work space */
165   ierr = PetscFree(smap);CHKERRQ(ierr);
166   ierr = PetscFree(lens);CHKERRQ(ierr);
167   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
168   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
169 
170   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
171   *B   = C;
172   PetscFunctionReturn(0);
173 }
174 
175 #undef __FUNCT__
176 #define __FUNCT__ "MatGetSubMatrix_SeqSBAIJ"
177 PetscErrorCode MatGetSubMatrix_SeqSBAIJ(Mat A,IS isrow,IS iscol,MatReuse scall,Mat *B)
178 {
179   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
180   IS             is1;
181   PetscErrorCode ierr;
182   PetscInt       *vary,*iary,nrows,i,bs=A->rmap->bs,count;
183   const PetscInt *irow;
184 
185   PetscFunctionBegin;
186   if (isrow != iscol) {
187     PetscBool isequal;
188     ierr = ISEqual(isrow,iscol,&isequal);CHKERRQ(ierr);
189     if (!isequal) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"For symmetric format, iscol must equal isrow");
190   }
191 
192   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
193   ierr = ISGetSize(isrow,&nrows);CHKERRQ(ierr);
194 
195   /* Verify if the indices corespond to each element in a block
196    and form the IS with compressed IS */
197   ierr = PetscMalloc2(a->mbs,&vary,a->mbs,&iary);CHKERRQ(ierr);
198   ierr = PetscMemzero(vary,a->mbs*sizeof(PetscInt));CHKERRQ(ierr);
199   for (i=0; i<nrows; i++) vary[irow[i]/bs]++;
200 
201   count = 0;
202   for (i=0; i<a->mbs; i++) {
203     if (vary[i]!=0 && vary[i]!=bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Index set does not match blocks");
204     if (vary[i]==bs) iary[count++] = i;
205   }
206   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
207   ierr = ISCreateGeneral(PETSC_COMM_SELF,count,iary,PETSC_COPY_VALUES,&is1);CHKERRQ(ierr);
208   ierr = PetscFree2(vary,iary);CHKERRQ(ierr);
209 
210   ierr = MatGetSubMatrix_SeqSBAIJ_Private(A,is1,scall,B);CHKERRQ(ierr);
211   ierr = ISDestroy(&is1);CHKERRQ(ierr);
212   PetscFunctionReturn(0);
213 }
214 
215 #undef __FUNCT__
216 #define __FUNCT__ "MatGetSubMatrices_SeqSBAIJ"
217 PetscErrorCode MatGetSubMatrices_SeqSBAIJ(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
218 {
219   PetscErrorCode ierr;
220   PetscInt       i;
221 
222   PetscFunctionBegin;
223   if (scall == MAT_INITIAL_MATRIX) {
224     ierr = PetscMalloc1(n+1,B);CHKERRQ(ierr);
225   }
226 
227   for (i=0; i<n; i++) {
228     ierr = MatGetSubMatrix_SeqSBAIJ(A,irow[i],icol[i],scall,&(*B)[i]);CHKERRQ(ierr);
229   }
230   PetscFunctionReturn(0);
231 }
232 
233 /* -------------------------------------------------------*/
234 /* Should check that shapes of vectors and matrices match */
235 /* -------------------------------------------------------*/
236 
237 #undef __FUNCT__
238 #define __FUNCT__ "MatMult_SeqSBAIJ_2"
239 PetscErrorCode MatMult_SeqSBAIJ_2(Mat A,Vec xx,Vec zz)
240 {
241   Mat_SeqSBAIJ      *a = (Mat_SeqSBAIJ*)A->data;
242   PetscScalar       *z,x1,x2,zero=0.0;
243   const PetscScalar *x,*xb;
244   const MatScalar   *v;
245   PetscErrorCode    ierr;
246   PetscInt          mbs = a->mbs,i,n,cval,j,jmin;
247   const PetscInt    *aj=a->j,*ai=a->i,*ib;
248   PetscInt          nonzerorow=0;
249 
250   PetscFunctionBegin;
251   ierr = VecSet(zz,zero);CHKERRQ(ierr);
252   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
253   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
254 
255   v  = a->a;
256   xb = x;
257 
258   for (i=0; i<mbs; i++) {
259     n           = ai[1] - ai[0]; /* length of i_th block row of A */
260     x1          = xb[0]; x2 = xb[1];
261     ib          = aj + *ai;
262     jmin        = 0;
263     nonzerorow += (n>0);
264     if (*ib == i) {     /* (diag of A)*x */
265       z[2*i]   += v[0]*x1 + v[2]*x2;
266       z[2*i+1] += v[2]*x1 + v[3]*x2;
267       v        += 4; jmin++;
268     }
269     PetscPrefetchBlock(ib+jmin+n,n,0,PETSC_PREFETCH_HINT_NTA); /* Indices for the next row (assumes same size as this one) */
270     PetscPrefetchBlock(v+4*n,4*n,0,PETSC_PREFETCH_HINT_NTA);   /* Entries for the next row */
271     for (j=jmin; j<n; j++) {
272       /* (strict lower triangular part of A)*x  */
273       cval       = ib[j]*2;
274       z[cval]   += v[0]*x1 + v[1]*x2;
275       z[cval+1] += v[2]*x1 + v[3]*x2;
276       /* (strict upper triangular part of A)*x  */
277       z[2*i]   += v[0]*x[cval] + v[2]*x[cval+1];
278       z[2*i+1] += v[1]*x[cval] + v[3]*x[cval+1];
279       v        += 4;
280     }
281     xb +=2; ai++;
282   }
283 
284   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
285   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
286   ierr = PetscLogFlops(8.0*(a->nz*2.0 - nonzerorow) - nonzerorow);CHKERRQ(ierr);
287   PetscFunctionReturn(0);
288 }
289 
290 #undef __FUNCT__
291 #define __FUNCT__ "MatMult_SeqSBAIJ_3"
292 PetscErrorCode MatMult_SeqSBAIJ_3(Mat A,Vec xx,Vec zz)
293 {
294   Mat_SeqSBAIJ      *a = (Mat_SeqSBAIJ*)A->data;
295   PetscScalar       *z,x1,x2,x3,zero=0.0;
296   const PetscScalar *x,*xb;
297   const MatScalar   *v;
298   PetscErrorCode    ierr;
299   PetscInt          mbs = a->mbs,i,n,cval,j,jmin;
300   const PetscInt    *aj = a->j,*ai = a->i,*ib;
301   PetscInt          nonzerorow=0;
302 
303   PetscFunctionBegin;
304   ierr = VecSet(zz,zero);CHKERRQ(ierr);
305   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
306   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
307 
308   v  = a->a;
309   xb = x;
310 
311   for (i=0; i<mbs; i++) {
312     n           = ai[1] - ai[0]; /* length of i_th block row of A */
313     x1          = xb[0]; x2 = xb[1]; x3 = xb[2];
314     ib          = aj + *ai;
315     jmin        = 0;
316     nonzerorow += (n>0);
317     if (*ib == i) {     /* (diag of A)*x */
318       z[3*i]   += v[0]*x1 + v[3]*x2 + v[6]*x3;
319       z[3*i+1] += v[3]*x1 + v[4]*x2 + v[7]*x3;
320       z[3*i+2] += v[6]*x1 + v[7]*x2 + v[8]*x3;
321       v        += 9; jmin++;
322     }
323     PetscPrefetchBlock(ib+jmin+n,n,0,PETSC_PREFETCH_HINT_NTA); /* Indices for the next row (assumes same size as this one) */
324     PetscPrefetchBlock(v+9*n,9*n,0,PETSC_PREFETCH_HINT_NTA);   /* Entries for the next row */
325     for (j=jmin; j<n; j++) {
326       /* (strict lower triangular part of A)*x  */
327       cval       = ib[j]*3;
328       z[cval]   += v[0]*x1 + v[1]*x2 + v[2]*x3;
329       z[cval+1] += v[3]*x1 + v[4]*x2 + v[5]*x3;
330       z[cval+2] += v[6]*x1 + v[7]*x2 + v[8]*x3;
331       /* (strict upper triangular part of A)*x  */
332       z[3*i]   += v[0]*x[cval] + v[3]*x[cval+1]+ v[6]*x[cval+2];
333       z[3*i+1] += v[1]*x[cval] + v[4]*x[cval+1]+ v[7]*x[cval+2];
334       z[3*i+2] += v[2]*x[cval] + v[5]*x[cval+1]+ v[8]*x[cval+2];
335       v        += 9;
336     }
337     xb +=3; ai++;
338   }
339 
340   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
341   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
342   ierr = PetscLogFlops(18.0*(a->nz*2.0 - nonzerorow) - nonzerorow);CHKERRQ(ierr);
343   PetscFunctionReturn(0);
344 }
345 
346 #undef __FUNCT__
347 #define __FUNCT__ "MatMult_SeqSBAIJ_4"
348 PetscErrorCode MatMult_SeqSBAIJ_4(Mat A,Vec xx,Vec zz)
349 {
350   Mat_SeqSBAIJ      *a = (Mat_SeqSBAIJ*)A->data;
351   PetscScalar       *z,x1,x2,x3,x4,zero=0.0;
352   const PetscScalar *x,*xb;
353   const MatScalar   *v;
354   PetscErrorCode    ierr;
355   PetscInt          mbs = a->mbs,i,n,cval,j,jmin;
356   const PetscInt    *aj = a->j,*ai = a->i,*ib;
357   PetscInt          nonzerorow = 0;
358 
359   PetscFunctionBegin;
360   ierr = VecSet(zz,zero);CHKERRQ(ierr);
361   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
362   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
363 
364   v  = a->a;
365   xb = x;
366 
367   for (i=0; i<mbs; i++) {
368     n           = ai[1] - ai[0]; /* length of i_th block row of A */
369     x1          = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
370     ib          = aj + *ai;
371     jmin        = 0;
372     nonzerorow += (n>0);
373     if (*ib == i) {     /* (diag of A)*x */
374       z[4*i]   += v[0]*x1 + v[4]*x2 +  v[8]*x3 + v[12]*x4;
375       z[4*i+1] += v[4]*x1 + v[5]*x2 +  v[9]*x3 + v[13]*x4;
376       z[4*i+2] += v[8]*x1 + v[9]*x2 + v[10]*x3 + v[14]*x4;
377       z[4*i+3] += v[12]*x1+ v[13]*x2+ v[14]*x3 + v[15]*x4;
378       v        += 16; jmin++;
379     }
380     PetscPrefetchBlock(ib+jmin+n,n,0,PETSC_PREFETCH_HINT_NTA); /* Indices for the next row (assumes same size as this one) */
381     PetscPrefetchBlock(v+16*n,16*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
382     for (j=jmin; j<n; j++) {
383       /* (strict lower triangular part of A)*x  */
384       cval       = ib[j]*4;
385       z[cval]   += v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4;
386       z[cval+1] += v[4]*x1 + v[5]*x2 + v[6]*x3 + v[7]*x4;
387       z[cval+2] += v[8]*x1 + v[9]*x2 + v[10]*x3 + v[11]*x4;
388       z[cval+3] += v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4;
389       /* (strict upper triangular part of A)*x  */
390       z[4*i]   += v[0]*x[cval] + v[4]*x[cval+1]+ v[8]*x[cval+2] + v[12]*x[cval+3];
391       z[4*i+1] += v[1]*x[cval] + v[5]*x[cval+1]+ v[9]*x[cval+2] + v[13]*x[cval+3];
392       z[4*i+2] += v[2]*x[cval] + v[6]*x[cval+1]+ v[10]*x[cval+2]+ v[14]*x[cval+3];
393       z[4*i+3] += v[3]*x[cval] + v[7]*x[cval+1]+ v[11]*x[cval+2]+ v[15]*x[cval+3];
394       v        += 16;
395     }
396     xb +=4; ai++;
397   }
398 
399   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
400   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
401   ierr = PetscLogFlops(32.0*(a->nz*2.0 - nonzerorow) - nonzerorow);CHKERRQ(ierr);
402   PetscFunctionReturn(0);
403 }
404 
405 #undef __FUNCT__
406 #define __FUNCT__ "MatMult_SeqSBAIJ_5"
407 PetscErrorCode MatMult_SeqSBAIJ_5(Mat A,Vec xx,Vec zz)
408 {
409   Mat_SeqSBAIJ      *a = (Mat_SeqSBAIJ*)A->data;
410   PetscScalar       *z,x1,x2,x3,x4,x5,zero=0.0;
411   const PetscScalar *x,*xb;
412   const MatScalar   *v;
413   PetscErrorCode    ierr;
414   PetscInt          mbs = a->mbs,i,n,cval,j,jmin;
415   const PetscInt    *aj = a->j,*ai = a->i,*ib;
416   PetscInt          nonzerorow=0;
417 
418   PetscFunctionBegin;
419   ierr = VecSet(zz,zero);CHKERRQ(ierr);
420   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
421   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
422 
423   v  = a->a;
424   xb = x;
425 
426   for (i=0; i<mbs; i++) {
427     n           = ai[1] - ai[0]; /* length of i_th block row of A */
428     x1          = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5=xb[4];
429     ib          = aj + *ai;
430     jmin        = 0;
431     nonzerorow += (n>0);
432     if (*ib == i) {      /* (diag of A)*x */
433       z[5*i]   += v[0]*x1  + v[5]*x2 + v[10]*x3 + v[15]*x4+ v[20]*x5;
434       z[5*i+1] += v[5]*x1  + v[6]*x2 + v[11]*x3 + v[16]*x4+ v[21]*x5;
435       z[5*i+2] += v[10]*x1 +v[11]*x2 + v[12]*x3 + v[17]*x4+ v[22]*x5;
436       z[5*i+3] += v[15]*x1 +v[16]*x2 + v[17]*x3 + v[18]*x4+ v[23]*x5;
437       z[5*i+4] += v[20]*x1 +v[21]*x2 + v[22]*x3 + v[23]*x4+ v[24]*x5;
438       v        += 25; jmin++;
439     }
440     PetscPrefetchBlock(ib+jmin+n,n,0,PETSC_PREFETCH_HINT_NTA); /* Indices for the next row (assumes same size as this one) */
441     PetscPrefetchBlock(v+25*n,25*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
442     for (j=jmin; j<n; j++) {
443       /* (strict lower triangular part of A)*x  */
444       cval       = ib[j]*5;
445       z[cval]   += v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4 + v[4]*x5;
446       z[cval+1] += v[5]*x1 + v[6]*x2 + v[7]*x3 + v[8]*x4 + v[9]*x5;
447       z[cval+2] += v[10]*x1 + v[11]*x2 + v[12]*x3 + v[13]*x4+ v[14]*x5;
448       z[cval+3] += v[15]*x1 + v[16]*x2 + v[17]*x3 + v[18]*x4+ v[19]*x5;
449       z[cval+4] += v[20]*x1 + v[21]*x2 + v[22]*x3 + v[23]*x4+ v[24]*x5;
450       /* (strict upper triangular part of A)*x  */
451       z[5*i]   +=v[0]*x[cval]+v[5]*x[cval+1]+v[10]*x[cval+2]+v[15]*x[cval+3]+v[20]*x[cval+4];
452       z[5*i+1] +=v[1]*x[cval]+v[6]*x[cval+1]+v[11]*x[cval+2]+v[16]*x[cval+3]+v[21]*x[cval+4];
453       z[5*i+2] +=v[2]*x[cval]+v[7]*x[cval+1]+v[12]*x[cval+2]+v[17]*x[cval+3]+v[22]*x[cval+4];
454       z[5*i+3] +=v[3]*x[cval]+v[8]*x[cval+1]+v[13]*x[cval+2]+v[18]*x[cval+3]+v[23]*x[cval+4];
455       z[5*i+4] +=v[4]*x[cval]+v[9]*x[cval+1]+v[14]*x[cval+2]+v[19]*x[cval+3]+v[24]*x[cval+4];
456       v        += 25;
457     }
458     xb +=5; ai++;
459   }
460 
461   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
462   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
463   ierr = PetscLogFlops(50.0*(a->nz*2.0 - nonzerorow) - nonzerorow);CHKERRQ(ierr);
464   PetscFunctionReturn(0);
465 }
466 
467 
468 #undef __FUNCT__
469 #define __FUNCT__ "MatMult_SeqSBAIJ_6"
470 PetscErrorCode MatMult_SeqSBAIJ_6(Mat A,Vec xx,Vec zz)
471 {
472   Mat_SeqSBAIJ      *a = (Mat_SeqSBAIJ*)A->data;
473   PetscScalar       *z,x1,x2,x3,x4,x5,x6,zero=0.0;
474   const PetscScalar *x,*xb;
475   const MatScalar   *v;
476   PetscErrorCode    ierr;
477   PetscInt          mbs = a->mbs,i,n,cval,j,jmin;
478   const PetscInt    *aj=a->j,*ai=a->i,*ib;
479   PetscInt          nonzerorow=0;
480 
481   PetscFunctionBegin;
482   ierr = VecSet(zz,zero);CHKERRQ(ierr);
483   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
484   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
485 
486   v  = a->a;
487   xb = x;
488 
489   for (i=0; i<mbs; i++) {
490     n           = ai[1] - ai[0]; /* length of i_th block row of A */
491     x1          = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5=xb[4]; x6=xb[5];
492     ib          = aj + *ai;
493     jmin        = 0;
494     nonzerorow += (n>0);
495     if (*ib == i) {      /* (diag of A)*x */
496       z[6*i]   += v[0]*x1  + v[6]*x2 + v[12]*x3 + v[18]*x4+ v[24]*x5 + v[30]*x6;
497       z[6*i+1] += v[6]*x1  + v[7]*x2 + v[13]*x3 + v[19]*x4+ v[25]*x5 + v[31]*x6;
498       z[6*i+2] += v[12]*x1 +v[13]*x2 + v[14]*x3 + v[20]*x4+ v[26]*x5 + v[32]*x6;
499       z[6*i+3] += v[18]*x1 +v[19]*x2 + v[20]*x3 + v[21]*x4+ v[27]*x5 + v[33]*x6;
500       z[6*i+4] += v[24]*x1 +v[25]*x2 + v[26]*x3 + v[27]*x4+ v[28]*x5 + v[34]*x6;
501       z[6*i+5] += v[30]*x1 +v[31]*x2 + v[32]*x3 + v[33]*x4+ v[34]*x5 + v[35]*x6;
502       v        += 36; jmin++;
503     }
504     PetscPrefetchBlock(ib+jmin+n,n,0,PETSC_PREFETCH_HINT_NTA); /* Indices for the next row (assumes same size as this one) */
505     PetscPrefetchBlock(v+36*n,36*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
506     for (j=jmin; j<n; j++) {
507       /* (strict lower triangular part of A)*x  */
508       cval       = ib[j]*6;
509       z[cval]   += v[0]*x1  + v[1]*x2 + v[2]*x3 + v[3]*x4+ v[4]*x5 + v[5]*x6;
510       z[cval+1] += v[6]*x1  + v[7]*x2 + v[8]*x3 + v[9]*x4+ v[10]*x5 + v[11]*x6;
511       z[cval+2] += v[12]*x1  + v[13]*x2 + v[14]*x3 + v[15]*x4+ v[16]*x5 + v[17]*x6;
512       z[cval+3] += v[18]*x1  + v[19]*x2 + v[20]*x3 + v[21]*x4+ v[22]*x5 + v[23]*x6;
513       z[cval+4] += v[24]*x1  + v[25]*x2 + v[26]*x3 + v[27]*x4+ v[28]*x5 + v[29]*x6;
514       z[cval+5] += v[30]*x1  + v[31]*x2 + v[32]*x3 + v[33]*x4+ v[34]*x5 + v[35]*x6;
515       /* (strict upper triangular part of A)*x  */
516       z[6*i]   +=v[0]*x[cval]+v[6]*x[cval+1]+v[12]*x[cval+2]+v[18]*x[cval+3]+v[24]*x[cval+4]+v[30]*x[cval+5];
517       z[6*i+1] +=v[1]*x[cval]+v[7]*x[cval+1]+v[13]*x[cval+2]+v[19]*x[cval+3]+v[25]*x[cval+4]+v[31]*x[cval+5];
518       z[6*i+2] +=v[2]*x[cval]+v[8]*x[cval+1]+v[14]*x[cval+2]+v[20]*x[cval+3]+v[26]*x[cval+4]+v[32]*x[cval+5];
519       z[6*i+3] +=v[3]*x[cval]+v[9]*x[cval+1]+v[15]*x[cval+2]+v[21]*x[cval+3]+v[27]*x[cval+4]+v[33]*x[cval+5];
520       z[6*i+4] +=v[4]*x[cval]+v[10]*x[cval+1]+v[16]*x[cval+2]+v[22]*x[cval+3]+v[28]*x[cval+4]+v[34]*x[cval+5];
521       z[6*i+5] +=v[5]*x[cval]+v[11]*x[cval+1]+v[17]*x[cval+2]+v[23]*x[cval+3]+v[29]*x[cval+4]+v[35]*x[cval+5];
522       v        += 36;
523     }
524     xb +=6; ai++;
525   }
526 
527   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
528   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
529   ierr = PetscLogFlops(72.0*(a->nz*2.0 - nonzerorow) - nonzerorow);CHKERRQ(ierr);
530   PetscFunctionReturn(0);
531 }
532 #undef __FUNCT__
533 #define __FUNCT__ "MatMult_SeqSBAIJ_7"
534 PetscErrorCode MatMult_SeqSBAIJ_7(Mat A,Vec xx,Vec zz)
535 {
536   Mat_SeqSBAIJ      *a = (Mat_SeqSBAIJ*)A->data;
537   PetscScalar       *z,x1,x2,x3,x4,x5,x6,x7,zero=0.0;
538   const PetscScalar *x,*xb;
539   const MatScalar   *v;
540   PetscErrorCode    ierr;
541   PetscInt          mbs = a->mbs,i,n,cval,j,jmin;
542   const PetscInt    *aj=a->j,*ai=a->i,*ib;
543   PetscInt          nonzerorow=0;
544 
545   PetscFunctionBegin;
546   ierr = VecSet(zz,zero);CHKERRQ(ierr);
547   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
548   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
549 
550   v  = a->a;
551   xb = x;
552 
553   for (i=0; i<mbs; i++) {
554     n           = ai[1] - ai[0]; /* length of i_th block row of A */
555     x1          = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5=xb[4]; x6=xb[5]; x7=xb[6];
556     ib          = aj + *ai;
557     jmin        = 0;
558     nonzerorow += (n>0);
559     if (*ib == i) {      /* (diag of A)*x */
560       z[7*i]   += v[0]*x1 + v[7]*x2 + v[14]*x3 + v[21]*x4+ v[28]*x5 + v[35]*x6+ v[42]*x7;
561       z[7*i+1] += v[7]*x1 + v[8]*x2 + v[15]*x3 + v[22]*x4+ v[29]*x5 + v[36]*x6+ v[43]*x7;
562       z[7*i+2] += v[14]*x1+ v[15]*x2 +v[16]*x3 + v[23]*x4+ v[30]*x5 + v[37]*x6+ v[44]*x7;
563       z[7*i+3] += v[21]*x1+ v[22]*x2 +v[23]*x3 + v[24]*x4+ v[31]*x5 + v[38]*x6+ v[45]*x7;
564       z[7*i+4] += v[28]*x1+ v[29]*x2 +v[30]*x3 + v[31]*x4+ v[32]*x5 + v[39]*x6+ v[46]*x7;
565       z[7*i+5] += v[35]*x1+ v[36]*x2 +v[37]*x3 + v[38]*x4+ v[39]*x5 + v[40]*x6+ v[47]*x7;
566       z[7*i+6] += v[42]*x1+ v[43]*x2 +v[44]*x3 + v[45]*x4+ v[46]*x5 + v[47]*x6+ v[48]*x7;
567       v        += 49; jmin++;
568     }
569     PetscPrefetchBlock(ib+jmin+n,n,0,PETSC_PREFETCH_HINT_NTA); /* Indices for the next row (assumes same size as this one) */
570     PetscPrefetchBlock(v+49*n,49*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
571     for (j=jmin; j<n; j++) {
572       /* (strict lower triangular part of A)*x  */
573       cval       = ib[j]*7;
574       z[cval]   += v[0]*x1  + v[1]*x2 + v[2]*x3 + v[3]*x4+ v[4]*x5 + v[5]*x6+ v[6]*x7;
575       z[cval+1] += v[7]*x1  + v[8]*x2 + v[9]*x3 + v[10]*x4+ v[11]*x5 + v[12]*x6+ v[13]*x7;
576       z[cval+2] += v[14]*x1  + v[15]*x2 + v[16]*x3 + v[17]*x4+ v[18]*x5 + v[19]*x6+ v[20]*x7;
577       z[cval+3] += v[21]*x1  + v[22]*x2 + v[23]*x3 + v[24]*x4+ v[25]*x5 + v[26]*x6+ v[27]*x7;
578       z[cval+4] += v[28]*x1  + v[29]*x2 + v[30]*x3 + v[31]*x4+ v[32]*x5 + v[33]*x6+ v[34]*x7;
579       z[cval+5] += v[35]*x1  + v[36]*x2 + v[37]*x3 + v[38]*x4+ v[39]*x5 + v[40]*x6+ v[41]*x7;
580       z[cval+6] += v[42]*x1  + v[43]*x2 + v[44]*x3 + v[45]*x4+ v[46]*x5 + v[47]*x6+ v[48]*x7;
581       /* (strict upper triangular part of A)*x  */
582       z[7*i]  +=v[0]*x[cval]+v[7]*x[cval+1]+v[14]*x[cval+2]+v[21]*x[cval+3]+v[28]*x[cval+4]+v[35]*x[cval+5]+v[42]*x[cval+6];
583       z[7*i+1]+=v[1]*x[cval]+v[8]*x[cval+1]+v[15]*x[cval+2]+v[22]*x[cval+3]+v[29]*x[cval+4]+v[36]*x[cval+5]+v[43]*x[cval+6];
584       z[7*i+2]+=v[2]*x[cval]+v[9]*x[cval+1]+v[16]*x[cval+2]+v[23]*x[cval+3]+v[30]*x[cval+4]+v[37]*x[cval+5]+v[44]*x[cval+6];
585       z[7*i+3]+=v[3]*x[cval]+v[10]*x[cval+1]+v[17]*x[cval+2]+v[24]*x[cval+3]+v[31]*x[cval+4]+v[38]*x[cval+5]+v[45]*x[cval+6];
586       z[7*i+4]+=v[4]*x[cval]+v[11]*x[cval+1]+v[18]*x[cval+2]+v[25]*x[cval+3]+v[32]*x[cval+4]+v[39]*x[cval+5]+v[46]*x[cval+6];
587       z[7*i+5]+=v[5]*x[cval]+v[12]*x[cval+1]+v[19]*x[cval+2]+v[26]*x[cval+3]+v[33]*x[cval+4]+v[40]*x[cval+5]+v[47]*x[cval+6];
588       z[7*i+6]+=v[6]*x[cval]+v[13]*x[cval+1]+v[20]*x[cval+2]+v[27]*x[cval+3]+v[34]*x[cval+4]+v[41]*x[cval+5]+v[48]*x[cval+6];
589       v       += 49;
590     }
591     xb +=7; ai++;
592   }
593   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
594   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
595   ierr = PetscLogFlops(98.0*(a->nz*2.0 - nonzerorow) - nonzerorow);CHKERRQ(ierr);
596   PetscFunctionReturn(0);
597 }
598 
599 /*
600     This will not work with MatScalar == float because it calls the BLAS
601 */
602 #undef __FUNCT__
603 #define __FUNCT__ "MatMult_SeqSBAIJ_N"
604 PetscErrorCode MatMult_SeqSBAIJ_N(Mat A,Vec xx,Vec zz)
605 {
606   Mat_SeqSBAIJ      *a = (Mat_SeqSBAIJ*)A->data;
607   PetscScalar       *z,*z_ptr,*zb,*work,*workt,zero=0.0;
608   const PetscScalar *x,*x_ptr,*xb;
609   const MatScalar   *v;
610   PetscErrorCode    ierr;
611   PetscInt          mbs =a->mbs,i,bs=A->rmap->bs,j,n,bs2=a->bs2,ncols,k;
612   const PetscInt    *idx,*aj,*ii;
613   PetscInt          nonzerorow=0;
614 
615   PetscFunctionBegin;
616   ierr = VecSet(zz,zero);CHKERRQ(ierr);
617   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);x_ptr = x;
618   ierr = VecGetArray(zz,&z);CHKERRQ(ierr); z_ptr=z;
619 
620   aj = a->j;
621   v  = a->a;
622   ii = a->i;
623 
624   if (!a->mult_work) {
625     ierr = PetscMalloc1(A->rmap->N+1,&a->mult_work);CHKERRQ(ierr);
626   }
627   work = a->mult_work;
628 
629   for (i=0; i<mbs; i++) {
630     n           = ii[1] - ii[0]; ncols = n*bs;
631     workt       = work; idx=aj+ii[0];
632     nonzerorow += (n>0);
633 
634     /* upper triangular part */
635     for (j=0; j<n; j++) {
636       xb = x_ptr + bs*(*idx++);
637       for (k=0; k<bs; k++) workt[k] = xb[k];
638       workt += bs;
639     }
640     /* z(i*bs:(i+1)*bs-1) += A(i,:)*x */
641     PetscKernel_w_gets_w_plus_Ar_times_v(bs,ncols,work,v,z);
642 
643     /* strict lower triangular part */
644     idx = aj+ii[0];
645     if (*idx == i) {
646       ncols -= bs; v += bs2; idx++; n--;
647     }
648 
649     if (ncols > 0) {
650       workt = work;
651       ierr  = PetscMemzero(workt,ncols*sizeof(PetscScalar));CHKERRQ(ierr);
652       PetscKernel_w_gets_w_plus_trans_Ar_times_v(bs,ncols,x,v,workt);
653       for (j=0; j<n; j++) {
654         zb = z_ptr + bs*(*idx++);
655         for (k=0; k<bs; k++) zb[k] += workt[k];
656         workt += bs;
657       }
658     }
659     x += bs; v += n*bs2; z += bs; ii++;
660   }
661 
662   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
663   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
664   ierr = PetscLogFlops(2.0*(a->nz*2.0 - nonzerorow)*bs2 - nonzerorow);CHKERRQ(ierr);
665   PetscFunctionReturn(0);
666 }
667 
668 #undef __FUNCT__
669 #define __FUNCT__ "MatMultAdd_SeqSBAIJ_1"
670 PetscErrorCode MatMultAdd_SeqSBAIJ_1(Mat A,Vec xx,Vec yy,Vec zz)
671 {
672   Mat_SeqSBAIJ      *a = (Mat_SeqSBAIJ*)A->data;
673   PetscScalar       *z,x1;
674   const PetscScalar *x,*xb;
675   const MatScalar   *v;
676   PetscErrorCode    ierr;
677   PetscInt          mbs =a->mbs,i,n,cval,j,jmin;
678   const PetscInt    *aj=a->j,*ai=a->i,*ib;
679   PetscInt          nonzerorow=0;
680 
681   PetscFunctionBegin;
682   ierr = VecCopy(yy,zz);CHKERRQ(ierr);
683   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
684   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
685   v    = a->a;
686   xb   = x;
687 
688   for (i=0; i<mbs; i++) {
689     n           = ai[1] - ai[0]; /* length of i_th row of A */
690     x1          = xb[0];
691     ib          = aj + *ai;
692     jmin        = 0;
693     nonzerorow += (n>0);
694     if (*ib == i) {            /* (diag of A)*x */
695       z[i] += *v++ * x[*ib++]; jmin++;
696     }
697     for (j=jmin; j<n; j++) {
698       cval    = *ib;
699       z[cval] += *v * x1;      /* (strict lower triangular part of A)*x  */
700       z[i] += *v++ * x[*ib++]; /* (strict upper triangular part of A)*x  */
701     }
702     xb++; ai++;
703   }
704 
705   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
706   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
707 
708   ierr = PetscLogFlops(2.0*(a->nz*2.0 - nonzerorow));CHKERRQ(ierr);
709   PetscFunctionReturn(0);
710 }
711 
712 #undef __FUNCT__
713 #define __FUNCT__ "MatMultAdd_SeqSBAIJ_2"
714 PetscErrorCode MatMultAdd_SeqSBAIJ_2(Mat A,Vec xx,Vec yy,Vec zz)
715 {
716   Mat_SeqSBAIJ      *a = (Mat_SeqSBAIJ*)A->data;
717   PetscScalar       *z,x1,x2;
718   const PetscScalar *x,*xb;
719   const MatScalar   *v;
720   PetscErrorCode    ierr;
721   PetscInt          mbs =a->mbs,i,n,cval,j,jmin;
722   const PetscInt    *aj=a->j,*ai=a->i,*ib;
723   PetscInt          nonzerorow=0;
724 
725   PetscFunctionBegin;
726   ierr = VecCopy(yy,zz);CHKERRQ(ierr);
727   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
728   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
729 
730   v  = a->a;
731   xb = x;
732 
733   for (i=0; i<mbs; i++) {
734     n           = ai[1] - ai[0]; /* length of i_th block row of A */
735     x1          = xb[0]; x2 = xb[1];
736     ib          = aj + *ai;
737     jmin        = 0;
738     nonzerorow += (n>0);
739     if (*ib == i) {      /* (diag of A)*x */
740       z[2*i]   += v[0]*x1 + v[2]*x2;
741       z[2*i+1] += v[2]*x1 + v[3]*x2;
742       v        += 4; jmin++;
743     }
744     PetscPrefetchBlock(ib+jmin+n,n,0,PETSC_PREFETCH_HINT_NTA); /* Indices for the next row (assumes same size as this one) */
745     PetscPrefetchBlock(v+4*n,4*n,0,PETSC_PREFETCH_HINT_NTA);   /* Entries for the next row */
746     for (j=jmin; j<n; j++) {
747       /* (strict lower triangular part of A)*x  */
748       cval       = ib[j]*2;
749       z[cval]   += v[0]*x1 + v[1]*x2;
750       z[cval+1] += v[2]*x1 + v[3]*x2;
751       /* (strict upper triangular part of A)*x  */
752       z[2*i]   += v[0]*x[cval] + v[2]*x[cval+1];
753       z[2*i+1] += v[1]*x[cval] + v[3]*x[cval+1];
754       v        += 4;
755     }
756     xb +=2; ai++;
757   }
758   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
759   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
760 
761   ierr = PetscLogFlops(4.0*(a->nz*2.0 - nonzerorow));CHKERRQ(ierr);
762   PetscFunctionReturn(0);
763 }
764 
765 #undef __FUNCT__
766 #define __FUNCT__ "MatMultAdd_SeqSBAIJ_3"
767 PetscErrorCode MatMultAdd_SeqSBAIJ_3(Mat A,Vec xx,Vec yy,Vec zz)
768 {
769   Mat_SeqSBAIJ      *a = (Mat_SeqSBAIJ*)A->data;
770   PetscScalar       *z,x1,x2,x3;
771   const PetscScalar *x,*xb;
772   const MatScalar   *v;
773   PetscErrorCode    ierr;
774   PetscInt          mbs = a->mbs,i,n,cval,j,jmin;
775   const PetscInt    *aj=a->j,*ai=a->i,*ib;
776   PetscInt          nonzerorow=0;
777 
778   PetscFunctionBegin;
779   ierr = VecCopy(yy,zz);CHKERRQ(ierr);
780   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
781   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
782 
783   v  = a->a;
784   xb = x;
785 
786   for (i=0; i<mbs; i++) {
787     n           = ai[1] - ai[0]; /* length of i_th block row of A */
788     x1          = xb[0]; x2 = xb[1]; x3 = xb[2];
789     ib          = aj + *ai;
790     jmin        = 0;
791     nonzerorow += (n>0);
792     if (*ib == i) {     /* (diag of A)*x */
793       z[3*i]   += v[0]*x1 + v[3]*x2 + v[6]*x3;
794       z[3*i+1] += v[3]*x1 + v[4]*x2 + v[7]*x3;
795       z[3*i+2] += v[6]*x1 + v[7]*x2 + v[8]*x3;
796       v        += 9; jmin++;
797     }
798     PetscPrefetchBlock(ib+jmin+n,n,0,PETSC_PREFETCH_HINT_NTA); /* Indices for the next row (assumes same size as this one) */
799     PetscPrefetchBlock(v+9*n,9*n,0,PETSC_PREFETCH_HINT_NTA);   /* Entries for the next row */
800     for (j=jmin; j<n; j++) {
801       /* (strict lower triangular part of A)*x  */
802       cval       = ib[j]*3;
803       z[cval]   += v[0]*x1 + v[1]*x2 + v[2]*x3;
804       z[cval+1] += v[3]*x1 + v[4]*x2 + v[5]*x3;
805       z[cval+2] += v[6]*x1 + v[7]*x2 + v[8]*x3;
806       /* (strict upper triangular part of A)*x  */
807       z[3*i]   += v[0]*x[cval] + v[3]*x[cval+1]+ v[6]*x[cval+2];
808       z[3*i+1] += v[1]*x[cval] + v[4]*x[cval+1]+ v[7]*x[cval+2];
809       z[3*i+2] += v[2]*x[cval] + v[5]*x[cval+1]+ v[8]*x[cval+2];
810       v        += 9;
811     }
812     xb +=3; ai++;
813   }
814 
815   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
816   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
817 
818   ierr = PetscLogFlops(18.0*(a->nz*2.0 - nonzerorow));CHKERRQ(ierr);
819   PetscFunctionReturn(0);
820 }
821 
822 #undef __FUNCT__
823 #define __FUNCT__ "MatMultAdd_SeqSBAIJ_4"
824 PetscErrorCode MatMultAdd_SeqSBAIJ_4(Mat A,Vec xx,Vec yy,Vec zz)
825 {
826   Mat_SeqSBAIJ      *a = (Mat_SeqSBAIJ*)A->data;
827   PetscScalar       *z,x1,x2,x3,x4;
828   const PetscScalar *x,*xb;
829   const MatScalar   *v;
830   PetscErrorCode    ierr;
831   PetscInt          mbs = a->mbs,i,n,cval,j,jmin;
832   const PetscInt    *aj=a->j,*ai=a->i,*ib;
833   PetscInt          nonzerorow=0;
834 
835   PetscFunctionBegin;
836   ierr = VecCopy(yy,zz);CHKERRQ(ierr);
837   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
838   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
839 
840   v  = a->a;
841   xb = x;
842 
843   for (i=0; i<mbs; i++) {
844     n           = ai[1] - ai[0]; /* length of i_th block row of A */
845     x1          = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
846     ib          = aj + *ai;
847     jmin        = 0;
848     nonzerorow += (n>0);
849     if (*ib == i) {      /* (diag of A)*x */
850       z[4*i]   += v[0]*x1 + v[4]*x2 +  v[8]*x3 + v[12]*x4;
851       z[4*i+1] += v[4]*x1 + v[5]*x2 +  v[9]*x3 + v[13]*x4;
852       z[4*i+2] += v[8]*x1 + v[9]*x2 + v[10]*x3 + v[14]*x4;
853       z[4*i+3] += v[12]*x1+ v[13]*x2+ v[14]*x3 + v[15]*x4;
854       v        += 16; jmin++;
855     }
856     PetscPrefetchBlock(ib+jmin+n,n,0,PETSC_PREFETCH_HINT_NTA); /* Indices for the next row (assumes same size as this one) */
857     PetscPrefetchBlock(v+16*n,16*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
858     for (j=jmin; j<n; j++) {
859       /* (strict lower triangular part of A)*x  */
860       cval       = ib[j]*4;
861       z[cval]   += v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4;
862       z[cval+1] += v[4]*x1 + v[5]*x2 + v[6]*x3 + v[7]*x4;
863       z[cval+2] += v[8]*x1 + v[9]*x2 + v[10]*x3 + v[11]*x4;
864       z[cval+3] += v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4;
865       /* (strict upper triangular part of A)*x  */
866       z[4*i]   += v[0]*x[cval] + v[4]*x[cval+1]+ v[8]*x[cval+2] + v[12]*x[cval+3];
867       z[4*i+1] += v[1]*x[cval] + v[5]*x[cval+1]+ v[9]*x[cval+2] + v[13]*x[cval+3];
868       z[4*i+2] += v[2]*x[cval] + v[6]*x[cval+1]+ v[10]*x[cval+2]+ v[14]*x[cval+3];
869       z[4*i+3] += v[3]*x[cval] + v[7]*x[cval+1]+ v[11]*x[cval+2]+ v[15]*x[cval+3];
870       v        += 16;
871     }
872     xb +=4; ai++;
873   }
874 
875   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
876   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
877 
878   ierr = PetscLogFlops(32.0*(a->nz*2.0 - nonzerorow));CHKERRQ(ierr);
879   PetscFunctionReturn(0);
880 }
881 
882 #undef __FUNCT__
883 #define __FUNCT__ "MatMultAdd_SeqSBAIJ_5"
884 PetscErrorCode MatMultAdd_SeqSBAIJ_5(Mat A,Vec xx,Vec yy,Vec zz)
885 {
886   Mat_SeqSBAIJ      *a = (Mat_SeqSBAIJ*)A->data;
887   PetscScalar       *z,x1,x2,x3,x4,x5;
888   const PetscScalar *x,*xb;
889   const MatScalar   *v;
890   PetscErrorCode    ierr;
891   PetscInt          mbs = a->mbs,i,n,cval,j,jmin;
892   const PetscInt    *aj=a->j,*ai=a->i,*ib;
893   PetscInt          nonzerorow=0;
894 
895   PetscFunctionBegin;
896   ierr = VecCopy(yy,zz);CHKERRQ(ierr);
897   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
898   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
899 
900   v  = a->a;
901   xb = x;
902 
903   for (i=0; i<mbs; i++) {
904     n           = ai[1] - ai[0]; /* length of i_th block row of A */
905     x1          = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5=xb[4];
906     ib          = aj + *ai;
907     jmin        = 0;
908     nonzerorow += (n>0);
909     if (*ib == i) {      /* (diag of A)*x */
910       z[5*i]   += v[0]*x1  + v[5]*x2 + v[10]*x3 + v[15]*x4+ v[20]*x5;
911       z[5*i+1] += v[5]*x1  + v[6]*x2 + v[11]*x3 + v[16]*x4+ v[21]*x5;
912       z[5*i+2] += v[10]*x1 +v[11]*x2 + v[12]*x3 + v[17]*x4+ v[22]*x5;
913       z[5*i+3] += v[15]*x1 +v[16]*x2 + v[17]*x3 + v[18]*x4+ v[23]*x5;
914       z[5*i+4] += v[20]*x1 +v[21]*x2 + v[22]*x3 + v[23]*x4+ v[24]*x5;
915       v        += 25; jmin++;
916     }
917     PetscPrefetchBlock(ib+jmin+n,n,0,PETSC_PREFETCH_HINT_NTA); /* Indices for the next row (assumes same size as this one) */
918     PetscPrefetchBlock(v+25*n,25*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
919     for (j=jmin; j<n; j++) {
920       /* (strict lower triangular part of A)*x  */
921       cval       = ib[j]*5;
922       z[cval]   += v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4 + v[4]*x5;
923       z[cval+1] += v[5]*x1 + v[6]*x2 + v[7]*x3 + v[8]*x4 + v[9]*x5;
924       z[cval+2] += v[10]*x1 + v[11]*x2 + v[12]*x3 + v[13]*x4+ v[14]*x5;
925       z[cval+3] += v[15]*x1 + v[16]*x2 + v[17]*x3 + v[18]*x4+ v[19]*x5;
926       z[cval+4] += v[20]*x1 + v[21]*x2 + v[22]*x3 + v[23]*x4+ v[24]*x5;
927       /* (strict upper triangular part of A)*x  */
928       z[5*i]   +=v[0]*x[cval]+v[5]*x[cval+1]+v[10]*x[cval+2]+v[15]*x[cval+3]+v[20]*x[cval+4];
929       z[5*i+1] +=v[1]*x[cval]+v[6]*x[cval+1]+v[11]*x[cval+2]+v[16]*x[cval+3]+v[21]*x[cval+4];
930       z[5*i+2] +=v[2]*x[cval]+v[7]*x[cval+1]+v[12]*x[cval+2]+v[17]*x[cval+3]+v[22]*x[cval+4];
931       z[5*i+3] +=v[3]*x[cval]+v[8]*x[cval+1]+v[13]*x[cval+2]+v[18]*x[cval+3]+v[23]*x[cval+4];
932       z[5*i+4] +=v[4]*x[cval]+v[9]*x[cval+1]+v[14]*x[cval+2]+v[19]*x[cval+3]+v[24]*x[cval+4];
933       v        += 25;
934     }
935     xb +=5; ai++;
936   }
937 
938   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
939   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
940 
941   ierr = PetscLogFlops(50.0*(a->nz*2.0 - nonzerorow));CHKERRQ(ierr);
942   PetscFunctionReturn(0);
943 }
944 #undef __FUNCT__
945 #define __FUNCT__ "MatMultAdd_SeqSBAIJ_6"
946 PetscErrorCode MatMultAdd_SeqSBAIJ_6(Mat A,Vec xx,Vec yy,Vec zz)
947 {
948   Mat_SeqSBAIJ      *a = (Mat_SeqSBAIJ*)A->data;
949   PetscScalar       *z,x1,x2,x3,x4,x5,x6;
950   const PetscScalar *x,*xb;
951   const MatScalar   *v;
952   PetscErrorCode    ierr;
953   PetscInt          mbs = a->mbs,i,n,cval,j,jmin;
954   const PetscInt    *aj=a->j,*ai=a->i,*ib;
955   PetscInt          nonzerorow=0;
956 
957   PetscFunctionBegin;
958   ierr = VecCopy(yy,zz);CHKERRQ(ierr);
959   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
960   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
961 
962   v  = a->a;
963   xb = x;
964 
965   for (i=0; i<mbs; i++) {
966     n           = ai[1] - ai[0]; /* length of i_th block row of A */
967     x1          = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5=xb[4]; x6=xb[5];
968     ib          = aj + *ai;
969     jmin        = 0;
970     nonzerorow += (n>0);
971     if (*ib == i) {     /* (diag of A)*x */
972       z[6*i]   += v[0]*x1  + v[6]*x2 + v[12]*x3 + v[18]*x4+ v[24]*x5 + v[30]*x6;
973       z[6*i+1] += v[6]*x1  + v[7]*x2 + v[13]*x3 + v[19]*x4+ v[25]*x5 + v[31]*x6;
974       z[6*i+2] += v[12]*x1 +v[13]*x2 + v[14]*x3 + v[20]*x4+ v[26]*x5 + v[32]*x6;
975       z[6*i+3] += v[18]*x1 +v[19]*x2 + v[20]*x3 + v[21]*x4+ v[27]*x5 + v[33]*x6;
976       z[6*i+4] += v[24]*x1 +v[25]*x2 + v[26]*x3 + v[27]*x4+ v[28]*x5 + v[34]*x6;
977       z[6*i+5] += v[30]*x1 +v[31]*x2 + v[32]*x3 + v[33]*x4+ v[34]*x5 + v[35]*x6;
978       v        += 36; jmin++;
979     }
980     PetscPrefetchBlock(ib+jmin+n,n,0,PETSC_PREFETCH_HINT_NTA); /* Indices for the next row (assumes same size as this one) */
981     PetscPrefetchBlock(v+36*n,36*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
982     for (j=jmin; j<n; j++) {
983       /* (strict lower triangular part of A)*x  */
984       cval       = ib[j]*6;
985       z[cval]   += v[0]*x1  + v[1]*x2 + v[2]*x3 + v[3]*x4+ v[4]*x5 + v[5]*x6;
986       z[cval+1] += v[6]*x1  + v[7]*x2 + v[8]*x3 + v[9]*x4+ v[10]*x5 + v[11]*x6;
987       z[cval+2] += v[12]*x1  + v[13]*x2 + v[14]*x3 + v[15]*x4+ v[16]*x5 + v[17]*x6;
988       z[cval+3] += v[18]*x1  + v[19]*x2 + v[20]*x3 + v[21]*x4+ v[22]*x5 + v[23]*x6;
989       z[cval+4] += v[24]*x1  + v[25]*x2 + v[26]*x3 + v[27]*x4+ v[28]*x5 + v[29]*x6;
990       z[cval+5] += v[30]*x1  + v[31]*x2 + v[32]*x3 + v[33]*x4+ v[34]*x5 + v[35]*x6;
991       /* (strict upper triangular part of A)*x  */
992       z[6*i]   +=v[0]*x[cval]+v[6]*x[cval+1]+v[12]*x[cval+2]+v[18]*x[cval+3]+v[24]*x[cval+4]+v[30]*x[cval+5];
993       z[6*i+1] +=v[1]*x[cval]+v[7]*x[cval+1]+v[13]*x[cval+2]+v[19]*x[cval+3]+v[25]*x[cval+4]+v[31]*x[cval+5];
994       z[6*i+2] +=v[2]*x[cval]+v[8]*x[cval+1]+v[14]*x[cval+2]+v[20]*x[cval+3]+v[26]*x[cval+4]+v[32]*x[cval+5];
995       z[6*i+3] +=v[3]*x[cval]+v[9]*x[cval+1]+v[15]*x[cval+2]+v[21]*x[cval+3]+v[27]*x[cval+4]+v[33]*x[cval+5];
996       z[6*i+4] +=v[4]*x[cval]+v[10]*x[cval+1]+v[16]*x[cval+2]+v[22]*x[cval+3]+v[28]*x[cval+4]+v[34]*x[cval+5];
997       z[6*i+5] +=v[5]*x[cval]+v[11]*x[cval+1]+v[17]*x[cval+2]+v[23]*x[cval+3]+v[29]*x[cval+4]+v[35]*x[cval+5];
998       v        += 36;
999     }
1000     xb +=6; ai++;
1001   }
1002 
1003   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
1004   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
1005 
1006   ierr = PetscLogFlops(72.0*(a->nz*2.0 - nonzerorow));CHKERRQ(ierr);
1007   PetscFunctionReturn(0);
1008 }
1009 
1010 #undef __FUNCT__
1011 #define __FUNCT__ "MatMultAdd_SeqSBAIJ_7"
1012 PetscErrorCode MatMultAdd_SeqSBAIJ_7(Mat A,Vec xx,Vec yy,Vec zz)
1013 {
1014   Mat_SeqSBAIJ      *a = (Mat_SeqSBAIJ*)A->data;
1015   PetscScalar       *z,x1,x2,x3,x4,x5,x6,x7;
1016   const PetscScalar *x,*xb;
1017   const MatScalar   *v;
1018   PetscErrorCode    ierr;
1019   PetscInt          mbs = a->mbs,i,n,cval,j,jmin;
1020   const PetscInt    *aj=a->j,*ai=a->i,*ib;
1021   PetscInt          nonzerorow=0;
1022 
1023   PetscFunctionBegin;
1024   ierr = VecCopy(yy,zz);CHKERRQ(ierr);
1025   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
1026   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
1027 
1028   v  = a->a;
1029   xb = x;
1030 
1031   for (i=0; i<mbs; i++) {
1032     n           = ai[1] - ai[0]; /* length of i_th block row of A */
1033     x1          = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5=xb[4]; x6=xb[5]; x7=xb[6];
1034     ib          = aj + *ai;
1035     jmin        = 0;
1036     nonzerorow += (n>0);
1037     if (*ib == i) {     /* (diag of A)*x */
1038       z[7*i]   += v[0]*x1 + v[7]*x2 + v[14]*x3 + v[21]*x4+ v[28]*x5 + v[35]*x6+ v[42]*x7;
1039       z[7*i+1] += v[7]*x1 + v[8]*x2 + v[15]*x3 + v[22]*x4+ v[29]*x5 + v[36]*x6+ v[43]*x7;
1040       z[7*i+2] += v[14]*x1+ v[15]*x2 +v[16]*x3 + v[23]*x4+ v[30]*x5 + v[37]*x6+ v[44]*x7;
1041       z[7*i+3] += v[21]*x1+ v[22]*x2 +v[23]*x3 + v[24]*x4+ v[31]*x5 + v[38]*x6+ v[45]*x7;
1042       z[7*i+4] += v[28]*x1+ v[29]*x2 +v[30]*x3 + v[31]*x4+ v[32]*x5 + v[39]*x6+ v[46]*x7;
1043       z[7*i+5] += v[35]*x1+ v[36]*x2 +v[37]*x3 + v[38]*x4+ v[39]*x5 + v[40]*x6+ v[47]*x7;
1044       z[7*i+6] += v[42]*x1+ v[43]*x2 +v[44]*x3 + v[45]*x4+ v[46]*x5 + v[47]*x6+ v[48]*x7;
1045       v        += 49; jmin++;
1046     }
1047     PetscPrefetchBlock(ib+jmin+n,n,0,PETSC_PREFETCH_HINT_NTA); /* Indices for the next row (assumes same size as this one) */
1048     PetscPrefetchBlock(v+49*n,49*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
1049     for (j=jmin; j<n; j++) {
1050       /* (strict lower triangular part of A)*x  */
1051       cval       = ib[j]*7;
1052       z[cval]   += v[0]*x1  + v[1]*x2 + v[2]*x3 + v[3]*x4+ v[4]*x5 + v[5]*x6+ v[6]*x7;
1053       z[cval+1] += v[7]*x1  + v[8]*x2 + v[9]*x3 + v[10]*x4+ v[11]*x5 + v[12]*x6+ v[13]*x7;
1054       z[cval+2] += v[14]*x1  + v[15]*x2 + v[16]*x3 + v[17]*x4+ v[18]*x5 + v[19]*x6+ v[20]*x7;
1055       z[cval+3] += v[21]*x1  + v[22]*x2 + v[23]*x3 + v[24]*x4+ v[25]*x5 + v[26]*x6+ v[27]*x7;
1056       z[cval+4] += v[28]*x1  + v[29]*x2 + v[30]*x3 + v[31]*x4+ v[32]*x5 + v[33]*x6+ v[34]*x7;
1057       z[cval+5] += v[35]*x1  + v[36]*x2 + v[37]*x3 + v[38]*x4+ v[39]*x5 + v[40]*x6+ v[41]*x7;
1058       z[cval+6] += v[42]*x1  + v[43]*x2 + v[44]*x3 + v[45]*x4+ v[46]*x5 + v[47]*x6+ v[48]*x7;
1059       /* (strict upper triangular part of A)*x  */
1060       z[7*i]  +=v[0]*x[cval]+v[7]*x[cval+1]+v[14]*x[cval+2]+v[21]*x[cval+3]+v[28]*x[cval+4]+v[35]*x[cval+5]+v[42]*x[cval+6];
1061       z[7*i+1]+=v[1]*x[cval]+v[8]*x[cval+1]+v[15]*x[cval+2]+v[22]*x[cval+3]+v[29]*x[cval+4]+v[36]*x[cval+5]+v[43]*x[cval+6];
1062       z[7*i+2]+=v[2]*x[cval]+v[9]*x[cval+1]+v[16]*x[cval+2]+v[23]*x[cval+3]+v[30]*x[cval+4]+v[37]*x[cval+5]+v[44]*x[cval+6];
1063       z[7*i+3]+=v[3]*x[cval]+v[10]*x[cval+1]+v[17]*x[cval+2]+v[24]*x[cval+3]+v[31]*x[cval+4]+v[38]*x[cval+5]+v[45]*x[cval+6];
1064       z[7*i+4]+=v[4]*x[cval]+v[11]*x[cval+1]+v[18]*x[cval+2]+v[25]*x[cval+3]+v[32]*x[cval+4]+v[39]*x[cval+5]+v[46]*x[cval+6];
1065       z[7*i+5]+=v[5]*x[cval]+v[12]*x[cval+1]+v[19]*x[cval+2]+v[26]*x[cval+3]+v[33]*x[cval+4]+v[40]*x[cval+5]+v[47]*x[cval+6];
1066       z[7*i+6]+=v[6]*x[cval]+v[13]*x[cval+1]+v[20]*x[cval+2]+v[27]*x[cval+3]+v[34]*x[cval+4]+v[41]*x[cval+5]+v[48]*x[cval+6];
1067       v       += 49;
1068     }
1069     xb +=7; ai++;
1070   }
1071 
1072   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
1073   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
1074 
1075   ierr = PetscLogFlops(98.0*(a->nz*2.0 - nonzerorow));CHKERRQ(ierr);
1076   PetscFunctionReturn(0);
1077 }
1078 
1079 #undef __FUNCT__
1080 #define __FUNCT__ "MatMultAdd_SeqSBAIJ_N"
1081 PetscErrorCode MatMultAdd_SeqSBAIJ_N(Mat A,Vec xx,Vec yy,Vec zz)
1082 {
1083   Mat_SeqSBAIJ      *a = (Mat_SeqSBAIJ*)A->data;
1084   PetscScalar       *z,*z_ptr=0,*zb,*work,*workt;
1085   const PetscScalar *x,*x_ptr,*xb;
1086   const MatScalar   *v;
1087   PetscErrorCode    ierr;
1088   PetscInt          mbs = a->mbs,i,bs=A->rmap->bs,j,n,bs2=a->bs2,ncols,k;
1089   const PetscInt    *idx,*aj,*ii;
1090   PetscInt          nonzerorow=0;
1091 
1092   PetscFunctionBegin;
1093   ierr = VecCopy(yy,zz);CHKERRQ(ierr);
1094   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr); x_ptr=x;
1095   ierr = VecGetArray(zz,&z);CHKERRQ(ierr); z_ptr=z;
1096 
1097   aj = a->j;
1098   v  = a->a;
1099   ii = a->i;
1100 
1101   if (!a->mult_work) {
1102     ierr = PetscMalloc1(A->rmap->n+1,&a->mult_work);CHKERRQ(ierr);
1103   }
1104   work = a->mult_work;
1105 
1106 
1107   for (i=0; i<mbs; i++) {
1108     n           = ii[1] - ii[0]; ncols = n*bs;
1109     workt       = work; idx=aj+ii[0];
1110     nonzerorow += (n>0);
1111 
1112     /* upper triangular part */
1113     for (j=0; j<n; j++) {
1114       xb = x_ptr + bs*(*idx++);
1115       for (k=0; k<bs; k++) workt[k] = xb[k];
1116       workt += bs;
1117     }
1118     /* z(i*bs:(i+1)*bs-1) += A(i,:)*x */
1119     PetscKernel_w_gets_w_plus_Ar_times_v(bs,ncols,work,v,z);
1120 
1121     /* strict lower triangular part */
1122     idx = aj+ii[0];
1123     if (*idx == i) {
1124       ncols -= bs; v += bs2; idx++; n--;
1125     }
1126     if (ncols > 0) {
1127       workt = work;
1128       ierr  = PetscMemzero(workt,ncols*sizeof(PetscScalar));CHKERRQ(ierr);
1129       PetscKernel_w_gets_w_plus_trans_Ar_times_v(bs,ncols,x,v,workt);
1130       for (j=0; j<n; j++) {
1131         zb = z_ptr + bs*(*idx++);
1132         for (k=0; k<bs; k++) zb[k] += workt[k];
1133         workt += bs;
1134       }
1135     }
1136 
1137     x += bs; v += n*bs2; z += bs; ii++;
1138   }
1139 
1140   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
1141   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
1142 
1143   ierr = PetscLogFlops(2.0*(a->nz*2.0 - nonzerorow));CHKERRQ(ierr);
1144   PetscFunctionReturn(0);
1145 }
1146 
1147 #undef __FUNCT__
1148 #define __FUNCT__ "MatScale_SeqSBAIJ"
1149 PetscErrorCode MatScale_SeqSBAIJ(Mat inA,PetscScalar alpha)
1150 {
1151   Mat_SeqSBAIJ   *a     = (Mat_SeqSBAIJ*)inA->data;
1152   PetscScalar    oalpha = alpha;
1153   PetscErrorCode ierr;
1154   PetscBLASInt   one = 1,totalnz;
1155 
1156   PetscFunctionBegin;
1157   ierr = PetscBLASIntCast(a->bs2*a->nz,&totalnz);CHKERRQ(ierr);
1158   PetscStackCallBLAS("BLASscal",BLASscal_(&totalnz,&oalpha,a->a,&one));
1159   ierr = PetscLogFlops(totalnz);CHKERRQ(ierr);
1160   PetscFunctionReturn(0);
1161 }
1162 
1163 #undef __FUNCT__
1164 #define __FUNCT__ "MatNorm_SeqSBAIJ"
1165 PetscErrorCode MatNorm_SeqSBAIJ(Mat A,NormType type,PetscReal *norm)
1166 {
1167   Mat_SeqSBAIJ    *a       = (Mat_SeqSBAIJ*)A->data;
1168   const MatScalar *v       = a->a;
1169   PetscReal       sum_diag = 0.0, sum_off = 0.0, *sum;
1170   PetscInt        i,j,k,bs = A->rmap->bs,bs2=a->bs2,k1,mbs=a->mbs,jmin,jmax,nexti,ik,*jl,*il;
1171   PetscErrorCode  ierr;
1172   const PetscInt  *aj=a->j,*col;
1173 
1174   PetscFunctionBegin;
1175   if (type == NORM_FROBENIUS) {
1176     for (k=0; k<mbs; k++) {
1177       jmin = a->i[k]; jmax = a->i[k+1];
1178       col  = aj + jmin;
1179       if (*col == k) {         /* diagonal block */
1180         for (i=0; i<bs2; i++) {
1181           sum_diag += PetscRealPart(PetscConj(*v)*(*v)); v++;
1182         }
1183         jmin++;
1184       }
1185       for (j=jmin; j<jmax; j++) {  /* off-diagonal blocks */
1186         for (i=0; i<bs2; i++) {
1187           sum_off += PetscRealPart(PetscConj(*v)*(*v)); v++;
1188         }
1189       }
1190     }
1191     *norm = PetscSqrtReal(sum_diag + 2*sum_off);
1192   } else if (type == NORM_INFINITY || type == NORM_1) {  /* maximum row/column sum */
1193     ierr = PetscMalloc3(bs,&sum,mbs,&il,mbs,&jl);CHKERRQ(ierr);
1194     for (i=0; i<mbs; i++) jl[i] = mbs;
1195     il[0] = 0;
1196 
1197     *norm = 0.0;
1198     for (k=0; k<mbs; k++) { /* k_th block row */
1199       for (j=0; j<bs; j++) sum[j]=0.0;
1200       /*-- col sum --*/
1201       i = jl[k]; /* first |A(i,k)| to be added */
1202       /* jl[k]=i: first nozero element in row i for submatrix A(1:k,k:n) (active window)
1203                   at step k */
1204       while (i<mbs) {
1205         nexti = jl[i];  /* next block row to be added */
1206         ik    = il[i];  /* block index of A(i,k) in the array a */
1207         for (j=0; j<bs; j++) {
1208           v = a->a + ik*bs2 + j*bs;
1209           for (k1=0; k1<bs; k1++) {
1210             sum[j] += PetscAbsScalar(*v); v++;
1211           }
1212         }
1213         /* update il, jl */
1214         jmin = ik + 1; /* block index of array a: points to the next nonzero of A in row i */
1215         jmax = a->i[i+1];
1216         if (jmin < jmax) {
1217           il[i] = jmin;
1218           j     = a->j[jmin];
1219           jl[i] = jl[j]; jl[j]=i;
1220         }
1221         i = nexti;
1222       }
1223       /*-- row sum --*/
1224       jmin = a->i[k]; jmax = a->i[k+1];
1225       for (i=jmin; i<jmax; i++) {
1226         for (j=0; j<bs; j++) {
1227           v = a->a + i*bs2 + j;
1228           for (k1=0; k1<bs; k1++) {
1229             sum[j] += PetscAbsScalar(*v); v += bs;
1230           }
1231         }
1232       }
1233       /* add k_th block row to il, jl */
1234       col = aj+jmin;
1235       if (*col == k) jmin++;
1236       if (jmin < jmax) {
1237         il[k] = jmin;
1238         j = a->j[jmin]; jl[k] = jl[j]; jl[j] = k;
1239       }
1240       for (j=0; j<bs; j++) {
1241         if (sum[j] > *norm) *norm = sum[j];
1242       }
1243     }
1244     ierr = PetscFree3(sum,il,jl);CHKERRQ(ierr);
1245   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for this norm yet");
1246   PetscFunctionReturn(0);
1247 }
1248 
1249 #undef __FUNCT__
1250 #define __FUNCT__ "MatEqual_SeqSBAIJ"
1251 PetscErrorCode MatEqual_SeqSBAIJ(Mat A,Mat B,PetscBool * flg)
1252 {
1253   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ*)B->data;
1254   PetscErrorCode ierr;
1255 
1256   PetscFunctionBegin;
1257   /* If the  matrix/block dimensions are not equal, or no of nonzeros or shift */
1258   if ((A->rmap->N != B->rmap->N) || (A->cmap->n != B->cmap->n) || (A->rmap->bs != B->rmap->bs)|| (a->nz != b->nz)) {
1259     *flg = PETSC_FALSE;
1260     PetscFunctionReturn(0);
1261   }
1262 
1263   /* if the a->i are the same */
1264   ierr = PetscMemcmp(a->i,b->i,(a->mbs+1)*sizeof(PetscInt),flg);CHKERRQ(ierr);
1265   if (!*flg) PetscFunctionReturn(0);
1266 
1267   /* if a->j are the same */
1268   ierr = PetscMemcmp(a->j,b->j,(a->nz)*sizeof(PetscInt),flg);CHKERRQ(ierr);
1269   if (!*flg) PetscFunctionReturn(0);
1270 
1271   /* if a->a are the same */
1272   ierr = PetscMemcmp(a->a,b->a,(a->nz)*(A->rmap->bs)*(A->rmap->bs)*sizeof(PetscScalar),flg);CHKERRQ(ierr);
1273   PetscFunctionReturn(0);
1274 }
1275 
1276 #undef __FUNCT__
1277 #define __FUNCT__ "MatGetDiagonal_SeqSBAIJ"
1278 PetscErrorCode MatGetDiagonal_SeqSBAIJ(Mat A,Vec v)
1279 {
1280   Mat_SeqSBAIJ    *a = (Mat_SeqSBAIJ*)A->data;
1281   PetscErrorCode  ierr;
1282   PetscInt        i,j,k,row,bs,ambs,bs2;
1283   const PetscInt  *ai,*aj;
1284   PetscScalar     *x,zero = 0.0;
1285   const MatScalar *aa,*aa_j;
1286 
1287   PetscFunctionBegin;
1288   bs = A->rmap->bs;
1289   if (A->factortype && bs>1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix with bs>1");
1290 
1291   aa   = a->a;
1292   ambs = a->mbs;
1293 
1294   if (A->factortype == MAT_FACTOR_CHOLESKY || A->factortype == MAT_FACTOR_ICC) {
1295     PetscInt *diag=a->diag;
1296     aa   = a->a;
1297     ambs = a->mbs;
1298     ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1299     for (i=0; i<ambs; i++) x[i] = 1.0/aa[diag[i]];
1300     ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1301     PetscFunctionReturn(0);
1302   }
1303 
1304   ai   = a->i;
1305   aj   = a->j;
1306   bs2  = a->bs2;
1307   ierr = VecSet(v,zero);CHKERRQ(ierr);
1308   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1309   for (i=0; i<ambs; i++) {
1310     j=ai[i];
1311     if (aj[j] == i) {    /* if this is a diagonal element */
1312       row  = i*bs;
1313       aa_j = aa + j*bs2;
1314       for (k=0; k<bs2; k+=(bs+1),row++) x[row] = aa_j[k];
1315     }
1316   }
1317   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1318   PetscFunctionReturn(0);
1319 }
1320 
1321 #undef __FUNCT__
1322 #define __FUNCT__ "MatDiagonalScale_SeqSBAIJ"
1323 PetscErrorCode MatDiagonalScale_SeqSBAIJ(Mat A,Vec ll,Vec rr)
1324 {
1325   Mat_SeqSBAIJ      *a = (Mat_SeqSBAIJ*)A->data;
1326   PetscScalar       x;
1327   const PetscScalar *l,*li,*ri;
1328   MatScalar         *aa,*v;
1329   PetscErrorCode    ierr;
1330   PetscInt          i,j,k,lm,M,m,*ai,*aj,mbs,tmp,bs,bs2;
1331   PetscBool         flg;
1332 
1333   PetscFunctionBegin;
1334   if (ll != rr) {
1335     ierr = VecEqual(ll,rr,&flg);CHKERRQ(ierr);
1336     if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"For symmetric format, left and right scaling vectors must be same\n");
1337   }
1338   if (!ll) PetscFunctionReturn(0);
1339   ai  = a->i;
1340   aj  = a->j;
1341   aa  = a->a;
1342   m   = A->rmap->N;
1343   bs  = A->rmap->bs;
1344   mbs = a->mbs;
1345   bs2 = a->bs2;
1346 
1347   ierr = VecGetArrayRead(ll,&l);CHKERRQ(ierr);
1348   ierr = VecGetLocalSize(ll,&lm);CHKERRQ(ierr);
1349   if (lm != m) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Left scaling vector wrong length");
1350   for (i=0; i<mbs; i++) { /* for each block row */
1351     M  = ai[i+1] - ai[i];
1352     li = l + i*bs;
1353     v  = aa + bs2*ai[i];
1354     for (j=0; j<M; j++) { /* for each block */
1355       ri = l + bs*aj[ai[i]+j];
1356       for (k=0; k<bs; k++) {
1357         x = ri[k];
1358         for (tmp=0; tmp<bs; tmp++) (*v++) *= li[tmp]*x;
1359       }
1360     }
1361   }
1362   ierr = VecRestoreArrayRead(ll,&l);CHKERRQ(ierr);
1363   ierr = PetscLogFlops(2.0*a->nz);CHKERRQ(ierr);
1364   PetscFunctionReturn(0);
1365 }
1366 
1367 #undef __FUNCT__
1368 #define __FUNCT__ "MatGetInfo_SeqSBAIJ"
1369 PetscErrorCode MatGetInfo_SeqSBAIJ(Mat A,MatInfoType flag,MatInfo *info)
1370 {
1371   Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data;
1372 
1373   PetscFunctionBegin;
1374   info->block_size   = a->bs2;
1375   info->nz_allocated = a->bs2*a->maxnz;   /*num. of nonzeros in upper triangular part */
1376   info->nz_used      = a->bs2*a->nz;   /*num. of nonzeros in upper triangular part */
1377   info->nz_unneeded  = (double)(info->nz_allocated - info->nz_used);
1378   info->assemblies   = A->num_ass;
1379   info->mallocs      = A->info.mallocs;
1380   info->memory       = ((PetscObject)A)->mem;
1381   if (A->factortype) {
1382     info->fill_ratio_given  = A->info.fill_ratio_given;
1383     info->fill_ratio_needed = A->info.fill_ratio_needed;
1384     info->factor_mallocs    = A->info.factor_mallocs;
1385   } else {
1386     info->fill_ratio_given  = 0;
1387     info->fill_ratio_needed = 0;
1388     info->factor_mallocs    = 0;
1389   }
1390   PetscFunctionReturn(0);
1391 }
1392 
1393 
1394 #undef __FUNCT__
1395 #define __FUNCT__ "MatZeroEntries_SeqSBAIJ"
1396 PetscErrorCode MatZeroEntries_SeqSBAIJ(Mat A)
1397 {
1398   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
1399   PetscErrorCode ierr;
1400 
1401   PetscFunctionBegin;
1402   ierr = PetscMemzero(a->a,a->bs2*a->i[a->mbs]*sizeof(MatScalar));CHKERRQ(ierr);
1403   PetscFunctionReturn(0);
1404 }
1405 
1406 #undef __FUNCT__
1407 #define __FUNCT__ "MatGetRowMaxAbs_SeqSBAIJ"
1408 /*
1409    This code does not work since it only checks the upper triangular part of
1410   the matrix. Hence it is not listed in the function table.
1411 */
1412 PetscErrorCode MatGetRowMaxAbs_SeqSBAIJ(Mat A,Vec v,PetscInt idx[])
1413 {
1414   Mat_SeqSBAIJ    *a = (Mat_SeqSBAIJ*)A->data;
1415   PetscErrorCode  ierr;
1416   PetscInt        i,j,n,row,col,bs,mbs;
1417   const PetscInt  *ai,*aj;
1418   PetscReal       atmp;
1419   const MatScalar *aa;
1420   PetscScalar     *x;
1421   PetscInt        ncols,brow,bcol,krow,kcol;
1422 
1423   PetscFunctionBegin;
1424   if (idx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Send email to petsc-maint@mcs.anl.gov");
1425   if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1426   bs  = A->rmap->bs;
1427   aa  = a->a;
1428   ai  = a->i;
1429   aj  = a->j;
1430   mbs = a->mbs;
1431 
1432   ierr = VecSet(v,0.0);CHKERRQ(ierr);
1433   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1434   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
1435   if (n != A->rmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
1436   for (i=0; i<mbs; i++) {
1437     ncols = ai[1] - ai[0]; ai++;
1438     brow  = bs*i;
1439     for (j=0; j<ncols; j++) {
1440       bcol = bs*(*aj);
1441       for (kcol=0; kcol<bs; kcol++) {
1442         col = bcol + kcol;      /* col index */
1443         for (krow=0; krow<bs; krow++) {
1444           atmp = PetscAbsScalar(*aa); aa++;
1445           row  = brow + krow;   /* row index */
1446           if (PetscRealPart(x[row]) < atmp) x[row] = atmp;
1447           if (*aj > i && PetscRealPart(x[col]) < atmp) x[col] = atmp;
1448         }
1449       }
1450       aj++;
1451     }
1452   }
1453   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
1454   PetscFunctionReturn(0);
1455 }
1456