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