xref: /petsc/src/mat/impls/baij/seq/baij2.c (revision 6a6a5d1dfdd5fcb6e31127fc1e7bbfb8f4c48c0a)
1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER
2*6a6a5d1dSBarry Smith static char vcid[] = "$Id: baij2.c,v 1.21 1997/12/05 20:01:54 balay Exp bsmith $";
3cac129eeSSatish Balay #endif
4cac129eeSSatish Balay 
52d61bbb3SSatish Balay #include "pinclude/pviewer.h"
62d61bbb3SSatish Balay #include "sys.h"
770f55243SBarry Smith #include "src/mat/impls/baij/seq/baij.h"
82d61bbb3SSatish Balay #include "src/vec/vecimpl.h"
92d61bbb3SSatish Balay #include "src/inline/spops.h"
10cac129eeSSatish Balay #include "petsc.h"
11736121d4SSatish Balay #include "src/inline/bitarray.h"
12cac129eeSSatish Balay 
132d61bbb3SSatish Balay 
145615d1e5SSatish Balay #undef __FUNC__
15d4bb536fSBarry Smith #define __FUNC__ "MatIncreaseOverlap_SeqBAIJ"
16736121d4SSatish Balay int MatIncreaseOverlap_SeqBAIJ(Mat A,int is_max,IS *is,int ov)
17a3192f15SSatish Balay {
18a3192f15SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
19218c64b6SSatish Balay   int         row, i,j,k,l,m,n, *idx,ierr, *nidx, isz, val, ival;
20218c64b6SSatish Balay   int         start, end, *ai, *aj,bs,*nidx2;
2142523b56SSatish Balay   BT          table;
22a3192f15SSatish Balay 
233a40ed3dSBarry Smith   PetscFunctionBegin;
24a3192f15SSatish Balay   m     = a->mbs;
25a3192f15SSatish Balay   ai    = a->i;
26a3192f15SSatish Balay   aj    = a->j;
27218c64b6SSatish Balay   bs    = a->bs;
28a3192f15SSatish Balay 
29a8c6a408SBarry Smith   if (ov < 0)  SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative overlap specified");
30a3192f15SSatish Balay 
3142523b56SSatish Balay   ierr  = BTCreate(m,table); CHKERRQ(ierr);
32a3192f15SSatish Balay   nidx  = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(nidx);
33218c64b6SSatish Balay   nidx2 = (int *)PetscMalloc((a->m+1)*sizeof(int)); CHKPTRQ(nidx2);
34a3192f15SSatish Balay 
35a3192f15SSatish Balay   for ( i=0; i<is_max; i++ ) {
36a3192f15SSatish Balay     /* Initialise the two local arrays */
37a3192f15SSatish Balay     isz  = 0;
3842523b56SSatish Balay     BTMemzero(m,table);
39a3192f15SSatish Balay 
40a3192f15SSatish Balay     /* Extract the indices, assume there can be duplicate entries */
41a3192f15SSatish Balay     ierr = ISGetIndices(is[i],&idx);  CHKERRQ(ierr);
42a3192f15SSatish Balay     ierr = ISGetSize(is[i],&n);  CHKERRQ(ierr);
43a3192f15SSatish Balay 
44a3192f15SSatish Balay     /* Enter these into the temp arrays i.e mark table[row], enter row into new index */
45a3192f15SSatish Balay     for ( j=0; j<n ; ++j){
46218c64b6SSatish Balay       ival = idx[j]/bs; /* convert the indices into block indices */
47a8c6a408SBarry Smith       if (ival>m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"index greater than mat-dim");
4842523b56SSatish Balay       if(!BTLookupSet(table, ival)) { nidx[isz++] = ival;}
49a3192f15SSatish Balay     }
50a3192f15SSatish Balay     ierr = ISRestoreIndices(is[i],&idx);  CHKERRQ(ierr);
51a3192f15SSatish Balay     ierr = ISDestroy(is[i]); CHKERRQ(ierr);
52a3192f15SSatish Balay 
53a3192f15SSatish Balay     k = 0;
54a3192f15SSatish Balay     for ( j=0; j<ov; j++){ /* for each overlap*/
55a3192f15SSatish Balay       n = isz;
56a3192f15SSatish Balay       for ( ; k<n ; k++){ /* do only those rows in nidx[k], which are not done yet */
57a3192f15SSatish Balay         row   = nidx[k];
58a3192f15SSatish Balay         start = ai[row];
59a3192f15SSatish Balay         end   = ai[row+1];
60a3192f15SSatish Balay         for ( l = start; l<end ; l++){
61a3192f15SSatish Balay           val = aj[l];
6242523b56SSatish Balay           if (!BTLookupSet(table,val)) {nidx[isz++] = val;}
63a3192f15SSatish Balay         }
64a3192f15SSatish Balay       }
65a3192f15SSatish Balay     }
66218c64b6SSatish Balay     /* expand the Index Set */
67218c64b6SSatish Balay     for (j=0; j<isz; j++ ) {
68218c64b6SSatish Balay       for (k=0; k<bs; k++ )
69218c64b6SSatish Balay         nidx2[j*bs+k] = nidx[j]*bs+k;
70218c64b6SSatish Balay     }
71029af93fSBarry Smith     ierr = ISCreateGeneral(PETSC_COMM_SELF, isz*bs, nidx2, (is+i)); CHKERRQ(ierr);
72a3192f15SSatish Balay   }
7342523b56SSatish Balay   BTDestroy(table);
74a3192f15SSatish Balay   PetscFree(nidx);
75218c64b6SSatish Balay   PetscFree(nidx2);
763a40ed3dSBarry Smith   PetscFunctionReturn(0);
77a3192f15SSatish Balay }
781c351548SSatish Balay 
795615d1e5SSatish Balay #undef __FUNC__
805615d1e5SSatish Balay #define __FUNC__ "MatGetSubMatrix_SeqBAIJ_Private"
81*6a6a5d1dSBarry Smith int MatGetSubMatrix_SeqBAIJ_Private(Mat A,IS isrow,IS iscol,int cs,MatGetSubMatrixCall scall,Mat *B)
82736121d4SSatish Balay {
83736121d4SSatish Balay   Mat_SeqBAIJ  *a = (Mat_SeqBAIJ *) A->data,*c;
84736121d4SSatish Balay   int          nznew, *smap, i, k, kstart, kend, ierr, oldcols = a->nbs,*lens;
85218c64b6SSatish Balay   int          row,mat_i,*mat_j,tcol,*mat_ilen;
86736121d4SSatish Balay   int          *irow, *icol, nrows, ncols,*ssmap,bs=a->bs, bs2=a->bs2;
87218c64b6SSatish Balay   int          *aj = a->j, *ai = a->i;
88218c64b6SSatish Balay   Scalar       *mat_a;
89736121d4SSatish Balay   Mat          C;
90736121d4SSatish Balay 
913a40ed3dSBarry Smith   PetscFunctionBegin;
92736121d4SSatish Balay   ierr = ISSorted(iscol,(PetscTruth*)&i);
93a8c6a408SBarry Smith   if (!i) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"IS is not sorted");
94736121d4SSatish Balay 
95736121d4SSatish Balay   ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr);
96218c64b6SSatish Balay   ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr);
97736121d4SSatish Balay   ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr);
98736121d4SSatish Balay   ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr);
99736121d4SSatish Balay 
100736121d4SSatish Balay   smap  = (int *) PetscMalloc((1+oldcols)*sizeof(int)); CHKPTRQ(smap);
101736121d4SSatish Balay   ssmap = smap;
102736121d4SSatish Balay   lens  = (int *) PetscMalloc((1+nrows)*sizeof(int)); CHKPTRQ(lens);
103736121d4SSatish Balay   PetscMemzero(smap,oldcols*sizeof(int));
104736121d4SSatish Balay   for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1;
105736121d4SSatish Balay   /* determine lens of each row */
106736121d4SSatish Balay   for (i=0; i<nrows; i++) {
107736121d4SSatish Balay     kstart  = ai[irow[i]];
108736121d4SSatish Balay     kend    = kstart + a->ilen[irow[i]];
109736121d4SSatish Balay     lens[i] = 0;
110736121d4SSatish Balay       for ( k=kstart; k<kend; k++ ) {
111736121d4SSatish Balay         if (ssmap[aj[k]]) {
112736121d4SSatish Balay           lens[i]++;
113736121d4SSatish Balay         }
114736121d4SSatish Balay       }
115736121d4SSatish Balay     }
116736121d4SSatish Balay   /* Create and fill new matrix */
117736121d4SSatish Balay   if (scall == MAT_REUSE_MATRIX) {
118736121d4SSatish Balay     c = (Mat_SeqBAIJ *)((*B)->data);
119736121d4SSatish Balay 
120a8c6a408SBarry Smith     if (c->mbs!=nrows || c->nbs!=ncols || c->bs!=bs) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Submatrix wrong size");
121736121d4SSatish Balay     if (PetscMemcmp(c->ilen,lens, c->mbs *sizeof(int))) {
122a8c6a408SBarry Smith       SETERRQ(PETSC_ERR_ARG_SIZ,0,"Cannot reuse matrix. wrong no of nonzeros");
123736121d4SSatish Balay     }
124736121d4SSatish Balay     PetscMemzero(c->ilen,c->mbs*sizeof(int));
125736121d4SSatish Balay     C = *B;
1263a40ed3dSBarry Smith   } else {
127736121d4SSatish Balay     ierr = MatCreateSeqBAIJ(A->comm,bs,nrows*bs,ncols*bs,0,lens,&C);CHKERRQ(ierr);
128736121d4SSatish Balay   }
129736121d4SSatish Balay   c = (Mat_SeqBAIJ *)(C->data);
130736121d4SSatish Balay   for (i=0; i<nrows; i++) {
131736121d4SSatish Balay     row    = irow[i];
132736121d4SSatish Balay     nznew  = 0;
133736121d4SSatish Balay     kstart = ai[row];
134736121d4SSatish Balay     kend   = kstart + a->ilen[row];
135736121d4SSatish Balay     mat_i  = c->i[i];
136736121d4SSatish Balay     mat_j  = c->j + mat_i;
137218c64b6SSatish Balay     mat_a  = c->a + mat_i*bs2;
138736121d4SSatish Balay     mat_ilen = c->ilen + i;
139736121d4SSatish Balay     for ( k=kstart; k<kend; k++ ) {
140736121d4SSatish Balay       if ((tcol=ssmap[a->j[k]])) {
141736121d4SSatish Balay         *mat_j++ = tcol - 1;
142736121d4SSatish Balay         PetscMemcpy(mat_a,a->a+k*bs2,bs2*sizeof(Scalar)); mat_a+=bs2;
143736121d4SSatish Balay         (*mat_ilen)++;
144736121d4SSatish Balay 
145736121d4SSatish Balay       }
146736121d4SSatish Balay     }
147736121d4SSatish Balay   }
148218c64b6SSatish Balay 
149736121d4SSatish Balay   /* Free work space */
150736121d4SSatish Balay   ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr);
151736121d4SSatish Balay   PetscFree(smap); PetscFree(lens);
1526d4a8577SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
1536d4a8577SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr);
154736121d4SSatish Balay 
155736121d4SSatish Balay   ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr);
156736121d4SSatish Balay   *B = C;
1573a40ed3dSBarry Smith   PetscFunctionReturn(0);
158736121d4SSatish Balay }
159736121d4SSatish Balay 
1605615d1e5SSatish Balay #undef __FUNC__
1615615d1e5SSatish Balay #define __FUNC__ "MatGetSubMatrix_SeqBAIJ"
162*6a6a5d1dSBarry Smith int MatGetSubMatrix_SeqBAIJ(Mat A,IS isrow,IS iscol,int cs,MatGetSubMatrixCall scall,Mat *B)
163218c64b6SSatish Balay {
164218c64b6SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
165218c64b6SSatish Balay   IS          is1,is2;
166218c64b6SSatish Balay   int         *vary,*iary,*irow,*icol,nrows,ncols,i,ierr,bs=a->bs,count;
167218c64b6SSatish Balay 
1683a40ed3dSBarry Smith   PetscFunctionBegin;
169218c64b6SSatish Balay   ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr);
170218c64b6SSatish Balay   ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr);
171218c64b6SSatish Balay   ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr);
172218c64b6SSatish Balay   ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr);
173218c64b6SSatish Balay 
174218c64b6SSatish Balay   /* Verify if the indices corespond to each element in a block
175218c64b6SSatish Balay    and form the IS with compressed IS */
176218c64b6SSatish Balay   vary = (int *) PetscMalloc(2*(a->mbs+1)*sizeof(int)); CHKPTRQ(vary);
177218c64b6SSatish Balay   iary = vary + a->mbs;
178218c64b6SSatish Balay   PetscMemzero(vary,(a->mbs)*sizeof(int));
179218c64b6SSatish Balay   for ( i=0; i<nrows; i++) vary[irow[i]/bs]++;
180218c64b6SSatish Balay   count = 0;
181218c64b6SSatish Balay   for (i=0; i<a->mbs; i++) {
182*6a6a5d1dSBarry Smith     if (vary[i]!=0 && vary[i]!=bs) SETERRA(1,0,"Index set does not match blocks");
183218c64b6SSatish Balay     if (vary[i]==bs) iary[count++] = i;
184218c64b6SSatish Balay   }
185029af93fSBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF, count, iary,&is1); CHKERRQ(ierr);
186218c64b6SSatish Balay 
187218c64b6SSatish Balay   PetscMemzero(vary,(a->mbs)*sizeof(int));
188218c64b6SSatish Balay   for ( i=0; i<ncols; i++) vary[icol[i]/bs]++;
189218c64b6SSatish Balay   count = 0;
190218c64b6SSatish Balay   for (i=0; i<a->mbs; i++) {
191e3372554SBarry Smith     if (vary[i]!=0 && vary[i]!=bs) SETERRA(1,0,"MatGetSubmatrices_SeqBAIJ:");
192218c64b6SSatish Balay     if (vary[i]==bs) iary[count++] = i;
193218c64b6SSatish Balay   }
194029af93fSBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF, count, iary,&is2); CHKERRQ(ierr);
195218c64b6SSatish Balay   ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr);
196218c64b6SSatish Balay   ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr);
197218c64b6SSatish Balay   PetscFree(vary);
198218c64b6SSatish Balay 
199*6a6a5d1dSBarry Smith   ierr = MatGetSubMatrix_SeqBAIJ_Private(A,is1,is2,cs,scall,B); CHKERRQ(ierr);
200218c64b6SSatish Balay   ISDestroy(is1);
201218c64b6SSatish Balay   ISDestroy(is2);
2023a40ed3dSBarry Smith   PetscFunctionReturn(0);
203218c64b6SSatish Balay }
204218c64b6SSatish Balay 
2055615d1e5SSatish Balay #undef __FUNC__
2065615d1e5SSatish Balay #define __FUNC__ "MatGetSubMatrices_SeqBAIJ"
207*6a6a5d1dSBarry Smith int MatGetSubMatrices_SeqBAIJ(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall,Mat **B)
208736121d4SSatish Balay {
209736121d4SSatish Balay   int ierr,i;
210736121d4SSatish Balay 
2113a40ed3dSBarry Smith   PetscFunctionBegin;
212736121d4SSatish Balay   if (scall == MAT_INITIAL_MATRIX) {
213736121d4SSatish Balay     *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B);
214736121d4SSatish Balay   }
215736121d4SSatish Balay 
216736121d4SSatish Balay   for ( i=0; i<n; i++ ) {
217*6a6a5d1dSBarry Smith     ierr = MatGetSubMatrix_SeqBAIJ(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr);
218736121d4SSatish Balay   }
2193a40ed3dSBarry Smith   PetscFunctionReturn(0);
220736121d4SSatish Balay }
221218c64b6SSatish Balay 
222218c64b6SSatish Balay 
2232d61bbb3SSatish Balay /* -------------------------------------------------------*/
2242d61bbb3SSatish Balay /* Should check that shapes of vectors and matrices match */
2252d61bbb3SSatish Balay /* -------------------------------------------------------*/
2262d61bbb3SSatish Balay #include "pinclude/plapack.h"
2272d61bbb3SSatish Balay 
2282d61bbb3SSatish Balay #undef __FUNC__
2292d61bbb3SSatish Balay #define __FUNC__ "MatMult_SeqBAIJ_1"
2302d61bbb3SSatish Balay int MatMult_SeqBAIJ_1(Mat A,Vec xx,Vec zz)
2312d61bbb3SSatish Balay {
2322d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
2332d61bbb3SSatish Balay   register Scalar *x,*z,*v,sum;
2342d61bbb3SSatish Balay   int             mbs=a->mbs,i,*idx,*ii,n;
2352d61bbb3SSatish Balay 
2362d61bbb3SSatish Balay   PetscFunctionBegin;
2372d61bbb3SSatish Balay   VecGetArray_Fast(xx,x);
2382d61bbb3SSatish Balay   VecGetArray_Fast(zz,z);
2392d61bbb3SSatish Balay 
2402d61bbb3SSatish Balay   idx   = a->j;
2412d61bbb3SSatish Balay   v     = a->a;
2422d61bbb3SSatish Balay   ii    = a->i;
2432d61bbb3SSatish Balay 
2442d61bbb3SSatish Balay   for ( i=0; i<mbs; i++ ) {
2452d61bbb3SSatish Balay     n    = ii[1] - ii[0]; ii++;
2462d61bbb3SSatish Balay     sum  = 0.0;
2472d61bbb3SSatish Balay     while (n--) sum += *v++ * x[*idx++];
2482d61bbb3SSatish Balay     z[i] = sum;
2492d61bbb3SSatish Balay   }
2502d61bbb3SSatish Balay   VecRestoreArray_Fast(xx,x);
2512d61bbb3SSatish Balay   VecRestoreArray_Fast(zz,z);
2522d61bbb3SSatish Balay   PLogFlops(2*a->nz - a->m);
2532d61bbb3SSatish Balay   PetscFunctionReturn(0);
2542d61bbb3SSatish Balay }
2552d61bbb3SSatish Balay 
2562d61bbb3SSatish Balay #undef __FUNC__
2572d61bbb3SSatish Balay #define __FUNC__ "MatMult_SeqBAIJ_2"
2582d61bbb3SSatish Balay int MatMult_SeqBAIJ_2(Mat A,Vec xx,Vec zz)
2592d61bbb3SSatish Balay {
2602d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
2612d61bbb3SSatish Balay   register Scalar *x,*z,*v,*xb,sum1,sum2;
2622d61bbb3SSatish Balay   register Scalar x1,x2;
2632d61bbb3SSatish Balay   int             mbs=a->mbs,i,*idx,*ii,j,n;
2642d61bbb3SSatish Balay 
2652d61bbb3SSatish Balay   PetscFunctionBegin;
2662d61bbb3SSatish Balay   VecGetArray_Fast(xx,x);
2672d61bbb3SSatish Balay   VecGetArray_Fast(zz,z);
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   }
2852d61bbb3SSatish Balay   VecRestoreArray_Fast(xx,x);
2862d61bbb3SSatish Balay   VecRestoreArray_Fast(zz,z);
2872d61bbb3SSatish Balay   PLogFlops(4*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;
2962d61bbb3SSatish Balay   register Scalar *x,*z,*v,*xb,sum1,sum2,sum3,x1,x2,x3;
2972d61bbb3SSatish Balay   int             mbs=a->mbs,i,*idx,*ii,j,n;
2982d61bbb3SSatish Balay 
2992d61bbb3SSatish Balay   PetscFunctionBegin;
3002d61bbb3SSatish Balay   VecGetArray_Fast(xx,x);
3012d61bbb3SSatish Balay   VecGetArray_Fast(zz,z);
3022d61bbb3SSatish Balay 
3032d61bbb3SSatish Balay   idx   = a->j;
3042d61bbb3SSatish Balay   v     = a->a;
3052d61bbb3SSatish Balay   ii    = a->i;
3062d61bbb3SSatish Balay 
3072d61bbb3SSatish Balay   for ( i=0; i<mbs; i++ ) {
3082d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
3092d61bbb3SSatish Balay     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0;
3102d61bbb3SSatish Balay     for ( j=0; j<n; j++ ) {
3112d61bbb3SSatish Balay       xb = x + 3*(*idx++); x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
3122d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[3]*x2 + v[6]*x3;
3132d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[4]*x2 + v[7]*x3;
3142d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[5]*x2 + v[8]*x3;
3152d61bbb3SSatish Balay       v += 9;
3162d61bbb3SSatish Balay     }
3172d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3;
3182d61bbb3SSatish Balay     z += 3;
3192d61bbb3SSatish Balay   }
3202d61bbb3SSatish Balay   VecRestoreArray_Fast(xx,x);
3212d61bbb3SSatish Balay   VecRestoreArray_Fast(zz,z);
3222d61bbb3SSatish Balay   PLogFlops(18*a->nz - a->m);
3232d61bbb3SSatish Balay   PetscFunctionReturn(0);
3242d61bbb3SSatish Balay }
3252d61bbb3SSatish Balay 
3262d61bbb3SSatish Balay #undef __FUNC__
3272d61bbb3SSatish Balay #define __FUNC__ "MatMult_SeqBAIJ_4"
3282d61bbb3SSatish Balay int MatMult_SeqBAIJ_4(Mat A,Vec xx,Vec zz)
3292d61bbb3SSatish Balay {
3302d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
3312d61bbb3SSatish Balay   register Scalar *x,*z,*v,*xb,sum1,sum2,sum3,sum4,x1,x2,x3,x4;
3322d61bbb3SSatish Balay   int             mbs=a->mbs,i,*idx,*ii,j,n;
3332d61bbb3SSatish Balay 
3342d61bbb3SSatish Balay   PetscFunctionBegin;
3352d61bbb3SSatish Balay   VecGetArray_Fast(xx,x);
3362d61bbb3SSatish Balay   VecGetArray_Fast(zz,z);
3372d61bbb3SSatish Balay 
3382d61bbb3SSatish Balay   idx   = a->j;
3392d61bbb3SSatish Balay   v     = a->a;
3402d61bbb3SSatish Balay   ii    = a->i;
3412d61bbb3SSatish Balay 
3422d61bbb3SSatish Balay   for ( i=0; i<mbs; i++ ) {
3432d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
3442d61bbb3SSatish Balay     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0;
3452d61bbb3SSatish Balay     for ( j=0; j<n; j++ ) {
3462d61bbb3SSatish Balay       xb = x + 4*(*idx++);
3472d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
3482d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[4]*x2 + v[8]*x3   + v[12]*x4;
3492d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[5]*x2 + v[9]*x3   + v[13]*x4;
3502d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[6]*x2 + v[10]*x3  + v[14]*x4;
3512d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[7]*x2 + v[11]*x3  + v[15]*x4;
3522d61bbb3SSatish Balay       v += 16;
3532d61bbb3SSatish Balay     }
3542d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4;
3552d61bbb3SSatish Balay     z += 4;
3562d61bbb3SSatish Balay   }
3572d61bbb3SSatish Balay   VecRestoreArray_Fast(xx,x);
3582d61bbb3SSatish Balay   VecRestoreArray_Fast(zz,z);
3592d61bbb3SSatish Balay   PLogFlops(32*a->nz - a->m);
3602d61bbb3SSatish Balay   PetscFunctionReturn(0);
3612d61bbb3SSatish Balay }
3622d61bbb3SSatish Balay 
3632d61bbb3SSatish Balay #undef __FUNC__
3642d61bbb3SSatish Balay #define __FUNC__ "MatMult_SeqBAIJ_5"
3652d61bbb3SSatish Balay int MatMult_SeqBAIJ_5(Mat A,Vec xx,Vec zz)
3662d61bbb3SSatish Balay {
3672d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
3682d61bbb3SSatish Balay   register Scalar sum1,sum2,sum3,sum4,sum5,x1,x2,x3,x4,x5;
3692d61bbb3SSatish Balay   register Scalar * restrict v,* restrict xb,* restrict z, * restrict x;
3702d61bbb3SSatish Balay   int             mbs=a->mbs,i,*idx,*ii,j,n;
3712d61bbb3SSatish Balay 
3722d61bbb3SSatish Balay   VecGetArray_Fast(xx,x);
3732d61bbb3SSatish Balay   VecGetArray_Fast(zz,z);
3742d61bbb3SSatish Balay 
3752d61bbb3SSatish Balay   idx   = a->j;
3762d61bbb3SSatish Balay   v     = a->a;
3772d61bbb3SSatish Balay   ii    = a->i;
3782d61bbb3SSatish Balay 
3792d61bbb3SSatish Balay   for ( i=0; i<mbs; i++ ) {
3802d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
3812d61bbb3SSatish Balay     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0;
3822d61bbb3SSatish Balay     for ( j=0; j<n; j++ ) {
3832d61bbb3SSatish Balay       xb = x + 5*(*idx++);
3842d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4];
3852d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[5]*x2 + v[10]*x3  + v[15]*x4 + v[20]*x5;
3862d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[6]*x2 + v[11]*x3  + v[16]*x4 + v[21]*x5;
3872d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[7]*x2 + v[12]*x3  + v[17]*x4 + v[22]*x5;
3882d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[8]*x2 + v[13]*x3  + v[18]*x4 + v[23]*x5;
3892d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[9]*x2 + v[14]*x3  + v[19]*x4 + v[24]*x5;
3902d61bbb3SSatish Balay       v += 25;
3912d61bbb3SSatish Balay     }
3922d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5;
3932d61bbb3SSatish Balay     z += 5;
3942d61bbb3SSatish Balay   }
3952d61bbb3SSatish Balay   VecRestoreArray_Fast(xx,x);
3962d61bbb3SSatish Balay   VecRestoreArray_Fast(zz,z);
3972d61bbb3SSatish Balay   PLogFlops(50*a->nz - a->m);
3982d61bbb3SSatish Balay   PetscFunctionReturn(0);
3992d61bbb3SSatish Balay }
4002d61bbb3SSatish Balay 
4012d61bbb3SSatish Balay #undef __FUNC__
4022d61bbb3SSatish Balay #define __FUNC__ "MatMult_SeqBAIJ_7"
4032d61bbb3SSatish Balay int MatMult_SeqBAIJ_7(Mat A,Vec xx,Vec zz)
4042d61bbb3SSatish Balay {
4052d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
4062d61bbb3SSatish Balay   register Scalar *x,*z,*v,*xb,sum1,sum2,sum3,sum4,sum5,sum6,sum7;
4072d61bbb3SSatish Balay   register Scalar x1,x2,x3,x4,x5,x6,x7;
4082d61bbb3SSatish Balay   int             mbs=a->mbs,i,*idx,*ii,j,n;
4092d61bbb3SSatish Balay 
4102d61bbb3SSatish Balay   VecGetArray_Fast(xx,x);
4112d61bbb3SSatish Balay   VecGetArray_Fast(zz,z);
4122d61bbb3SSatish Balay 
4132d61bbb3SSatish Balay   idx   = a->j;
4142d61bbb3SSatish Balay   v     = a->a;
4152d61bbb3SSatish Balay   ii    = a->i;
4162d61bbb3SSatish Balay 
4172d61bbb3SSatish Balay   for ( i=0; i<mbs; i++ ) {
4182d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
4192d61bbb3SSatish Balay     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0; sum7 = 0.0;
4202d61bbb3SSatish Balay     for ( j=0; j<n; j++ ) {
4212d61bbb3SSatish Balay       xb = x + 7*(*idx++);
4222d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6];
4232d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[7]*x2  + v[14]*x3  + v[21]*x4 + v[28]*x5 + v[35]*x6 + v[42]*x7;
4242d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[8]*x2  + v[15]*x3  + v[22]*x4 + v[29]*x5 + v[36]*x6 + v[43]*x7;
4252d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[9]*x2  + v[16]*x3  + v[23]*x4 + v[30]*x5 + v[37]*x6 + v[44]*x7;
4262d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[10]*x2 + v[17]*x3  + v[24]*x4 + v[31]*x5 + v[38]*x6 + v[45]*x7;
4272d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[11]*x2 + v[18]*x3  + v[25]*x4 + v[32]*x5 + v[39]*x6 + v[46]*x7;
4282d61bbb3SSatish Balay       sum6 += v[5]*x1 + v[12]*x2 + v[19]*x3  + v[26]*x4 + v[33]*x5 + v[40]*x6 + v[47]*x7;
4292d61bbb3SSatish Balay       sum7 += v[6]*x1 + v[13]*x2 + v[20]*x3  + v[27]*x4 + v[34]*x5 + v[41]*x6 + v[48]*x7;
4302d61bbb3SSatish Balay       v += 49;
4312d61bbb3SSatish Balay     }
4322d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
4332d61bbb3SSatish Balay     z += 7;
4342d61bbb3SSatish Balay   }
4352d61bbb3SSatish Balay 
4362d61bbb3SSatish Balay   VecRestoreArray_Fast(xx,x);
4372d61bbb3SSatish Balay   VecRestoreArray_Fast(zz,z);
4382d61bbb3SSatish Balay   PLogFlops(98*a->nz - a->m);
4392d61bbb3SSatish Balay   PetscFunctionReturn(0);
4402d61bbb3SSatish Balay }
4412d61bbb3SSatish Balay 
4422d61bbb3SSatish Balay #undef __FUNC__
4432d61bbb3SSatish Balay #define __FUNC__ "MatMult_SeqBAIJ_N"
4442d61bbb3SSatish Balay int MatMult_SeqBAIJ_N(Mat A,Vec xx,Vec zz)
4452d61bbb3SSatish Balay {
4462d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
4472d61bbb3SSatish Balay   register Scalar *x,*z,*v,*xb;
4482d61bbb3SSatish Balay   Scalar          _DOne = 1.0,*work,*workt,_DZero = 0.0;
4492d61bbb3SSatish Balay   int             mbs=a->mbs,i,*idx,*ii,bs=a->bs,j,n,bs2=a->bs2;
4502d61bbb3SSatish Balay   int             _One = 1,ncols,k;
4512d61bbb3SSatish Balay 
4522d61bbb3SSatish Balay   PetscFunctionBegin;
4532d61bbb3SSatish Balay   VecGetArray_Fast(xx,x);
4542d61bbb3SSatish Balay   VecGetArray_Fast(zz,z);
4552d61bbb3SSatish Balay 
4562d61bbb3SSatish Balay   idx   = a->j;
4572d61bbb3SSatish Balay   v     = a->a;
4582d61bbb3SSatish Balay   ii    = a->i;
459218c64b6SSatish Balay 
460218c64b6SSatish Balay 
4612d61bbb3SSatish Balay   if (!a->mult_work) {
4622d61bbb3SSatish Balay     k = PetscMax(a->m,a->n);
4632d61bbb3SSatish Balay     a->mult_work = (Scalar *) PetscMalloc((k+1)*sizeof(Scalar));CHKPTRQ(a->mult_work);
4642d61bbb3SSatish Balay   }
4652d61bbb3SSatish Balay   work = a->mult_work;
4662d61bbb3SSatish Balay   for ( i=0; i<mbs; i++ ) {
4672d61bbb3SSatish Balay     n     = ii[1] - ii[0]; ii++;
4682d61bbb3SSatish Balay     ncols = n*bs;
4692d61bbb3SSatish Balay     workt = work;
4702d61bbb3SSatish Balay     for ( j=0; j<n; j++ ) {
4712d61bbb3SSatish Balay       xb = x + bs*(*idx++);
4722d61bbb3SSatish Balay       for ( k=0; k<bs; k++ ) workt[k] = xb[k];
4732d61bbb3SSatish Balay       workt += bs;
4742d61bbb3SSatish Balay     }
4752d61bbb3SSatish Balay     LAgemv_("N",&bs,&ncols,&_DOne,v,&bs,work,&_One,&_DZero,z,&_One);
4762d61bbb3SSatish Balay     v += n*bs2;
4772d61bbb3SSatish Balay     z += bs;
4782d61bbb3SSatish Balay   }
4792d61bbb3SSatish Balay   VecRestoreArray_Fast(xx,x);
4802d61bbb3SSatish Balay   VecRestoreArray_Fast(zz,z);
4812d61bbb3SSatish Balay   PLogFlops(2*a->nz*bs2 - a->m);
4822d61bbb3SSatish Balay   PetscFunctionReturn(0);
4832d61bbb3SSatish Balay }
4842d61bbb3SSatish Balay 
4852d61bbb3SSatish Balay #undef __FUNC__
4862d61bbb3SSatish Balay #define __FUNC__ "MatMultAdd_SeqBAIJ_1"
4872d61bbb3SSatish Balay int MatMultAdd_SeqBAIJ_1(Mat A,Vec xx,Vec yy,Vec zz)
4882d61bbb3SSatish Balay {
4892d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
4902d61bbb3SSatish Balay   register Scalar *x,*y,*z,*v,sum;
4912d61bbb3SSatish Balay   int             mbs=a->mbs,i,*idx,*ii,n;
4922d61bbb3SSatish Balay 
4932d61bbb3SSatish Balay   PetscFunctionBegin;
4942d61bbb3SSatish Balay   VecGetArray_Fast(xx,x);
4952d61bbb3SSatish Balay   VecGetArray_Fast(yy,y);
4962d61bbb3SSatish Balay   VecGetArray_Fast(zz,z);
4972d61bbb3SSatish Balay 
4982d61bbb3SSatish Balay   idx   = a->j;
4992d61bbb3SSatish Balay   v     = a->a;
5002d61bbb3SSatish Balay   ii    = a->i;
5012d61bbb3SSatish Balay 
5022d61bbb3SSatish Balay   for ( i=0; i<mbs; i++ ) {
5032d61bbb3SSatish Balay     n    = ii[1] - ii[0]; ii++;
5042d61bbb3SSatish Balay     sum  = y[i];
5052d61bbb3SSatish Balay     while (n--) sum += *v++ * x[*idx++];
5062d61bbb3SSatish Balay     z[i] = sum;
5072d61bbb3SSatish Balay   }
5082d61bbb3SSatish Balay   VecRestoreArray_Fast(xx,x);
5092d61bbb3SSatish Balay   VecRestoreArray_Fast(yy,y);
5102d61bbb3SSatish Balay   VecRestoreArray_Fast(zz,z);
5112d61bbb3SSatish Balay   PLogFlops(2*a->nz);
5122d61bbb3SSatish Balay   PetscFunctionReturn(0);
5132d61bbb3SSatish Balay }
5142d61bbb3SSatish Balay 
5152d61bbb3SSatish Balay #undef __FUNC__
5162d61bbb3SSatish Balay #define __FUNC__ "MatMultAdd_SeqBAIJ_2"
5172d61bbb3SSatish Balay int MatMultAdd_SeqBAIJ_2(Mat A,Vec xx,Vec yy,Vec zz)
5182d61bbb3SSatish Balay {
5192d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
5202d61bbb3SSatish Balay   register Scalar *x,*y,*z,*v,*xb,sum1,sum2;
5212d61bbb3SSatish Balay   register Scalar x1,x2;
5222d61bbb3SSatish Balay   int             mbs=a->mbs,i,*idx,*ii,j,n;
5232d61bbb3SSatish Balay 
5242d61bbb3SSatish Balay   PetscFunctionBegin;
5252d61bbb3SSatish Balay   VecGetArray_Fast(xx,x);
5262d61bbb3SSatish Balay   VecGetArray_Fast(yy,y);
5272d61bbb3SSatish Balay   VecGetArray_Fast(zz,z);
5282d61bbb3SSatish Balay 
5292d61bbb3SSatish Balay   idx   = a->j;
5302d61bbb3SSatish Balay   v     = a->a;
5312d61bbb3SSatish Balay   ii    = a->i;
5322d61bbb3SSatish Balay 
5332d61bbb3SSatish Balay   for ( i=0; i<mbs; i++ ) {
5342d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
5352d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1];
5362d61bbb3SSatish Balay     for ( j=0; j<n; j++ ) {
5372d61bbb3SSatish Balay       xb = x + 2*(*idx++); x1 = xb[0]; x2 = xb[1];
5382d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[2]*x2;
5392d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[3]*x2;
5402d61bbb3SSatish Balay       v += 4;
5412d61bbb3SSatish Balay     }
5422d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2;
5432d61bbb3SSatish Balay     z += 2; y += 2;
5442d61bbb3SSatish Balay   }
5452d61bbb3SSatish Balay   VecRestoreArray_Fast(xx,x);
5462d61bbb3SSatish Balay   VecRestoreArray_Fast(yy,y);
5472d61bbb3SSatish Balay   VecRestoreArray_Fast(zz,z);
5482d61bbb3SSatish Balay   PLogFlops(4*a->nz);
5492d61bbb3SSatish Balay   PetscFunctionReturn(0);
5502d61bbb3SSatish Balay }
5512d61bbb3SSatish Balay 
5522d61bbb3SSatish Balay #undef __FUNC__
5532d61bbb3SSatish Balay #define __FUNC__ "MatMultAdd_SeqBAIJ_3"
5542d61bbb3SSatish Balay int MatMultAdd_SeqBAIJ_3(Mat A,Vec xx,Vec yy,Vec zz)
5552d61bbb3SSatish Balay {
5562d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
5572d61bbb3SSatish Balay   register Scalar *x,*y,*z,*v,*xb,sum1,sum2,sum3,x1,x2,x3;
5582d61bbb3SSatish Balay   int             mbs=a->mbs,i,*idx,*ii,j,n;
5592d61bbb3SSatish Balay 
5602d61bbb3SSatish Balay   PetscFunctionBegin;
5612d61bbb3SSatish Balay   VecGetArray_Fast(xx,x);
5622d61bbb3SSatish Balay   VecGetArray_Fast(yy,y);
5632d61bbb3SSatish Balay   VecGetArray_Fast(zz,z);
5642d61bbb3SSatish Balay 
5652d61bbb3SSatish Balay   idx   = a->j;
5662d61bbb3SSatish Balay   v     = a->a;
5672d61bbb3SSatish Balay   ii    = a->i;
5682d61bbb3SSatish Balay 
5692d61bbb3SSatish Balay   for ( i=0; i<mbs; i++ ) {
5702d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
5712d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2];
5722d61bbb3SSatish Balay     for ( j=0; j<n; j++ ) {
5732d61bbb3SSatish Balay       xb = x + 3*(*idx++); x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
5742d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[3]*x2 + v[6]*x3;
5752d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[4]*x2 + v[7]*x3;
5762d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[5]*x2 + v[8]*x3;
5772d61bbb3SSatish Balay       v += 9;
5782d61bbb3SSatish Balay     }
5792d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3;
5802d61bbb3SSatish Balay     z += 3; y += 3;
5812d61bbb3SSatish Balay   }
5822d61bbb3SSatish Balay   VecRestoreArray_Fast(xx,x);
5832d61bbb3SSatish Balay   VecRestoreArray_Fast(yy,y);
5842d61bbb3SSatish Balay   VecRestoreArray_Fast(zz,z);
5852d61bbb3SSatish Balay   PLogFlops(18*a->nz);
5862d61bbb3SSatish Balay   PetscFunctionReturn(0);
5872d61bbb3SSatish Balay }
5882d61bbb3SSatish Balay 
5892d61bbb3SSatish Balay #undef __FUNC__
5902d61bbb3SSatish Balay #define __FUNC__ "MatMultAdd_SeqBAIJ_4"
5912d61bbb3SSatish Balay int MatMultAdd_SeqBAIJ_4(Mat A,Vec xx,Vec yy,Vec zz)
5922d61bbb3SSatish Balay {
5932d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
5942d61bbb3SSatish Balay   register Scalar *x,*y,*z,*v,*xb,sum1,sum2,sum3,sum4,x1,x2,x3,x4;
5952d61bbb3SSatish Balay   int             mbs=a->mbs,i,*idx,*ii;
5962d61bbb3SSatish Balay   int             j,n;
5972d61bbb3SSatish Balay 
5982d61bbb3SSatish Balay   PetscFunctionBegin;
5992d61bbb3SSatish Balay   VecGetArray_Fast(xx,x);
6002d61bbb3SSatish Balay   VecGetArray_Fast(yy,y);
6012d61bbb3SSatish Balay   VecGetArray_Fast(zz,z);
6022d61bbb3SSatish Balay 
6032d61bbb3SSatish Balay   idx   = a->j;
6042d61bbb3SSatish Balay   v     = a->a;
6052d61bbb3SSatish Balay   ii    = a->i;
6062d61bbb3SSatish Balay 
6072d61bbb3SSatish Balay   for ( i=0; i<mbs; i++ ) {
6082d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
6092d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3];
6102d61bbb3SSatish Balay     for ( j=0; j<n; j++ ) {
6112d61bbb3SSatish Balay       xb = x + 4*(*idx++);
6122d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
6132d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[4]*x2 + v[8]*x3   + v[12]*x4;
6142d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[5]*x2 + v[9]*x3   + v[13]*x4;
6152d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[6]*x2 + v[10]*x3  + v[14]*x4;
6162d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[7]*x2 + v[11]*x3  + v[15]*x4;
6172d61bbb3SSatish Balay       v += 16;
6182d61bbb3SSatish Balay     }
6192d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4;
6202d61bbb3SSatish Balay     z += 4; y += 4;
6212d61bbb3SSatish Balay   }
6222d61bbb3SSatish Balay   VecRestoreArray_Fast(xx,x);
6232d61bbb3SSatish Balay   VecRestoreArray_Fast(yy,y);
6242d61bbb3SSatish Balay   VecRestoreArray_Fast(zz,z);
6252d61bbb3SSatish Balay   PLogFlops(32*a->nz);
6262d61bbb3SSatish Balay   PetscFunctionReturn(0);
6272d61bbb3SSatish Balay }
6282d61bbb3SSatish Balay 
6292d61bbb3SSatish Balay #undef __FUNC__
6302d61bbb3SSatish Balay #define __FUNC__ "MatMultAdd_SeqBAIJ_5"
6312d61bbb3SSatish Balay int MatMultAdd_SeqBAIJ_5(Mat A,Vec xx,Vec yy,Vec zz)
6322d61bbb3SSatish Balay {
6332d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
6342d61bbb3SSatish Balay   register Scalar *x,*y,*z,*v,*xb,sum1,sum2,sum3,sum4,sum5,x1,x2,x3,x4,x5;
6352d61bbb3SSatish Balay   int             mbs=a->mbs,i,*idx,*ii,j,n;
6362d61bbb3SSatish Balay 
6372d61bbb3SSatish Balay   PetscFunctionBegin;
6382d61bbb3SSatish Balay   VecGetArray_Fast(xx,x);
6392d61bbb3SSatish Balay   VecGetArray_Fast(yy,y);
6402d61bbb3SSatish Balay   VecGetArray_Fast(zz,z);
6412d61bbb3SSatish Balay 
6422d61bbb3SSatish Balay   idx   = a->j;
6432d61bbb3SSatish Balay   v     = a->a;
6442d61bbb3SSatish Balay   ii    = a->i;
6452d61bbb3SSatish Balay 
6462d61bbb3SSatish Balay   for ( i=0; i<mbs; i++ ) {
6472d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
6482d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4];
6492d61bbb3SSatish Balay     for ( j=0; j<n; j++ ) {
6502d61bbb3SSatish Balay       xb = x + 5*(*idx++);
6512d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4];
6522d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[5]*x2 + v[10]*x3  + v[15]*x4 + v[20]*x5;
6532d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[6]*x2 + v[11]*x3  + v[16]*x4 + v[21]*x5;
6542d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[7]*x2 + v[12]*x3  + v[17]*x4 + v[22]*x5;
6552d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[8]*x2 + v[13]*x3  + v[18]*x4 + v[23]*x5;
6562d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[9]*x2 + v[14]*x3  + v[19]*x4 + v[24]*x5;
6572d61bbb3SSatish Balay       v += 25;
6582d61bbb3SSatish Balay     }
6592d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5;
6602d61bbb3SSatish Balay     z += 5; y += 5;
6612d61bbb3SSatish Balay   }
6622d61bbb3SSatish Balay   VecRestoreArray_Fast(xx,x);
6632d61bbb3SSatish Balay   VecRestoreArray_Fast(yy,y);
6642d61bbb3SSatish Balay   VecRestoreArray_Fast(zz,z);
6652d61bbb3SSatish Balay   PLogFlops(50*a->nz);
6662d61bbb3SSatish Balay   PetscFunctionReturn(0);
6672d61bbb3SSatish Balay }
6682d61bbb3SSatish Balay 
6692d61bbb3SSatish Balay #undef __FUNC__
6702d61bbb3SSatish Balay #define __FUNC__ "MatMultAdd_SeqBAIJ_7"
6712d61bbb3SSatish Balay int MatMultAdd_SeqBAIJ_7(Mat A,Vec xx,Vec yy,Vec zz)
6722d61bbb3SSatish Balay {
6732d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
6742d61bbb3SSatish Balay   register Scalar *x,*y,*z,*v,*xb,sum1,sum2,sum3,sum4,sum5,sum6,sum7;
6752d61bbb3SSatish Balay   register Scalar x1,x2,x3,x4,x5,x6,x7;
6762d61bbb3SSatish Balay   int             mbs=a->mbs,i,*idx,*ii,j,n;
6772d61bbb3SSatish Balay 
6782d61bbb3SSatish Balay   PetscFunctionBegin;
6792d61bbb3SSatish Balay   VecGetArray_Fast(xx,x);
6802d61bbb3SSatish Balay   VecGetArray_Fast(yy,y);
6812d61bbb3SSatish Balay   VecGetArray_Fast(zz,z);
6822d61bbb3SSatish Balay 
6832d61bbb3SSatish Balay   idx   = a->j;
6842d61bbb3SSatish Balay   v     = a->a;
6852d61bbb3SSatish Balay   ii    = a->i;
6862d61bbb3SSatish Balay 
6872d61bbb3SSatish Balay   for ( i=0; i<mbs; i++ ) {
6882d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
6892d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4]; sum6 = y[5]; sum7 = y[6];
6902d61bbb3SSatish Balay     for ( j=0; j<n; j++ ) {
6912d61bbb3SSatish Balay       xb = x + 7*(*idx++);
6922d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6];
6932d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[7]*x2  + v[14]*x3  + v[21]*x4 + v[28]*x5 + v[35]*x6 + v[42]*x7;
6942d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[8]*x2  + v[15]*x3  + v[22]*x4 + v[29]*x5 + v[36]*x6 + v[43]*x7;
6952d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[9]*x2  + v[16]*x3  + v[23]*x4 + v[30]*x5 + v[37]*x6 + v[44]*x7;
6962d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[10]*x2 + v[17]*x3  + v[24]*x4 + v[31]*x5 + v[38]*x6 + v[45]*x7;
6972d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[11]*x2 + v[18]*x3  + v[25]*x4 + v[32]*x5 + v[39]*x6 + v[46]*x7;
6982d61bbb3SSatish Balay       sum6 += v[5]*x1 + v[12]*x2 + v[19]*x3  + v[26]*x4 + v[33]*x5 + v[40]*x6 + v[47]*x7;
6992d61bbb3SSatish Balay       sum7 += v[6]*x1 + v[13]*x2 + v[20]*x3  + v[27]*x4 + v[34]*x5 + v[41]*x6 + v[48]*x7;
7002d61bbb3SSatish Balay       v += 49;
7012d61bbb3SSatish Balay     }
7022d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
7032d61bbb3SSatish Balay     z += 7; y += 7;
7042d61bbb3SSatish Balay   }
7052d61bbb3SSatish Balay   VecRestoreArray_Fast(xx,x);
7062d61bbb3SSatish Balay   VecRestoreArray_Fast(yy,y);
7072d61bbb3SSatish Balay   VecRestoreArray_Fast(zz,z);
7082d61bbb3SSatish Balay   PLogFlops(98*a->nz);
7092d61bbb3SSatish Balay   PetscFunctionReturn(0);
7102d61bbb3SSatish Balay }
711218c64b6SSatish Balay 
712218c64b6SSatish Balay 
7132d61bbb3SSatish Balay #undef __FUNC__
7142d61bbb3SSatish Balay #define __FUNC__ "MatMultAdd_SeqBAIJ_N"
7152d61bbb3SSatish Balay int MatMultAdd_SeqBAIJ_N(Mat A,Vec xx,Vec yy,Vec zz)
7162d61bbb3SSatish Balay {
7172d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
7182d61bbb3SSatish Balay   register Scalar *x,*z,*v,*xb;
7192d61bbb3SSatish Balay   int             mbs=a->mbs,i,*idx,*ii,bs=a->bs,j,n,bs2=a->bs2,ierr;
7202d61bbb3SSatish Balay   int             _One = 1,ncols,k; Scalar _DOne = 1.0, *work,*workt;
721218c64b6SSatish Balay 
7222d61bbb3SSatish Balay   PetscFunctionBegin;
7232d61bbb3SSatish Balay   if ( xx != yy) { ierr = VecCopy(yy,zz); CHKERRQ(ierr); }
7242d61bbb3SSatish Balay 
7252d61bbb3SSatish Balay   VecGetArray_Fast(xx,x);
7262d61bbb3SSatish Balay   VecGetArray_Fast(zz,z);
7272d61bbb3SSatish Balay 
7282d61bbb3SSatish Balay   idx   = a->j;
7292d61bbb3SSatish Balay   v     = a->a;
7302d61bbb3SSatish Balay   ii    = a->i;
7312d61bbb3SSatish Balay 
7322d61bbb3SSatish Balay 
7332d61bbb3SSatish Balay   if (!a->mult_work) {
7342d61bbb3SSatish Balay     k = PetscMax(a->m,a->n);
7352d61bbb3SSatish Balay     a->mult_work = (Scalar *) PetscMalloc(k*sizeof(Scalar));CHKPTRQ(a->mult_work);
7362d61bbb3SSatish Balay   }
7372d61bbb3SSatish Balay   work = a->mult_work;
7382d61bbb3SSatish Balay   for ( i=0; i<mbs; i++ ) {
7392d61bbb3SSatish Balay     n     = ii[1] - ii[0]; ii++;
7402d61bbb3SSatish Balay     ncols = n*bs;
7412d61bbb3SSatish Balay     workt = work;
7422d61bbb3SSatish Balay     for ( j=0; j<n; j++ ) {
7432d61bbb3SSatish Balay       xb = x + bs*(*idx++);
7442d61bbb3SSatish Balay       for ( k=0; k<bs; k++ ) workt[k] = xb[k];
7452d61bbb3SSatish Balay       workt += bs;
7462d61bbb3SSatish Balay     }
7472d61bbb3SSatish Balay     LAgemv_("N",&bs,&ncols,&_DOne,v,&bs,work,&_One,&_DOne,z,&_One);
7482d61bbb3SSatish Balay     v += n*bs2;
7492d61bbb3SSatish Balay     z += bs;
7502d61bbb3SSatish Balay   }
7512d61bbb3SSatish Balay   VecRestoreArray_Fast(xx,x);
7522d61bbb3SSatish Balay   VecRestoreArray_Fast(zz,z);
7532d61bbb3SSatish Balay   PLogFlops(2*a->nz*bs2 );
7542d61bbb3SSatish Balay   PetscFunctionReturn(0);
7552d61bbb3SSatish Balay }
7562d61bbb3SSatish Balay 
7572d61bbb3SSatish Balay #undef __FUNC__
7582d61bbb3SSatish Balay #define __FUNC__ "MatMultTrans_SeqBAIJ"
7592d61bbb3SSatish Balay int MatMultTrans_SeqBAIJ(Mat A,Vec xx,Vec zz)
7602d61bbb3SSatish Balay {
7612d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
7622d61bbb3SSatish Balay   Scalar          *xg,*zg,*zb;
7632d61bbb3SSatish Balay   register Scalar *x,*z,*v,*xb,x1,x2,x3,x4,x5;
7642d61bbb3SSatish Balay   int             mbs=a->mbs,i,*idx,*ii,*ai=a->i,rval,N=a->n;
7652d61bbb3SSatish Balay   int             bs=a->bs,j,n,bs2=a->bs2,*ib,ierr;
7662d61bbb3SSatish Balay 
7672d61bbb3SSatish Balay 
7682d61bbb3SSatish Balay   PetscFunctionBegin;
7692d61bbb3SSatish Balay   VecGetArray_Fast(xx,xg); x = xg;
7702d61bbb3SSatish Balay   VecGetArray_Fast(zz,zg); z = zg;
7712d61bbb3SSatish Balay   PetscMemzero(z,N*sizeof(Scalar));
7722d61bbb3SSatish Balay 
7732d61bbb3SSatish Balay   idx   = a->j;
7742d61bbb3SSatish Balay   v     = a->a;
7752d61bbb3SSatish Balay   ii    = a->i;
7762d61bbb3SSatish Balay 
7772d61bbb3SSatish Balay   switch (bs) {
7782d61bbb3SSatish Balay   case 1:
7792d61bbb3SSatish Balay     for ( i=0; i<mbs; i++ ) {
7802d61bbb3SSatish Balay       n  = ii[1] - ii[0]; ii++;
7812d61bbb3SSatish Balay       xb = x + i; x1 = xb[0];
7822d61bbb3SSatish Balay       ib = idx + ai[i];
7832d61bbb3SSatish Balay       for ( j=0; j<n; j++ ) {
7842d61bbb3SSatish Balay         rval    = ib[j];
7852d61bbb3SSatish Balay         z[rval] += *v++ * x1;
7862d61bbb3SSatish Balay       }
7872d61bbb3SSatish Balay     }
7882d61bbb3SSatish Balay     break;
7892d61bbb3SSatish Balay   case 2:
7902d61bbb3SSatish Balay     for ( i=0; i<mbs; i++ ) {
7912d61bbb3SSatish Balay       n  = ii[1] - ii[0]; ii++;
7922d61bbb3SSatish Balay       xb = x + 2*i; x1 = xb[0]; x2 = xb[1];
7932d61bbb3SSatish Balay       ib = idx + ai[i];
7942d61bbb3SSatish Balay       for ( j=0; j<n; j++ ) {
7952d61bbb3SSatish Balay         rval      = ib[j]*2;
7962d61bbb3SSatish Balay         z[rval++] += v[0]*x1 + v[1]*x2;
7972d61bbb3SSatish Balay         z[rval++] += v[2]*x1 + v[3]*x2;
7982d61bbb3SSatish Balay         v += 4;
7992d61bbb3SSatish Balay       }
8002d61bbb3SSatish Balay     }
8012d61bbb3SSatish Balay     break;
8022d61bbb3SSatish Balay   case 3:
8032d61bbb3SSatish Balay     for ( i=0; i<mbs; i++ ) {
8042d61bbb3SSatish Balay       n  = ii[1] - ii[0]; ii++;
8052d61bbb3SSatish Balay       xb = x + 3*i; x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
8062d61bbb3SSatish Balay       ib = idx + ai[i];
8072d61bbb3SSatish Balay       for ( j=0; j<n; j++ ) {
8082d61bbb3SSatish Balay         rval      = ib[j]*3;
8092d61bbb3SSatish Balay         z[rval++] += v[0]*x1 + v[1]*x2 + v[2]*x3;
8102d61bbb3SSatish Balay         z[rval++] += v[3]*x1 + v[4]*x2 + v[5]*x3;
8112d61bbb3SSatish Balay         z[rval++] += v[6]*x1 + v[7]*x2 + v[8]*x3;
8122d61bbb3SSatish Balay         v += 9;
8132d61bbb3SSatish Balay       }
8142d61bbb3SSatish Balay     }
8152d61bbb3SSatish Balay     break;
8162d61bbb3SSatish Balay   case 4:
8172d61bbb3SSatish Balay     for ( i=0; i<mbs; i++ ) {
8182d61bbb3SSatish Balay       n  = ii[1] - ii[0]; ii++;
8192d61bbb3SSatish Balay       xb = x + 4*i; x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
8202d61bbb3SSatish Balay       ib = idx + ai[i];
8212d61bbb3SSatish Balay       for ( j=0; j<n; j++ ) {
8222d61bbb3SSatish Balay         rval      = ib[j]*4;
8232d61bbb3SSatish Balay         z[rval++] +=  v[0]*x1 +  v[1]*x2 +  v[2]*x3 +  v[3]*x4;
8242d61bbb3SSatish Balay         z[rval++] +=  v[4]*x1 +  v[5]*x2 +  v[6]*x3 +  v[7]*x4;
8252d61bbb3SSatish Balay         z[rval++] +=  v[8]*x1 +  v[9]*x2 + v[10]*x3 + v[11]*x4;
8262d61bbb3SSatish Balay         z[rval++] += v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4;
8272d61bbb3SSatish Balay         v += 16;
8282d61bbb3SSatish Balay       }
8292d61bbb3SSatish Balay     }
8302d61bbb3SSatish Balay     break;
8312d61bbb3SSatish Balay   case 5:
8322d61bbb3SSatish Balay     for ( i=0; i<mbs; i++ ) {
8332d61bbb3SSatish Balay       n  = ii[1] - ii[0]; ii++;
8342d61bbb3SSatish Balay       xb = x + 5*i; x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
8352d61bbb3SSatish Balay       x4 = xb[3];   x5 = xb[4];
8362d61bbb3SSatish Balay       ib = idx + ai[i];
8372d61bbb3SSatish Balay       for ( j=0; j<n; j++ ) {
8382d61bbb3SSatish Balay         rval      = ib[j]*5;
8392d61bbb3SSatish Balay         z[rval++] +=  v[0]*x1 +  v[1]*x2 +  v[2]*x3 +  v[3]*x4 +  v[4]*x5;
8402d61bbb3SSatish Balay         z[rval++] +=  v[5]*x1 +  v[6]*x2 +  v[7]*x3 +  v[8]*x4 +  v[9]*x5;
8412d61bbb3SSatish Balay         z[rval++] += v[10]*x1 + v[11]*x2 + v[12]*x3 + v[13]*x4 + v[14]*x5;
8422d61bbb3SSatish Balay         z[rval++] += v[15]*x1 + v[16]*x2 + v[17]*x3 + v[18]*x4 + v[19]*x5;
8432d61bbb3SSatish Balay         z[rval++] += v[20]*x1 + v[21]*x2 + v[22]*x3 + v[23]*x4 + v[24]*x5;
8442d61bbb3SSatish Balay         v += 25;
8452d61bbb3SSatish Balay       }
8462d61bbb3SSatish Balay     }
8472d61bbb3SSatish Balay     break;
8482d61bbb3SSatish Balay       /* block sizes larger then 5 by 5 are handled by BLAS */
8492d61bbb3SSatish Balay     default: {
8502d61bbb3SSatish Balay       int  _One = 1,ncols,k; Scalar _DOne = 1.0, *work,*workt;
8512d61bbb3SSatish Balay       if (!a->mult_work) {
8522d61bbb3SSatish Balay         k = PetscMax(a->m,a->n);
8532d61bbb3SSatish Balay         a->mult_work = (Scalar *) PetscMalloc(k*sizeof(Scalar));
8542d61bbb3SSatish Balay         CHKPTRQ(a->mult_work);
8552d61bbb3SSatish Balay       }
8562d61bbb3SSatish Balay       work = a->mult_work;
8572d61bbb3SSatish Balay       for ( i=0; i<mbs; i++ ) {
8582d61bbb3SSatish Balay         n     = ii[1] - ii[0]; ii++;
8592d61bbb3SSatish Balay         ncols = n*bs;
8602d61bbb3SSatish Balay         PetscMemzero(work,ncols*sizeof(Scalar));
8612d61bbb3SSatish Balay         LAgemv_("T",&bs,&ncols,&_DOne,v,&bs,x,&_One,&_DOne,work,&_One);
8622d61bbb3SSatish Balay         v += n*bs2;
8632d61bbb3SSatish Balay         x += bs;
8642d61bbb3SSatish Balay         workt = work;
8652d61bbb3SSatish Balay         for ( j=0; j<n; j++ ) {
8662d61bbb3SSatish Balay           zb = z + bs*(*idx++);
8672d61bbb3SSatish Balay           for ( k=0; k<bs; k++ ) zb[k] += workt[k] ;
8682d61bbb3SSatish Balay           workt += bs;
8692d61bbb3SSatish Balay         }
8702d61bbb3SSatish Balay       }
8712d61bbb3SSatish Balay     }
8722d61bbb3SSatish Balay   }
8732d61bbb3SSatish Balay   ierr = VecRestoreArray(xx,&xg); CHKERRQ(ierr);
8742d61bbb3SSatish Balay   ierr = VecRestoreArray(zz,&zg); CHKERRQ(ierr);
8752d61bbb3SSatish Balay   PLogFlops(2*a->nz*a->bs2 - a->n);
8762d61bbb3SSatish Balay   PetscFunctionReturn(0);
8772d61bbb3SSatish Balay }
8782d61bbb3SSatish Balay 
8792d61bbb3SSatish Balay #undef __FUNC__
8802d61bbb3SSatish Balay #define __FUNC__ "MatMultTransAdd_SeqBAIJ"
8812d61bbb3SSatish Balay int MatMultTransAdd_SeqBAIJ(Mat A,Vec xx,Vec yy,Vec zz)
8822d61bbb3SSatish Balay 
8832d61bbb3SSatish Balay {
8842d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *) A->data;
8852d61bbb3SSatish Balay   Scalar          *xg,*zg,*zb;
8862d61bbb3SSatish Balay   register Scalar *x,*z,*v,*xb,x1,x2,x3,x4,x5;
8872d61bbb3SSatish Balay   int             mbs=a->mbs,i,*idx,*ii,*ai=a->i,rval,N=a->n;
8882d61bbb3SSatish Balay   int             bs=a->bs,j,n,bs2=a->bs2,*ib,ierr;
8892d61bbb3SSatish Balay 
8902d61bbb3SSatish Balay   PetscFunctionBegin;
8912d61bbb3SSatish Balay   VecGetArray_Fast(xx,xg); x = xg;
8922d61bbb3SSatish Balay   VecGetArray_Fast(zz,zg); z = zg;
8932d61bbb3SSatish Balay 
8942d61bbb3SSatish Balay   if ( yy != zz ) { ierr = VecCopy(yy,zz); CHKERRQ(ierr); }
8952d61bbb3SSatish Balay   else PetscMemzero(z,N*sizeof(Scalar));
8962d61bbb3SSatish Balay 
8972d61bbb3SSatish Balay   idx   = a->j;
8982d61bbb3SSatish Balay   v     = a->a;
8992d61bbb3SSatish Balay   ii    = a->i;
9002d61bbb3SSatish Balay 
9012d61bbb3SSatish Balay   switch (bs) {
9022d61bbb3SSatish Balay   case 1:
9032d61bbb3SSatish Balay     for ( i=0; i<mbs; i++ ) {
9042d61bbb3SSatish Balay       n  = ii[1] - ii[0]; ii++;
9052d61bbb3SSatish Balay       xb = x + i; x1 = xb[0];
9062d61bbb3SSatish Balay       ib = idx + ai[i];
9072d61bbb3SSatish Balay       for ( j=0; j<n; j++ ) {
9082d61bbb3SSatish Balay         rval    = ib[j];
9092d61bbb3SSatish Balay         z[rval] += *v++ * x1;
9102d61bbb3SSatish Balay       }
9112d61bbb3SSatish Balay     }
9122d61bbb3SSatish Balay     break;
9132d61bbb3SSatish Balay   case 2:
9142d61bbb3SSatish Balay     for ( i=0; i<mbs; i++ ) {
9152d61bbb3SSatish Balay       n  = ii[1] - ii[0]; ii++;
9162d61bbb3SSatish Balay       xb = x + 2*i; x1 = xb[0]; x2 = xb[1];
9172d61bbb3SSatish Balay       ib = idx + ai[i];
9182d61bbb3SSatish Balay       for ( j=0; j<n; j++ ) {
9192d61bbb3SSatish Balay         rval      = ib[j]*2;
9202d61bbb3SSatish Balay         z[rval++] += v[0]*x1 + v[1]*x2;
9212d61bbb3SSatish Balay         z[rval++] += v[2]*x1 + v[3]*x2;
9222d61bbb3SSatish Balay         v += 4;
9232d61bbb3SSatish Balay       }
9242d61bbb3SSatish Balay     }
9252d61bbb3SSatish Balay     break;
9262d61bbb3SSatish Balay   case 3:
9272d61bbb3SSatish Balay     for ( i=0; i<mbs; i++ ) {
9282d61bbb3SSatish Balay       n  = ii[1] - ii[0]; ii++;
9292d61bbb3SSatish Balay       xb = x + 3*i; x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
9302d61bbb3SSatish Balay       ib = idx + ai[i];
9312d61bbb3SSatish Balay       for ( j=0; j<n; j++ ) {
9322d61bbb3SSatish Balay         rval      = ib[j]*3;
9332d61bbb3SSatish Balay         z[rval++] += v[0]*x1 + v[1]*x2 + v[2]*x3;
9342d61bbb3SSatish Balay         z[rval++] += v[3]*x1 + v[4]*x2 + v[5]*x3;
9352d61bbb3SSatish Balay         z[rval++] += v[6]*x1 + v[7]*x2 + v[8]*x3;
9362d61bbb3SSatish Balay         v += 9;
9372d61bbb3SSatish Balay       }
9382d61bbb3SSatish Balay     }
9392d61bbb3SSatish Balay     break;
9402d61bbb3SSatish Balay   case 4:
9412d61bbb3SSatish Balay     for ( i=0; i<mbs; i++ ) {
9422d61bbb3SSatish Balay       n  = ii[1] - ii[0]; ii++;
9432d61bbb3SSatish Balay       xb = x + 4*i; x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
9442d61bbb3SSatish Balay       ib = idx + ai[i];
9452d61bbb3SSatish Balay       for ( j=0; j<n; j++ ) {
9462d61bbb3SSatish Balay         rval      = ib[j]*4;
9472d61bbb3SSatish Balay         z[rval++] +=  v[0]*x1 +  v[1]*x2 +  v[2]*x3 +  v[3]*x4;
9482d61bbb3SSatish Balay         z[rval++] +=  v[4]*x1 +  v[5]*x2 +  v[6]*x3 +  v[7]*x4;
9492d61bbb3SSatish Balay         z[rval++] +=  v[8]*x1 +  v[9]*x2 + v[10]*x3 + v[11]*x4;
9502d61bbb3SSatish Balay         z[rval++] += v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4;
9512d61bbb3SSatish Balay         v += 16;
9522d61bbb3SSatish Balay       }
9532d61bbb3SSatish Balay     }
9542d61bbb3SSatish Balay     break;
9552d61bbb3SSatish Balay   case 5:
9562d61bbb3SSatish Balay     for ( i=0; i<mbs; i++ ) {
9572d61bbb3SSatish Balay       n  = ii[1] - ii[0]; ii++;
9582d61bbb3SSatish Balay       xb = x + 5*i; x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
9592d61bbb3SSatish Balay       x4 = xb[3];   x5 = xb[4];
9602d61bbb3SSatish Balay       ib = idx + ai[i];
9612d61bbb3SSatish Balay       for ( j=0; j<n; j++ ) {
9622d61bbb3SSatish Balay         rval      = ib[j]*5;
9632d61bbb3SSatish Balay         z[rval++] +=  v[0]*x1 +  v[1]*x2 +  v[2]*x3 +  v[3]*x4 +  v[4]*x5;
9642d61bbb3SSatish Balay         z[rval++] +=  v[5]*x1 +  v[6]*x2 +  v[7]*x3 +  v[8]*x4 +  v[9]*x5;
9652d61bbb3SSatish Balay         z[rval++] += v[10]*x1 + v[11]*x2 + v[12]*x3 + v[13]*x4 + v[14]*x5;
9662d61bbb3SSatish Balay         z[rval++] += v[15]*x1 + v[16]*x2 + v[17]*x3 + v[18]*x4 + v[19]*x5;
9672d61bbb3SSatish Balay         z[rval++] += v[20]*x1 + v[21]*x2 + v[22]*x3 + v[23]*x4 + v[24]*x5;
9682d61bbb3SSatish Balay         v += 25;
9692d61bbb3SSatish Balay       }
9702d61bbb3SSatish Balay     }
9712d61bbb3SSatish Balay     break;
9722d61bbb3SSatish Balay       /* block sizes larger then 5 by 5 are handled by BLAS */
9732d61bbb3SSatish Balay     default: {
9742d61bbb3SSatish Balay       int  _One = 1,ncols,k; Scalar _DOne = 1.0, *work,*workt;
9752d61bbb3SSatish Balay       if (!a->mult_work) {
9762d61bbb3SSatish Balay         k = PetscMax(a->m,a->n);
9772d61bbb3SSatish Balay         a->mult_work = (Scalar *) PetscMalloc(k*sizeof(Scalar));
9782d61bbb3SSatish Balay         CHKPTRQ(a->mult_work);
9792d61bbb3SSatish Balay       }
9802d61bbb3SSatish Balay       work = a->mult_work;
9812d61bbb3SSatish Balay       for ( i=0; i<mbs; i++ ) {
9822d61bbb3SSatish Balay         n     = ii[1] - ii[0]; ii++;
9832d61bbb3SSatish Balay         ncols = n*bs;
9842d61bbb3SSatish Balay         PetscMemzero(work,ncols*sizeof(Scalar));
9852d61bbb3SSatish Balay         LAgemv_("T",&bs,&ncols,&_DOne,v,&bs,x,&_One,&_DOne,work,&_One);
9862d61bbb3SSatish Balay         v += n*bs2;
9872d61bbb3SSatish Balay         x += bs;
9882d61bbb3SSatish Balay         workt = work;
9892d61bbb3SSatish Balay         for ( j=0; j<n; j++ ) {
9902d61bbb3SSatish Balay           zb = z + bs*(*idx++);
9912d61bbb3SSatish Balay           for ( k=0; k<bs; k++ ) zb[k] += workt[k] ;
9922d61bbb3SSatish Balay           workt += bs;
9932d61bbb3SSatish Balay         }
9942d61bbb3SSatish Balay       }
9952d61bbb3SSatish Balay     }
9962d61bbb3SSatish Balay   }
9972d61bbb3SSatish Balay   ierr = VecRestoreArray(xx,&xg); CHKERRQ(ierr);
9982d61bbb3SSatish Balay   ierr = VecRestoreArray(zz,&zg); CHKERRQ(ierr);
9992d61bbb3SSatish Balay   PLogFlops(2*a->nz*a->bs2);
10002d61bbb3SSatish Balay   PetscFunctionReturn(0);
10012d61bbb3SSatish Balay }
10022d61bbb3SSatish Balay 
10032d61bbb3SSatish Balay #undef __FUNC__
10042d61bbb3SSatish Balay #define __FUNC__ "MatScale_SeqBAIJ"
10052d61bbb3SSatish Balay int MatScale_SeqBAIJ(Scalar *alpha,Mat inA)
10062d61bbb3SSatish Balay {
10072d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) inA->data;
10082d61bbb3SSatish Balay   int         one = 1, totalnz = a->bs2*a->nz;
10092d61bbb3SSatish Balay 
10102d61bbb3SSatish Balay   PetscFunctionBegin;
10112d61bbb3SSatish Balay   BLscal_( &totalnz, alpha, a->a, &one );
10122d61bbb3SSatish Balay   PLogFlops(totalnz);
10132d61bbb3SSatish Balay   PetscFunctionReturn(0);
10142d61bbb3SSatish Balay }
10152d61bbb3SSatish Balay 
10162d61bbb3SSatish Balay #undef __FUNC__
10172d61bbb3SSatish Balay #define __FUNC__ "MatNorm_SeqBAIJ"
10182d61bbb3SSatish Balay int MatNorm_SeqBAIJ(Mat A,NormType type,double *norm)
10192d61bbb3SSatish Balay {
10202d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
10212d61bbb3SSatish Balay   Scalar      *v = a->a;
10222d61bbb3SSatish Balay   double      sum = 0.0;
10232d61bbb3SSatish Balay   int         i,nz=a->nz,bs2=a->bs2;
10242d61bbb3SSatish Balay 
10252d61bbb3SSatish Balay   PetscFunctionBegin;
10262d61bbb3SSatish Balay   if (type == NORM_FROBENIUS) {
10272d61bbb3SSatish Balay     for (i=0; i< bs2*nz; i++ ) {
10282d61bbb3SSatish Balay #if defined(USE_PETSC_COMPLEX)
10292d61bbb3SSatish Balay       sum += real(conj(*v)*(*v)); v++;
10302d61bbb3SSatish Balay #else
10312d61bbb3SSatish Balay       sum += (*v)*(*v); v++;
10322d61bbb3SSatish Balay #endif
10332d61bbb3SSatish Balay     }
10342d61bbb3SSatish Balay     *norm = sqrt(sum);
10352d61bbb3SSatish Balay   }
10362d61bbb3SSatish Balay   else {
10372d61bbb3SSatish Balay     SETERRQ(PETSC_ERR_SUP,0,"No support for this norm yet");
10382d61bbb3SSatish Balay   }
10392d61bbb3SSatish Balay   PetscFunctionReturn(0);
10402d61bbb3SSatish Balay }
10412d61bbb3SSatish Balay 
10422d61bbb3SSatish Balay 
10432d61bbb3SSatish Balay #undef __FUNC__
10442d61bbb3SSatish Balay #define __FUNC__ "MatEqual_SeqBAIJ"
10452d61bbb3SSatish Balay int MatEqual_SeqBAIJ(Mat A,Mat B, PetscTruth* flg)
10462d61bbb3SSatish Balay {
10472d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)B->data;
10482d61bbb3SSatish Balay 
10492d61bbb3SSatish Balay   PetscFunctionBegin;
10502d61bbb3SSatish Balay   if (B->type !=MATSEQBAIJ) SETERRQ(PETSC_ERR_ARG_INCOMP,0,"Matrices must be same type");
10512d61bbb3SSatish Balay 
10522d61bbb3SSatish Balay   /* If the  matrix/block dimensions are not equal, or no of nonzeros or shift */
10532d61bbb3SSatish Balay   if ((a->m != b->m) || (a->n !=b->n) || (a->bs != b->bs)|| (a->nz != b->nz)) {
10542d61bbb3SSatish Balay     *flg = PETSC_FALSE; PetscFunctionReturn(0);
10552d61bbb3SSatish Balay   }
10562d61bbb3SSatish Balay 
10572d61bbb3SSatish Balay   /* if the a->i are the same */
10582d61bbb3SSatish Balay   if (PetscMemcmp(a->i,b->i, (a->mbs+1)*sizeof(int))) {
10592d61bbb3SSatish Balay     *flg = PETSC_FALSE; PetscFunctionReturn(0);
10602d61bbb3SSatish Balay   }
10612d61bbb3SSatish Balay 
10622d61bbb3SSatish Balay   /* if a->j are the same */
10632d61bbb3SSatish Balay   if (PetscMemcmp(a->j,b->j,(a->nz)*sizeof(int))) {
10642d61bbb3SSatish Balay     *flg = PETSC_FALSE; PetscFunctionReturn(0);
10652d61bbb3SSatish Balay   }
10662d61bbb3SSatish Balay 
10672d61bbb3SSatish Balay   /* if a->a are the same */
10682d61bbb3SSatish Balay   if (PetscMemcmp(a->a, b->a,(a->nz)*(a->bs)*(a->bs)*sizeof(Scalar))) {
10692d61bbb3SSatish Balay     *flg = PETSC_FALSE; PetscFunctionReturn(0);
10702d61bbb3SSatish Balay   }
10712d61bbb3SSatish Balay   *flg = PETSC_TRUE;
10722d61bbb3SSatish Balay   PetscFunctionReturn(0);
10732d61bbb3SSatish Balay 
10742d61bbb3SSatish Balay }
10752d61bbb3SSatish Balay 
10762d61bbb3SSatish Balay #undef __FUNC__
10772d61bbb3SSatish Balay #define __FUNC__ "MatGetDiagonal_SeqBAIJ"
10782d61bbb3SSatish Balay int MatGetDiagonal_SeqBAIJ(Mat A,Vec v)
10792d61bbb3SSatish Balay {
10802d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
10812d61bbb3SSatish Balay   int         i,j,k,n,row,bs,*ai,*aj,ambs,bs2;
10822d61bbb3SSatish Balay   Scalar      *x, zero = 0.0,*aa,*aa_j;
10832d61bbb3SSatish Balay 
10842d61bbb3SSatish Balay   PetscFunctionBegin;
10852d61bbb3SSatish Balay   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix");
10862d61bbb3SSatish Balay   bs   = a->bs;
10872d61bbb3SSatish Balay   aa   = a->a;
10882d61bbb3SSatish Balay   ai   = a->i;
10892d61bbb3SSatish Balay   aj   = a->j;
10902d61bbb3SSatish Balay   ambs = a->mbs;
10912d61bbb3SSatish Balay   bs2  = a->bs2;
10922d61bbb3SSatish Balay 
10932d61bbb3SSatish Balay   VecSet(&zero,v);
10942d61bbb3SSatish Balay   VecGetArray_Fast(v,x); VecGetLocalSize_Fast(v,n);
10952d61bbb3SSatish Balay   if (n != a->m) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Nonconforming matrix and vector");
10962d61bbb3SSatish Balay   for ( i=0; i<ambs; i++ ) {
10972d61bbb3SSatish Balay     for ( j=ai[i]; j<ai[i+1]; j++ ) {
10982d61bbb3SSatish Balay       if (aj[j] == i) {
10992d61bbb3SSatish Balay         row  = i*bs;
11002d61bbb3SSatish Balay         aa_j = aa+j*bs2;
11012d61bbb3SSatish Balay         for (k=0; k<bs2; k+=(bs+1),row++) x[row] = aa_j[k];
11022d61bbb3SSatish Balay         break;
11032d61bbb3SSatish Balay       }
11042d61bbb3SSatish Balay     }
11052d61bbb3SSatish Balay   }
11062d61bbb3SSatish Balay   PetscFunctionReturn(0);
11072d61bbb3SSatish Balay }
11082d61bbb3SSatish Balay 
11092d61bbb3SSatish Balay #undef __FUNC__
11102d61bbb3SSatish Balay #define __FUNC__ "MatDiagonalScale_SeqBAIJ"
11112d61bbb3SSatish Balay int MatDiagonalScale_SeqBAIJ(Mat A,Vec ll,Vec rr)
11122d61bbb3SSatish Balay {
11132d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
11142d61bbb3SSatish Balay   Scalar      *l,*r,x,*v,*aa,*li,*ri;
11152d61bbb3SSatish Balay   int         i,j,k,lm,rn,M,m,n,*ai,*aj,mbs,tmp,bs,bs2;
11162d61bbb3SSatish Balay 
11172d61bbb3SSatish Balay   PetscFunctionBegin;
11182d61bbb3SSatish Balay   ai  = a->i;
11192d61bbb3SSatish Balay   aj  = a->j;
11202d61bbb3SSatish Balay   aa  = a->a;
11212d61bbb3SSatish Balay   m   = a->m;
11222d61bbb3SSatish Balay   n   = a->n;
11232d61bbb3SSatish Balay   bs  = a->bs;
11242d61bbb3SSatish Balay   mbs = a->mbs;
11252d61bbb3SSatish Balay   bs2 = a->bs2;
11262d61bbb3SSatish Balay   if (ll) {
11272d61bbb3SSatish Balay     VecGetArray_Fast(ll,l); VecGetLocalSize_Fast(ll,lm);
11282d61bbb3SSatish Balay     if (lm != m) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Left scaling vector wrong length");
11292d61bbb3SSatish Balay     for ( i=0; i<mbs; i++ ) { /* for each block row */
11302d61bbb3SSatish Balay       M  = ai[i+1] - ai[i];
11312d61bbb3SSatish Balay       li = l + i*bs;
11322d61bbb3SSatish Balay       v  = aa + bs2*ai[i];
11332d61bbb3SSatish Balay       for ( j=0; j<M; j++ ) { /* for each block */
11342d61bbb3SSatish Balay         for ( k=0; k<bs2; k++ ) {
11352d61bbb3SSatish Balay           (*v++) *= li[k%bs];
11362d61bbb3SSatish Balay         }
11372d61bbb3SSatish Balay       }
11382d61bbb3SSatish Balay     }
11392d61bbb3SSatish Balay   }
11402d61bbb3SSatish Balay 
11412d61bbb3SSatish Balay   if (rr) {
11422d61bbb3SSatish Balay     VecGetArray_Fast(rr,r); VecGetLocalSize_Fast(rr,rn);
11432d61bbb3SSatish Balay     if (rn != n) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Right scaling vector wrong length");
11442d61bbb3SSatish Balay     for ( i=0; i<mbs; i++ ) { /* for each block row */
11452d61bbb3SSatish Balay       M  = ai[i+1] - ai[i];
11462d61bbb3SSatish Balay       v  = aa + bs2*ai[i];
11472d61bbb3SSatish Balay       for ( j=0; j<M; j++ ) { /* for each block */
11482d61bbb3SSatish Balay         ri = r + bs*aj[ai[i]+j];
11492d61bbb3SSatish Balay         for ( k=0; k<bs; k++ ) {
11502d61bbb3SSatish Balay           x = ri[k];
11512d61bbb3SSatish Balay           for ( tmp=0; tmp<bs; tmp++ ) (*v++) *= x;
11522d61bbb3SSatish Balay         }
11532d61bbb3SSatish Balay       }
11542d61bbb3SSatish Balay     }
11552d61bbb3SSatish Balay   }
11562d61bbb3SSatish Balay   PetscFunctionReturn(0);
11572d61bbb3SSatish Balay }
11582d61bbb3SSatish Balay 
11592d61bbb3SSatish Balay 
11602d61bbb3SSatish Balay #undef __FUNC__
11612d61bbb3SSatish Balay #define __FUNC__ "MatGetInfo_SeqBAIJ"
11622d61bbb3SSatish Balay int MatGetInfo_SeqBAIJ(Mat A,MatInfoType flag,MatInfo *info)
11632d61bbb3SSatish Balay {
11642d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
11652d61bbb3SSatish Balay 
11662d61bbb3SSatish Balay   PetscFunctionBegin;
11672d61bbb3SSatish Balay   info->rows_global    = (double)a->m;
11682d61bbb3SSatish Balay   info->columns_global = (double)a->n;
11692d61bbb3SSatish Balay   info->rows_local     = (double)a->m;
11702d61bbb3SSatish Balay   info->columns_local  = (double)a->n;
11712d61bbb3SSatish Balay   info->block_size     = a->bs2;
11722d61bbb3SSatish Balay   info->nz_allocated   = a->maxnz;
11732d61bbb3SSatish Balay   info->nz_used        = a->bs2*a->nz;
11742d61bbb3SSatish Balay   info->nz_unneeded    = (double)(info->nz_allocated - info->nz_used);
11752d61bbb3SSatish Balay   /*  if (info->nz_unneeded != A->info.nz_unneeded)
11762d61bbb3SSatish Balay     printf("space descrepancy: maxnz-nz = %d, nz_unneeded = %d\n",(int)info->nz_unneeded,(int)A->info.nz_unneeded); */
11772d61bbb3SSatish Balay   info->assemblies   = A->num_ass;
11782d61bbb3SSatish Balay   info->mallocs      = a->reallocs;
11792d61bbb3SSatish Balay   info->memory       = A->mem;
11802d61bbb3SSatish Balay   if (A->factor) {
11812d61bbb3SSatish Balay     info->fill_ratio_given  = A->info.fill_ratio_given;
11822d61bbb3SSatish Balay     info->fill_ratio_needed = A->info.fill_ratio_needed;
11832d61bbb3SSatish Balay     info->factor_mallocs    = A->info.factor_mallocs;
11842d61bbb3SSatish Balay   } else {
11852d61bbb3SSatish Balay     info->fill_ratio_given  = 0;
11862d61bbb3SSatish Balay     info->fill_ratio_needed = 0;
11872d61bbb3SSatish Balay     info->factor_mallocs    = 0;
11882d61bbb3SSatish Balay   }
11892d61bbb3SSatish Balay   PetscFunctionReturn(0);
11902d61bbb3SSatish Balay }
11912d61bbb3SSatish Balay 
11922d61bbb3SSatish Balay 
11932d61bbb3SSatish Balay #undef __FUNC__
11942d61bbb3SSatish Balay #define __FUNC__ "MatZeroEntries_SeqBAIJ"
11952d61bbb3SSatish Balay int MatZeroEntries_SeqBAIJ(Mat A)
11962d61bbb3SSatish Balay {
11972d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data;
11982d61bbb3SSatish Balay 
11992d61bbb3SSatish Balay   PetscFunctionBegin;
12002d61bbb3SSatish Balay   PetscMemzero(a->a,a->bs2*a->i[a->mbs]*sizeof(Scalar));
12012d61bbb3SSatish Balay   PetscFunctionReturn(0);
12022d61bbb3SSatish Balay }
1203