xref: /petsc/src/mat/impls/baij/seq/baij2.c (revision 606d414c6e034e02e67059b83ebaefc3ebe99698)
1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER
2*606d414cSSatish Balay static char vcid[] = "$Id: baij2.c,v 1.47 1999/06/08 22:56:09 balay Exp balay $";
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;
1942523b56SSatish Balay   BT          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 
2942523b56SSatish Balay   ierr  = BTCreate(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;
3642523b56SSatish Balay     BTMemzero(m,table);
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");
4642523b56SSatish Balay       if(!BTLookupSet(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];
6042523b56SSatish Balay           if (!BTLookupSet(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   }
7142523b56SSatish Balay   BTDestroy(table);
72*606d414cSSatish Balay   ierr = PetscFree(nidx);CHKERRQ(ierr);
73*606d414cSSatish 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;
88736121d4SSatish Balay 
893a40ed3dSBarry Smith   PetscFunctionBegin;
90888f2ed8SSatish Balay   ierr = ISSorted(iscol,(PetscTruth*)&i);CHKERRQ(ierr);
91a8c6a408SBarry Smith   if (!i) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"IS is not sorted");
92736121d4SSatish Balay 
93736121d4SSatish Balay   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
94218c64b6SSatish Balay   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
95736121d4SSatish Balay   ierr = ISGetSize(isrow,&nrows);CHKERRQ(ierr);
96736121d4SSatish Balay   ierr = ISGetSize(iscol,&ncols);CHKERRQ(ierr);
97736121d4SSatish Balay 
98736121d4SSatish Balay   smap  = (int *) PetscMalloc((1+oldcols)*sizeof(int));CHKPTRQ(smap);
99736121d4SSatish Balay   ssmap = smap;
100736121d4SSatish Balay   lens  = (int *) PetscMalloc((1+nrows)*sizeof(int));CHKPTRQ(lens);
101549d3d68SSatish Balay   ierr  = PetscMemzero(smap,oldcols*sizeof(int));CHKERRQ(ierr);
102736121d4SSatish Balay   for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1;
103736121d4SSatish Balay   /* determine lens of each row */
104736121d4SSatish Balay   for (i=0; i<nrows; i++) {
105736121d4SSatish Balay     kstart  = ai[irow[i]];
106736121d4SSatish Balay     kend    = kstart + a->ilen[irow[i]];
107736121d4SSatish Balay     lens[i] = 0;
108736121d4SSatish Balay       for ( k=kstart; k<kend; k++ ) {
109736121d4SSatish Balay         if (ssmap[aj[k]]) {
110736121d4SSatish Balay           lens[i]++;
111736121d4SSatish Balay         }
112736121d4SSatish Balay       }
113736121d4SSatish Balay     }
114736121d4SSatish Balay   /* Create and fill new matrix */
115736121d4SSatish Balay   if (scall == MAT_REUSE_MATRIX) {
116736121d4SSatish Balay     c = (Mat_SeqBAIJ *)((*B)->data);
117736121d4SSatish Balay 
118a8c6a408SBarry Smith     if (c->mbs!=nrows || c->nbs!=ncols || c->bs!=bs) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Submatrix wrong size");
119736121d4SSatish Balay     if (PetscMemcmp(c->ilen,lens, c->mbs *sizeof(int))) {
120a8c6a408SBarry Smith       SETERRQ(PETSC_ERR_ARG_SIZ,0,"Cannot reuse matrix. wrong no of nonzeros");
121736121d4SSatish Balay     }
122549d3d68SSatish Balay     ierr = PetscMemzero(c->ilen,c->mbs*sizeof(int));CHKERRQ(ierr);
123736121d4SSatish Balay     C = *B;
1243a40ed3dSBarry Smith   } else {
125736121d4SSatish Balay     ierr = MatCreateSeqBAIJ(A->comm,bs,nrows*bs,ncols*bs,0,lens,&C);CHKERRQ(ierr);
126736121d4SSatish Balay   }
127736121d4SSatish Balay   c = (Mat_SeqBAIJ *)(C->data);
128736121d4SSatish Balay   for (i=0; i<nrows; i++) {
129736121d4SSatish Balay     row    = irow[i];
130736121d4SSatish Balay     kstart = ai[row];
131736121d4SSatish Balay     kend   = kstart + a->ilen[row];
132736121d4SSatish Balay     mat_i  = c->i[i];
133736121d4SSatish Balay     mat_j  = c->j + mat_i;
134218c64b6SSatish Balay     mat_a  = c->a + mat_i*bs2;
135736121d4SSatish Balay     mat_ilen = c->ilen + i;
136736121d4SSatish Balay     for ( k=kstart; k<kend; k++ ) {
137736121d4SSatish Balay       if ((tcol=ssmap[a->j[k]])) {
138736121d4SSatish Balay         *mat_j++ = tcol - 1;
139549d3d68SSatish Balay         ierr     = PetscMemcpy(mat_a,a->a+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
140549d3d68SSatish Balay         mat_a   += bs2;
141736121d4SSatish Balay         (*mat_ilen)++;
142736121d4SSatish Balay       }
143736121d4SSatish Balay     }
144736121d4SSatish Balay   }
145218c64b6SSatish Balay 
146736121d4SSatish Balay   /* Free work space */
147736121d4SSatish Balay   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
148*606d414cSSatish Balay   ierr = PetscFree(smap);CHKERRQ(ierr);
149*606d414cSSatish Balay   ierr = PetscFree(lens);CHKERRQ(ierr);
1506d4a8577SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1516d4a8577SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
152736121d4SSatish Balay 
153736121d4SSatish Balay   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
154736121d4SSatish Balay   *B = C;
1553a40ed3dSBarry Smith   PetscFunctionReturn(0);
156736121d4SSatish Balay }
157736121d4SSatish Balay 
1585615d1e5SSatish Balay #undef __FUNC__
1595615d1e5SSatish Balay #define __FUNC__ "MatGetSubMatrix_SeqBAIJ"
1607b2a1423SBarry Smith int MatGetSubMatrix_SeqBAIJ(Mat A,IS isrow,IS iscol,int cs,MatReuse scall,Mat *B)
161218c64b6SSatish Balay {
162218c64b6SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
163218c64b6SSatish Balay   IS          is1,is2;
164218c64b6SSatish Balay   int         *vary,*iary,*irow,*icol,nrows,ncols,i,ierr,bs=a->bs,count;
165218c64b6SSatish Balay 
1663a40ed3dSBarry Smith   PetscFunctionBegin;
167218c64b6SSatish Balay   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
168218c64b6SSatish Balay   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
169218c64b6SSatish Balay   ierr = ISGetSize(isrow,&nrows);CHKERRQ(ierr);
170218c64b6SSatish Balay   ierr = ISGetSize(iscol,&ncols);CHKERRQ(ierr);
171218c64b6SSatish Balay 
172218c64b6SSatish Balay   /* Verify if the indices corespond to each element in a block
173218c64b6SSatish Balay    and form the IS with compressed IS */
174218c64b6SSatish Balay   vary = (int *) PetscMalloc(2*(a->mbs+1)*sizeof(int));CHKPTRQ(vary);
175218c64b6SSatish Balay   iary = vary + a->mbs;
176549d3d68SSatish Balay   ierr = PetscMemzero(vary,(a->mbs)*sizeof(int));CHKERRQ(ierr);
177218c64b6SSatish Balay   for ( i=0; i<nrows; i++) vary[irow[i]/bs]++;
178218c64b6SSatish Balay   count = 0;
179218c64b6SSatish Balay   for (i=0; i<a->mbs; i++) {
1806a6a5d1dSBarry Smith     if (vary[i]!=0 && vary[i]!=bs) SETERRA(1,0,"Index set does not match blocks");
181218c64b6SSatish Balay     if (vary[i]==bs) iary[count++] = i;
182218c64b6SSatish Balay   }
183029af93fSBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF, count, iary,&is1);CHKERRQ(ierr);
184218c64b6SSatish Balay 
185549d3d68SSatish Balay   ierr = PetscMemzero(vary,(a->mbs)*sizeof(int));CHKERRQ(ierr);
186218c64b6SSatish Balay   for ( i=0; i<ncols; i++) vary[icol[i]/bs]++;
187218c64b6SSatish Balay   count = 0;
188218c64b6SSatish Balay   for (i=0; i<a->mbs; i++) {
189e3372554SBarry Smith     if (vary[i]!=0 && vary[i]!=bs) SETERRA(1,0,"MatGetSubmatrices_SeqBAIJ:");
190218c64b6SSatish Balay     if (vary[i]==bs) iary[count++] = i;
191218c64b6SSatish Balay   }
192029af93fSBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF, count, iary,&is2);CHKERRQ(ierr);
193218c64b6SSatish Balay   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
194218c64b6SSatish Balay   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
195*606d414cSSatish Balay   ierr = PetscFree(vary);CHKERRQ(ierr);
196218c64b6SSatish Balay 
1976a6a5d1dSBarry Smith   ierr = MatGetSubMatrix_SeqBAIJ_Private(A,is1,is2,cs,scall,B);CHKERRQ(ierr);
198218c64b6SSatish Balay   ISDestroy(is1);
199218c64b6SSatish Balay   ISDestroy(is2);
2003a40ed3dSBarry Smith   PetscFunctionReturn(0);
201218c64b6SSatish Balay }
202218c64b6SSatish Balay 
2035615d1e5SSatish Balay #undef __FUNC__
2045615d1e5SSatish Balay #define __FUNC__ "MatGetSubMatrices_SeqBAIJ"
2057b2a1423SBarry Smith int MatGetSubMatrices_SeqBAIJ(Mat A,int n, IS *irow,IS *icol,MatReuse scall,Mat **B)
206736121d4SSatish Balay {
207736121d4SSatish Balay   int ierr,i;
208736121d4SSatish Balay 
2093a40ed3dSBarry Smith   PetscFunctionBegin;
210736121d4SSatish Balay   if (scall == MAT_INITIAL_MATRIX) {
211736121d4SSatish Balay     *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) );CHKPTRQ(*B);
212736121d4SSatish Balay   }
213736121d4SSatish Balay 
214736121d4SSatish Balay   for ( i=0; i<n; i++ ) {
2156a6a5d1dSBarry Smith     ierr = MatGetSubMatrix_SeqBAIJ(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr);
216736121d4SSatish Balay   }
2173a40ed3dSBarry Smith   PetscFunctionReturn(0);
218736121d4SSatish Balay }
219218c64b6SSatish Balay 
220218c64b6SSatish Balay 
2212d61bbb3SSatish Balay /* -------------------------------------------------------*/
2222d61bbb3SSatish Balay /* Should check that shapes of vectors and matrices match */
2232d61bbb3SSatish Balay /* -------------------------------------------------------*/
224eadb2fb4SBarry Smith #include "pinclude/blaslapack.h"
2252d61bbb3SSatish Balay 
2262d61bbb3SSatish Balay #undef __FUNC__
2272d61bbb3SSatish Balay #define __FUNC__ "MatMult_SeqBAIJ_1"
2282d61bbb3SSatish Balay int MatMult_SeqBAIJ_1(Mat A,Vec xx,Vec zz)
2292d61bbb3SSatish Balay {
2302d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
2313f1db9ecSBarry Smith   Scalar          *x,*z,sum;
2323f1db9ecSBarry Smith   MatScalar       *v;
233e1311b90SBarry Smith   int             mbs=a->mbs,i,*idx,*ii,n,ierr;
2342d61bbb3SSatish Balay 
2352d61bbb3SSatish Balay   PetscFunctionBegin;
236e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
237e1311b90SBarry Smith   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
2382d61bbb3SSatish Balay 
2392d61bbb3SSatish Balay   idx   = a->j;
2402d61bbb3SSatish Balay   v     = a->a;
2412d61bbb3SSatish Balay   ii    = a->i;
2422d61bbb3SSatish Balay 
2432d61bbb3SSatish Balay   for ( i=0; i<mbs; i++ ) {
2442d61bbb3SSatish Balay     n    = ii[1] - ii[0]; ii++;
2452d61bbb3SSatish Balay     sum  = 0.0;
2462d61bbb3SSatish Balay     while (n--) sum += *v++ * x[*idx++];
2472d61bbb3SSatish Balay     z[i] = sum;
2482d61bbb3SSatish Balay   }
249e1311b90SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
250e1311b90SBarry Smith   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
2512d61bbb3SSatish Balay   PLogFlops(2*a->nz - a->m);
2522d61bbb3SSatish Balay   PetscFunctionReturn(0);
2532d61bbb3SSatish Balay }
2542d61bbb3SSatish Balay 
2552d61bbb3SSatish Balay #undef __FUNC__
2562d61bbb3SSatish Balay #define __FUNC__ "MatMult_SeqBAIJ_2"
2572d61bbb3SSatish Balay int MatMult_SeqBAIJ_2(Mat A,Vec xx,Vec zz)
2582d61bbb3SSatish Balay {
2592d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
2603f1db9ecSBarry Smith   Scalar          *x,*z,*xb,sum1,sum2;
261e1311b90SBarry Smith   Scalar          x1,x2;
2623f1db9ecSBarry Smith   MatScalar       *v;
263e1311b90SBarry Smith   int             ierr,mbs=a->mbs,i,*idx,*ii,j,n;
2642d61bbb3SSatish Balay 
2652d61bbb3SSatish Balay   PetscFunctionBegin;
266e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
267e1311b90SBarry Smith   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
2682d61bbb3SSatish Balay 
2692d61bbb3SSatish Balay   idx   = a->j;
2702d61bbb3SSatish Balay   v     = a->a;
2712d61bbb3SSatish Balay   ii    = a->i;
2722d61bbb3SSatish Balay 
2732d61bbb3SSatish Balay   for ( i=0; i<mbs; i++ ) {
2742d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
2752d61bbb3SSatish Balay     sum1 = 0.0; sum2 = 0.0;
2762d61bbb3SSatish Balay     for ( j=0; j<n; j++ ) {
2772d61bbb3SSatish Balay       xb = x + 2*(*idx++); x1 = xb[0]; x2 = xb[1];
2782d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[2]*x2;
2792d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[3]*x2;
2802d61bbb3SSatish Balay       v += 4;
2812d61bbb3SSatish Balay     }
2822d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2;
2832d61bbb3SSatish Balay     z += 2;
2842d61bbb3SSatish Balay   }
285e1311b90SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
286e1311b90SBarry Smith   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
28715091d37SBarry Smith   PLogFlops(8*a->nz - a->m);
2882d61bbb3SSatish Balay   PetscFunctionReturn(0);
2892d61bbb3SSatish Balay }
2902d61bbb3SSatish Balay 
2912d61bbb3SSatish Balay #undef __FUNC__
2922d61bbb3SSatish Balay #define __FUNC__ "MatMult_SeqBAIJ_3"
2932d61bbb3SSatish Balay int MatMult_SeqBAIJ_3(Mat A,Vec xx,Vec zz)
2942d61bbb3SSatish Balay {
2952d61bbb3SSatish Balay   Mat_SeqBAIJ  *a = (Mat_SeqBAIJ *) A->data;
2963f1db9ecSBarry Smith   Scalar       *x,*z,*xb,sum1,sum2,sum3,x1,x2,x3;
2973f1db9ecSBarry Smith   MatScalar    *v;
298e1311b90SBarry Smith   int          ierr,mbs=a->mbs,i,*idx,*ii,j,n;
2992d61bbb3SSatish Balay 
300aa482453SBarry Smith #if defined(PETSC_HAVE_PRAGMA_DISJOINT)
301fee21e36SBarry Smith #pragma disjoint(*v,*z,*xb)
302fee21e36SBarry Smith #endif
303fee21e36SBarry Smith 
3042d61bbb3SSatish Balay   PetscFunctionBegin;
305e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
306e1311b90SBarry Smith   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
3072d61bbb3SSatish Balay 
3082d61bbb3SSatish Balay   idx   = a->j;
3092d61bbb3SSatish Balay   v     = a->a;
3102d61bbb3SSatish Balay   ii    = a->i;
3112d61bbb3SSatish Balay 
3122d61bbb3SSatish Balay   for ( i=0; i<mbs; i++ ) {
3132d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
3142d61bbb3SSatish Balay     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0;
3152d61bbb3SSatish Balay     for ( j=0; j<n; j++ ) {
3162d61bbb3SSatish Balay       xb = x + 3*(*idx++); x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
3172d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[3]*x2 + v[6]*x3;
3182d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[4]*x2 + v[7]*x3;
3192d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[5]*x2 + v[8]*x3;
3202d61bbb3SSatish Balay       v += 9;
3212d61bbb3SSatish Balay     }
3222d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3;
3232d61bbb3SSatish Balay     z += 3;
3242d61bbb3SSatish Balay   }
325e1311b90SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
326e1311b90SBarry Smith   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
3272d61bbb3SSatish Balay   PLogFlops(18*a->nz - a->m);
3282d61bbb3SSatish Balay   PetscFunctionReturn(0);
3292d61bbb3SSatish Balay }
3302d61bbb3SSatish Balay 
3312d61bbb3SSatish Balay #undef __FUNC__
3322d61bbb3SSatish Balay #define __FUNC__ "MatMult_SeqBAIJ_4"
3332d61bbb3SSatish Balay int MatMult_SeqBAIJ_4(Mat A,Vec xx,Vec zz)
3342d61bbb3SSatish Balay {
3352d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
3363f1db9ecSBarry Smith   Scalar          *x,*z,*xb,sum1,sum2,sum3,sum4,x1,x2,x3,x4;
3373f1db9ecSBarry Smith   MatScalar       *v;
338e1311b90SBarry Smith   int             ierr,mbs=a->mbs,i,*idx,*ii,j,n;
3392d61bbb3SSatish Balay 
3402d61bbb3SSatish Balay   PetscFunctionBegin;
341e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
342e1311b90SBarry Smith   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
3432d61bbb3SSatish Balay 
3442d61bbb3SSatish Balay   idx   = a->j;
3452d61bbb3SSatish Balay   v     = a->a;
3462d61bbb3SSatish Balay   ii    = a->i;
3472d61bbb3SSatish Balay 
3482d61bbb3SSatish Balay   for ( i=0; i<mbs; i++ ) {
3492d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
3502d61bbb3SSatish Balay     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0;
3512d61bbb3SSatish Balay     for ( j=0; j<n; j++ ) {
3522d61bbb3SSatish Balay       xb = x + 4*(*idx++);
3532d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
3542d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[4]*x2 + v[8]*x3   + v[12]*x4;
3552d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[5]*x2 + v[9]*x3   + v[13]*x4;
3562d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[6]*x2 + v[10]*x3  + v[14]*x4;
3572d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[7]*x2 + v[11]*x3  + v[15]*x4;
3582d61bbb3SSatish Balay       v += 16;
3592d61bbb3SSatish Balay     }
3602d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4;
3612d61bbb3SSatish Balay     z += 4;
3622d61bbb3SSatish Balay   }
363e1311b90SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
364e1311b90SBarry Smith   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
3652d61bbb3SSatish Balay   PLogFlops(32*a->nz - a->m);
3662d61bbb3SSatish Balay   PetscFunctionReturn(0);
3672d61bbb3SSatish Balay }
3682d61bbb3SSatish Balay 
3692d61bbb3SSatish Balay #undef __FUNC__
3702d61bbb3SSatish Balay #define __FUNC__ "MatMult_SeqBAIJ_5"
3712d61bbb3SSatish Balay int MatMult_SeqBAIJ_5(Mat A,Vec xx,Vec zz)
3722d61bbb3SSatish Balay {
3732d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
3743f1db9ecSBarry Smith   Scalar          sum1,sum2,sum3,sum4,sum5,x1,x2,x3,x4,x5,*xb,*z,*x;
3753f1db9ecSBarry Smith   MatScalar       *v;
376e1311b90SBarry Smith   int             ierr,mbs=a->mbs,i,*idx,*ii,j,n;
3772d61bbb3SSatish Balay 
378e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
379e1311b90SBarry Smith   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
3802d61bbb3SSatish Balay 
3812d61bbb3SSatish Balay   idx   = a->j;
3822d61bbb3SSatish Balay   v     = a->a;
3832d61bbb3SSatish Balay   ii    = a->i;
3842d61bbb3SSatish Balay 
3852d61bbb3SSatish Balay   for ( i=0; i<mbs; i++ ) {
3862d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
3872d61bbb3SSatish Balay     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0;
3882d61bbb3SSatish Balay     for ( j=0; j<n; j++ ) {
3892d61bbb3SSatish Balay       xb = x + 5*(*idx++);
3902d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4];
3912d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[5]*x2 + v[10]*x3  + v[15]*x4 + v[20]*x5;
3922d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[6]*x2 + v[11]*x3  + v[16]*x4 + v[21]*x5;
3932d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[7]*x2 + v[12]*x3  + v[17]*x4 + v[22]*x5;
3942d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[8]*x2 + v[13]*x3  + v[18]*x4 + v[23]*x5;
3952d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[9]*x2 + v[14]*x3  + v[19]*x4 + v[24]*x5;
3962d61bbb3SSatish Balay       v += 25;
3972d61bbb3SSatish Balay     }
3982d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5;
3992d61bbb3SSatish Balay     z += 5;
4002d61bbb3SSatish Balay   }
401e1311b90SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
402e1311b90SBarry Smith   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
4032d61bbb3SSatish Balay   PLogFlops(50*a->nz - a->m);
4042d61bbb3SSatish Balay   PetscFunctionReturn(0);
4052d61bbb3SSatish Balay }
4062d61bbb3SSatish Balay 
40715091d37SBarry Smith 
40815091d37SBarry Smith #undef __FUNC__
40915091d37SBarry Smith #define __FUNC__ "MatMult_SeqBAIJ_6"
41015091d37SBarry Smith int MatMult_SeqBAIJ_6(Mat A,Vec xx,Vec zz)
41115091d37SBarry Smith {
41215091d37SBarry Smith   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
41315091d37SBarry Smith   Scalar          *x,*z,*xb,sum1,sum2,sum3,sum4,sum5,sum6;
41415091d37SBarry Smith   Scalar          x1,x2,x3,x4,x5,x6;
41515091d37SBarry Smith   MatScalar       *v;
41615091d37SBarry Smith   int             ierr,mbs=a->mbs,i,*idx,*ii,j,n;
41715091d37SBarry Smith 
41815091d37SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
41915091d37SBarry Smith   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
42015091d37SBarry Smith 
42115091d37SBarry Smith   idx   = a->j;
42215091d37SBarry Smith   v     = a->a;
42315091d37SBarry Smith   ii    = a->i;
42415091d37SBarry Smith 
42515091d37SBarry Smith   for ( i=0; i<mbs; i++ ) {
42615091d37SBarry Smith     n  = ii[1] - ii[0]; ii++;
42715091d37SBarry Smith     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0;
42815091d37SBarry Smith     for ( j=0; j<n; j++ ) {
42915091d37SBarry Smith       xb = x + 6*(*idx++);
43015091d37SBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5];
43115091d37SBarry Smith       sum1 += v[0]*x1 + v[6]*x2  + v[12]*x3  + v[18]*x4 + v[24]*x5 + v[30]*x6;
43215091d37SBarry Smith       sum2 += v[1]*x1 + v[7]*x2  + v[13]*x3  + v[19]*x4 + v[25]*x5 + v[31]*x6;
43315091d37SBarry Smith       sum3 += v[2]*x1 + v[8]*x2  + v[14]*x3  + v[20]*x4 + v[26]*x5 + v[32]*x6;
43415091d37SBarry Smith       sum4 += v[3]*x1 + v[9]*x2  + v[15]*x3  + v[21]*x4 + v[27]*x5 + v[33]*x6;
43515091d37SBarry Smith       sum5 += v[4]*x1 + v[10]*x2 + v[16]*x3  + v[22]*x4 + v[28]*x5 + v[34]*x6;
43615091d37SBarry Smith       sum6 += v[5]*x1 + v[11]*x2 + v[17]*x3  + v[23]*x4 + v[29]*x5 + v[35]*x6;
43715091d37SBarry Smith       v += 36;
43815091d37SBarry Smith     }
43915091d37SBarry Smith     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6;
44015091d37SBarry Smith     z += 6;
44115091d37SBarry Smith   }
44215091d37SBarry Smith 
44315091d37SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
44415091d37SBarry Smith   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
44515091d37SBarry Smith   PLogFlops(72*a->nz - a->m);
44615091d37SBarry Smith   PetscFunctionReturn(0);
44715091d37SBarry Smith }
4482d61bbb3SSatish Balay #undef __FUNC__
4492d61bbb3SSatish Balay #define __FUNC__ "MatMult_SeqBAIJ_7"
4502d61bbb3SSatish Balay int MatMult_SeqBAIJ_7(Mat A,Vec xx,Vec zz)
4512d61bbb3SSatish Balay {
4522d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
4533f1db9ecSBarry Smith   Scalar          *x,*z,*xb,sum1,sum2,sum3,sum4,sum5,sum6,sum7;
454e1311b90SBarry Smith   Scalar          x1,x2,x3,x4,x5,x6,x7;
4553f1db9ecSBarry Smith   MatScalar       *v;
456e1311b90SBarry Smith   int             ierr,mbs=a->mbs,i,*idx,*ii,j,n;
4572d61bbb3SSatish Balay 
458e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
459e1311b90SBarry Smith   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
4602d61bbb3SSatish Balay 
4612d61bbb3SSatish Balay   idx   = a->j;
4622d61bbb3SSatish Balay   v     = a->a;
4632d61bbb3SSatish Balay   ii    = a->i;
4642d61bbb3SSatish Balay 
4652d61bbb3SSatish Balay   for ( i=0; i<mbs; i++ ) {
4662d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
4672d61bbb3SSatish Balay     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0; sum7 = 0.0;
4682d61bbb3SSatish Balay     for ( j=0; j<n; j++ ) {
4692d61bbb3SSatish Balay       xb = x + 7*(*idx++);
4702d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6];
4712d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[7]*x2  + v[14]*x3  + v[21]*x4 + v[28]*x5 + v[35]*x6 + v[42]*x7;
4722d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[8]*x2  + v[15]*x3  + v[22]*x4 + v[29]*x5 + v[36]*x6 + v[43]*x7;
4732d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[9]*x2  + v[16]*x3  + v[23]*x4 + v[30]*x5 + v[37]*x6 + v[44]*x7;
4742d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[10]*x2 + v[17]*x3  + v[24]*x4 + v[31]*x5 + v[38]*x6 + v[45]*x7;
4752d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[11]*x2 + v[18]*x3  + v[25]*x4 + v[32]*x5 + v[39]*x6 + v[46]*x7;
4762d61bbb3SSatish Balay       sum6 += v[5]*x1 + v[12]*x2 + v[19]*x3  + v[26]*x4 + v[33]*x5 + v[40]*x6 + v[47]*x7;
4772d61bbb3SSatish Balay       sum7 += v[6]*x1 + v[13]*x2 + v[20]*x3  + v[27]*x4 + v[34]*x5 + v[41]*x6 + v[48]*x7;
4782d61bbb3SSatish Balay       v += 49;
4792d61bbb3SSatish Balay     }
4802d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
4812d61bbb3SSatish Balay     z += 7;
4822d61bbb3SSatish Balay   }
4832d61bbb3SSatish Balay 
484e1311b90SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
485e1311b90SBarry Smith   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
4862d61bbb3SSatish Balay   PLogFlops(98*a->nz - a->m);
4872d61bbb3SSatish Balay   PetscFunctionReturn(0);
4882d61bbb3SSatish Balay }
4892d61bbb3SSatish Balay 
4903f1db9ecSBarry Smith /*
4913f1db9ecSBarry Smith     This will not work with MatScalar == float because it calls the BLAS
4923f1db9ecSBarry Smith */
4932d61bbb3SSatish Balay #undef __FUNC__
4942d61bbb3SSatish Balay #define __FUNC__ "MatMult_SeqBAIJ_N"
4952d61bbb3SSatish Balay int MatMult_SeqBAIJ_N(Mat A,Vec xx,Vec zz)
4962d61bbb3SSatish Balay {
4972d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
4983f1db9ecSBarry Smith   Scalar          *x,*z,*xb, *work,*workt;
4993f1db9ecSBarry Smith   MatScalar       *v;
500e1311b90SBarry Smith   int             ierr,mbs=a->mbs,i,*idx,*ii,bs=a->bs,j,n,bs2=a->bs2;
5013f1db9ecSBarry Smith   int             ncols,k;
5022d61bbb3SSatish Balay 
5032d61bbb3SSatish Balay   PetscFunctionBegin;
504e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
505e1311b90SBarry Smith   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
5062d61bbb3SSatish Balay 
5072d61bbb3SSatish Balay   idx   = a->j;
5082d61bbb3SSatish Balay   v     = a->a;
5092d61bbb3SSatish Balay   ii    = a->i;
510218c64b6SSatish Balay 
511218c64b6SSatish Balay 
5122d61bbb3SSatish Balay   if (!a->mult_work) {
5132d61bbb3SSatish Balay     k = PetscMax(a->m,a->n);
5142d61bbb3SSatish Balay     a->mult_work = (Scalar *) PetscMalloc((k+1)*sizeof(Scalar));CHKPTRQ(a->mult_work);
5152d61bbb3SSatish Balay   }
5162d61bbb3SSatish Balay   work = a->mult_work;
5172d61bbb3SSatish Balay   for ( i=0; i<mbs; i++ ) {
5182d61bbb3SSatish Balay     n     = ii[1] - ii[0]; ii++;
5192d61bbb3SSatish Balay     ncols = n*bs;
5202d61bbb3SSatish Balay     workt = work;
5212d61bbb3SSatish Balay     for ( j=0; j<n; j++ ) {
5222d61bbb3SSatish Balay       xb = x + bs*(*idx++);
5232d61bbb3SSatish Balay       for ( k=0; k<bs; k++ ) workt[k] = xb[k];
5242d61bbb3SSatish Balay       workt += bs;
5252d61bbb3SSatish Balay     }
5263f1db9ecSBarry Smith     Kernel_w_gets_Ar_times_v(bs,ncols,work,v,z);
5273f1db9ecSBarry Smith     /* LAgemv_("N",&bs,&ncols,&_DOne,v,&bs,work,&_One,&_DZero,z,&_One); */
5282d61bbb3SSatish Balay     v += n*bs2;
5292d61bbb3SSatish Balay     z += bs;
5302d61bbb3SSatish Balay   }
531e1311b90SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
532e1311b90SBarry Smith   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
5332d61bbb3SSatish Balay   PLogFlops(2*a->nz*bs2 - a->m);
5342d61bbb3SSatish Balay   PetscFunctionReturn(0);
5352d61bbb3SSatish Balay }
5362d61bbb3SSatish Balay 
5372d61bbb3SSatish Balay #undef __FUNC__
5382d61bbb3SSatish Balay #define __FUNC__ "MatMultAdd_SeqBAIJ_1"
5392d61bbb3SSatish Balay int MatMultAdd_SeqBAIJ_1(Mat A,Vec xx,Vec yy,Vec zz)
5402d61bbb3SSatish Balay {
5412d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
5423f1db9ecSBarry Smith   Scalar          *x,*y,*z,sum;
5433f1db9ecSBarry Smith   MatScalar       *v;
544e1311b90SBarry Smith   int             ierr,mbs=a->mbs,i,*idx,*ii,n;
5452d61bbb3SSatish Balay 
5462d61bbb3SSatish Balay   PetscFunctionBegin;
547e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
548e1311b90SBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
5492e8a6d31SBarry Smith   if (zz != yy) {
550e1311b90SBarry Smith     ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
5512e8a6d31SBarry Smith   } else {
5522e8a6d31SBarry Smith     z = y;
5532e8a6d31SBarry Smith   }
5542d61bbb3SSatish Balay 
5552d61bbb3SSatish Balay   idx   = a->j;
5562d61bbb3SSatish Balay   v     = a->a;
5572d61bbb3SSatish Balay   ii    = a->i;
5582d61bbb3SSatish Balay 
5592d61bbb3SSatish Balay   for ( i=0; i<mbs; i++ ) {
5602d61bbb3SSatish Balay     n    = ii[1] - ii[0]; ii++;
5612d61bbb3SSatish Balay     sum  = y[i];
5622d61bbb3SSatish Balay     while (n--) sum += *v++ * x[*idx++];
5632d61bbb3SSatish Balay     z[i] = sum;
5642d61bbb3SSatish Balay   }
565e1311b90SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
566e1311b90SBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
5672e8a6d31SBarry Smith   if (zz != yy) {
568e1311b90SBarry Smith     ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
5692e8a6d31SBarry Smith   }
5702d61bbb3SSatish Balay   PLogFlops(2*a->nz);
5712d61bbb3SSatish Balay   PetscFunctionReturn(0);
5722d61bbb3SSatish Balay }
5732d61bbb3SSatish Balay 
5742d61bbb3SSatish Balay #undef __FUNC__
5752d61bbb3SSatish Balay #define __FUNC__ "MatMultAdd_SeqBAIJ_2"
5762d61bbb3SSatish Balay int MatMultAdd_SeqBAIJ_2(Mat A,Vec xx,Vec yy,Vec zz)
5772d61bbb3SSatish Balay {
5782d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
5793f1db9ecSBarry Smith   Scalar          *x,*y,*z,*xb,sum1,sum2;
580e1311b90SBarry Smith   Scalar          x1,x2;
5813f1db9ecSBarry Smith   MatScalar       *v;
582e1311b90SBarry Smith   int             ierr,mbs=a->mbs,i,*idx,*ii,j,n;
5832d61bbb3SSatish Balay 
5842d61bbb3SSatish Balay   PetscFunctionBegin;
585e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
586e1311b90SBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
5872e8a6d31SBarry Smith   if (zz != yy) {
588e1311b90SBarry Smith     ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
5892e8a6d31SBarry Smith   } else {
5902e8a6d31SBarry Smith     z = y;
5912e8a6d31SBarry Smith   }
5922d61bbb3SSatish Balay 
5932d61bbb3SSatish Balay   idx   = a->j;
5942d61bbb3SSatish Balay   v     = a->a;
5952d61bbb3SSatish Balay   ii    = a->i;
5962d61bbb3SSatish Balay 
5972d61bbb3SSatish Balay   for ( i=0; i<mbs; i++ ) {
5982d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
5992d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1];
6002d61bbb3SSatish Balay     for ( j=0; j<n; j++ ) {
6012d61bbb3SSatish Balay       xb = x + 2*(*idx++); x1 = xb[0]; x2 = xb[1];
6022d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[2]*x2;
6032d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[3]*x2;
6042d61bbb3SSatish Balay       v += 4;
6052d61bbb3SSatish Balay     }
6062d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2;
6072d61bbb3SSatish Balay     z += 2; y += 2;
6082d61bbb3SSatish Balay   }
609e1311b90SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
610e1311b90SBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
6112e8a6d31SBarry Smith   if (zz != yy) {
612e1311b90SBarry Smith     ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
6132e8a6d31SBarry Smith   }
6142d61bbb3SSatish Balay   PLogFlops(4*a->nz);
6152d61bbb3SSatish Balay   PetscFunctionReturn(0);
6162d61bbb3SSatish Balay }
6172d61bbb3SSatish Balay 
6182d61bbb3SSatish Balay #undef __FUNC__
6192d61bbb3SSatish Balay #define __FUNC__ "MatMultAdd_SeqBAIJ_3"
6202d61bbb3SSatish Balay int MatMultAdd_SeqBAIJ_3(Mat A,Vec xx,Vec yy,Vec zz)
6212d61bbb3SSatish Balay {
6222d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
6233f1db9ecSBarry Smith   Scalar          *x,*y,*z,*xb,sum1,sum2,sum3,x1,x2,x3;
6243f1db9ecSBarry Smith   MatScalar       *v;
625e1311b90SBarry Smith   int             ierr,mbs=a->mbs,i,*idx,*ii,j,n;
6262d61bbb3SSatish Balay 
6272d61bbb3SSatish Balay   PetscFunctionBegin;
628e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
629e1311b90SBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
6302e8a6d31SBarry Smith   if (zz != yy) {
631e1311b90SBarry Smith     ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
6322e8a6d31SBarry Smith   } else {
6332e8a6d31SBarry Smith     z = y;
6342e8a6d31SBarry Smith   }
6352d61bbb3SSatish Balay 
6362d61bbb3SSatish Balay   idx   = a->j;
6372d61bbb3SSatish Balay   v     = a->a;
6382d61bbb3SSatish Balay   ii    = a->i;
6392d61bbb3SSatish Balay 
6402d61bbb3SSatish Balay   for ( i=0; i<mbs; i++ ) {
6412d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
6422d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2];
6432d61bbb3SSatish Balay     for ( j=0; j<n; j++ ) {
6442d61bbb3SSatish Balay       xb = x + 3*(*idx++); x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
6452d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[3]*x2 + v[6]*x3;
6462d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[4]*x2 + v[7]*x3;
6472d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[5]*x2 + v[8]*x3;
6482d61bbb3SSatish Balay       v += 9;
6492d61bbb3SSatish Balay     }
6502d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3;
6512d61bbb3SSatish Balay     z += 3; y += 3;
6522d61bbb3SSatish Balay   }
653e1311b90SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
654e1311b90SBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
6552e8a6d31SBarry Smith   if (zz != yy) {
656e1311b90SBarry Smith     ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
6572e8a6d31SBarry Smith   }
6582d61bbb3SSatish Balay   PLogFlops(18*a->nz);
6592d61bbb3SSatish Balay   PetscFunctionReturn(0);
6602d61bbb3SSatish Balay }
6612d61bbb3SSatish Balay 
6622d61bbb3SSatish Balay #undef __FUNC__
6632d61bbb3SSatish Balay #define __FUNC__ "MatMultAdd_SeqBAIJ_4"
6642d61bbb3SSatish Balay int MatMultAdd_SeqBAIJ_4(Mat A,Vec xx,Vec yy,Vec zz)
6652d61bbb3SSatish Balay {
6662d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
6673f1db9ecSBarry Smith   Scalar          *x,*y,*z,*xb,sum1,sum2,sum3,sum4,x1,x2,x3,x4;
6683f1db9ecSBarry Smith   MatScalar       *v;
669e1311b90SBarry Smith   int             ierr,mbs=a->mbs,i,*idx,*ii;
6702d61bbb3SSatish Balay   int             j,n;
6712d61bbb3SSatish Balay 
6722d61bbb3SSatish Balay   PetscFunctionBegin;
673e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
674e1311b90SBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
6752e8a6d31SBarry Smith   if (zz != yy) {
676e1311b90SBarry Smith     ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
6772e8a6d31SBarry Smith   } else {
6782e8a6d31SBarry Smith     z = y;
6792e8a6d31SBarry Smith   }
6802d61bbb3SSatish Balay 
6812d61bbb3SSatish Balay   idx   = a->j;
6822d61bbb3SSatish Balay   v     = a->a;
6832d61bbb3SSatish Balay   ii    = a->i;
6842d61bbb3SSatish Balay 
6852d61bbb3SSatish Balay   for ( i=0; i<mbs; i++ ) {
6862d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
6872d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3];
6882d61bbb3SSatish Balay     for ( j=0; j<n; j++ ) {
6892d61bbb3SSatish Balay       xb = x + 4*(*idx++);
6902d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
6912d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[4]*x2 + v[8]*x3   + v[12]*x4;
6922d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[5]*x2 + v[9]*x3   + v[13]*x4;
6932d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[6]*x2 + v[10]*x3  + v[14]*x4;
6942d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[7]*x2 + v[11]*x3  + v[15]*x4;
6952d61bbb3SSatish Balay       v += 16;
6962d61bbb3SSatish Balay     }
6972d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4;
6982d61bbb3SSatish Balay     z += 4; y += 4;
6992d61bbb3SSatish Balay   }
700e1311b90SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
701e1311b90SBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
7022e8a6d31SBarry Smith   if (zz != yy) {
703e1311b90SBarry Smith     ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
7042e8a6d31SBarry Smith   }
7052d61bbb3SSatish Balay   PLogFlops(32*a->nz);
7062d61bbb3SSatish Balay   PetscFunctionReturn(0);
7072d61bbb3SSatish Balay }
7082d61bbb3SSatish Balay 
7092d61bbb3SSatish Balay #undef __FUNC__
7102d61bbb3SSatish Balay #define __FUNC__ "MatMultAdd_SeqBAIJ_5"
7112d61bbb3SSatish Balay int MatMultAdd_SeqBAIJ_5(Mat A,Vec xx,Vec yy,Vec zz)
7122d61bbb3SSatish Balay {
7132d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
7143f1db9ecSBarry Smith   Scalar          *x,*y,*z,*xb,sum1,sum2,sum3,sum4,sum5,x1,x2,x3,x4,x5;
7153f1db9ecSBarry Smith   MatScalar       *v;
716e1311b90SBarry Smith   int             ierr,mbs=a->mbs,i,*idx,*ii,j,n;
7172d61bbb3SSatish Balay 
7182d61bbb3SSatish Balay   PetscFunctionBegin;
719e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
720e1311b90SBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
7212e8a6d31SBarry Smith   if (zz != yy) {
722e1311b90SBarry Smith     ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
7232e8a6d31SBarry Smith   } else {
7242e8a6d31SBarry Smith     z = y;
7252e8a6d31SBarry Smith   }
7262d61bbb3SSatish Balay 
7272d61bbb3SSatish Balay   idx   = a->j;
7282d61bbb3SSatish Balay   v     = a->a;
7292d61bbb3SSatish Balay   ii    = a->i;
7302d61bbb3SSatish Balay 
7312d61bbb3SSatish Balay   for ( i=0; i<mbs; i++ ) {
7322d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
7332d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4];
7342d61bbb3SSatish Balay     for ( j=0; j<n; j++ ) {
7352d61bbb3SSatish Balay       xb = x + 5*(*idx++);
7362d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4];
7372d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[5]*x2 + v[10]*x3  + v[15]*x4 + v[20]*x5;
7382d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[6]*x2 + v[11]*x3  + v[16]*x4 + v[21]*x5;
7392d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[7]*x2 + v[12]*x3  + v[17]*x4 + v[22]*x5;
7402d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[8]*x2 + v[13]*x3  + v[18]*x4 + v[23]*x5;
7412d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[9]*x2 + v[14]*x3  + v[19]*x4 + v[24]*x5;
7422d61bbb3SSatish Balay       v += 25;
7432d61bbb3SSatish Balay     }
7442d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5;
7452d61bbb3SSatish Balay     z += 5; y += 5;
7462d61bbb3SSatish Balay   }
747e1311b90SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
748e1311b90SBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
7492e8a6d31SBarry Smith   if (zz != yy) {
750e1311b90SBarry Smith     ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
7512e8a6d31SBarry Smith   }
7522d61bbb3SSatish Balay   PLogFlops(50*a->nz);
7532d61bbb3SSatish Balay   PetscFunctionReturn(0);
7542d61bbb3SSatish Balay }
75515091d37SBarry Smith #undef __FUNC__
75615091d37SBarry Smith #define __FUNC__ "MatMultAdd_SeqBAIJ_6"
75715091d37SBarry Smith int MatMultAdd_SeqBAIJ_6(Mat A,Vec xx,Vec yy,Vec zz)
75815091d37SBarry Smith {
75915091d37SBarry Smith   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
76015091d37SBarry Smith   Scalar          *x,*y,*z,*xb,sum1,sum2,sum3,sum4,sum5,sum6;
76115091d37SBarry Smith   Scalar          x1,x2,x3,x4,x5,x6;
76215091d37SBarry Smith   MatScalar       *v;
76315091d37SBarry Smith   int             ierr,mbs=a->mbs,i,*idx,*ii,j,n;
76415091d37SBarry Smith 
76515091d37SBarry Smith   PetscFunctionBegin;
76615091d37SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
76715091d37SBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
76815091d37SBarry Smith   if (zz != yy) {
76915091d37SBarry Smith     ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
77015091d37SBarry Smith   } else {
77115091d37SBarry Smith     z = y;
77215091d37SBarry Smith   }
77315091d37SBarry Smith 
77415091d37SBarry Smith   idx   = a->j;
77515091d37SBarry Smith   v     = a->a;
77615091d37SBarry Smith   ii    = a->i;
77715091d37SBarry Smith 
77815091d37SBarry Smith   for ( i=0; i<mbs; i++ ) {
77915091d37SBarry Smith     n  = ii[1] - ii[0]; ii++;
78015091d37SBarry Smith     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4]; sum6 = y[5];
78115091d37SBarry Smith     for ( j=0; j<n; j++ ) {
7823b95cb0eSSatish Balay       xb = x + 6*(*idx++);
78315091d37SBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5];
78415091d37SBarry Smith       sum1 += v[0]*x1 + v[6]*x2  + v[12]*x3  + v[18]*x4 + v[24]*x5 + v[30]*x6;
78515091d37SBarry Smith       sum2 += v[1]*x1 + v[7]*x2  + v[13]*x3  + v[19]*x4 + v[25]*x5 + v[31]*x6;
78615091d37SBarry Smith       sum3 += v[2]*x1 + v[8]*x2  + v[14]*x3  + v[20]*x4 + v[26]*x5 + v[32]*x6;
78715091d37SBarry Smith       sum4 += v[3]*x1 + v[9]*x2  + v[15]*x3  + v[21]*x4 + v[27]*x5 + v[33]*x6;
78815091d37SBarry Smith       sum5 += v[4]*x1 + v[10]*x2 + v[16]*x3  + v[22]*x4 + v[28]*x5 + v[34]*x6;
78915091d37SBarry Smith       sum6 += v[5]*x1 + v[11]*x2 + v[17]*x3  + v[23]*x4 + v[29]*x5 + v[35]*x6;
79015091d37SBarry Smith       v += 36;
79115091d37SBarry Smith     }
79215091d37SBarry Smith     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6;
79315091d37SBarry Smith     z += 6; y += 6;
79415091d37SBarry Smith   }
79515091d37SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
79615091d37SBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
79715091d37SBarry Smith   if (zz != yy) {
79815091d37SBarry Smith     ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
79915091d37SBarry Smith   }
80015091d37SBarry Smith   PLogFlops(72*a->nz);
80115091d37SBarry Smith   PetscFunctionReturn(0);
80215091d37SBarry Smith }
8032d61bbb3SSatish Balay 
8042d61bbb3SSatish Balay #undef __FUNC__
8052d61bbb3SSatish Balay #define __FUNC__ "MatMultAdd_SeqBAIJ_7"
8062d61bbb3SSatish Balay int MatMultAdd_SeqBAIJ_7(Mat A,Vec xx,Vec yy,Vec zz)
8072d61bbb3SSatish Balay {
8082d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
8093f1db9ecSBarry Smith   Scalar          *x,*y,*z,*xb,sum1,sum2,sum3,sum4,sum5,sum6,sum7;
810e1311b90SBarry Smith   Scalar          x1,x2,x3,x4,x5,x6,x7;
8113f1db9ecSBarry Smith   MatScalar       *v;
812e1311b90SBarry Smith   int             ierr,mbs=a->mbs,i,*idx,*ii,j,n;
8132d61bbb3SSatish Balay 
8142d61bbb3SSatish Balay   PetscFunctionBegin;
815e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
816e1311b90SBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
8172e8a6d31SBarry Smith   if (zz != yy) {
818e1311b90SBarry Smith     ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
8192e8a6d31SBarry Smith   } else {
8202e8a6d31SBarry Smith     z = y;
8212e8a6d31SBarry Smith   }
8222d61bbb3SSatish Balay 
8232d61bbb3SSatish Balay   idx   = a->j;
8242d61bbb3SSatish Balay   v     = a->a;
8252d61bbb3SSatish Balay   ii    = a->i;
8262d61bbb3SSatish Balay 
8272d61bbb3SSatish Balay   for ( i=0; i<mbs; i++ ) {
8282d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
8292d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4]; sum6 = y[5]; sum7 = y[6];
8302d61bbb3SSatish Balay     for ( j=0; j<n; j++ ) {
8312d61bbb3SSatish Balay       xb = x + 7*(*idx++);
8322d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6];
8332d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[7]*x2  + v[14]*x3  + v[21]*x4 + v[28]*x5 + v[35]*x6 + v[42]*x7;
8342d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[8]*x2  + v[15]*x3  + v[22]*x4 + v[29]*x5 + v[36]*x6 + v[43]*x7;
8352d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[9]*x2  + v[16]*x3  + v[23]*x4 + v[30]*x5 + v[37]*x6 + v[44]*x7;
8362d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[10]*x2 + v[17]*x3  + v[24]*x4 + v[31]*x5 + v[38]*x6 + v[45]*x7;
8372d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[11]*x2 + v[18]*x3  + v[25]*x4 + v[32]*x5 + v[39]*x6 + v[46]*x7;
8382d61bbb3SSatish Balay       sum6 += v[5]*x1 + v[12]*x2 + v[19]*x3  + v[26]*x4 + v[33]*x5 + v[40]*x6 + v[47]*x7;
8392d61bbb3SSatish Balay       sum7 += v[6]*x1 + v[13]*x2 + v[20]*x3  + v[27]*x4 + v[34]*x5 + v[41]*x6 + v[48]*x7;
8402d61bbb3SSatish Balay       v += 49;
8412d61bbb3SSatish Balay     }
8422d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
8432d61bbb3SSatish Balay     z += 7; y += 7;
8442d61bbb3SSatish Balay   }
845e1311b90SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
846e1311b90SBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
8472e8a6d31SBarry Smith   if (zz != yy) {
848e1311b90SBarry Smith     ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
8492e8a6d31SBarry Smith   }
8502d61bbb3SSatish Balay   PLogFlops(98*a->nz);
8512d61bbb3SSatish Balay   PetscFunctionReturn(0);
8522d61bbb3SSatish Balay }
853218c64b6SSatish Balay 
8542d61bbb3SSatish Balay #undef __FUNC__
8552d61bbb3SSatish Balay #define __FUNC__ "MatMultAdd_SeqBAIJ_N"
8562d61bbb3SSatish Balay int MatMultAdd_SeqBAIJ_N(Mat A,Vec xx,Vec yy,Vec zz)
8572d61bbb3SSatish Balay {
8582d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ *) A->data;
8593f1db9ecSBarry Smith   Scalar         *x,*z,*xb,*work,*workt;
8603f1db9ecSBarry Smith   MatScalar      *v;
8612d61bbb3SSatish Balay   int            mbs=a->mbs,i,*idx,*ii,bs=a->bs,j,n,bs2=a->bs2,ierr;
8623f1db9ecSBarry Smith   int            ncols,k;
863218c64b6SSatish Balay 
8642d61bbb3SSatish Balay   PetscFunctionBegin;
8652d61bbb3SSatish Balay   if ( xx != yy) { ierr = VecCopy(yy,zz);CHKERRQ(ierr); }
8662d61bbb3SSatish Balay 
867e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
868e1311b90SBarry Smith   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
8692d61bbb3SSatish Balay 
8702d61bbb3SSatish Balay   idx   = a->j;
8712d61bbb3SSatish Balay   v     = a->a;
8722d61bbb3SSatish Balay   ii    = a->i;
8732d61bbb3SSatish Balay 
8742d61bbb3SSatish Balay 
8752d61bbb3SSatish Balay   if (!a->mult_work) {
8762d61bbb3SSatish Balay     k = PetscMax(a->m,a->n);
877544c85edSSatish Balay     a->mult_work = (Scalar *) PetscMalloc((k+1)*sizeof(Scalar));CHKPTRQ(a->mult_work);
8782d61bbb3SSatish Balay   }
8792d61bbb3SSatish Balay   work = a->mult_work;
8802d61bbb3SSatish Balay   for ( i=0; i<mbs; i++ ) {
8812d61bbb3SSatish Balay     n     = ii[1] - ii[0]; ii++;
8822d61bbb3SSatish Balay     ncols = n*bs;
8832d61bbb3SSatish Balay     workt = work;
8842d61bbb3SSatish Balay     for ( j=0; j<n; j++ ) {
8852d61bbb3SSatish Balay       xb = x + bs*(*idx++);
8862d61bbb3SSatish Balay       for ( k=0; k<bs; k++ ) workt[k] = xb[k];
8872d61bbb3SSatish Balay       workt += bs;
8882d61bbb3SSatish Balay     }
8893f1db9ecSBarry Smith     Kernel_w_gets_w_plus_Ar_times_v(bs,ncols,work,v,z);
8903f1db9ecSBarry Smith     /* LAgemv_("N",&bs,&ncols,&_DOne,v,&bs,work,&_One,&_DOne,z,&_One); */
8912d61bbb3SSatish Balay     v += n*bs2;
8922d61bbb3SSatish Balay     z += bs;
8932d61bbb3SSatish Balay   }
894e1311b90SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
895e1311b90SBarry Smith   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
8962d61bbb3SSatish Balay   PLogFlops(2*a->nz*bs2 );
8972d61bbb3SSatish Balay   PetscFunctionReturn(0);
8982d61bbb3SSatish Balay }
8992d61bbb3SSatish Balay 
9002d61bbb3SSatish Balay #undef __FUNC__
9012d61bbb3SSatish Balay #define __FUNC__ "MatMultTrans_SeqBAIJ"
9022d61bbb3SSatish Balay int MatMultTrans_SeqBAIJ(Mat A,Vec xx,Vec zz)
9032d61bbb3SSatish Balay {
9042d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
9052d61bbb3SSatish Balay   Scalar          *xg,*zg,*zb;
9063f1db9ecSBarry Smith   Scalar          *x,*z,*xb,x1,x2,x3,x4,x5;
9073f1db9ecSBarry Smith   MatScalar       *v;
9082d61bbb3SSatish Balay   int             mbs=a->mbs,i,*idx,*ii,*ai=a->i,rval,N=a->n;
9092d61bbb3SSatish Balay   int             bs=a->bs,j,n,bs2=a->bs2,*ib,ierr;
9102d61bbb3SSatish Balay 
9112d61bbb3SSatish Balay 
9122d61bbb3SSatish Balay   PetscFunctionBegin;
913e1311b90SBarry Smith   ierr = VecGetArray(xx,&xg);CHKERRQ(ierr); x = xg;
914e1311b90SBarry Smith   ierr = VecGetArray(zz,&zg);CHKERRQ(ierr); z = zg;
915549d3d68SSatish Balay   ierr = PetscMemzero(z,N*sizeof(Scalar));CHKERRQ(ierr);
9162d61bbb3SSatish Balay 
9172d61bbb3SSatish Balay   idx   = a->j;
9182d61bbb3SSatish Balay   v     = a->a;
9192d61bbb3SSatish Balay   ii    = a->i;
9202d61bbb3SSatish Balay 
9212d61bbb3SSatish Balay   switch (bs) {
9222d61bbb3SSatish Balay   case 1:
9232d61bbb3SSatish Balay     for ( i=0; i<mbs; i++ ) {
9242d61bbb3SSatish Balay       n  = ii[1] - ii[0]; ii++;
9252d61bbb3SSatish Balay       xb = x + i; x1 = xb[0];
9262d61bbb3SSatish Balay       ib = idx + ai[i];
9272d61bbb3SSatish Balay       for ( j=0; j<n; j++ ) {
9282d61bbb3SSatish Balay         rval    = ib[j];
9292d61bbb3SSatish Balay         z[rval] += *v++ * x1;
9302d61bbb3SSatish Balay       }
9312d61bbb3SSatish Balay     }
9322d61bbb3SSatish Balay     break;
9332d61bbb3SSatish Balay   case 2:
9342d61bbb3SSatish Balay     for ( i=0; i<mbs; i++ ) {
9352d61bbb3SSatish Balay       n  = ii[1] - ii[0]; ii++;
9362d61bbb3SSatish Balay       xb = x + 2*i; x1 = xb[0]; x2 = xb[1];
9372d61bbb3SSatish Balay       ib = idx + ai[i];
9382d61bbb3SSatish Balay       for ( j=0; j<n; j++ ) {
9392d61bbb3SSatish Balay         rval      = ib[j]*2;
9402d61bbb3SSatish Balay         z[rval++] += v[0]*x1 + v[1]*x2;
9412d61bbb3SSatish Balay         z[rval++] += v[2]*x1 + v[3]*x2;
9422d61bbb3SSatish Balay         v += 4;
9432d61bbb3SSatish Balay       }
9442d61bbb3SSatish Balay     }
9452d61bbb3SSatish Balay     break;
9462d61bbb3SSatish Balay   case 3:
9472d61bbb3SSatish Balay     for ( i=0; i<mbs; i++ ) {
9482d61bbb3SSatish Balay       n  = ii[1] - ii[0]; ii++;
9492d61bbb3SSatish Balay       xb = x + 3*i; x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
9502d61bbb3SSatish Balay       ib = idx + ai[i];
9512d61bbb3SSatish Balay       for ( j=0; j<n; j++ ) {
9522d61bbb3SSatish Balay         rval      = ib[j]*3;
9532d61bbb3SSatish Balay         z[rval++] += v[0]*x1 + v[1]*x2 + v[2]*x3;
9542d61bbb3SSatish Balay         z[rval++] += v[3]*x1 + v[4]*x2 + v[5]*x3;
9552d61bbb3SSatish Balay         z[rval++] += v[6]*x1 + v[7]*x2 + v[8]*x3;
9562d61bbb3SSatish Balay         v += 9;
9572d61bbb3SSatish Balay       }
9582d61bbb3SSatish Balay     }
9592d61bbb3SSatish Balay     break;
9602d61bbb3SSatish Balay   case 4:
9612d61bbb3SSatish Balay     for ( i=0; i<mbs; i++ ) {
9622d61bbb3SSatish Balay       n  = ii[1] - ii[0]; ii++;
9632d61bbb3SSatish Balay       xb = x + 4*i; x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
9642d61bbb3SSatish Balay       ib = idx + ai[i];
9652d61bbb3SSatish Balay       for ( j=0; j<n; j++ ) {
9662d61bbb3SSatish Balay         rval      = ib[j]*4;
9672d61bbb3SSatish Balay         z[rval++] +=  v[0]*x1 +  v[1]*x2 +  v[2]*x3 +  v[3]*x4;
9682d61bbb3SSatish Balay         z[rval++] +=  v[4]*x1 +  v[5]*x2 +  v[6]*x3 +  v[7]*x4;
9692d61bbb3SSatish Balay         z[rval++] +=  v[8]*x1 +  v[9]*x2 + v[10]*x3 + v[11]*x4;
9702d61bbb3SSatish Balay         z[rval++] += v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4;
9712d61bbb3SSatish Balay         v += 16;
9722d61bbb3SSatish Balay       }
9732d61bbb3SSatish Balay     }
9742d61bbb3SSatish Balay     break;
9752d61bbb3SSatish Balay   case 5:
9762d61bbb3SSatish Balay     for ( i=0; i<mbs; i++ ) {
9772d61bbb3SSatish Balay       n  = ii[1] - ii[0]; ii++;
9782d61bbb3SSatish Balay       xb = x + 5*i; x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
9792d61bbb3SSatish Balay       x4 = xb[3];   x5 = xb[4];
9802d61bbb3SSatish Balay       ib = idx + ai[i];
9812d61bbb3SSatish Balay       for ( j=0; j<n; j++ ) {
9822d61bbb3SSatish Balay         rval      = ib[j]*5;
9832d61bbb3SSatish Balay         z[rval++] +=  v[0]*x1 +  v[1]*x2 +  v[2]*x3 +  v[3]*x4 +  v[4]*x5;
9842d61bbb3SSatish Balay         z[rval++] +=  v[5]*x1 +  v[6]*x2 +  v[7]*x3 +  v[8]*x4 +  v[9]*x5;
9852d61bbb3SSatish Balay         z[rval++] += v[10]*x1 + v[11]*x2 + v[12]*x3 + v[13]*x4 + v[14]*x5;
9862d61bbb3SSatish Balay         z[rval++] += v[15]*x1 + v[16]*x2 + v[17]*x3 + v[18]*x4 + v[19]*x5;
9872d61bbb3SSatish Balay         z[rval++] += v[20]*x1 + v[21]*x2 + v[22]*x3 + v[23]*x4 + v[24]*x5;
9882d61bbb3SSatish Balay         v += 25;
9892d61bbb3SSatish Balay       }
9902d61bbb3SSatish Balay     }
9912d61bbb3SSatish Balay     break;
9922d61bbb3SSatish Balay       /* block sizes larger then 5 by 5 are handled by BLAS */
9932d61bbb3SSatish Balay     default: {
9943f1db9ecSBarry Smith       int       ncols,k;
9953f1db9ecSBarry Smith       Scalar    *work,*workt;
9963f1db9ecSBarry Smith 
9972d61bbb3SSatish Balay       if (!a->mult_work) {
9982d61bbb3SSatish Balay         k = PetscMax(a->m,a->n);
999544c85edSSatish Balay         a->mult_work = (Scalar *) PetscMalloc((k+1)*sizeof(Scalar));CHKPTRQ(a->mult_work);
10002d61bbb3SSatish Balay       }
10012d61bbb3SSatish Balay       work = a->mult_work;
10022d61bbb3SSatish Balay       for ( i=0; i<mbs; i++ ) {
10032d61bbb3SSatish Balay         n     = ii[1] - ii[0]; ii++;
10042d61bbb3SSatish Balay         ncols = n*bs;
1005549d3d68SSatish Balay         ierr  = PetscMemzero(work,ncols*sizeof(Scalar));CHKERRQ(ierr);
1006d824769bSBarry Smith         Kernel_w_gets_w_plus_trans_Ar_times_v(bs,ncols,x,v,work);
10073f1db9ecSBarry Smith         /* LAgemv_("T",&bs,&ncols,&_DOne,v,&bs,x,&_One,&_DOne,work,&_One); */
10082d61bbb3SSatish Balay         v += n*bs2;
10092d61bbb3SSatish Balay         x += bs;
10102d61bbb3SSatish Balay         workt = work;
10112d61bbb3SSatish Balay         for ( j=0; j<n; j++ ) {
10122d61bbb3SSatish Balay           zb = z + bs*(*idx++);
10132d61bbb3SSatish Balay           for ( k=0; k<bs; k++ ) zb[k] += workt[k] ;
10142d61bbb3SSatish Balay           workt += bs;
10152d61bbb3SSatish Balay         }
10162d61bbb3SSatish Balay       }
10172d61bbb3SSatish Balay     }
10182d61bbb3SSatish Balay   }
10192d61bbb3SSatish Balay   ierr = VecRestoreArray(xx,&xg);CHKERRQ(ierr);
10202d61bbb3SSatish Balay   ierr = VecRestoreArray(zz,&zg);CHKERRQ(ierr);
10212d61bbb3SSatish Balay   PLogFlops(2*a->nz*a->bs2 - a->n);
10222d61bbb3SSatish Balay   PetscFunctionReturn(0);
10232d61bbb3SSatish Balay }
10242d61bbb3SSatish Balay 
10252d61bbb3SSatish Balay #undef __FUNC__
10262d61bbb3SSatish Balay #define __FUNC__ "MatMultTransAdd_SeqBAIJ"
10272d61bbb3SSatish Balay int MatMultTransAdd_SeqBAIJ(Mat A,Vec xx,Vec yy,Vec zz)
10282d61bbb3SSatish Balay 
10292d61bbb3SSatish Balay {
10302d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
10312d61bbb3SSatish Balay   Scalar          *xg,*zg,*zb;
10323f1db9ecSBarry Smith   Scalar          *x,*z,*xb,x1,x2,x3,x4,x5;
10333f1db9ecSBarry Smith   MatScalar       *v;
10342d61bbb3SSatish Balay   int             mbs=a->mbs,i,*idx,*ii,*ai=a->i,rval,N=a->n;
10352d61bbb3SSatish Balay   int             bs=a->bs,j,n,bs2=a->bs2,*ib,ierr;
10362d61bbb3SSatish Balay 
10372d61bbb3SSatish Balay   PetscFunctionBegin;
1038e1311b90SBarry Smith   ierr = VecGetArray(xx,&xg);CHKERRQ(ierr); x = xg;
1039e1311b90SBarry Smith   ierr = VecGetArray(zz,&zg);CHKERRQ(ierr); z = zg;
10402d61bbb3SSatish Balay 
10412d61bbb3SSatish Balay   if ( yy != zz ) { ierr = VecCopy(yy,zz);CHKERRQ(ierr); }
1042549d3d68SSatish Balay   else {
1043549d3d68SSatish Balay     ierr = PetscMemzero(z,N*sizeof(Scalar));CHKERRQ(ierr);
1044549d3d68SSatish Balay   }
10452d61bbb3SSatish Balay 
10462d61bbb3SSatish Balay   idx   = a->j;
10472d61bbb3SSatish Balay   v     = a->a;
10482d61bbb3SSatish Balay   ii    = a->i;
10492d61bbb3SSatish Balay 
10502d61bbb3SSatish Balay   switch (bs) {
10512d61bbb3SSatish Balay   case 1:
10522d61bbb3SSatish Balay     for ( i=0; i<mbs; i++ ) {
10532d61bbb3SSatish Balay       n  = ii[1] - ii[0]; ii++;
10542d61bbb3SSatish Balay       xb = x + i; x1 = xb[0];
10552d61bbb3SSatish Balay       ib = idx + ai[i];
10562d61bbb3SSatish Balay       for ( j=0; j<n; j++ ) {
10572d61bbb3SSatish Balay         rval    = ib[j];
10582d61bbb3SSatish Balay         z[rval] += *v++ * x1;
10592d61bbb3SSatish Balay       }
10602d61bbb3SSatish Balay     }
10612d61bbb3SSatish Balay     break;
10622d61bbb3SSatish Balay   case 2:
10632d61bbb3SSatish Balay     for ( i=0; i<mbs; i++ ) {
10642d61bbb3SSatish Balay       n  = ii[1] - ii[0]; ii++;
10652d61bbb3SSatish Balay       xb = x + 2*i; x1 = xb[0]; x2 = xb[1];
10662d61bbb3SSatish Balay       ib = idx + ai[i];
10672d61bbb3SSatish Balay       for ( j=0; j<n; j++ ) {
10682d61bbb3SSatish Balay         rval      = ib[j]*2;
10692d61bbb3SSatish Balay         z[rval++] += v[0]*x1 + v[1]*x2;
10702d61bbb3SSatish Balay         z[rval++] += v[2]*x1 + v[3]*x2;
10712d61bbb3SSatish Balay         v += 4;
10722d61bbb3SSatish Balay       }
10732d61bbb3SSatish Balay     }
10742d61bbb3SSatish Balay     break;
10752d61bbb3SSatish Balay   case 3:
10762d61bbb3SSatish Balay     for ( i=0; i<mbs; i++ ) {
10772d61bbb3SSatish Balay       n  = ii[1] - ii[0]; ii++;
10782d61bbb3SSatish Balay       xb = x + 3*i; x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
10792d61bbb3SSatish Balay       ib = idx + ai[i];
10802d61bbb3SSatish Balay       for ( j=0; j<n; j++ ) {
10812d61bbb3SSatish Balay         rval      = ib[j]*3;
10822d61bbb3SSatish Balay         z[rval++] += v[0]*x1 + v[1]*x2 + v[2]*x3;
10832d61bbb3SSatish Balay         z[rval++] += v[3]*x1 + v[4]*x2 + v[5]*x3;
10842d61bbb3SSatish Balay         z[rval++] += v[6]*x1 + v[7]*x2 + v[8]*x3;
10852d61bbb3SSatish Balay         v += 9;
10862d61bbb3SSatish Balay       }
10872d61bbb3SSatish Balay     }
10882d61bbb3SSatish Balay     break;
10892d61bbb3SSatish Balay   case 4:
10902d61bbb3SSatish Balay     for ( i=0; i<mbs; i++ ) {
10912d61bbb3SSatish Balay       n  = ii[1] - ii[0]; ii++;
10922d61bbb3SSatish Balay       xb = x + 4*i; x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
10932d61bbb3SSatish Balay       ib = idx + ai[i];
10942d61bbb3SSatish Balay       for ( j=0; j<n; j++ ) {
10952d61bbb3SSatish Balay         rval      = ib[j]*4;
10962d61bbb3SSatish Balay         z[rval++] +=  v[0]*x1 +  v[1]*x2 +  v[2]*x3 +  v[3]*x4;
10972d61bbb3SSatish Balay         z[rval++] +=  v[4]*x1 +  v[5]*x2 +  v[6]*x3 +  v[7]*x4;
10982d61bbb3SSatish Balay         z[rval++] +=  v[8]*x1 +  v[9]*x2 + v[10]*x3 + v[11]*x4;
10992d61bbb3SSatish Balay         z[rval++] += v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4;
11002d61bbb3SSatish Balay         v += 16;
11012d61bbb3SSatish Balay       }
11022d61bbb3SSatish Balay     }
11032d61bbb3SSatish Balay     break;
11042d61bbb3SSatish Balay   case 5:
11052d61bbb3SSatish Balay     for ( i=0; i<mbs; i++ ) {
11062d61bbb3SSatish Balay       n  = ii[1] - ii[0]; ii++;
11072d61bbb3SSatish Balay       xb = x + 5*i; x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
11082d61bbb3SSatish Balay       x4 = xb[3];   x5 = xb[4];
11092d61bbb3SSatish Balay       ib = idx + ai[i];
11102d61bbb3SSatish Balay       for ( j=0; j<n; j++ ) {
11112d61bbb3SSatish Balay         rval      = ib[j]*5;
11122d61bbb3SSatish Balay         z[rval++] +=  v[0]*x1 +  v[1]*x2 +  v[2]*x3 +  v[3]*x4 +  v[4]*x5;
11132d61bbb3SSatish Balay         z[rval++] +=  v[5]*x1 +  v[6]*x2 +  v[7]*x3 +  v[8]*x4 +  v[9]*x5;
11142d61bbb3SSatish Balay         z[rval++] += v[10]*x1 + v[11]*x2 + v[12]*x3 + v[13]*x4 + v[14]*x5;
11152d61bbb3SSatish Balay         z[rval++] += v[15]*x1 + v[16]*x2 + v[17]*x3 + v[18]*x4 + v[19]*x5;
11162d61bbb3SSatish Balay         z[rval++] += v[20]*x1 + v[21]*x2 + v[22]*x3 + v[23]*x4 + v[24]*x5;
11172d61bbb3SSatish Balay         v += 25;
11182d61bbb3SSatish Balay       }
11192d61bbb3SSatish Balay     }
11202d61bbb3SSatish Balay     break;
11212d61bbb3SSatish Balay       /* block sizes larger then 5 by 5 are handled by BLAS */
11222d61bbb3SSatish Balay     default: {
11233f1db9ecSBarry Smith       int       ncols,k;
11243f1db9ecSBarry Smith       Scalar    *work,*workt;
11253f1db9ecSBarry Smith 
11262d61bbb3SSatish Balay       if (!a->mult_work) {
11272d61bbb3SSatish Balay         k = PetscMax(a->m,a->n);
1128544c85edSSatish Balay         a->mult_work = (Scalar *) PetscMalloc((k+1)*sizeof(Scalar));CHKPTRQ(a->mult_work);
11292d61bbb3SSatish Balay       }
11302d61bbb3SSatish Balay       work = a->mult_work;
11312d61bbb3SSatish Balay       for ( i=0; i<mbs; i++ ) {
11322d61bbb3SSatish Balay         n     = ii[1] - ii[0]; ii++;
11332d61bbb3SSatish Balay         ncols = n*bs;
1134549d3d68SSatish Balay         ierr  = PetscMemzero(work,ncols*sizeof(Scalar));CHKERRQ(ierr);
11353f1db9ecSBarry Smith         Kernel_w_gets_w_plus_trans_Ar_times_v(bs,ncols,x,v,work);
11363f1db9ecSBarry Smith         /* LAgemv_("T",&bs,&ncols,&_DOne,v,&bs,x,&_One,&_DOne,work,&_One); */
11372d61bbb3SSatish Balay         v += n*bs2;
11382d61bbb3SSatish Balay         x += bs;
11392d61bbb3SSatish Balay         workt = work;
11402d61bbb3SSatish Balay         for ( j=0; j<n; j++ ) {
11412d61bbb3SSatish Balay           zb = z + bs*(*idx++);
11422d61bbb3SSatish Balay           for ( k=0; k<bs; k++ ) zb[k] += workt[k] ;
11432d61bbb3SSatish Balay           workt += bs;
11442d61bbb3SSatish Balay         }
11452d61bbb3SSatish Balay       }
11462d61bbb3SSatish Balay     }
11472d61bbb3SSatish Balay   }
11482d61bbb3SSatish Balay   ierr = VecRestoreArray(xx,&xg);CHKERRQ(ierr);
11492d61bbb3SSatish Balay   ierr = VecRestoreArray(zz,&zg);CHKERRQ(ierr);
11502d61bbb3SSatish Balay   PLogFlops(2*a->nz*a->bs2);
11512d61bbb3SSatish Balay   PetscFunctionReturn(0);
11522d61bbb3SSatish Balay }
11532d61bbb3SSatish Balay 
11542d61bbb3SSatish Balay #undef __FUNC__
11552d61bbb3SSatish Balay #define __FUNC__ "MatScale_SeqBAIJ"
11562d61bbb3SSatish Balay int MatScale_SeqBAIJ(Scalar *alpha,Mat inA)
11572d61bbb3SSatish Balay {
11582d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) inA->data;
11592d61bbb3SSatish Balay   int         one = 1, totalnz = a->bs2*a->nz;
11602d61bbb3SSatish Balay 
11612d61bbb3SSatish Balay   PetscFunctionBegin;
11622d61bbb3SSatish Balay   BLscal_( &totalnz, alpha, a->a, &one );
11632d61bbb3SSatish Balay   PLogFlops(totalnz);
11642d61bbb3SSatish Balay   PetscFunctionReturn(0);
11652d61bbb3SSatish Balay }
11662d61bbb3SSatish Balay 
11672d61bbb3SSatish Balay #undef __FUNC__
11682d61bbb3SSatish Balay #define __FUNC__ "MatNorm_SeqBAIJ"
11692d61bbb3SSatish Balay int MatNorm_SeqBAIJ(Mat A,NormType type,double *norm)
11702d61bbb3SSatish Balay {
11712d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
11723f1db9ecSBarry Smith   MatScalar   *v = a->a;
11732d61bbb3SSatish Balay   double      sum = 0.0;
1174253ffaf7SBarry Smith   int         i,j,k,bs = a->bs, nz=a->nz,bs2=a->bs2;
11752d61bbb3SSatish Balay 
11762d61bbb3SSatish Balay   PetscFunctionBegin;
11772d61bbb3SSatish Balay   if (type == NORM_FROBENIUS) {
11782d61bbb3SSatish Balay     for (i=0; i< bs2*nz; i++ ) {
1179aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
1180e20fef11SSatish Balay       sum += PetscReal(PetscConj(*v)*(*v)); v++;
11812d61bbb3SSatish Balay #else
11822d61bbb3SSatish Balay       sum += (*v)*(*v); v++;
11832d61bbb3SSatish Balay #endif
11842d61bbb3SSatish Balay     }
11852d61bbb3SSatish Balay     *norm = sqrt(sum);
1186596552b5SBarry Smith   }  else if (type == NORM_INFINITY) { /* maximum row sum */
1187596552b5SBarry Smith     *norm = 0.0;
1188596552b5SBarry Smith     for ( k=0; k<bs; k++ ) {
1189596552b5SBarry Smith       for ( j=0; j<a->m; j++ ) {
1190596552b5SBarry Smith         v = a->a + bs2*a->i[j] + k;
1191596552b5SBarry Smith         sum = 0.0;
1192596552b5SBarry Smith         for ( i=0; i<a->i[j+1]-a->i[j]; i++ ) {
1193596552b5SBarry Smith           sum += PetscAbsScalar(*v);
1194596552b5SBarry Smith           v   += bs;
11952d61bbb3SSatish Balay         }
1196596552b5SBarry Smith         if (sum > *norm) *norm = sum;
1197596552b5SBarry Smith       }
1198596552b5SBarry Smith     }
1199596552b5SBarry Smith   } else {
12002d61bbb3SSatish Balay     SETERRQ(PETSC_ERR_SUP,0,"No support for this norm yet");
12012d61bbb3SSatish Balay   }
12022d61bbb3SSatish Balay   PetscFunctionReturn(0);
12032d61bbb3SSatish Balay }
12042d61bbb3SSatish Balay 
12052d61bbb3SSatish Balay 
12062d61bbb3SSatish Balay #undef __FUNC__
12072d61bbb3SSatish Balay #define __FUNC__ "MatEqual_SeqBAIJ"
12082d61bbb3SSatish Balay int MatEqual_SeqBAIJ(Mat A,Mat B, PetscTruth* flg)
12092d61bbb3SSatish Balay {
12102d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)B->data;
12112d61bbb3SSatish Balay 
12122d61bbb3SSatish Balay   PetscFunctionBegin;
12132d61bbb3SSatish Balay   if (B->type !=MATSEQBAIJ) SETERRQ(PETSC_ERR_ARG_INCOMP,0,"Matrices must be same type");
12142d61bbb3SSatish Balay 
12152d61bbb3SSatish Balay   /* If the  matrix/block dimensions are not equal, or no of nonzeros or shift */
12162d61bbb3SSatish Balay   if ((a->m != b->m) || (a->n !=b->n) || (a->bs != b->bs)|| (a->nz != b->nz)) {
12172d61bbb3SSatish Balay     *flg = PETSC_FALSE; PetscFunctionReturn(0);
12182d61bbb3SSatish Balay   }
12192d61bbb3SSatish Balay 
12202d61bbb3SSatish Balay   /* if the a->i are the same */
12212d61bbb3SSatish Balay   if (PetscMemcmp(a->i,b->i, (a->mbs+1)*sizeof(int))) {
12222d61bbb3SSatish Balay     *flg = PETSC_FALSE; PetscFunctionReturn(0);
12232d61bbb3SSatish Balay   }
12242d61bbb3SSatish Balay 
12252d61bbb3SSatish Balay   /* if a->j are the same */
12262d61bbb3SSatish Balay   if (PetscMemcmp(a->j,b->j,(a->nz)*sizeof(int))) {
12272d61bbb3SSatish Balay     *flg = PETSC_FALSE; PetscFunctionReturn(0);
12282d61bbb3SSatish Balay   }
12292d61bbb3SSatish Balay 
12302d61bbb3SSatish Balay   /* if a->a are the same */
12312d61bbb3SSatish Balay   if (PetscMemcmp(a->a, b->a,(a->nz)*(a->bs)*(a->bs)*sizeof(Scalar))) {
12322d61bbb3SSatish Balay     *flg = PETSC_FALSE; PetscFunctionReturn(0);
12332d61bbb3SSatish Balay   }
12342d61bbb3SSatish Balay   *flg = PETSC_TRUE;
12352d61bbb3SSatish Balay   PetscFunctionReturn(0);
12362d61bbb3SSatish Balay 
12372d61bbb3SSatish Balay }
12382d61bbb3SSatish Balay 
12392d61bbb3SSatish Balay #undef __FUNC__
12402d61bbb3SSatish Balay #define __FUNC__ "MatGetDiagonal_SeqBAIJ"
12412d61bbb3SSatish Balay int MatGetDiagonal_SeqBAIJ(Mat A,Vec v)
12422d61bbb3SSatish Balay {
12432d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
1244e1311b90SBarry Smith   int         ierr,i,j,k,n,row,bs,*ai,*aj,ambs,bs2;
12453f1db9ecSBarry Smith   Scalar      *x, zero = 0.0;
12463f1db9ecSBarry Smith   MatScalar   *aa,*aa_j;
12472d61bbb3SSatish Balay 
12482d61bbb3SSatish Balay   PetscFunctionBegin;
12492d61bbb3SSatish Balay   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix");
12502d61bbb3SSatish Balay   bs   = a->bs;
12512d61bbb3SSatish Balay   aa   = a->a;
12522d61bbb3SSatish Balay   ai   = a->i;
12532d61bbb3SSatish Balay   aj   = a->j;
12542d61bbb3SSatish Balay   ambs = a->mbs;
12552d61bbb3SSatish Balay   bs2  = a->bs2;
12562d61bbb3SSatish Balay 
1257e1311b90SBarry Smith   ierr = VecSet(&zero,v);CHKERRQ(ierr);
1258e1311b90SBarry Smith   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1259e1311b90SBarry Smith   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
12602d61bbb3SSatish Balay   if (n != a->m) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Nonconforming matrix and vector");
12612d61bbb3SSatish Balay   for ( i=0; i<ambs; i++ ) {
12622d61bbb3SSatish Balay     for ( j=ai[i]; j<ai[i+1]; j++ ) {
12632d61bbb3SSatish Balay       if (aj[j] == i) {
12642d61bbb3SSatish Balay         row  = i*bs;
12652d61bbb3SSatish Balay         aa_j = aa+j*bs2;
12662d61bbb3SSatish Balay         for (k=0; k<bs2; k+=(bs+1),row++) x[row] = aa_j[k];
12672d61bbb3SSatish Balay         break;
12682d61bbb3SSatish Balay       }
12692d61bbb3SSatish Balay     }
12702d61bbb3SSatish Balay   }
127115091d37SBarry Smith   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
12722d61bbb3SSatish Balay   PetscFunctionReturn(0);
12732d61bbb3SSatish Balay }
12742d61bbb3SSatish Balay 
12752d61bbb3SSatish Balay #undef __FUNC__
12762d61bbb3SSatish Balay #define __FUNC__ "MatDiagonalScale_SeqBAIJ"
12772d61bbb3SSatish Balay int MatDiagonalScale_SeqBAIJ(Mat A,Vec ll,Vec rr)
12782d61bbb3SSatish Balay {
12792d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
12803f1db9ecSBarry Smith   Scalar      *l,*r,x,*li,*ri;
12813f1db9ecSBarry Smith   MatScalar   *aa,*v;
1282e1311b90SBarry Smith   int         ierr,i,j,k,lm,rn,M,m,n,*ai,*aj,mbs,tmp,bs,bs2;
12832d61bbb3SSatish Balay 
12842d61bbb3SSatish Balay   PetscFunctionBegin;
12852d61bbb3SSatish Balay   ai  = a->i;
12862d61bbb3SSatish Balay   aj  = a->j;
12872d61bbb3SSatish Balay   aa  = a->a;
12882d61bbb3SSatish Balay   m   = a->m;
12892d61bbb3SSatish Balay   n   = a->n;
12902d61bbb3SSatish Balay   bs  = a->bs;
12912d61bbb3SSatish Balay   mbs = a->mbs;
12922d61bbb3SSatish Balay   bs2 = a->bs2;
12932d61bbb3SSatish Balay   if (ll) {
1294e1311b90SBarry Smith     ierr = VecGetArray(ll,&l);CHKERRQ(ierr);
1295e1311b90SBarry Smith     ierr = VecGetLocalSize(ll,&lm);CHKERRQ(ierr);
12962d61bbb3SSatish Balay     if (lm != m) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Left scaling vector wrong length");
12972d61bbb3SSatish Balay     for ( i=0; i<mbs; i++ ) { /* for each block row */
12982d61bbb3SSatish Balay       M  = ai[i+1] - ai[i];
12992d61bbb3SSatish Balay       li = l + i*bs;
13002d61bbb3SSatish Balay       v  = aa + bs2*ai[i];
13012d61bbb3SSatish Balay       for ( j=0; j<M; j++ ) { /* for each block */
13022d61bbb3SSatish Balay         for ( k=0; k<bs2; k++ ) {
13032d61bbb3SSatish Balay           (*v++) *= li[k%bs];
13042d61bbb3SSatish Balay         }
13052d61bbb3SSatish Balay       }
13062d61bbb3SSatish Balay     }
130760d550acSSatish Balay     ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr);
130896896c26SBarry Smith     PLogFlops(a->nz);
13092d61bbb3SSatish Balay   }
13102d61bbb3SSatish Balay 
13112d61bbb3SSatish Balay   if (rr) {
1312e1311b90SBarry Smith     ierr = VecGetArray(rr,&r);CHKERRQ(ierr);
1313e1311b90SBarry Smith     ierr = VecGetLocalSize(rr,&rn);CHKERRQ(ierr);
13142d61bbb3SSatish Balay     if (rn != n) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Right scaling vector wrong length");
13152d61bbb3SSatish Balay     for ( i=0; i<mbs; i++ ) { /* for each block row */
13162d61bbb3SSatish Balay       M  = ai[i+1] - ai[i];
13172d61bbb3SSatish Balay       v  = aa + bs2*ai[i];
13182d61bbb3SSatish Balay       for ( j=0; j<M; j++ ) { /* for each block */
13192d61bbb3SSatish Balay         ri = r + bs*aj[ai[i]+j];
13202d61bbb3SSatish Balay         for ( k=0; k<bs; k++ ) {
13212d61bbb3SSatish Balay           x = ri[k];
13222d61bbb3SSatish Balay           for ( tmp=0; tmp<bs; tmp++ ) (*v++) *= x;
13232d61bbb3SSatish Balay         }
13242d61bbb3SSatish Balay       }
13252d61bbb3SSatish Balay     }
132660d550acSSatish Balay     ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr);
132796896c26SBarry Smith     PLogFlops(a->nz);
13282d61bbb3SSatish Balay   }
13292d61bbb3SSatish Balay   PetscFunctionReturn(0);
13302d61bbb3SSatish Balay }
13312d61bbb3SSatish Balay 
13322d61bbb3SSatish Balay 
13332d61bbb3SSatish Balay #undef __FUNC__
13342d61bbb3SSatish Balay #define __FUNC__ "MatGetInfo_SeqBAIJ"
13352d61bbb3SSatish Balay int MatGetInfo_SeqBAIJ(Mat A,MatInfoType flag,MatInfo *info)
13362d61bbb3SSatish Balay {
13372d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
13382d61bbb3SSatish Balay 
13392d61bbb3SSatish Balay   PetscFunctionBegin;
13402d61bbb3SSatish Balay   info->rows_global    = (double)a->m;
13412d61bbb3SSatish Balay   info->columns_global = (double)a->n;
13422d61bbb3SSatish Balay   info->rows_local     = (double)a->m;
13432d61bbb3SSatish Balay   info->columns_local  = (double)a->n;
13442d61bbb3SSatish Balay   info->block_size     = a->bs2;
13452d61bbb3SSatish Balay   info->nz_allocated   = a->maxnz;
13462d61bbb3SSatish Balay   info->nz_used        = a->bs2*a->nz;
13472d61bbb3SSatish Balay   info->nz_unneeded    = (double)(info->nz_allocated - info->nz_used);
13482d61bbb3SSatish Balay   /*  if (info->nz_unneeded != A->info.nz_unneeded)
13492d61bbb3SSatish Balay     printf("space descrepancy: maxnz-nz = %d, nz_unneeded = %d\n",(int)info->nz_unneeded,(int)A->info.nz_unneeded); */
13502d61bbb3SSatish Balay   info->assemblies   = A->num_ass;
13512d61bbb3SSatish Balay   info->mallocs      = a->reallocs;
13522d61bbb3SSatish Balay   info->memory       = A->mem;
13532d61bbb3SSatish Balay   if (A->factor) {
13542d61bbb3SSatish Balay     info->fill_ratio_given  = A->info.fill_ratio_given;
13552d61bbb3SSatish Balay     info->fill_ratio_needed = A->info.fill_ratio_needed;
13562d61bbb3SSatish Balay     info->factor_mallocs    = A->info.factor_mallocs;
13572d61bbb3SSatish Balay   } else {
13582d61bbb3SSatish Balay     info->fill_ratio_given  = 0;
13592d61bbb3SSatish Balay     info->fill_ratio_needed = 0;
13602d61bbb3SSatish Balay     info->factor_mallocs    = 0;
13612d61bbb3SSatish Balay   }
13622d61bbb3SSatish Balay   PetscFunctionReturn(0);
13632d61bbb3SSatish Balay }
13642d61bbb3SSatish Balay 
13652d61bbb3SSatish Balay 
13662d61bbb3SSatish Balay #undef __FUNC__
13672d61bbb3SSatish Balay #define __FUNC__ "MatZeroEntries_SeqBAIJ"
13682d61bbb3SSatish Balay int MatZeroEntries_SeqBAIJ(Mat A)
13692d61bbb3SSatish Balay {
13702d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
1371549d3d68SSatish Balay   int         ierr;
13722d61bbb3SSatish Balay 
13732d61bbb3SSatish Balay   PetscFunctionBegin;
1374549d3d68SSatish Balay   ierr = PetscMemzero(a->a,a->bs2*a->i[a->mbs]*sizeof(MatScalar));CHKERRQ(ierr);
13752d61bbb3SSatish Balay   PetscFunctionReturn(0);
13762d61bbb3SSatish Balay }
1377