xref: /petsc/src/mat/impls/baij/seq/baij2.c (revision 6831982abb6453c2f3e25efecb5d0951d333371e)
1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER
2*6831982aSBarry Smith static char vcid[] = "$Id: baij2.c,v 1.50 1999/10/13 20:37:28 bsmith Exp bsmith $";
3cac129eeSSatish Balay #endif
4cac129eeSSatish Balay 
52d61bbb3SSatish Balay #include "sys.h"
670f55243SBarry Smith #include "src/mat/impls/baij/seq/baij.h"
72d61bbb3SSatish Balay #include "src/vec/vecimpl.h"
82d61bbb3SSatish Balay #include "src/inline/spops.h"
93f1db9ecSBarry Smith #include "src/inline/ilu.h"
10eeb667a2SSatish Balay #include "bitarray.h"
11cac129eeSSatish Balay 
125615d1e5SSatish Balay #undef __FUNC__
13d4bb536fSBarry Smith #define __FUNC__ "MatIncreaseOverlap_SeqBAIJ"
14736121d4SSatish Balay int MatIncreaseOverlap_SeqBAIJ(Mat A,int is_max,IS *is,int ov)
15a3192f15SSatish Balay {
16a3192f15SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
17218c64b6SSatish Balay   int         row, i,j,k,l,m,n, *idx,ierr, *nidx, isz, val, ival;
18218c64b6SSatish Balay   int         start, end, *ai, *aj,bs,*nidx2;
19*6831982aSBarry Smith   BTPetsc     table;
20a3192f15SSatish Balay 
213a40ed3dSBarry Smith   PetscFunctionBegin;
22a3192f15SSatish Balay   m     = a->mbs;
23a3192f15SSatish Balay   ai    = a->i;
24a3192f15SSatish Balay   aj    = a->j;
25218c64b6SSatish Balay   bs    = a->bs;
26a3192f15SSatish Balay 
27a8c6a408SBarry Smith   if (ov < 0)  SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative overlap specified");
28a3192f15SSatish Balay 
29*6831982aSBarry Smith   ierr  = PetscBTCreate(m,table);CHKERRQ(ierr);
30a3192f15SSatish Balay   nidx  = (int *) PetscMalloc((m+1)*sizeof(int));CHKPTRQ(nidx);
31218c64b6SSatish Balay   nidx2 = (int *)PetscMalloc((a->m+1)*sizeof(int));CHKPTRQ(nidx2);
32a3192f15SSatish Balay 
33a3192f15SSatish Balay   for ( i=0; i<is_max; i++ ) {
34a3192f15SSatish Balay     /* Initialise the two local arrays */
35a3192f15SSatish Balay     isz  = 0;
36*6831982aSBarry Smith     ierr = PetscBTMemzero(m,table);CHKERRQ(ierr);
37a3192f15SSatish Balay 
38a3192f15SSatish Balay     /* Extract the indices, assume there can be duplicate entries */
39a3192f15SSatish Balay     ierr = ISGetIndices(is[i],&idx);CHKERRQ(ierr);
40a3192f15SSatish Balay     ierr = ISGetSize(is[i],&n);CHKERRQ(ierr);
41a3192f15SSatish Balay 
42a3192f15SSatish Balay     /* Enter these into the temp arrays i.e mark table[row], enter row into new index */
43a3192f15SSatish Balay     for ( j=0; j<n ; ++j){
44218c64b6SSatish Balay       ival = idx[j]/bs; /* convert the indices into block indices */
45a8c6a408SBarry Smith       if (ival>m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"index greater than mat-dim");
46*6831982aSBarry Smith       if(!PetscBTLoopupSet(table, ival)) { nidx[isz++] = ival;}
47a3192f15SSatish Balay     }
48a3192f15SSatish Balay     ierr = ISRestoreIndices(is[i],&idx);CHKERRQ(ierr);
49a3192f15SSatish Balay     ierr = ISDestroy(is[i]);CHKERRQ(ierr);
50a3192f15SSatish Balay 
51a3192f15SSatish Balay     k = 0;
52a3192f15SSatish Balay     for ( j=0; j<ov; j++){ /* for each overlap*/
53a3192f15SSatish Balay       n = isz;
54a3192f15SSatish Balay       for ( ; k<n ; k++){ /* do only those rows in nidx[k], which are not done yet */
55a3192f15SSatish Balay         row   = nidx[k];
56a3192f15SSatish Balay         start = ai[row];
57a3192f15SSatish Balay         end   = ai[row+1];
58a3192f15SSatish Balay         for ( l = start; l<end ; l++){
59a3192f15SSatish Balay           val = aj[l];
60*6831982aSBarry Smith           if (!PetscBTLoopupSet(table,val)) {nidx[isz++] = val;}
61a3192f15SSatish Balay         }
62a3192f15SSatish Balay       }
63a3192f15SSatish Balay     }
64218c64b6SSatish Balay     /* expand the Index Set */
65218c64b6SSatish Balay     for (j=0; j<isz; j++ ) {
66218c64b6SSatish Balay       for (k=0; k<bs; k++ )
67218c64b6SSatish Balay         nidx2[j*bs+k] = nidx[j]*bs+k;
68218c64b6SSatish Balay     }
69029af93fSBarry Smith     ierr = ISCreateGeneral(PETSC_COMM_SELF, isz*bs, nidx2, (is+i));CHKERRQ(ierr);
70a3192f15SSatish Balay   }
71*6831982aSBarry Smith   ierr = PetscBTDestroy(table);CHKERRQ(ierr);
72606d414cSSatish Balay   ierr = PetscFree(nidx);CHKERRQ(ierr);
73606d414cSSatish Balay   ierr = PetscFree(nidx2);CHKERRQ(ierr);
743a40ed3dSBarry Smith   PetscFunctionReturn(0);
75a3192f15SSatish Balay }
761c351548SSatish Balay 
775615d1e5SSatish Balay #undef __FUNC__
785615d1e5SSatish Balay #define __FUNC__ "MatGetSubMatrix_SeqBAIJ_Private"
797b2a1423SBarry Smith int MatGetSubMatrix_SeqBAIJ_Private(Mat A,IS isrow,IS iscol,int cs,MatReuse scall,Mat *B)
80736121d4SSatish Balay {
81736121d4SSatish Balay   Mat_SeqBAIJ  *a = (Mat_SeqBAIJ *) A->data,*c;
8240011551SBarry Smith   int          *smap, i, k, kstart, kend, ierr, oldcols = a->nbs,*lens;
83218c64b6SSatish Balay   int          row,mat_i,*mat_j,tcol,*mat_ilen;
84736121d4SSatish Balay   int          *irow, *icol, nrows, ncols,*ssmap,bs=a->bs, bs2=a->bs2;
85218c64b6SSatish Balay   int          *aj = a->j, *ai = a->i;
863f1db9ecSBarry Smith   MatScalar    *mat_a;
87736121d4SSatish Balay   Mat          C;
880f5bd95cSBarry Smith   PetscTruth   flag;
89736121d4SSatish Balay 
903a40ed3dSBarry Smith   PetscFunctionBegin;
91888f2ed8SSatish Balay   ierr = ISSorted(iscol,(PetscTruth*)&i);CHKERRQ(ierr);
92a8c6a408SBarry Smith   if (!i) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"IS is not sorted");
93736121d4SSatish Balay 
94736121d4SSatish Balay   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
95218c64b6SSatish Balay   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
96736121d4SSatish Balay   ierr = ISGetSize(isrow,&nrows);CHKERRQ(ierr);
97736121d4SSatish Balay   ierr = ISGetSize(iscol,&ncols);CHKERRQ(ierr);
98736121d4SSatish Balay 
99736121d4SSatish Balay   smap  = (int *) PetscMalloc((1+oldcols)*sizeof(int));CHKPTRQ(smap);
100736121d4SSatish Balay   ssmap = smap;
101736121d4SSatish Balay   lens  = (int *) PetscMalloc((1+nrows)*sizeof(int));CHKPTRQ(lens);
102549d3d68SSatish Balay   ierr  = PetscMemzero(smap,oldcols*sizeof(int));CHKERRQ(ierr);
103736121d4SSatish Balay   for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1;
104736121d4SSatish Balay   /* determine lens of each row */
105736121d4SSatish Balay   for (i=0; i<nrows; i++) {
106736121d4SSatish Balay     kstart  = ai[irow[i]];
107736121d4SSatish Balay     kend    = kstart + a->ilen[irow[i]];
108736121d4SSatish Balay     lens[i] = 0;
109736121d4SSatish Balay       for ( k=kstart; k<kend; k++ ) {
110736121d4SSatish Balay         if (ssmap[aj[k]]) {
111736121d4SSatish Balay           lens[i]++;
112736121d4SSatish Balay         }
113736121d4SSatish Balay       }
114736121d4SSatish Balay     }
115736121d4SSatish Balay   /* Create and fill new matrix */
116736121d4SSatish Balay   if (scall == MAT_REUSE_MATRIX) {
117736121d4SSatish Balay     c = (Mat_SeqBAIJ *)((*B)->data);
118736121d4SSatish Balay 
119a8c6a408SBarry Smith     if (c->mbs!=nrows || c->nbs!=ncols || c->bs!=bs) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Submatrix wrong size");
1200f5bd95cSBarry Smith     ierr = PetscMemcmp(c->ilen,lens, c->mbs *sizeof(int),&flag);CHKERRQ(ierr);
1210f5bd95cSBarry Smith     if (flag == PETSC_FALSE) {
122a8c6a408SBarry Smith       SETERRQ(PETSC_ERR_ARG_SIZ,0,"Cannot reuse matrix. wrong no of nonzeros");
123736121d4SSatish Balay     }
124549d3d68SSatish Balay     ierr = PetscMemzero(c->ilen,c->mbs*sizeof(int));CHKERRQ(ierr);
125736121d4SSatish Balay     C = *B;
1263a40ed3dSBarry Smith   } else {
127736121d4SSatish Balay     ierr = MatCreateSeqBAIJ(A->comm,bs,nrows*bs,ncols*bs,0,lens,&C);CHKERRQ(ierr);
128736121d4SSatish Balay   }
129736121d4SSatish Balay   c = (Mat_SeqBAIJ *)(C->data);
130736121d4SSatish Balay   for (i=0; i<nrows; i++) {
131736121d4SSatish Balay     row    = irow[i];
132736121d4SSatish Balay     kstart = ai[row];
133736121d4SSatish Balay     kend   = kstart + a->ilen[row];
134736121d4SSatish Balay     mat_i  = c->i[i];
135736121d4SSatish Balay     mat_j  = c->j + mat_i;
136218c64b6SSatish Balay     mat_a  = c->a + mat_i*bs2;
137736121d4SSatish Balay     mat_ilen = c->ilen + i;
138736121d4SSatish Balay     for ( k=kstart; k<kend; k++ ) {
139736121d4SSatish Balay       if ((tcol=ssmap[a->j[k]])) {
140736121d4SSatish Balay         *mat_j++ = tcol - 1;
141549d3d68SSatish Balay         ierr     = PetscMemcpy(mat_a,a->a+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
142549d3d68SSatish Balay         mat_a   += bs2;
143736121d4SSatish Balay         (*mat_ilen)++;
144736121d4SSatish Balay       }
145736121d4SSatish Balay     }
146736121d4SSatish Balay   }
147218c64b6SSatish Balay 
148736121d4SSatish Balay   /* Free work space */
149736121d4SSatish Balay   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
150606d414cSSatish Balay   ierr = PetscFree(smap);CHKERRQ(ierr);
151606d414cSSatish Balay   ierr = PetscFree(lens);CHKERRQ(ierr);
1526d4a8577SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1536d4a8577SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
154736121d4SSatish Balay 
155736121d4SSatish Balay   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
156736121d4SSatish Balay   *B = C;
1573a40ed3dSBarry Smith   PetscFunctionReturn(0);
158736121d4SSatish Balay }
159736121d4SSatish Balay 
1605615d1e5SSatish Balay #undef __FUNC__
1615615d1e5SSatish Balay #define __FUNC__ "MatGetSubMatrix_SeqBAIJ"
1627b2a1423SBarry Smith int MatGetSubMatrix_SeqBAIJ(Mat A,IS isrow,IS iscol,int cs,MatReuse scall,Mat *B)
163218c64b6SSatish Balay {
164218c64b6SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
165218c64b6SSatish Balay   IS          is1,is2;
166218c64b6SSatish Balay   int         *vary,*iary,*irow,*icol,nrows,ncols,i,ierr,bs=a->bs,count;
167218c64b6SSatish Balay 
1683a40ed3dSBarry Smith   PetscFunctionBegin;
169218c64b6SSatish Balay   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
170218c64b6SSatish Balay   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
171218c64b6SSatish Balay   ierr = ISGetSize(isrow,&nrows);CHKERRQ(ierr);
172218c64b6SSatish Balay   ierr = ISGetSize(iscol,&ncols);CHKERRQ(ierr);
173218c64b6SSatish Balay 
174218c64b6SSatish Balay   /* Verify if the indices corespond to each element in a block
175218c64b6SSatish Balay    and form the IS with compressed IS */
176218c64b6SSatish Balay   vary = (int *) PetscMalloc(2*(a->mbs+1)*sizeof(int));CHKPTRQ(vary);
177218c64b6SSatish Balay   iary = vary + a->mbs;
178549d3d68SSatish Balay   ierr = PetscMemzero(vary,(a->mbs)*sizeof(int));CHKERRQ(ierr);
179218c64b6SSatish Balay   for ( i=0; i<nrows; i++) vary[irow[i]/bs]++;
180218c64b6SSatish Balay   count = 0;
181218c64b6SSatish Balay   for (i=0; i<a->mbs; i++) {
1826a6a5d1dSBarry Smith     if (vary[i]!=0 && vary[i]!=bs) SETERRA(1,0,"Index set does not match blocks");
183218c64b6SSatish Balay     if (vary[i]==bs) iary[count++] = i;
184218c64b6SSatish Balay   }
185029af93fSBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF, count, iary,&is1);CHKERRQ(ierr);
186218c64b6SSatish Balay 
187549d3d68SSatish Balay   ierr = PetscMemzero(vary,(a->mbs)*sizeof(int));CHKERRQ(ierr);
188218c64b6SSatish Balay   for ( i=0; i<ncols; i++) vary[icol[i]/bs]++;
189218c64b6SSatish Balay   count = 0;
190218c64b6SSatish Balay   for (i=0; i<a->mbs; i++) {
191e3372554SBarry Smith     if (vary[i]!=0 && vary[i]!=bs) SETERRA(1,0,"MatGetSubmatrices_SeqBAIJ:");
192218c64b6SSatish Balay     if (vary[i]==bs) iary[count++] = i;
193218c64b6SSatish Balay   }
194029af93fSBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF, count, iary,&is2);CHKERRQ(ierr);
195218c64b6SSatish Balay   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
196218c64b6SSatish Balay   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
197606d414cSSatish Balay   ierr = PetscFree(vary);CHKERRQ(ierr);
198218c64b6SSatish Balay 
1996a6a5d1dSBarry Smith   ierr = MatGetSubMatrix_SeqBAIJ_Private(A,is1,is2,cs,scall,B);CHKERRQ(ierr);
200218c64b6SSatish Balay   ISDestroy(is1);
201218c64b6SSatish Balay   ISDestroy(is2);
2023a40ed3dSBarry Smith   PetscFunctionReturn(0);
203218c64b6SSatish Balay }
204218c64b6SSatish Balay 
2055615d1e5SSatish Balay #undef __FUNC__
2065615d1e5SSatish Balay #define __FUNC__ "MatGetSubMatrices_SeqBAIJ"
2077b2a1423SBarry Smith int MatGetSubMatrices_SeqBAIJ(Mat A,int n, IS *irow,IS *icol,MatReuse scall,Mat **B)
208736121d4SSatish Balay {
209736121d4SSatish Balay   int ierr,i;
210736121d4SSatish Balay 
2113a40ed3dSBarry Smith   PetscFunctionBegin;
212736121d4SSatish Balay   if (scall == MAT_INITIAL_MATRIX) {
213736121d4SSatish Balay     *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) );CHKPTRQ(*B);
214736121d4SSatish Balay   }
215736121d4SSatish Balay 
216736121d4SSatish Balay   for ( i=0; i<n; i++ ) {
2176a6a5d1dSBarry Smith     ierr = MatGetSubMatrix_SeqBAIJ(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr);
218736121d4SSatish Balay   }
2193a40ed3dSBarry Smith   PetscFunctionReturn(0);
220736121d4SSatish Balay }
221218c64b6SSatish Balay 
222218c64b6SSatish Balay 
2232d61bbb3SSatish Balay /* -------------------------------------------------------*/
2242d61bbb3SSatish Balay /* Should check that shapes of vectors and matrices match */
2252d61bbb3SSatish Balay /* -------------------------------------------------------*/
226eadb2fb4SBarry Smith #include "pinclude/blaslapack.h"
2272d61bbb3SSatish Balay 
2282d61bbb3SSatish Balay #undef __FUNC__
2292d61bbb3SSatish Balay #define __FUNC__ "MatMult_SeqBAIJ_1"
2302d61bbb3SSatish Balay int MatMult_SeqBAIJ_1(Mat A,Vec xx,Vec zz)
2312d61bbb3SSatish Balay {
2322d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
2333f1db9ecSBarry Smith   Scalar          *x,*z,sum;
2343f1db9ecSBarry Smith   MatScalar       *v;
235e1311b90SBarry Smith   int             mbs=a->mbs,i,*idx,*ii,n,ierr;
2362d61bbb3SSatish Balay 
2372d61bbb3SSatish Balay   PetscFunctionBegin;
238e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
239e1311b90SBarry Smith   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
2402d61bbb3SSatish Balay 
2412d61bbb3SSatish Balay   idx   = a->j;
2422d61bbb3SSatish Balay   v     = a->a;
2432d61bbb3SSatish Balay   ii    = a->i;
2442d61bbb3SSatish Balay 
2452d61bbb3SSatish Balay   for ( i=0; i<mbs; i++ ) {
2462d61bbb3SSatish Balay     n    = ii[1] - ii[0]; ii++;
2472d61bbb3SSatish Balay     sum  = 0.0;
2482d61bbb3SSatish Balay     while (n--) sum += *v++ * x[*idx++];
2492d61bbb3SSatish Balay     z[i] = sum;
2502d61bbb3SSatish Balay   }
251e1311b90SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
252e1311b90SBarry Smith   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
2532d61bbb3SSatish Balay   PLogFlops(2*a->nz - a->m);
2542d61bbb3SSatish Balay   PetscFunctionReturn(0);
2552d61bbb3SSatish Balay }
2562d61bbb3SSatish Balay 
2572d61bbb3SSatish Balay #undef __FUNC__
2582d61bbb3SSatish Balay #define __FUNC__ "MatMult_SeqBAIJ_2"
2592d61bbb3SSatish Balay int MatMult_SeqBAIJ_2(Mat A,Vec xx,Vec zz)
2602d61bbb3SSatish Balay {
2612d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
2623f1db9ecSBarry Smith   Scalar          *x,*z,*xb,sum1,sum2;
263e1311b90SBarry Smith   Scalar          x1,x2;
2643f1db9ecSBarry Smith   MatScalar       *v;
265e1311b90SBarry Smith   int             ierr,mbs=a->mbs,i,*idx,*ii,j,n;
2662d61bbb3SSatish Balay 
2672d61bbb3SSatish Balay   PetscFunctionBegin;
268e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
269e1311b90SBarry Smith   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
2702d61bbb3SSatish Balay 
2712d61bbb3SSatish Balay   idx   = a->j;
2722d61bbb3SSatish Balay   v     = a->a;
2732d61bbb3SSatish Balay   ii    = a->i;
2742d61bbb3SSatish Balay 
2752d61bbb3SSatish Balay   for ( i=0; i<mbs; i++ ) {
2762d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
2772d61bbb3SSatish Balay     sum1 = 0.0; sum2 = 0.0;
2782d61bbb3SSatish Balay     for ( j=0; j<n; j++ ) {
2792d61bbb3SSatish Balay       xb = x + 2*(*idx++); x1 = xb[0]; x2 = xb[1];
2802d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[2]*x2;
2812d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[3]*x2;
2822d61bbb3SSatish Balay       v += 4;
2832d61bbb3SSatish Balay     }
2842d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2;
2852d61bbb3SSatish Balay     z += 2;
2862d61bbb3SSatish Balay   }
287e1311b90SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
288e1311b90SBarry Smith   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
28915091d37SBarry Smith   PLogFlops(8*a->nz - a->m);
2902d61bbb3SSatish Balay   PetscFunctionReturn(0);
2912d61bbb3SSatish Balay }
2922d61bbb3SSatish Balay 
2932d61bbb3SSatish Balay #undef __FUNC__
2942d61bbb3SSatish Balay #define __FUNC__ "MatMult_SeqBAIJ_3"
2952d61bbb3SSatish Balay int MatMult_SeqBAIJ_3(Mat A,Vec xx,Vec zz)
2962d61bbb3SSatish Balay {
2972d61bbb3SSatish Balay   Mat_SeqBAIJ  *a = (Mat_SeqBAIJ *) A->data;
2983f1db9ecSBarry Smith   Scalar       *x,*z,*xb,sum1,sum2,sum3,x1,x2,x3;
2993f1db9ecSBarry Smith   MatScalar    *v;
300e1311b90SBarry Smith   int          ierr,mbs=a->mbs,i,*idx,*ii,j,n;
3012d61bbb3SSatish Balay 
302aa482453SBarry Smith #if defined(PETSC_HAVE_PRAGMA_DISJOINT)
303fee21e36SBarry Smith #pragma disjoint(*v,*z,*xb)
304fee21e36SBarry Smith #endif
305fee21e36SBarry Smith 
3062d61bbb3SSatish Balay   PetscFunctionBegin;
307e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
308e1311b90SBarry Smith   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
3092d61bbb3SSatish Balay 
3102d61bbb3SSatish Balay   idx   = a->j;
3112d61bbb3SSatish Balay   v     = a->a;
3122d61bbb3SSatish Balay   ii    = a->i;
3132d61bbb3SSatish Balay 
3142d61bbb3SSatish Balay   for ( i=0; i<mbs; i++ ) {
3152d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
3162d61bbb3SSatish Balay     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0;
3172d61bbb3SSatish Balay     for ( j=0; j<n; j++ ) {
3182d61bbb3SSatish Balay       xb = x + 3*(*idx++); x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
3192d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[3]*x2 + v[6]*x3;
3202d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[4]*x2 + v[7]*x3;
3212d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[5]*x2 + v[8]*x3;
3222d61bbb3SSatish Balay       v += 9;
3232d61bbb3SSatish Balay     }
3242d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3;
3252d61bbb3SSatish Balay     z += 3;
3262d61bbb3SSatish Balay   }
327e1311b90SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
328e1311b90SBarry Smith   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
3292d61bbb3SSatish Balay   PLogFlops(18*a->nz - a->m);
3302d61bbb3SSatish Balay   PetscFunctionReturn(0);
3312d61bbb3SSatish Balay }
3322d61bbb3SSatish Balay 
3332d61bbb3SSatish Balay #undef __FUNC__
3342d61bbb3SSatish Balay #define __FUNC__ "MatMult_SeqBAIJ_4"
3352d61bbb3SSatish Balay int MatMult_SeqBAIJ_4(Mat A,Vec xx,Vec zz)
3362d61bbb3SSatish Balay {
3372d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
3383f1db9ecSBarry Smith   Scalar          *x,*z,*xb,sum1,sum2,sum3,sum4,x1,x2,x3,x4;
3393f1db9ecSBarry Smith   MatScalar       *v;
340e1311b90SBarry Smith   int             ierr,mbs=a->mbs,i,*idx,*ii,j,n;
3412d61bbb3SSatish Balay 
3422d61bbb3SSatish Balay   PetscFunctionBegin;
343e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
344e1311b90SBarry Smith   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
3452d61bbb3SSatish Balay 
3462d61bbb3SSatish Balay   idx   = a->j;
3472d61bbb3SSatish Balay   v     = a->a;
3482d61bbb3SSatish Balay   ii    = a->i;
3492d61bbb3SSatish Balay 
3502d61bbb3SSatish Balay   for ( i=0; i<mbs; i++ ) {
3512d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
3522d61bbb3SSatish Balay     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0;
3532d61bbb3SSatish Balay     for ( j=0; j<n; j++ ) {
3542d61bbb3SSatish Balay       xb = x + 4*(*idx++);
3552d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
3562d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[4]*x2 + v[8]*x3   + v[12]*x4;
3572d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[5]*x2 + v[9]*x3   + v[13]*x4;
3582d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[6]*x2 + v[10]*x3  + v[14]*x4;
3592d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[7]*x2 + v[11]*x3  + v[15]*x4;
3602d61bbb3SSatish Balay       v += 16;
3612d61bbb3SSatish Balay     }
3622d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4;
3632d61bbb3SSatish Balay     z += 4;
3642d61bbb3SSatish Balay   }
365e1311b90SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
366e1311b90SBarry Smith   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
3672d61bbb3SSatish Balay   PLogFlops(32*a->nz - a->m);
3682d61bbb3SSatish Balay   PetscFunctionReturn(0);
3692d61bbb3SSatish Balay }
3702d61bbb3SSatish Balay 
3712d61bbb3SSatish Balay #undef __FUNC__
3722d61bbb3SSatish Balay #define __FUNC__ "MatMult_SeqBAIJ_5"
3732d61bbb3SSatish Balay int MatMult_SeqBAIJ_5(Mat A,Vec xx,Vec zz)
3742d61bbb3SSatish Balay {
3752d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
3763f1db9ecSBarry Smith   Scalar          sum1,sum2,sum3,sum4,sum5,x1,x2,x3,x4,x5,*xb,*z,*x;
3773f1db9ecSBarry Smith   MatScalar       *v;
378e1311b90SBarry Smith   int             ierr,mbs=a->mbs,i,*idx,*ii,j,n;
3792d61bbb3SSatish Balay 
380433994e6SBarry Smith   PetscFunctionBegin;
381e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
382e1311b90SBarry Smith   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
3832d61bbb3SSatish Balay 
3842d61bbb3SSatish Balay   idx   = a->j;
3852d61bbb3SSatish Balay   v     = a->a;
3862d61bbb3SSatish Balay   ii    = a->i;
3872d61bbb3SSatish Balay 
3882d61bbb3SSatish Balay   for ( i=0; i<mbs; i++ ) {
3892d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
3902d61bbb3SSatish Balay     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0;
3912d61bbb3SSatish Balay     for ( j=0; j<n; j++ ) {
3922d61bbb3SSatish Balay       xb = x + 5*(*idx++);
3932d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4];
3942d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[5]*x2 + v[10]*x3  + v[15]*x4 + v[20]*x5;
3952d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[6]*x2 + v[11]*x3  + v[16]*x4 + v[21]*x5;
3962d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[7]*x2 + v[12]*x3  + v[17]*x4 + v[22]*x5;
3972d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[8]*x2 + v[13]*x3  + v[18]*x4 + v[23]*x5;
3982d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[9]*x2 + v[14]*x3  + v[19]*x4 + v[24]*x5;
3992d61bbb3SSatish Balay       v += 25;
4002d61bbb3SSatish Balay     }
4012d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5;
4022d61bbb3SSatish Balay     z += 5;
4032d61bbb3SSatish Balay   }
404e1311b90SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
405e1311b90SBarry Smith   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
4062d61bbb3SSatish Balay   PLogFlops(50*a->nz - a->m);
4072d61bbb3SSatish Balay   PetscFunctionReturn(0);
4082d61bbb3SSatish Balay }
4092d61bbb3SSatish Balay 
41015091d37SBarry Smith 
41115091d37SBarry Smith #undef __FUNC__
41215091d37SBarry Smith #define __FUNC__ "MatMult_SeqBAIJ_6"
41315091d37SBarry Smith int MatMult_SeqBAIJ_6(Mat A,Vec xx,Vec zz)
41415091d37SBarry Smith {
41515091d37SBarry Smith   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
41615091d37SBarry Smith   Scalar          *x,*z,*xb,sum1,sum2,sum3,sum4,sum5,sum6;
41715091d37SBarry Smith   Scalar          x1,x2,x3,x4,x5,x6;
41815091d37SBarry Smith   MatScalar       *v;
41915091d37SBarry Smith   int             ierr,mbs=a->mbs,i,*idx,*ii,j,n;
42015091d37SBarry Smith 
421433994e6SBarry Smith   PetscFunctionBegin;
42215091d37SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
42315091d37SBarry Smith   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
42415091d37SBarry Smith 
42515091d37SBarry Smith   idx   = a->j;
42615091d37SBarry Smith   v     = a->a;
42715091d37SBarry Smith   ii    = a->i;
42815091d37SBarry Smith 
42915091d37SBarry Smith   for ( i=0; i<mbs; i++ ) {
43015091d37SBarry Smith     n  = ii[1] - ii[0]; ii++;
43115091d37SBarry Smith     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0;
43215091d37SBarry Smith     for ( j=0; j<n; j++ ) {
43315091d37SBarry Smith       xb = x + 6*(*idx++);
43415091d37SBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5];
43515091d37SBarry Smith       sum1 += v[0]*x1 + v[6]*x2  + v[12]*x3  + v[18]*x4 + v[24]*x5 + v[30]*x6;
43615091d37SBarry Smith       sum2 += v[1]*x1 + v[7]*x2  + v[13]*x3  + v[19]*x4 + v[25]*x5 + v[31]*x6;
43715091d37SBarry Smith       sum3 += v[2]*x1 + v[8]*x2  + v[14]*x3  + v[20]*x4 + v[26]*x5 + v[32]*x6;
43815091d37SBarry Smith       sum4 += v[3]*x1 + v[9]*x2  + v[15]*x3  + v[21]*x4 + v[27]*x5 + v[33]*x6;
43915091d37SBarry Smith       sum5 += v[4]*x1 + v[10]*x2 + v[16]*x3  + v[22]*x4 + v[28]*x5 + v[34]*x6;
44015091d37SBarry Smith       sum6 += v[5]*x1 + v[11]*x2 + v[17]*x3  + v[23]*x4 + v[29]*x5 + v[35]*x6;
44115091d37SBarry Smith       v += 36;
44215091d37SBarry Smith     }
44315091d37SBarry Smith     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6;
44415091d37SBarry Smith     z += 6;
44515091d37SBarry Smith   }
44615091d37SBarry Smith 
44715091d37SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
44815091d37SBarry Smith   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
44915091d37SBarry Smith   PLogFlops(72*a->nz - a->m);
45015091d37SBarry Smith   PetscFunctionReturn(0);
45115091d37SBarry Smith }
4522d61bbb3SSatish Balay #undef __FUNC__
4532d61bbb3SSatish Balay #define __FUNC__ "MatMult_SeqBAIJ_7"
4542d61bbb3SSatish Balay int MatMult_SeqBAIJ_7(Mat A,Vec xx,Vec zz)
4552d61bbb3SSatish Balay {
4562d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
4573f1db9ecSBarry Smith   Scalar          *x,*z,*xb,sum1,sum2,sum3,sum4,sum5,sum6,sum7;
458e1311b90SBarry Smith   Scalar          x1,x2,x3,x4,x5,x6,x7;
4593f1db9ecSBarry Smith   MatScalar       *v;
460e1311b90SBarry Smith   int             ierr,mbs=a->mbs,i,*idx,*ii,j,n;
4612d61bbb3SSatish Balay 
462433994e6SBarry Smith   PetscFunctionBegin;
463e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
464e1311b90SBarry Smith   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
4652d61bbb3SSatish Balay 
4662d61bbb3SSatish Balay   idx   = a->j;
4672d61bbb3SSatish Balay   v     = a->a;
4682d61bbb3SSatish Balay   ii    = a->i;
4692d61bbb3SSatish Balay 
4702d61bbb3SSatish Balay   for ( i=0; i<mbs; i++ ) {
4712d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
4722d61bbb3SSatish Balay     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0; sum7 = 0.0;
4732d61bbb3SSatish Balay     for ( j=0; j<n; j++ ) {
4742d61bbb3SSatish Balay       xb = x + 7*(*idx++);
4752d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6];
4762d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[7]*x2  + v[14]*x3  + v[21]*x4 + v[28]*x5 + v[35]*x6 + v[42]*x7;
4772d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[8]*x2  + v[15]*x3  + v[22]*x4 + v[29]*x5 + v[36]*x6 + v[43]*x7;
4782d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[9]*x2  + v[16]*x3  + v[23]*x4 + v[30]*x5 + v[37]*x6 + v[44]*x7;
4792d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[10]*x2 + v[17]*x3  + v[24]*x4 + v[31]*x5 + v[38]*x6 + v[45]*x7;
4802d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[11]*x2 + v[18]*x3  + v[25]*x4 + v[32]*x5 + v[39]*x6 + v[46]*x7;
4812d61bbb3SSatish Balay       sum6 += v[5]*x1 + v[12]*x2 + v[19]*x3  + v[26]*x4 + v[33]*x5 + v[40]*x6 + v[47]*x7;
4822d61bbb3SSatish Balay       sum7 += v[6]*x1 + v[13]*x2 + v[20]*x3  + v[27]*x4 + v[34]*x5 + v[41]*x6 + v[48]*x7;
4832d61bbb3SSatish Balay       v += 49;
4842d61bbb3SSatish Balay     }
4852d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
4862d61bbb3SSatish Balay     z += 7;
4872d61bbb3SSatish Balay   }
4882d61bbb3SSatish Balay 
489e1311b90SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
490e1311b90SBarry Smith   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
4912d61bbb3SSatish Balay   PLogFlops(98*a->nz - a->m);
4922d61bbb3SSatish Balay   PetscFunctionReturn(0);
4932d61bbb3SSatish Balay }
4942d61bbb3SSatish Balay 
4953f1db9ecSBarry Smith /*
4963f1db9ecSBarry Smith     This will not work with MatScalar == float because it calls the BLAS
4973f1db9ecSBarry Smith */
4982d61bbb3SSatish Balay #undef __FUNC__
4992d61bbb3SSatish Balay #define __FUNC__ "MatMult_SeqBAIJ_N"
5002d61bbb3SSatish Balay int MatMult_SeqBAIJ_N(Mat A,Vec xx,Vec zz)
5012d61bbb3SSatish Balay {
5022d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
5033f1db9ecSBarry Smith   Scalar          *x,*z,*xb, *work,*workt;
5043f1db9ecSBarry Smith   MatScalar       *v;
505e1311b90SBarry Smith   int             ierr,mbs=a->mbs,i,*idx,*ii,bs=a->bs,j,n,bs2=a->bs2;
5063f1db9ecSBarry Smith   int             ncols,k;
5072d61bbb3SSatish Balay 
5082d61bbb3SSatish Balay   PetscFunctionBegin;
509e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
510e1311b90SBarry Smith   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
5112d61bbb3SSatish Balay 
5122d61bbb3SSatish Balay   idx   = a->j;
5132d61bbb3SSatish Balay   v     = a->a;
5142d61bbb3SSatish Balay   ii    = a->i;
515218c64b6SSatish Balay 
516218c64b6SSatish Balay 
5172d61bbb3SSatish Balay   if (!a->mult_work) {
5182d61bbb3SSatish Balay     k = PetscMax(a->m,a->n);
5192d61bbb3SSatish Balay     a->mult_work = (Scalar *) PetscMalloc((k+1)*sizeof(Scalar));CHKPTRQ(a->mult_work);
5202d61bbb3SSatish Balay   }
5212d61bbb3SSatish Balay   work = a->mult_work;
5222d61bbb3SSatish Balay   for ( i=0; i<mbs; i++ ) {
5232d61bbb3SSatish Balay     n     = ii[1] - ii[0]; ii++;
5242d61bbb3SSatish Balay     ncols = n*bs;
5252d61bbb3SSatish Balay     workt = work;
5262d61bbb3SSatish Balay     for ( j=0; j<n; j++ ) {
5272d61bbb3SSatish Balay       xb = x + bs*(*idx++);
5282d61bbb3SSatish Balay       for ( k=0; k<bs; k++ ) workt[k] = xb[k];
5292d61bbb3SSatish Balay       workt += bs;
5302d61bbb3SSatish Balay     }
5313f1db9ecSBarry Smith     Kernel_w_gets_Ar_times_v(bs,ncols,work,v,z);
5323f1db9ecSBarry Smith     /* LAgemv_("N",&bs,&ncols,&_DOne,v,&bs,work,&_One,&_DZero,z,&_One); */
5332d61bbb3SSatish Balay     v += n*bs2;
5342d61bbb3SSatish Balay     z += bs;
5352d61bbb3SSatish Balay   }
536e1311b90SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
537e1311b90SBarry Smith   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
5382d61bbb3SSatish Balay   PLogFlops(2*a->nz*bs2 - a->m);
5392d61bbb3SSatish Balay   PetscFunctionReturn(0);
5402d61bbb3SSatish Balay }
5412d61bbb3SSatish Balay 
5422d61bbb3SSatish Balay #undef __FUNC__
5432d61bbb3SSatish Balay #define __FUNC__ "MatMultAdd_SeqBAIJ_1"
5442d61bbb3SSatish Balay int MatMultAdd_SeqBAIJ_1(Mat A,Vec xx,Vec yy,Vec zz)
5452d61bbb3SSatish Balay {
5462d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
5473f1db9ecSBarry Smith   Scalar          *x,*y,*z,sum;
5483f1db9ecSBarry Smith   MatScalar       *v;
549e1311b90SBarry Smith   int             ierr,mbs=a->mbs,i,*idx,*ii,n;
5502d61bbb3SSatish Balay 
5512d61bbb3SSatish Balay   PetscFunctionBegin;
552e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
553e1311b90SBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
5542e8a6d31SBarry Smith   if (zz != yy) {
555e1311b90SBarry Smith     ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
5562e8a6d31SBarry Smith   } else {
5572e8a6d31SBarry Smith     z = y;
5582e8a6d31SBarry Smith   }
5592d61bbb3SSatish Balay 
5602d61bbb3SSatish Balay   idx   = a->j;
5612d61bbb3SSatish Balay   v     = a->a;
5622d61bbb3SSatish Balay   ii    = a->i;
5632d61bbb3SSatish Balay 
5642d61bbb3SSatish Balay   for ( i=0; i<mbs; i++ ) {
5652d61bbb3SSatish Balay     n    = ii[1] - ii[0]; ii++;
5662d61bbb3SSatish Balay     sum  = y[i];
5672d61bbb3SSatish Balay     while (n--) sum += *v++ * x[*idx++];
5682d61bbb3SSatish Balay     z[i] = sum;
5692d61bbb3SSatish Balay   }
570e1311b90SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
571e1311b90SBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
5722e8a6d31SBarry Smith   if (zz != yy) {
573e1311b90SBarry Smith     ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
5742e8a6d31SBarry Smith   }
5752d61bbb3SSatish Balay   PLogFlops(2*a->nz);
5762d61bbb3SSatish Balay   PetscFunctionReturn(0);
5772d61bbb3SSatish Balay }
5782d61bbb3SSatish Balay 
5792d61bbb3SSatish Balay #undef __FUNC__
5802d61bbb3SSatish Balay #define __FUNC__ "MatMultAdd_SeqBAIJ_2"
5812d61bbb3SSatish Balay int MatMultAdd_SeqBAIJ_2(Mat A,Vec xx,Vec yy,Vec zz)
5822d61bbb3SSatish Balay {
5832d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
5843f1db9ecSBarry Smith   Scalar          *x,*y,*z,*xb,sum1,sum2;
585e1311b90SBarry Smith   Scalar          x1,x2;
5863f1db9ecSBarry Smith   MatScalar       *v;
587e1311b90SBarry Smith   int             ierr,mbs=a->mbs,i,*idx,*ii,j,n;
5882d61bbb3SSatish Balay 
5892d61bbb3SSatish Balay   PetscFunctionBegin;
590e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
591e1311b90SBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
5922e8a6d31SBarry Smith   if (zz != yy) {
593e1311b90SBarry Smith     ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
5942e8a6d31SBarry Smith   } else {
5952e8a6d31SBarry Smith     z = y;
5962e8a6d31SBarry Smith   }
5972d61bbb3SSatish Balay 
5982d61bbb3SSatish Balay   idx   = a->j;
5992d61bbb3SSatish Balay   v     = a->a;
6002d61bbb3SSatish Balay   ii    = a->i;
6012d61bbb3SSatish Balay 
6022d61bbb3SSatish Balay   for ( i=0; i<mbs; i++ ) {
6032d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
6042d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1];
6052d61bbb3SSatish Balay     for ( j=0; j<n; j++ ) {
6062d61bbb3SSatish Balay       xb = x + 2*(*idx++); x1 = xb[0]; x2 = xb[1];
6072d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[2]*x2;
6082d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[3]*x2;
6092d61bbb3SSatish Balay       v += 4;
6102d61bbb3SSatish Balay     }
6112d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2;
6122d61bbb3SSatish Balay     z += 2; y += 2;
6132d61bbb3SSatish Balay   }
614e1311b90SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
615e1311b90SBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
6162e8a6d31SBarry Smith   if (zz != yy) {
617e1311b90SBarry Smith     ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
6182e8a6d31SBarry Smith   }
6192d61bbb3SSatish Balay   PLogFlops(4*a->nz);
6202d61bbb3SSatish Balay   PetscFunctionReturn(0);
6212d61bbb3SSatish Balay }
6222d61bbb3SSatish Balay 
6232d61bbb3SSatish Balay #undef __FUNC__
6242d61bbb3SSatish Balay #define __FUNC__ "MatMultAdd_SeqBAIJ_3"
6252d61bbb3SSatish Balay int MatMultAdd_SeqBAIJ_3(Mat A,Vec xx,Vec yy,Vec zz)
6262d61bbb3SSatish Balay {
6272d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
6283f1db9ecSBarry Smith   Scalar          *x,*y,*z,*xb,sum1,sum2,sum3,x1,x2,x3;
6293f1db9ecSBarry Smith   MatScalar       *v;
630e1311b90SBarry Smith   int             ierr,mbs=a->mbs,i,*idx,*ii,j,n;
6312d61bbb3SSatish Balay 
6322d61bbb3SSatish Balay   PetscFunctionBegin;
633e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
634e1311b90SBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
6352e8a6d31SBarry Smith   if (zz != yy) {
636e1311b90SBarry Smith     ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
6372e8a6d31SBarry Smith   } else {
6382e8a6d31SBarry Smith     z = y;
6392e8a6d31SBarry Smith   }
6402d61bbb3SSatish Balay 
6412d61bbb3SSatish Balay   idx   = a->j;
6422d61bbb3SSatish Balay   v     = a->a;
6432d61bbb3SSatish Balay   ii    = a->i;
6442d61bbb3SSatish Balay 
6452d61bbb3SSatish Balay   for ( i=0; i<mbs; i++ ) {
6462d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
6472d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2];
6482d61bbb3SSatish Balay     for ( j=0; j<n; j++ ) {
6492d61bbb3SSatish Balay       xb = x + 3*(*idx++); x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
6502d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[3]*x2 + v[6]*x3;
6512d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[4]*x2 + v[7]*x3;
6522d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[5]*x2 + v[8]*x3;
6532d61bbb3SSatish Balay       v += 9;
6542d61bbb3SSatish Balay     }
6552d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3;
6562d61bbb3SSatish Balay     z += 3; y += 3;
6572d61bbb3SSatish Balay   }
658e1311b90SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
659e1311b90SBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
6602e8a6d31SBarry Smith   if (zz != yy) {
661e1311b90SBarry Smith     ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
6622e8a6d31SBarry Smith   }
6632d61bbb3SSatish Balay   PLogFlops(18*a->nz);
6642d61bbb3SSatish Balay   PetscFunctionReturn(0);
6652d61bbb3SSatish Balay }
6662d61bbb3SSatish Balay 
6672d61bbb3SSatish Balay #undef __FUNC__
6682d61bbb3SSatish Balay #define __FUNC__ "MatMultAdd_SeqBAIJ_4"
6692d61bbb3SSatish Balay int MatMultAdd_SeqBAIJ_4(Mat A,Vec xx,Vec yy,Vec zz)
6702d61bbb3SSatish Balay {
6712d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
6723f1db9ecSBarry Smith   Scalar          *x,*y,*z,*xb,sum1,sum2,sum3,sum4,x1,x2,x3,x4;
6733f1db9ecSBarry Smith   MatScalar       *v;
674e1311b90SBarry Smith   int             ierr,mbs=a->mbs,i,*idx,*ii;
6752d61bbb3SSatish Balay   int             j,n;
6762d61bbb3SSatish Balay 
6772d61bbb3SSatish Balay   PetscFunctionBegin;
678e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
679e1311b90SBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
6802e8a6d31SBarry Smith   if (zz != yy) {
681e1311b90SBarry Smith     ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
6822e8a6d31SBarry Smith   } else {
6832e8a6d31SBarry Smith     z = y;
6842e8a6d31SBarry Smith   }
6852d61bbb3SSatish Balay 
6862d61bbb3SSatish Balay   idx   = a->j;
6872d61bbb3SSatish Balay   v     = a->a;
6882d61bbb3SSatish Balay   ii    = a->i;
6892d61bbb3SSatish Balay 
6902d61bbb3SSatish Balay   for ( i=0; i<mbs; i++ ) {
6912d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
6922d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3];
6932d61bbb3SSatish Balay     for ( j=0; j<n; j++ ) {
6942d61bbb3SSatish Balay       xb = x + 4*(*idx++);
6952d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
6962d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[4]*x2 + v[8]*x3   + v[12]*x4;
6972d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[5]*x2 + v[9]*x3   + v[13]*x4;
6982d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[6]*x2 + v[10]*x3  + v[14]*x4;
6992d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[7]*x2 + v[11]*x3  + v[15]*x4;
7002d61bbb3SSatish Balay       v += 16;
7012d61bbb3SSatish Balay     }
7022d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4;
7032d61bbb3SSatish Balay     z += 4; y += 4;
7042d61bbb3SSatish Balay   }
705e1311b90SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
706e1311b90SBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
7072e8a6d31SBarry Smith   if (zz != yy) {
708e1311b90SBarry Smith     ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
7092e8a6d31SBarry Smith   }
7102d61bbb3SSatish Balay   PLogFlops(32*a->nz);
7112d61bbb3SSatish Balay   PetscFunctionReturn(0);
7122d61bbb3SSatish Balay }
7132d61bbb3SSatish Balay 
7142d61bbb3SSatish Balay #undef __FUNC__
7152d61bbb3SSatish Balay #define __FUNC__ "MatMultAdd_SeqBAIJ_5"
7162d61bbb3SSatish Balay int MatMultAdd_SeqBAIJ_5(Mat A,Vec xx,Vec yy,Vec zz)
7172d61bbb3SSatish Balay {
7182d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
7193f1db9ecSBarry Smith   Scalar          *x,*y,*z,*xb,sum1,sum2,sum3,sum4,sum5,x1,x2,x3,x4,x5;
7203f1db9ecSBarry Smith   MatScalar       *v;
721e1311b90SBarry Smith   int             ierr,mbs=a->mbs,i,*idx,*ii,j,n;
7222d61bbb3SSatish Balay 
7232d61bbb3SSatish Balay   PetscFunctionBegin;
724e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
725e1311b90SBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
7262e8a6d31SBarry Smith   if (zz != yy) {
727e1311b90SBarry Smith     ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
7282e8a6d31SBarry Smith   } else {
7292e8a6d31SBarry Smith     z = y;
7302e8a6d31SBarry Smith   }
7312d61bbb3SSatish Balay 
7322d61bbb3SSatish Balay   idx   = a->j;
7332d61bbb3SSatish Balay   v     = a->a;
7342d61bbb3SSatish Balay   ii    = a->i;
7352d61bbb3SSatish Balay 
7362d61bbb3SSatish Balay   for ( i=0; i<mbs; i++ ) {
7372d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
7382d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4];
7392d61bbb3SSatish Balay     for ( j=0; j<n; j++ ) {
7402d61bbb3SSatish Balay       xb = x + 5*(*idx++);
7412d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4];
7422d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[5]*x2 + v[10]*x3  + v[15]*x4 + v[20]*x5;
7432d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[6]*x2 + v[11]*x3  + v[16]*x4 + v[21]*x5;
7442d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[7]*x2 + v[12]*x3  + v[17]*x4 + v[22]*x5;
7452d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[8]*x2 + v[13]*x3  + v[18]*x4 + v[23]*x5;
7462d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[9]*x2 + v[14]*x3  + v[19]*x4 + v[24]*x5;
7472d61bbb3SSatish Balay       v += 25;
7482d61bbb3SSatish Balay     }
7492d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5;
7502d61bbb3SSatish Balay     z += 5; y += 5;
7512d61bbb3SSatish Balay   }
752e1311b90SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
753e1311b90SBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
7542e8a6d31SBarry Smith   if (zz != yy) {
755e1311b90SBarry Smith     ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
7562e8a6d31SBarry Smith   }
7572d61bbb3SSatish Balay   PLogFlops(50*a->nz);
7582d61bbb3SSatish Balay   PetscFunctionReturn(0);
7592d61bbb3SSatish Balay }
76015091d37SBarry Smith #undef __FUNC__
76115091d37SBarry Smith #define __FUNC__ "MatMultAdd_SeqBAIJ_6"
76215091d37SBarry Smith int MatMultAdd_SeqBAIJ_6(Mat A,Vec xx,Vec yy,Vec zz)
76315091d37SBarry Smith {
76415091d37SBarry Smith   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
76515091d37SBarry Smith   Scalar          *x,*y,*z,*xb,sum1,sum2,sum3,sum4,sum5,sum6;
76615091d37SBarry Smith   Scalar          x1,x2,x3,x4,x5,x6;
76715091d37SBarry Smith   MatScalar       *v;
76815091d37SBarry Smith   int             ierr,mbs=a->mbs,i,*idx,*ii,j,n;
76915091d37SBarry Smith 
77015091d37SBarry Smith   PetscFunctionBegin;
77115091d37SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
77215091d37SBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
77315091d37SBarry Smith   if (zz != yy) {
77415091d37SBarry Smith     ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
77515091d37SBarry Smith   } else {
77615091d37SBarry Smith     z = y;
77715091d37SBarry Smith   }
77815091d37SBarry Smith 
77915091d37SBarry Smith   idx   = a->j;
78015091d37SBarry Smith   v     = a->a;
78115091d37SBarry Smith   ii    = a->i;
78215091d37SBarry Smith 
78315091d37SBarry Smith   for ( i=0; i<mbs; i++ ) {
78415091d37SBarry Smith     n  = ii[1] - ii[0]; ii++;
78515091d37SBarry Smith     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4]; sum6 = y[5];
78615091d37SBarry Smith     for ( j=0; j<n; j++ ) {
7873b95cb0eSSatish Balay       xb = x + 6*(*idx++);
78815091d37SBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5];
78915091d37SBarry Smith       sum1 += v[0]*x1 + v[6]*x2  + v[12]*x3  + v[18]*x4 + v[24]*x5 + v[30]*x6;
79015091d37SBarry Smith       sum2 += v[1]*x1 + v[7]*x2  + v[13]*x3  + v[19]*x4 + v[25]*x5 + v[31]*x6;
79115091d37SBarry Smith       sum3 += v[2]*x1 + v[8]*x2  + v[14]*x3  + v[20]*x4 + v[26]*x5 + v[32]*x6;
79215091d37SBarry Smith       sum4 += v[3]*x1 + v[9]*x2  + v[15]*x3  + v[21]*x4 + v[27]*x5 + v[33]*x6;
79315091d37SBarry Smith       sum5 += v[4]*x1 + v[10]*x2 + v[16]*x3  + v[22]*x4 + v[28]*x5 + v[34]*x6;
79415091d37SBarry Smith       sum6 += v[5]*x1 + v[11]*x2 + v[17]*x3  + v[23]*x4 + v[29]*x5 + v[35]*x6;
79515091d37SBarry Smith       v += 36;
79615091d37SBarry Smith     }
79715091d37SBarry Smith     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6;
79815091d37SBarry Smith     z += 6; y += 6;
79915091d37SBarry Smith   }
80015091d37SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
80115091d37SBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
80215091d37SBarry Smith   if (zz != yy) {
80315091d37SBarry Smith     ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
80415091d37SBarry Smith   }
80515091d37SBarry Smith   PLogFlops(72*a->nz);
80615091d37SBarry Smith   PetscFunctionReturn(0);
80715091d37SBarry Smith }
8082d61bbb3SSatish Balay 
8092d61bbb3SSatish Balay #undef __FUNC__
8102d61bbb3SSatish Balay #define __FUNC__ "MatMultAdd_SeqBAIJ_7"
8112d61bbb3SSatish Balay int MatMultAdd_SeqBAIJ_7(Mat A,Vec xx,Vec yy,Vec zz)
8122d61bbb3SSatish Balay {
8132d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
8143f1db9ecSBarry Smith   Scalar          *x,*y,*z,*xb,sum1,sum2,sum3,sum4,sum5,sum6,sum7;
815e1311b90SBarry Smith   Scalar          x1,x2,x3,x4,x5,x6,x7;
8163f1db9ecSBarry Smith   MatScalar       *v;
817e1311b90SBarry Smith   int             ierr,mbs=a->mbs,i,*idx,*ii,j,n;
8182d61bbb3SSatish Balay 
8192d61bbb3SSatish Balay   PetscFunctionBegin;
820e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
821e1311b90SBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
8222e8a6d31SBarry Smith   if (zz != yy) {
823e1311b90SBarry Smith     ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
8242e8a6d31SBarry Smith   } else {
8252e8a6d31SBarry Smith     z = y;
8262e8a6d31SBarry Smith   }
8272d61bbb3SSatish Balay 
8282d61bbb3SSatish Balay   idx   = a->j;
8292d61bbb3SSatish Balay   v     = a->a;
8302d61bbb3SSatish Balay   ii    = a->i;
8312d61bbb3SSatish Balay 
8322d61bbb3SSatish Balay   for ( i=0; i<mbs; i++ ) {
8332d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
8342d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4]; sum6 = y[5]; sum7 = y[6];
8352d61bbb3SSatish Balay     for ( j=0; j<n; j++ ) {
8362d61bbb3SSatish Balay       xb = x + 7*(*idx++);
8372d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6];
8382d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[7]*x2  + v[14]*x3  + v[21]*x4 + v[28]*x5 + v[35]*x6 + v[42]*x7;
8392d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[8]*x2  + v[15]*x3  + v[22]*x4 + v[29]*x5 + v[36]*x6 + v[43]*x7;
8402d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[9]*x2  + v[16]*x3  + v[23]*x4 + v[30]*x5 + v[37]*x6 + v[44]*x7;
8412d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[10]*x2 + v[17]*x3  + v[24]*x4 + v[31]*x5 + v[38]*x6 + v[45]*x7;
8422d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[11]*x2 + v[18]*x3  + v[25]*x4 + v[32]*x5 + v[39]*x6 + v[46]*x7;
8432d61bbb3SSatish Balay       sum6 += v[5]*x1 + v[12]*x2 + v[19]*x3  + v[26]*x4 + v[33]*x5 + v[40]*x6 + v[47]*x7;
8442d61bbb3SSatish Balay       sum7 += v[6]*x1 + v[13]*x2 + v[20]*x3  + v[27]*x4 + v[34]*x5 + v[41]*x6 + v[48]*x7;
8452d61bbb3SSatish Balay       v += 49;
8462d61bbb3SSatish Balay     }
8472d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
8482d61bbb3SSatish Balay     z += 7; y += 7;
8492d61bbb3SSatish Balay   }
850e1311b90SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
851e1311b90SBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
8522e8a6d31SBarry Smith   if (zz != yy) {
853e1311b90SBarry Smith     ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
8542e8a6d31SBarry Smith   }
8552d61bbb3SSatish Balay   PLogFlops(98*a->nz);
8562d61bbb3SSatish Balay   PetscFunctionReturn(0);
8572d61bbb3SSatish Balay }
858218c64b6SSatish Balay 
8592d61bbb3SSatish Balay #undef __FUNC__
8602d61bbb3SSatish Balay #define __FUNC__ "MatMultAdd_SeqBAIJ_N"
8612d61bbb3SSatish Balay int MatMultAdd_SeqBAIJ_N(Mat A,Vec xx,Vec yy,Vec zz)
8622d61bbb3SSatish Balay {
8632d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ *) A->data;
8643f1db9ecSBarry Smith   Scalar         *x,*z,*xb,*work,*workt;
8653f1db9ecSBarry Smith   MatScalar      *v;
8662d61bbb3SSatish Balay   int            mbs=a->mbs,i,*idx,*ii,bs=a->bs,j,n,bs2=a->bs2,ierr;
8673f1db9ecSBarry Smith   int            ncols,k;
868218c64b6SSatish Balay 
8692d61bbb3SSatish Balay   PetscFunctionBegin;
8702d61bbb3SSatish Balay   if ( xx != yy) { ierr = VecCopy(yy,zz);CHKERRQ(ierr); }
8712d61bbb3SSatish Balay 
872e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
873e1311b90SBarry Smith   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
8742d61bbb3SSatish Balay 
8752d61bbb3SSatish Balay   idx   = a->j;
8762d61bbb3SSatish Balay   v     = a->a;
8772d61bbb3SSatish Balay   ii    = a->i;
8782d61bbb3SSatish Balay 
8792d61bbb3SSatish Balay 
8802d61bbb3SSatish Balay   if (!a->mult_work) {
8812d61bbb3SSatish Balay     k = PetscMax(a->m,a->n);
882544c85edSSatish Balay     a->mult_work = (Scalar *) PetscMalloc((k+1)*sizeof(Scalar));CHKPTRQ(a->mult_work);
8832d61bbb3SSatish Balay   }
8842d61bbb3SSatish Balay   work = a->mult_work;
8852d61bbb3SSatish Balay   for ( i=0; i<mbs; i++ ) {
8862d61bbb3SSatish Balay     n     = ii[1] - ii[0]; ii++;
8872d61bbb3SSatish Balay     ncols = n*bs;
8882d61bbb3SSatish Balay     workt = work;
8892d61bbb3SSatish Balay     for ( j=0; j<n; j++ ) {
8902d61bbb3SSatish Balay       xb = x + bs*(*idx++);
8912d61bbb3SSatish Balay       for ( k=0; k<bs; k++ ) workt[k] = xb[k];
8922d61bbb3SSatish Balay       workt += bs;
8932d61bbb3SSatish Balay     }
8943f1db9ecSBarry Smith     Kernel_w_gets_w_plus_Ar_times_v(bs,ncols,work,v,z);
8953f1db9ecSBarry Smith     /* LAgemv_("N",&bs,&ncols,&_DOne,v,&bs,work,&_One,&_DOne,z,&_One); */
8962d61bbb3SSatish Balay     v += n*bs2;
8972d61bbb3SSatish Balay     z += bs;
8982d61bbb3SSatish Balay   }
899e1311b90SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
900e1311b90SBarry Smith   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
9012d61bbb3SSatish Balay   PLogFlops(2*a->nz*bs2 );
9022d61bbb3SSatish Balay   PetscFunctionReturn(0);
9032d61bbb3SSatish Balay }
9042d61bbb3SSatish Balay 
9052d61bbb3SSatish Balay #undef __FUNC__
9062d61bbb3SSatish Balay #define __FUNC__ "MatMultTrans_SeqBAIJ"
9072d61bbb3SSatish Balay int MatMultTrans_SeqBAIJ(Mat A,Vec xx,Vec zz)
9082d61bbb3SSatish Balay {
9092d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
9102d61bbb3SSatish Balay   Scalar          *xg,*zg,*zb;
9113f1db9ecSBarry Smith   Scalar          *x,*z,*xb,x1,x2,x3,x4,x5;
9123f1db9ecSBarry Smith   MatScalar       *v;
9132d61bbb3SSatish Balay   int             mbs=a->mbs,i,*idx,*ii,*ai=a->i,rval,N=a->n;
9142d61bbb3SSatish Balay   int             bs=a->bs,j,n,bs2=a->bs2,*ib,ierr;
9152d61bbb3SSatish Balay 
9162d61bbb3SSatish Balay 
9172d61bbb3SSatish Balay   PetscFunctionBegin;
918e1311b90SBarry Smith   ierr = VecGetArray(xx,&xg);CHKERRQ(ierr); x = xg;
919e1311b90SBarry Smith   ierr = VecGetArray(zz,&zg);CHKERRQ(ierr); z = zg;
920549d3d68SSatish Balay   ierr = PetscMemzero(z,N*sizeof(Scalar));CHKERRQ(ierr);
9212d61bbb3SSatish Balay 
9222d61bbb3SSatish Balay   idx   = a->j;
9232d61bbb3SSatish Balay   v     = a->a;
9242d61bbb3SSatish Balay   ii    = a->i;
9252d61bbb3SSatish Balay 
9262d61bbb3SSatish Balay   switch (bs) {
9272d61bbb3SSatish Balay   case 1:
9282d61bbb3SSatish Balay     for ( i=0; i<mbs; i++ ) {
9292d61bbb3SSatish Balay       n  = ii[1] - ii[0]; ii++;
9302d61bbb3SSatish Balay       xb = x + i; x1 = xb[0];
9312d61bbb3SSatish Balay       ib = idx + ai[i];
9322d61bbb3SSatish Balay       for ( j=0; j<n; j++ ) {
9332d61bbb3SSatish Balay         rval    = ib[j];
9342d61bbb3SSatish Balay         z[rval] += *v++ * x1;
9352d61bbb3SSatish Balay       }
9362d61bbb3SSatish Balay     }
9372d61bbb3SSatish Balay     break;
9382d61bbb3SSatish Balay   case 2:
9392d61bbb3SSatish Balay     for ( i=0; i<mbs; i++ ) {
9402d61bbb3SSatish Balay       n  = ii[1] - ii[0]; ii++;
9412d61bbb3SSatish Balay       xb = x + 2*i; x1 = xb[0]; x2 = xb[1];
9422d61bbb3SSatish Balay       ib = idx + ai[i];
9432d61bbb3SSatish Balay       for ( j=0; j<n; j++ ) {
9442d61bbb3SSatish Balay         rval      = ib[j]*2;
9452d61bbb3SSatish Balay         z[rval++] += v[0]*x1 + v[1]*x2;
9462d61bbb3SSatish Balay         z[rval++] += v[2]*x1 + v[3]*x2;
9472d61bbb3SSatish Balay         v += 4;
9482d61bbb3SSatish Balay       }
9492d61bbb3SSatish Balay     }
9502d61bbb3SSatish Balay     break;
9512d61bbb3SSatish Balay   case 3:
9522d61bbb3SSatish Balay     for ( i=0; i<mbs; i++ ) {
9532d61bbb3SSatish Balay       n  = ii[1] - ii[0]; ii++;
9542d61bbb3SSatish Balay       xb = x + 3*i; x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
9552d61bbb3SSatish Balay       ib = idx + ai[i];
9562d61bbb3SSatish Balay       for ( j=0; j<n; j++ ) {
9572d61bbb3SSatish Balay         rval      = ib[j]*3;
9582d61bbb3SSatish Balay         z[rval++] += v[0]*x1 + v[1]*x2 + v[2]*x3;
9592d61bbb3SSatish Balay         z[rval++] += v[3]*x1 + v[4]*x2 + v[5]*x3;
9602d61bbb3SSatish Balay         z[rval++] += v[6]*x1 + v[7]*x2 + v[8]*x3;
9612d61bbb3SSatish Balay         v += 9;
9622d61bbb3SSatish Balay       }
9632d61bbb3SSatish Balay     }
9642d61bbb3SSatish Balay     break;
9652d61bbb3SSatish Balay   case 4:
9662d61bbb3SSatish Balay     for ( i=0; i<mbs; i++ ) {
9672d61bbb3SSatish Balay       n  = ii[1] - ii[0]; ii++;
9682d61bbb3SSatish Balay       xb = x + 4*i; x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
9692d61bbb3SSatish Balay       ib = idx + ai[i];
9702d61bbb3SSatish Balay       for ( j=0; j<n; j++ ) {
9712d61bbb3SSatish Balay         rval      = ib[j]*4;
9722d61bbb3SSatish Balay         z[rval++] +=  v[0]*x1 +  v[1]*x2 +  v[2]*x3 +  v[3]*x4;
9732d61bbb3SSatish Balay         z[rval++] +=  v[4]*x1 +  v[5]*x2 +  v[6]*x3 +  v[7]*x4;
9742d61bbb3SSatish Balay         z[rval++] +=  v[8]*x1 +  v[9]*x2 + v[10]*x3 + v[11]*x4;
9752d61bbb3SSatish Balay         z[rval++] += v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4;
9762d61bbb3SSatish Balay         v += 16;
9772d61bbb3SSatish Balay       }
9782d61bbb3SSatish Balay     }
9792d61bbb3SSatish Balay     break;
9802d61bbb3SSatish Balay   case 5:
9812d61bbb3SSatish Balay     for ( i=0; i<mbs; i++ ) {
9822d61bbb3SSatish Balay       n  = ii[1] - ii[0]; ii++;
9832d61bbb3SSatish Balay       xb = x + 5*i; x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
9842d61bbb3SSatish Balay       x4 = xb[3];   x5 = xb[4];
9852d61bbb3SSatish Balay       ib = idx + ai[i];
9862d61bbb3SSatish Balay       for ( j=0; j<n; j++ ) {
9872d61bbb3SSatish Balay         rval      = ib[j]*5;
9882d61bbb3SSatish Balay         z[rval++] +=  v[0]*x1 +  v[1]*x2 +  v[2]*x3 +  v[3]*x4 +  v[4]*x5;
9892d61bbb3SSatish Balay         z[rval++] +=  v[5]*x1 +  v[6]*x2 +  v[7]*x3 +  v[8]*x4 +  v[9]*x5;
9902d61bbb3SSatish Balay         z[rval++] += v[10]*x1 + v[11]*x2 + v[12]*x3 + v[13]*x4 + v[14]*x5;
9912d61bbb3SSatish Balay         z[rval++] += v[15]*x1 + v[16]*x2 + v[17]*x3 + v[18]*x4 + v[19]*x5;
9922d61bbb3SSatish Balay         z[rval++] += v[20]*x1 + v[21]*x2 + v[22]*x3 + v[23]*x4 + v[24]*x5;
9932d61bbb3SSatish Balay         v += 25;
9942d61bbb3SSatish Balay       }
9952d61bbb3SSatish Balay     }
9962d61bbb3SSatish Balay     break;
9972d61bbb3SSatish Balay       /* block sizes larger then 5 by 5 are handled by BLAS */
9982d61bbb3SSatish Balay     default: {
9993f1db9ecSBarry Smith       int       ncols,k;
10003f1db9ecSBarry Smith       Scalar    *work,*workt;
10013f1db9ecSBarry Smith 
10022d61bbb3SSatish Balay       if (!a->mult_work) {
10032d61bbb3SSatish Balay         k = PetscMax(a->m,a->n);
1004544c85edSSatish Balay         a->mult_work = (Scalar *) PetscMalloc((k+1)*sizeof(Scalar));CHKPTRQ(a->mult_work);
10052d61bbb3SSatish Balay       }
10062d61bbb3SSatish Balay       work = a->mult_work;
10072d61bbb3SSatish Balay       for ( i=0; i<mbs; i++ ) {
10082d61bbb3SSatish Balay         n     = ii[1] - ii[0]; ii++;
10092d61bbb3SSatish Balay         ncols = n*bs;
1010549d3d68SSatish Balay         ierr  = PetscMemzero(work,ncols*sizeof(Scalar));CHKERRQ(ierr);
1011d824769bSBarry Smith         Kernel_w_gets_w_plus_trans_Ar_times_v(bs,ncols,x,v,work);
10123f1db9ecSBarry Smith         /* LAgemv_("T",&bs,&ncols,&_DOne,v,&bs,x,&_One,&_DOne,work,&_One); */
10132d61bbb3SSatish Balay         v += n*bs2;
10142d61bbb3SSatish Balay         x += bs;
10152d61bbb3SSatish Balay         workt = work;
10162d61bbb3SSatish Balay         for ( j=0; j<n; j++ ) {
10172d61bbb3SSatish Balay           zb = z + bs*(*idx++);
10182d61bbb3SSatish Balay           for ( k=0; k<bs; k++ ) zb[k] += workt[k] ;
10192d61bbb3SSatish Balay           workt += bs;
10202d61bbb3SSatish Balay         }
10212d61bbb3SSatish Balay       }
10222d61bbb3SSatish Balay     }
10232d61bbb3SSatish Balay   }
10242d61bbb3SSatish Balay   ierr = VecRestoreArray(xx,&xg);CHKERRQ(ierr);
10252d61bbb3SSatish Balay   ierr = VecRestoreArray(zz,&zg);CHKERRQ(ierr);
10262d61bbb3SSatish Balay   PLogFlops(2*a->nz*a->bs2 - a->n);
10272d61bbb3SSatish Balay   PetscFunctionReturn(0);
10282d61bbb3SSatish Balay }
10292d61bbb3SSatish Balay 
10302d61bbb3SSatish Balay #undef __FUNC__
10312d61bbb3SSatish Balay #define __FUNC__ "MatMultTransAdd_SeqBAIJ"
10322d61bbb3SSatish Balay int MatMultTransAdd_SeqBAIJ(Mat A,Vec xx,Vec yy,Vec zz)
10332d61bbb3SSatish Balay 
10342d61bbb3SSatish Balay {
10352d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
10362d61bbb3SSatish Balay   Scalar          *xg,*zg,*zb;
10373f1db9ecSBarry Smith   Scalar          *x,*z,*xb,x1,x2,x3,x4,x5;
10383f1db9ecSBarry Smith   MatScalar       *v;
10392d61bbb3SSatish Balay   int             mbs=a->mbs,i,*idx,*ii,*ai=a->i,rval,N=a->n;
10402d61bbb3SSatish Balay   int             bs=a->bs,j,n,bs2=a->bs2,*ib,ierr;
10412d61bbb3SSatish Balay 
10422d61bbb3SSatish Balay   PetscFunctionBegin;
1043e1311b90SBarry Smith   ierr = VecGetArray(xx,&xg);CHKERRQ(ierr); x = xg;
1044e1311b90SBarry Smith   ierr = VecGetArray(zz,&zg);CHKERRQ(ierr); z = zg;
10452d61bbb3SSatish Balay 
10462d61bbb3SSatish Balay   if ( yy != zz ) { ierr = VecCopy(yy,zz);CHKERRQ(ierr); }
1047549d3d68SSatish Balay   else {
1048549d3d68SSatish Balay     ierr = PetscMemzero(z,N*sizeof(Scalar));CHKERRQ(ierr);
1049549d3d68SSatish Balay   }
10502d61bbb3SSatish Balay 
10512d61bbb3SSatish Balay   idx   = a->j;
10522d61bbb3SSatish Balay   v     = a->a;
10532d61bbb3SSatish Balay   ii    = a->i;
10542d61bbb3SSatish Balay 
10552d61bbb3SSatish Balay   switch (bs) {
10562d61bbb3SSatish Balay   case 1:
10572d61bbb3SSatish Balay     for ( i=0; i<mbs; i++ ) {
10582d61bbb3SSatish Balay       n  = ii[1] - ii[0]; ii++;
10592d61bbb3SSatish Balay       xb = x + i; x1 = xb[0];
10602d61bbb3SSatish Balay       ib = idx + ai[i];
10612d61bbb3SSatish Balay       for ( j=0; j<n; j++ ) {
10622d61bbb3SSatish Balay         rval    = ib[j];
10632d61bbb3SSatish Balay         z[rval] += *v++ * x1;
10642d61bbb3SSatish Balay       }
10652d61bbb3SSatish Balay     }
10662d61bbb3SSatish Balay     break;
10672d61bbb3SSatish Balay   case 2:
10682d61bbb3SSatish Balay     for ( i=0; i<mbs; i++ ) {
10692d61bbb3SSatish Balay       n  = ii[1] - ii[0]; ii++;
10702d61bbb3SSatish Balay       xb = x + 2*i; x1 = xb[0]; x2 = xb[1];
10712d61bbb3SSatish Balay       ib = idx + ai[i];
10722d61bbb3SSatish Balay       for ( j=0; j<n; j++ ) {
10732d61bbb3SSatish Balay         rval      = ib[j]*2;
10742d61bbb3SSatish Balay         z[rval++] += v[0]*x1 + v[1]*x2;
10752d61bbb3SSatish Balay         z[rval++] += v[2]*x1 + v[3]*x2;
10762d61bbb3SSatish Balay         v += 4;
10772d61bbb3SSatish Balay       }
10782d61bbb3SSatish Balay     }
10792d61bbb3SSatish Balay     break;
10802d61bbb3SSatish Balay   case 3:
10812d61bbb3SSatish Balay     for ( i=0; i<mbs; i++ ) {
10822d61bbb3SSatish Balay       n  = ii[1] - ii[0]; ii++;
10832d61bbb3SSatish Balay       xb = x + 3*i; x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
10842d61bbb3SSatish Balay       ib = idx + ai[i];
10852d61bbb3SSatish Balay       for ( j=0; j<n; j++ ) {
10862d61bbb3SSatish Balay         rval      = ib[j]*3;
10872d61bbb3SSatish Balay         z[rval++] += v[0]*x1 + v[1]*x2 + v[2]*x3;
10882d61bbb3SSatish Balay         z[rval++] += v[3]*x1 + v[4]*x2 + v[5]*x3;
10892d61bbb3SSatish Balay         z[rval++] += v[6]*x1 + v[7]*x2 + v[8]*x3;
10902d61bbb3SSatish Balay         v += 9;
10912d61bbb3SSatish Balay       }
10922d61bbb3SSatish Balay     }
10932d61bbb3SSatish Balay     break;
10942d61bbb3SSatish Balay   case 4:
10952d61bbb3SSatish Balay     for ( i=0; i<mbs; i++ ) {
10962d61bbb3SSatish Balay       n  = ii[1] - ii[0]; ii++;
10972d61bbb3SSatish Balay       xb = x + 4*i; x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
10982d61bbb3SSatish Balay       ib = idx + ai[i];
10992d61bbb3SSatish Balay       for ( j=0; j<n; j++ ) {
11002d61bbb3SSatish Balay         rval      = ib[j]*4;
11012d61bbb3SSatish Balay         z[rval++] +=  v[0]*x1 +  v[1]*x2 +  v[2]*x3 +  v[3]*x4;
11022d61bbb3SSatish Balay         z[rval++] +=  v[4]*x1 +  v[5]*x2 +  v[6]*x3 +  v[7]*x4;
11032d61bbb3SSatish Balay         z[rval++] +=  v[8]*x1 +  v[9]*x2 + v[10]*x3 + v[11]*x4;
11042d61bbb3SSatish Balay         z[rval++] += v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4;
11052d61bbb3SSatish Balay         v += 16;
11062d61bbb3SSatish Balay       }
11072d61bbb3SSatish Balay     }
11082d61bbb3SSatish Balay     break;
11092d61bbb3SSatish Balay   case 5:
11102d61bbb3SSatish Balay     for ( i=0; i<mbs; i++ ) {
11112d61bbb3SSatish Balay       n  = ii[1] - ii[0]; ii++;
11122d61bbb3SSatish Balay       xb = x + 5*i; x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
11132d61bbb3SSatish Balay       x4 = xb[3];   x5 = xb[4];
11142d61bbb3SSatish Balay       ib = idx + ai[i];
11152d61bbb3SSatish Balay       for ( j=0; j<n; j++ ) {
11162d61bbb3SSatish Balay         rval      = ib[j]*5;
11172d61bbb3SSatish Balay         z[rval++] +=  v[0]*x1 +  v[1]*x2 +  v[2]*x3 +  v[3]*x4 +  v[4]*x5;
11182d61bbb3SSatish Balay         z[rval++] +=  v[5]*x1 +  v[6]*x2 +  v[7]*x3 +  v[8]*x4 +  v[9]*x5;
11192d61bbb3SSatish Balay         z[rval++] += v[10]*x1 + v[11]*x2 + v[12]*x3 + v[13]*x4 + v[14]*x5;
11202d61bbb3SSatish Balay         z[rval++] += v[15]*x1 + v[16]*x2 + v[17]*x3 + v[18]*x4 + v[19]*x5;
11212d61bbb3SSatish Balay         z[rval++] += v[20]*x1 + v[21]*x2 + v[22]*x3 + v[23]*x4 + v[24]*x5;
11222d61bbb3SSatish Balay         v += 25;
11232d61bbb3SSatish Balay       }
11242d61bbb3SSatish Balay     }
11252d61bbb3SSatish Balay     break;
11262d61bbb3SSatish Balay       /* block sizes larger then 5 by 5 are handled by BLAS */
11272d61bbb3SSatish Balay     default: {
11283f1db9ecSBarry Smith       int       ncols,k;
11293f1db9ecSBarry Smith       Scalar    *work,*workt;
11303f1db9ecSBarry Smith 
11312d61bbb3SSatish Balay       if (!a->mult_work) {
11322d61bbb3SSatish Balay         k = PetscMax(a->m,a->n);
1133544c85edSSatish Balay         a->mult_work = (Scalar *) PetscMalloc((k+1)*sizeof(Scalar));CHKPTRQ(a->mult_work);
11342d61bbb3SSatish Balay       }
11352d61bbb3SSatish Balay       work = a->mult_work;
11362d61bbb3SSatish Balay       for ( i=0; i<mbs; i++ ) {
11372d61bbb3SSatish Balay         n     = ii[1] - ii[0]; ii++;
11382d61bbb3SSatish Balay         ncols = n*bs;
1139549d3d68SSatish Balay         ierr  = PetscMemzero(work,ncols*sizeof(Scalar));CHKERRQ(ierr);
11403f1db9ecSBarry Smith         Kernel_w_gets_w_plus_trans_Ar_times_v(bs,ncols,x,v,work);
11413f1db9ecSBarry Smith         /* LAgemv_("T",&bs,&ncols,&_DOne,v,&bs,x,&_One,&_DOne,work,&_One); */
11422d61bbb3SSatish Balay         v += n*bs2;
11432d61bbb3SSatish Balay         x += bs;
11442d61bbb3SSatish Balay         workt = work;
11452d61bbb3SSatish Balay         for ( j=0; j<n; j++ ) {
11462d61bbb3SSatish Balay           zb = z + bs*(*idx++);
11472d61bbb3SSatish Balay           for ( k=0; k<bs; k++ ) zb[k] += workt[k] ;
11482d61bbb3SSatish Balay           workt += bs;
11492d61bbb3SSatish Balay         }
11502d61bbb3SSatish Balay       }
11512d61bbb3SSatish Balay     }
11522d61bbb3SSatish Balay   }
11532d61bbb3SSatish Balay   ierr = VecRestoreArray(xx,&xg);CHKERRQ(ierr);
11542d61bbb3SSatish Balay   ierr = VecRestoreArray(zz,&zg);CHKERRQ(ierr);
11552d61bbb3SSatish Balay   PLogFlops(2*a->nz*a->bs2);
11562d61bbb3SSatish Balay   PetscFunctionReturn(0);
11572d61bbb3SSatish Balay }
11582d61bbb3SSatish Balay 
11592d61bbb3SSatish Balay #undef __FUNC__
11602d61bbb3SSatish Balay #define __FUNC__ "MatScale_SeqBAIJ"
11612d61bbb3SSatish Balay int MatScale_SeqBAIJ(Scalar *alpha,Mat inA)
11622d61bbb3SSatish Balay {
11632d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) inA->data;
11642d61bbb3SSatish Balay   int         one = 1, totalnz = a->bs2*a->nz;
11652d61bbb3SSatish Balay 
11662d61bbb3SSatish Balay   PetscFunctionBegin;
11672d61bbb3SSatish Balay   BLscal_( &totalnz, alpha, a->a, &one );
11682d61bbb3SSatish Balay   PLogFlops(totalnz);
11692d61bbb3SSatish Balay   PetscFunctionReturn(0);
11702d61bbb3SSatish Balay }
11712d61bbb3SSatish Balay 
11722d61bbb3SSatish Balay #undef __FUNC__
11732d61bbb3SSatish Balay #define __FUNC__ "MatNorm_SeqBAIJ"
11742d61bbb3SSatish Balay int MatNorm_SeqBAIJ(Mat A,NormType type,double *norm)
11752d61bbb3SSatish Balay {
11762d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
11773f1db9ecSBarry Smith   MatScalar   *v = a->a;
11782d61bbb3SSatish Balay   double      sum = 0.0;
1179253ffaf7SBarry Smith   int         i,j,k,bs = a->bs, nz=a->nz,bs2=a->bs2;
11802d61bbb3SSatish Balay 
11812d61bbb3SSatish Balay   PetscFunctionBegin;
11822d61bbb3SSatish Balay   if (type == NORM_FROBENIUS) {
11832d61bbb3SSatish Balay     for (i=0; i< bs2*nz; i++ ) {
1184aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
1185e20fef11SSatish Balay       sum += PetscReal(PetscConj(*v)*(*v)); v++;
11862d61bbb3SSatish Balay #else
11872d61bbb3SSatish Balay       sum += (*v)*(*v); v++;
11882d61bbb3SSatish Balay #endif
11892d61bbb3SSatish Balay     }
11902d61bbb3SSatish Balay     *norm = sqrt(sum);
1191596552b5SBarry Smith   }  else if (type == NORM_INFINITY) { /* maximum row sum */
1192596552b5SBarry Smith     *norm = 0.0;
1193596552b5SBarry Smith     for ( k=0; k<bs; k++ ) {
1194596552b5SBarry Smith       for ( j=0; j<a->m; j++ ) {
1195596552b5SBarry Smith         v = a->a + bs2*a->i[j] + k;
1196596552b5SBarry Smith         sum = 0.0;
1197596552b5SBarry Smith         for ( i=0; i<a->i[j+1]-a->i[j]; i++ ) {
1198596552b5SBarry Smith           sum += PetscAbsScalar(*v);
1199596552b5SBarry Smith           v   += bs;
12002d61bbb3SSatish Balay         }
1201596552b5SBarry Smith         if (sum > *norm) *norm = sum;
1202596552b5SBarry Smith       }
1203596552b5SBarry Smith     }
1204596552b5SBarry Smith   } else {
12052d61bbb3SSatish Balay     SETERRQ(PETSC_ERR_SUP,0,"No support for this norm yet");
12062d61bbb3SSatish Balay   }
12072d61bbb3SSatish Balay   PetscFunctionReturn(0);
12082d61bbb3SSatish Balay }
12092d61bbb3SSatish Balay 
12102d61bbb3SSatish Balay 
12112d61bbb3SSatish Balay #undef __FUNC__
12122d61bbb3SSatish Balay #define __FUNC__ "MatEqual_SeqBAIJ"
12132d61bbb3SSatish Balay int MatEqual_SeqBAIJ(Mat A,Mat B, PetscTruth* flg)
12142d61bbb3SSatish Balay {
12152d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)B->data;
12160f5bd95cSBarry Smith   int         ierr;
12172d61bbb3SSatish Balay 
12182d61bbb3SSatish Balay   PetscFunctionBegin;
12192d61bbb3SSatish Balay   if (B->type !=MATSEQBAIJ) SETERRQ(PETSC_ERR_ARG_INCOMP,0,"Matrices must be same type");
12202d61bbb3SSatish Balay 
12212d61bbb3SSatish Balay   /* If the  matrix/block dimensions are not equal, or no of nonzeros or shift */
12222d61bbb3SSatish Balay   if ((a->m != b->m) || (a->n !=b->n) || (a->bs != b->bs)|| (a->nz != b->nz)) {
12232d61bbb3SSatish Balay     *flg = PETSC_FALSE; PetscFunctionReturn(0);
12242d61bbb3SSatish Balay   }
12252d61bbb3SSatish Balay 
12262d61bbb3SSatish Balay   /* if the a->i are the same */
12270f5bd95cSBarry Smith   ierr = PetscMemcmp(a->i,b->i, (a->mbs+1)*sizeof(int),flg);CHKERRQ(ierr);
12280f5bd95cSBarry Smith   if (*flg == PETSC_FALSE) {
12290f5bd95cSBarry Smith     PetscFunctionReturn(0);
12302d61bbb3SSatish Balay   }
12312d61bbb3SSatish Balay 
12322d61bbb3SSatish Balay   /* if a->j are the same */
12330f5bd95cSBarry Smith   ierr = PetscMemcmp(a->j,b->j,(a->nz)*sizeof(int),flg);CHKERRQ(ierr);
12340f5bd95cSBarry Smith   if (*flg == PETSC_FALSE) {
12350f5bd95cSBarry Smith     PetscFunctionReturn(0);
12362d61bbb3SSatish Balay   }
12372d61bbb3SSatish Balay   /* if a->a are the same */
12380f5bd95cSBarry Smith   ierr = PetscMemcmp(a->a, b->a,(a->nz)*(a->bs)*(a->bs)*sizeof(Scalar),flg);CHKERRQ(ierr);
12392d61bbb3SSatish Balay   PetscFunctionReturn(0);
12402d61bbb3SSatish Balay 
12412d61bbb3SSatish Balay }
12422d61bbb3SSatish Balay 
12432d61bbb3SSatish Balay #undef __FUNC__
12442d61bbb3SSatish Balay #define __FUNC__ "MatGetDiagonal_SeqBAIJ"
12452d61bbb3SSatish Balay int MatGetDiagonal_SeqBAIJ(Mat A,Vec v)
12462d61bbb3SSatish Balay {
12472d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
1248e1311b90SBarry Smith   int         ierr,i,j,k,n,row,bs,*ai,*aj,ambs,bs2;
12493f1db9ecSBarry Smith   Scalar      *x, zero = 0.0;
12503f1db9ecSBarry Smith   MatScalar   *aa,*aa_j;
12512d61bbb3SSatish Balay 
12522d61bbb3SSatish Balay   PetscFunctionBegin;
12532d61bbb3SSatish Balay   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix");
12542d61bbb3SSatish Balay   bs   = a->bs;
12552d61bbb3SSatish Balay   aa   = a->a;
12562d61bbb3SSatish Balay   ai   = a->i;
12572d61bbb3SSatish Balay   aj   = a->j;
12582d61bbb3SSatish Balay   ambs = a->mbs;
12592d61bbb3SSatish Balay   bs2  = a->bs2;
12602d61bbb3SSatish Balay 
1261e1311b90SBarry Smith   ierr = VecSet(&zero,v);CHKERRQ(ierr);
1262e1311b90SBarry Smith   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1263e1311b90SBarry Smith   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
12642d61bbb3SSatish Balay   if (n != a->m) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Nonconforming matrix and vector");
12652d61bbb3SSatish Balay   for ( i=0; i<ambs; i++ ) {
12662d61bbb3SSatish Balay     for ( j=ai[i]; j<ai[i+1]; j++ ) {
12672d61bbb3SSatish Balay       if (aj[j] == i) {
12682d61bbb3SSatish Balay         row  = i*bs;
12692d61bbb3SSatish Balay         aa_j = aa+j*bs2;
12702d61bbb3SSatish Balay         for (k=0; k<bs2; k+=(bs+1),row++) x[row] = aa_j[k];
12712d61bbb3SSatish Balay         break;
12722d61bbb3SSatish Balay       }
12732d61bbb3SSatish Balay     }
12742d61bbb3SSatish Balay   }
127515091d37SBarry Smith   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
12762d61bbb3SSatish Balay   PetscFunctionReturn(0);
12772d61bbb3SSatish Balay }
12782d61bbb3SSatish Balay 
12792d61bbb3SSatish Balay #undef __FUNC__
12802d61bbb3SSatish Balay #define __FUNC__ "MatDiagonalScale_SeqBAIJ"
12812d61bbb3SSatish Balay int MatDiagonalScale_SeqBAIJ(Mat A,Vec ll,Vec rr)
12822d61bbb3SSatish Balay {
12832d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
12843f1db9ecSBarry Smith   Scalar      *l,*r,x,*li,*ri;
12853f1db9ecSBarry Smith   MatScalar   *aa,*v;
1286e1311b90SBarry Smith   int         ierr,i,j,k,lm,rn,M,m,n,*ai,*aj,mbs,tmp,bs,bs2;
12872d61bbb3SSatish Balay 
12882d61bbb3SSatish Balay   PetscFunctionBegin;
12892d61bbb3SSatish Balay   ai  = a->i;
12902d61bbb3SSatish Balay   aj  = a->j;
12912d61bbb3SSatish Balay   aa  = a->a;
12922d61bbb3SSatish Balay   m   = a->m;
12932d61bbb3SSatish Balay   n   = a->n;
12942d61bbb3SSatish Balay   bs  = a->bs;
12952d61bbb3SSatish Balay   mbs = a->mbs;
12962d61bbb3SSatish Balay   bs2 = a->bs2;
12972d61bbb3SSatish Balay   if (ll) {
1298e1311b90SBarry Smith     ierr = VecGetArray(ll,&l);CHKERRQ(ierr);
1299e1311b90SBarry Smith     ierr = VecGetLocalSize(ll,&lm);CHKERRQ(ierr);
13002d61bbb3SSatish Balay     if (lm != m) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Left scaling vector wrong length");
13012d61bbb3SSatish Balay     for ( i=0; i<mbs; i++ ) { /* for each block row */
13022d61bbb3SSatish Balay       M  = ai[i+1] - ai[i];
13032d61bbb3SSatish Balay       li = l + i*bs;
13042d61bbb3SSatish Balay       v  = aa + bs2*ai[i];
13052d61bbb3SSatish Balay       for ( j=0; j<M; j++ ) { /* for each block */
13062d61bbb3SSatish Balay         for ( k=0; k<bs2; k++ ) {
13072d61bbb3SSatish Balay           (*v++) *= li[k%bs];
13082d61bbb3SSatish Balay         }
13092d61bbb3SSatish Balay       }
13102d61bbb3SSatish Balay     }
131160d550acSSatish Balay     ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr);
131296896c26SBarry Smith     PLogFlops(a->nz);
13132d61bbb3SSatish Balay   }
13142d61bbb3SSatish Balay 
13152d61bbb3SSatish Balay   if (rr) {
1316e1311b90SBarry Smith     ierr = VecGetArray(rr,&r);CHKERRQ(ierr);
1317e1311b90SBarry Smith     ierr = VecGetLocalSize(rr,&rn);CHKERRQ(ierr);
13182d61bbb3SSatish Balay     if (rn != n) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Right scaling vector wrong length");
13192d61bbb3SSatish Balay     for ( i=0; i<mbs; i++ ) { /* for each block row */
13202d61bbb3SSatish Balay       M  = ai[i+1] - ai[i];
13212d61bbb3SSatish Balay       v  = aa + bs2*ai[i];
13222d61bbb3SSatish Balay       for ( j=0; j<M; j++ ) { /* for each block */
13232d61bbb3SSatish Balay         ri = r + bs*aj[ai[i]+j];
13242d61bbb3SSatish Balay         for ( k=0; k<bs; k++ ) {
13252d61bbb3SSatish Balay           x = ri[k];
13262d61bbb3SSatish Balay           for ( tmp=0; tmp<bs; tmp++ ) (*v++) *= x;
13272d61bbb3SSatish Balay         }
13282d61bbb3SSatish Balay       }
13292d61bbb3SSatish Balay     }
133060d550acSSatish Balay     ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr);
133196896c26SBarry Smith     PLogFlops(a->nz);
13322d61bbb3SSatish Balay   }
13332d61bbb3SSatish Balay   PetscFunctionReturn(0);
13342d61bbb3SSatish Balay }
13352d61bbb3SSatish Balay 
13362d61bbb3SSatish Balay 
13372d61bbb3SSatish Balay #undef __FUNC__
13382d61bbb3SSatish Balay #define __FUNC__ "MatGetInfo_SeqBAIJ"
13392d61bbb3SSatish Balay int MatGetInfo_SeqBAIJ(Mat A,MatInfoType flag,MatInfo *info)
13402d61bbb3SSatish Balay {
13412d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
13422d61bbb3SSatish Balay 
13432d61bbb3SSatish Balay   PetscFunctionBegin;
13442d61bbb3SSatish Balay   info->rows_global    = (double)a->m;
13452d61bbb3SSatish Balay   info->columns_global = (double)a->n;
13462d61bbb3SSatish Balay   info->rows_local     = (double)a->m;
13472d61bbb3SSatish Balay   info->columns_local  = (double)a->n;
13482d61bbb3SSatish Balay   info->block_size     = a->bs2;
13492d61bbb3SSatish Balay   info->nz_allocated   = a->maxnz;
13502d61bbb3SSatish Balay   info->nz_used        = a->bs2*a->nz;
13512d61bbb3SSatish Balay   info->nz_unneeded    = (double)(info->nz_allocated - info->nz_used);
13522d61bbb3SSatish Balay   info->assemblies   = A->num_ass;
13532d61bbb3SSatish Balay   info->mallocs      = a->reallocs;
13542d61bbb3SSatish Balay   info->memory       = A->mem;
13552d61bbb3SSatish Balay   if (A->factor) {
13562d61bbb3SSatish Balay     info->fill_ratio_given  = A->info.fill_ratio_given;
13572d61bbb3SSatish Balay     info->fill_ratio_needed = A->info.fill_ratio_needed;
13582d61bbb3SSatish Balay     info->factor_mallocs    = A->info.factor_mallocs;
13592d61bbb3SSatish Balay   } else {
13602d61bbb3SSatish Balay     info->fill_ratio_given  = 0;
13612d61bbb3SSatish Balay     info->fill_ratio_needed = 0;
13622d61bbb3SSatish Balay     info->factor_mallocs    = 0;
13632d61bbb3SSatish Balay   }
13642d61bbb3SSatish Balay   PetscFunctionReturn(0);
13652d61bbb3SSatish Balay }
13662d61bbb3SSatish Balay 
13672d61bbb3SSatish Balay 
13682d61bbb3SSatish Balay #undef __FUNC__
13692d61bbb3SSatish Balay #define __FUNC__ "MatZeroEntries_SeqBAIJ"
13702d61bbb3SSatish Balay int MatZeroEntries_SeqBAIJ(Mat A)
13712d61bbb3SSatish Balay {
13722d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
1373549d3d68SSatish Balay   int         ierr;
13742d61bbb3SSatish Balay 
13752d61bbb3SSatish Balay   PetscFunctionBegin;
1376549d3d68SSatish Balay   ierr = PetscMemzero(a->a,a->bs2*a->i[a->mbs]*sizeof(MatScalar));CHKERRQ(ierr);
13772d61bbb3SSatish Balay   PetscFunctionReturn(0);
13782d61bbb3SSatish Balay }
1379