xref: /petsc/src/mat/impls/baij/seq/baij2.c (revision 3b95cb0efbec1b084a4bcf5ff0914a41873f2466)
1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER
2*3b95cb0eSSatish Balay static char vcid[] = "$Id: baij2.c,v 1.40 1999/03/17 23:23:13 bsmith 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);
72a3192f15SSatish Balay   PetscFree(nidx);
73218c64b6SSatish Balay   PetscFree(nidx2);
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;
90736121d4SSatish Balay   ierr = ISSorted(iscol,(PetscTruth*)&i);
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);
101736121d4SSatish Balay   PetscMemzero(smap,oldcols*sizeof(int));
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     }
122736121d4SSatish Balay     PetscMemzero(c->ilen,c->mbs*sizeof(int));
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;
1393f1db9ecSBarry Smith         PetscMemcpy(mat_a,a->a+k*bs2,bs2*sizeof(MatScalar)); mat_a+=bs2;
140736121d4SSatish Balay         (*mat_ilen)++;
141736121d4SSatish Balay       }
142736121d4SSatish Balay     }
143736121d4SSatish Balay   }
144218c64b6SSatish Balay 
145736121d4SSatish Balay   /* Free work space */
146736121d4SSatish Balay   ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr);
147736121d4SSatish Balay   PetscFree(smap); PetscFree(lens);
1486d4a8577SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
1496d4a8577SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
150736121d4SSatish Balay 
151736121d4SSatish Balay   ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr);
152736121d4SSatish Balay   *B = C;
1533a40ed3dSBarry Smith   PetscFunctionReturn(0);
154736121d4SSatish Balay }
155736121d4SSatish Balay 
1565615d1e5SSatish Balay #undef __FUNC__
1575615d1e5SSatish Balay #define __FUNC__ "MatGetSubMatrix_SeqBAIJ"
1587b2a1423SBarry Smith int MatGetSubMatrix_SeqBAIJ(Mat A,IS isrow,IS iscol,int cs,MatReuse scall,Mat *B)
159218c64b6SSatish Balay {
160218c64b6SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
161218c64b6SSatish Balay   IS          is1,is2;
162218c64b6SSatish Balay   int         *vary,*iary,*irow,*icol,nrows,ncols,i,ierr,bs=a->bs,count;
163218c64b6SSatish Balay 
1643a40ed3dSBarry Smith   PetscFunctionBegin;
165218c64b6SSatish Balay   ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr);
166218c64b6SSatish Balay   ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr);
167218c64b6SSatish Balay   ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr);
168218c64b6SSatish Balay   ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr);
169218c64b6SSatish Balay 
170218c64b6SSatish Balay   /* Verify if the indices corespond to each element in a block
171218c64b6SSatish Balay    and form the IS with compressed IS */
172218c64b6SSatish Balay   vary = (int *) PetscMalloc(2*(a->mbs+1)*sizeof(int)); CHKPTRQ(vary);
173218c64b6SSatish Balay   iary = vary + a->mbs;
174218c64b6SSatish Balay   PetscMemzero(vary,(a->mbs)*sizeof(int));
175218c64b6SSatish Balay   for ( i=0; i<nrows; i++) vary[irow[i]/bs]++;
176218c64b6SSatish Balay   count = 0;
177218c64b6SSatish Balay   for (i=0; i<a->mbs; i++) {
1786a6a5d1dSBarry Smith     if (vary[i]!=0 && vary[i]!=bs) SETERRA(1,0,"Index set does not match blocks");
179218c64b6SSatish Balay     if (vary[i]==bs) iary[count++] = i;
180218c64b6SSatish Balay   }
181029af93fSBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF, count, iary,&is1); CHKERRQ(ierr);
182218c64b6SSatish Balay 
183218c64b6SSatish Balay   PetscMemzero(vary,(a->mbs)*sizeof(int));
184218c64b6SSatish Balay   for ( i=0; i<ncols; i++) vary[icol[i]/bs]++;
185218c64b6SSatish Balay   count = 0;
186218c64b6SSatish Balay   for (i=0; i<a->mbs; i++) {
187e3372554SBarry Smith     if (vary[i]!=0 && vary[i]!=bs) SETERRA(1,0,"MatGetSubmatrices_SeqBAIJ:");
188218c64b6SSatish Balay     if (vary[i]==bs) iary[count++] = i;
189218c64b6SSatish Balay   }
190029af93fSBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF, count, iary,&is2); CHKERRQ(ierr);
191218c64b6SSatish Balay   ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr);
192218c64b6SSatish Balay   ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr);
193218c64b6SSatish Balay   PetscFree(vary);
194218c64b6SSatish Balay 
1956a6a5d1dSBarry Smith   ierr = MatGetSubMatrix_SeqBAIJ_Private(A,is1,is2,cs,scall,B); CHKERRQ(ierr);
196218c64b6SSatish Balay   ISDestroy(is1);
197218c64b6SSatish Balay   ISDestroy(is2);
1983a40ed3dSBarry Smith   PetscFunctionReturn(0);
199218c64b6SSatish Balay }
200218c64b6SSatish Balay 
2015615d1e5SSatish Balay #undef __FUNC__
2025615d1e5SSatish Balay #define __FUNC__ "MatGetSubMatrices_SeqBAIJ"
2037b2a1423SBarry Smith int MatGetSubMatrices_SeqBAIJ(Mat A,int n, IS *irow,IS *icol,MatReuse scall,Mat **B)
204736121d4SSatish Balay {
205736121d4SSatish Balay   int ierr,i;
206736121d4SSatish Balay 
2073a40ed3dSBarry Smith   PetscFunctionBegin;
208736121d4SSatish Balay   if (scall == MAT_INITIAL_MATRIX) {
209736121d4SSatish Balay     *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B);
210736121d4SSatish Balay   }
211736121d4SSatish Balay 
212736121d4SSatish Balay   for ( i=0; i<n; i++ ) {
2136a6a5d1dSBarry Smith     ierr = MatGetSubMatrix_SeqBAIJ(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr);
214736121d4SSatish Balay   }
2153a40ed3dSBarry Smith   PetscFunctionReturn(0);
216736121d4SSatish Balay }
217218c64b6SSatish Balay 
218218c64b6SSatish Balay 
2192d61bbb3SSatish Balay /* -------------------------------------------------------*/
2202d61bbb3SSatish Balay /* Should check that shapes of vectors and matrices match */
2212d61bbb3SSatish Balay /* -------------------------------------------------------*/
222eadb2fb4SBarry Smith #include "pinclude/blaslapack.h"
2232d61bbb3SSatish Balay 
2242d61bbb3SSatish Balay #undef __FUNC__
2252d61bbb3SSatish Balay #define __FUNC__ "MatMult_SeqBAIJ_1"
2262d61bbb3SSatish Balay int MatMult_SeqBAIJ_1(Mat A,Vec xx,Vec zz)
2272d61bbb3SSatish Balay {
2282d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
2293f1db9ecSBarry Smith   Scalar          *x,*z,sum;
2303f1db9ecSBarry Smith   MatScalar       *v;
231e1311b90SBarry Smith   int             mbs=a->mbs,i,*idx,*ii,n,ierr;
2322d61bbb3SSatish Balay 
2332d61bbb3SSatish Balay   PetscFunctionBegin;
234e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
235e1311b90SBarry Smith   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
2362d61bbb3SSatish Balay 
2372d61bbb3SSatish Balay   idx   = a->j;
2382d61bbb3SSatish Balay   v     = a->a;
2392d61bbb3SSatish Balay   ii    = a->i;
2402d61bbb3SSatish Balay 
2412d61bbb3SSatish Balay   for ( i=0; i<mbs; i++ ) {
2422d61bbb3SSatish Balay     n    = ii[1] - ii[0]; ii++;
2432d61bbb3SSatish Balay     sum  = 0.0;
2442d61bbb3SSatish Balay     while (n--) sum += *v++ * x[*idx++];
2452d61bbb3SSatish Balay     z[i] = sum;
2462d61bbb3SSatish Balay   }
247e1311b90SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
248e1311b90SBarry Smith   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
2492d61bbb3SSatish Balay   PLogFlops(2*a->nz - a->m);
2502d61bbb3SSatish Balay   PetscFunctionReturn(0);
2512d61bbb3SSatish Balay }
2522d61bbb3SSatish Balay 
2532d61bbb3SSatish Balay #undef __FUNC__
2542d61bbb3SSatish Balay #define __FUNC__ "MatMult_SeqBAIJ_2"
2552d61bbb3SSatish Balay int MatMult_SeqBAIJ_2(Mat A,Vec xx,Vec zz)
2562d61bbb3SSatish Balay {
2572d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
2583f1db9ecSBarry Smith   Scalar          *x,*z,*xb,sum1,sum2;
259e1311b90SBarry Smith   Scalar          x1,x2;
2603f1db9ecSBarry Smith   MatScalar       *v;
261e1311b90SBarry Smith   int             ierr,mbs=a->mbs,i,*idx,*ii,j,n;
2622d61bbb3SSatish Balay 
2632d61bbb3SSatish Balay   PetscFunctionBegin;
264e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
265e1311b90SBarry Smith   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
2662d61bbb3SSatish Balay 
2672d61bbb3SSatish Balay   idx   = a->j;
2682d61bbb3SSatish Balay   v     = a->a;
2692d61bbb3SSatish Balay   ii    = a->i;
2702d61bbb3SSatish Balay 
2712d61bbb3SSatish Balay   for ( i=0; i<mbs; i++ ) {
2722d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
2732d61bbb3SSatish Balay     sum1 = 0.0; sum2 = 0.0;
2742d61bbb3SSatish Balay     for ( j=0; j<n; j++ ) {
2752d61bbb3SSatish Balay       xb = x + 2*(*idx++); x1 = xb[0]; x2 = xb[1];
2762d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[2]*x2;
2772d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[3]*x2;
2782d61bbb3SSatish Balay       v += 4;
2792d61bbb3SSatish Balay     }
2802d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2;
2812d61bbb3SSatish Balay     z += 2;
2822d61bbb3SSatish Balay   }
283e1311b90SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
284e1311b90SBarry Smith   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
28515091d37SBarry Smith   PLogFlops(8*a->nz - a->m);
2862d61bbb3SSatish Balay   PetscFunctionReturn(0);
2872d61bbb3SSatish Balay }
2882d61bbb3SSatish Balay 
2892d61bbb3SSatish Balay #undef __FUNC__
2902d61bbb3SSatish Balay #define __FUNC__ "MatMult_SeqBAIJ_3"
2912d61bbb3SSatish Balay int MatMult_SeqBAIJ_3(Mat A,Vec xx,Vec zz)
2922d61bbb3SSatish Balay {
2932d61bbb3SSatish Balay   Mat_SeqBAIJ  *a = (Mat_SeqBAIJ *) A->data;
2943f1db9ecSBarry Smith   Scalar       *x,*z,*xb,sum1,sum2,sum3,x1,x2,x3;
2953f1db9ecSBarry Smith   MatScalar    *v;
296e1311b90SBarry Smith   int          ierr,mbs=a->mbs,i,*idx,*ii,j,n;
2972d61bbb3SSatish Balay 
298fee21e36SBarry Smith #if defined(HAVE_PRAGMA_DISJOINT)
299fee21e36SBarry Smith #pragma disjoint(*v,*z,*xb)
300fee21e36SBarry Smith #endif
301fee21e36SBarry Smith 
3022d61bbb3SSatish Balay   PetscFunctionBegin;
303e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
304e1311b90SBarry Smith   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
3052d61bbb3SSatish Balay 
3062d61bbb3SSatish Balay   idx   = a->j;
3072d61bbb3SSatish Balay   v     = a->a;
3082d61bbb3SSatish Balay   ii    = a->i;
3092d61bbb3SSatish Balay 
3102d61bbb3SSatish Balay   for ( i=0; i<mbs; i++ ) {
3112d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
3122d61bbb3SSatish Balay     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0;
3132d61bbb3SSatish Balay     for ( j=0; j<n; j++ ) {
3142d61bbb3SSatish Balay       xb = x + 3*(*idx++); x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
3152d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[3]*x2 + v[6]*x3;
3162d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[4]*x2 + v[7]*x3;
3172d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[5]*x2 + v[8]*x3;
3182d61bbb3SSatish Balay       v += 9;
3192d61bbb3SSatish Balay     }
3202d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3;
3212d61bbb3SSatish Balay     z += 3;
3222d61bbb3SSatish Balay   }
323e1311b90SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
324e1311b90SBarry Smith   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
3252d61bbb3SSatish Balay   PLogFlops(18*a->nz - a->m);
3262d61bbb3SSatish Balay   PetscFunctionReturn(0);
3272d61bbb3SSatish Balay }
3282d61bbb3SSatish Balay 
3292d61bbb3SSatish Balay #undef __FUNC__
3302d61bbb3SSatish Balay #define __FUNC__ "MatMult_SeqBAIJ_4"
3312d61bbb3SSatish Balay int MatMult_SeqBAIJ_4(Mat A,Vec xx,Vec zz)
3322d61bbb3SSatish Balay {
3332d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
3343f1db9ecSBarry Smith   Scalar          *x,*z,*xb,sum1,sum2,sum3,sum4,x1,x2,x3,x4;
3353f1db9ecSBarry Smith   MatScalar       *v;
336e1311b90SBarry Smith   int             ierr,mbs=a->mbs,i,*idx,*ii,j,n;
3372d61bbb3SSatish Balay 
3382d61bbb3SSatish Balay   PetscFunctionBegin;
339e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
340e1311b90SBarry Smith   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
3412d61bbb3SSatish Balay 
3422d61bbb3SSatish Balay   idx   = a->j;
3432d61bbb3SSatish Balay   v     = a->a;
3442d61bbb3SSatish Balay   ii    = a->i;
3452d61bbb3SSatish Balay 
3462d61bbb3SSatish Balay   for ( i=0; i<mbs; i++ ) {
3472d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
3482d61bbb3SSatish Balay     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0;
3492d61bbb3SSatish Balay     for ( j=0; j<n; j++ ) {
3502d61bbb3SSatish Balay       xb = x + 4*(*idx++);
3512d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
3522d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[4]*x2 + v[8]*x3   + v[12]*x4;
3532d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[5]*x2 + v[9]*x3   + v[13]*x4;
3542d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[6]*x2 + v[10]*x3  + v[14]*x4;
3552d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[7]*x2 + v[11]*x3  + v[15]*x4;
3562d61bbb3SSatish Balay       v += 16;
3572d61bbb3SSatish Balay     }
3582d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4;
3592d61bbb3SSatish Balay     z += 4;
3602d61bbb3SSatish Balay   }
361e1311b90SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
362e1311b90SBarry Smith   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
3632d61bbb3SSatish Balay   PLogFlops(32*a->nz - a->m);
3642d61bbb3SSatish Balay   PetscFunctionReturn(0);
3652d61bbb3SSatish Balay }
3662d61bbb3SSatish Balay 
3672d61bbb3SSatish Balay #undef __FUNC__
3682d61bbb3SSatish Balay #define __FUNC__ "MatMult_SeqBAIJ_5"
3692d61bbb3SSatish Balay int MatMult_SeqBAIJ_5(Mat A,Vec xx,Vec zz)
3702d61bbb3SSatish Balay {
3712d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
3723f1db9ecSBarry Smith   Scalar          sum1,sum2,sum3,sum4,sum5,x1,x2,x3,x4,x5,*xb,*z,*x;
3733f1db9ecSBarry Smith   MatScalar       *v;
374e1311b90SBarry Smith   int             ierr,mbs=a->mbs,i,*idx,*ii,j,n;
3752d61bbb3SSatish Balay 
376e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
377e1311b90SBarry Smith   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
3782d61bbb3SSatish Balay 
3792d61bbb3SSatish Balay   idx   = a->j;
3802d61bbb3SSatish Balay   v     = a->a;
3812d61bbb3SSatish Balay   ii    = a->i;
3822d61bbb3SSatish Balay 
3832d61bbb3SSatish Balay   for ( i=0; i<mbs; i++ ) {
3842d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
3852d61bbb3SSatish Balay     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0;
3862d61bbb3SSatish Balay     for ( j=0; j<n; j++ ) {
3872d61bbb3SSatish Balay       xb = x + 5*(*idx++);
3882d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4];
3892d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[5]*x2 + v[10]*x3  + v[15]*x4 + v[20]*x5;
3902d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[6]*x2 + v[11]*x3  + v[16]*x4 + v[21]*x5;
3912d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[7]*x2 + v[12]*x3  + v[17]*x4 + v[22]*x5;
3922d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[8]*x2 + v[13]*x3  + v[18]*x4 + v[23]*x5;
3932d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[9]*x2 + v[14]*x3  + v[19]*x4 + v[24]*x5;
3942d61bbb3SSatish Balay       v += 25;
3952d61bbb3SSatish Balay     }
3962d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5;
3972d61bbb3SSatish Balay     z += 5;
3982d61bbb3SSatish Balay   }
399e1311b90SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
400e1311b90SBarry Smith   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
4012d61bbb3SSatish Balay   PLogFlops(50*a->nz - a->m);
4022d61bbb3SSatish Balay   PetscFunctionReturn(0);
4032d61bbb3SSatish Balay }
4042d61bbb3SSatish Balay 
40515091d37SBarry Smith 
40615091d37SBarry Smith #undef __FUNC__
40715091d37SBarry Smith #define __FUNC__ "MatMult_SeqBAIJ_6"
40815091d37SBarry Smith int MatMult_SeqBAIJ_6(Mat A,Vec xx,Vec zz)
40915091d37SBarry Smith {
41015091d37SBarry Smith   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
41115091d37SBarry Smith   Scalar          *x,*z,*xb,sum1,sum2,sum3,sum4,sum5,sum6;
41215091d37SBarry Smith   Scalar          x1,x2,x3,x4,x5,x6;
41315091d37SBarry Smith   MatScalar       *v;
41415091d37SBarry Smith   int             ierr,mbs=a->mbs,i,*idx,*ii,j,n;
41515091d37SBarry Smith 
41615091d37SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
41715091d37SBarry Smith   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
41815091d37SBarry Smith 
41915091d37SBarry Smith   idx   = a->j;
42015091d37SBarry Smith   v     = a->a;
42115091d37SBarry Smith   ii    = a->i;
42215091d37SBarry Smith 
42315091d37SBarry Smith   for ( i=0; i<mbs; i++ ) {
42415091d37SBarry Smith     n  = ii[1] - ii[0]; ii++;
42515091d37SBarry Smith     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0;
42615091d37SBarry Smith     for ( j=0; j<n; j++ ) {
42715091d37SBarry Smith       xb = x + 6*(*idx++);
42815091d37SBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5];
42915091d37SBarry Smith       sum1 += v[0]*x1 + v[6]*x2  + v[12]*x3  + v[18]*x4 + v[24]*x5 + v[30]*x6;
43015091d37SBarry Smith       sum2 += v[1]*x1 + v[7]*x2  + v[13]*x3  + v[19]*x4 + v[25]*x5 + v[31]*x6;
43115091d37SBarry Smith       sum3 += v[2]*x1 + v[8]*x2  + v[14]*x3  + v[20]*x4 + v[26]*x5 + v[32]*x6;
43215091d37SBarry Smith       sum4 += v[3]*x1 + v[9]*x2  + v[15]*x3  + v[21]*x4 + v[27]*x5 + v[33]*x6;
43315091d37SBarry Smith       sum5 += v[4]*x1 + v[10]*x2 + v[16]*x3  + v[22]*x4 + v[28]*x5 + v[34]*x6;
43415091d37SBarry Smith       sum6 += v[5]*x1 + v[11]*x2 + v[17]*x3  + v[23]*x4 + v[29]*x5 + v[35]*x6;
43515091d37SBarry Smith       v += 36;
43615091d37SBarry Smith     }
43715091d37SBarry Smith     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6;
43815091d37SBarry Smith     z += 6;
43915091d37SBarry Smith   }
44015091d37SBarry Smith 
44115091d37SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
44215091d37SBarry Smith   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
44315091d37SBarry Smith   PLogFlops(72*a->nz - a->m);
44415091d37SBarry Smith   PetscFunctionReturn(0);
44515091d37SBarry Smith }
4462d61bbb3SSatish Balay #undef __FUNC__
4472d61bbb3SSatish Balay #define __FUNC__ "MatMult_SeqBAIJ_7"
4482d61bbb3SSatish Balay int MatMult_SeqBAIJ_7(Mat A,Vec xx,Vec zz)
4492d61bbb3SSatish Balay {
4502d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
4513f1db9ecSBarry Smith   Scalar          *x,*z,*xb,sum1,sum2,sum3,sum4,sum5,sum6,sum7;
452e1311b90SBarry Smith   Scalar          x1,x2,x3,x4,x5,x6,x7;
4533f1db9ecSBarry Smith   MatScalar       *v;
454e1311b90SBarry Smith   int             ierr,mbs=a->mbs,i,*idx,*ii,j,n;
4552d61bbb3SSatish Balay 
456e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
457e1311b90SBarry Smith   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
4582d61bbb3SSatish Balay 
4592d61bbb3SSatish Balay   idx   = a->j;
4602d61bbb3SSatish Balay   v     = a->a;
4612d61bbb3SSatish Balay   ii    = a->i;
4622d61bbb3SSatish Balay 
4632d61bbb3SSatish Balay   for ( i=0; i<mbs; i++ ) {
4642d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
4652d61bbb3SSatish Balay     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0; sum7 = 0.0;
4662d61bbb3SSatish Balay     for ( j=0; j<n; j++ ) {
4672d61bbb3SSatish Balay       xb = x + 7*(*idx++);
4682d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6];
4692d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[7]*x2  + v[14]*x3  + v[21]*x4 + v[28]*x5 + v[35]*x6 + v[42]*x7;
4702d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[8]*x2  + v[15]*x3  + v[22]*x4 + v[29]*x5 + v[36]*x6 + v[43]*x7;
4712d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[9]*x2  + v[16]*x3  + v[23]*x4 + v[30]*x5 + v[37]*x6 + v[44]*x7;
4722d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[10]*x2 + v[17]*x3  + v[24]*x4 + v[31]*x5 + v[38]*x6 + v[45]*x7;
4732d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[11]*x2 + v[18]*x3  + v[25]*x4 + v[32]*x5 + v[39]*x6 + v[46]*x7;
4742d61bbb3SSatish Balay       sum6 += v[5]*x1 + v[12]*x2 + v[19]*x3  + v[26]*x4 + v[33]*x5 + v[40]*x6 + v[47]*x7;
4752d61bbb3SSatish Balay       sum7 += v[6]*x1 + v[13]*x2 + v[20]*x3  + v[27]*x4 + v[34]*x5 + v[41]*x6 + v[48]*x7;
4762d61bbb3SSatish Balay       v += 49;
4772d61bbb3SSatish Balay     }
4782d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
4792d61bbb3SSatish Balay     z += 7;
4802d61bbb3SSatish Balay   }
4812d61bbb3SSatish Balay 
482e1311b90SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
483e1311b90SBarry Smith   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
4842d61bbb3SSatish Balay   PLogFlops(98*a->nz - a->m);
4852d61bbb3SSatish Balay   PetscFunctionReturn(0);
4862d61bbb3SSatish Balay }
4872d61bbb3SSatish Balay 
4883f1db9ecSBarry Smith /*
4893f1db9ecSBarry Smith     This will not work with MatScalar == float because it calls the BLAS
4903f1db9ecSBarry Smith */
4912d61bbb3SSatish Balay #undef __FUNC__
4922d61bbb3SSatish Balay #define __FUNC__ "MatMult_SeqBAIJ_N"
4932d61bbb3SSatish Balay int MatMult_SeqBAIJ_N(Mat A,Vec xx,Vec zz)
4942d61bbb3SSatish Balay {
4952d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
4963f1db9ecSBarry Smith   Scalar          *x,*z,*xb, *work,*workt;
4973f1db9ecSBarry Smith   MatScalar       *v;
498e1311b90SBarry Smith   int             ierr,mbs=a->mbs,i,*idx,*ii,bs=a->bs,j,n,bs2=a->bs2;
4993f1db9ecSBarry Smith   int             ncols,k;
5002d61bbb3SSatish Balay 
5012d61bbb3SSatish Balay   PetscFunctionBegin;
502e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
503e1311b90SBarry Smith   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
5042d61bbb3SSatish Balay 
5052d61bbb3SSatish Balay   idx   = a->j;
5062d61bbb3SSatish Balay   v     = a->a;
5072d61bbb3SSatish Balay   ii    = a->i;
508218c64b6SSatish Balay 
509218c64b6SSatish Balay 
5102d61bbb3SSatish Balay   if (!a->mult_work) {
5112d61bbb3SSatish Balay     k = PetscMax(a->m,a->n);
5122d61bbb3SSatish Balay     a->mult_work = (Scalar *) PetscMalloc((k+1)*sizeof(Scalar));CHKPTRQ(a->mult_work);
5132d61bbb3SSatish Balay   }
5142d61bbb3SSatish Balay   work = a->mult_work;
5152d61bbb3SSatish Balay   for ( i=0; i<mbs; i++ ) {
5162d61bbb3SSatish Balay     n     = ii[1] - ii[0]; ii++;
5172d61bbb3SSatish Balay     ncols = n*bs;
5182d61bbb3SSatish Balay     workt = work;
5192d61bbb3SSatish Balay     for ( j=0; j<n; j++ ) {
5202d61bbb3SSatish Balay       xb = x + bs*(*idx++);
5212d61bbb3SSatish Balay       for ( k=0; k<bs; k++ ) workt[k] = xb[k];
5222d61bbb3SSatish Balay       workt += bs;
5232d61bbb3SSatish Balay     }
5243f1db9ecSBarry Smith     Kernel_w_gets_Ar_times_v(bs,ncols,work,v,z);
5253f1db9ecSBarry Smith     /* LAgemv_("N",&bs,&ncols,&_DOne,v,&bs,work,&_One,&_DZero,z,&_One); */
5262d61bbb3SSatish Balay     v += n*bs2;
5272d61bbb3SSatish Balay     z += bs;
5282d61bbb3SSatish Balay   }
529e1311b90SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
530e1311b90SBarry Smith   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
5312d61bbb3SSatish Balay   PLogFlops(2*a->nz*bs2 - a->m);
5322d61bbb3SSatish Balay   PetscFunctionReturn(0);
5332d61bbb3SSatish Balay }
5342d61bbb3SSatish Balay 
5352d61bbb3SSatish Balay #undef __FUNC__
5362d61bbb3SSatish Balay #define __FUNC__ "MatMultAdd_SeqBAIJ_1"
5372d61bbb3SSatish Balay int MatMultAdd_SeqBAIJ_1(Mat A,Vec xx,Vec yy,Vec zz)
5382d61bbb3SSatish Balay {
5392d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
5403f1db9ecSBarry Smith   Scalar          *x,*y,*z,sum;
5413f1db9ecSBarry Smith   MatScalar       *v;
542e1311b90SBarry Smith   int             ierr,mbs=a->mbs,i,*idx,*ii,n;
5432d61bbb3SSatish Balay 
5442d61bbb3SSatish Balay   PetscFunctionBegin;
545e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
546e1311b90SBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
5472e8a6d31SBarry Smith   if (zz != yy) {
548e1311b90SBarry Smith     ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
5492e8a6d31SBarry Smith   } else {
5502e8a6d31SBarry Smith     z = y;
5512e8a6d31SBarry Smith   }
5522d61bbb3SSatish Balay 
5532d61bbb3SSatish Balay   idx   = a->j;
5542d61bbb3SSatish Balay   v     = a->a;
5552d61bbb3SSatish Balay   ii    = a->i;
5562d61bbb3SSatish Balay 
5572d61bbb3SSatish Balay   for ( i=0; i<mbs; i++ ) {
5582d61bbb3SSatish Balay     n    = ii[1] - ii[0]; ii++;
5592d61bbb3SSatish Balay     sum  = y[i];
5602d61bbb3SSatish Balay     while (n--) sum += *v++ * x[*idx++];
5612d61bbb3SSatish Balay     z[i] = sum;
5622d61bbb3SSatish Balay   }
563e1311b90SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
564e1311b90SBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
5652e8a6d31SBarry Smith   if (zz != yy) {
566e1311b90SBarry Smith     ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
5672e8a6d31SBarry Smith   }
5682d61bbb3SSatish Balay   PLogFlops(2*a->nz);
5692d61bbb3SSatish Balay   PetscFunctionReturn(0);
5702d61bbb3SSatish Balay }
5712d61bbb3SSatish Balay 
5722d61bbb3SSatish Balay #undef __FUNC__
5732d61bbb3SSatish Balay #define __FUNC__ "MatMultAdd_SeqBAIJ_2"
5742d61bbb3SSatish Balay int MatMultAdd_SeqBAIJ_2(Mat A,Vec xx,Vec yy,Vec zz)
5752d61bbb3SSatish Balay {
5762d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
5773f1db9ecSBarry Smith   Scalar          *x,*y,*z,*xb,sum1,sum2;
578e1311b90SBarry Smith   Scalar          x1,x2;
5793f1db9ecSBarry Smith   MatScalar       *v;
580e1311b90SBarry Smith   int             ierr,mbs=a->mbs,i,*idx,*ii,j,n;
5812d61bbb3SSatish Balay 
5822d61bbb3SSatish Balay   PetscFunctionBegin;
583e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
584e1311b90SBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
5852e8a6d31SBarry Smith   if (zz != yy) {
586e1311b90SBarry Smith     ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
5872e8a6d31SBarry Smith   } else {
5882e8a6d31SBarry Smith     z = y;
5892e8a6d31SBarry Smith   }
5902d61bbb3SSatish Balay 
5912d61bbb3SSatish Balay   idx   = a->j;
5922d61bbb3SSatish Balay   v     = a->a;
5932d61bbb3SSatish Balay   ii    = a->i;
5942d61bbb3SSatish Balay 
5952d61bbb3SSatish Balay   for ( i=0; i<mbs; i++ ) {
5962d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
5972d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1];
5982d61bbb3SSatish Balay     for ( j=0; j<n; j++ ) {
5992d61bbb3SSatish Balay       xb = x + 2*(*idx++); x1 = xb[0]; x2 = xb[1];
6002d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[2]*x2;
6012d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[3]*x2;
6022d61bbb3SSatish Balay       v += 4;
6032d61bbb3SSatish Balay     }
6042d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2;
6052d61bbb3SSatish Balay     z += 2; y += 2;
6062d61bbb3SSatish Balay   }
607e1311b90SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
608e1311b90SBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
6092e8a6d31SBarry Smith   if (zz != yy) {
610e1311b90SBarry Smith     ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
6112e8a6d31SBarry Smith   }
6122d61bbb3SSatish Balay   PLogFlops(4*a->nz);
6132d61bbb3SSatish Balay   PetscFunctionReturn(0);
6142d61bbb3SSatish Balay }
6152d61bbb3SSatish Balay 
6162d61bbb3SSatish Balay #undef __FUNC__
6172d61bbb3SSatish Balay #define __FUNC__ "MatMultAdd_SeqBAIJ_3"
6182d61bbb3SSatish Balay int MatMultAdd_SeqBAIJ_3(Mat A,Vec xx,Vec yy,Vec zz)
6192d61bbb3SSatish Balay {
6202d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
6213f1db9ecSBarry Smith   Scalar          *x,*y,*z,*xb,sum1,sum2,sum3,x1,x2,x3;
6223f1db9ecSBarry Smith   MatScalar       *v;
623e1311b90SBarry Smith   int             ierr,mbs=a->mbs,i,*idx,*ii,j,n;
6242d61bbb3SSatish Balay 
6252d61bbb3SSatish Balay   PetscFunctionBegin;
626e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
627e1311b90SBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
6282e8a6d31SBarry Smith   if (zz != yy) {
629e1311b90SBarry Smith     ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
6302e8a6d31SBarry Smith   } else {
6312e8a6d31SBarry Smith     z = y;
6322e8a6d31SBarry Smith   }
6332d61bbb3SSatish Balay 
6342d61bbb3SSatish Balay   idx   = a->j;
6352d61bbb3SSatish Balay   v     = a->a;
6362d61bbb3SSatish Balay   ii    = a->i;
6372d61bbb3SSatish Balay 
6382d61bbb3SSatish Balay   for ( i=0; i<mbs; i++ ) {
6392d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
6402d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2];
6412d61bbb3SSatish Balay     for ( j=0; j<n; j++ ) {
6422d61bbb3SSatish Balay       xb = x + 3*(*idx++); x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
6432d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[3]*x2 + v[6]*x3;
6442d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[4]*x2 + v[7]*x3;
6452d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[5]*x2 + v[8]*x3;
6462d61bbb3SSatish Balay       v += 9;
6472d61bbb3SSatish Balay     }
6482d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3;
6492d61bbb3SSatish Balay     z += 3; y += 3;
6502d61bbb3SSatish Balay   }
651e1311b90SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
652e1311b90SBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
6532e8a6d31SBarry Smith   if (zz != yy) {
654e1311b90SBarry Smith     ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
6552e8a6d31SBarry Smith   }
6562d61bbb3SSatish Balay   PLogFlops(18*a->nz);
6572d61bbb3SSatish Balay   PetscFunctionReturn(0);
6582d61bbb3SSatish Balay }
6592d61bbb3SSatish Balay 
6602d61bbb3SSatish Balay #undef __FUNC__
6612d61bbb3SSatish Balay #define __FUNC__ "MatMultAdd_SeqBAIJ_4"
6622d61bbb3SSatish Balay int MatMultAdd_SeqBAIJ_4(Mat A,Vec xx,Vec yy,Vec zz)
6632d61bbb3SSatish Balay {
6642d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
6653f1db9ecSBarry Smith   Scalar          *x,*y,*z,*xb,sum1,sum2,sum3,sum4,x1,x2,x3,x4;
6663f1db9ecSBarry Smith   MatScalar       *v;
667e1311b90SBarry Smith   int             ierr,mbs=a->mbs,i,*idx,*ii;
6682d61bbb3SSatish Balay   int             j,n;
6692d61bbb3SSatish Balay 
6702d61bbb3SSatish Balay   PetscFunctionBegin;
671e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
672e1311b90SBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
6732e8a6d31SBarry Smith   if (zz != yy) {
674e1311b90SBarry Smith     ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
6752e8a6d31SBarry Smith   } else {
6762e8a6d31SBarry Smith     z = y;
6772e8a6d31SBarry Smith   }
6782d61bbb3SSatish Balay 
6792d61bbb3SSatish Balay   idx   = a->j;
6802d61bbb3SSatish Balay   v     = a->a;
6812d61bbb3SSatish Balay   ii    = a->i;
6822d61bbb3SSatish Balay 
6832d61bbb3SSatish Balay   for ( i=0; i<mbs; i++ ) {
6842d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
6852d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3];
6862d61bbb3SSatish Balay     for ( j=0; j<n; j++ ) {
6872d61bbb3SSatish Balay       xb = x + 4*(*idx++);
6882d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
6892d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[4]*x2 + v[8]*x3   + v[12]*x4;
6902d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[5]*x2 + v[9]*x3   + v[13]*x4;
6912d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[6]*x2 + v[10]*x3  + v[14]*x4;
6922d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[7]*x2 + v[11]*x3  + v[15]*x4;
6932d61bbb3SSatish Balay       v += 16;
6942d61bbb3SSatish Balay     }
6952d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4;
6962d61bbb3SSatish Balay     z += 4; y += 4;
6972d61bbb3SSatish Balay   }
698e1311b90SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
699e1311b90SBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
7002e8a6d31SBarry Smith   if (zz != yy) {
701e1311b90SBarry Smith     ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
7022e8a6d31SBarry Smith   }
7032d61bbb3SSatish Balay   PLogFlops(32*a->nz);
7042d61bbb3SSatish Balay   PetscFunctionReturn(0);
7052d61bbb3SSatish Balay }
7062d61bbb3SSatish Balay 
7072d61bbb3SSatish Balay #undef __FUNC__
7082d61bbb3SSatish Balay #define __FUNC__ "MatMultAdd_SeqBAIJ_5"
7092d61bbb3SSatish Balay int MatMultAdd_SeqBAIJ_5(Mat A,Vec xx,Vec yy,Vec zz)
7102d61bbb3SSatish Balay {
7112d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
7123f1db9ecSBarry Smith   Scalar          *x,*y,*z,*xb,sum1,sum2,sum3,sum4,sum5,x1,x2,x3,x4,x5;
7133f1db9ecSBarry Smith   MatScalar       *v;
714e1311b90SBarry Smith   int             ierr,mbs=a->mbs,i,*idx,*ii,j,n;
7152d61bbb3SSatish Balay 
7162d61bbb3SSatish Balay   PetscFunctionBegin;
717e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
718e1311b90SBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
7192e8a6d31SBarry Smith   if (zz != yy) {
720e1311b90SBarry Smith     ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
7212e8a6d31SBarry Smith   } else {
7222e8a6d31SBarry Smith     z = y;
7232e8a6d31SBarry Smith   }
7242d61bbb3SSatish Balay 
7252d61bbb3SSatish Balay   idx   = a->j;
7262d61bbb3SSatish Balay   v     = a->a;
7272d61bbb3SSatish Balay   ii    = a->i;
7282d61bbb3SSatish Balay 
7292d61bbb3SSatish Balay   for ( i=0; i<mbs; i++ ) {
7302d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
7312d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4];
7322d61bbb3SSatish Balay     for ( j=0; j<n; j++ ) {
7332d61bbb3SSatish Balay       xb = x + 5*(*idx++);
7342d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4];
7352d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[5]*x2 + v[10]*x3  + v[15]*x4 + v[20]*x5;
7362d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[6]*x2 + v[11]*x3  + v[16]*x4 + v[21]*x5;
7372d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[7]*x2 + v[12]*x3  + v[17]*x4 + v[22]*x5;
7382d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[8]*x2 + v[13]*x3  + v[18]*x4 + v[23]*x5;
7392d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[9]*x2 + v[14]*x3  + v[19]*x4 + v[24]*x5;
7402d61bbb3SSatish Balay       v += 25;
7412d61bbb3SSatish Balay     }
7422d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5;
7432d61bbb3SSatish Balay     z += 5; y += 5;
7442d61bbb3SSatish Balay   }
745e1311b90SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
746e1311b90SBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
7472e8a6d31SBarry Smith   if (zz != yy) {
748e1311b90SBarry Smith     ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
7492e8a6d31SBarry Smith   }
7502d61bbb3SSatish Balay   PLogFlops(50*a->nz);
7512d61bbb3SSatish Balay   PetscFunctionReturn(0);
7522d61bbb3SSatish Balay }
75315091d37SBarry Smith #undef __FUNC__
75415091d37SBarry Smith #define __FUNC__ "MatMultAdd_SeqBAIJ_6"
75515091d37SBarry Smith int MatMultAdd_SeqBAIJ_6(Mat A,Vec xx,Vec yy,Vec zz)
75615091d37SBarry Smith {
75715091d37SBarry Smith   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
75815091d37SBarry Smith   Scalar          *x,*y,*z,*xb,sum1,sum2,sum3,sum4,sum5,sum6;
75915091d37SBarry Smith   Scalar          x1,x2,x3,x4,x5,x6;
76015091d37SBarry Smith   MatScalar       *v;
76115091d37SBarry Smith   int             ierr,mbs=a->mbs,i,*idx,*ii,j,n;
76215091d37SBarry Smith 
76315091d37SBarry Smith   PetscFunctionBegin;
76415091d37SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
76515091d37SBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
76615091d37SBarry Smith   if (zz != yy) {
76715091d37SBarry Smith     ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
76815091d37SBarry Smith   } else {
76915091d37SBarry Smith     z = y;
77015091d37SBarry Smith   }
77115091d37SBarry Smith 
77215091d37SBarry Smith   idx   = a->j;
77315091d37SBarry Smith   v     = a->a;
77415091d37SBarry Smith   ii    = a->i;
77515091d37SBarry Smith 
77615091d37SBarry Smith   for ( i=0; i<mbs; i++ ) {
77715091d37SBarry Smith     n  = ii[1] - ii[0]; ii++;
77815091d37SBarry Smith     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4]; sum6 = y[5];
77915091d37SBarry Smith     for ( j=0; j<n; j++ ) {
780*3b95cb0eSSatish Balay       xb = x + 6*(*idx++);
78115091d37SBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5];
78215091d37SBarry Smith       sum1 += v[0]*x1 + v[6]*x2  + v[12]*x3  + v[18]*x4 + v[24]*x5 + v[30]*x6;
78315091d37SBarry Smith       sum2 += v[1]*x1 + v[7]*x2  + v[13]*x3  + v[19]*x4 + v[25]*x5 + v[31]*x6;
78415091d37SBarry Smith       sum3 += v[2]*x1 + v[8]*x2  + v[14]*x3  + v[20]*x4 + v[26]*x5 + v[32]*x6;
78515091d37SBarry Smith       sum4 += v[3]*x1 + v[9]*x2  + v[15]*x3  + v[21]*x4 + v[27]*x5 + v[33]*x6;
78615091d37SBarry Smith       sum5 += v[4]*x1 + v[10]*x2 + v[16]*x3  + v[22]*x4 + v[28]*x5 + v[34]*x6;
78715091d37SBarry Smith       sum6 += v[5]*x1 + v[11]*x2 + v[17]*x3  + v[23]*x4 + v[29]*x5 + v[35]*x6;
78815091d37SBarry Smith       v += 36;
78915091d37SBarry Smith     }
79015091d37SBarry Smith     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6;
79115091d37SBarry Smith     z += 6; y += 6;
79215091d37SBarry Smith   }
79315091d37SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
79415091d37SBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
79515091d37SBarry Smith   if (zz != yy) {
79615091d37SBarry Smith     ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
79715091d37SBarry Smith   }
79815091d37SBarry Smith   PLogFlops(72*a->nz);
79915091d37SBarry Smith   PetscFunctionReturn(0);
80015091d37SBarry Smith }
8012d61bbb3SSatish Balay 
8022d61bbb3SSatish Balay #undef __FUNC__
8032d61bbb3SSatish Balay #define __FUNC__ "MatMultAdd_SeqBAIJ_7"
8042d61bbb3SSatish Balay int MatMultAdd_SeqBAIJ_7(Mat A,Vec xx,Vec yy,Vec zz)
8052d61bbb3SSatish Balay {
8062d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
8073f1db9ecSBarry Smith   Scalar          *x,*y,*z,*xb,sum1,sum2,sum3,sum4,sum5,sum6,sum7;
808e1311b90SBarry Smith   Scalar          x1,x2,x3,x4,x5,x6,x7;
8093f1db9ecSBarry Smith   MatScalar       *v;
810e1311b90SBarry Smith   int             ierr,mbs=a->mbs,i,*idx,*ii,j,n;
8112d61bbb3SSatish Balay 
8122d61bbb3SSatish Balay   PetscFunctionBegin;
813e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
814e1311b90SBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
8152e8a6d31SBarry Smith   if (zz != yy) {
816e1311b90SBarry Smith     ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
8172e8a6d31SBarry Smith   } else {
8182e8a6d31SBarry Smith     z = y;
8192e8a6d31SBarry Smith   }
8202d61bbb3SSatish Balay 
8212d61bbb3SSatish Balay   idx   = a->j;
8222d61bbb3SSatish Balay   v     = a->a;
8232d61bbb3SSatish Balay   ii    = a->i;
8242d61bbb3SSatish Balay 
8252d61bbb3SSatish Balay   for ( i=0; i<mbs; i++ ) {
8262d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
8272d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4]; sum6 = y[5]; sum7 = y[6];
8282d61bbb3SSatish Balay     for ( j=0; j<n; j++ ) {
8292d61bbb3SSatish Balay       xb = x + 7*(*idx++);
8302d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6];
8312d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[7]*x2  + v[14]*x3  + v[21]*x4 + v[28]*x5 + v[35]*x6 + v[42]*x7;
8322d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[8]*x2  + v[15]*x3  + v[22]*x4 + v[29]*x5 + v[36]*x6 + v[43]*x7;
8332d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[9]*x2  + v[16]*x3  + v[23]*x4 + v[30]*x5 + v[37]*x6 + v[44]*x7;
8342d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[10]*x2 + v[17]*x3  + v[24]*x4 + v[31]*x5 + v[38]*x6 + v[45]*x7;
8352d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[11]*x2 + v[18]*x3  + v[25]*x4 + v[32]*x5 + v[39]*x6 + v[46]*x7;
8362d61bbb3SSatish Balay       sum6 += v[5]*x1 + v[12]*x2 + v[19]*x3  + v[26]*x4 + v[33]*x5 + v[40]*x6 + v[47]*x7;
8372d61bbb3SSatish Balay       sum7 += v[6]*x1 + v[13]*x2 + v[20]*x3  + v[27]*x4 + v[34]*x5 + v[41]*x6 + v[48]*x7;
8382d61bbb3SSatish Balay       v += 49;
8392d61bbb3SSatish Balay     }
8402d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
8412d61bbb3SSatish Balay     z += 7; y += 7;
8422d61bbb3SSatish Balay   }
843e1311b90SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
844e1311b90SBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
8452e8a6d31SBarry Smith   if (zz != yy) {
846e1311b90SBarry Smith     ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
8472e8a6d31SBarry Smith   }
8482d61bbb3SSatish Balay   PLogFlops(98*a->nz);
8492d61bbb3SSatish Balay   PetscFunctionReturn(0);
8502d61bbb3SSatish Balay }
851218c64b6SSatish Balay 
8522d61bbb3SSatish Balay #undef __FUNC__
8532d61bbb3SSatish Balay #define __FUNC__ "MatMultAdd_SeqBAIJ_N"
8542d61bbb3SSatish Balay int MatMultAdd_SeqBAIJ_N(Mat A,Vec xx,Vec yy,Vec zz)
8552d61bbb3SSatish Balay {
8562d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ *) A->data;
8573f1db9ecSBarry Smith   Scalar         *x,*z,*xb,*work,*workt;
8583f1db9ecSBarry Smith   MatScalar      *v;
8592d61bbb3SSatish Balay   int            mbs=a->mbs,i,*idx,*ii,bs=a->bs,j,n,bs2=a->bs2,ierr;
8603f1db9ecSBarry Smith   int            ncols,k;
861218c64b6SSatish Balay 
8622d61bbb3SSatish Balay   PetscFunctionBegin;
8632d61bbb3SSatish Balay   if ( xx != yy) { ierr = VecCopy(yy,zz); CHKERRQ(ierr); }
8642d61bbb3SSatish Balay 
865e1311b90SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
866e1311b90SBarry Smith   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
8672d61bbb3SSatish Balay 
8682d61bbb3SSatish Balay   idx   = a->j;
8692d61bbb3SSatish Balay   v     = a->a;
8702d61bbb3SSatish Balay   ii    = a->i;
8712d61bbb3SSatish Balay 
8722d61bbb3SSatish Balay 
8732d61bbb3SSatish Balay   if (!a->mult_work) {
8742d61bbb3SSatish Balay     k = PetscMax(a->m,a->n);
875544c85edSSatish Balay     a->mult_work = (Scalar *) PetscMalloc((k+1)*sizeof(Scalar));CHKPTRQ(a->mult_work);
8762d61bbb3SSatish Balay   }
8772d61bbb3SSatish Balay   work = a->mult_work;
8782d61bbb3SSatish Balay   for ( i=0; i<mbs; i++ ) {
8792d61bbb3SSatish Balay     n     = ii[1] - ii[0]; ii++;
8802d61bbb3SSatish Balay     ncols = n*bs;
8812d61bbb3SSatish Balay     workt = work;
8822d61bbb3SSatish Balay     for ( j=0; j<n; j++ ) {
8832d61bbb3SSatish Balay       xb = x + bs*(*idx++);
8842d61bbb3SSatish Balay       for ( k=0; k<bs; k++ ) workt[k] = xb[k];
8852d61bbb3SSatish Balay       workt += bs;
8862d61bbb3SSatish Balay     }
8873f1db9ecSBarry Smith     Kernel_w_gets_w_plus_Ar_times_v(bs,ncols,work,v,z);
8883f1db9ecSBarry Smith     /* LAgemv_("N",&bs,&ncols,&_DOne,v,&bs,work,&_One,&_DOne,z,&_One); */
8892d61bbb3SSatish Balay     v += n*bs2;
8902d61bbb3SSatish Balay     z += bs;
8912d61bbb3SSatish Balay   }
892e1311b90SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
893e1311b90SBarry Smith   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
8942d61bbb3SSatish Balay   PLogFlops(2*a->nz*bs2 );
8952d61bbb3SSatish Balay   PetscFunctionReturn(0);
8962d61bbb3SSatish Balay }
8972d61bbb3SSatish Balay 
8982d61bbb3SSatish Balay #undef __FUNC__
8992d61bbb3SSatish Balay #define __FUNC__ "MatMultTrans_SeqBAIJ"
9002d61bbb3SSatish Balay int MatMultTrans_SeqBAIJ(Mat A,Vec xx,Vec zz)
9012d61bbb3SSatish Balay {
9022d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
9032d61bbb3SSatish Balay   Scalar          *xg,*zg,*zb;
9043f1db9ecSBarry Smith   Scalar          *x,*z,*xb,x1,x2,x3,x4,x5;
9053f1db9ecSBarry Smith   MatScalar       *v;
9062d61bbb3SSatish Balay   int             mbs=a->mbs,i,*idx,*ii,*ai=a->i,rval,N=a->n;
9072d61bbb3SSatish Balay   int             bs=a->bs,j,n,bs2=a->bs2,*ib,ierr;
9082d61bbb3SSatish Balay 
9092d61bbb3SSatish Balay 
9102d61bbb3SSatish Balay   PetscFunctionBegin;
911e1311b90SBarry Smith   ierr = VecGetArray(xx,&xg);CHKERRQ(ierr); x = xg;
912e1311b90SBarry Smith   ierr = VecGetArray(zz,&zg);CHKERRQ(ierr); z = zg;
9132d61bbb3SSatish Balay   PetscMemzero(z,N*sizeof(Scalar));
9142d61bbb3SSatish Balay 
9152d61bbb3SSatish Balay   idx   = a->j;
9162d61bbb3SSatish Balay   v     = a->a;
9172d61bbb3SSatish Balay   ii    = a->i;
9182d61bbb3SSatish Balay 
9192d61bbb3SSatish Balay   switch (bs) {
9202d61bbb3SSatish Balay   case 1:
9212d61bbb3SSatish Balay     for ( i=0; i<mbs; i++ ) {
9222d61bbb3SSatish Balay       n  = ii[1] - ii[0]; ii++;
9232d61bbb3SSatish Balay       xb = x + i; x1 = xb[0];
9242d61bbb3SSatish Balay       ib = idx + ai[i];
9252d61bbb3SSatish Balay       for ( j=0; j<n; j++ ) {
9262d61bbb3SSatish Balay         rval    = ib[j];
9272d61bbb3SSatish Balay         z[rval] += *v++ * x1;
9282d61bbb3SSatish Balay       }
9292d61bbb3SSatish Balay     }
9302d61bbb3SSatish Balay     break;
9312d61bbb3SSatish Balay   case 2:
9322d61bbb3SSatish Balay     for ( i=0; i<mbs; i++ ) {
9332d61bbb3SSatish Balay       n  = ii[1] - ii[0]; ii++;
9342d61bbb3SSatish Balay       xb = x + 2*i; x1 = xb[0]; x2 = xb[1];
9352d61bbb3SSatish Balay       ib = idx + ai[i];
9362d61bbb3SSatish Balay       for ( j=0; j<n; j++ ) {
9372d61bbb3SSatish Balay         rval      = ib[j]*2;
9382d61bbb3SSatish Balay         z[rval++] += v[0]*x1 + v[1]*x2;
9392d61bbb3SSatish Balay         z[rval++] += v[2]*x1 + v[3]*x2;
9402d61bbb3SSatish Balay         v += 4;
9412d61bbb3SSatish Balay       }
9422d61bbb3SSatish Balay     }
9432d61bbb3SSatish Balay     break;
9442d61bbb3SSatish Balay   case 3:
9452d61bbb3SSatish Balay     for ( i=0; i<mbs; i++ ) {
9462d61bbb3SSatish Balay       n  = ii[1] - ii[0]; ii++;
9472d61bbb3SSatish Balay       xb = x + 3*i; x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
9482d61bbb3SSatish Balay       ib = idx + ai[i];
9492d61bbb3SSatish Balay       for ( j=0; j<n; j++ ) {
9502d61bbb3SSatish Balay         rval      = ib[j]*3;
9512d61bbb3SSatish Balay         z[rval++] += v[0]*x1 + v[1]*x2 + v[2]*x3;
9522d61bbb3SSatish Balay         z[rval++] += v[3]*x1 + v[4]*x2 + v[5]*x3;
9532d61bbb3SSatish Balay         z[rval++] += v[6]*x1 + v[7]*x2 + v[8]*x3;
9542d61bbb3SSatish Balay         v += 9;
9552d61bbb3SSatish Balay       }
9562d61bbb3SSatish Balay     }
9572d61bbb3SSatish Balay     break;
9582d61bbb3SSatish Balay   case 4:
9592d61bbb3SSatish Balay     for ( i=0; i<mbs; i++ ) {
9602d61bbb3SSatish Balay       n  = ii[1] - ii[0]; ii++;
9612d61bbb3SSatish Balay       xb = x + 4*i; x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
9622d61bbb3SSatish Balay       ib = idx + ai[i];
9632d61bbb3SSatish Balay       for ( j=0; j<n; j++ ) {
9642d61bbb3SSatish Balay         rval      = ib[j]*4;
9652d61bbb3SSatish Balay         z[rval++] +=  v[0]*x1 +  v[1]*x2 +  v[2]*x3 +  v[3]*x4;
9662d61bbb3SSatish Balay         z[rval++] +=  v[4]*x1 +  v[5]*x2 +  v[6]*x3 +  v[7]*x4;
9672d61bbb3SSatish Balay         z[rval++] +=  v[8]*x1 +  v[9]*x2 + v[10]*x3 + v[11]*x4;
9682d61bbb3SSatish Balay         z[rval++] += v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4;
9692d61bbb3SSatish Balay         v += 16;
9702d61bbb3SSatish Balay       }
9712d61bbb3SSatish Balay     }
9722d61bbb3SSatish Balay     break;
9732d61bbb3SSatish Balay   case 5:
9742d61bbb3SSatish Balay     for ( i=0; i<mbs; i++ ) {
9752d61bbb3SSatish Balay       n  = ii[1] - ii[0]; ii++;
9762d61bbb3SSatish Balay       xb = x + 5*i; x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
9772d61bbb3SSatish Balay       x4 = xb[3];   x5 = xb[4];
9782d61bbb3SSatish Balay       ib = idx + ai[i];
9792d61bbb3SSatish Balay       for ( j=0; j<n; j++ ) {
9802d61bbb3SSatish Balay         rval      = ib[j]*5;
9812d61bbb3SSatish Balay         z[rval++] +=  v[0]*x1 +  v[1]*x2 +  v[2]*x3 +  v[3]*x4 +  v[4]*x5;
9822d61bbb3SSatish Balay         z[rval++] +=  v[5]*x1 +  v[6]*x2 +  v[7]*x3 +  v[8]*x4 +  v[9]*x5;
9832d61bbb3SSatish Balay         z[rval++] += v[10]*x1 + v[11]*x2 + v[12]*x3 + v[13]*x4 + v[14]*x5;
9842d61bbb3SSatish Balay         z[rval++] += v[15]*x1 + v[16]*x2 + v[17]*x3 + v[18]*x4 + v[19]*x5;
9852d61bbb3SSatish Balay         z[rval++] += v[20]*x1 + v[21]*x2 + v[22]*x3 + v[23]*x4 + v[24]*x5;
9862d61bbb3SSatish Balay         v += 25;
9872d61bbb3SSatish Balay       }
9882d61bbb3SSatish Balay     }
9892d61bbb3SSatish Balay     break;
9902d61bbb3SSatish Balay       /* block sizes larger then 5 by 5 are handled by BLAS */
9912d61bbb3SSatish Balay     default: {
9923f1db9ecSBarry Smith       int       ncols,k;
9933f1db9ecSBarry Smith       Scalar    *work,*workt;
9943f1db9ecSBarry Smith 
9952d61bbb3SSatish Balay       if (!a->mult_work) {
9962d61bbb3SSatish Balay         k = PetscMax(a->m,a->n);
997544c85edSSatish Balay         a->mult_work = (Scalar *) PetscMalloc((k+1)*sizeof(Scalar));CHKPTRQ(a->mult_work);
9982d61bbb3SSatish Balay       }
9992d61bbb3SSatish Balay       work = a->mult_work;
10002d61bbb3SSatish Balay       for ( i=0; i<mbs; i++ ) {
10012d61bbb3SSatish Balay         n     = ii[1] - ii[0]; ii++;
10022d61bbb3SSatish Balay         ncols = n*bs;
10038184e3edSSatish Balay         PetscMemzero(work,ncols*sizeof(Scalar));
1004d824769bSBarry Smith         Kernel_w_gets_w_plus_trans_Ar_times_v(bs,ncols,x,v,work);
10053f1db9ecSBarry Smith         /* LAgemv_("T",&bs,&ncols,&_DOne,v,&bs,x,&_One,&_DOne,work,&_One); */
10062d61bbb3SSatish Balay         v += n*bs2;
10072d61bbb3SSatish Balay         x += bs;
10082d61bbb3SSatish Balay         workt = work;
10092d61bbb3SSatish Balay         for ( j=0; j<n; j++ ) {
10102d61bbb3SSatish Balay           zb = z + bs*(*idx++);
10112d61bbb3SSatish Balay           for ( k=0; k<bs; k++ ) zb[k] += workt[k] ;
10122d61bbb3SSatish Balay           workt += bs;
10132d61bbb3SSatish Balay         }
10142d61bbb3SSatish Balay       }
10152d61bbb3SSatish Balay     }
10162d61bbb3SSatish Balay   }
10172d61bbb3SSatish Balay   ierr = VecRestoreArray(xx,&xg); CHKERRQ(ierr);
10182d61bbb3SSatish Balay   ierr = VecRestoreArray(zz,&zg); CHKERRQ(ierr);
10192d61bbb3SSatish Balay   PLogFlops(2*a->nz*a->bs2 - a->n);
10202d61bbb3SSatish Balay   PetscFunctionReturn(0);
10212d61bbb3SSatish Balay }
10222d61bbb3SSatish Balay 
10232d61bbb3SSatish Balay #undef __FUNC__
10242d61bbb3SSatish Balay #define __FUNC__ "MatMultTransAdd_SeqBAIJ"
10252d61bbb3SSatish Balay int MatMultTransAdd_SeqBAIJ(Mat A,Vec xx,Vec yy,Vec zz)
10262d61bbb3SSatish Balay 
10272d61bbb3SSatish Balay {
10282d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
10292d61bbb3SSatish Balay   Scalar          *xg,*zg,*zb;
10303f1db9ecSBarry Smith   Scalar          *x,*z,*xb,x1,x2,x3,x4,x5;
10313f1db9ecSBarry Smith   MatScalar       *v;
10322d61bbb3SSatish Balay   int             mbs=a->mbs,i,*idx,*ii,*ai=a->i,rval,N=a->n;
10332d61bbb3SSatish Balay   int             bs=a->bs,j,n,bs2=a->bs2,*ib,ierr;
10342d61bbb3SSatish Balay 
10352d61bbb3SSatish Balay   PetscFunctionBegin;
1036e1311b90SBarry Smith   ierr = VecGetArray(xx,&xg);CHKERRQ(ierr); x = xg;
1037e1311b90SBarry Smith   ierr = VecGetArray(zz,&zg);CHKERRQ(ierr); z = zg;
10382d61bbb3SSatish Balay 
10392d61bbb3SSatish Balay   if ( yy != zz ) { ierr = VecCopy(yy,zz); CHKERRQ(ierr); }
10402d61bbb3SSatish Balay   else PetscMemzero(z,N*sizeof(Scalar));
10412d61bbb3SSatish Balay 
10422d61bbb3SSatish Balay   idx   = a->j;
10432d61bbb3SSatish Balay   v     = a->a;
10442d61bbb3SSatish Balay   ii    = a->i;
10452d61bbb3SSatish Balay 
10462d61bbb3SSatish Balay   switch (bs) {
10472d61bbb3SSatish Balay   case 1:
10482d61bbb3SSatish Balay     for ( i=0; i<mbs; i++ ) {
10492d61bbb3SSatish Balay       n  = ii[1] - ii[0]; ii++;
10502d61bbb3SSatish Balay       xb = x + i; x1 = xb[0];
10512d61bbb3SSatish Balay       ib = idx + ai[i];
10522d61bbb3SSatish Balay       for ( j=0; j<n; j++ ) {
10532d61bbb3SSatish Balay         rval    = ib[j];
10542d61bbb3SSatish Balay         z[rval] += *v++ * x1;
10552d61bbb3SSatish Balay       }
10562d61bbb3SSatish Balay     }
10572d61bbb3SSatish Balay     break;
10582d61bbb3SSatish Balay   case 2:
10592d61bbb3SSatish Balay     for ( i=0; i<mbs; i++ ) {
10602d61bbb3SSatish Balay       n  = ii[1] - ii[0]; ii++;
10612d61bbb3SSatish Balay       xb = x + 2*i; x1 = xb[0]; x2 = xb[1];
10622d61bbb3SSatish Balay       ib = idx + ai[i];
10632d61bbb3SSatish Balay       for ( j=0; j<n; j++ ) {
10642d61bbb3SSatish Balay         rval      = ib[j]*2;
10652d61bbb3SSatish Balay         z[rval++] += v[0]*x1 + v[1]*x2;
10662d61bbb3SSatish Balay         z[rval++] += v[2]*x1 + v[3]*x2;
10672d61bbb3SSatish Balay         v += 4;
10682d61bbb3SSatish Balay       }
10692d61bbb3SSatish Balay     }
10702d61bbb3SSatish Balay     break;
10712d61bbb3SSatish Balay   case 3:
10722d61bbb3SSatish Balay     for ( i=0; i<mbs; i++ ) {
10732d61bbb3SSatish Balay       n  = ii[1] - ii[0]; ii++;
10742d61bbb3SSatish Balay       xb = x + 3*i; x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
10752d61bbb3SSatish Balay       ib = idx + ai[i];
10762d61bbb3SSatish Balay       for ( j=0; j<n; j++ ) {
10772d61bbb3SSatish Balay         rval      = ib[j]*3;
10782d61bbb3SSatish Balay         z[rval++] += v[0]*x1 + v[1]*x2 + v[2]*x3;
10792d61bbb3SSatish Balay         z[rval++] += v[3]*x1 + v[4]*x2 + v[5]*x3;
10802d61bbb3SSatish Balay         z[rval++] += v[6]*x1 + v[7]*x2 + v[8]*x3;
10812d61bbb3SSatish Balay         v += 9;
10822d61bbb3SSatish Balay       }
10832d61bbb3SSatish Balay     }
10842d61bbb3SSatish Balay     break;
10852d61bbb3SSatish Balay   case 4:
10862d61bbb3SSatish Balay     for ( i=0; i<mbs; i++ ) {
10872d61bbb3SSatish Balay       n  = ii[1] - ii[0]; ii++;
10882d61bbb3SSatish Balay       xb = x + 4*i; x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
10892d61bbb3SSatish Balay       ib = idx + ai[i];
10902d61bbb3SSatish Balay       for ( j=0; j<n; j++ ) {
10912d61bbb3SSatish Balay         rval      = ib[j]*4;
10922d61bbb3SSatish Balay         z[rval++] +=  v[0]*x1 +  v[1]*x2 +  v[2]*x3 +  v[3]*x4;
10932d61bbb3SSatish Balay         z[rval++] +=  v[4]*x1 +  v[5]*x2 +  v[6]*x3 +  v[7]*x4;
10942d61bbb3SSatish Balay         z[rval++] +=  v[8]*x1 +  v[9]*x2 + v[10]*x3 + v[11]*x4;
10952d61bbb3SSatish Balay         z[rval++] += v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4;
10962d61bbb3SSatish Balay         v += 16;
10972d61bbb3SSatish Balay       }
10982d61bbb3SSatish Balay     }
10992d61bbb3SSatish Balay     break;
11002d61bbb3SSatish Balay   case 5:
11012d61bbb3SSatish Balay     for ( i=0; i<mbs; i++ ) {
11022d61bbb3SSatish Balay       n  = ii[1] - ii[0]; ii++;
11032d61bbb3SSatish Balay       xb = x + 5*i; x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
11042d61bbb3SSatish Balay       x4 = xb[3];   x5 = xb[4];
11052d61bbb3SSatish Balay       ib = idx + ai[i];
11062d61bbb3SSatish Balay       for ( j=0; j<n; j++ ) {
11072d61bbb3SSatish Balay         rval      = ib[j]*5;
11082d61bbb3SSatish Balay         z[rval++] +=  v[0]*x1 +  v[1]*x2 +  v[2]*x3 +  v[3]*x4 +  v[4]*x5;
11092d61bbb3SSatish Balay         z[rval++] +=  v[5]*x1 +  v[6]*x2 +  v[7]*x3 +  v[8]*x4 +  v[9]*x5;
11102d61bbb3SSatish Balay         z[rval++] += v[10]*x1 + v[11]*x2 + v[12]*x3 + v[13]*x4 + v[14]*x5;
11112d61bbb3SSatish Balay         z[rval++] += v[15]*x1 + v[16]*x2 + v[17]*x3 + v[18]*x4 + v[19]*x5;
11122d61bbb3SSatish Balay         z[rval++] += v[20]*x1 + v[21]*x2 + v[22]*x3 + v[23]*x4 + v[24]*x5;
11132d61bbb3SSatish Balay         v += 25;
11142d61bbb3SSatish Balay       }
11152d61bbb3SSatish Balay     }
11162d61bbb3SSatish Balay     break;
11172d61bbb3SSatish Balay       /* block sizes larger then 5 by 5 are handled by BLAS */
11182d61bbb3SSatish Balay     default: {
11193f1db9ecSBarry Smith       int       ncols,k;
11203f1db9ecSBarry Smith       Scalar    *work,*workt;
11213f1db9ecSBarry Smith 
11222d61bbb3SSatish Balay       if (!a->mult_work) {
11232d61bbb3SSatish Balay         k = PetscMax(a->m,a->n);
1124544c85edSSatish Balay         a->mult_work = (Scalar *) PetscMalloc((k+1)*sizeof(Scalar));CHKPTRQ(a->mult_work);
11252d61bbb3SSatish Balay       }
11262d61bbb3SSatish Balay       work = a->mult_work;
11272d61bbb3SSatish Balay       for ( i=0; i<mbs; i++ ) {
11282d61bbb3SSatish Balay         n     = ii[1] - ii[0]; ii++;
11292d61bbb3SSatish Balay         ncols = n*bs;
11302d61bbb3SSatish Balay         PetscMemzero(work,ncols*sizeof(Scalar));
11313f1db9ecSBarry Smith         Kernel_w_gets_w_plus_trans_Ar_times_v(bs,ncols,x,v,work);
11323f1db9ecSBarry Smith         /* LAgemv_("T",&bs,&ncols,&_DOne,v,&bs,x,&_One,&_DOne,work,&_One); */
11332d61bbb3SSatish Balay         v += n*bs2;
11342d61bbb3SSatish Balay         x += bs;
11352d61bbb3SSatish Balay         workt = work;
11362d61bbb3SSatish Balay         for ( j=0; j<n; j++ ) {
11372d61bbb3SSatish Balay           zb = z + bs*(*idx++);
11382d61bbb3SSatish Balay           for ( k=0; k<bs; k++ ) zb[k] += workt[k] ;
11392d61bbb3SSatish Balay           workt += bs;
11402d61bbb3SSatish Balay         }
11412d61bbb3SSatish Balay       }
11422d61bbb3SSatish Balay     }
11432d61bbb3SSatish Balay   }
11442d61bbb3SSatish Balay   ierr = VecRestoreArray(xx,&xg); CHKERRQ(ierr);
11452d61bbb3SSatish Balay   ierr = VecRestoreArray(zz,&zg); CHKERRQ(ierr);
11462d61bbb3SSatish Balay   PLogFlops(2*a->nz*a->bs2);
11472d61bbb3SSatish Balay   PetscFunctionReturn(0);
11482d61bbb3SSatish Balay }
11492d61bbb3SSatish Balay 
11502d61bbb3SSatish Balay #undef __FUNC__
11512d61bbb3SSatish Balay #define __FUNC__ "MatScale_SeqBAIJ"
11522d61bbb3SSatish Balay int MatScale_SeqBAIJ(Scalar *alpha,Mat inA)
11532d61bbb3SSatish Balay {
11542d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) inA->data;
11552d61bbb3SSatish Balay   int         one = 1, totalnz = a->bs2*a->nz;
11562d61bbb3SSatish Balay 
11572d61bbb3SSatish Balay   PetscFunctionBegin;
11582d61bbb3SSatish Balay   BLscal_( &totalnz, alpha, a->a, &one );
11592d61bbb3SSatish Balay   PLogFlops(totalnz);
11602d61bbb3SSatish Balay   PetscFunctionReturn(0);
11612d61bbb3SSatish Balay }
11622d61bbb3SSatish Balay 
11632d61bbb3SSatish Balay #undef __FUNC__
11642d61bbb3SSatish Balay #define __FUNC__ "MatNorm_SeqBAIJ"
11652d61bbb3SSatish Balay int MatNorm_SeqBAIJ(Mat A,NormType type,double *norm)
11662d61bbb3SSatish Balay {
11672d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
11683f1db9ecSBarry Smith   MatScalar   *v = a->a;
11692d61bbb3SSatish Balay   double      sum = 0.0;
1170253ffaf7SBarry Smith   int         i,j,k,bs = a->bs, nz=a->nz,bs2=a->bs2;
11712d61bbb3SSatish Balay 
11722d61bbb3SSatish Balay   PetscFunctionBegin;
11732d61bbb3SSatish Balay   if (type == NORM_FROBENIUS) {
11742d61bbb3SSatish Balay     for (i=0; i< bs2*nz; i++ ) {
11752d61bbb3SSatish Balay #if defined(USE_PETSC_COMPLEX)
1176e20fef11SSatish Balay       sum += PetscReal(PetscConj(*v)*(*v)); v++;
11772d61bbb3SSatish Balay #else
11782d61bbb3SSatish Balay       sum += (*v)*(*v); v++;
11792d61bbb3SSatish Balay #endif
11802d61bbb3SSatish Balay     }
11812d61bbb3SSatish Balay     *norm = sqrt(sum);
1182596552b5SBarry Smith   }  else if (type == NORM_INFINITY) { /* maximum row sum */
1183596552b5SBarry Smith     *norm = 0.0;
1184596552b5SBarry Smith     for ( k=0; k<bs; k++ ) {
1185596552b5SBarry Smith       for ( j=0; j<a->m; j++ ) {
1186596552b5SBarry Smith         v = a->a + bs2*a->i[j] + k;
1187596552b5SBarry Smith         sum = 0.0;
1188596552b5SBarry Smith         for ( i=0; i<a->i[j+1]-a->i[j]; i++ ) {
1189596552b5SBarry Smith           sum += PetscAbsScalar(*v);
1190596552b5SBarry Smith           v   += bs;
11912d61bbb3SSatish Balay         }
1192596552b5SBarry Smith         if (sum > *norm) *norm = sum;
1193596552b5SBarry Smith       }
1194596552b5SBarry Smith     }
1195596552b5SBarry Smith   } else {
11962d61bbb3SSatish Balay     SETERRQ(PETSC_ERR_SUP,0,"No support for this norm yet");
11972d61bbb3SSatish Balay   }
11982d61bbb3SSatish Balay   PetscFunctionReturn(0);
11992d61bbb3SSatish Balay }
12002d61bbb3SSatish Balay 
12012d61bbb3SSatish Balay 
12022d61bbb3SSatish Balay #undef __FUNC__
12032d61bbb3SSatish Balay #define __FUNC__ "MatEqual_SeqBAIJ"
12042d61bbb3SSatish Balay int MatEqual_SeqBAIJ(Mat A,Mat B, PetscTruth* flg)
12052d61bbb3SSatish Balay {
12062d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)B->data;
12072d61bbb3SSatish Balay 
12082d61bbb3SSatish Balay   PetscFunctionBegin;
12092d61bbb3SSatish Balay   if (B->type !=MATSEQBAIJ) SETERRQ(PETSC_ERR_ARG_INCOMP,0,"Matrices must be same type");
12102d61bbb3SSatish Balay 
12112d61bbb3SSatish Balay   /* If the  matrix/block dimensions are not equal, or no of nonzeros or shift */
12122d61bbb3SSatish Balay   if ((a->m != b->m) || (a->n !=b->n) || (a->bs != b->bs)|| (a->nz != b->nz)) {
12132d61bbb3SSatish Balay     *flg = PETSC_FALSE; PetscFunctionReturn(0);
12142d61bbb3SSatish Balay   }
12152d61bbb3SSatish Balay 
12162d61bbb3SSatish Balay   /* if the a->i are the same */
12172d61bbb3SSatish Balay   if (PetscMemcmp(a->i,b->i, (a->mbs+1)*sizeof(int))) {
12182d61bbb3SSatish Balay     *flg = PETSC_FALSE; PetscFunctionReturn(0);
12192d61bbb3SSatish Balay   }
12202d61bbb3SSatish Balay 
12212d61bbb3SSatish Balay   /* if a->j are the same */
12222d61bbb3SSatish Balay   if (PetscMemcmp(a->j,b->j,(a->nz)*sizeof(int))) {
12232d61bbb3SSatish Balay     *flg = PETSC_FALSE; PetscFunctionReturn(0);
12242d61bbb3SSatish Balay   }
12252d61bbb3SSatish Balay 
12262d61bbb3SSatish Balay   /* if a->a are the same */
12272d61bbb3SSatish Balay   if (PetscMemcmp(a->a, b->a,(a->nz)*(a->bs)*(a->bs)*sizeof(Scalar))) {
12282d61bbb3SSatish Balay     *flg = PETSC_FALSE; PetscFunctionReturn(0);
12292d61bbb3SSatish Balay   }
12302d61bbb3SSatish Balay   *flg = PETSC_TRUE;
12312d61bbb3SSatish Balay   PetscFunctionReturn(0);
12322d61bbb3SSatish Balay 
12332d61bbb3SSatish Balay }
12342d61bbb3SSatish Balay 
12352d61bbb3SSatish Balay #undef __FUNC__
12362d61bbb3SSatish Balay #define __FUNC__ "MatGetDiagonal_SeqBAIJ"
12372d61bbb3SSatish Balay int MatGetDiagonal_SeqBAIJ(Mat A,Vec v)
12382d61bbb3SSatish Balay {
12392d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
1240e1311b90SBarry Smith   int         ierr,i,j,k,n,row,bs,*ai,*aj,ambs,bs2;
12413f1db9ecSBarry Smith   Scalar      *x, zero = 0.0;
12423f1db9ecSBarry Smith   MatScalar   *aa,*aa_j;
12432d61bbb3SSatish Balay 
12442d61bbb3SSatish Balay   PetscFunctionBegin;
12452d61bbb3SSatish Balay   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix");
12462d61bbb3SSatish Balay   bs   = a->bs;
12472d61bbb3SSatish Balay   aa   = a->a;
12482d61bbb3SSatish Balay   ai   = a->i;
12492d61bbb3SSatish Balay   aj   = a->j;
12502d61bbb3SSatish Balay   ambs = a->mbs;
12512d61bbb3SSatish Balay   bs2  = a->bs2;
12522d61bbb3SSatish Balay 
1253e1311b90SBarry Smith   ierr = VecSet(&zero,v);CHKERRQ(ierr);
1254e1311b90SBarry Smith   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1255e1311b90SBarry Smith   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
12562d61bbb3SSatish Balay   if (n != a->m) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Nonconforming matrix and vector");
12572d61bbb3SSatish Balay   for ( i=0; i<ambs; i++ ) {
12582d61bbb3SSatish Balay     for ( j=ai[i]; j<ai[i+1]; j++ ) {
12592d61bbb3SSatish Balay       if (aj[j] == i) {
12602d61bbb3SSatish Balay         row  = i*bs;
12612d61bbb3SSatish Balay         aa_j = aa+j*bs2;
12622d61bbb3SSatish Balay         for (k=0; k<bs2; k+=(bs+1),row++) x[row] = aa_j[k];
12632d61bbb3SSatish Balay         break;
12642d61bbb3SSatish Balay       }
12652d61bbb3SSatish Balay     }
12662d61bbb3SSatish Balay   }
126715091d37SBarry Smith   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
12682d61bbb3SSatish Balay   PetscFunctionReturn(0);
12692d61bbb3SSatish Balay }
12702d61bbb3SSatish Balay 
12712d61bbb3SSatish Balay #undef __FUNC__
12722d61bbb3SSatish Balay #define __FUNC__ "MatDiagonalScale_SeqBAIJ"
12732d61bbb3SSatish Balay int MatDiagonalScale_SeqBAIJ(Mat A,Vec ll,Vec rr)
12742d61bbb3SSatish Balay {
12752d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
12763f1db9ecSBarry Smith   Scalar      *l,*r,x,*li,*ri;
12773f1db9ecSBarry Smith   MatScalar   *aa,*v;
1278e1311b90SBarry Smith   int         ierr,i,j,k,lm,rn,M,m,n,*ai,*aj,mbs,tmp,bs,bs2;
12792d61bbb3SSatish Balay 
12802d61bbb3SSatish Balay   PetscFunctionBegin;
12812d61bbb3SSatish Balay   ai  = a->i;
12822d61bbb3SSatish Balay   aj  = a->j;
12832d61bbb3SSatish Balay   aa  = a->a;
12842d61bbb3SSatish Balay   m   = a->m;
12852d61bbb3SSatish Balay   n   = a->n;
12862d61bbb3SSatish Balay   bs  = a->bs;
12872d61bbb3SSatish Balay   mbs = a->mbs;
12882d61bbb3SSatish Balay   bs2 = a->bs2;
12892d61bbb3SSatish Balay   if (ll) {
1290e1311b90SBarry Smith     ierr = VecGetArray(ll,&l); CHKERRQ(ierr);
1291e1311b90SBarry Smith     ierr = VecGetLocalSize(ll,&lm);CHKERRQ(ierr);
12922d61bbb3SSatish Balay     if (lm != m) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Left scaling vector wrong length");
12932d61bbb3SSatish Balay     for ( i=0; i<mbs; i++ ) { /* for each block row */
12942d61bbb3SSatish Balay       M  = ai[i+1] - ai[i];
12952d61bbb3SSatish Balay       li = l + i*bs;
12962d61bbb3SSatish Balay       v  = aa + bs2*ai[i];
12972d61bbb3SSatish Balay       for ( j=0; j<M; j++ ) { /* for each block */
12982d61bbb3SSatish Balay         for ( k=0; k<bs2; k++ ) {
12992d61bbb3SSatish Balay           (*v++) *= li[k%bs];
13002d61bbb3SSatish Balay         }
13012d61bbb3SSatish Balay       }
13022d61bbb3SSatish Balay     }
13032d61bbb3SSatish Balay   }
13042d61bbb3SSatish Balay 
13052d61bbb3SSatish Balay   if (rr) {
1306e1311b90SBarry Smith     ierr = VecGetArray(rr,&r); CHKERRQ(ierr);
1307e1311b90SBarry Smith     ierr = VecGetLocalSize(rr,&rn);CHKERRQ(ierr);
13082d61bbb3SSatish Balay     if (rn != n) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Right scaling vector wrong length");
13092d61bbb3SSatish Balay     for ( i=0; i<mbs; i++ ) { /* for each block row */
13102d61bbb3SSatish Balay       M  = ai[i+1] - ai[i];
13112d61bbb3SSatish Balay       v  = aa + bs2*ai[i];
13122d61bbb3SSatish Balay       for ( j=0; j<M; j++ ) { /* for each block */
13132d61bbb3SSatish Balay         ri = r + bs*aj[ai[i]+j];
13142d61bbb3SSatish Balay         for ( k=0; k<bs; k++ ) {
13152d61bbb3SSatish Balay           x = ri[k];
13162d61bbb3SSatish Balay           for ( tmp=0; tmp<bs; tmp++ ) (*v++) *= x;
13172d61bbb3SSatish Balay         }
13182d61bbb3SSatish Balay       }
13192d61bbb3SSatish Balay     }
13202d61bbb3SSatish Balay   }
13212d61bbb3SSatish Balay   PetscFunctionReturn(0);
13222d61bbb3SSatish Balay }
13232d61bbb3SSatish Balay 
13242d61bbb3SSatish Balay 
13252d61bbb3SSatish Balay #undef __FUNC__
13262d61bbb3SSatish Balay #define __FUNC__ "MatGetInfo_SeqBAIJ"
13272d61bbb3SSatish Balay int MatGetInfo_SeqBAIJ(Mat A,MatInfoType flag,MatInfo *info)
13282d61bbb3SSatish Balay {
13292d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
13302d61bbb3SSatish Balay 
13312d61bbb3SSatish Balay   PetscFunctionBegin;
13322d61bbb3SSatish Balay   info->rows_global    = (double)a->m;
13332d61bbb3SSatish Balay   info->columns_global = (double)a->n;
13342d61bbb3SSatish Balay   info->rows_local     = (double)a->m;
13352d61bbb3SSatish Balay   info->columns_local  = (double)a->n;
13362d61bbb3SSatish Balay   info->block_size     = a->bs2;
13372d61bbb3SSatish Balay   info->nz_allocated   = a->maxnz;
13382d61bbb3SSatish Balay   info->nz_used        = a->bs2*a->nz;
13392d61bbb3SSatish Balay   info->nz_unneeded    = (double)(info->nz_allocated - info->nz_used);
13402d61bbb3SSatish Balay   /*  if (info->nz_unneeded != A->info.nz_unneeded)
13412d61bbb3SSatish Balay     printf("space descrepancy: maxnz-nz = %d, nz_unneeded = %d\n",(int)info->nz_unneeded,(int)A->info.nz_unneeded); */
13422d61bbb3SSatish Balay   info->assemblies   = A->num_ass;
13432d61bbb3SSatish Balay   info->mallocs      = a->reallocs;
13442d61bbb3SSatish Balay   info->memory       = A->mem;
13452d61bbb3SSatish Balay   if (A->factor) {
13462d61bbb3SSatish Balay     info->fill_ratio_given  = A->info.fill_ratio_given;
13472d61bbb3SSatish Balay     info->fill_ratio_needed = A->info.fill_ratio_needed;
13482d61bbb3SSatish Balay     info->factor_mallocs    = A->info.factor_mallocs;
13492d61bbb3SSatish Balay   } else {
13502d61bbb3SSatish Balay     info->fill_ratio_given  = 0;
13512d61bbb3SSatish Balay     info->fill_ratio_needed = 0;
13522d61bbb3SSatish Balay     info->factor_mallocs    = 0;
13532d61bbb3SSatish Balay   }
13542d61bbb3SSatish Balay   PetscFunctionReturn(0);
13552d61bbb3SSatish Balay }
13562d61bbb3SSatish Balay 
13572d61bbb3SSatish Balay 
13582d61bbb3SSatish Balay #undef __FUNC__
13592d61bbb3SSatish Balay #define __FUNC__ "MatZeroEntries_SeqBAIJ"
13602d61bbb3SSatish Balay int MatZeroEntries_SeqBAIJ(Mat A)
13612d61bbb3SSatish Balay {
13622d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
13632d61bbb3SSatish Balay 
13642d61bbb3SSatish Balay   PetscFunctionBegin;
13653f1db9ecSBarry Smith   PetscMemzero(a->a,a->bs2*a->i[a->mbs]*sizeof(MatScalar));
13662d61bbb3SSatish Balay   PetscFunctionReturn(0);
13672d61bbb3SSatish Balay }
1368