xref: /petsc/src/mat/impls/baij/seq/baij2.c (revision 6849ba73f22fecb8f92ef896a42e4e8bd4cd6965)
1cac129eeSSatish Balay 
270f55243SBarry Smith #include "src/mat/impls/baij/seq/baij.h"
32d61bbb3SSatish Balay #include "src/inline/spops.h"
43f1db9ecSBarry Smith #include "src/inline/ilu.h"
50a835dfdSSatish Balay #include "petscbt.h"
6cac129eeSSatish Balay 
74a2ae208SSatish Balay #undef __FUNCT__
84a2ae208SSatish Balay #define __FUNCT__ "MatIncreaseOverlap_SeqBAIJ"
9dfbe8321SBarry Smith PetscErrorCode MatIncreaseOverlap_SeqBAIJ(Mat A,int is_max,IS is[],int ov)
10a3192f15SSatish Balay {
11a3192f15SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
12*6849ba73SBarry Smith   PetscErrorCode ierr;
13*6849ba73SBarry Smith   int         row,i,j,k,l,m,n,*idx,*nidx,isz,val,ival;
14218c64b6SSatish Balay   int         start,end,*ai,*aj,bs,*nidx2;
15f1af5d2fSBarry Smith   PetscBT     table;
16a3192f15SSatish Balay 
173a40ed3dSBarry Smith   PetscFunctionBegin;
18a3192f15SSatish Balay   m     = a->mbs;
19a3192f15SSatish Balay   ai    = a->i;
20a3192f15SSatish Balay   aj    = a->j;
21218c64b6SSatish Balay   bs    = a->bs;
22a3192f15SSatish Balay 
2329bbc08cSBarry Smith   if (ov < 0)  SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap specified");
24a3192f15SSatish Balay 
256831982aSBarry Smith   ierr = PetscBTCreate(m,table);CHKERRQ(ierr);
26b0a32e0cSBarry Smith   ierr = PetscMalloc((m+1)*sizeof(int),&nidx);CHKERRQ(ierr);
27b0a32e0cSBarry Smith   ierr = PetscMalloc((A->m+1)*sizeof(int),&nidx2);CHKERRQ(ierr);
28a3192f15SSatish Balay 
29a3192f15SSatish Balay   for (i=0; i<is_max; i++) {
30a3192f15SSatish Balay     /* Initialise the two local arrays */
31a3192f15SSatish Balay     isz  = 0;
326831982aSBarry Smith     ierr = PetscBTMemzero(m,table);CHKERRQ(ierr);
33a3192f15SSatish Balay 
34a3192f15SSatish Balay     /* Extract the indices, assume there can be duplicate entries */
35a3192f15SSatish Balay     ierr = ISGetIndices(is[i],&idx);CHKERRQ(ierr);
36b9b97703SBarry Smith     ierr = ISGetLocalSize(is[i],&n);CHKERRQ(ierr);
37a3192f15SSatish Balay 
38a3192f15SSatish Balay     /* Enter these into the temp arrays i.e mark table[row], enter row into new index */
39a3192f15SSatish Balay     for (j=0; j<n ; ++j){
40218c64b6SSatish Balay       ival = idx[j]/bs; /* convert the indices into block indices */
415eee224dSHong Zhang       if (ival>=m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"index greater than mat-dim");
42f1af5d2fSBarry Smith       if(!PetscBTLookupSet(table,ival)) { nidx[isz++] = ival;}
43a3192f15SSatish Balay     }
44a3192f15SSatish Balay     ierr = ISRestoreIndices(is[i],&idx);CHKERRQ(ierr);
45a3192f15SSatish Balay     ierr = ISDestroy(is[i]);CHKERRQ(ierr);
46a3192f15SSatish Balay 
47a3192f15SSatish Balay     k = 0;
48a3192f15SSatish Balay     for (j=0; j<ov; j++){ /* for each overlap*/
49a3192f15SSatish Balay       n = isz;
50a3192f15SSatish Balay       for (; k<n ; k++){ /* do only those rows in nidx[k], which are not done yet */
51a3192f15SSatish Balay         row   = nidx[k];
52a3192f15SSatish Balay         start = ai[row];
53a3192f15SSatish Balay         end   = ai[row+1];
54a3192f15SSatish Balay         for (l = start; l<end ; l++){
55a3192f15SSatish Balay           val = aj[l];
56f1af5d2fSBarry Smith           if (!PetscBTLookupSet(table,val)) {nidx[isz++] = val;}
57a3192f15SSatish Balay         }
58a3192f15SSatish Balay       }
59a3192f15SSatish Balay     }
60218c64b6SSatish Balay     /* expand the Index Set */
61218c64b6SSatish Balay     for (j=0; j<isz; j++) {
62218c64b6SSatish Balay       for (k=0; k<bs; k++)
63218c64b6SSatish Balay         nidx2[j*bs+k] = nidx[j]*bs+k;
64218c64b6SSatish Balay     }
65329f5518SBarry Smith     ierr = ISCreateGeneral(PETSC_COMM_SELF,isz*bs,nidx2,is+i);CHKERRQ(ierr);
66a3192f15SSatish Balay   }
676831982aSBarry Smith   ierr = PetscBTDestroy(table);CHKERRQ(ierr);
68606d414cSSatish Balay   ierr = PetscFree(nidx);CHKERRQ(ierr);
69606d414cSSatish Balay   ierr = PetscFree(nidx2);CHKERRQ(ierr);
703a40ed3dSBarry Smith   PetscFunctionReturn(0);
71a3192f15SSatish Balay }
721c351548SSatish Balay 
734a2ae208SSatish Balay #undef __FUNCT__
744a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqBAIJ_Private"
75dfbe8321SBarry Smith PetscErrorCode MatGetSubMatrix_SeqBAIJ_Private(Mat A,IS isrow,IS iscol,int cs,MatReuse scall,Mat *B)
76736121d4SSatish Balay {
77736121d4SSatish Balay   Mat_SeqBAIJ  *a = (Mat_SeqBAIJ*)A->data,*c;
78*6849ba73SBarry Smith   PetscErrorCode ierr;
79*6849ba73SBarry Smith   int          *smap,i,k,kstart,kend,oldcols = a->nbs,*lens;
80218c64b6SSatish Balay   int          row,mat_i,*mat_j,tcol,*mat_ilen;
81736121d4SSatish Balay   int          *irow,*icol,nrows,ncols,*ssmap,bs=a->bs,bs2=a->bs2;
82218c64b6SSatish Balay   int          *aj = a->j,*ai = a->i;
833f1db9ecSBarry Smith   MatScalar    *mat_a;
84736121d4SSatish Balay   Mat          C;
850f5bd95cSBarry Smith   PetscTruth   flag;
86736121d4SSatish Balay 
873a40ed3dSBarry Smith   PetscFunctionBegin;
88888f2ed8SSatish Balay   ierr = ISSorted(iscol,(PetscTruth*)&i);CHKERRQ(ierr);
8929bbc08cSBarry Smith   if (!i) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"IS is not sorted");
90736121d4SSatish Balay 
91736121d4SSatish Balay   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
92218c64b6SSatish Balay   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
93b9b97703SBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
94b9b97703SBarry Smith   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
95736121d4SSatish Balay 
96b0a32e0cSBarry Smith   ierr = PetscMalloc((1+oldcols)*sizeof(int),&smap);CHKERRQ(ierr);
97736121d4SSatish Balay   ssmap = smap;
98b0a32e0cSBarry Smith   ierr = PetscMalloc((1+nrows)*sizeof(int),&lens);CHKERRQ(ierr);
99549d3d68SSatish Balay   ierr  = PetscMemzero(smap,oldcols*sizeof(int));CHKERRQ(ierr);
100736121d4SSatish Balay   for (i=0; i<ncols; i++) smap[icol[i]] = i+1;
101736121d4SSatish Balay   /* determine lens of each row */
102736121d4SSatish Balay   for (i=0; i<nrows; i++) {
103736121d4SSatish Balay     kstart  = ai[irow[i]];
104736121d4SSatish Balay     kend    = kstart + a->ilen[irow[i]];
105736121d4SSatish Balay     lens[i] = 0;
106736121d4SSatish Balay       for (k=kstart; k<kend; k++) {
107736121d4SSatish Balay         if (ssmap[aj[k]]) {
108736121d4SSatish Balay           lens[i]++;
109736121d4SSatish Balay         }
110736121d4SSatish Balay       }
111736121d4SSatish Balay     }
112736121d4SSatish Balay   /* Create and fill new matrix */
113736121d4SSatish Balay   if (scall == MAT_REUSE_MATRIX) {
114736121d4SSatish Balay     c = (Mat_SeqBAIJ *)((*B)->data);
115736121d4SSatish Balay 
11629bbc08cSBarry Smith     if (c->mbs!=nrows || c->nbs!=ncols || c->bs!=bs) SETERRQ(PETSC_ERR_ARG_SIZ,"Submatrix wrong size");
1170f5bd95cSBarry Smith     ierr = PetscMemcmp(c->ilen,lens,c->mbs *sizeof(int),&flag);CHKERRQ(ierr);
1180f5bd95cSBarry Smith     if (flag == PETSC_FALSE) {
11929bbc08cSBarry Smith       SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong no of nonzeros");
120736121d4SSatish Balay     }
121549d3d68SSatish Balay     ierr = PetscMemzero(c->ilen,c->mbs*sizeof(int));CHKERRQ(ierr);
122736121d4SSatish Balay     C = *B;
1233a40ed3dSBarry Smith   } else {
124e2d9671bSKris Buschelman     ierr = MatCreate(A->comm,nrows*bs,ncols*bs,PETSC_DETERMINE,PETSC_DETERMINE,&C);CHKERRQ(ierr);
125e2d9671bSKris Buschelman     ierr = MatSetType(C,A->type_name);CHKERRQ(ierr);
126e2d9671bSKris Buschelman     ierr = MatSeqBAIJSetPreallocation(C,bs,0,lens);CHKERRQ(ierr);
127736121d4SSatish Balay   }
128736121d4SSatish Balay   c = (Mat_SeqBAIJ *)(C->data);
129736121d4SSatish Balay   for (i=0; i<nrows; i++) {
130736121d4SSatish Balay     row    = irow[i];
131736121d4SSatish Balay     kstart = ai[row];
132736121d4SSatish Balay     kend   = kstart + a->ilen[row];
133736121d4SSatish Balay     mat_i  = c->i[i];
134736121d4SSatish Balay     mat_j  = c->j + mat_i;
135218c64b6SSatish Balay     mat_a  = c->a + mat_i*bs2;
136736121d4SSatish Balay     mat_ilen = c->ilen + i;
137736121d4SSatish Balay     for (k=kstart; k<kend; k++) {
138736121d4SSatish Balay       if ((tcol=ssmap[a->j[k]])) {
139736121d4SSatish Balay         *mat_j++ = tcol - 1;
140549d3d68SSatish Balay         ierr     = PetscMemcpy(mat_a,a->a+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
141549d3d68SSatish Balay         mat_a   += bs2;
142736121d4SSatish Balay         (*mat_ilen)++;
143736121d4SSatish Balay       }
144736121d4SSatish Balay     }
145736121d4SSatish Balay   }
146218c64b6SSatish Balay 
147736121d4SSatish Balay   /* Free work space */
148736121d4SSatish Balay   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
149606d414cSSatish Balay   ierr = PetscFree(smap);CHKERRQ(ierr);
150606d414cSSatish Balay   ierr = PetscFree(lens);CHKERRQ(ierr);
1516d4a8577SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1526d4a8577SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
153736121d4SSatish Balay 
154736121d4SSatish Balay   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
155736121d4SSatish Balay   *B = C;
1563a40ed3dSBarry Smith   PetscFunctionReturn(0);
157736121d4SSatish Balay }
158736121d4SSatish Balay 
1594a2ae208SSatish Balay #undef __FUNCT__
1604a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqBAIJ"
161dfbe8321SBarry Smith PetscErrorCode MatGetSubMatrix_SeqBAIJ(Mat A,IS isrow,IS iscol,int cs,MatReuse scall,Mat *B)
162218c64b6SSatish Balay {
163218c64b6SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
164218c64b6SSatish Balay   IS          is1,is2;
165*6849ba73SBarry Smith   PetscErrorCode ierr;
166*6849ba73SBarry Smith   int         *vary,*iary,*irow,*icol,nrows,ncols,i,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);
171b9b97703SBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
172b9b97703SBarry Smith   ierr = ISGetLocalSize(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 */
176b0a32e0cSBarry Smith   ierr = PetscMalloc(2*(a->mbs+1)*sizeof(int),&vary);CHKERRQ(ierr);
177218c64b6SSatish Balay   iary = vary + a->mbs;
178549d3d68SSatish Balay   ierr = PetscMemzero(vary,(a->mbs)*sizeof(int));CHKERRQ(ierr);
179218c64b6SSatish Balay   for (i=0; i<nrows; i++) vary[irow[i]/bs]++;
180218c64b6SSatish Balay   count = 0;
181218c64b6SSatish Balay   for (i=0; i<a->mbs; i++) {
182634064b4SBarry Smith     if (vary[i]!=0 && vary[i]!=bs) SETERRQ(PETSC_ERR_ARG_SIZ,"Index set does not match blocks");
183218c64b6SSatish Balay     if (vary[i]==bs) iary[count++] = i;
184218c64b6SSatish Balay   }
185029af93fSBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,count,iary,&is1);CHKERRQ(ierr);
186218c64b6SSatish Balay 
187549d3d68SSatish Balay   ierr = PetscMemzero(vary,(a->mbs)*sizeof(int));CHKERRQ(ierr);
188218c64b6SSatish Balay   for (i=0; i<ncols; i++) vary[icol[i]/bs]++;
189218c64b6SSatish Balay   count = 0;
190218c64b6SSatish Balay   for (i=0; i<a->mbs; i++) {
191634064b4SBarry Smith     if (vary[i]!=0 && vary[i]!=bs) SETERRQ(PETSC_ERR_PLIB,"Internal error in PETSc");
192218c64b6SSatish Balay     if (vary[i]==bs) iary[count++] = i;
193218c64b6SSatish Balay   }
194029af93fSBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,count,iary,&is2);CHKERRQ(ierr);
195218c64b6SSatish Balay   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
196218c64b6SSatish Balay   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
197606d414cSSatish Balay   ierr = PetscFree(vary);CHKERRQ(ierr);
198218c64b6SSatish Balay 
1996a6a5d1dSBarry Smith   ierr = MatGetSubMatrix_SeqBAIJ_Private(A,is1,is2,cs,scall,B);CHKERRQ(ierr);
200218c64b6SSatish Balay   ISDestroy(is1);
201218c64b6SSatish Balay   ISDestroy(is2);
2023a40ed3dSBarry Smith   PetscFunctionReturn(0);
203218c64b6SSatish Balay }
204218c64b6SSatish Balay 
2054a2ae208SSatish Balay #undef __FUNCT__
2064a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_SeqBAIJ"
207dfbe8321SBarry Smith PetscErrorCode MatGetSubMatrices_SeqBAIJ(Mat A,int n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
208736121d4SSatish Balay {
209*6849ba73SBarry Smith   PetscErrorCode ierr;
210*6849ba73SBarry Smith   int i;
211736121d4SSatish Balay 
2123a40ed3dSBarry Smith   PetscFunctionBegin;
213736121d4SSatish Balay   if (scall == MAT_INITIAL_MATRIX) {
21482502324SSatish Balay     ierr = PetscMalloc((n+1)*sizeof(Mat),B);CHKERRQ(ierr);
215736121d4SSatish Balay   }
216736121d4SSatish Balay 
217736121d4SSatish Balay   for (i=0; i<n; i++) {
2186a6a5d1dSBarry Smith     ierr = MatGetSubMatrix_SeqBAIJ(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr);
219736121d4SSatish Balay   }
2203a40ed3dSBarry Smith   PetscFunctionReturn(0);
221736121d4SSatish Balay }
222218c64b6SSatish Balay 
223218c64b6SSatish Balay 
2242d61bbb3SSatish Balay /* -------------------------------------------------------*/
2252d61bbb3SSatish Balay /* Should check that shapes of vectors and matrices match */
2262d61bbb3SSatish Balay /* -------------------------------------------------------*/
227d9eff348SSatish Balay #include "petscblaslapack.h"
2282d61bbb3SSatish Balay 
2294a2ae208SSatish Balay #undef __FUNCT__
2304a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_1"
231dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_1(Mat A,Vec xx,Vec zz)
2322d61bbb3SSatish Balay {
2332d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ*)A->data;
23487828ca2SBarry Smith   PetscScalar     *x,*z,sum;
2353f1db9ecSBarry Smith   MatScalar       *v;
236*6849ba73SBarry Smith   PetscErrorCode ierr;
237*6849ba73SBarry Smith   int             mbs=a->mbs,i,*idx,*ii,n;
2382d61bbb3SSatish Balay 
2392d61bbb3SSatish Balay   PetscFunctionBegin;
2401ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2411ebc52fbSHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
2422d61bbb3SSatish Balay 
2432d61bbb3SSatish Balay   idx   = a->j;
2442d61bbb3SSatish Balay   v     = a->a;
2452d61bbb3SSatish Balay   ii    = a->i;
2462d61bbb3SSatish Balay 
2472d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
2482d61bbb3SSatish Balay     n    = ii[1] - ii[0]; ii++;
2492d61bbb3SSatish Balay     sum  = 0.0;
2502d61bbb3SSatish Balay     while (n--) sum += *v++ * x[*idx++];
2512d61bbb3SSatish Balay     z[i] = sum;
2522d61bbb3SSatish Balay   }
2531ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2541ebc52fbSHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
255b0a32e0cSBarry Smith   PetscLogFlops(2*a->nz - A->m);
2562d61bbb3SSatish Balay   PetscFunctionReturn(0);
2572d61bbb3SSatish Balay }
2582d61bbb3SSatish Balay 
2594a2ae208SSatish Balay #undef __FUNCT__
2604a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_2"
261dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_2(Mat A,Vec xx,Vec zz)
2622d61bbb3SSatish Balay {
2632d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ*)A->data;
26487828ca2SBarry Smith   PetscScalar     *x,*z,*xb,sum1,sum2;
26587828ca2SBarry Smith   PetscScalar     x1,x2;
2663f1db9ecSBarry Smith   MatScalar       *v;
267dfbe8321SBarry Smith   PetscErrorCode ierr;
268dfbe8321SBarry Smith   int mbs=a->mbs,i,*idx,*ii,j,n;
2692d61bbb3SSatish Balay 
2702d61bbb3SSatish Balay   PetscFunctionBegin;
2711ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2721ebc52fbSHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
2732d61bbb3SSatish Balay 
2742d61bbb3SSatish Balay   idx   = a->j;
2752d61bbb3SSatish Balay   v     = a->a;
2762d61bbb3SSatish Balay   ii    = a->i;
2772d61bbb3SSatish Balay 
2782d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
2792d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
2802d61bbb3SSatish Balay     sum1 = 0.0; sum2 = 0.0;
2812d61bbb3SSatish Balay     for (j=0; j<n; j++) {
2822d61bbb3SSatish Balay       xb = x + 2*(*idx++); x1 = xb[0]; x2 = xb[1];
2832d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[2]*x2;
2842d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[3]*x2;
2852d61bbb3SSatish Balay       v += 4;
2862d61bbb3SSatish Balay     }
2872d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2;
2882d61bbb3SSatish Balay     z += 2;
2892d61bbb3SSatish Balay   }
2901ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2911ebc52fbSHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
292b0a32e0cSBarry Smith   PetscLogFlops(8*a->nz - A->m);
2932d61bbb3SSatish Balay   PetscFunctionReturn(0);
2942d61bbb3SSatish Balay }
2952d61bbb3SSatish Balay 
2964a2ae208SSatish Balay #undef __FUNCT__
2974a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_3"
298dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_3(Mat A,Vec xx,Vec zz)
2992d61bbb3SSatish Balay {
3002d61bbb3SSatish Balay   Mat_SeqBAIJ  *a = (Mat_SeqBAIJ*)A->data;
30187828ca2SBarry Smith   PetscScalar  *x,*z,*xb,sum1,sum2,sum3,x1,x2,x3;
3023f1db9ecSBarry Smith   MatScalar    *v;
303dfbe8321SBarry Smith   PetscErrorCode ierr;
304dfbe8321SBarry Smith   int  mbs=a->mbs,i,*idx,*ii,j,n;
3052d61bbb3SSatish Balay 
306b6410449SSatish Balay #if defined(PETSC_HAVE_PRAGMA_DISJOINT)
307fee21e36SBarry Smith #pragma disjoint(*v,*z,*xb)
308fee21e36SBarry Smith #endif
309fee21e36SBarry Smith 
3102d61bbb3SSatish Balay   PetscFunctionBegin;
3111ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
3121ebc52fbSHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
3132d61bbb3SSatish Balay 
3142d61bbb3SSatish Balay   idx   = a->j;
3152d61bbb3SSatish Balay   v     = a->a;
3162d61bbb3SSatish Balay   ii    = a->i;
3172d61bbb3SSatish Balay 
3182d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
3192d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
3202d61bbb3SSatish Balay     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0;
3212d61bbb3SSatish Balay     for (j=0; j<n; j++) {
3222d61bbb3SSatish Balay       xb = x + 3*(*idx++); x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
3232d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[3]*x2 + v[6]*x3;
3242d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[4]*x2 + v[7]*x3;
3252d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[5]*x2 + v[8]*x3;
3262d61bbb3SSatish Balay       v += 9;
3272d61bbb3SSatish Balay     }
3282d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3;
3292d61bbb3SSatish Balay     z += 3;
3302d61bbb3SSatish Balay   }
3311ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
3321ebc52fbSHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
333b0a32e0cSBarry Smith   PetscLogFlops(18*a->nz - A->m);
3342d61bbb3SSatish Balay   PetscFunctionReturn(0);
3352d61bbb3SSatish Balay }
3362d61bbb3SSatish Balay 
3374a2ae208SSatish Balay #undef __FUNCT__
3384a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_4"
339dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_4(Mat A,Vec xx,Vec zz)
3402d61bbb3SSatish Balay {
3412d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ*)A->data;
34287828ca2SBarry Smith   PetscScalar     *x,*z,*xb,sum1,sum2,sum3,sum4,x1,x2,x3,x4;
3433f1db9ecSBarry Smith   MatScalar       *v;
344dfbe8321SBarry Smith   PetscErrorCode ierr;
345dfbe8321SBarry Smith   int mbs=a->mbs,i,*idx,*ii,j,n;
3462d61bbb3SSatish Balay 
3472d61bbb3SSatish Balay   PetscFunctionBegin;
3481ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
3491ebc52fbSHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
3502d61bbb3SSatish Balay 
3512d61bbb3SSatish Balay   idx   = a->j;
3522d61bbb3SSatish Balay   v     = a->a;
3532d61bbb3SSatish Balay   ii    = a->i;
3542d61bbb3SSatish Balay 
3552d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
3562d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
3572d61bbb3SSatish Balay     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0;
3582d61bbb3SSatish Balay     for (j=0; j<n; j++) {
3592d61bbb3SSatish Balay       xb = x + 4*(*idx++);
3602d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
3612d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[4]*x2 + v[8]*x3   + v[12]*x4;
3622d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[5]*x2 + v[9]*x3   + v[13]*x4;
3632d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[6]*x2 + v[10]*x3  + v[14]*x4;
3642d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[7]*x2 + v[11]*x3  + v[15]*x4;
3652d61bbb3SSatish Balay       v += 16;
3662d61bbb3SSatish Balay     }
3672d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4;
3682d61bbb3SSatish Balay     z += 4;
3692d61bbb3SSatish Balay   }
3701ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
3711ebc52fbSHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
372b0a32e0cSBarry Smith   PetscLogFlops(32*a->nz - A->m);
3732d61bbb3SSatish Balay   PetscFunctionReturn(0);
3742d61bbb3SSatish Balay }
3752d61bbb3SSatish Balay 
3764a2ae208SSatish Balay #undef __FUNCT__
3774a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_5"
378dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_5(Mat A,Vec xx,Vec zz)
3792d61bbb3SSatish Balay {
3802d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ*)A->data;
38187828ca2SBarry Smith   PetscScalar     sum1,sum2,sum3,sum4,sum5,x1,x2,x3,x4,x5,*xb,*z,*x;
3823f1db9ecSBarry Smith   MatScalar       *v;
383dfbe8321SBarry Smith   PetscErrorCode ierr;
384dfbe8321SBarry Smith   int mbs=a->mbs,i,*idx,*ii,j,n;
3852d61bbb3SSatish Balay 
386433994e6SBarry Smith   PetscFunctionBegin;
3871ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
3881ebc52fbSHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
3892d61bbb3SSatish Balay 
3902d61bbb3SSatish Balay   idx   = a->j;
3912d61bbb3SSatish Balay   v     = a->a;
3922d61bbb3SSatish Balay   ii    = a->i;
3932d61bbb3SSatish Balay 
3942d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
3952d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
3962d61bbb3SSatish Balay     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0;
3972d61bbb3SSatish Balay     for (j=0; j<n; j++) {
3982d61bbb3SSatish Balay       xb = x + 5*(*idx++);
3992d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4];
4002d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[5]*x2 + v[10]*x3  + v[15]*x4 + v[20]*x5;
4012d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[6]*x2 + v[11]*x3  + v[16]*x4 + v[21]*x5;
4022d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[7]*x2 + v[12]*x3  + v[17]*x4 + v[22]*x5;
4032d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[8]*x2 + v[13]*x3  + v[18]*x4 + v[23]*x5;
4042d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[9]*x2 + v[14]*x3  + v[19]*x4 + v[24]*x5;
4052d61bbb3SSatish Balay       v += 25;
4062d61bbb3SSatish Balay     }
4072d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5;
4082d61bbb3SSatish Balay     z += 5;
4092d61bbb3SSatish Balay   }
4101ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
4111ebc52fbSHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
412b0a32e0cSBarry Smith   PetscLogFlops(50*a->nz - A->m);
4132d61bbb3SSatish Balay   PetscFunctionReturn(0);
4142d61bbb3SSatish Balay }
4152d61bbb3SSatish Balay 
41615091d37SBarry Smith 
4174a2ae208SSatish Balay #undef __FUNCT__
4184a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_6"
419dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_6(Mat A,Vec xx,Vec zz)
42015091d37SBarry Smith {
42115091d37SBarry Smith   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ*)A->data;
42287828ca2SBarry Smith   PetscScalar     *x,*z,*xb,sum1,sum2,sum3,sum4,sum5,sum6;
42387828ca2SBarry Smith   PetscScalar     x1,x2,x3,x4,x5,x6;
42415091d37SBarry Smith   MatScalar       *v;
425dfbe8321SBarry Smith   PetscErrorCode ierr;
426dfbe8321SBarry Smith   int   mbs=a->mbs,i,*idx,*ii,j,n;
42715091d37SBarry Smith 
428433994e6SBarry Smith   PetscFunctionBegin;
4291ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
4301ebc52fbSHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
43115091d37SBarry Smith 
43215091d37SBarry Smith   idx   = a->j;
43315091d37SBarry Smith   v     = a->a;
43415091d37SBarry Smith   ii    = a->i;
43515091d37SBarry Smith 
43615091d37SBarry Smith   for (i=0; i<mbs; i++) {
43715091d37SBarry Smith     n  = ii[1] - ii[0]; ii++;
43815091d37SBarry Smith     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0;
43915091d37SBarry Smith     for (j=0; j<n; j++) {
44015091d37SBarry Smith       xb = x + 6*(*idx++);
44115091d37SBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5];
44215091d37SBarry Smith       sum1 += v[0]*x1 + v[6]*x2  + v[12]*x3  + v[18]*x4 + v[24]*x5 + v[30]*x6;
44315091d37SBarry Smith       sum2 += v[1]*x1 + v[7]*x2  + v[13]*x3  + v[19]*x4 + v[25]*x5 + v[31]*x6;
44415091d37SBarry Smith       sum3 += v[2]*x1 + v[8]*x2  + v[14]*x3  + v[20]*x4 + v[26]*x5 + v[32]*x6;
44515091d37SBarry Smith       sum4 += v[3]*x1 + v[9]*x2  + v[15]*x3  + v[21]*x4 + v[27]*x5 + v[33]*x6;
44615091d37SBarry Smith       sum5 += v[4]*x1 + v[10]*x2 + v[16]*x3  + v[22]*x4 + v[28]*x5 + v[34]*x6;
44715091d37SBarry Smith       sum6 += v[5]*x1 + v[11]*x2 + v[17]*x3  + v[23]*x4 + v[29]*x5 + v[35]*x6;
44815091d37SBarry Smith       v += 36;
44915091d37SBarry Smith     }
45015091d37SBarry Smith     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6;
45115091d37SBarry Smith     z += 6;
45215091d37SBarry Smith   }
45315091d37SBarry Smith 
4541ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
4551ebc52fbSHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
456b0a32e0cSBarry Smith   PetscLogFlops(72*a->nz - A->m);
45715091d37SBarry Smith   PetscFunctionReturn(0);
45815091d37SBarry Smith }
4594a2ae208SSatish Balay #undef __FUNCT__
4604a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_7"
461dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_7(Mat A,Vec xx,Vec zz)
4622d61bbb3SSatish Balay {
4632d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ*)A->data;
46487828ca2SBarry Smith   PetscScalar     *x,*z,*xb,sum1,sum2,sum3,sum4,sum5,sum6,sum7;
46587828ca2SBarry Smith   PetscScalar     x1,x2,x3,x4,x5,x6,x7;
4663f1db9ecSBarry Smith   MatScalar       *v;
467dfbe8321SBarry Smith   PetscErrorCode ierr;
468dfbe8321SBarry Smith   int mbs=a->mbs,i,*idx,*ii,j,n;
4692d61bbb3SSatish Balay 
470433994e6SBarry Smith   PetscFunctionBegin;
4711ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
4721ebc52fbSHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
4732d61bbb3SSatish Balay 
4742d61bbb3SSatish Balay   idx   = a->j;
4752d61bbb3SSatish Balay   v     = a->a;
4762d61bbb3SSatish Balay   ii    = a->i;
4772d61bbb3SSatish Balay 
4782d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
4792d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
4802d61bbb3SSatish Balay     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0; sum7 = 0.0;
4812d61bbb3SSatish Balay     for (j=0; j<n; j++) {
4822d61bbb3SSatish Balay       xb = x + 7*(*idx++);
4832d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6];
4842d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[7]*x2  + v[14]*x3  + v[21]*x4 + v[28]*x5 + v[35]*x6 + v[42]*x7;
4852d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[8]*x2  + v[15]*x3  + v[22]*x4 + v[29]*x5 + v[36]*x6 + v[43]*x7;
4862d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[9]*x2  + v[16]*x3  + v[23]*x4 + v[30]*x5 + v[37]*x6 + v[44]*x7;
4872d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[10]*x2 + v[17]*x3  + v[24]*x4 + v[31]*x5 + v[38]*x6 + v[45]*x7;
4882d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[11]*x2 + v[18]*x3  + v[25]*x4 + v[32]*x5 + v[39]*x6 + v[46]*x7;
4892d61bbb3SSatish Balay       sum6 += v[5]*x1 + v[12]*x2 + v[19]*x3  + v[26]*x4 + v[33]*x5 + v[40]*x6 + v[47]*x7;
4902d61bbb3SSatish Balay       sum7 += v[6]*x1 + v[13]*x2 + v[20]*x3  + v[27]*x4 + v[34]*x5 + v[41]*x6 + v[48]*x7;
4912d61bbb3SSatish Balay       v += 49;
4922d61bbb3SSatish Balay     }
4932d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
4942d61bbb3SSatish Balay     z += 7;
4952d61bbb3SSatish Balay   }
4962d61bbb3SSatish Balay 
4971ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
4981ebc52fbSHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
499b0a32e0cSBarry Smith   PetscLogFlops(98*a->nz - A->m);
5002d61bbb3SSatish Balay   PetscFunctionReturn(0);
5012d61bbb3SSatish Balay }
5022d61bbb3SSatish Balay 
5033f1db9ecSBarry Smith /*
5043f1db9ecSBarry Smith     This will not work with MatScalar == float because it calls the BLAS
5053f1db9ecSBarry Smith */
5064a2ae208SSatish Balay #undef __FUNCT__
5074a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_N"
508dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_N(Mat A,Vec xx,Vec zz)
5092d61bbb3SSatish Balay {
5102d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ*)A->data;
51187828ca2SBarry Smith   PetscScalar     *x,*z,*xb,*work,*workt;
5123f1db9ecSBarry Smith   MatScalar       *v;
513dfbe8321SBarry Smith   PetscErrorCode ierr;
514dfbe8321SBarry Smith   int mbs=a->mbs,i,*idx,*ii,bs=a->bs,j,n,bs2=a->bs2;
5153f1db9ecSBarry Smith   int             ncols,k;
5162d61bbb3SSatish Balay 
5172d61bbb3SSatish Balay   PetscFunctionBegin;
5181ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
5191ebc52fbSHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
5202d61bbb3SSatish Balay 
5212d61bbb3SSatish Balay   idx   = a->j;
5222d61bbb3SSatish Balay   v     = a->a;
5232d61bbb3SSatish Balay   ii    = a->i;
524218c64b6SSatish Balay 
525218c64b6SSatish Balay 
5262d61bbb3SSatish Balay   if (!a->mult_work) {
527273d9f13SBarry Smith     k    = PetscMax(A->m,A->n);
52887828ca2SBarry Smith     ierr = PetscMalloc((k+1)*sizeof(PetscScalar),&a->mult_work);CHKERRQ(ierr);
5292d61bbb3SSatish Balay   }
5302d61bbb3SSatish Balay   work = a->mult_work;
5312d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
5322d61bbb3SSatish Balay     n     = ii[1] - ii[0]; ii++;
5332d61bbb3SSatish Balay     ncols = n*bs;
5342d61bbb3SSatish Balay     workt = work;
5352d61bbb3SSatish Balay     for (j=0; j<n; j++) {
5362d61bbb3SSatish Balay       xb = x + bs*(*idx++);
5372d61bbb3SSatish Balay       for (k=0; k<bs; k++) workt[k] = xb[k];
5382d61bbb3SSatish Balay       workt += bs;
5392d61bbb3SSatish Balay     }
5403f1db9ecSBarry Smith     Kernel_w_gets_Ar_times_v(bs,ncols,work,v,z);
5413f1db9ecSBarry Smith     /* LAgemv_("N",&bs,&ncols,&_DOne,v,&bs,work,&_One,&_DZero,z,&_One); */
5422d61bbb3SSatish Balay     v += n*bs2;
5432d61bbb3SSatish Balay     z += bs;
5442d61bbb3SSatish Balay   }
5451ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
5461ebc52fbSHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
547b0a32e0cSBarry Smith   PetscLogFlops(2*a->nz*bs2 - A->m);
5482d61bbb3SSatish Balay   PetscFunctionReturn(0);
5492d61bbb3SSatish Balay }
5502d61bbb3SSatish Balay 
5514a2ae208SSatish Balay #undef __FUNCT__
5524a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_1"
553dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_1(Mat A,Vec xx,Vec yy,Vec zz)
5542d61bbb3SSatish Balay {
5552d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ*)A->data;
55687828ca2SBarry Smith   PetscScalar     *x,*y,*z,sum;
5573f1db9ecSBarry Smith   MatScalar       *v;
558dfbe8321SBarry Smith   PetscErrorCode ierr;
559dfbe8321SBarry Smith   int mbs=a->mbs,i,*idx,*ii,n;
5602d61bbb3SSatish Balay 
5612d61bbb3SSatish Balay   PetscFunctionBegin;
5621ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
5631ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
5642e8a6d31SBarry Smith   if (zz != yy) {
5651ebc52fbSHong Zhang     ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
5662e8a6d31SBarry Smith   } else {
5672e8a6d31SBarry Smith     z = y;
5682e8a6d31SBarry Smith   }
5692d61bbb3SSatish Balay 
5702d61bbb3SSatish Balay   idx   = a->j;
5712d61bbb3SSatish Balay   v     = a->a;
5722d61bbb3SSatish Balay   ii    = a->i;
5732d61bbb3SSatish Balay 
5742d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
5752d61bbb3SSatish Balay     n    = ii[1] - ii[0]; ii++;
5762d61bbb3SSatish Balay     sum  = y[i];
5772d61bbb3SSatish Balay     while (n--) sum += *v++ * x[*idx++];
5782d61bbb3SSatish Balay     z[i] = sum;
5792d61bbb3SSatish Balay   }
5801ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
5811ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
5822e8a6d31SBarry Smith   if (zz != yy) {
5831ebc52fbSHong Zhang     ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
5842e8a6d31SBarry Smith   }
585b0a32e0cSBarry Smith   PetscLogFlops(2*a->nz);
5862d61bbb3SSatish Balay   PetscFunctionReturn(0);
5872d61bbb3SSatish Balay }
5882d61bbb3SSatish Balay 
5894a2ae208SSatish Balay #undef __FUNCT__
5904a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_2"
591dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_2(Mat A,Vec xx,Vec yy,Vec zz)
5922d61bbb3SSatish Balay {
5932d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ*)A->data;
59487828ca2SBarry Smith   PetscScalar     *x,*y,*z,*xb,sum1,sum2;
59587828ca2SBarry Smith   PetscScalar     x1,x2;
5963f1db9ecSBarry Smith   MatScalar       *v;
597dfbe8321SBarry Smith   PetscErrorCode ierr;
598dfbe8321SBarry Smith   int mbs=a->mbs,i,*idx,*ii,j,n;
5992d61bbb3SSatish Balay 
6002d61bbb3SSatish Balay   PetscFunctionBegin;
6011ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
6021ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
6032e8a6d31SBarry Smith   if (zz != yy) {
6041ebc52fbSHong Zhang     ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
6052e8a6d31SBarry Smith   } else {
6062e8a6d31SBarry Smith     z = y;
6072e8a6d31SBarry Smith   }
6082d61bbb3SSatish Balay 
6092d61bbb3SSatish Balay   idx   = a->j;
6102d61bbb3SSatish Balay   v     = a->a;
6112d61bbb3SSatish Balay   ii    = a->i;
6122d61bbb3SSatish Balay 
6132d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
6142d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
6152d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1];
6162d61bbb3SSatish Balay     for (j=0; j<n; j++) {
6172d61bbb3SSatish Balay       xb = x + 2*(*idx++); x1 = xb[0]; x2 = xb[1];
6182d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[2]*x2;
6192d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[3]*x2;
6202d61bbb3SSatish Balay       v += 4;
6212d61bbb3SSatish Balay     }
6222d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2;
6232d61bbb3SSatish Balay     z += 2; y += 2;
6242d61bbb3SSatish Balay   }
6251ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
6261ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
6272e8a6d31SBarry Smith   if (zz != yy) {
6281ebc52fbSHong Zhang     ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
6292e8a6d31SBarry Smith   }
630b0a32e0cSBarry Smith   PetscLogFlops(4*a->nz);
6312d61bbb3SSatish Balay   PetscFunctionReturn(0);
6322d61bbb3SSatish Balay }
6332d61bbb3SSatish Balay 
6344a2ae208SSatish Balay #undef __FUNCT__
6354a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_3"
636dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_3(Mat A,Vec xx,Vec yy,Vec zz)
6372d61bbb3SSatish Balay {
6382d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ*)A->data;
63987828ca2SBarry Smith   PetscScalar     *x,*y,*z,*xb,sum1,sum2,sum3,x1,x2,x3;
6403f1db9ecSBarry Smith   MatScalar       *v;
641dfbe8321SBarry Smith   PetscErrorCode ierr;
642dfbe8321SBarry Smith   int mbs=a->mbs,i,*idx,*ii,j,n;
6432d61bbb3SSatish Balay 
6442d61bbb3SSatish Balay   PetscFunctionBegin;
6451ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
6461ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
6472e8a6d31SBarry Smith   if (zz != yy) {
6481ebc52fbSHong Zhang     ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
6492e8a6d31SBarry Smith   } else {
6502e8a6d31SBarry Smith     z = y;
6512e8a6d31SBarry Smith   }
6522d61bbb3SSatish Balay 
6532d61bbb3SSatish Balay   idx   = a->j;
6542d61bbb3SSatish Balay   v     = a->a;
6552d61bbb3SSatish Balay   ii    = a->i;
6562d61bbb3SSatish Balay 
6572d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
6582d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
6592d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2];
6602d61bbb3SSatish Balay     for (j=0; j<n; j++) {
6612d61bbb3SSatish Balay       xb = x + 3*(*idx++); x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
6622d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[3]*x2 + v[6]*x3;
6632d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[4]*x2 + v[7]*x3;
6642d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[5]*x2 + v[8]*x3;
6652d61bbb3SSatish Balay       v += 9;
6662d61bbb3SSatish Balay     }
6672d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3;
6682d61bbb3SSatish Balay     z += 3; y += 3;
6692d61bbb3SSatish Balay   }
6701ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
6711ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
6722e8a6d31SBarry Smith   if (zz != yy) {
6731ebc52fbSHong Zhang     ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
6742e8a6d31SBarry Smith   }
675b0a32e0cSBarry Smith   PetscLogFlops(18*a->nz);
6762d61bbb3SSatish Balay   PetscFunctionReturn(0);
6772d61bbb3SSatish Balay }
6782d61bbb3SSatish Balay 
6794a2ae208SSatish Balay #undef __FUNCT__
6804a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_4"
681dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_4(Mat A,Vec xx,Vec yy,Vec zz)
6822d61bbb3SSatish Balay {
6832d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ*)A->data;
68487828ca2SBarry Smith   PetscScalar     *x,*y,*z,*xb,sum1,sum2,sum3,sum4,x1,x2,x3,x4;
6853f1db9ecSBarry Smith   MatScalar       *v;
686dfbe8321SBarry Smith   PetscErrorCode ierr;
687dfbe8321SBarry Smith   int  mbs=a->mbs,i,*idx,*ii;
6882d61bbb3SSatish Balay   int             j,n;
6892d61bbb3SSatish Balay 
6902d61bbb3SSatish Balay   PetscFunctionBegin;
6911ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
6921ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
6932e8a6d31SBarry Smith   if (zz != yy) {
6941ebc52fbSHong Zhang     ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
6952e8a6d31SBarry Smith   } else {
6962e8a6d31SBarry Smith     z = y;
6972e8a6d31SBarry Smith   }
6982d61bbb3SSatish Balay 
6992d61bbb3SSatish Balay   idx   = a->j;
7002d61bbb3SSatish Balay   v     = a->a;
7012d61bbb3SSatish Balay   ii    = a->i;
7022d61bbb3SSatish Balay 
7032d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
7042d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
7052d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3];
7062d61bbb3SSatish Balay     for (j=0; j<n; j++) {
7072d61bbb3SSatish Balay       xb = x + 4*(*idx++);
7082d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
7092d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[4]*x2 + v[8]*x3   + v[12]*x4;
7102d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[5]*x2 + v[9]*x3   + v[13]*x4;
7112d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[6]*x2 + v[10]*x3  + v[14]*x4;
7122d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[7]*x2 + v[11]*x3  + v[15]*x4;
7132d61bbb3SSatish Balay       v += 16;
7142d61bbb3SSatish Balay     }
7152d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4;
7162d61bbb3SSatish Balay     z += 4; y += 4;
7172d61bbb3SSatish Balay   }
7181ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
7191ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
7202e8a6d31SBarry Smith   if (zz != yy) {
7211ebc52fbSHong Zhang     ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
7222e8a6d31SBarry Smith   }
723b0a32e0cSBarry Smith   PetscLogFlops(32*a->nz);
7242d61bbb3SSatish Balay   PetscFunctionReturn(0);
7252d61bbb3SSatish Balay }
7262d61bbb3SSatish Balay 
7274a2ae208SSatish Balay #undef __FUNCT__
7284a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_5"
729dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_5(Mat A,Vec xx,Vec yy,Vec zz)
7302d61bbb3SSatish Balay {
7312d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ*)A->data;
73287828ca2SBarry Smith   PetscScalar     *x,*y,*z,*xb,sum1,sum2,sum3,sum4,sum5,x1,x2,x3,x4,x5;
7333f1db9ecSBarry Smith   MatScalar       *v;
734dfbe8321SBarry Smith   PetscErrorCode ierr;
735dfbe8321SBarry Smith   int mbs=a->mbs,i,*idx,*ii,j,n;
7362d61bbb3SSatish Balay 
7372d61bbb3SSatish Balay   PetscFunctionBegin;
7381ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
7391ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
7402e8a6d31SBarry Smith   if (zz != yy) {
7411ebc52fbSHong Zhang     ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
7422e8a6d31SBarry Smith   } else {
7432e8a6d31SBarry Smith     z = y;
7442e8a6d31SBarry Smith   }
7452d61bbb3SSatish Balay 
7462d61bbb3SSatish Balay   idx   = a->j;
7472d61bbb3SSatish Balay   v     = a->a;
7482d61bbb3SSatish Balay   ii    = a->i;
7492d61bbb3SSatish Balay 
7502d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
7512d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
7522d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4];
7532d61bbb3SSatish Balay     for (j=0; j<n; j++) {
7542d61bbb3SSatish Balay       xb = x + 5*(*idx++);
7552d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4];
7562d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[5]*x2 + v[10]*x3  + v[15]*x4 + v[20]*x5;
7572d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[6]*x2 + v[11]*x3  + v[16]*x4 + v[21]*x5;
7582d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[7]*x2 + v[12]*x3  + v[17]*x4 + v[22]*x5;
7592d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[8]*x2 + v[13]*x3  + v[18]*x4 + v[23]*x5;
7602d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[9]*x2 + v[14]*x3  + v[19]*x4 + v[24]*x5;
7612d61bbb3SSatish Balay       v += 25;
7622d61bbb3SSatish Balay     }
7632d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5;
7642d61bbb3SSatish Balay     z += 5; y += 5;
7652d61bbb3SSatish Balay   }
7661ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
7671ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
7682e8a6d31SBarry Smith   if (zz != yy) {
7691ebc52fbSHong Zhang     ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
7702e8a6d31SBarry Smith   }
771b0a32e0cSBarry Smith   PetscLogFlops(50*a->nz);
7722d61bbb3SSatish Balay   PetscFunctionReturn(0);
7732d61bbb3SSatish Balay }
7744a2ae208SSatish Balay #undef __FUNCT__
7754a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_6"
776dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_6(Mat A,Vec xx,Vec yy,Vec zz)
77715091d37SBarry Smith {
77815091d37SBarry Smith   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ*)A->data;
77987828ca2SBarry Smith   PetscScalar     *x,*y,*z,*xb,sum1,sum2,sum3,sum4,sum5,sum6;
78087828ca2SBarry Smith   PetscScalar     x1,x2,x3,x4,x5,x6;
78115091d37SBarry Smith   MatScalar       *v;
782dfbe8321SBarry Smith   PetscErrorCode ierr;
783dfbe8321SBarry Smith   int mbs=a->mbs,i,*idx,*ii,j,n;
78415091d37SBarry Smith 
78515091d37SBarry Smith   PetscFunctionBegin;
7861ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
7871ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
78815091d37SBarry Smith   if (zz != yy) {
7891ebc52fbSHong Zhang     ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
79015091d37SBarry Smith   } else {
79115091d37SBarry Smith     z = y;
79215091d37SBarry Smith   }
79315091d37SBarry Smith 
79415091d37SBarry Smith   idx   = a->j;
79515091d37SBarry Smith   v     = a->a;
79615091d37SBarry Smith   ii    = a->i;
79715091d37SBarry Smith 
79815091d37SBarry Smith   for (i=0; i<mbs; i++) {
79915091d37SBarry Smith     n  = ii[1] - ii[0]; ii++;
80015091d37SBarry Smith     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4]; sum6 = y[5];
80115091d37SBarry Smith     for (j=0; j<n; j++) {
8023b95cb0eSSatish Balay       xb = x + 6*(*idx++);
80315091d37SBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5];
80415091d37SBarry Smith       sum1 += v[0]*x1 + v[6]*x2  + v[12]*x3  + v[18]*x4 + v[24]*x5 + v[30]*x6;
80515091d37SBarry Smith       sum2 += v[1]*x1 + v[7]*x2  + v[13]*x3  + v[19]*x4 + v[25]*x5 + v[31]*x6;
80615091d37SBarry Smith       sum3 += v[2]*x1 + v[8]*x2  + v[14]*x3  + v[20]*x4 + v[26]*x5 + v[32]*x6;
80715091d37SBarry Smith       sum4 += v[3]*x1 + v[9]*x2  + v[15]*x3  + v[21]*x4 + v[27]*x5 + v[33]*x6;
80815091d37SBarry Smith       sum5 += v[4]*x1 + v[10]*x2 + v[16]*x3  + v[22]*x4 + v[28]*x5 + v[34]*x6;
80915091d37SBarry Smith       sum6 += v[5]*x1 + v[11]*x2 + v[17]*x3  + v[23]*x4 + v[29]*x5 + v[35]*x6;
81015091d37SBarry Smith       v += 36;
81115091d37SBarry Smith     }
81215091d37SBarry Smith     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6;
81315091d37SBarry Smith     z += 6; y += 6;
81415091d37SBarry Smith   }
8151ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
8161ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
81715091d37SBarry Smith   if (zz != yy) {
8181ebc52fbSHong Zhang     ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
81915091d37SBarry Smith   }
820b0a32e0cSBarry Smith   PetscLogFlops(72*a->nz);
82115091d37SBarry Smith   PetscFunctionReturn(0);
82215091d37SBarry Smith }
8232d61bbb3SSatish Balay 
8244a2ae208SSatish Balay #undef __FUNCT__
8254a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_7"
826dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_7(Mat A,Vec xx,Vec yy,Vec zz)
8272d61bbb3SSatish Balay {
8282d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ*)A->data;
82987828ca2SBarry Smith   PetscScalar     *x,*y,*z,*xb,sum1,sum2,sum3,sum4,sum5,sum6,sum7;
83087828ca2SBarry Smith   PetscScalar     x1,x2,x3,x4,x5,x6,x7;
8313f1db9ecSBarry Smith   MatScalar       *v;
832dfbe8321SBarry Smith   PetscErrorCode ierr;
833dfbe8321SBarry Smith   int mbs=a->mbs,i,*idx,*ii,j,n;
8342d61bbb3SSatish Balay 
8352d61bbb3SSatish Balay   PetscFunctionBegin;
8361ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
8371ebc52fbSHong Zhang   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
8382e8a6d31SBarry Smith   if (zz != yy) {
8391ebc52fbSHong Zhang     ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
8402e8a6d31SBarry Smith   } else {
8412e8a6d31SBarry Smith     z = y;
8422e8a6d31SBarry Smith   }
8432d61bbb3SSatish Balay 
8442d61bbb3SSatish Balay   idx   = a->j;
8452d61bbb3SSatish Balay   v     = a->a;
8462d61bbb3SSatish Balay   ii    = a->i;
8472d61bbb3SSatish Balay 
8482d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
8492d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
8502d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4]; sum6 = y[5]; sum7 = y[6];
8512d61bbb3SSatish Balay     for (j=0; j<n; j++) {
8522d61bbb3SSatish Balay       xb = x + 7*(*idx++);
8532d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6];
8542d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[7]*x2  + v[14]*x3  + v[21]*x4 + v[28]*x5 + v[35]*x6 + v[42]*x7;
8552d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[8]*x2  + v[15]*x3  + v[22]*x4 + v[29]*x5 + v[36]*x6 + v[43]*x7;
8562d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[9]*x2  + v[16]*x3  + v[23]*x4 + v[30]*x5 + v[37]*x6 + v[44]*x7;
8572d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[10]*x2 + v[17]*x3  + v[24]*x4 + v[31]*x5 + v[38]*x6 + v[45]*x7;
8582d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[11]*x2 + v[18]*x3  + v[25]*x4 + v[32]*x5 + v[39]*x6 + v[46]*x7;
8592d61bbb3SSatish Balay       sum6 += v[5]*x1 + v[12]*x2 + v[19]*x3  + v[26]*x4 + v[33]*x5 + v[40]*x6 + v[47]*x7;
8602d61bbb3SSatish Balay       sum7 += v[6]*x1 + v[13]*x2 + v[20]*x3  + v[27]*x4 + v[34]*x5 + v[41]*x6 + v[48]*x7;
8612d61bbb3SSatish Balay       v += 49;
8622d61bbb3SSatish Balay     }
8632d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
8642d61bbb3SSatish Balay     z += 7; y += 7;
8652d61bbb3SSatish Balay   }
8661ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
8671ebc52fbSHong Zhang   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
8682e8a6d31SBarry Smith   if (zz != yy) {
8691ebc52fbSHong Zhang     ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
8702e8a6d31SBarry Smith   }
871b0a32e0cSBarry Smith   PetscLogFlops(98*a->nz);
8722d61bbb3SSatish Balay   PetscFunctionReturn(0);
8732d61bbb3SSatish Balay }
874218c64b6SSatish Balay 
8754a2ae208SSatish Balay #undef __FUNCT__
8764a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_N"
877dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_N(Mat A,Vec xx,Vec yy,Vec zz)
8782d61bbb3SSatish Balay {
8792d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
880b2e7d493SHong Zhang   PetscScalar    *x,*z,*xb,*work,*workt,*y;
8813f1db9ecSBarry Smith   MatScalar      *v;
882*6849ba73SBarry Smith   PetscErrorCode ierr;
883*6849ba73SBarry Smith   int            mbs=a->mbs,i,*idx,*ii,bs=a->bs,j,n,bs2=a->bs2;
8843f1db9ecSBarry Smith   int            ncols,k;
885218c64b6SSatish Balay 
8862d61bbb3SSatish Balay   PetscFunctionBegin;
8871ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
8881ebc52fbSHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
889b2e7d493SHong Zhang   if (zz != yy) {
8901ebc52fbSHong Zhang     ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
891b2e7d493SHong Zhang     ierr = PetscMemcpy(z,y,yy->n*sizeof(PetscScalar));CHKERRQ(ierr);
8921ebc52fbSHong Zhang     ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
893b2e7d493SHong Zhang   }
8942d61bbb3SSatish Balay 
8952d61bbb3SSatish Balay   idx   = a->j;
8962d61bbb3SSatish Balay   v     = a->a;
8972d61bbb3SSatish Balay   ii    = a->i;
8982d61bbb3SSatish Balay 
8992d61bbb3SSatish Balay 
9002d61bbb3SSatish Balay   if (!a->mult_work) {
901273d9f13SBarry Smith     k    = PetscMax(A->m,A->n);
90287828ca2SBarry Smith     ierr = PetscMalloc((k+1)*sizeof(PetscScalar),&a->mult_work);CHKERRQ(ierr);
9032d61bbb3SSatish Balay   }
9042d61bbb3SSatish Balay   work = a->mult_work;
9052d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
9062d61bbb3SSatish Balay     n     = ii[1] - ii[0]; ii++;
9072d61bbb3SSatish Balay     ncols = n*bs;
9082d61bbb3SSatish Balay     workt = work;
9092d61bbb3SSatish Balay     for (j=0; j<n; j++) {
9102d61bbb3SSatish Balay       xb = x + bs*(*idx++);
9112d61bbb3SSatish Balay       for (k=0; k<bs; k++) workt[k] = xb[k];
9122d61bbb3SSatish Balay       workt += bs;
9132d61bbb3SSatish Balay     }
9143f1db9ecSBarry Smith     Kernel_w_gets_w_plus_Ar_times_v(bs,ncols,work,v,z);
9153f1db9ecSBarry Smith     /* LAgemv_("N",&bs,&ncols,&_DOne,v,&bs,work,&_One,&_DOne,z,&_One); */
9162d61bbb3SSatish Balay     v += n*bs2;
9172d61bbb3SSatish Balay     z += bs;
9182d61bbb3SSatish Balay   }
9191ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
9201ebc52fbSHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
921b0a32e0cSBarry Smith   PetscLogFlops(2*a->nz*bs2);
9222d61bbb3SSatish Balay   PetscFunctionReturn(0);
9232d61bbb3SSatish Balay }
9242d61bbb3SSatish Balay 
9254a2ae208SSatish Balay #undef __FUNCT__
9264a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_SeqBAIJ"
927dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_SeqBAIJ(Mat A,Vec xx,Vec zz)
9282d61bbb3SSatish Balay {
9292d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ*)A->data;
93087828ca2SBarry Smith   PetscScalar     *xg,*zg,*zb,zero = 0.0;
93187828ca2SBarry Smith   PetscScalar     *x,*z,*xb,x1,x2,x3,x4,x5,x6,x7;
9323f1db9ecSBarry Smith   MatScalar       *v;
933f1af5d2fSBarry Smith   int             mbs=a->mbs,i,*idx,*ii,*ai=a->i,rval;
934*6849ba73SBarry Smith   PetscErrorCode ierr;
935*6849ba73SBarry Smith   int             bs=a->bs,j,n,bs2=a->bs2,*ib;
9362d61bbb3SSatish Balay 
9372d61bbb3SSatish Balay   PetscFunctionBegin;
938f1af5d2fSBarry Smith   ierr = VecSet(&zero,zz);CHKERRQ(ierr);
9391ebc52fbSHong Zhang   ierr = VecGetArray(xx,&xg);CHKERRQ(ierr); x = xg;
9401ebc52fbSHong Zhang   ierr = VecGetArray(zz,&zg);CHKERRQ(ierr); z = zg;
9412d61bbb3SSatish Balay 
9422d61bbb3SSatish Balay   idx   = a->j;
9432d61bbb3SSatish Balay   v     = a->a;
9442d61bbb3SSatish Balay   ii    = a->i;
945f1af5d2fSBarry Smith   xb    = x;
9462d61bbb3SSatish Balay   switch (bs) {
9472d61bbb3SSatish Balay   case 1:
9482d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
9492d61bbb3SSatish Balay       n  = ii[1] - ii[0]; ii++;
950f1af5d2fSBarry Smith       x1 = xb[0];
9512d61bbb3SSatish Balay       ib = idx + ai[i];
9522d61bbb3SSatish Balay       for (j=0; j<n; j++) {
9532d61bbb3SSatish Balay         rval    = ib[j];
954f1af5d2fSBarry Smith         z[rval] += *v * x1;
955f1af5d2fSBarry Smith         v++;
9562d61bbb3SSatish Balay       }
957f1af5d2fSBarry Smith       xb++;
9582d61bbb3SSatish Balay     }
9592d61bbb3SSatish Balay     break;
9602d61bbb3SSatish Balay   case 2:
9612d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
9622d61bbb3SSatish Balay       n  = ii[1] - ii[0]; ii++;
963f1af5d2fSBarry Smith       x1 = xb[0]; x2 = xb[1];
9642d61bbb3SSatish Balay       ib = idx + ai[i];
9652d61bbb3SSatish Balay       for (j=0; j<n; j++) {
9662d61bbb3SSatish Balay         rval      = ib[j]*2;
9672d61bbb3SSatish Balay         z[rval++] += v[0]*x1 + v[1]*x2;
968f1af5d2fSBarry Smith         z[rval]   += v[2]*x1 + v[3]*x2;
9692d61bbb3SSatish Balay         v  += 4;
9702d61bbb3SSatish Balay       }
971f1af5d2fSBarry Smith       xb += 2;
9722d61bbb3SSatish Balay     }
9732d61bbb3SSatish Balay     break;
9742d61bbb3SSatish Balay   case 3:
9752d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
9762d61bbb3SSatish Balay       n  = ii[1] - ii[0]; ii++;
977f1af5d2fSBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
9782d61bbb3SSatish Balay       ib = idx + ai[i];
9792d61bbb3SSatish Balay       for (j=0; j<n; j++) {
9802d61bbb3SSatish Balay         rval      = ib[j]*3;
9812d61bbb3SSatish Balay         z[rval++] += v[0]*x1 + v[1]*x2 + v[2]*x3;
9822d61bbb3SSatish Balay         z[rval++] += v[3]*x1 + v[4]*x2 + v[5]*x3;
983f1af5d2fSBarry Smith         z[rval]   += v[6]*x1 + v[7]*x2 + v[8]*x3;
9842d61bbb3SSatish Balay         v  += 9;
9852d61bbb3SSatish Balay       }
986f1af5d2fSBarry Smith       xb += 3;
9872d61bbb3SSatish Balay     }
9882d61bbb3SSatish Balay     break;
9892d61bbb3SSatish Balay   case 4:
9902d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
9912d61bbb3SSatish Balay       n  = ii[1] - ii[0]; ii++;
992f1af5d2fSBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
9932d61bbb3SSatish Balay       ib = idx + ai[i];
9942d61bbb3SSatish Balay       for (j=0; j<n; j++) {
9952d61bbb3SSatish Balay         rval      = ib[j]*4;
9962d61bbb3SSatish Balay         z[rval++] +=  v[0]*x1 +  v[1]*x2 +  v[2]*x3 +  v[3]*x4;
9972d61bbb3SSatish Balay         z[rval++] +=  v[4]*x1 +  v[5]*x2 +  v[6]*x3 +  v[7]*x4;
9982d61bbb3SSatish Balay         z[rval++] +=  v[8]*x1 +  v[9]*x2 + v[10]*x3 + v[11]*x4;
999f1af5d2fSBarry Smith         z[rval]   += v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4;
10002d61bbb3SSatish Balay         v  += 16;
10012d61bbb3SSatish Balay       }
1002f1af5d2fSBarry Smith       xb += 4;
10032d61bbb3SSatish Balay     }
10042d61bbb3SSatish Balay     break;
10052d61bbb3SSatish Balay   case 5:
10062d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
10072d61bbb3SSatish Balay       n  = ii[1] - ii[0]; ii++;
1008f1af5d2fSBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
10092d61bbb3SSatish Balay       x4 = xb[3]; x5 = xb[4];
10102d61bbb3SSatish Balay       ib = idx + ai[i];
10112d61bbb3SSatish Balay       for (j=0; j<n; j++) {
10122d61bbb3SSatish Balay         rval      = ib[j]*5;
10132d61bbb3SSatish Balay         z[rval++] +=  v[0]*x1 +  v[1]*x2 +  v[2]*x3 +  v[3]*x4 +  v[4]*x5;
10142d61bbb3SSatish Balay         z[rval++] +=  v[5]*x1 +  v[6]*x2 +  v[7]*x3 +  v[8]*x4 +  v[9]*x5;
10152d61bbb3SSatish Balay         z[rval++] += v[10]*x1 + v[11]*x2 + v[12]*x3 + v[13]*x4 + v[14]*x5;
10162d61bbb3SSatish Balay         z[rval++] += v[15]*x1 + v[16]*x2 + v[17]*x3 + v[18]*x4 + v[19]*x5;
1017f1af5d2fSBarry Smith         z[rval]   += v[20]*x1 + v[21]*x2 + v[22]*x3 + v[23]*x4 + v[24]*x5;
10182d61bbb3SSatish Balay         v  += 25;
10192d61bbb3SSatish Balay       }
1020f1af5d2fSBarry Smith       xb += 5;
10212d61bbb3SSatish Balay     }
10222d61bbb3SSatish Balay     break;
1023f1af5d2fSBarry Smith   case 6:
1024f1af5d2fSBarry Smith     for (i=0; i<mbs; i++) {
1025f1af5d2fSBarry Smith       n  = ii[1] - ii[0]; ii++;
1026f1af5d2fSBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
1027f1af5d2fSBarry Smith       x4 = xb[3]; x5 = xb[4]; x6 = xb[5];
1028f1af5d2fSBarry Smith       ib = idx + ai[i];
1029f1af5d2fSBarry Smith       for (j=0; j<n; j++) {
1030f1af5d2fSBarry Smith         rval      = ib[j]*6;
1031f1af5d2fSBarry Smith         z[rval++] +=  v[0]*x1 +  v[1]*x2 +  v[2]*x3 +  v[3]*x4 + v[4]*x5 + v[5]*x6;
1032f1af5d2fSBarry Smith         z[rval++] +=  v[6]*x1 +  v[7]*x2 +  v[8]*x3 +  v[9]*x4 + v[10]*x5 + v[11]*x6;
1033f1af5d2fSBarry Smith         z[rval++] += v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4 + v[16]*x5 + v[17]*x6;
1034f1af5d2fSBarry Smith         z[rval++] += v[18]*x1 + v[19]*x2 + v[20]*x3 + v[21]*x4 + v[22]*x5 + v[23]*x6;
1035f1af5d2fSBarry Smith         z[rval++] += v[24]*x1 + v[25]*x2 + v[26]*x3 + v[27]*x4 + v[28]*x5 + v[29]*x6;
1036f1af5d2fSBarry Smith         z[rval]   += v[30]*x1 + v[31]*x2 + v[32]*x3 + v[33]*x4 + v[34]*x5 + v[35]*x6;
1037f1af5d2fSBarry Smith         v  += 36;
1038f1af5d2fSBarry Smith       }
1039f1af5d2fSBarry Smith       xb += 6;
1040f1af5d2fSBarry Smith     }
1041f1af5d2fSBarry Smith     break;
1042f1af5d2fSBarry Smith   case 7:
1043f1af5d2fSBarry Smith     for (i=0; i<mbs; i++) {
1044f1af5d2fSBarry Smith       n  = ii[1] - ii[0]; ii++;
1045f1af5d2fSBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
1046f1af5d2fSBarry Smith       x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6];
1047f1af5d2fSBarry Smith       ib = idx + ai[i];
1048f1af5d2fSBarry Smith       for (j=0; j<n; j++) {
1049f1af5d2fSBarry Smith         rval      = ib[j]*7;
1050f1af5d2fSBarry Smith         z[rval++] +=  v[0]*x1 +  v[1]*x2 +  v[2]*x3 +  v[3]*x4 + v[4]*x5 + v[5]*x6 + v[6]*x7;
1051f1af5d2fSBarry Smith         z[rval++] +=  v[7]*x1 +  v[8]*x2 +  v[9]*x3 + v[10]*x4 + v[11]*x5 + v[12]*x6 + v[13]*x7;
1052f1af5d2fSBarry Smith         z[rval++] += v[14]*x1 + v[15]*x2 + v[16]*x3 + v[17]*x4 + v[18]*x5 + v[19]*x6 + v[20]*x7;
1053f1af5d2fSBarry Smith         z[rval++] += v[21]*x1 + v[22]*x2 + v[23]*x3 + v[24]*x4 + v[25]*x5 + v[26]*x6 + v[27]*x7;
1054f1af5d2fSBarry Smith         z[rval++] += v[28]*x1 + v[29]*x2 + v[30]*x3 + v[31]*x4 + v[32]*x5 + v[33]*x6 + v[34]*x7;
1055f1af5d2fSBarry Smith         z[rval++] += v[35]*x1 + v[36]*x2 + v[37]*x3 + v[38]*x4 + v[39]*x5 + v[40]*x6 + v[41]*x7;
1056f1af5d2fSBarry Smith         z[rval]   += v[42]*x1 + v[43]*x2 + v[44]*x3 + v[45]*x4 + v[46]*x5 + v[47]*x6 + v[48]*x7;
1057f1af5d2fSBarry Smith         v  += 49;
1058f1af5d2fSBarry Smith       }
1059f1af5d2fSBarry Smith       xb += 7;
1060f1af5d2fSBarry Smith     }
1061f1af5d2fSBarry Smith     break;
1062f1af5d2fSBarry Smith   default: {       /* block sizes larger then 7 by 7 are handled by BLAS */
10633f1db9ecSBarry Smith       int          ncols,k;
106487828ca2SBarry Smith       PetscScalar  *work,*workt;
10653f1db9ecSBarry Smith 
10662d61bbb3SSatish Balay       if (!a->mult_work) {
1067273d9f13SBarry Smith         k = PetscMax(A->m,A->n);
106887828ca2SBarry Smith         ierr = PetscMalloc((k+1)*sizeof(PetscScalar),&a->mult_work);CHKERRQ(ierr);
10692d61bbb3SSatish Balay       }
10702d61bbb3SSatish Balay       work = a->mult_work;
10712d61bbb3SSatish Balay       for (i=0; i<mbs; i++) {
10722d61bbb3SSatish Balay         n     = ii[1] - ii[0]; ii++;
10732d61bbb3SSatish Balay         ncols = n*bs;
107487828ca2SBarry Smith         ierr  = PetscMemzero(work,ncols*sizeof(PetscScalar));CHKERRQ(ierr);
1075d824769bSBarry Smith         Kernel_w_gets_w_plus_trans_Ar_times_v(bs,ncols,x,v,work);
10763f1db9ecSBarry Smith         /* LAgemv_("T",&bs,&ncols,&_DOne,v,&bs,x,&_One,&_DOne,work,&_One); */
10772d61bbb3SSatish Balay         v += n*bs2;
10782d61bbb3SSatish Balay         x += bs;
10792d61bbb3SSatish Balay         workt = work;
10802d61bbb3SSatish Balay         for (j=0; j<n; j++) {
10812d61bbb3SSatish Balay           zb = z + bs*(*idx++);
10822d61bbb3SSatish Balay           for (k=0; k<bs; k++) zb[k] += workt[k] ;
10832d61bbb3SSatish Balay           workt += bs;
10842d61bbb3SSatish Balay         }
10852d61bbb3SSatish Balay       }
10862d61bbb3SSatish Balay     }
10872d61bbb3SSatish Balay   }
10881ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&xg);CHKERRQ(ierr);
10891ebc52fbSHong Zhang   ierr = VecRestoreArray(zz,&zg);CHKERRQ(ierr);
1090b0a32e0cSBarry Smith   PetscLogFlops(2*a->nz*a->bs2 - A->n);
10912d61bbb3SSatish Balay   PetscFunctionReturn(0);
10922d61bbb3SSatish Balay }
10932d61bbb3SSatish Balay 
10944a2ae208SSatish Balay #undef __FUNCT__
10954a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_SeqBAIJ"
1096dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_SeqBAIJ(Mat A,Vec xx,Vec yy,Vec zz)
10972d61bbb3SSatish Balay 
10982d61bbb3SSatish Balay {
10992d61bbb3SSatish Balay   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ*)A->data;
110087828ca2SBarry Smith   PetscScalar     *xg,*zg,*zb,*x,*z,*xb,x1,x2,x3,x4,x5;
11013f1db9ecSBarry Smith   MatScalar       *v;
1102*6849ba73SBarry Smith   PetscErrorCode ierr;
1103*6849ba73SBarry Smith   int             mbs=a->mbs,i,*idx,*ii,*ai=a->i,rval,bs=a->bs,j,n,bs2=a->bs2,*ib;
11042d61bbb3SSatish Balay 
11052d61bbb3SSatish Balay   PetscFunctionBegin;
1106ef66eb69SBarry Smith   if (yy != zz) { ierr = VecCopy(yy,zz);CHKERRQ(ierr); }
11071ebc52fbSHong Zhang   ierr = VecGetArray(xx,&xg);CHKERRQ(ierr); x = xg;
11081ebc52fbSHong Zhang   ierr = VecGetArray(zz,&zg);CHKERRQ(ierr); z = zg;
11092d61bbb3SSatish Balay 
11102d61bbb3SSatish Balay 
11112d61bbb3SSatish Balay   idx   = a->j;
11122d61bbb3SSatish Balay   v     = a->a;
11132d61bbb3SSatish Balay   ii    = a->i;
1114f1af5d2fSBarry Smith   xb    = x;
11152d61bbb3SSatish Balay 
11162d61bbb3SSatish Balay   switch (bs) {
11172d61bbb3SSatish Balay   case 1:
11182d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
11192d61bbb3SSatish Balay       n  = ii[1] - ii[0]; ii++;
1120f1af5d2fSBarry Smith       x1 = xb[0];
11212d61bbb3SSatish Balay       ib = idx + ai[i];
11222d61bbb3SSatish Balay       for (j=0; j<n; j++) {
11232d61bbb3SSatish Balay         rval    = ib[j];
1124f1af5d2fSBarry Smith         z[rval] += *v * x1;
1125f1af5d2fSBarry Smith         v++;
11262d61bbb3SSatish Balay       }
1127f1af5d2fSBarry Smith       xb++;
11282d61bbb3SSatish Balay     }
11292d61bbb3SSatish Balay     break;
11302d61bbb3SSatish Balay   case 2:
11312d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
11322d61bbb3SSatish Balay       n  = ii[1] - ii[0]; ii++;
1133f1af5d2fSBarry Smith       x1 = xb[0]; x2 = xb[1];
11342d61bbb3SSatish Balay       ib = idx + ai[i];
11352d61bbb3SSatish Balay       for (j=0; j<n; j++) {
11362d61bbb3SSatish Balay         rval      = ib[j]*2;
11372d61bbb3SSatish Balay         z[rval++] += v[0]*x1 + v[1]*x2;
11382d61bbb3SSatish Balay         z[rval++] += v[2]*x1 + v[3]*x2;
11392d61bbb3SSatish Balay         v  += 4;
11402d61bbb3SSatish Balay       }
1141f1af5d2fSBarry Smith       xb += 2;
11422d61bbb3SSatish Balay     }
11432d61bbb3SSatish Balay     break;
11442d61bbb3SSatish Balay   case 3:
11452d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
11462d61bbb3SSatish Balay       n  = ii[1] - ii[0]; ii++;
1147f1af5d2fSBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
11482d61bbb3SSatish Balay       ib = idx + ai[i];
11492d61bbb3SSatish Balay       for (j=0; j<n; j++) {
11502d61bbb3SSatish Balay         rval      = ib[j]*3;
11512d61bbb3SSatish Balay         z[rval++] += v[0]*x1 + v[1]*x2 + v[2]*x3;
11522d61bbb3SSatish Balay         z[rval++] += v[3]*x1 + v[4]*x2 + v[5]*x3;
11532d61bbb3SSatish Balay         z[rval++] += v[6]*x1 + v[7]*x2 + v[8]*x3;
11542d61bbb3SSatish Balay         v  += 9;
11552d61bbb3SSatish Balay       }
1156f1af5d2fSBarry Smith       xb += 3;
11572d61bbb3SSatish Balay     }
11582d61bbb3SSatish Balay     break;
11592d61bbb3SSatish Balay   case 4:
11602d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
11612d61bbb3SSatish Balay       n  = ii[1] - ii[0]; ii++;
1162f1af5d2fSBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
11632d61bbb3SSatish Balay       ib = idx + ai[i];
11642d61bbb3SSatish Balay       for (j=0; j<n; j++) {
11652d61bbb3SSatish Balay         rval      = ib[j]*4;
11662d61bbb3SSatish Balay         z[rval++] +=  v[0]*x1 +  v[1]*x2 +  v[2]*x3 +  v[3]*x4;
11672d61bbb3SSatish Balay         z[rval++] +=  v[4]*x1 +  v[5]*x2 +  v[6]*x3 +  v[7]*x4;
11682d61bbb3SSatish Balay         z[rval++] +=  v[8]*x1 +  v[9]*x2 + v[10]*x3 + v[11]*x4;
11692d61bbb3SSatish Balay         z[rval++] += v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4;
11702d61bbb3SSatish Balay         v  += 16;
11712d61bbb3SSatish Balay       }
1172f1af5d2fSBarry Smith       xb += 4;
11732d61bbb3SSatish Balay     }
11742d61bbb3SSatish Balay     break;
11752d61bbb3SSatish Balay   case 5:
11762d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
11772d61bbb3SSatish Balay       n  = ii[1] - ii[0]; ii++;
1178f1af5d2fSBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
11792d61bbb3SSatish Balay       x4 = xb[3]; x5 = xb[4];
11802d61bbb3SSatish Balay       ib = idx + ai[i];
11812d61bbb3SSatish Balay       for (j=0; j<n; j++) {
11822d61bbb3SSatish Balay         rval      = ib[j]*5;
11832d61bbb3SSatish Balay         z[rval++] +=  v[0]*x1 +  v[1]*x2 +  v[2]*x3 +  v[3]*x4 +  v[4]*x5;
11842d61bbb3SSatish Balay         z[rval++] +=  v[5]*x1 +  v[6]*x2 +  v[7]*x3 +  v[8]*x4 +  v[9]*x5;
11852d61bbb3SSatish Balay         z[rval++] += v[10]*x1 + v[11]*x2 + v[12]*x3 + v[13]*x4 + v[14]*x5;
11862d61bbb3SSatish Balay         z[rval++] += v[15]*x1 + v[16]*x2 + v[17]*x3 + v[18]*x4 + v[19]*x5;
11872d61bbb3SSatish Balay         z[rval++] += v[20]*x1 + v[21]*x2 + v[22]*x3 + v[23]*x4 + v[24]*x5;
11882d61bbb3SSatish Balay         v  += 25;
11892d61bbb3SSatish Balay       }
1190f1af5d2fSBarry Smith       xb += 5;
11912d61bbb3SSatish Balay     }
11922d61bbb3SSatish Balay     break;
1193f1af5d2fSBarry Smith   default: {      /* block sizes larger then 5 by 5 are handled by BLAS */
11943f1db9ecSBarry Smith       int          ncols,k;
119587828ca2SBarry Smith       PetscScalar  *work,*workt;
11963f1db9ecSBarry Smith 
11972d61bbb3SSatish Balay       if (!a->mult_work) {
1198273d9f13SBarry Smith         k = PetscMax(A->m,A->n);
119987828ca2SBarry Smith         ierr = PetscMalloc((k+1)*sizeof(PetscScalar),&a->mult_work);CHKERRQ(ierr);
12002d61bbb3SSatish Balay       }
12012d61bbb3SSatish Balay       work = a->mult_work;
12022d61bbb3SSatish Balay       for (i=0; i<mbs; i++) {
12032d61bbb3SSatish Balay         n     = ii[1] - ii[0]; ii++;
12042d61bbb3SSatish Balay         ncols = n*bs;
120587828ca2SBarry Smith         ierr  = PetscMemzero(work,ncols*sizeof(PetscScalar));CHKERRQ(ierr);
12063f1db9ecSBarry Smith         Kernel_w_gets_w_plus_trans_Ar_times_v(bs,ncols,x,v,work);
12073f1db9ecSBarry Smith         /* LAgemv_("T",&bs,&ncols,&_DOne,v,&bs,x,&_One,&_DOne,work,&_One); */
12082d61bbb3SSatish Balay         v += n*bs2;
12092d61bbb3SSatish Balay         x += bs;
12102d61bbb3SSatish Balay         workt = work;
12112d61bbb3SSatish Balay         for (j=0; j<n; j++) {
12122d61bbb3SSatish Balay           zb = z + bs*(*idx++);
12132d61bbb3SSatish Balay           for (k=0; k<bs; k++) zb[k] += workt[k] ;
12142d61bbb3SSatish Balay           workt += bs;
12152d61bbb3SSatish Balay         }
12162d61bbb3SSatish Balay       }
12172d61bbb3SSatish Balay     }
12182d61bbb3SSatish Balay   }
12191ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&xg);CHKERRQ(ierr);
12201ebc52fbSHong Zhang   ierr = VecRestoreArray(zz,&zg);CHKERRQ(ierr);
1221b0a32e0cSBarry Smith   PetscLogFlops(2*a->nz*a->bs2);
12222d61bbb3SSatish Balay   PetscFunctionReturn(0);
12232d61bbb3SSatish Balay }
12242d61bbb3SSatish Balay 
12254a2ae208SSatish Balay #undef __FUNCT__
12264a2ae208SSatish Balay #define __FUNCT__ "MatScale_SeqBAIJ"
1227dfbe8321SBarry Smith PetscErrorCode MatScale_SeqBAIJ(const PetscScalar *alpha,Mat inA)
12282d61bbb3SSatish Balay {
12292d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)inA->data;
123087052302SDinesh Kaushik   int         totalnz = a->bs2*a->nz;
12313eda8832SBarry Smith #if defined(PETSC_USE_MAT_SINGLE)
12323eda8832SBarry Smith   int         i;
123387052302SDinesh Kaushik #else
123487052302SDinesh Kaushik   int         one = 1;
12353eda8832SBarry Smith #endif
12362d61bbb3SSatish Balay 
12372d61bbb3SSatish Balay   PetscFunctionBegin;
12383eda8832SBarry Smith #if defined(PETSC_USE_MAT_SINGLE)
12393eda8832SBarry Smith   for (i=0; i<totalnz; i++) a->a[i] *= *alpha;
12403eda8832SBarry Smith #else
1241268466fbSBarry Smith   BLscal_(&totalnz,(PetscScalar*)alpha,a->a,&one);
12423eda8832SBarry Smith #endif
1243b0a32e0cSBarry Smith   PetscLogFlops(totalnz);
12442d61bbb3SSatish Balay   PetscFunctionReturn(0);
12452d61bbb3SSatish Balay }
12462d61bbb3SSatish Balay 
12474a2ae208SSatish Balay #undef __FUNCT__
12484a2ae208SSatish Balay #define __FUNCT__ "MatNorm_SeqBAIJ"
1249dfbe8321SBarry Smith PetscErrorCode MatNorm_SeqBAIJ(Mat A,NormType type,PetscReal *norm)
12502d61bbb3SSatish Balay {
12512d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
12523f1db9ecSBarry Smith   MatScalar   *v = a->a;
1253329f5518SBarry Smith   PetscReal   sum = 0.0;
12540e90e235SBarry Smith   int         i,j,k,bs = a->bs,nz=a->nz,bs2=a->bs2,k1;
12552d61bbb3SSatish Balay 
12562d61bbb3SSatish Balay   PetscFunctionBegin;
12572d61bbb3SSatish Balay   if (type == NORM_FROBENIUS) {
12582d61bbb3SSatish Balay     for (i=0; i< bs2*nz; i++) {
1259aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
1260329f5518SBarry Smith       sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
12612d61bbb3SSatish Balay #else
12622d61bbb3SSatish Balay       sum += (*v)*(*v); v++;
12632d61bbb3SSatish Balay #endif
12642d61bbb3SSatish Balay     }
12652d61bbb3SSatish Balay     *norm = sqrt(sum);
1266596552b5SBarry Smith   }  else if (type == NORM_INFINITY) { /* maximum row sum */
1267596552b5SBarry Smith     *norm = 0.0;
1268596552b5SBarry Smith     for (k=0; k<bs; k++) {
126974f84c7bSSatish Balay       for (j=0; j<a->mbs; j++) {
1270596552b5SBarry Smith         v = a->a + bs2*a->i[j] + k;
1271596552b5SBarry Smith         sum = 0.0;
1272596552b5SBarry Smith         for (i=0; i<a->i[j+1]-a->i[j]; i++) {
12730e90e235SBarry Smith           for (k1=0; k1<bs; k1++){
1274596552b5SBarry Smith             sum += PetscAbsScalar(*v);
1275596552b5SBarry Smith             v   += bs;
12762d61bbb3SSatish Balay           }
12770e90e235SBarry Smith         }
1278596552b5SBarry Smith         if (sum > *norm) *norm = sum;
1279596552b5SBarry Smith       }
1280596552b5SBarry Smith     }
1281596552b5SBarry Smith   } else {
128229bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"No support for this norm yet");
12832d61bbb3SSatish Balay   }
12842d61bbb3SSatish Balay   PetscFunctionReturn(0);
12852d61bbb3SSatish Balay }
12862d61bbb3SSatish Balay 
12872d61bbb3SSatish Balay 
12884a2ae208SSatish Balay #undef __FUNCT__
12894a2ae208SSatish Balay #define __FUNCT__ "MatEqual_SeqBAIJ"
1290dfbe8321SBarry Smith PetscErrorCode MatEqual_SeqBAIJ(Mat A,Mat B,PetscTruth* flg)
12912d61bbb3SSatish Balay {
12922d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data,*b = (Mat_SeqBAIJ *)B->data;
1293dfbe8321SBarry Smith   PetscErrorCode ierr;
12942d61bbb3SSatish Balay 
12952d61bbb3SSatish Balay   PetscFunctionBegin;
12962d61bbb3SSatish Balay   /* If the  matrix/block dimensions are not equal, or no of nonzeros or shift */
1297273d9f13SBarry Smith   if ((A->m != B->m) || (A->n != B->n) || (a->bs != b->bs)|| (a->nz != b->nz)) {
1298273d9f13SBarry Smith     *flg = PETSC_FALSE;
1299273d9f13SBarry Smith     PetscFunctionReturn(0);
13002d61bbb3SSatish Balay   }
13012d61bbb3SSatish Balay 
13022d61bbb3SSatish Balay   /* if the a->i are the same */
13030f5bd95cSBarry Smith   ierr = PetscMemcmp(a->i,b->i,(a->mbs+1)*sizeof(int),flg);CHKERRQ(ierr);
13040f5bd95cSBarry Smith   if (*flg == PETSC_FALSE) {
13050f5bd95cSBarry Smith     PetscFunctionReturn(0);
13062d61bbb3SSatish Balay   }
13072d61bbb3SSatish Balay 
13082d61bbb3SSatish Balay   /* if a->j are the same */
13090f5bd95cSBarry Smith   ierr = PetscMemcmp(a->j,b->j,(a->nz)*sizeof(int),flg);CHKERRQ(ierr);
13100f5bd95cSBarry Smith   if (*flg == PETSC_FALSE) {
13110f5bd95cSBarry Smith     PetscFunctionReturn(0);
13122d61bbb3SSatish Balay   }
13132d61bbb3SSatish Balay   /* if a->a are the same */
131487828ca2SBarry Smith   ierr = PetscMemcmp(a->a,b->a,(a->nz)*(a->bs)*(a->bs)*sizeof(PetscScalar),flg);CHKERRQ(ierr);
13152d61bbb3SSatish Balay   PetscFunctionReturn(0);
13162d61bbb3SSatish Balay 
13172d61bbb3SSatish Balay }
13182d61bbb3SSatish Balay 
13194a2ae208SSatish Balay #undef __FUNCT__
13204a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_SeqBAIJ"
1321dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_SeqBAIJ(Mat A,Vec v)
13222d61bbb3SSatish Balay {
13232d61bbb3SSatish Balay   Mat_SeqBAIJ  *a = (Mat_SeqBAIJ*)A->data;
1324dfbe8321SBarry Smith   PetscErrorCode ierr;
1325dfbe8321SBarry Smith   int  i,j,k,n,row,bs,*ai,*aj,ambs,bs2;
132687828ca2SBarry Smith   PetscScalar  *x,zero = 0.0;
13273f1db9ecSBarry Smith   MatScalar    *aa,*aa_j;
13282d61bbb3SSatish Balay 
13292d61bbb3SSatish Balay   PetscFunctionBegin;
133029bbc08cSBarry Smith   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
13312d61bbb3SSatish Balay   bs   = a->bs;
13322d61bbb3SSatish Balay   aa   = a->a;
13332d61bbb3SSatish Balay   ai   = a->i;
13342d61bbb3SSatish Balay   aj   = a->j;
13352d61bbb3SSatish Balay   ambs = a->mbs;
13362d61bbb3SSatish Balay   bs2  = a->bs2;
13372d61bbb3SSatish Balay 
1338e1311b90SBarry Smith   ierr = VecSet(&zero,v);CHKERRQ(ierr);
13391ebc52fbSHong Zhang   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1340e1311b90SBarry Smith   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
1341273d9f13SBarry Smith   if (n != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
13422d61bbb3SSatish Balay   for (i=0; i<ambs; i++) {
13432d61bbb3SSatish Balay     for (j=ai[i]; j<ai[i+1]; j++) {
13442d61bbb3SSatish Balay       if (aj[j] == i) {
13452d61bbb3SSatish Balay         row  = i*bs;
13462d61bbb3SSatish Balay         aa_j = aa+j*bs2;
13472d61bbb3SSatish Balay         for (k=0; k<bs2; k+=(bs+1),row++) x[row] = aa_j[k];
13482d61bbb3SSatish Balay         break;
13492d61bbb3SSatish Balay       }
13502d61bbb3SSatish Balay     }
13512d61bbb3SSatish Balay   }
13521ebc52fbSHong Zhang   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
13532d61bbb3SSatish Balay   PetscFunctionReturn(0);
13542d61bbb3SSatish Balay }
13552d61bbb3SSatish Balay 
13564a2ae208SSatish Balay #undef __FUNCT__
13574a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_SeqBAIJ"
1358dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_SeqBAIJ(Mat A,Vec ll,Vec rr)
13592d61bbb3SSatish Balay {
13602d61bbb3SSatish Balay   Mat_SeqBAIJ  *a = (Mat_SeqBAIJ*)A->data;
136187828ca2SBarry Smith   PetscScalar  *l,*r,x,*li,*ri;
13623f1db9ecSBarry Smith   MatScalar    *aa,*v;
1363dfbe8321SBarry Smith   PetscErrorCode ierr;
1364dfbe8321SBarry Smith   int i,j,k,lm,rn,M,m,n,*ai,*aj,mbs,tmp,bs,bs2;
13652d61bbb3SSatish Balay 
13662d61bbb3SSatish Balay   PetscFunctionBegin;
13672d61bbb3SSatish Balay   ai  = a->i;
13682d61bbb3SSatish Balay   aj  = a->j;
13692d61bbb3SSatish Balay   aa  = a->a;
1370273d9f13SBarry Smith   m   = A->m;
1371273d9f13SBarry Smith   n   = A->n;
13722d61bbb3SSatish Balay   bs  = a->bs;
13732d61bbb3SSatish Balay   mbs = a->mbs;
13742d61bbb3SSatish Balay   bs2 = a->bs2;
13752d61bbb3SSatish Balay   if (ll) {
13761ebc52fbSHong Zhang     ierr = VecGetArray(ll,&l);CHKERRQ(ierr);
1377e1311b90SBarry Smith     ierr = VecGetLocalSize(ll,&lm);CHKERRQ(ierr);
137829bbc08cSBarry Smith     if (lm != m) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vector wrong length");
13792d61bbb3SSatish Balay     for (i=0; i<mbs; i++) { /* for each block row */
13802d61bbb3SSatish Balay       M  = ai[i+1] - ai[i];
13812d61bbb3SSatish Balay       li = l + i*bs;
13822d61bbb3SSatish Balay       v  = aa + bs2*ai[i];
13832d61bbb3SSatish Balay       for (j=0; j<M; j++) { /* for each block */
13842d61bbb3SSatish Balay         for (k=0; k<bs2; k++) {
13852d61bbb3SSatish Balay           (*v++) *= li[k%bs];
13862d61bbb3SSatish Balay         }
13872d61bbb3SSatish Balay       }
13882d61bbb3SSatish Balay     }
13891ebc52fbSHong Zhang     ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr);
1390b0a32e0cSBarry Smith     PetscLogFlops(a->nz);
13912d61bbb3SSatish Balay   }
13922d61bbb3SSatish Balay 
13932d61bbb3SSatish Balay   if (rr) {
13941ebc52fbSHong Zhang     ierr = VecGetArray(rr,&r);CHKERRQ(ierr);
1395e1311b90SBarry Smith     ierr = VecGetLocalSize(rr,&rn);CHKERRQ(ierr);
139629bbc08cSBarry Smith     if (rn != n) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vector wrong length");
13972d61bbb3SSatish Balay     for (i=0; i<mbs; i++) { /* for each block row */
13982d61bbb3SSatish Balay       M  = ai[i+1] - ai[i];
13992d61bbb3SSatish Balay       v  = aa + bs2*ai[i];
14002d61bbb3SSatish Balay       for (j=0; j<M; j++) { /* for each block */
14012d61bbb3SSatish Balay         ri = r + bs*aj[ai[i]+j];
14022d61bbb3SSatish Balay         for (k=0; k<bs; k++) {
14032d61bbb3SSatish Balay           x = ri[k];
14042d61bbb3SSatish Balay           for (tmp=0; tmp<bs; tmp++) (*v++) *= x;
14052d61bbb3SSatish Balay         }
14062d61bbb3SSatish Balay       }
14072d61bbb3SSatish Balay     }
14081ebc52fbSHong Zhang     ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr);
1409b0a32e0cSBarry Smith     PetscLogFlops(a->nz);
14102d61bbb3SSatish Balay   }
14112d61bbb3SSatish Balay   PetscFunctionReturn(0);
14122d61bbb3SSatish Balay }
14132d61bbb3SSatish Balay 
14142d61bbb3SSatish Balay 
14154a2ae208SSatish Balay #undef __FUNCT__
14164a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_SeqBAIJ"
1417dfbe8321SBarry Smith PetscErrorCode MatGetInfo_SeqBAIJ(Mat A,MatInfoType flag,MatInfo *info)
14182d61bbb3SSatish Balay {
14192d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
14202d61bbb3SSatish Balay 
14212d61bbb3SSatish Balay   PetscFunctionBegin;
1422273d9f13SBarry Smith   info->rows_global    = (double)A->m;
1423273d9f13SBarry Smith   info->columns_global = (double)A->n;
1424273d9f13SBarry Smith   info->rows_local     = (double)A->m;
1425273d9f13SBarry Smith   info->columns_local  = (double)A->n;
14262d61bbb3SSatish Balay   info->block_size     = a->bs2;
14272d61bbb3SSatish Balay   info->nz_allocated   = a->maxnz;
14282d61bbb3SSatish Balay   info->nz_used        = a->bs2*a->nz;
14292d61bbb3SSatish Balay   info->nz_unneeded    = (double)(info->nz_allocated - info->nz_used);
14302d61bbb3SSatish Balay   info->assemblies   = A->num_ass;
14312d61bbb3SSatish Balay   info->mallocs      = a->reallocs;
14322d61bbb3SSatish Balay   info->memory       = A->mem;
14332d61bbb3SSatish Balay   if (A->factor) {
14342d61bbb3SSatish Balay     info->fill_ratio_given  = A->info.fill_ratio_given;
14352d61bbb3SSatish Balay     info->fill_ratio_needed = A->info.fill_ratio_needed;
14362d61bbb3SSatish Balay     info->factor_mallocs    = A->info.factor_mallocs;
14372d61bbb3SSatish Balay   } else {
14382d61bbb3SSatish Balay     info->fill_ratio_given  = 0;
14392d61bbb3SSatish Balay     info->fill_ratio_needed = 0;
14402d61bbb3SSatish Balay     info->factor_mallocs    = 0;
14412d61bbb3SSatish Balay   }
14422d61bbb3SSatish Balay   PetscFunctionReturn(0);
14432d61bbb3SSatish Balay }
14442d61bbb3SSatish Balay 
14452d61bbb3SSatish Balay 
14464a2ae208SSatish Balay #undef __FUNCT__
14474a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_SeqBAIJ"
1448dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_SeqBAIJ(Mat A)
14492d61bbb3SSatish Balay {
14502d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
1451dfbe8321SBarry Smith   PetscErrorCode ierr;
14522d61bbb3SSatish Balay 
14532d61bbb3SSatish Balay   PetscFunctionBegin;
1454549d3d68SSatish Balay   ierr = PetscMemzero(a->a,a->bs2*a->i[a->mbs]*sizeof(MatScalar));CHKERRQ(ierr);
14552d61bbb3SSatish Balay   PetscFunctionReturn(0);
14562d61bbb3SSatish Balay }
1457