xref: /petsc/src/mat/impls/baij/seq/baij2.c (revision dcf5cc722893849f9d6a49faec7bb9e7955fb956)
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"
9690b6cddSBarry Smith PetscErrorCode MatIncreaseOverlap_SeqBAIJ(Mat A,PetscInt is_max,IS is[],PetscInt ov)
10a3192f15SSatish Balay {
11a3192f15SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
126849ba73SBarry Smith   PetscErrorCode ierr;
13690b6cddSBarry Smith   PetscInt       row,i,j,k,l,m,n,*idx,*nidx,isz,val,ival;
14690b6cddSBarry Smith   PetscInt       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;
21521d7252SBarry Smith   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);
26690b6cddSBarry Smith   ierr = PetscMalloc((m+1)*sizeof(PetscInt),&nidx);CHKERRQ(ierr);
27690b6cddSBarry Smith   ierr = PetscMalloc((A->m+1)*sizeof(PetscInt),&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"
75690b6cddSBarry Smith PetscErrorCode MatGetSubMatrix_SeqBAIJ_Private(Mat A,IS isrow,IS iscol,PetscInt cs,MatReuse scall,Mat *B)
76736121d4SSatish Balay {
77736121d4SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*c;
786849ba73SBarry Smith   PetscErrorCode ierr;
79690b6cddSBarry Smith   PetscInt       *smap,i,k,kstart,kend,oldcols = a->nbs,*lens;
80690b6cddSBarry Smith   PetscInt       row,mat_i,*mat_j,tcol,*mat_ilen;
81521d7252SBarry Smith   PetscInt       *irow,*icol,nrows,ncols,*ssmap,bs=A->bs,bs2=a->bs2;
82690b6cddSBarry Smith   PetscInt       *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 
96690b6cddSBarry Smith   ierr = PetscMalloc((1+oldcols)*sizeof(PetscInt),&smap);CHKERRQ(ierr);
97736121d4SSatish Balay   ssmap = smap;
98690b6cddSBarry Smith   ierr = PetscMalloc((1+nrows)*sizeof(PetscInt),&lens);CHKERRQ(ierr);
99690b6cddSBarry Smith   ierr  = PetscMemzero(smap,oldcols*sizeof(PetscInt));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 
116521d7252SBarry Smith     if (c->mbs!=nrows || c->nbs!=ncols || (*B)->bs!=bs) SETERRQ(PETSC_ERR_ARG_SIZ,"Submatrix wrong size");
117690b6cddSBarry Smith     ierr = PetscMemcmp(c->ilen,lens,c->mbs *sizeof(PetscInt),&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     }
121690b6cddSBarry Smith     ierr = PetscMemzero(c->ilen,c->mbs*sizeof(PetscInt));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"
161690b6cddSBarry Smith PetscErrorCode MatGetSubMatrix_SeqBAIJ(Mat A,IS isrow,IS iscol,PetscInt cs,MatReuse scall,Mat *B)
162218c64b6SSatish Balay {
163218c64b6SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
164218c64b6SSatish Balay   IS             is1,is2;
1656849ba73SBarry Smith   PetscErrorCode ierr;
166521d7252SBarry Smith   PetscInt       *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 */
176690b6cddSBarry Smith   ierr = PetscMalloc(2*(a->mbs+1)*sizeof(PetscInt),&vary);CHKERRQ(ierr);
177218c64b6SSatish Balay   iary = vary + a->mbs;
178690b6cddSBarry Smith   ierr = PetscMemzero(vary,(a->mbs)*sizeof(PetscInt));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 
187690b6cddSBarry Smith   ierr = PetscMemzero(vary,(a->mbs)*sizeof(PetscInt));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"
207690b6cddSBarry Smith PetscErrorCode MatGetSubMatrices_SeqBAIJ(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
208736121d4SSatish Balay {
2096849ba73SBarry Smith   PetscErrorCode ierr;
210690b6cddSBarry Smith   PetscInt       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;
2366849ba73SBarry Smith   PetscErrorCode ierr;
2377b2bb3b9SHong Zhang   PetscInt       mbs,i,*idx,*ii,n,*ridx=PETSC_NULL;
23826e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
2392d61bbb3SSatish Balay 
2402d61bbb3SSatish Balay   PetscFunctionBegin;
2411ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
2421ebc52fbSHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
2432d61bbb3SSatish Balay 
2442d61bbb3SSatish Balay   idx = a->j;
2452d61bbb3SSatish Balay   v   = a->a;
24626e093fcSHong Zhang   if (usecprow){
24726e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
24826e093fcSHong Zhang     ii   = a->compressedrow.i;
2497b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
25026e093fcSHong Zhang   } else {
25126e093fcSHong Zhang     mbs = a->mbs;
2522d61bbb3SSatish Balay     ii  = a->i;
25326e093fcSHong Zhang   }
2542d61bbb3SSatish Balay 
2552d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
2562d61bbb3SSatish Balay     n    = ii[1] - ii[0]; ii++;
2572d61bbb3SSatish Balay     sum  = 0.0;
2582d61bbb3SSatish Balay     while (n--) sum += *v++ * x[*idx++];
25926e093fcSHong Zhang     if (usecprow){
2607b2bb3b9SHong Zhang       z[ridx[i]] = sum;
26126e093fcSHong Zhang     } else {
2622d61bbb3SSatish Balay       z[i] = sum;
2632d61bbb3SSatish Balay     }
26426e093fcSHong Zhang   }
2651ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
2661ebc52fbSHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
26726e093fcSHong Zhang   PetscLogFlops(2*a->nz - mbs);
2682d61bbb3SSatish Balay   PetscFunctionReturn(0);
2692d61bbb3SSatish Balay }
2702d61bbb3SSatish Balay 
2714a2ae208SSatish Balay #undef __FUNCT__
2724a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_2"
273dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_2(Mat A,Vec xx,Vec zz)
2742d61bbb3SSatish Balay {
2752d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
276*dcf5cc72SBarry Smith   PetscScalar    *x,*z = 0,*xb,sum1,sum2,*zarray;
27787828ca2SBarry Smith   PetscScalar    x1,x2;
2783f1db9ecSBarry Smith   MatScalar      *v;
279dfbe8321SBarry Smith   PetscErrorCode ierr;
2807b2bb3b9SHong Zhang   PetscInt       mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL;
28126e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
2822d61bbb3SSatish Balay 
2832d61bbb3SSatish Balay   PetscFunctionBegin;
2841ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
28526e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
2862d61bbb3SSatish Balay 
2872d61bbb3SSatish Balay   idx = a->j;
2882d61bbb3SSatish Balay   v   = a->a;
28926e093fcSHong Zhang   if (usecprow){
29026e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
29126e093fcSHong Zhang     ii   = a->compressedrow.i;
2927b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
29326e093fcSHong Zhang   } else {
29426e093fcSHong Zhang     mbs = a->mbs;
2952d61bbb3SSatish Balay     ii  = a->i;
29626e093fcSHong Zhang     z   = zarray;
29726e093fcSHong Zhang   }
2982d61bbb3SSatish Balay 
2992d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
3002d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
3012d61bbb3SSatish Balay     sum1 = 0.0; sum2 = 0.0;
3022d61bbb3SSatish Balay     for (j=0; j<n; j++) {
3032d61bbb3SSatish Balay       xb = x + 2*(*idx++); x1 = xb[0]; x2 = xb[1];
3042d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[2]*x2;
3052d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[3]*x2;
3062d61bbb3SSatish Balay       v += 4;
3072d61bbb3SSatish Balay     }
3087b2bb3b9SHong Zhang     if (usecprow) z = zarray + 2*ridx[i];
3092d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2;
31026e093fcSHong Zhang     if (!usecprow) z += 2;
3112d61bbb3SSatish Balay   }
3121ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
31326e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
31426e093fcSHong Zhang   PetscLogFlops(8*a->nz - 2*mbs);
3152d61bbb3SSatish Balay   PetscFunctionReturn(0);
3162d61bbb3SSatish Balay }
3172d61bbb3SSatish Balay 
3184a2ae208SSatish Balay #undef __FUNCT__
3194a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_3"
320dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_3(Mat A,Vec xx,Vec zz)
3212d61bbb3SSatish Balay {
3222d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
323*dcf5cc72SBarry Smith   PetscScalar    *x,*z = 0,*xb,sum1,sum2,sum3,x1,x2,x3,*zarray;
3243f1db9ecSBarry Smith   MatScalar      *v;
325dfbe8321SBarry Smith   PetscErrorCode ierr;
3267b2bb3b9SHong Zhang   PetscInt       mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL;
32726e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
32826e093fcSHong Zhang 
3292d61bbb3SSatish Balay 
330b6410449SSatish Balay #if defined(PETSC_HAVE_PRAGMA_DISJOINT)
331fee21e36SBarry Smith #pragma disjoint(*v,*z,*xb)
332fee21e36SBarry Smith #endif
333fee21e36SBarry Smith 
3342d61bbb3SSatish Balay   PetscFunctionBegin;
3351ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
33626e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
3372d61bbb3SSatish Balay 
3382d61bbb3SSatish Balay   idx = a->j;
3392d61bbb3SSatish Balay   v   = a->a;
34026e093fcSHong Zhang   if (usecprow){
34126e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
34226e093fcSHong Zhang     ii   = a->compressedrow.i;
3437b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
34426e093fcSHong Zhang   } else {
34526e093fcSHong Zhang     mbs = a->mbs;
3462d61bbb3SSatish Balay     ii  = a->i;
34726e093fcSHong Zhang     z   = zarray;
34826e093fcSHong Zhang   }
3492d61bbb3SSatish Balay 
3502d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
3512d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
3522d61bbb3SSatish Balay     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0;
3532d61bbb3SSatish Balay     for (j=0; j<n; j++) {
3542d61bbb3SSatish Balay       xb = x + 3*(*idx++); x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
3552d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[3]*x2 + v[6]*x3;
3562d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[4]*x2 + v[7]*x3;
3572d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[5]*x2 + v[8]*x3;
3582d61bbb3SSatish Balay       v += 9;
3592d61bbb3SSatish Balay     }
3607b2bb3b9SHong Zhang     if (usecprow) z = zarray + 3*ridx[i];
3612d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3;
36226e093fcSHong Zhang     if (!usecprow) z += 3;
3632d61bbb3SSatish Balay   }
3641ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
36526e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
36626e093fcSHong Zhang   PetscLogFlops(18*a->nz - 3*mbs);
3672d61bbb3SSatish Balay   PetscFunctionReturn(0);
3682d61bbb3SSatish Balay }
3692d61bbb3SSatish Balay 
3704a2ae208SSatish Balay #undef __FUNCT__
3714a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_4"
372dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_4(Mat A,Vec xx,Vec zz)
3732d61bbb3SSatish Balay {
3742d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
375*dcf5cc72SBarry Smith   PetscScalar    *x,*z = 0,*xb,sum1,sum2,sum3,sum4,x1,x2,x3,x4,*zarray;
3763f1db9ecSBarry Smith   MatScalar      *v;
377dfbe8321SBarry Smith   PetscErrorCode ierr;
3787b2bb3b9SHong Zhang   PetscInt       mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL;
37926e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
3802d61bbb3SSatish Balay 
3812d61bbb3SSatish Balay   PetscFunctionBegin;
3821ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
38326e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
3842d61bbb3SSatish Balay 
3852d61bbb3SSatish Balay   idx = a->j;
3862d61bbb3SSatish Balay   v   = a->a;
38726e093fcSHong Zhang   if (usecprow){
38826e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
38926e093fcSHong Zhang     ii   = a->compressedrow.i;
3907b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
39126e093fcSHong Zhang   } else {
39226e093fcSHong Zhang     mbs = a->mbs;
3932d61bbb3SSatish Balay     ii  = a->i;
39426e093fcSHong Zhang     z   = zarray;
39526e093fcSHong Zhang   }
3962d61bbb3SSatish Balay 
3972d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
3982d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
3992d61bbb3SSatish Balay     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0;
4002d61bbb3SSatish Balay     for (j=0; j<n; j++) {
4012d61bbb3SSatish Balay       xb = x + 4*(*idx++);
4022d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
4032d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[4]*x2 + v[8]*x3   + v[12]*x4;
4042d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[5]*x2 + v[9]*x3   + v[13]*x4;
4052d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[6]*x2 + v[10]*x3  + v[14]*x4;
4062d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[7]*x2 + v[11]*x3  + v[15]*x4;
4072d61bbb3SSatish Balay       v += 16;
4082d61bbb3SSatish Balay     }
4097b2bb3b9SHong Zhang     if (usecprow) z = zarray + 4*ridx[i];
4102d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4;
41126e093fcSHong Zhang     if (!usecprow) z += 4;
4122d61bbb3SSatish Balay   }
4131ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
41426e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
41526e093fcSHong Zhang   PetscLogFlops(32*a->nz - 4*mbs);
4162d61bbb3SSatish Balay   PetscFunctionReturn(0);
4172d61bbb3SSatish Balay }
4182d61bbb3SSatish Balay 
4194a2ae208SSatish Balay #undef __FUNCT__
4204a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_5"
421dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_5(Mat A,Vec xx,Vec zz)
4222d61bbb3SSatish Balay {
4232d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
424*dcf5cc72SBarry Smith   PetscScalar    sum1,sum2,sum3,sum4,sum5,x1,x2,x3,x4,x5,*xb,*z = 0,*x,*zarray;
4253f1db9ecSBarry Smith   MatScalar      *v;
426dfbe8321SBarry Smith   PetscErrorCode ierr;
4277b2bb3b9SHong Zhang   PetscInt       mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL;
42826e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
4292d61bbb3SSatish Balay 
430433994e6SBarry Smith   PetscFunctionBegin;
4311ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
43226e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
4332d61bbb3SSatish Balay 
4342d61bbb3SSatish Balay   idx = a->j;
4352d61bbb3SSatish Balay   v   = a->a;
43626e093fcSHong Zhang   if (usecprow){
43726e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
43826e093fcSHong Zhang     ii   = a->compressedrow.i;
4397b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
44026e093fcSHong Zhang   } else {
44126e093fcSHong Zhang     mbs = a->mbs;
4422d61bbb3SSatish Balay     ii  = a->i;
44326e093fcSHong Zhang     z   = zarray;
44426e093fcSHong Zhang   }
4452d61bbb3SSatish Balay 
4462d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
4472d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
4482d61bbb3SSatish Balay     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0;
4492d61bbb3SSatish Balay     for (j=0; j<n; j++) {
4502d61bbb3SSatish Balay       xb = x + 5*(*idx++);
4512d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4];
4522d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[5]*x2 + v[10]*x3  + v[15]*x4 + v[20]*x5;
4532d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[6]*x2 + v[11]*x3  + v[16]*x4 + v[21]*x5;
4542d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[7]*x2 + v[12]*x3  + v[17]*x4 + v[22]*x5;
4552d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[8]*x2 + v[13]*x3  + v[18]*x4 + v[23]*x5;
4562d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[9]*x2 + v[14]*x3  + v[19]*x4 + v[24]*x5;
4572d61bbb3SSatish Balay       v += 25;
4582d61bbb3SSatish Balay     }
4597b2bb3b9SHong Zhang     if (usecprow) z = zarray + 5*ridx[i];
4602d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5;
46126e093fcSHong Zhang     if (!usecprow) z += 5;
4622d61bbb3SSatish Balay   }
4631ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
46426e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
46526e093fcSHong Zhang   PetscLogFlops(50*a->nz - 5*mbs);
4662d61bbb3SSatish Balay   PetscFunctionReturn(0);
4672d61bbb3SSatish Balay }
4682d61bbb3SSatish Balay 
46915091d37SBarry Smith 
4704a2ae208SSatish Balay #undef __FUNCT__
4714a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_6"
472dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_6(Mat A,Vec xx,Vec zz)
47315091d37SBarry Smith {
47415091d37SBarry Smith   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
475*dcf5cc72SBarry Smith   PetscScalar    *x,*z = 0,*xb,sum1,sum2,sum3,sum4,sum5,sum6;
47626e093fcSHong Zhang   PetscScalar    x1,x2,x3,x4,x5,x6,*zarray;
47715091d37SBarry Smith   MatScalar      *v;
478dfbe8321SBarry Smith   PetscErrorCode ierr;
4797b2bb3b9SHong Zhang   PetscInt       mbs=a->mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL;
48026e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
48115091d37SBarry Smith 
482433994e6SBarry Smith   PetscFunctionBegin;
4831ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
48426e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
48515091d37SBarry Smith 
48615091d37SBarry Smith   idx = a->j;
48715091d37SBarry Smith   v   = a->a;
48826e093fcSHong Zhang   if (usecprow){
48926e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
49026e093fcSHong Zhang     ii   = a->compressedrow.i;
4917b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
49226e093fcSHong Zhang   } else {
49326e093fcSHong Zhang     mbs = a->mbs;
49415091d37SBarry Smith     ii  = a->i;
49526e093fcSHong Zhang     z   = zarray;
49626e093fcSHong Zhang   }
49715091d37SBarry Smith 
49815091d37SBarry Smith   for (i=0; i<mbs; i++) {
49915091d37SBarry Smith     n  = ii[1] - ii[0]; ii++;
50015091d37SBarry Smith     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0;
50115091d37SBarry Smith     for (j=0; j<n; j++) {
50215091d37SBarry Smith       xb = x + 6*(*idx++);
50315091d37SBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5];
50415091d37SBarry Smith       sum1 += v[0]*x1 + v[6]*x2  + v[12]*x3  + v[18]*x4 + v[24]*x5 + v[30]*x6;
50515091d37SBarry Smith       sum2 += v[1]*x1 + v[7]*x2  + v[13]*x3  + v[19]*x4 + v[25]*x5 + v[31]*x6;
50615091d37SBarry Smith       sum3 += v[2]*x1 + v[8]*x2  + v[14]*x3  + v[20]*x4 + v[26]*x5 + v[32]*x6;
50715091d37SBarry Smith       sum4 += v[3]*x1 + v[9]*x2  + v[15]*x3  + v[21]*x4 + v[27]*x5 + v[33]*x6;
50815091d37SBarry Smith       sum5 += v[4]*x1 + v[10]*x2 + v[16]*x3  + v[22]*x4 + v[28]*x5 + v[34]*x6;
50915091d37SBarry Smith       sum6 += v[5]*x1 + v[11]*x2 + v[17]*x3  + v[23]*x4 + v[29]*x5 + v[35]*x6;
51015091d37SBarry Smith       v += 36;
51115091d37SBarry Smith     }
5127b2bb3b9SHong Zhang     if (usecprow) z = zarray + 6*ridx[i];
51315091d37SBarry Smith     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6;
51426e093fcSHong Zhang     if (!usecprow) z += 6;
51515091d37SBarry Smith   }
51615091d37SBarry Smith 
5171ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
51826e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
51926e093fcSHong Zhang   PetscLogFlops(72*a->nz - 6*mbs);
52015091d37SBarry Smith   PetscFunctionReturn(0);
52115091d37SBarry Smith }
5224a2ae208SSatish Balay #undef __FUNCT__
5234a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_7"
524dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_7(Mat A,Vec xx,Vec zz)
5252d61bbb3SSatish Balay {
5262d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
527*dcf5cc72SBarry Smith   PetscScalar    *x,*z = 0,*xb,sum1,sum2,sum3,sum4,sum5,sum6,sum7;
52826e093fcSHong Zhang   PetscScalar    x1,x2,x3,x4,x5,x6,x7,*zarray;
5293f1db9ecSBarry Smith   MatScalar      *v;
530dfbe8321SBarry Smith   PetscErrorCode ierr;
5317b2bb3b9SHong Zhang   PetscInt       mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL;
53226e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
5332d61bbb3SSatish Balay 
534433994e6SBarry Smith   PetscFunctionBegin;
5351ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
53626e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
5372d61bbb3SSatish Balay 
5382d61bbb3SSatish Balay   idx = a->j;
5392d61bbb3SSatish Balay   v   = a->a;
54026e093fcSHong Zhang   if (usecprow){
54126e093fcSHong Zhang     mbs    = a->compressedrow.nrows;
54226e093fcSHong Zhang     ii     = a->compressedrow.i;
5437b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
54426e093fcSHong Zhang   } else {
54526e093fcSHong Zhang     mbs = a->mbs;
5462d61bbb3SSatish Balay     ii  = a->i;
54726e093fcSHong Zhang     z   = zarray;
54826e093fcSHong Zhang   }
5492d61bbb3SSatish Balay 
5502d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
5512d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
5522d61bbb3SSatish Balay     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0; sum7 = 0.0;
5532d61bbb3SSatish Balay     for (j=0; j<n; j++) {
5542d61bbb3SSatish Balay       xb = x + 7*(*idx++);
5552d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6];
5562d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[7]*x2  + v[14]*x3  + v[21]*x4 + v[28]*x5 + v[35]*x6 + v[42]*x7;
5572d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[8]*x2  + v[15]*x3  + v[22]*x4 + v[29]*x5 + v[36]*x6 + v[43]*x7;
5582d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[9]*x2  + v[16]*x3  + v[23]*x4 + v[30]*x5 + v[37]*x6 + v[44]*x7;
5592d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[10]*x2 + v[17]*x3  + v[24]*x4 + v[31]*x5 + v[38]*x6 + v[45]*x7;
5602d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[11]*x2 + v[18]*x3  + v[25]*x4 + v[32]*x5 + v[39]*x6 + v[46]*x7;
5612d61bbb3SSatish Balay       sum6 += v[5]*x1 + v[12]*x2 + v[19]*x3  + v[26]*x4 + v[33]*x5 + v[40]*x6 + v[47]*x7;
5622d61bbb3SSatish Balay       sum7 += v[6]*x1 + v[13]*x2 + v[20]*x3  + v[27]*x4 + v[34]*x5 + v[41]*x6 + v[48]*x7;
5632d61bbb3SSatish Balay       v += 49;
5642d61bbb3SSatish Balay     }
5657b2bb3b9SHong Zhang     if (usecprow) z = zarray + 7*ridx[i];
5662d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
56726e093fcSHong Zhang     if (!usecprow) z += 7;
5682d61bbb3SSatish Balay   }
5692d61bbb3SSatish Balay 
5701ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
57126e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
57226e093fcSHong Zhang   PetscLogFlops(98*a->nz - 7*mbs);
5732d61bbb3SSatish Balay   PetscFunctionReturn(0);
5742d61bbb3SSatish Balay }
5752d61bbb3SSatish Balay 
5763f1db9ecSBarry Smith /*
5773f1db9ecSBarry Smith     This will not work with MatScalar == float because it calls the BLAS
5783f1db9ecSBarry Smith */
5794a2ae208SSatish Balay #undef __FUNCT__
5804a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_N"
581dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_N(Mat A,Vec xx,Vec zz)
5822d61bbb3SSatish Balay {
5832d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
584*dcf5cc72SBarry Smith   PetscScalar    *x,*z = 0,*xb,*work,*workt,*zarray;
5853f1db9ecSBarry Smith   MatScalar      *v;
586dfbe8321SBarry Smith   PetscErrorCode ierr;
587521d7252SBarry Smith   PetscInt       mbs=a->mbs,i,*idx,*ii,bs=A->bs,j,n,bs2=a->bs2;
5887b2bb3b9SHong Zhang   PetscInt       ncols,k,*ridx=PETSC_NULL;
58926e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
5902d61bbb3SSatish Balay 
5912d61bbb3SSatish Balay   PetscFunctionBegin;
5921ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
59326e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
5942d61bbb3SSatish Balay 
5952d61bbb3SSatish Balay   idx = a->j;
5962d61bbb3SSatish Balay   v   = a->a;
59726e093fcSHong Zhang   if (usecprow){
59826e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
59926e093fcSHong Zhang     ii   = a->compressedrow.i;
6007b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
60126e093fcSHong Zhang   } else {
60226e093fcSHong Zhang     mbs = a->mbs;
6032d61bbb3SSatish Balay     ii  = a->i;
60426e093fcSHong Zhang     z   = zarray;
60526e093fcSHong Zhang   }
606218c64b6SSatish Balay 
6072d61bbb3SSatish Balay   if (!a->mult_work) {
608273d9f13SBarry Smith     k    = PetscMax(A->m,A->n);
60987828ca2SBarry Smith     ierr = PetscMalloc((k+1)*sizeof(PetscScalar),&a->mult_work);CHKERRQ(ierr);
6102d61bbb3SSatish Balay   }
6112d61bbb3SSatish Balay   work = a->mult_work;
6122d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
6132d61bbb3SSatish Balay     n     = ii[1] - ii[0]; ii++;
6142d61bbb3SSatish Balay     ncols = n*bs;
6152d61bbb3SSatish Balay     workt = work;
6162d61bbb3SSatish Balay     for (j=0; j<n; j++) {
6172d61bbb3SSatish Balay       xb = x + bs*(*idx++);
6182d61bbb3SSatish Balay       for (k=0; k<bs; k++) workt[k] = xb[k];
6192d61bbb3SSatish Balay       workt += bs;
6202d61bbb3SSatish Balay     }
6217b2bb3b9SHong Zhang     if (usecprow) z = zarray + bs*ridx[i];
622737d121aSSatish Balay     Kernel_w_gets_Ar_times_v(bs,ncols,work,v,z);
6233f1db9ecSBarry Smith     /* LAgemv_("N",&bs,&ncols,&_DOne,v,&bs,work,&_One,&_DZero,z,&_One); */
6242d61bbb3SSatish Balay     v += n*bs2;
62526e093fcSHong Zhang     if (!usecprow) z += bs;
6262d61bbb3SSatish Balay   }
6271ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
62826e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
62926e093fcSHong Zhang   PetscLogFlops(2*a->nz*bs2 - bs*mbs);
6302d61bbb3SSatish Balay   PetscFunctionReturn(0);
6312d61bbb3SSatish Balay }
6322d61bbb3SSatish Balay 
6334a2ae208SSatish Balay #undef __FUNCT__
6344a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_1"
635dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_1(Mat A,Vec xx,Vec yy,Vec zz)
6362d61bbb3SSatish Balay {
6372d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
638*dcf5cc72SBarry Smith   PetscScalar    *x,*y = 0,*z = 0,sum,*yarray,*zarray;
6393f1db9ecSBarry Smith   MatScalar      *v;
640dfbe8321SBarry Smith   PetscErrorCode ierr;
6417b2bb3b9SHong Zhang   PetscInt       mbs,i,*idx,*ii,n,*ridx=PETSC_NULL;
64226e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
6432d61bbb3SSatish Balay 
6442d61bbb3SSatish Balay   PetscFunctionBegin;
6451ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
64626e093fcSHong Zhang   ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr);
6472e8a6d31SBarry Smith   if (zz != yy) {
64826e093fcSHong Zhang     ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
6492e8a6d31SBarry Smith   } else {
65026e093fcSHong Zhang     zarray = yarray;
6512e8a6d31SBarry Smith   }
6522d61bbb3SSatish Balay 
6532d61bbb3SSatish Balay   idx = a->j;
6542d61bbb3SSatish Balay   v   = a->a;
65526e093fcSHong Zhang   if (usecprow){
65626e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
65726e093fcSHong Zhang     ii   = a->compressedrow.i;
6587b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
65926e093fcSHong Zhang   } else {
66026e093fcSHong Zhang     mbs = a->mbs;
6612d61bbb3SSatish Balay     ii  = a->i;
66226e093fcSHong Zhang     y   = yarray;
66326e093fcSHong Zhang     z   = zarray;
66426e093fcSHong Zhang   }
6652d61bbb3SSatish Balay 
6662d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
6672d61bbb3SSatish Balay     n    = ii[1] - ii[0]; ii++;
66826e093fcSHong Zhang     if (usecprow){
6697b2bb3b9SHong Zhang       z = zarray + ridx[i];
6707b2bb3b9SHong Zhang       y = yarray + ridx[i];
67126e093fcSHong Zhang     }
67226e093fcSHong Zhang     sum = y[0];
6732d61bbb3SSatish Balay     while (n--) sum += *v++ * x[*idx++];
67426e093fcSHong Zhang     z[0] = sum;
67526e093fcSHong Zhang     if (!usecprow){
67626e093fcSHong Zhang       z++; y++;
67726e093fcSHong Zhang     }
6782d61bbb3SSatish Balay   }
6791ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
68026e093fcSHong Zhang   ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr);
6812e8a6d31SBarry Smith   if (zz != yy) {
68226e093fcSHong Zhang     ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
6832e8a6d31SBarry Smith   }
684b0a32e0cSBarry Smith   PetscLogFlops(2*a->nz);
6852d61bbb3SSatish Balay   PetscFunctionReturn(0);
6862d61bbb3SSatish Balay }
6872d61bbb3SSatish Balay 
6884a2ae208SSatish Balay #undef __FUNCT__
6894a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_2"
690dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_2(Mat A,Vec xx,Vec yy,Vec zz)
6912d61bbb3SSatish Balay {
6922d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
693*dcf5cc72SBarry Smith   PetscScalar    *x,*y = 0,*z = 0,*xb,sum1,sum2;
69426e093fcSHong Zhang   PetscScalar    x1,x2,*yarray,*zarray;
6953f1db9ecSBarry Smith   MatScalar      *v;
696dfbe8321SBarry Smith   PetscErrorCode ierr;
6977b2bb3b9SHong Zhang   PetscInt       mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL;
69826e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
6992d61bbb3SSatish Balay 
7002d61bbb3SSatish Balay   PetscFunctionBegin;
7011ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
70226e093fcSHong Zhang   ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr);
7032e8a6d31SBarry Smith   if (zz != yy) {
70426e093fcSHong Zhang     ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
7052e8a6d31SBarry Smith   } else {
70626e093fcSHong Zhang     zarray = yarray;
7072e8a6d31SBarry Smith   }
7082d61bbb3SSatish Balay 
7092d61bbb3SSatish Balay   idx = a->j;
7102d61bbb3SSatish Balay   v   = a->a;
71126e093fcSHong Zhang   if (usecprow){
71226e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
71326e093fcSHong Zhang     ii   = a->compressedrow.i;
7147b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
71526e093fcSHong Zhang   } else {
71626e093fcSHong Zhang     mbs = a->mbs;
7172d61bbb3SSatish Balay     ii  = a->i;
71826e093fcSHong Zhang     y   = yarray;
71926e093fcSHong Zhang     z   = zarray;
72026e093fcSHong Zhang   }
7212d61bbb3SSatish Balay 
7222d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
7232d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
72426e093fcSHong Zhang     if (usecprow){
7257b2bb3b9SHong Zhang       z = zarray + 2*ridx[i];
7267b2bb3b9SHong Zhang       y = yarray + 2*ridx[i];
72726e093fcSHong Zhang     }
7282d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1];
7292d61bbb3SSatish Balay     for (j=0; j<n; j++) {
7302d61bbb3SSatish Balay       xb = x + 2*(*idx++); x1 = xb[0]; x2 = xb[1];
7312d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[2]*x2;
7322d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[3]*x2;
7332d61bbb3SSatish Balay       v += 4;
7342d61bbb3SSatish Balay     }
7352d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2;
73626e093fcSHong Zhang     if (!usecprow){
7372d61bbb3SSatish Balay       z += 2; y += 2;
7382d61bbb3SSatish Balay     }
73926e093fcSHong Zhang   }
7401ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
74126e093fcSHong Zhang   ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr);
7422e8a6d31SBarry Smith   if (zz != yy) {
74326e093fcSHong Zhang     ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
7442e8a6d31SBarry Smith   }
745b0a32e0cSBarry Smith   PetscLogFlops(4*a->nz);
7462d61bbb3SSatish Balay   PetscFunctionReturn(0);
7472d61bbb3SSatish Balay }
7482d61bbb3SSatish Balay 
7494a2ae208SSatish Balay #undef __FUNCT__
7504a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_3"
751dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_3(Mat A,Vec xx,Vec yy,Vec zz)
7522d61bbb3SSatish Balay {
7532d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
754*dcf5cc72SBarry Smith   PetscScalar    *x,*y = 0,*z = 0,*xb,sum1,sum2,sum3,x1,x2,x3,*yarray,*zarray;
7553f1db9ecSBarry Smith   MatScalar      *v;
756dfbe8321SBarry Smith   PetscErrorCode ierr;
7577b2bb3b9SHong Zhang   PetscInt       mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL;
75826e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
7592d61bbb3SSatish Balay 
7602d61bbb3SSatish Balay   PetscFunctionBegin;
7611ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
76226e093fcSHong Zhang   ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr);
7632e8a6d31SBarry Smith   if (zz != yy) {
76426e093fcSHong Zhang     ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
7652e8a6d31SBarry Smith   } else {
76626e093fcSHong Zhang     zarray = yarray;
7672e8a6d31SBarry Smith   }
7682d61bbb3SSatish Balay 
7692d61bbb3SSatish Balay   idx = a->j;
7702d61bbb3SSatish Balay   v   = a->a;
77126e093fcSHong Zhang   if (usecprow){
77226e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
77326e093fcSHong Zhang     ii   = a->compressedrow.i;
7747b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
77526e093fcSHong Zhang   } else {
77626e093fcSHong Zhang     mbs = a->mbs;
7772d61bbb3SSatish Balay     ii  = a->i;
77826e093fcSHong Zhang     y   = yarray;
77926e093fcSHong Zhang     z   = zarray;
78026e093fcSHong Zhang   }
7812d61bbb3SSatish Balay 
7822d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
7832d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
78426e093fcSHong Zhang     if (usecprow){
7857b2bb3b9SHong Zhang       z = zarray + 3*ridx[i];
7867b2bb3b9SHong Zhang       y = yarray + 3*ridx[i];
78726e093fcSHong Zhang     }
7882d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2];
7892d61bbb3SSatish Balay     for (j=0; j<n; j++) {
7902d61bbb3SSatish Balay       xb = x + 3*(*idx++); x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
7912d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[3]*x2 + v[6]*x3;
7922d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[4]*x2 + v[7]*x3;
7932d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[5]*x2 + v[8]*x3;
7942d61bbb3SSatish Balay       v += 9;
7952d61bbb3SSatish Balay     }
7962d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3;
79726e093fcSHong Zhang     if (!usecprow){
7982d61bbb3SSatish Balay       z += 3; y += 3;
7992d61bbb3SSatish Balay     }
80026e093fcSHong Zhang   }
8011ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
80226e093fcSHong Zhang   ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr);
8032e8a6d31SBarry Smith   if (zz != yy) {
80426e093fcSHong Zhang     ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
8052e8a6d31SBarry Smith   }
806b0a32e0cSBarry Smith   PetscLogFlops(18*a->nz);
8072d61bbb3SSatish Balay   PetscFunctionReturn(0);
8082d61bbb3SSatish Balay }
8092d61bbb3SSatish Balay 
8104a2ae208SSatish Balay #undef __FUNCT__
8114a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_4"
812dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_4(Mat A,Vec xx,Vec yy,Vec zz)
8132d61bbb3SSatish Balay {
8142d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
815*dcf5cc72SBarry Smith   PetscScalar    *x,*y = 0,*z = 0,*xb,sum1,sum2,sum3,sum4,x1,x2,x3,x4,*yarray,*zarray;
8163f1db9ecSBarry Smith   MatScalar      *v;
817dfbe8321SBarry Smith   PetscErrorCode ierr;
8187b2bb3b9SHong Zhang   PetscInt       mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL;
81926e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
8202d61bbb3SSatish Balay 
8212d61bbb3SSatish Balay   PetscFunctionBegin;
8221ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
82326e093fcSHong Zhang   ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr);
8242e8a6d31SBarry Smith   if (zz != yy) {
82526e093fcSHong Zhang     ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
8262e8a6d31SBarry Smith   } else {
82726e093fcSHong Zhang     zarray = yarray;
8282e8a6d31SBarry Smith   }
8292d61bbb3SSatish Balay 
8302d61bbb3SSatish Balay   idx   = a->j;
8312d61bbb3SSatish Balay   v     = a->a;
83226e093fcSHong Zhang   if (usecprow){
83326e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
83426e093fcSHong Zhang     ii   = a->compressedrow.i;
8357b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
83626e093fcSHong Zhang   } else {
83726e093fcSHong Zhang     mbs = a->mbs;
8382d61bbb3SSatish Balay     ii  = a->i;
83926e093fcSHong Zhang     y   = yarray;
84026e093fcSHong Zhang     z   = zarray;
84126e093fcSHong Zhang   }
8422d61bbb3SSatish Balay 
8432d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
8442d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
84526e093fcSHong Zhang     if (usecprow){
8467b2bb3b9SHong Zhang       z = zarray + 4*ridx[i];
8477b2bb3b9SHong Zhang       y = yarray + 4*ridx[i];
84826e093fcSHong Zhang     }
8492d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3];
8502d61bbb3SSatish Balay     for (j=0; j<n; j++) {
8512d61bbb3SSatish Balay       xb = x + 4*(*idx++);
8522d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
8532d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[4]*x2 + v[8]*x3   + v[12]*x4;
8542d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[5]*x2 + v[9]*x3   + v[13]*x4;
8552d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[6]*x2 + v[10]*x3  + v[14]*x4;
8562d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[7]*x2 + v[11]*x3  + v[15]*x4;
8572d61bbb3SSatish Balay       v += 16;
8582d61bbb3SSatish Balay     }
8592d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4;
86026e093fcSHong Zhang     if (!usecprow){
8612d61bbb3SSatish Balay       z += 4; y += 4;
8622d61bbb3SSatish Balay     }
86326e093fcSHong Zhang   }
8641ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
86526e093fcSHong Zhang   ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr);
8662e8a6d31SBarry Smith   if (zz != yy) {
86726e093fcSHong Zhang     ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
8682e8a6d31SBarry Smith   }
869b0a32e0cSBarry Smith   PetscLogFlops(32*a->nz);
8702d61bbb3SSatish Balay   PetscFunctionReturn(0);
8712d61bbb3SSatish Balay }
8722d61bbb3SSatish Balay 
8734a2ae208SSatish Balay #undef __FUNCT__
8744a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_5"
875dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_5(Mat A,Vec xx,Vec yy,Vec zz)
8762d61bbb3SSatish Balay {
8772d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
878*dcf5cc72SBarry Smith   PetscScalar    *x,*y = 0,*z = 0,*xb,sum1,sum2,sum3,sum4,sum5,x1,x2,x3,x4,x5;
87926e093fcSHong Zhang   PetscScalar    *yarray,*zarray;
8803f1db9ecSBarry Smith   MatScalar      *v;
881dfbe8321SBarry Smith   PetscErrorCode ierr;
8827b2bb3b9SHong Zhang   PetscInt       mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL;
88326e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
8842d61bbb3SSatish Balay 
8852d61bbb3SSatish Balay   PetscFunctionBegin;
8861ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
88726e093fcSHong Zhang   ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr);
8882e8a6d31SBarry Smith   if (zz != yy) {
88926e093fcSHong Zhang     ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
8902e8a6d31SBarry Smith   } else {
89126e093fcSHong Zhang     zarray = yarray;
8922e8a6d31SBarry Smith   }
8932d61bbb3SSatish Balay 
8942d61bbb3SSatish Balay   idx = a->j;
8952d61bbb3SSatish Balay   v   = a->a;
89626e093fcSHong Zhang   if (usecprow){
89726e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
89826e093fcSHong Zhang     ii   = a->compressedrow.i;
8997b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
90026e093fcSHong Zhang   } else {
90126e093fcSHong Zhang     mbs = a->mbs;
9022d61bbb3SSatish Balay     ii  = a->i;
90326e093fcSHong Zhang     y   = yarray;
90426e093fcSHong Zhang     z   = zarray;
90526e093fcSHong Zhang   }
9062d61bbb3SSatish Balay 
9072d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
9082d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
90926e093fcSHong Zhang     if (usecprow){
9107b2bb3b9SHong Zhang       z = zarray + 5*ridx[i];
9117b2bb3b9SHong Zhang       y = yarray + 5*ridx[i];
91226e093fcSHong Zhang     }
9132d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4];
9142d61bbb3SSatish Balay     for (j=0; j<n; j++) {
9152d61bbb3SSatish Balay       xb = x + 5*(*idx++);
9162d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4];
9172d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[5]*x2 + v[10]*x3  + v[15]*x4 + v[20]*x5;
9182d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[6]*x2 + v[11]*x3  + v[16]*x4 + v[21]*x5;
9192d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[7]*x2 + v[12]*x3  + v[17]*x4 + v[22]*x5;
9202d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[8]*x2 + v[13]*x3  + v[18]*x4 + v[23]*x5;
9212d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[9]*x2 + v[14]*x3  + v[19]*x4 + v[24]*x5;
9222d61bbb3SSatish Balay       v += 25;
9232d61bbb3SSatish Balay     }
9242d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5;
92526e093fcSHong Zhang     if (!usecprow){
9262d61bbb3SSatish Balay       z += 5; y += 5;
9272d61bbb3SSatish Balay     }
92826e093fcSHong Zhang   }
9291ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
93026e093fcSHong Zhang   ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr);
9312e8a6d31SBarry Smith   if (zz != yy) {
93226e093fcSHong Zhang     ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
9332e8a6d31SBarry Smith   }
934b0a32e0cSBarry Smith   PetscLogFlops(50*a->nz);
9352d61bbb3SSatish Balay   PetscFunctionReturn(0);
9362d61bbb3SSatish Balay }
9374a2ae208SSatish Balay #undef __FUNCT__
9384a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_6"
939dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_6(Mat A,Vec xx,Vec yy,Vec zz)
94015091d37SBarry Smith {
94115091d37SBarry Smith   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
942*dcf5cc72SBarry Smith   PetscScalar    *x,*y = 0,*z = 0,*xb,sum1,sum2,sum3,sum4,sum5,sum6;
94326e093fcSHong Zhang   PetscScalar    x1,x2,x3,x4,x5,x6,*yarray,*zarray;
94415091d37SBarry Smith   MatScalar      *v;
945dfbe8321SBarry Smith   PetscErrorCode ierr;
9467b2bb3b9SHong Zhang   PetscInt       mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL;
94726e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
94815091d37SBarry Smith 
94915091d37SBarry Smith   PetscFunctionBegin;
9501ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
95126e093fcSHong Zhang   ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr);
95215091d37SBarry Smith   if (zz != yy) {
95326e093fcSHong Zhang     ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
95415091d37SBarry Smith   } else {
95526e093fcSHong Zhang     zarray = yarray;
95615091d37SBarry Smith   }
95715091d37SBarry Smith 
95815091d37SBarry Smith   idx = a->j;
95915091d37SBarry Smith   v   = a->a;
96026e093fcSHong Zhang   if (usecprow){
96126e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
96226e093fcSHong Zhang     ii   = a->compressedrow.i;
9637b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
96426e093fcSHong Zhang   } else {
96526e093fcSHong Zhang     mbs = a->mbs;
96615091d37SBarry Smith     ii  = a->i;
96726e093fcSHong Zhang     y   = yarray;
96826e093fcSHong Zhang     z   = zarray;
96926e093fcSHong Zhang   }
97015091d37SBarry Smith 
97115091d37SBarry Smith   for (i=0; i<mbs; i++) {
97215091d37SBarry Smith     n  = ii[1] - ii[0]; ii++;
97326e093fcSHong Zhang     if (usecprow){
9747b2bb3b9SHong Zhang       z = zarray + 6*ridx[i];
9757b2bb3b9SHong Zhang       y = yarray + 6*ridx[i];
97626e093fcSHong Zhang     }
97715091d37SBarry Smith     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4]; sum6 = y[5];
97815091d37SBarry Smith     for (j=0; j<n; j++) {
9793b95cb0eSSatish Balay       xb = x + 6*(*idx++);
98015091d37SBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5];
98115091d37SBarry Smith       sum1 += v[0]*x1 + v[6]*x2  + v[12]*x3  + v[18]*x4 + v[24]*x5 + v[30]*x6;
98215091d37SBarry Smith       sum2 += v[1]*x1 + v[7]*x2  + v[13]*x3  + v[19]*x4 + v[25]*x5 + v[31]*x6;
98315091d37SBarry Smith       sum3 += v[2]*x1 + v[8]*x2  + v[14]*x3  + v[20]*x4 + v[26]*x5 + v[32]*x6;
98415091d37SBarry Smith       sum4 += v[3]*x1 + v[9]*x2  + v[15]*x3  + v[21]*x4 + v[27]*x5 + v[33]*x6;
98515091d37SBarry Smith       sum5 += v[4]*x1 + v[10]*x2 + v[16]*x3  + v[22]*x4 + v[28]*x5 + v[34]*x6;
98615091d37SBarry Smith       sum6 += v[5]*x1 + v[11]*x2 + v[17]*x3  + v[23]*x4 + v[29]*x5 + v[35]*x6;
98715091d37SBarry Smith       v += 36;
98815091d37SBarry Smith     }
98915091d37SBarry Smith     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6;
99026e093fcSHong Zhang     if (!usecprow){
99115091d37SBarry Smith       z += 6; y += 6;
99215091d37SBarry Smith     }
99326e093fcSHong Zhang   }
9941ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
99526e093fcSHong Zhang   ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr);
99615091d37SBarry Smith   if (zz != yy) {
99726e093fcSHong Zhang     ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
99815091d37SBarry Smith   }
999b0a32e0cSBarry Smith   PetscLogFlops(72*a->nz);
100015091d37SBarry Smith   PetscFunctionReturn(0);
100115091d37SBarry Smith }
10022d61bbb3SSatish Balay 
10034a2ae208SSatish Balay #undef __FUNCT__
10044a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_7"
1005dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_7(Mat A,Vec xx,Vec yy,Vec zz)
10062d61bbb3SSatish Balay {
10072d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
1008*dcf5cc72SBarry Smith   PetscScalar    *x,*y = 0,*z = 0,*xb,sum1,sum2,sum3,sum4,sum5,sum6,sum7;
100926e093fcSHong Zhang   PetscScalar    x1,x2,x3,x4,x5,x6,x7,*yarray,*zarray;
10103f1db9ecSBarry Smith   MatScalar      *v;
1011dfbe8321SBarry Smith   PetscErrorCode ierr;
10127b2bb3b9SHong Zhang   PetscInt       mbs=a->mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL;
101326e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
10142d61bbb3SSatish Balay 
10152d61bbb3SSatish Balay   PetscFunctionBegin;
10161ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
101726e093fcSHong Zhang   ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr);
10182e8a6d31SBarry Smith   if (zz != yy) {
101926e093fcSHong Zhang     ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
10202e8a6d31SBarry Smith   } else {
102126e093fcSHong Zhang     zarray = yarray;
10222e8a6d31SBarry Smith   }
10232d61bbb3SSatish Balay 
10242d61bbb3SSatish Balay   idx = a->j;
10252d61bbb3SSatish Balay   v   = a->a;
102626e093fcSHong Zhang   if (usecprow){
102726e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
102826e093fcSHong Zhang     ii   = a->compressedrow.i;
10297b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
103026e093fcSHong Zhang   } else {
103126e093fcSHong Zhang     mbs = a->mbs;
10322d61bbb3SSatish Balay     ii  = a->i;
103326e093fcSHong Zhang     y   = yarray;
103426e093fcSHong Zhang     z   = zarray;
103526e093fcSHong Zhang   }
10362d61bbb3SSatish Balay 
10372d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
10382d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
103926e093fcSHong Zhang     if (usecprow){
10407b2bb3b9SHong Zhang       z = zarray + 7*ridx[i];
10417b2bb3b9SHong Zhang       y = yarray + 7*ridx[i];
104226e093fcSHong Zhang     }
10432d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4]; sum6 = y[5]; sum7 = y[6];
10442d61bbb3SSatish Balay     for (j=0; j<n; j++) {
10452d61bbb3SSatish Balay       xb = x + 7*(*idx++);
10462d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6];
10472d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[7]*x2  + v[14]*x3  + v[21]*x4 + v[28]*x5 + v[35]*x6 + v[42]*x7;
10482d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[8]*x2  + v[15]*x3  + v[22]*x4 + v[29]*x5 + v[36]*x6 + v[43]*x7;
10492d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[9]*x2  + v[16]*x3  + v[23]*x4 + v[30]*x5 + v[37]*x6 + v[44]*x7;
10502d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[10]*x2 + v[17]*x3  + v[24]*x4 + v[31]*x5 + v[38]*x6 + v[45]*x7;
10512d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[11]*x2 + v[18]*x3  + v[25]*x4 + v[32]*x5 + v[39]*x6 + v[46]*x7;
10522d61bbb3SSatish Balay       sum6 += v[5]*x1 + v[12]*x2 + v[19]*x3  + v[26]*x4 + v[33]*x5 + v[40]*x6 + v[47]*x7;
10532d61bbb3SSatish Balay       sum7 += v[6]*x1 + v[13]*x2 + v[20]*x3  + v[27]*x4 + v[34]*x5 + v[41]*x6 + v[48]*x7;
10542d61bbb3SSatish Balay       v += 49;
10552d61bbb3SSatish Balay     }
10562d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
105726e093fcSHong Zhang     if (!usecprow){
10582d61bbb3SSatish Balay       z += 7; y += 7;
10592d61bbb3SSatish Balay     }
106026e093fcSHong Zhang   }
10611ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
106226e093fcSHong Zhang   ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr);
10632e8a6d31SBarry Smith   if (zz != yy) {
106426e093fcSHong Zhang     ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
10652e8a6d31SBarry Smith   }
1066b0a32e0cSBarry Smith   PetscLogFlops(98*a->nz);
10672d61bbb3SSatish Balay   PetscFunctionReturn(0);
10682d61bbb3SSatish Balay }
1069218c64b6SSatish Balay 
10704a2ae208SSatish Balay #undef __FUNCT__
10714a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_N"
1072dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_N(Mat A,Vec xx,Vec yy,Vec zz)
10732d61bbb3SSatish Balay {
10742d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
1075*dcf5cc72SBarry Smith   PetscScalar    *x,*z = 0,*xb,*work,*workt,*y,*zarray;
10763f1db9ecSBarry Smith   MatScalar      *v;
10776849ba73SBarry Smith   PetscErrorCode ierr;
107826e093fcSHong Zhang   PetscInt       mbs,i,*idx,*ii,bs=A->bs,j,n,bs2=a->bs2;
10797b2bb3b9SHong Zhang   PetscInt       ncols,k,*ridx=PETSC_NULL;
108026e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
1081218c64b6SSatish Balay 
10822d61bbb3SSatish Balay   PetscFunctionBegin;
10831ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
108426e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
1085b2e7d493SHong Zhang   if (zz != yy) {
10861ebc52fbSHong Zhang     ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
108726e093fcSHong Zhang     ierr = PetscMemcpy(zarray,y,yy->n*sizeof(PetscScalar));CHKERRQ(ierr);
10881ebc52fbSHong Zhang     ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
1089b2e7d493SHong Zhang   }
10902d61bbb3SSatish Balay 
10912d61bbb3SSatish Balay   idx = a->j;
10922d61bbb3SSatish Balay   v   = a->a;
109326e093fcSHong Zhang   if (usecprow){
109426e093fcSHong Zhang     mbs    = a->compressedrow.nrows;
109526e093fcSHong Zhang     ii     = a->compressedrow.i;
10967b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
109726e093fcSHong Zhang   } else {
109826e093fcSHong Zhang     mbs = a->mbs;
10992d61bbb3SSatish Balay     ii  = a->i;
110026e093fcSHong Zhang     z   = zarray;
110126e093fcSHong Zhang   }
11022d61bbb3SSatish Balay 
11032d61bbb3SSatish Balay   if (!a->mult_work) {
1104273d9f13SBarry Smith     k    = PetscMax(A->m,A->n);
110587828ca2SBarry Smith     ierr = PetscMalloc((k+1)*sizeof(PetscScalar),&a->mult_work);CHKERRQ(ierr);
11062d61bbb3SSatish Balay   }
11072d61bbb3SSatish Balay   work = a->mult_work;
11082d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
11092d61bbb3SSatish Balay     n     = ii[1] - ii[0]; ii++;
11102d61bbb3SSatish Balay     ncols = n*bs;
11112d61bbb3SSatish Balay     workt = work;
11122d61bbb3SSatish Balay     for (j=0; j<n; j++) {
11132d61bbb3SSatish Balay       xb = x + bs*(*idx++);
11142d61bbb3SSatish Balay       for (k=0; k<bs; k++) workt[k] = xb[k];
11152d61bbb3SSatish Balay       workt += bs;
11162d61bbb3SSatish Balay     }
11177b2bb3b9SHong Zhang     if (usecprow) z = zarray + bs*ridx[i];
11183f1db9ecSBarry Smith     Kernel_w_gets_w_plus_Ar_times_v(bs,ncols,work,v,z);
11193f1db9ecSBarry Smith     /* LAgemv_("N",&bs,&ncols,&_DOne,v,&bs,work,&_One,&_DOne,z,&_One); */
11202d61bbb3SSatish Balay     v += n*bs2;
112126e093fcSHong Zhang     if (!usecprow){
11222d61bbb3SSatish Balay       z += bs;
11232d61bbb3SSatish Balay     }
112426e093fcSHong Zhang   }
11251ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
112626e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
1127b0a32e0cSBarry Smith   PetscLogFlops(2*a->nz*bs2);
11282d61bbb3SSatish Balay   PetscFunctionReturn(0);
11292d61bbb3SSatish Balay }
11302d61bbb3SSatish Balay 
11314a2ae208SSatish Balay #undef __FUNCT__
11324a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_SeqBAIJ"
1133dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_SeqBAIJ(Mat A,Vec xx,Vec zz)
11342d61bbb3SSatish Balay {
11352d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
11363447b6efSHong Zhang   PetscScalar    zero = 0.0;
11376849ba73SBarry Smith   PetscErrorCode ierr;
11382d61bbb3SSatish Balay 
11392d61bbb3SSatish Balay   PetscFunctionBegin;
1140f1af5d2fSBarry Smith   ierr = VecSet(&zero,zz);CHKERRQ(ierr);
11413447b6efSHong Zhang   ierr = MatMultTransposeAdd_SeqBAIJ(A,xx,zz,zz);CHKERRQ(ierr);
11422d61bbb3SSatish Balay 
1143b0a32e0cSBarry Smith   PetscLogFlops(2*a->nz*a->bs2 - A->n);
11442d61bbb3SSatish Balay   PetscFunctionReturn(0);
11452d61bbb3SSatish Balay }
11462d61bbb3SSatish Balay 
11474a2ae208SSatish Balay #undef __FUNCT__
11484a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_SeqBAIJ"
1149dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_SeqBAIJ(Mat A,Vec xx,Vec yy,Vec zz)
11502d61bbb3SSatish Balay 
11512d61bbb3SSatish Balay {
11522d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
1153*dcf5cc72SBarry Smith   PetscScalar    *zb,*x,*z,*xb = 0,x1,x2,x3,x4,x5;
11543f1db9ecSBarry Smith   MatScalar      *v;
11556849ba73SBarry Smith   PetscErrorCode ierr;
11567b2bb3b9SHong Zhang   PetscInt       mbs,i,*idx,*ii,rval,bs=A->bs,j,n,bs2=a->bs2,*ib,*ridx=PETSC_NULL;
11573447b6efSHong Zhang   Mat_CompressedRow cprow = a->compressedrow;
11583447b6efSHong Zhang   PetscTruth        usecprow=cprow.use;
11592d61bbb3SSatish Balay 
11602d61bbb3SSatish Balay   PetscFunctionBegin;
1161ef66eb69SBarry Smith   if (yy != zz) { ierr = VecCopy(yy,zz);CHKERRQ(ierr); }
11623447b6efSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
11633447b6efSHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
11642d61bbb3SSatish Balay 
11652d61bbb3SSatish Balay   idx = a->j;
11662d61bbb3SSatish Balay   v   = a->a;
11673447b6efSHong Zhang   if (usecprow){
11683447b6efSHong Zhang     mbs  = cprow.nrows;
11693447b6efSHong Zhang     ii   = cprow.i;
11707b2bb3b9SHong Zhang     ridx = cprow.rindex;
11713447b6efSHong Zhang   } else {
11723447b6efSHong Zhang     mbs=a->mbs;
11732d61bbb3SSatish Balay     ii = a->i;
1174f1af5d2fSBarry Smith     xb = x;
11753447b6efSHong Zhang   }
11762d61bbb3SSatish Balay 
11772d61bbb3SSatish Balay   switch (bs) {
11782d61bbb3SSatish Balay   case 1:
11792d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
11807b2bb3b9SHong Zhang       if (usecprow) xb = x + ridx[i];
1181f1af5d2fSBarry Smith       x1 = xb[0];
11823447b6efSHong Zhang       ib = idx + ii[0];
11833447b6efSHong Zhang       n  = ii[1] - ii[0]; ii++;
11842d61bbb3SSatish Balay       for (j=0; j<n; j++) {
11852d61bbb3SSatish Balay         rval    = ib[j];
1186f1af5d2fSBarry Smith         z[rval] += *v * x1;
1187f1af5d2fSBarry Smith         v++;
11882d61bbb3SSatish Balay       }
11893447b6efSHong Zhang       if (!usecprow) xb++;
11902d61bbb3SSatish Balay     }
11912d61bbb3SSatish Balay     break;
11922d61bbb3SSatish Balay   case 2:
11932d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
11947b2bb3b9SHong Zhang       if (usecprow) xb = x + 2*ridx[i];
1195f1af5d2fSBarry Smith       x1 = xb[0]; x2 = xb[1];
11963447b6efSHong Zhang       ib = idx + ii[0];
11973447b6efSHong Zhang       n  = ii[1] - ii[0]; ii++;
11982d61bbb3SSatish Balay       for (j=0; j<n; j++) {
11992d61bbb3SSatish Balay         rval      = ib[j]*2;
12002d61bbb3SSatish Balay         z[rval++] += v[0]*x1 + v[1]*x2;
12012d61bbb3SSatish Balay         z[rval++] += v[2]*x1 + v[3]*x2;
12022d61bbb3SSatish Balay         v  += 4;
12032d61bbb3SSatish Balay       }
12043447b6efSHong Zhang       if (!usecprow) xb += 2;
12052d61bbb3SSatish Balay     }
12062d61bbb3SSatish Balay     break;
12072d61bbb3SSatish Balay   case 3:
12082d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
12097b2bb3b9SHong Zhang       if (usecprow) xb = x + 3*ridx[i];
1210f1af5d2fSBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
12113447b6efSHong Zhang       ib = idx + ii[0];
12123447b6efSHong Zhang       n  = ii[1] - ii[0]; ii++;
12132d61bbb3SSatish Balay       for (j=0; j<n; j++) {
12142d61bbb3SSatish Balay         rval      = ib[j]*3;
12152d61bbb3SSatish Balay         z[rval++] += v[0]*x1 + v[1]*x2 + v[2]*x3;
12162d61bbb3SSatish Balay         z[rval++] += v[3]*x1 + v[4]*x2 + v[5]*x3;
12172d61bbb3SSatish Balay         z[rval++] += v[6]*x1 + v[7]*x2 + v[8]*x3;
12182d61bbb3SSatish Balay         v  += 9;
12192d61bbb3SSatish Balay       }
12203447b6efSHong Zhang       if (!usecprow) xb += 3;
12212d61bbb3SSatish Balay     }
12222d61bbb3SSatish Balay     break;
12232d61bbb3SSatish Balay   case 4:
12242d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
12257b2bb3b9SHong Zhang       if (usecprow) xb = x + 4*ridx[i];
1226f1af5d2fSBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
12273447b6efSHong Zhang       ib = idx + ii[0];
12283447b6efSHong Zhang       n  = ii[1] - ii[0]; ii++;
12292d61bbb3SSatish Balay       for (j=0; j<n; j++) {
12302d61bbb3SSatish Balay         rval      = ib[j]*4;
12312d61bbb3SSatish Balay         z[rval++] +=  v[0]*x1 +  v[1]*x2 +  v[2]*x3 +  v[3]*x4;
12322d61bbb3SSatish Balay         z[rval++] +=  v[4]*x1 +  v[5]*x2 +  v[6]*x3 +  v[7]*x4;
12332d61bbb3SSatish Balay         z[rval++] +=  v[8]*x1 +  v[9]*x2 + v[10]*x3 + v[11]*x4;
12342d61bbb3SSatish Balay         z[rval++] += v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4;
12352d61bbb3SSatish Balay         v  += 16;
12362d61bbb3SSatish Balay       }
12373447b6efSHong Zhang       if (!usecprow) xb += 4;
12382d61bbb3SSatish Balay     }
12392d61bbb3SSatish Balay     break;
12402d61bbb3SSatish Balay   case 5:
12412d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
12427b2bb3b9SHong Zhang       if (usecprow) xb = x + 5*ridx[i];
1243f1af5d2fSBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
12442d61bbb3SSatish Balay       x4 = xb[3]; x5 = xb[4];
12453447b6efSHong Zhang       ib = idx + ii[0];
12463447b6efSHong Zhang       n  = ii[1] - ii[0]; ii++;
12472d61bbb3SSatish Balay       for (j=0; j<n; j++) {
12482d61bbb3SSatish Balay         rval      = ib[j]*5;
12492d61bbb3SSatish Balay         z[rval++] +=  v[0]*x1 +  v[1]*x2 +  v[2]*x3 +  v[3]*x4 +  v[4]*x5;
12502d61bbb3SSatish Balay         z[rval++] +=  v[5]*x1 +  v[6]*x2 +  v[7]*x3 +  v[8]*x4 +  v[9]*x5;
12512d61bbb3SSatish Balay         z[rval++] += v[10]*x1 + v[11]*x2 + v[12]*x3 + v[13]*x4 + v[14]*x5;
12522d61bbb3SSatish Balay         z[rval++] += v[15]*x1 + v[16]*x2 + v[17]*x3 + v[18]*x4 + v[19]*x5;
12532d61bbb3SSatish Balay         z[rval++] += v[20]*x1 + v[21]*x2 + v[22]*x3 + v[23]*x4 + v[24]*x5;
12542d61bbb3SSatish Balay         v  += 25;
12552d61bbb3SSatish Balay       }
12563447b6efSHong Zhang       if (!usecprow) xb += 5;
12572d61bbb3SSatish Balay     }
12582d61bbb3SSatish Balay     break;
1259f1af5d2fSBarry Smith   default: {      /* block sizes larger then 5 by 5 are handled by BLAS */
1260690b6cddSBarry Smith       PetscInt     ncols,k;
12613447b6efSHong Zhang       PetscScalar  *work,*workt,*xtmp;
12623f1db9ecSBarry Smith 
12632d61bbb3SSatish Balay       if (!a->mult_work) {
1264273d9f13SBarry Smith         k = PetscMax(A->m,A->n);
126587828ca2SBarry Smith         ierr = PetscMalloc((k+1)*sizeof(PetscScalar),&a->mult_work);CHKERRQ(ierr);
12662d61bbb3SSatish Balay       }
12672d61bbb3SSatish Balay       work = a->mult_work;
12683447b6efSHong Zhang       xtmp = x;
12692d61bbb3SSatish Balay       for (i=0; i<mbs; i++) {
12702d61bbb3SSatish Balay         n     = ii[1] - ii[0]; ii++;
12712d61bbb3SSatish Balay         ncols = n*bs;
127287828ca2SBarry Smith         ierr  = PetscMemzero(work,ncols*sizeof(PetscScalar));CHKERRQ(ierr);
12733447b6efSHong Zhang         if (usecprow) {
12747b2bb3b9SHong Zhang           xtmp = x + bs*ridx[i];
12753447b6efSHong Zhang         }
12763447b6efSHong Zhang         Kernel_w_gets_w_plus_trans_Ar_times_v(bs,ncols,xtmp,v,work);
12773447b6efSHong Zhang         /* LAgemv_("T",&bs,&ncols,&_DOne,v,&bs,xtmp,&_One,&_DOne,work,&_One); */
12782d61bbb3SSatish Balay         v += n*bs2;
12793447b6efSHong Zhang         if (!usecprow) xtmp += bs;
12802d61bbb3SSatish Balay         workt = work;
12812d61bbb3SSatish Balay         for (j=0; j<n; j++) {
12822d61bbb3SSatish Balay           zb = z + bs*(*idx++);
12832d61bbb3SSatish Balay           for (k=0; k<bs; k++) zb[k] += workt[k] ;
12842d61bbb3SSatish Balay           workt += bs;
12852d61bbb3SSatish Balay         }
12862d61bbb3SSatish Balay       }
12872d61bbb3SSatish Balay     }
12882d61bbb3SSatish Balay   }
12893447b6efSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
12903447b6efSHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
1291b0a32e0cSBarry Smith   PetscLogFlops(2*a->nz*a->bs2);
12922d61bbb3SSatish Balay   PetscFunctionReturn(0);
12932d61bbb3SSatish Balay }
12942d61bbb3SSatish Balay 
12954a2ae208SSatish Balay #undef __FUNCT__
12964a2ae208SSatish Balay #define __FUNCT__ "MatScale_SeqBAIJ"
1297dfbe8321SBarry Smith PetscErrorCode MatScale_SeqBAIJ(const PetscScalar *alpha,Mat inA)
12982d61bbb3SSatish Balay {
12992d61bbb3SSatish Balay   Mat_SeqBAIJ  *a = (Mat_SeqBAIJ*)inA->data;
1300690b6cddSBarry Smith   PetscInt     totalnz = a->bs2*a->nz;
13013eda8832SBarry Smith #if defined(PETSC_USE_MAT_SINGLE)
1302690b6cddSBarry Smith   PetscInt     i;
130387052302SDinesh Kaushik #else
13044ce68768SBarry Smith   PetscBLASInt tnz = (PetscBLASInt) totalnz,one = 1;
13053eda8832SBarry Smith #endif
13062d61bbb3SSatish Balay 
13072d61bbb3SSatish Balay   PetscFunctionBegin;
13083eda8832SBarry Smith #if defined(PETSC_USE_MAT_SINGLE)
13093eda8832SBarry Smith   for (i=0; i<totalnz; i++) a->a[i] *= *alpha;
13103eda8832SBarry Smith #else
13114ce68768SBarry Smith   BLscal_(&tnz,(PetscScalar*)alpha,a->a,&one);
13123eda8832SBarry Smith #endif
1313b0a32e0cSBarry Smith   PetscLogFlops(totalnz);
13142d61bbb3SSatish Balay   PetscFunctionReturn(0);
13152d61bbb3SSatish Balay }
13162d61bbb3SSatish Balay 
13174a2ae208SSatish Balay #undef __FUNCT__
13184a2ae208SSatish Balay #define __FUNCT__ "MatNorm_SeqBAIJ"
1319dfbe8321SBarry Smith PetscErrorCode MatNorm_SeqBAIJ(Mat A,NormType type,PetscReal *norm)
13202d61bbb3SSatish Balay {
13212d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
13223f1db9ecSBarry Smith   MatScalar   *v = a->a;
1323329f5518SBarry Smith   PetscReal   sum = 0.0;
1324521d7252SBarry Smith   PetscInt    i,j,k,bs = A->bs,nz=a->nz,bs2=a->bs2,k1;
13252d61bbb3SSatish Balay 
13262d61bbb3SSatish Balay   PetscFunctionBegin;
13272d61bbb3SSatish Balay   if (type == NORM_FROBENIUS) {
13282d61bbb3SSatish Balay     for (i=0; i< bs2*nz; i++) {
1329aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
1330329f5518SBarry Smith       sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
13312d61bbb3SSatish Balay #else
13322d61bbb3SSatish Balay       sum += (*v)*(*v); v++;
13332d61bbb3SSatish Balay #endif
13342d61bbb3SSatish Balay     }
13352d61bbb3SSatish Balay     *norm = sqrt(sum);
1336596552b5SBarry Smith   }  else if (type == NORM_INFINITY) { /* maximum row sum */
1337596552b5SBarry Smith     *norm = 0.0;
1338596552b5SBarry Smith     for (k=0; k<bs; k++) {
133974f84c7bSSatish Balay       for (j=0; j<a->mbs; j++) {
1340596552b5SBarry Smith         v = a->a + bs2*a->i[j] + k;
1341596552b5SBarry Smith         sum = 0.0;
1342596552b5SBarry Smith         for (i=0; i<a->i[j+1]-a->i[j]; i++) {
13430e90e235SBarry Smith           for (k1=0; k1<bs; k1++){
1344596552b5SBarry Smith             sum += PetscAbsScalar(*v);
1345596552b5SBarry Smith             v   += bs;
13462d61bbb3SSatish Balay           }
13470e90e235SBarry Smith         }
1348596552b5SBarry Smith         if (sum > *norm) *norm = sum;
1349596552b5SBarry Smith       }
1350596552b5SBarry Smith     }
1351596552b5SBarry Smith   } else {
135229bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"No support for this norm yet");
13532d61bbb3SSatish Balay   }
13542d61bbb3SSatish Balay   PetscFunctionReturn(0);
13552d61bbb3SSatish Balay }
13562d61bbb3SSatish Balay 
13572d61bbb3SSatish Balay 
13584a2ae208SSatish Balay #undef __FUNCT__
13594a2ae208SSatish Balay #define __FUNCT__ "MatEqual_SeqBAIJ"
1360dfbe8321SBarry Smith PetscErrorCode MatEqual_SeqBAIJ(Mat A,Mat B,PetscTruth* flg)
13612d61bbb3SSatish Balay {
13622d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ *)A->data,*b = (Mat_SeqBAIJ *)B->data;
1363dfbe8321SBarry Smith   PetscErrorCode ierr;
13642d61bbb3SSatish Balay 
13652d61bbb3SSatish Balay   PetscFunctionBegin;
13662d61bbb3SSatish Balay   /* If the  matrix/block dimensions are not equal, or no of nonzeros or shift */
1367521d7252SBarry Smith   if ((A->m != B->m) || (A->n != B->n) || (A->bs != B->bs)|| (a->nz != b->nz)) {
1368273d9f13SBarry Smith     *flg = PETSC_FALSE;
1369273d9f13SBarry Smith     PetscFunctionReturn(0);
13702d61bbb3SSatish Balay   }
13712d61bbb3SSatish Balay 
13722d61bbb3SSatish Balay   /* if the a->i are the same */
1373690b6cddSBarry Smith   ierr = PetscMemcmp(a->i,b->i,(a->mbs+1)*sizeof(PetscInt),flg);CHKERRQ(ierr);
13740f5bd95cSBarry Smith   if (*flg == PETSC_FALSE) {
13750f5bd95cSBarry Smith     PetscFunctionReturn(0);
13762d61bbb3SSatish Balay   }
13772d61bbb3SSatish Balay 
13782d61bbb3SSatish Balay   /* if a->j are the same */
1379690b6cddSBarry Smith   ierr = PetscMemcmp(a->j,b->j,(a->nz)*sizeof(PetscInt),flg);CHKERRQ(ierr);
13800f5bd95cSBarry Smith   if (*flg == PETSC_FALSE) {
13810f5bd95cSBarry Smith     PetscFunctionReturn(0);
13822d61bbb3SSatish Balay   }
13832d61bbb3SSatish Balay   /* if a->a are the same */
1384521d7252SBarry Smith   ierr = PetscMemcmp(a->a,b->a,(a->nz)*(A->bs)*(B->bs)*sizeof(PetscScalar),flg);CHKERRQ(ierr);
13852d61bbb3SSatish Balay   PetscFunctionReturn(0);
13862d61bbb3SSatish Balay 
13872d61bbb3SSatish Balay }
13882d61bbb3SSatish Balay 
13894a2ae208SSatish Balay #undef __FUNCT__
13904a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_SeqBAIJ"
1391dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_SeqBAIJ(Mat A,Vec v)
13922d61bbb3SSatish Balay {
13932d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
1394dfbe8321SBarry Smith   PetscErrorCode ierr;
1395690b6cddSBarry Smith   PetscInt       i,j,k,n,row,bs,*ai,*aj,ambs,bs2;
139687828ca2SBarry Smith   PetscScalar    *x,zero = 0.0;
13973f1db9ecSBarry Smith   MatScalar      *aa,*aa_j;
13982d61bbb3SSatish Balay 
13992d61bbb3SSatish Balay   PetscFunctionBegin;
140029bbc08cSBarry Smith   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1401521d7252SBarry Smith   bs   = A->bs;
14022d61bbb3SSatish Balay   aa   = a->a;
14032d61bbb3SSatish Balay   ai   = a->i;
14042d61bbb3SSatish Balay   aj   = a->j;
14052d61bbb3SSatish Balay   ambs = a->mbs;
14062d61bbb3SSatish Balay   bs2  = a->bs2;
14072d61bbb3SSatish Balay 
1408e1311b90SBarry Smith   ierr = VecSet(&zero,v);CHKERRQ(ierr);
14091ebc52fbSHong Zhang   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1410e1311b90SBarry Smith   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
1411273d9f13SBarry Smith   if (n != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
14122d61bbb3SSatish Balay   for (i=0; i<ambs; i++) {
14132d61bbb3SSatish Balay     for (j=ai[i]; j<ai[i+1]; j++) {
14142d61bbb3SSatish Balay       if (aj[j] == i) {
14152d61bbb3SSatish Balay         row  = i*bs;
14162d61bbb3SSatish Balay         aa_j = aa+j*bs2;
14172d61bbb3SSatish Balay         for (k=0; k<bs2; k+=(bs+1),row++) x[row] = aa_j[k];
14182d61bbb3SSatish Balay         break;
14192d61bbb3SSatish Balay       }
14202d61bbb3SSatish Balay     }
14212d61bbb3SSatish Balay   }
14221ebc52fbSHong Zhang   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
14232d61bbb3SSatish Balay   PetscFunctionReturn(0);
14242d61bbb3SSatish Balay }
14252d61bbb3SSatish Balay 
14264a2ae208SSatish Balay #undef __FUNCT__
14274a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_SeqBAIJ"
1428dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_SeqBAIJ(Mat A,Vec ll,Vec rr)
14292d61bbb3SSatish Balay {
14302d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
143187828ca2SBarry Smith   PetscScalar    *l,*r,x,*li,*ri;
14323f1db9ecSBarry Smith   MatScalar      *aa,*v;
1433dfbe8321SBarry Smith   PetscErrorCode ierr;
1434690b6cddSBarry Smith   PetscInt       i,j,k,lm,rn,M,m,n,*ai,*aj,mbs,tmp,bs,bs2;
14352d61bbb3SSatish Balay 
14362d61bbb3SSatish Balay   PetscFunctionBegin;
14372d61bbb3SSatish Balay   ai  = a->i;
14382d61bbb3SSatish Balay   aj  = a->j;
14392d61bbb3SSatish Balay   aa  = a->a;
1440273d9f13SBarry Smith   m   = A->m;
1441273d9f13SBarry Smith   n   = A->n;
1442521d7252SBarry Smith   bs  = A->bs;
14432d61bbb3SSatish Balay   mbs = a->mbs;
14442d61bbb3SSatish Balay   bs2 = a->bs2;
14452d61bbb3SSatish Balay   if (ll) {
14461ebc52fbSHong Zhang     ierr = VecGetArray(ll,&l);CHKERRQ(ierr);
1447e1311b90SBarry Smith     ierr = VecGetLocalSize(ll,&lm);CHKERRQ(ierr);
144829bbc08cSBarry Smith     if (lm != m) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vector wrong length");
14492d61bbb3SSatish Balay     for (i=0; i<mbs; i++) { /* for each block row */
14502d61bbb3SSatish Balay       M  = ai[i+1] - ai[i];
14512d61bbb3SSatish Balay       li = l + i*bs;
14522d61bbb3SSatish Balay       v  = aa + bs2*ai[i];
14532d61bbb3SSatish Balay       for (j=0; j<M; j++) { /* for each block */
14542d61bbb3SSatish Balay         for (k=0; k<bs2; k++) {
14552d61bbb3SSatish Balay           (*v++) *= li[k%bs];
14562d61bbb3SSatish Balay         }
14572d61bbb3SSatish Balay       }
14582d61bbb3SSatish Balay     }
14591ebc52fbSHong Zhang     ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr);
1460b0a32e0cSBarry Smith     PetscLogFlops(a->nz);
14612d61bbb3SSatish Balay   }
14622d61bbb3SSatish Balay 
14632d61bbb3SSatish Balay   if (rr) {
14641ebc52fbSHong Zhang     ierr = VecGetArray(rr,&r);CHKERRQ(ierr);
1465e1311b90SBarry Smith     ierr = VecGetLocalSize(rr,&rn);CHKERRQ(ierr);
146629bbc08cSBarry Smith     if (rn != n) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vector wrong length");
14672d61bbb3SSatish Balay     for (i=0; i<mbs; i++) { /* for each block row */
14682d61bbb3SSatish Balay       M  = ai[i+1] - ai[i];
14692d61bbb3SSatish Balay       v  = aa + bs2*ai[i];
14702d61bbb3SSatish Balay       for (j=0; j<M; j++) { /* for each block */
14712d61bbb3SSatish Balay         ri = r + bs*aj[ai[i]+j];
14722d61bbb3SSatish Balay         for (k=0; k<bs; k++) {
14732d61bbb3SSatish Balay           x = ri[k];
14742d61bbb3SSatish Balay           for (tmp=0; tmp<bs; tmp++) (*v++) *= x;
14752d61bbb3SSatish Balay         }
14762d61bbb3SSatish Balay       }
14772d61bbb3SSatish Balay     }
14781ebc52fbSHong Zhang     ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr);
1479b0a32e0cSBarry Smith     PetscLogFlops(a->nz);
14802d61bbb3SSatish Balay   }
14812d61bbb3SSatish Balay   PetscFunctionReturn(0);
14822d61bbb3SSatish Balay }
14832d61bbb3SSatish Balay 
14842d61bbb3SSatish Balay 
14854a2ae208SSatish Balay #undef __FUNCT__
14864a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_SeqBAIJ"
1487dfbe8321SBarry Smith PetscErrorCode MatGetInfo_SeqBAIJ(Mat A,MatInfoType flag,MatInfo *info)
14882d61bbb3SSatish Balay {
14892d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
14902d61bbb3SSatish Balay 
14912d61bbb3SSatish Balay   PetscFunctionBegin;
1492273d9f13SBarry Smith   info->rows_global    = (double)A->m;
1493273d9f13SBarry Smith   info->columns_global = (double)A->n;
1494273d9f13SBarry Smith   info->rows_local     = (double)A->m;
1495273d9f13SBarry Smith   info->columns_local  = (double)A->n;
14962d61bbb3SSatish Balay   info->block_size     = a->bs2;
14972d61bbb3SSatish Balay   info->nz_allocated   = a->maxnz;
14982d61bbb3SSatish Balay   info->nz_used        = a->bs2*a->nz;
14992d61bbb3SSatish Balay   info->nz_unneeded    = (double)(info->nz_allocated - info->nz_used);
15002d61bbb3SSatish Balay   info->assemblies   = A->num_ass;
15012d61bbb3SSatish Balay   info->mallocs      = a->reallocs;
15022d61bbb3SSatish Balay   info->memory       = A->mem;
15032d61bbb3SSatish Balay   if (A->factor) {
15042d61bbb3SSatish Balay     info->fill_ratio_given  = A->info.fill_ratio_given;
15052d61bbb3SSatish Balay     info->fill_ratio_needed = A->info.fill_ratio_needed;
15062d61bbb3SSatish Balay     info->factor_mallocs    = A->info.factor_mallocs;
15072d61bbb3SSatish Balay   } else {
15082d61bbb3SSatish Balay     info->fill_ratio_given  = 0;
15092d61bbb3SSatish Balay     info->fill_ratio_needed = 0;
15102d61bbb3SSatish Balay     info->factor_mallocs    = 0;
15112d61bbb3SSatish Balay   }
15122d61bbb3SSatish Balay   PetscFunctionReturn(0);
15132d61bbb3SSatish Balay }
15142d61bbb3SSatish Balay 
15152d61bbb3SSatish Balay 
15164a2ae208SSatish Balay #undef __FUNCT__
15174a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_SeqBAIJ"
1518dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_SeqBAIJ(Mat A)
15192d61bbb3SSatish Balay {
15202d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
1521dfbe8321SBarry Smith   PetscErrorCode ierr;
15222d61bbb3SSatish Balay 
15232d61bbb3SSatish Balay   PetscFunctionBegin;
1524549d3d68SSatish Balay   ierr = PetscMemzero(a->a,a->bs2*a->i[a->mbs]*sizeof(MatScalar));CHKERRQ(ierr);
15252d61bbb3SSatish Balay   PetscFunctionReturn(0);
15262d61bbb3SSatish Balay }
1527