xref: /petsc/src/mat/impls/baij/seq/baij2.c (revision fca92195c093e8a5faf0c844190ed016dcbb215a)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
2cac129eeSSatish Balay 
37c4f633dSBarry Smith #include "../src/mat/impls/baij/seq/baij.h"
4c60f0209SBarry Smith #include "../src/mat/blockinvert.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;
135d0c19d7SBarry Smith   PetscInt       row,i,j,k,l,m,n,*nidx,isz,val,ival;
145d0c19d7SBarry Smith   const PetscInt *idx;
15690b6cddSBarry Smith   PetscInt       start,end,*ai,*aj,bs,*nidx2;
16f1af5d2fSBarry Smith   PetscBT        table;
17a3192f15SSatish Balay 
183a40ed3dSBarry Smith   PetscFunctionBegin;
19a3192f15SSatish Balay   m     = a->mbs;
20a3192f15SSatish Balay   ai    = a->i;
21a3192f15SSatish Balay   aj    = a->j;
22d0f46423SBarry Smith   bs    = A->rmap->bs;
23a3192f15SSatish Balay 
2429bbc08cSBarry Smith   if (ov < 0)  SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap specified");
25a3192f15SSatish Balay 
266831982aSBarry Smith   ierr = PetscBTCreate(m,table);CHKERRQ(ierr);
27690b6cddSBarry Smith   ierr = PetscMalloc((m+1)*sizeof(PetscInt),&nidx);CHKERRQ(ierr);
28d0f46423SBarry Smith   ierr = PetscMalloc((A->rmap->N+1)*sizeof(PetscInt),&nidx2);CHKERRQ(ierr);
29a3192f15SSatish Balay 
30a3192f15SSatish Balay   for (i=0; i<is_max; i++) {
31a3192f15SSatish Balay     /* Initialise the two local arrays */
32a3192f15SSatish Balay     isz  = 0;
336831982aSBarry Smith     ierr = PetscBTMemzero(m,table);CHKERRQ(ierr);
34a3192f15SSatish Balay 
35a3192f15SSatish Balay     /* Extract the indices, assume there can be duplicate entries */
36a3192f15SSatish Balay     ierr = ISGetIndices(is[i],&idx);CHKERRQ(ierr);
37b9b97703SBarry Smith     ierr = ISGetLocalSize(is[i],&n);CHKERRQ(ierr);
38a3192f15SSatish Balay 
39a3192f15SSatish Balay     /* Enter these into the temp arrays i.e mark table[row], enter row into new index */
40a3192f15SSatish Balay     for (j=0; j<n ; ++j){
41218c64b6SSatish Balay       ival = idx[j]/bs; /* convert the indices into block indices */
425eee224dSHong Zhang       if (ival>=m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"index greater than mat-dim");
43f1af5d2fSBarry Smith       if(!PetscBTLookupSet(table,ival)) { nidx[isz++] = ival;}
44a3192f15SSatish Balay     }
45a3192f15SSatish Balay     ierr = ISRestoreIndices(is[i],&idx);CHKERRQ(ierr);
46a3192f15SSatish Balay     ierr = ISDestroy(is[i]);CHKERRQ(ierr);
47a3192f15SSatish Balay 
48a3192f15SSatish Balay     k = 0;
49a3192f15SSatish Balay     for (j=0; j<ov; j++){ /* for each overlap*/
50a3192f15SSatish Balay       n = isz;
51a3192f15SSatish Balay       for (; k<n ; k++){ /* do only those rows in nidx[k], which are not done yet */
52a3192f15SSatish Balay         row   = nidx[k];
53a3192f15SSatish Balay         start = ai[row];
54a3192f15SSatish Balay         end   = ai[row+1];
55a3192f15SSatish Balay         for (l = start; l<end ; l++){
56a3192f15SSatish Balay           val = aj[l];
57f1af5d2fSBarry Smith           if (!PetscBTLookupSet(table,val)) {nidx[isz++] = val;}
58a3192f15SSatish Balay         }
59a3192f15SSatish Balay       }
60a3192f15SSatish Balay     }
61218c64b6SSatish Balay     /* expand the Index Set */
62218c64b6SSatish Balay     for (j=0; j<isz; j++) {
63218c64b6SSatish Balay       for (k=0; k<bs; k++)
64218c64b6SSatish Balay         nidx2[j*bs+k] = nidx[j]*bs+k;
65218c64b6SSatish Balay     }
66329f5518SBarry Smith     ierr = ISCreateGeneral(PETSC_COMM_SELF,isz*bs,nidx2,is+i);CHKERRQ(ierr);
67a3192f15SSatish Balay   }
686831982aSBarry Smith   ierr = PetscBTDestroy(table);CHKERRQ(ierr);
69606d414cSSatish Balay   ierr = PetscFree(nidx);CHKERRQ(ierr);
70606d414cSSatish Balay   ierr = PetscFree(nidx2);CHKERRQ(ierr);
713a40ed3dSBarry Smith   PetscFunctionReturn(0);
72a3192f15SSatish Balay }
731c351548SSatish Balay 
744a2ae208SSatish Balay #undef __FUNCT__
754a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqBAIJ_Private"
764aa3045dSJed Brown PetscErrorCode MatGetSubMatrix_SeqBAIJ_Private(Mat A,IS isrow,IS iscol,MatReuse scall,Mat *B)
77736121d4SSatish Balay {
78736121d4SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*c;
796849ba73SBarry Smith   PetscErrorCode ierr;
80690b6cddSBarry Smith   PetscInt       *smap,i,k,kstart,kend,oldcols = a->nbs,*lens;
81690b6cddSBarry Smith   PetscInt       row,mat_i,*mat_j,tcol,*mat_ilen;
825d0c19d7SBarry Smith   const PetscInt *irow,*icol;
835d0c19d7SBarry Smith   PetscInt       nrows,ncols,*ssmap,bs=A->rmap->bs,bs2=a->bs2;
84690b6cddSBarry Smith   PetscInt       *aj = a->j,*ai = a->i;
853f1db9ecSBarry Smith   MatScalar      *mat_a;
86736121d4SSatish Balay   Mat            C;
8714ca34e6SBarry Smith   PetscTruth     flag,sorted;
88736121d4SSatish Balay 
893a40ed3dSBarry Smith   PetscFunctionBegin;
9014ca34e6SBarry Smith   ierr = ISSorted(iscol,&sorted);CHKERRQ(ierr);
9114ca34e6SBarry Smith   if (!sorted) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"IS is not sorted");
92736121d4SSatish Balay 
93736121d4SSatish Balay   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
94218c64b6SSatish Balay   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
95b9b97703SBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
96b9b97703SBarry Smith   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
97736121d4SSatish Balay 
98690b6cddSBarry Smith   ierr = PetscMalloc((1+oldcols)*sizeof(PetscInt),&smap);CHKERRQ(ierr);
99736121d4SSatish Balay   ssmap = smap;
100690b6cddSBarry Smith   ierr = PetscMalloc((1+nrows)*sizeof(PetscInt),&lens);CHKERRQ(ierr);
101690b6cddSBarry Smith   ierr  = PetscMemzero(smap,oldcols*sizeof(PetscInt));CHKERRQ(ierr);
102736121d4SSatish Balay   for (i=0; i<ncols; i++) smap[icol[i]] = i+1;
103736121d4SSatish Balay   /* determine lens of each row */
104736121d4SSatish Balay   for (i=0; i<nrows; i++) {
105736121d4SSatish Balay     kstart  = ai[irow[i]];
106736121d4SSatish Balay     kend    = kstart + a->ilen[irow[i]];
107736121d4SSatish Balay     lens[i] = 0;
108736121d4SSatish Balay       for (k=kstart; k<kend; k++) {
109736121d4SSatish Balay         if (ssmap[aj[k]]) {
110736121d4SSatish Balay           lens[i]++;
111736121d4SSatish Balay         }
112736121d4SSatish Balay       }
113736121d4SSatish Balay     }
114736121d4SSatish Balay   /* Create and fill new matrix */
115736121d4SSatish Balay   if (scall == MAT_REUSE_MATRIX) {
116736121d4SSatish Balay     c = (Mat_SeqBAIJ *)((*B)->data);
117736121d4SSatish Balay 
118d0f46423SBarry Smith     if (c->mbs!=nrows || c->nbs!=ncols || (*B)->rmap->bs!=bs) SETERRQ(PETSC_ERR_ARG_SIZ,"Submatrix wrong size");
119690b6cddSBarry Smith     ierr = PetscMemcmp(c->ilen,lens,c->mbs *sizeof(PetscInt),&flag);CHKERRQ(ierr);
120abc0a331SBarry Smith     if (!flag) {
12129bbc08cSBarry Smith       SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong no of nonzeros");
122736121d4SSatish Balay     }
123690b6cddSBarry Smith     ierr = PetscMemzero(c->ilen,c->mbs*sizeof(PetscInt));CHKERRQ(ierr);
124736121d4SSatish Balay     C = *B;
1253a40ed3dSBarry Smith   } else {
1267adad957SLisandro Dalcin     ierr = MatCreate(((PetscObject)A)->comm,&C);CHKERRQ(ierr);
127f69a0ea3SMatthew Knepley     ierr = MatSetSizes(C,nrows*bs,ncols*bs,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1287adad957SLisandro Dalcin     ierr = MatSetType(C,((PetscObject)A)->type_name);CHKERRQ(ierr);
129ab93d7beSBarry Smith     ierr = MatSeqBAIJSetPreallocation_SeqBAIJ(C,bs,0,lens);CHKERRQ(ierr);
130736121d4SSatish Balay   }
131736121d4SSatish Balay   c = (Mat_SeqBAIJ *)(C->data);
132736121d4SSatish Balay   for (i=0; i<nrows; i++) {
133736121d4SSatish Balay     row    = irow[i];
134736121d4SSatish Balay     kstart = ai[row];
135736121d4SSatish Balay     kend   = kstart + a->ilen[row];
136736121d4SSatish Balay     mat_i  = c->i[i];
137736121d4SSatish Balay     mat_j  = c->j + mat_i;
138218c64b6SSatish Balay     mat_a  = c->a + mat_i*bs2;
139736121d4SSatish Balay     mat_ilen = c->ilen + i;
140736121d4SSatish Balay     for (k=kstart; k<kend; k++) {
141736121d4SSatish Balay       if ((tcol=ssmap[a->j[k]])) {
142736121d4SSatish Balay         *mat_j++ = tcol - 1;
143549d3d68SSatish Balay         ierr     = PetscMemcpy(mat_a,a->a+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
144549d3d68SSatish Balay         mat_a   += bs2;
145736121d4SSatish Balay         (*mat_ilen)++;
146736121d4SSatish Balay       }
147736121d4SSatish Balay     }
148736121d4SSatish Balay   }
149218c64b6SSatish Balay 
150736121d4SSatish Balay   /* Free work space */
151736121d4SSatish Balay   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
152606d414cSSatish Balay   ierr = PetscFree(smap);CHKERRQ(ierr);
153606d414cSSatish Balay   ierr = PetscFree(lens);CHKERRQ(ierr);
1546d4a8577SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1556d4a8577SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
156736121d4SSatish Balay 
157736121d4SSatish Balay   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
158736121d4SSatish Balay   *B = C;
1593a40ed3dSBarry Smith   PetscFunctionReturn(0);
160736121d4SSatish Balay }
161736121d4SSatish Balay 
1624a2ae208SSatish Balay #undef __FUNCT__
1634a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqBAIJ"
1644aa3045dSJed Brown PetscErrorCode MatGetSubMatrix_SeqBAIJ(Mat A,IS isrow,IS iscol,MatReuse scall,Mat *B)
165218c64b6SSatish Balay {
166218c64b6SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
167218c64b6SSatish Balay   IS             is1,is2;
1686849ba73SBarry Smith   PetscErrorCode ierr;
1695d0c19d7SBarry Smith   PetscInt       *vary,*iary,nrows,ncols,i,bs=A->rmap->bs,count;
1705d0c19d7SBarry Smith   const PetscInt *irow,*icol;
171218c64b6SSatish Balay 
1723a40ed3dSBarry Smith   PetscFunctionBegin;
173218c64b6SSatish Balay   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
174218c64b6SSatish Balay   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
175b9b97703SBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
176b9b97703SBarry Smith   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
177218c64b6SSatish Balay 
178218c64b6SSatish Balay   /* Verify if the indices corespond to each element in a block
179218c64b6SSatish Balay    and form the IS with compressed IS */
180*fca92195SBarry Smith   ierr = PetscMalloc2(a->mbs,PetscInt,&vary,a->mbs,PetscInt,&iary);CHKERRQ(ierr);
181*fca92195SBarry Smith   ierr = PetscMemzero(vary,a->mbs*sizeof(PetscInt));CHKERRQ(ierr);
182218c64b6SSatish Balay   for (i=0; i<nrows; i++) vary[irow[i]/bs]++;
183218c64b6SSatish Balay   count = 0;
184218c64b6SSatish Balay   for (i=0; i<a->mbs; i++) {
185634064b4SBarry Smith     if (vary[i]!=0 && vary[i]!=bs) SETERRQ(PETSC_ERR_ARG_SIZ,"Index set does not match blocks");
186218c64b6SSatish Balay     if (vary[i]==bs) iary[count++] = i;
187218c64b6SSatish Balay   }
188029af93fSBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,count,iary,&is1);CHKERRQ(ierr);
189218c64b6SSatish Balay 
190690b6cddSBarry Smith   ierr = PetscMemzero(vary,(a->mbs)*sizeof(PetscInt));CHKERRQ(ierr);
191218c64b6SSatish Balay   for (i=0; i<ncols; i++) vary[icol[i]/bs]++;
192218c64b6SSatish Balay   count = 0;
193218c64b6SSatish Balay   for (i=0; i<a->mbs; i++) {
194634064b4SBarry Smith     if (vary[i]!=0 && vary[i]!=bs) SETERRQ(PETSC_ERR_PLIB,"Internal error in PETSc");
195218c64b6SSatish Balay     if (vary[i]==bs) iary[count++] = i;
196218c64b6SSatish Balay   }
197029af93fSBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,count,iary,&is2);CHKERRQ(ierr);
198218c64b6SSatish Balay   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
199218c64b6SSatish Balay   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
200*fca92195SBarry Smith   ierr = PetscFree2(vary,iary);CHKERRQ(ierr);
201218c64b6SSatish Balay 
2024aa3045dSJed Brown   ierr = MatGetSubMatrix_SeqBAIJ_Private(A,is1,is2,scall,B);CHKERRQ(ierr);
203218c64b6SSatish Balay   ISDestroy(is1);
204218c64b6SSatish Balay   ISDestroy(is2);
2053a40ed3dSBarry Smith   PetscFunctionReturn(0);
206218c64b6SSatish Balay }
207218c64b6SSatish Balay 
2084a2ae208SSatish Balay #undef __FUNCT__
2094a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_SeqBAIJ"
210690b6cddSBarry Smith PetscErrorCode MatGetSubMatrices_SeqBAIJ(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
211736121d4SSatish Balay {
2126849ba73SBarry Smith   PetscErrorCode ierr;
213690b6cddSBarry Smith   PetscInt       i;
214736121d4SSatish Balay 
2153a40ed3dSBarry Smith   PetscFunctionBegin;
216736121d4SSatish Balay   if (scall == MAT_INITIAL_MATRIX) {
21782502324SSatish Balay     ierr = PetscMalloc((n+1)*sizeof(Mat),B);CHKERRQ(ierr);
218736121d4SSatish Balay   }
219736121d4SSatish Balay 
220736121d4SSatish Balay   for (i=0; i<n; i++) {
2214aa3045dSJed Brown     ierr = MatGetSubMatrix_SeqBAIJ(A,irow[i],icol[i],scall,&(*B)[i]);CHKERRQ(ierr);
222736121d4SSatish Balay   }
2233a40ed3dSBarry Smith   PetscFunctionReturn(0);
224736121d4SSatish Balay }
225218c64b6SSatish Balay 
226218c64b6SSatish Balay 
2272d61bbb3SSatish Balay /* -------------------------------------------------------*/
2282d61bbb3SSatish Balay /* Should check that shapes of vectors and matrices match */
2292d61bbb3SSatish Balay /* -------------------------------------------------------*/
230d9eff348SSatish Balay #include "petscblaslapack.h"
2312d61bbb3SSatish Balay 
2324a2ae208SSatish Balay #undef __FUNCT__
2334a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_1"
234dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_1(Mat A,Vec xx,Vec zz)
2352d61bbb3SSatish Balay {
2362d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
237d9fead3dSBarry Smith   PetscScalar       *z,sum;
238d9fead3dSBarry Smith   const PetscScalar *x;
239d9fead3dSBarry Smith   const MatScalar   *v;
2406849ba73SBarry Smith   PetscErrorCode    ierr;
2412162cab8SBarry Smith   PetscInt          mbs,i,n,nonzerorow=0;
2422162cab8SBarry Smith   const PetscInt    *idx,*ii,*ridx=PETSC_NULL;
24326e093fcSHong Zhang   PetscTruth        usecprow=a->compressedrow.use;
2442d61bbb3SSatish Balay 
2452d61bbb3SSatish Balay   PetscFunctionBegin;
246d9fead3dSBarry Smith   ierr = VecGetArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
2471ebc52fbSHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
2482d61bbb3SSatish Balay 
24926e093fcSHong Zhang   if (usecprow){
25026e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
25126e093fcSHong Zhang     ii   = a->compressedrow.i;
2527b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
253a1c3900fSBarry Smith     ierr = PetscMemzero(z,mbs*sizeof(PetscScalar));CHKERRQ(ierr);
25426e093fcSHong Zhang   } else {
25526e093fcSHong Zhang     mbs = a->mbs;
2562d61bbb3SSatish Balay     ii  = a->i;
25726e093fcSHong Zhang   }
2582d61bbb3SSatish Balay 
2592d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
260ee54c7eeSHong Zhang     n    = ii[1] - ii[0];
261ee54c7eeSHong Zhang     v    = a->a + ii[0];
262ee54c7eeSHong Zhang     idx  = a->j + ii[0];
263ee54c7eeSHong Zhang     ii++;
2642d61bbb3SSatish Balay     sum  = 0.0;
2652162cab8SBarry Smith     PetscSparseDensePlusDot(sum,x,v,idx,n);
26626e093fcSHong Zhang     if (usecprow){
2677b2bb3b9SHong Zhang       z[ridx[i]] = sum;
26826e093fcSHong Zhang     } else {
269122f12eaSBarry Smith       nonzerorow += (n>0);
2702d61bbb3SSatish Balay       z[i] = sum;
2712d61bbb3SSatish Balay     }
27226e093fcSHong Zhang   }
273d9fead3dSBarry Smith   ierr = VecRestoreArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
2741ebc52fbSHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
275dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*a->nz - nonzerorow);CHKERRQ(ierr);
2762d61bbb3SSatish Balay   PetscFunctionReturn(0);
2772d61bbb3SSatish Balay }
2782d61bbb3SSatish Balay 
2794a2ae208SSatish Balay #undef __FUNCT__
2804a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_2"
281dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_2(Mat A,Vec xx,Vec zz)
2822d61bbb3SSatish Balay {
2832d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
284d9fead3dSBarry Smith   PetscScalar       *z = 0,sum1,sum2,*zarray;
285d9fead3dSBarry Smith   const PetscScalar *x,*xb;
28687828ca2SBarry Smith   PetscScalar       x1,x2;
287d9fead3dSBarry Smith   const MatScalar   *v;
288dfbe8321SBarry Smith   PetscErrorCode    ierr;
28998c9bda7SSatish Balay   PetscInt          mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL,nonzerorow=0;
29026e093fcSHong Zhang   PetscTruth        usecprow=a->compressedrow.use;
2912d61bbb3SSatish Balay 
2922d61bbb3SSatish Balay   PetscFunctionBegin;
293d9fead3dSBarry Smith   ierr = VecGetArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
29426e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
2952d61bbb3SSatish Balay 
2962d61bbb3SSatish Balay   idx = a->j;
2972d61bbb3SSatish Balay   v   = a->a;
29826e093fcSHong Zhang   if (usecprow){
29926e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
30026e093fcSHong Zhang     ii   = a->compressedrow.i;
3017b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
30226e093fcSHong Zhang   } else {
30326e093fcSHong Zhang     mbs = a->mbs;
3042d61bbb3SSatish Balay     ii  = a->i;
30526e093fcSHong Zhang     z   = zarray;
30626e093fcSHong Zhang   }
3072d61bbb3SSatish Balay 
3082d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
3092d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
3102d61bbb3SSatish Balay     sum1 = 0.0; sum2 = 0.0;
31198c9bda7SSatish Balay     nonzerorow += (n>0);
3122d61bbb3SSatish Balay     for (j=0; j<n; j++) {
3132d61bbb3SSatish Balay       xb = x + 2*(*idx++); x1 = xb[0]; x2 = xb[1];
3142d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[2]*x2;
3152d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[3]*x2;
3162d61bbb3SSatish Balay       v += 4;
3172d61bbb3SSatish Balay     }
3187b2bb3b9SHong Zhang     if (usecprow) z = zarray + 2*ridx[i];
3192d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2;
32026e093fcSHong Zhang     if (!usecprow) z += 2;
3212d61bbb3SSatish Balay   }
322d9fead3dSBarry Smith   ierr = VecRestoreArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
32326e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
324dc0b31edSSatish Balay   ierr = PetscLogFlops(8.0*a->nz - 2.0*nonzerorow);CHKERRQ(ierr);
3252d61bbb3SSatish Balay   PetscFunctionReturn(0);
3262d61bbb3SSatish Balay }
3272d61bbb3SSatish Balay 
3284a2ae208SSatish Balay #undef __FUNCT__
3294a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_3"
330dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_3(Mat A,Vec xx,Vec zz)
3312d61bbb3SSatish Balay {
3322d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
333d9fead3dSBarry Smith   PetscScalar       *z = 0,sum1,sum2,sum3,x1,x2,x3,*zarray;
334d9fead3dSBarry Smith   const PetscScalar *x,*xb;
335d9fead3dSBarry Smith   const MatScalar   *v;
336dfbe8321SBarry Smith   PetscErrorCode    ierr;
337028cd4eaSSatish Balay   PetscInt          mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL,nonzerorow=0;
338028cd4eaSSatish Balay   PetscTruth        usecprow=a->compressedrow.use;
33926e093fcSHong Zhang 
3402d61bbb3SSatish Balay 
341b6410449SSatish Balay #if defined(PETSC_HAVE_PRAGMA_DISJOINT)
342fee21e36SBarry Smith #pragma disjoint(*v,*z,*xb)
343fee21e36SBarry Smith #endif
344fee21e36SBarry Smith 
3452d61bbb3SSatish Balay   PetscFunctionBegin;
346d9fead3dSBarry Smith   ierr = VecGetArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
34726e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
3482d61bbb3SSatish Balay 
3492d61bbb3SSatish Balay   idx = a->j;
3502d61bbb3SSatish Balay   v   = a->a;
35126e093fcSHong Zhang   if (usecprow){
35226e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
35326e093fcSHong Zhang     ii   = a->compressedrow.i;
3547b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
35526e093fcSHong Zhang   } else {
35626e093fcSHong Zhang     mbs = a->mbs;
3572d61bbb3SSatish Balay     ii  = a->i;
35826e093fcSHong Zhang     z   = zarray;
35926e093fcSHong Zhang   }
3602d61bbb3SSatish Balay 
3612d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
3622d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
3632d61bbb3SSatish Balay     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0;
36498c9bda7SSatish Balay     nonzerorow += (n>0);
3652d61bbb3SSatish Balay     for (j=0; j<n; j++) {
3662d61bbb3SSatish Balay       xb = x + 3*(*idx++); x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
3672d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[3]*x2 + v[6]*x3;
3682d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[4]*x2 + v[7]*x3;
3692d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[5]*x2 + v[8]*x3;
3702d61bbb3SSatish Balay       v += 9;
3712d61bbb3SSatish Balay     }
3727b2bb3b9SHong Zhang     if (usecprow) z = zarray + 3*ridx[i];
3732d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3;
37426e093fcSHong Zhang     if (!usecprow) z += 3;
3752d61bbb3SSatish Balay   }
376d9fead3dSBarry Smith   ierr = VecRestoreArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
37726e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
378dc0b31edSSatish Balay   ierr = PetscLogFlops(18.0*a->nz - 3.0*nonzerorow);CHKERRQ(ierr);
3792d61bbb3SSatish Balay   PetscFunctionReturn(0);
3802d61bbb3SSatish Balay }
3812d61bbb3SSatish Balay 
3824a2ae208SSatish Balay #undef __FUNCT__
3834a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_4"
384dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_4(Mat A,Vec xx,Vec zz)
3852d61bbb3SSatish Balay {
3862d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
387d9fead3dSBarry Smith   PetscScalar       *z = 0,sum1,sum2,sum3,sum4,x1,x2,x3,x4,*zarray;
388d9fead3dSBarry Smith   const PetscScalar *x,*xb;
389d9fead3dSBarry Smith   const MatScalar   *v;
390dfbe8321SBarry Smith   PetscErrorCode    ierr;
39198c9bda7SSatish Balay   PetscInt          mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL,nonzerorow=0;
39226e093fcSHong Zhang   PetscTruth        usecprow=a->compressedrow.use;
3932d61bbb3SSatish Balay 
3942d61bbb3SSatish Balay   PetscFunctionBegin;
395d9fead3dSBarry Smith   ierr = VecGetArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
39626e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
3972d61bbb3SSatish Balay 
3982d61bbb3SSatish Balay   idx = a->j;
3992d61bbb3SSatish Balay   v   = a->a;
40026e093fcSHong Zhang   if (usecprow){
40126e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
40226e093fcSHong Zhang     ii   = a->compressedrow.i;
4037b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
40426e093fcSHong Zhang   } else {
40526e093fcSHong Zhang     mbs = a->mbs;
4062d61bbb3SSatish Balay     ii  = a->i;
40726e093fcSHong Zhang     z   = zarray;
40826e093fcSHong Zhang   }
4092d61bbb3SSatish Balay 
4102d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
4112d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
4122d61bbb3SSatish Balay     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0;
41398c9bda7SSatish Balay     nonzerorow += (n>0);
4142d61bbb3SSatish Balay     for (j=0; j<n; j++) {
4152d61bbb3SSatish Balay       xb = x + 4*(*idx++);
4162d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
4172d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[4]*x2 + v[8]*x3   + v[12]*x4;
4182d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[5]*x2 + v[9]*x3   + v[13]*x4;
4192d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[6]*x2 + v[10]*x3  + v[14]*x4;
4202d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[7]*x2 + v[11]*x3  + v[15]*x4;
4212d61bbb3SSatish Balay       v += 16;
4222d61bbb3SSatish Balay     }
4237b2bb3b9SHong Zhang     if (usecprow) z = zarray + 4*ridx[i];
4242d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4;
42526e093fcSHong Zhang     if (!usecprow) z += 4;
4262d61bbb3SSatish Balay   }
427d9fead3dSBarry Smith   ierr = VecRestoreArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
42826e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
429dc0b31edSSatish Balay   ierr = PetscLogFlops(32.0*a->nz - 4.0*nonzerorow);CHKERRQ(ierr);
4302d61bbb3SSatish Balay   PetscFunctionReturn(0);
4312d61bbb3SSatish Balay }
4322d61bbb3SSatish Balay 
4334a2ae208SSatish Balay #undef __FUNCT__
4344a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_5"
435dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_5(Mat A,Vec xx,Vec zz)
4362d61bbb3SSatish Balay {
4372d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
438d9fead3dSBarry Smith   PetscScalar       sum1,sum2,sum3,sum4,sum5,x1,x2,x3,x4,x5,*z = 0,*zarray;
439d9fead3dSBarry Smith   const PetscScalar *xb,*x;
440d9fead3dSBarry Smith   const MatScalar   *v;
441dfbe8321SBarry Smith   PetscErrorCode    ierr;
44298c9bda7SSatish Balay   PetscInt          mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL,nonzerorow=0;
44326e093fcSHong Zhang   PetscTruth        usecprow=a->compressedrow.use;
4442d61bbb3SSatish Balay 
445433994e6SBarry Smith   PetscFunctionBegin;
446d9fead3dSBarry Smith   ierr = VecGetArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
44726e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
4482d61bbb3SSatish Balay 
4492d61bbb3SSatish Balay   idx = a->j;
4502d61bbb3SSatish Balay   v   = a->a;
45126e093fcSHong Zhang   if (usecprow){
45226e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
45326e093fcSHong Zhang     ii   = a->compressedrow.i;
4547b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
45526e093fcSHong Zhang   } else {
45626e093fcSHong Zhang     mbs = a->mbs;
4572d61bbb3SSatish Balay     ii  = a->i;
45826e093fcSHong Zhang     z   = zarray;
45926e093fcSHong Zhang   }
4602d61bbb3SSatish Balay 
4612d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
4622d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
4632d61bbb3SSatish Balay     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0;
46498c9bda7SSatish Balay     nonzerorow += (n>0);
4652d61bbb3SSatish Balay     for (j=0; j<n; j++) {
4662d61bbb3SSatish Balay       xb = x + 5*(*idx++);
4672d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4];
4682d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[5]*x2 + v[10]*x3  + v[15]*x4 + v[20]*x5;
4692d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[6]*x2 + v[11]*x3  + v[16]*x4 + v[21]*x5;
4702d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[7]*x2 + v[12]*x3  + v[17]*x4 + v[22]*x5;
4712d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[8]*x2 + v[13]*x3  + v[18]*x4 + v[23]*x5;
4722d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[9]*x2 + v[14]*x3  + v[19]*x4 + v[24]*x5;
4732d61bbb3SSatish Balay       v += 25;
4742d61bbb3SSatish Balay     }
4757b2bb3b9SHong Zhang     if (usecprow) z = zarray + 5*ridx[i];
4762d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5;
47726e093fcSHong Zhang     if (!usecprow) z += 5;
4782d61bbb3SSatish Balay   }
479d9fead3dSBarry Smith   ierr = VecRestoreArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
48026e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
481dc0b31edSSatish Balay   ierr = PetscLogFlops(50.0*a->nz - 5.0*nonzerorow);CHKERRQ(ierr);
4822d61bbb3SSatish Balay   PetscFunctionReturn(0);
4832d61bbb3SSatish Balay }
4842d61bbb3SSatish Balay 
48515091d37SBarry Smith 
4864a2ae208SSatish Balay #undef __FUNCT__
4874a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_6"
488dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_6(Mat A,Vec xx,Vec zz)
48915091d37SBarry Smith {
49015091d37SBarry Smith   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
491d9fead3dSBarry Smith   PetscScalar       *z = 0,sum1,sum2,sum3,sum4,sum5,sum6;
492d9fead3dSBarry Smith   const PetscScalar *x,*xb;
49326e093fcSHong Zhang   PetscScalar       x1,x2,x3,x4,x5,x6,*zarray;
494d9fead3dSBarry Smith   const MatScalar   *v;
495dfbe8321SBarry Smith   PetscErrorCode    ierr;
49698c9bda7SSatish Balay   PetscInt          mbs=a->mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL,nonzerorow=0;
49726e093fcSHong Zhang   PetscTruth        usecprow=a->compressedrow.use;
49815091d37SBarry Smith 
499433994e6SBarry Smith   PetscFunctionBegin;
500d9fead3dSBarry Smith   ierr = VecGetArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
50126e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
50215091d37SBarry Smith 
50315091d37SBarry Smith   idx = a->j;
50415091d37SBarry Smith   v   = a->a;
50526e093fcSHong Zhang   if (usecprow){
50626e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
50726e093fcSHong Zhang     ii   = a->compressedrow.i;
5087b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
50926e093fcSHong Zhang   } else {
51026e093fcSHong Zhang     mbs = a->mbs;
51115091d37SBarry Smith     ii  = a->i;
51226e093fcSHong Zhang     z   = zarray;
51326e093fcSHong Zhang   }
51415091d37SBarry Smith 
51515091d37SBarry Smith   for (i=0; i<mbs; i++) {
51615091d37SBarry Smith     n  = ii[1] - ii[0]; ii++;
51715091d37SBarry Smith     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0;
51898c9bda7SSatish Balay     nonzerorow += (n>0);
51915091d37SBarry Smith     for (j=0; j<n; j++) {
52015091d37SBarry Smith       xb = x + 6*(*idx++);
52115091d37SBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5];
52215091d37SBarry Smith       sum1 += v[0]*x1 + v[6]*x2  + v[12]*x3  + v[18]*x4 + v[24]*x5 + v[30]*x6;
52315091d37SBarry Smith       sum2 += v[1]*x1 + v[7]*x2  + v[13]*x3  + v[19]*x4 + v[25]*x5 + v[31]*x6;
52415091d37SBarry Smith       sum3 += v[2]*x1 + v[8]*x2  + v[14]*x3  + v[20]*x4 + v[26]*x5 + v[32]*x6;
52515091d37SBarry Smith       sum4 += v[3]*x1 + v[9]*x2  + v[15]*x3  + v[21]*x4 + v[27]*x5 + v[33]*x6;
52615091d37SBarry Smith       sum5 += v[4]*x1 + v[10]*x2 + v[16]*x3  + v[22]*x4 + v[28]*x5 + v[34]*x6;
52715091d37SBarry Smith       sum6 += v[5]*x1 + v[11]*x2 + v[17]*x3  + v[23]*x4 + v[29]*x5 + v[35]*x6;
52815091d37SBarry Smith       v += 36;
52915091d37SBarry Smith     }
5307b2bb3b9SHong Zhang     if (usecprow) z = zarray + 6*ridx[i];
53115091d37SBarry Smith     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6;
53226e093fcSHong Zhang     if (!usecprow) z += 6;
53315091d37SBarry Smith   }
53415091d37SBarry Smith 
535d9fead3dSBarry Smith   ierr = VecRestoreArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
53626e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
537dc0b31edSSatish Balay   ierr = PetscLogFlops(72.0*a->nz - 6.0*nonzerorow);CHKERRQ(ierr);
53815091d37SBarry Smith   PetscFunctionReturn(0);
53915091d37SBarry Smith }
5404a2ae208SSatish Balay #undef __FUNCT__
5414a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_7"
542dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_7(Mat A,Vec xx,Vec zz)
5432d61bbb3SSatish Balay {
5442d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
545d9fead3dSBarry Smith   PetscScalar       *z = 0,sum1,sum2,sum3,sum4,sum5,sum6,sum7;
546d9fead3dSBarry Smith   const PetscScalar *x,*xb;
54726e093fcSHong Zhang   PetscScalar       x1,x2,x3,x4,x5,x6,x7,*zarray;
548d9fead3dSBarry Smith   const MatScalar   *v;
549dfbe8321SBarry Smith   PetscErrorCode    ierr;
55098c9bda7SSatish Balay   PetscInt          mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL,nonzerorow=0;
55126e093fcSHong Zhang   PetscTruth        usecprow=a->compressedrow.use;
5522d61bbb3SSatish Balay 
553433994e6SBarry Smith   PetscFunctionBegin;
554d9fead3dSBarry Smith   ierr = VecGetArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
55526e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
5562d61bbb3SSatish Balay 
5572d61bbb3SSatish Balay   idx = a->j;
5582d61bbb3SSatish Balay   v   = a->a;
55926e093fcSHong Zhang   if (usecprow){
56026e093fcSHong Zhang     mbs    = a->compressedrow.nrows;
56126e093fcSHong Zhang     ii     = a->compressedrow.i;
5627b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
56326e093fcSHong Zhang   } else {
56426e093fcSHong Zhang     mbs = a->mbs;
5652d61bbb3SSatish Balay     ii  = a->i;
56626e093fcSHong Zhang     z   = zarray;
56726e093fcSHong Zhang   }
5682d61bbb3SSatish Balay 
5692d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
5702d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
5712d61bbb3SSatish Balay     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0; sum7 = 0.0;
57298c9bda7SSatish Balay     nonzerorow += (n>0);
5732d61bbb3SSatish Balay     for (j=0; j<n; j++) {
5742d61bbb3SSatish Balay       xb = x + 7*(*idx++);
5752d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6];
5762d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[7]*x2  + v[14]*x3  + v[21]*x4 + v[28]*x5 + v[35]*x6 + v[42]*x7;
5772d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[8]*x2  + v[15]*x3  + v[22]*x4 + v[29]*x5 + v[36]*x6 + v[43]*x7;
5782d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[9]*x2  + v[16]*x3  + v[23]*x4 + v[30]*x5 + v[37]*x6 + v[44]*x7;
5792d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[10]*x2 + v[17]*x3  + v[24]*x4 + v[31]*x5 + v[38]*x6 + v[45]*x7;
5802d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[11]*x2 + v[18]*x3  + v[25]*x4 + v[32]*x5 + v[39]*x6 + v[46]*x7;
5812d61bbb3SSatish Balay       sum6 += v[5]*x1 + v[12]*x2 + v[19]*x3  + v[26]*x4 + v[33]*x5 + v[40]*x6 + v[47]*x7;
5822d61bbb3SSatish Balay       sum7 += v[6]*x1 + v[13]*x2 + v[20]*x3  + v[27]*x4 + v[34]*x5 + v[41]*x6 + v[48]*x7;
5832d61bbb3SSatish Balay       v += 49;
5842d61bbb3SSatish Balay     }
5857b2bb3b9SHong Zhang     if (usecprow) z = zarray + 7*ridx[i];
5862d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
58726e093fcSHong Zhang     if (!usecprow) z += 7;
5882d61bbb3SSatish Balay   }
5892d61bbb3SSatish Balay 
590d9fead3dSBarry Smith   ierr = VecRestoreArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
59126e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
592dc0b31edSSatish Balay   ierr = PetscLogFlops(98.0*a->nz - 7.0*nonzerorow);CHKERRQ(ierr);
5932d61bbb3SSatish Balay   PetscFunctionReturn(0);
5942d61bbb3SSatish Balay }
5952d61bbb3SSatish Balay 
5963f1db9ecSBarry Smith /*
5973f1db9ecSBarry Smith     This will not work with MatScalar == float because it calls the BLAS
5983f1db9ecSBarry Smith */
5994a2ae208SSatish Balay #undef __FUNCT__
6004a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_N"
601dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_N(Mat A,Vec xx,Vec zz)
6022d61bbb3SSatish Balay {
6032d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
604dcf5cc72SBarry Smith   PetscScalar    *x,*z = 0,*xb,*work,*workt,*zarray;
6053f1db9ecSBarry Smith   MatScalar      *v;
606dfbe8321SBarry Smith   PetscErrorCode ierr;
607d0f46423SBarry Smith   PetscInt       mbs=a->mbs,i,*idx,*ii,bs=A->rmap->bs,j,n,bs2=a->bs2;
60898c9bda7SSatish Balay   PetscInt       ncols,k,*ridx=PETSC_NULL,nonzerorow=0;
60926e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
6102d61bbb3SSatish Balay 
6112d61bbb3SSatish Balay   PetscFunctionBegin;
6121ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
61326e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
6142d61bbb3SSatish Balay 
6152d61bbb3SSatish Balay   idx = a->j;
6162d61bbb3SSatish Balay   v   = a->a;
61726e093fcSHong Zhang   if (usecprow){
61826e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
61926e093fcSHong Zhang     ii   = a->compressedrow.i;
6207b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
62126e093fcSHong Zhang   } else {
62226e093fcSHong Zhang     mbs = a->mbs;
6232d61bbb3SSatish Balay     ii  = a->i;
62426e093fcSHong Zhang     z   = zarray;
62526e093fcSHong Zhang   }
626218c64b6SSatish Balay 
6272d61bbb3SSatish Balay   if (!a->mult_work) {
628d0f46423SBarry Smith     k    = PetscMax(A->rmap->n,A->cmap->n);
62987828ca2SBarry Smith     ierr = PetscMalloc((k+1)*sizeof(PetscScalar),&a->mult_work);CHKERRQ(ierr);
6302d61bbb3SSatish Balay   }
6312d61bbb3SSatish Balay   work = a->mult_work;
6322d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
6332d61bbb3SSatish Balay     n     = ii[1] - ii[0]; ii++;
6342d61bbb3SSatish Balay     ncols = n*bs;
6352d61bbb3SSatish Balay     workt = work;
63698c9bda7SSatish Balay     nonzerorow += (n>0);
6372d61bbb3SSatish Balay     for (j=0; j<n; j++) {
6382d61bbb3SSatish Balay       xb = x + bs*(*idx++);
6392d61bbb3SSatish Balay       for (k=0; k<bs; k++) workt[k] = xb[k];
6402d61bbb3SSatish Balay       workt += bs;
6412d61bbb3SSatish Balay     }
6427b2bb3b9SHong Zhang     if (usecprow) z = zarray + bs*ridx[i];
643737d121aSSatish Balay     Kernel_w_gets_Ar_times_v(bs,ncols,work,v,z);
64471044d3cSBarry Smith     /* BLASgemv_("N",&bs,&ncols,&_DOne,v,&bs,work,&_One,&_DZero,z,&_One); */
6452d61bbb3SSatish Balay     v += n*bs2;
64626e093fcSHong Zhang     if (!usecprow) z += bs;
6472d61bbb3SSatish Balay   }
6481ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
64926e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
650dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*a->nz*bs2 - bs*nonzerorow);CHKERRQ(ierr);
6512d61bbb3SSatish Balay   PetscFunctionReturn(0);
6522d61bbb3SSatish Balay }
6532d61bbb3SSatish Balay 
6546a65c766SHong Zhang extern PetscErrorCode VecCopy_Seq(Vec,Vec);
6554a2ae208SSatish Balay #undef __FUNCT__
6564a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_1"
657dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_1(Mat A,Vec xx,Vec yy,Vec zz)
6582d61bbb3SSatish Balay {
6592d61bbb3SSatish Balay   Mat_SeqBAIJ        *a = (Mat_SeqBAIJ*)A->data;
660122f12eaSBarry Smith   const PetscScalar  *x;
661122f12eaSBarry Smith   PetscScalar        *y,*z,sum;
662122f12eaSBarry Smith   const MatScalar    *v;
663dfbe8321SBarry Smith   PetscErrorCode     ierr;
664122f12eaSBarry Smith   PetscInt           mbs=a->mbs,i,n,*ridx=PETSC_NULL,nonzerorow=0;
665122f12eaSBarry Smith   const PetscInt     *idx,*ii;
66626e093fcSHong Zhang   PetscTruth         usecprow=a->compressedrow.use;
6672d61bbb3SSatish Balay 
6682d61bbb3SSatish Balay   PetscFunctionBegin;
669122f12eaSBarry Smith   ierr = VecGetArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
670122f12eaSBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
6712e8a6d31SBarry Smith   if (zz != yy) {
672122f12eaSBarry Smith     ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
6732e8a6d31SBarry Smith   } else {
674122f12eaSBarry Smith     z = y;
6752e8a6d31SBarry Smith   }
6762d61bbb3SSatish Balay 
6772d61bbb3SSatish Balay   idx = a->j;
6782d61bbb3SSatish Balay   v   = a->a;
67926e093fcSHong Zhang   if (usecprow){
6804eb6d288SHong Zhang     if (zz != yy){
681122f12eaSBarry Smith       ierr = PetscMemcpy(z,y,mbs*sizeof(PetscScalar));CHKERRQ(ierr);
6824eb6d288SHong Zhang     }
68326e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
68426e093fcSHong Zhang     ii   = a->compressedrow.i;
6857b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
68626e093fcSHong Zhang   } else {
6872d61bbb3SSatish Balay     ii  = a->i;
68826e093fcSHong Zhang   }
6892d61bbb3SSatish Balay 
6902d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
691122f12eaSBarry Smith     n    = ii[1] - ii[0];
692122f12eaSBarry Smith     ii++;
69326e093fcSHong Zhang     if (!usecprow){
694122f12eaSBarry Smith       nonzerorow += (n>0);
695122f12eaSBarry Smith       sum = y[i];
696122f12eaSBarry Smith     } else {
697122f12eaSBarry Smith       sum = y[ridx[i]];
698122f12eaSBarry Smith     }
699122f12eaSBarry Smith     PetscSparseDensePlusDot(sum,x,v,idx,n);
700122f12eaSBarry Smith     v += n;
701122f12eaSBarry Smith     idx += n;
702122f12eaSBarry Smith     if (usecprow){
703122f12eaSBarry Smith       z[ridx[i]] = sum;
704122f12eaSBarry Smith     } else {
705122f12eaSBarry Smith       z[i] = sum;
70626e093fcSHong Zhang     }
7072d61bbb3SSatish Balay   }
708122f12eaSBarry Smith   ierr = VecRestoreArray(xx,(PetscScalar**)&x);CHKERRQ(ierr);
709122f12eaSBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
7102e8a6d31SBarry Smith   if (zz != yy) {
711122f12eaSBarry Smith     ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
7122e8a6d31SBarry Smith   }
713122f12eaSBarry Smith   ierr = PetscLogFlops(2.0*a->nz - nonzerorow);CHKERRQ(ierr);
7142d61bbb3SSatish Balay   PetscFunctionReturn(0);
7152d61bbb3SSatish Balay }
7162d61bbb3SSatish Balay 
7174a2ae208SSatish Balay #undef __FUNCT__
7184a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_2"
719dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_2(Mat A,Vec xx,Vec yy,Vec zz)
7202d61bbb3SSatish Balay {
7212d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
722dcf5cc72SBarry Smith   PetscScalar    *x,*y = 0,*z = 0,*xb,sum1,sum2;
72326e093fcSHong Zhang   PetscScalar    x1,x2,*yarray,*zarray;
7243f1db9ecSBarry Smith   MatScalar      *v;
725dfbe8321SBarry Smith   PetscErrorCode ierr;
7264eb6d288SHong Zhang   PetscInt       mbs=a->mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL;
72726e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
7282d61bbb3SSatish Balay 
7292d61bbb3SSatish Balay   PetscFunctionBegin;
7301ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
73126e093fcSHong Zhang   ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr);
7322e8a6d31SBarry Smith   if (zz != yy) {
73326e093fcSHong Zhang     ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
7342e8a6d31SBarry Smith   } else {
73526e093fcSHong Zhang     zarray = yarray;
7362e8a6d31SBarry Smith   }
7372d61bbb3SSatish Balay 
7382d61bbb3SSatish Balay   idx = a->j;
7392d61bbb3SSatish Balay   v   = a->a;
74026e093fcSHong Zhang   if (usecprow){
7414eb6d288SHong Zhang     if (zz != yy){
7424eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,2*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
7434eb6d288SHong Zhang     }
74426e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
74526e093fcSHong Zhang     ii   = a->compressedrow.i;
7467b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
7474eb6d288SHong Zhang     if (zz != yy){
7484eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,a->mbs*sizeof(PetscScalar));CHKERRQ(ierr);
7494eb6d288SHong Zhang     }
75026e093fcSHong Zhang   } else {
7512d61bbb3SSatish Balay     ii  = a->i;
75226e093fcSHong Zhang     y   = yarray;
75326e093fcSHong Zhang     z   = zarray;
75426e093fcSHong Zhang   }
7552d61bbb3SSatish Balay 
7562d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
7572d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
75826e093fcSHong Zhang     if (usecprow){
7597b2bb3b9SHong Zhang       z = zarray + 2*ridx[i];
7607b2bb3b9SHong Zhang       y = yarray + 2*ridx[i];
76126e093fcSHong Zhang     }
7622d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1];
7632d61bbb3SSatish Balay     for (j=0; j<n; j++) {
7642d61bbb3SSatish Balay       xb = x + 2*(*idx++); x1 = xb[0]; x2 = xb[1];
7652d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[2]*x2;
7662d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[3]*x2;
7672d61bbb3SSatish Balay       v += 4;
7682d61bbb3SSatish Balay     }
7692d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2;
77026e093fcSHong Zhang     if (!usecprow){
7712d61bbb3SSatish Balay       z += 2; y += 2;
7722d61bbb3SSatish Balay     }
77326e093fcSHong Zhang   }
7741ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
77526e093fcSHong Zhang   ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr);
7762e8a6d31SBarry Smith   if (zz != yy) {
77726e093fcSHong Zhang     ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
7782e8a6d31SBarry Smith   }
779dc0b31edSSatish Balay   ierr = PetscLogFlops(4.0*a->nz);CHKERRQ(ierr);
7802d61bbb3SSatish Balay   PetscFunctionReturn(0);
7812d61bbb3SSatish Balay }
7822d61bbb3SSatish Balay 
7834a2ae208SSatish Balay #undef __FUNCT__
7844a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_3"
785dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_3(Mat A,Vec xx,Vec yy,Vec zz)
7862d61bbb3SSatish Balay {
7872d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
788dcf5cc72SBarry Smith   PetscScalar    *x,*y = 0,*z = 0,*xb,sum1,sum2,sum3,x1,x2,x3,*yarray,*zarray;
7893f1db9ecSBarry Smith   MatScalar      *v;
790dfbe8321SBarry Smith   PetscErrorCode ierr;
7914eb6d288SHong Zhang   PetscInt       mbs=a->mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL;
79226e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
7932d61bbb3SSatish Balay 
7942d61bbb3SSatish Balay   PetscFunctionBegin;
7951ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
79626e093fcSHong Zhang   ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr);
7972e8a6d31SBarry Smith   if (zz != yy) {
79826e093fcSHong Zhang     ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
7992e8a6d31SBarry Smith   } else {
80026e093fcSHong Zhang     zarray = yarray;
8012e8a6d31SBarry Smith   }
8022d61bbb3SSatish Balay 
8032d61bbb3SSatish Balay   idx = a->j;
8042d61bbb3SSatish Balay   v   = a->a;
80526e093fcSHong Zhang   if (usecprow){
8064eb6d288SHong Zhang     if (zz != yy){
8074eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,3*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
8084eb6d288SHong Zhang     }
80926e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
81026e093fcSHong Zhang     ii   = a->compressedrow.i;
8117b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
81226e093fcSHong Zhang   } else {
8132d61bbb3SSatish Balay     ii  = a->i;
81426e093fcSHong Zhang     y   = yarray;
81526e093fcSHong Zhang     z   = zarray;
81626e093fcSHong Zhang   }
8172d61bbb3SSatish Balay 
8182d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
8192d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
82026e093fcSHong Zhang     if (usecprow){
8217b2bb3b9SHong Zhang       z = zarray + 3*ridx[i];
8227b2bb3b9SHong Zhang       y = yarray + 3*ridx[i];
82326e093fcSHong Zhang     }
8242d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2];
8252d61bbb3SSatish Balay     for (j=0; j<n; j++) {
8262d61bbb3SSatish Balay       xb = x + 3*(*idx++); x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
8272d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[3]*x2 + v[6]*x3;
8282d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[4]*x2 + v[7]*x3;
8292d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[5]*x2 + v[8]*x3;
8302d61bbb3SSatish Balay       v += 9;
8312d61bbb3SSatish Balay     }
8322d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3;
83326e093fcSHong Zhang     if (!usecprow){
8342d61bbb3SSatish Balay       z += 3; y += 3;
8352d61bbb3SSatish Balay     }
83626e093fcSHong Zhang   }
8371ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
83826e093fcSHong Zhang   ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr);
8392e8a6d31SBarry Smith   if (zz != yy) {
84026e093fcSHong Zhang     ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
8412e8a6d31SBarry Smith   }
842dc0b31edSSatish Balay   ierr = PetscLogFlops(18.0*a->nz);CHKERRQ(ierr);
8432d61bbb3SSatish Balay   PetscFunctionReturn(0);
8442d61bbb3SSatish Balay }
8452d61bbb3SSatish Balay 
8464a2ae208SSatish Balay #undef __FUNCT__
8474a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_4"
848dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_4(Mat A,Vec xx,Vec yy,Vec zz)
8492d61bbb3SSatish Balay {
8502d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
851dcf5cc72SBarry Smith   PetscScalar    *x,*y = 0,*z = 0,*xb,sum1,sum2,sum3,sum4,x1,x2,x3,x4,*yarray,*zarray;
8523f1db9ecSBarry Smith   MatScalar      *v;
853dfbe8321SBarry Smith   PetscErrorCode ierr;
8544eb6d288SHong Zhang   PetscInt       mbs=a->mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL;
85526e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
8562d61bbb3SSatish Balay 
8572d61bbb3SSatish Balay   PetscFunctionBegin;
8581ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
85926e093fcSHong Zhang   ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr);
8602e8a6d31SBarry Smith   if (zz != yy) {
86126e093fcSHong Zhang     ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
8622e8a6d31SBarry Smith   } else {
86326e093fcSHong Zhang     zarray = yarray;
8642e8a6d31SBarry Smith   }
8652d61bbb3SSatish Balay 
8662d61bbb3SSatish Balay   idx   = a->j;
8672d61bbb3SSatish Balay   v     = a->a;
86826e093fcSHong Zhang   if (usecprow){
8694eb6d288SHong Zhang     if (zz != yy){
8704eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,4*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
8714eb6d288SHong Zhang     }
87226e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
87326e093fcSHong Zhang     ii   = a->compressedrow.i;
8747b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
87526e093fcSHong Zhang   } else {
8762d61bbb3SSatish Balay     ii  = a->i;
87726e093fcSHong Zhang     y   = yarray;
87826e093fcSHong Zhang     z   = zarray;
87926e093fcSHong Zhang   }
8802d61bbb3SSatish Balay 
8812d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
8822d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
88326e093fcSHong Zhang     if (usecprow){
8847b2bb3b9SHong Zhang       z = zarray + 4*ridx[i];
8857b2bb3b9SHong Zhang       y = yarray + 4*ridx[i];
88626e093fcSHong Zhang     }
8872d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3];
8882d61bbb3SSatish Balay     for (j=0; j<n; j++) {
8892d61bbb3SSatish Balay       xb = x + 4*(*idx++);
8902d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
8912d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[4]*x2 + v[8]*x3   + v[12]*x4;
8922d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[5]*x2 + v[9]*x3   + v[13]*x4;
8932d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[6]*x2 + v[10]*x3  + v[14]*x4;
8942d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[7]*x2 + v[11]*x3  + v[15]*x4;
8952d61bbb3SSatish Balay       v += 16;
8962d61bbb3SSatish Balay     }
8972d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4;
89826e093fcSHong Zhang     if (!usecprow){
8992d61bbb3SSatish Balay       z += 4; y += 4;
9002d61bbb3SSatish Balay     }
90126e093fcSHong Zhang   }
9021ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
90326e093fcSHong Zhang   ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr);
9042e8a6d31SBarry Smith   if (zz != yy) {
90526e093fcSHong Zhang     ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
9062e8a6d31SBarry Smith   }
907dc0b31edSSatish Balay   ierr = PetscLogFlops(32.0*a->nz);CHKERRQ(ierr);
9082d61bbb3SSatish Balay   PetscFunctionReturn(0);
9092d61bbb3SSatish Balay }
9102d61bbb3SSatish Balay 
9114a2ae208SSatish Balay #undef __FUNCT__
9124a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_5"
913dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_5(Mat A,Vec xx,Vec yy,Vec zz)
9142d61bbb3SSatish Balay {
9152d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
916dcf5cc72SBarry Smith   PetscScalar    *x,*y = 0,*z = 0,*xb,sum1,sum2,sum3,sum4,sum5,x1,x2,x3,x4,x5;
91726e093fcSHong Zhang   PetscScalar    *yarray,*zarray;
9183f1db9ecSBarry Smith   MatScalar      *v;
919dfbe8321SBarry Smith   PetscErrorCode ierr;
9204eb6d288SHong Zhang   PetscInt       mbs=a->mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL;
92126e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
9222d61bbb3SSatish Balay 
9232d61bbb3SSatish Balay   PetscFunctionBegin;
9241ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
92526e093fcSHong Zhang   ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr);
9262e8a6d31SBarry Smith   if (zz != yy) {
92726e093fcSHong Zhang     ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
9282e8a6d31SBarry Smith   } else {
92926e093fcSHong Zhang     zarray = yarray;
9302e8a6d31SBarry Smith   }
9312d61bbb3SSatish Balay 
9322d61bbb3SSatish Balay   idx = a->j;
9332d61bbb3SSatish Balay   v   = a->a;
93426e093fcSHong Zhang   if (usecprow){
9354eb6d288SHong Zhang     if (zz != yy){
9364eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,5*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
9374eb6d288SHong Zhang     }
93826e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
93926e093fcSHong Zhang     ii   = a->compressedrow.i;
9407b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
94126e093fcSHong Zhang   } else {
9422d61bbb3SSatish Balay     ii  = a->i;
94326e093fcSHong Zhang     y   = yarray;
94426e093fcSHong Zhang     z   = zarray;
94526e093fcSHong Zhang   }
9462d61bbb3SSatish Balay 
9472d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
9482d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
94926e093fcSHong Zhang     if (usecprow){
9507b2bb3b9SHong Zhang       z = zarray + 5*ridx[i];
9517b2bb3b9SHong Zhang       y = yarray + 5*ridx[i];
95226e093fcSHong Zhang     }
9532d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4];
9542d61bbb3SSatish Balay     for (j=0; j<n; j++) {
9552d61bbb3SSatish Balay       xb = x + 5*(*idx++);
9562d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4];
9572d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[5]*x2 + v[10]*x3  + v[15]*x4 + v[20]*x5;
9582d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[6]*x2 + v[11]*x3  + v[16]*x4 + v[21]*x5;
9592d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[7]*x2 + v[12]*x3  + v[17]*x4 + v[22]*x5;
9602d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[8]*x2 + v[13]*x3  + v[18]*x4 + v[23]*x5;
9612d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[9]*x2 + v[14]*x3  + v[19]*x4 + v[24]*x5;
9622d61bbb3SSatish Balay       v += 25;
9632d61bbb3SSatish Balay     }
9642d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5;
96526e093fcSHong Zhang     if (!usecprow){
9662d61bbb3SSatish Balay       z += 5; y += 5;
9672d61bbb3SSatish Balay     }
96826e093fcSHong Zhang   }
9691ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
97026e093fcSHong Zhang   ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr);
9712e8a6d31SBarry Smith   if (zz != yy) {
97226e093fcSHong Zhang     ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
9732e8a6d31SBarry Smith   }
974dc0b31edSSatish Balay   ierr = PetscLogFlops(50.0*a->nz);CHKERRQ(ierr);
9752d61bbb3SSatish Balay   PetscFunctionReturn(0);
9762d61bbb3SSatish Balay }
9774a2ae208SSatish Balay #undef __FUNCT__
9784a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_6"
979dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_6(Mat A,Vec xx,Vec yy,Vec zz)
98015091d37SBarry Smith {
98115091d37SBarry Smith   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
982dcf5cc72SBarry Smith   PetscScalar    *x,*y = 0,*z = 0,*xb,sum1,sum2,sum3,sum4,sum5,sum6;
98326e093fcSHong Zhang   PetscScalar    x1,x2,x3,x4,x5,x6,*yarray,*zarray;
98415091d37SBarry Smith   MatScalar      *v;
985dfbe8321SBarry Smith   PetscErrorCode ierr;
9864eb6d288SHong Zhang   PetscInt       mbs=a->mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL;
98726e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
98815091d37SBarry Smith 
98915091d37SBarry Smith   PetscFunctionBegin;
9901ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
99126e093fcSHong Zhang   ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr);
99215091d37SBarry Smith   if (zz != yy) {
99326e093fcSHong Zhang     ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
99415091d37SBarry Smith   } else {
99526e093fcSHong Zhang     zarray = yarray;
99615091d37SBarry Smith   }
99715091d37SBarry Smith 
99815091d37SBarry Smith   idx = a->j;
99915091d37SBarry Smith   v   = a->a;
100026e093fcSHong Zhang   if (usecprow){
10014eb6d288SHong Zhang     if (zz != yy){
10024eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,6*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
10034eb6d288SHong Zhang     }
100426e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
100526e093fcSHong Zhang     ii   = a->compressedrow.i;
10067b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
100726e093fcSHong Zhang   } else {
100815091d37SBarry Smith     ii  = a->i;
100926e093fcSHong Zhang     y   = yarray;
101026e093fcSHong Zhang     z   = zarray;
101126e093fcSHong Zhang   }
101215091d37SBarry Smith 
101315091d37SBarry Smith   for (i=0; i<mbs; i++) {
101415091d37SBarry Smith     n  = ii[1] - ii[0]; ii++;
101526e093fcSHong Zhang     if (usecprow){
10167b2bb3b9SHong Zhang       z = zarray + 6*ridx[i];
10177b2bb3b9SHong Zhang       y = yarray + 6*ridx[i];
101826e093fcSHong Zhang     }
101915091d37SBarry Smith     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4]; sum6 = y[5];
102015091d37SBarry Smith     for (j=0; j<n; j++) {
10213b95cb0eSSatish Balay       xb = x + 6*(*idx++);
102215091d37SBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5];
102315091d37SBarry Smith       sum1 += v[0]*x1 + v[6]*x2  + v[12]*x3  + v[18]*x4 + v[24]*x5 + v[30]*x6;
102415091d37SBarry Smith       sum2 += v[1]*x1 + v[7]*x2  + v[13]*x3  + v[19]*x4 + v[25]*x5 + v[31]*x6;
102515091d37SBarry Smith       sum3 += v[2]*x1 + v[8]*x2  + v[14]*x3  + v[20]*x4 + v[26]*x5 + v[32]*x6;
102615091d37SBarry Smith       sum4 += v[3]*x1 + v[9]*x2  + v[15]*x3  + v[21]*x4 + v[27]*x5 + v[33]*x6;
102715091d37SBarry Smith       sum5 += v[4]*x1 + v[10]*x2 + v[16]*x3  + v[22]*x4 + v[28]*x5 + v[34]*x6;
102815091d37SBarry Smith       sum6 += v[5]*x1 + v[11]*x2 + v[17]*x3  + v[23]*x4 + v[29]*x5 + v[35]*x6;
102915091d37SBarry Smith       v += 36;
103015091d37SBarry Smith     }
103115091d37SBarry Smith     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6;
103226e093fcSHong Zhang     if (!usecprow){
103315091d37SBarry Smith       z += 6; y += 6;
103415091d37SBarry Smith     }
103526e093fcSHong Zhang   }
10361ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
103726e093fcSHong Zhang   ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr);
103815091d37SBarry Smith   if (zz != yy) {
103926e093fcSHong Zhang     ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
104015091d37SBarry Smith   }
1041dc0b31edSSatish Balay   ierr = PetscLogFlops(72.0*a->nz);CHKERRQ(ierr);
104215091d37SBarry Smith   PetscFunctionReturn(0);
104315091d37SBarry Smith }
10442d61bbb3SSatish Balay 
10454a2ae208SSatish Balay #undef __FUNCT__
10464a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_7"
1047dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_7(Mat A,Vec xx,Vec yy,Vec zz)
10482d61bbb3SSatish Balay {
10492d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
1050dcf5cc72SBarry Smith   PetscScalar    *x,*y = 0,*z = 0,*xb,sum1,sum2,sum3,sum4,sum5,sum6,sum7;
105126e093fcSHong Zhang   PetscScalar    x1,x2,x3,x4,x5,x6,x7,*yarray,*zarray;
10523f1db9ecSBarry Smith   MatScalar      *v;
1053dfbe8321SBarry Smith   PetscErrorCode ierr;
10547b2bb3b9SHong Zhang   PetscInt       mbs=a->mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL;
105526e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
10562d61bbb3SSatish Balay 
10572d61bbb3SSatish Balay   PetscFunctionBegin;
10581ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
105926e093fcSHong Zhang   ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr);
10602e8a6d31SBarry Smith   if (zz != yy) {
106126e093fcSHong Zhang     ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
10622e8a6d31SBarry Smith   } else {
106326e093fcSHong Zhang     zarray = yarray;
10642e8a6d31SBarry Smith   }
10652d61bbb3SSatish Balay 
10662d61bbb3SSatish Balay   idx = a->j;
10672d61bbb3SSatish Balay   v   = a->a;
106826e093fcSHong Zhang   if (usecprow){
10694eb6d288SHong Zhang     if (zz != yy){
10704eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,7*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
10714eb6d288SHong Zhang     }
107226e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
107326e093fcSHong Zhang     ii   = a->compressedrow.i;
10747b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
107526e093fcSHong Zhang   } else {
10762d61bbb3SSatish Balay     ii  = a->i;
107726e093fcSHong Zhang     y   = yarray;
107826e093fcSHong Zhang     z   = zarray;
107926e093fcSHong Zhang   }
10802d61bbb3SSatish Balay 
10812d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
10822d61bbb3SSatish Balay     n  = ii[1] - ii[0]; ii++;
108326e093fcSHong Zhang     if (usecprow){
10847b2bb3b9SHong Zhang       z = zarray + 7*ridx[i];
10857b2bb3b9SHong Zhang       y = yarray + 7*ridx[i];
108626e093fcSHong Zhang     }
10872d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4]; sum6 = y[5]; sum7 = y[6];
10882d61bbb3SSatish Balay     for (j=0; j<n; j++) {
10892d61bbb3SSatish Balay       xb = x + 7*(*idx++);
10902d61bbb3SSatish Balay       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6];
10912d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[7]*x2  + v[14]*x3  + v[21]*x4 + v[28]*x5 + v[35]*x6 + v[42]*x7;
10922d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[8]*x2  + v[15]*x3  + v[22]*x4 + v[29]*x5 + v[36]*x6 + v[43]*x7;
10932d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[9]*x2  + v[16]*x3  + v[23]*x4 + v[30]*x5 + v[37]*x6 + v[44]*x7;
10942d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[10]*x2 + v[17]*x3  + v[24]*x4 + v[31]*x5 + v[38]*x6 + v[45]*x7;
10952d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[11]*x2 + v[18]*x3  + v[25]*x4 + v[32]*x5 + v[39]*x6 + v[46]*x7;
10962d61bbb3SSatish Balay       sum6 += v[5]*x1 + v[12]*x2 + v[19]*x3  + v[26]*x4 + v[33]*x5 + v[40]*x6 + v[47]*x7;
10972d61bbb3SSatish Balay       sum7 += v[6]*x1 + v[13]*x2 + v[20]*x3  + v[27]*x4 + v[34]*x5 + v[41]*x6 + v[48]*x7;
10982d61bbb3SSatish Balay       v += 49;
10992d61bbb3SSatish Balay     }
11002d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
110126e093fcSHong Zhang     if (!usecprow){
11022d61bbb3SSatish Balay       z += 7; y += 7;
11032d61bbb3SSatish Balay     }
110426e093fcSHong Zhang   }
11051ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
110626e093fcSHong Zhang   ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr);
11072e8a6d31SBarry Smith   if (zz != yy) {
110826e093fcSHong Zhang     ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
11092e8a6d31SBarry Smith   }
1110dc0b31edSSatish Balay   ierr = PetscLogFlops(98.0*a->nz);CHKERRQ(ierr);
11112d61bbb3SSatish Balay   PetscFunctionReturn(0);
11122d61bbb3SSatish Balay }
1113218c64b6SSatish Balay 
11144a2ae208SSatish Balay #undef __FUNCT__
11154a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_N"
1116dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_N(Mat A,Vec xx,Vec yy,Vec zz)
11172d61bbb3SSatish Balay {
11182d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
1119bba805e6SBarry Smith   PetscScalar    *x,*z = 0,*xb,*work,*workt,*zarray;
11203f1db9ecSBarry Smith   MatScalar      *v;
11216849ba73SBarry Smith   PetscErrorCode ierr;
1122d0f46423SBarry Smith   PetscInt       mbs,i,*idx,*ii,bs=A->rmap->bs,j,n,bs2=a->bs2;
11237b2bb3b9SHong Zhang   PetscInt       ncols,k,*ridx=PETSC_NULL;
112426e093fcSHong Zhang   PetscTruth     usecprow=a->compressedrow.use;
1125218c64b6SSatish Balay 
11262d61bbb3SSatish Balay   PetscFunctionBegin;
11276a65c766SHong Zhang   ierr = VecCopy_Seq(yy,zz);CHKERRQ(ierr);
11281ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
112926e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
11302d61bbb3SSatish Balay 
11312d61bbb3SSatish Balay   idx = a->j;
11322d61bbb3SSatish Balay   v   = a->a;
113326e093fcSHong Zhang   if (usecprow){
113426e093fcSHong Zhang     mbs    = a->compressedrow.nrows;
113526e093fcSHong Zhang     ii     = a->compressedrow.i;
11367b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
113726e093fcSHong Zhang   } else {
113826e093fcSHong Zhang     mbs = a->mbs;
11392d61bbb3SSatish Balay     ii  = a->i;
114026e093fcSHong Zhang     z   = zarray;
114126e093fcSHong Zhang   }
11422d61bbb3SSatish Balay 
11432d61bbb3SSatish Balay   if (!a->mult_work) {
1144d0f46423SBarry Smith     k    = PetscMax(A->rmap->n,A->cmap->n);
114587828ca2SBarry Smith     ierr = PetscMalloc((k+1)*sizeof(PetscScalar),&a->mult_work);CHKERRQ(ierr);
11462d61bbb3SSatish Balay   }
11472d61bbb3SSatish Balay   work = a->mult_work;
11482d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
11492d61bbb3SSatish Balay     n     = ii[1] - ii[0]; ii++;
11502d61bbb3SSatish Balay     ncols = n*bs;
11512d61bbb3SSatish Balay     workt = work;
11522d61bbb3SSatish Balay     for (j=0; j<n; j++) {
11532d61bbb3SSatish Balay       xb = x + bs*(*idx++);
11542d61bbb3SSatish Balay       for (k=0; k<bs; k++) workt[k] = xb[k];
11552d61bbb3SSatish Balay       workt += bs;
11562d61bbb3SSatish Balay     }
11577b2bb3b9SHong Zhang     if (usecprow) z = zarray + bs*ridx[i];
11583f1db9ecSBarry Smith     Kernel_w_gets_w_plus_Ar_times_v(bs,ncols,work,v,z);
115971044d3cSBarry Smith     /* BLASgemv_("N",&bs,&ncols,&_DOne,v,&bs,work,&_One,&_DOne,z,&_One); */
11602d61bbb3SSatish Balay     v += n*bs2;
116126e093fcSHong Zhang     if (!usecprow){
11622d61bbb3SSatish Balay       z += bs;
11632d61bbb3SSatish Balay     }
116426e093fcSHong Zhang   }
11651ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
116626e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
1167dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*a->nz*bs2);CHKERRQ(ierr);
11682d61bbb3SSatish Balay   PetscFunctionReturn(0);
11692d61bbb3SSatish Balay }
11702d61bbb3SSatish Balay 
11714a2ae208SSatish Balay #undef __FUNCT__
1172547795f9SHong Zhang #define __FUNCT__ "MatMultHermitianTranspose_SeqBAIJ"
1173547795f9SHong Zhang PetscErrorCode MatMultHermitianTranspose_SeqBAIJ(Mat A,Vec xx,Vec zz)
1174547795f9SHong Zhang {
1175547795f9SHong Zhang   PetscScalar    zero = 0.0;
1176547795f9SHong Zhang   PetscErrorCode ierr;
1177547795f9SHong Zhang 
1178547795f9SHong Zhang   PetscFunctionBegin;
1179547795f9SHong Zhang   ierr = VecSet(zz,zero);CHKERRQ(ierr);
1180547795f9SHong Zhang   ierr = MatMultHermitianTransposeAdd_SeqBAIJ(A,xx,zz,zz);CHKERRQ(ierr);
1181547795f9SHong Zhang   PetscFunctionReturn(0);
1182547795f9SHong Zhang }
1183547795f9SHong Zhang 
1184547795f9SHong Zhang #undef __FUNCT__
11854a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_SeqBAIJ"
1186dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_SeqBAIJ(Mat A,Vec xx,Vec zz)
11872d61bbb3SSatish Balay {
11883447b6efSHong Zhang   PetscScalar    zero = 0.0;
11896849ba73SBarry Smith   PetscErrorCode ierr;
11902d61bbb3SSatish Balay 
11912d61bbb3SSatish Balay   PetscFunctionBegin;
11922dcb1b2aSMatthew Knepley   ierr = VecSet(zz,zero);CHKERRQ(ierr);
11933447b6efSHong Zhang   ierr = MatMultTransposeAdd_SeqBAIJ(A,xx,zz,zz);CHKERRQ(ierr);
11942d61bbb3SSatish Balay   PetscFunctionReturn(0);
11952d61bbb3SSatish Balay }
11962d61bbb3SSatish Balay 
11974a2ae208SSatish Balay #undef __FUNCT__
1198547795f9SHong Zhang #define __FUNCT__ "MatMultHermitianTransposeAdd_SeqBAIJ"
1199547795f9SHong Zhang PetscErrorCode MatMultHermitianTransposeAdd_SeqBAIJ(Mat A,Vec xx,Vec yy,Vec zz)
1200547795f9SHong Zhang 
1201547795f9SHong Zhang {
1202547795f9SHong Zhang   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
1203547795f9SHong Zhang   PetscScalar       *zb,*x,*z,*xb = 0,x1,x2,x3,x4,x5;
1204547795f9SHong Zhang   MatScalar         *v;
1205547795f9SHong Zhang   PetscErrorCode    ierr;
1206547795f9SHong Zhang   PetscInt          mbs,i,*idx,*ii,rval,bs=A->rmap->bs,j,n,bs2=a->bs2,*ib,*ridx=PETSC_NULL;
1207547795f9SHong Zhang   Mat_CompressedRow cprow = a->compressedrow;
1208547795f9SHong Zhang   PetscTruth        usecprow=cprow.use;
1209547795f9SHong Zhang 
1210547795f9SHong Zhang   PetscFunctionBegin;
1211547795f9SHong Zhang   if (yy != zz) { ierr = VecCopy_Seq(yy,zz);CHKERRQ(ierr); }
1212547795f9SHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1213547795f9SHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
1214547795f9SHong Zhang 
1215547795f9SHong Zhang   idx = a->j;
1216547795f9SHong Zhang   v   = a->a;
1217547795f9SHong Zhang   if (usecprow){
1218547795f9SHong Zhang     mbs  = cprow.nrows;
1219547795f9SHong Zhang     ii   = cprow.i;
1220547795f9SHong Zhang     ridx = cprow.rindex;
1221547795f9SHong Zhang   } else {
1222547795f9SHong Zhang     mbs=a->mbs;
1223547795f9SHong Zhang     ii = a->i;
1224547795f9SHong Zhang     xb = x;
1225547795f9SHong Zhang   }
1226547795f9SHong Zhang 
1227547795f9SHong Zhang   switch (bs) {
1228547795f9SHong Zhang   case 1:
1229547795f9SHong Zhang     for (i=0; i<mbs; i++) {
1230547795f9SHong Zhang       if (usecprow) xb = x + ridx[i];
1231547795f9SHong Zhang       x1 = xb[0];
1232547795f9SHong Zhang       ib = idx + ii[0];
1233547795f9SHong Zhang       n  = ii[1] - ii[0]; ii++;
1234547795f9SHong Zhang       for (j=0; j<n; j++) {
1235547795f9SHong Zhang         rval    = ib[j];
1236547795f9SHong Zhang         z[rval] += PetscConj(*v) * x1;
1237547795f9SHong Zhang         v++;
1238547795f9SHong Zhang       }
1239547795f9SHong Zhang       if (!usecprow) xb++;
1240547795f9SHong Zhang     }
1241547795f9SHong Zhang     break;
1242547795f9SHong Zhang   case 2:
1243547795f9SHong Zhang     for (i=0; i<mbs; i++) {
1244547795f9SHong Zhang       if (usecprow) xb = x + 2*ridx[i];
1245547795f9SHong Zhang       x1 = xb[0]; x2 = xb[1];
1246547795f9SHong Zhang       ib = idx + ii[0];
1247547795f9SHong Zhang       n  = ii[1] - ii[0]; ii++;
1248547795f9SHong Zhang       for (j=0; j<n; j++) {
1249547795f9SHong Zhang         rval      = ib[j]*2;
1250547795f9SHong Zhang         z[rval++] += PetscConj(v[0])*x1 + PetscConj(v[1])*x2;
1251547795f9SHong Zhang         z[rval++] += PetscConj(v[2])*x1 + PetscConj(v[3])*x2;
1252547795f9SHong Zhang         v  += 4;
1253547795f9SHong Zhang       }
1254547795f9SHong Zhang       if (!usecprow) xb += 2;
1255547795f9SHong Zhang     }
1256547795f9SHong Zhang     break;
1257547795f9SHong Zhang   case 3:
1258547795f9SHong Zhang     for (i=0; i<mbs; i++) {
1259547795f9SHong Zhang       if (usecprow) xb = x + 3*ridx[i];
1260547795f9SHong Zhang       x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
1261547795f9SHong Zhang       ib = idx + ii[0];
1262547795f9SHong Zhang       n  = ii[1] - ii[0]; ii++;
1263547795f9SHong Zhang       for (j=0; j<n; j++) {
1264547795f9SHong Zhang         rval      = ib[j]*3;
1265547795f9SHong Zhang         z[rval++] += PetscConj(v[0])*x1 + PetscConj(v[1])*x2 + PetscConj(v[2])*x3;
1266547795f9SHong Zhang         z[rval++] += PetscConj(v[3])*x1 + PetscConj(v[4])*x2 + PetscConj(v[5])*x3;
1267547795f9SHong Zhang         z[rval++] += PetscConj(v[6])*x1 + PetscConj(v[7])*x2 + PetscConj(v[8])*x3;
1268547795f9SHong Zhang         v  += 9;
1269547795f9SHong Zhang       }
1270547795f9SHong Zhang       if (!usecprow) xb += 3;
1271547795f9SHong Zhang     }
1272547795f9SHong Zhang     break;
1273547795f9SHong Zhang   case 4:
1274547795f9SHong Zhang     for (i=0; i<mbs; i++) {
1275547795f9SHong Zhang       if (usecprow) xb = x + 4*ridx[i];
1276547795f9SHong Zhang       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
1277547795f9SHong Zhang       ib = idx + ii[0];
1278547795f9SHong Zhang       n  = ii[1] - ii[0]; ii++;
1279547795f9SHong Zhang       for (j=0; j<n; j++) {
1280547795f9SHong Zhang         rval      = ib[j]*4;
1281547795f9SHong Zhang         z[rval++] +=  PetscConj(v[0])*x1 + PetscConj(v[1])*x2  + PetscConj(v[2])*x3  + PetscConj(v[3])*x4;
1282547795f9SHong Zhang         z[rval++] +=  PetscConj(v[4])*x1 + PetscConj(v[5])*x2  + PetscConj(v[6])*x3  + PetscConj(v[7])*x4;
1283547795f9SHong Zhang         z[rval++] +=  PetscConj(v[8])*x1 + PetscConj(v[9])*x2  + PetscConj(v[10])*x3 + PetscConj(v[11])*x4;
1284547795f9SHong Zhang         z[rval++] += PetscConj(v[12])*x1 + PetscConj(v[13])*x2 + PetscConj(v[14])*x3 + PetscConj(v[15])*x4;
1285547795f9SHong Zhang         v  += 16;
1286547795f9SHong Zhang       }
1287547795f9SHong Zhang       if (!usecprow) xb += 4;
1288547795f9SHong Zhang     }
1289547795f9SHong Zhang     break;
1290547795f9SHong Zhang   case 5:
1291547795f9SHong Zhang     for (i=0; i<mbs; i++) {
1292547795f9SHong Zhang       if (usecprow) xb = x + 5*ridx[i];
1293547795f9SHong Zhang       x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
1294547795f9SHong Zhang       x4 = xb[3]; x5 = xb[4];
1295547795f9SHong Zhang       ib = idx + ii[0];
1296547795f9SHong Zhang       n  = ii[1] - ii[0]; ii++;
1297547795f9SHong Zhang       for (j=0; j<n; j++) {
1298547795f9SHong Zhang         rval      = ib[j]*5;
1299547795f9SHong Zhang         z[rval++] +=  PetscConj(v[0])*x1 +  PetscConj(v[1])*x2 +  PetscConj(v[2])*x3 +  PetscConj(v[3])*x4 +  PetscConj(v[4])*x5;
1300547795f9SHong Zhang         z[rval++] +=  PetscConj(v[5])*x1 +  PetscConj(v[6])*x2 +  PetscConj(v[7])*x3 +  PetscConj(v[8])*x4 +  PetscConj(v[9])*x5;
1301547795f9SHong Zhang         z[rval++] += PetscConj(v[10])*x1 + PetscConj(v[11])*x2 + PetscConj(v[12])*x3 + PetscConj(v[13])*x4 + PetscConj(v[14])*x5;
1302547795f9SHong Zhang         z[rval++] += PetscConj(v[15])*x1 + PetscConj(v[16])*x2 + PetscConj(v[17])*x3 + PetscConj(v[18])*x4 + PetscConj(v[19])*x5;
1303547795f9SHong Zhang         z[rval++] += PetscConj(v[20])*x1 + PetscConj(v[21])*x2 + PetscConj(v[22])*x3 + PetscConj(v[23])*x4 + PetscConj(v[24])*x5;
1304547795f9SHong Zhang         v  += 25;
1305547795f9SHong Zhang       }
1306547795f9SHong Zhang       if (!usecprow) xb += 5;
1307547795f9SHong Zhang     }
1308547795f9SHong Zhang     break;
1309547795f9SHong Zhang   default: {      /* block sizes larger than 5 by 5 are handled by BLAS */
1310547795f9SHong Zhang       PetscInt     ncols,k;
1311547795f9SHong Zhang       PetscScalar  *work,*workt,*xtmp;
1312547795f9SHong Zhang 
1313547795f9SHong Zhang       SETERRQ(PETSC_ERR_SUP,"block size larger than 5 is not supported yet");
1314547795f9SHong Zhang       if (!a->mult_work) {
1315547795f9SHong Zhang         k = PetscMax(A->rmap->n,A->cmap->n);
1316547795f9SHong Zhang         ierr = PetscMalloc((k+1)*sizeof(PetscScalar),&a->mult_work);CHKERRQ(ierr);
1317547795f9SHong Zhang       }
1318547795f9SHong Zhang       work = a->mult_work;
1319547795f9SHong Zhang       xtmp = x;
1320547795f9SHong Zhang       for (i=0; i<mbs; i++) {
1321547795f9SHong Zhang         n     = ii[1] - ii[0]; ii++;
1322547795f9SHong Zhang         ncols = n*bs;
1323547795f9SHong Zhang         ierr  = PetscMemzero(work,ncols*sizeof(PetscScalar));CHKERRQ(ierr);
1324547795f9SHong Zhang         if (usecprow) {
1325547795f9SHong Zhang           xtmp = x + bs*ridx[i];
1326547795f9SHong Zhang         }
1327547795f9SHong Zhang         Kernel_w_gets_w_plus_trans_Ar_times_v(bs,ncols,xtmp,v,work);
1328547795f9SHong Zhang         /* BLASgemv_("T",&bs,&ncols,&_DOne,v,&bs,xtmp,&_One,&_DOne,work,&_One); */
1329547795f9SHong Zhang         v += n*bs2;
1330547795f9SHong Zhang         if (!usecprow) xtmp += bs;
1331547795f9SHong Zhang         workt = work;
1332547795f9SHong Zhang         for (j=0; j<n; j++) {
1333547795f9SHong Zhang           zb = z + bs*(*idx++);
1334547795f9SHong Zhang           for (k=0; k<bs; k++) zb[k] += workt[k] ;
1335547795f9SHong Zhang           workt += bs;
1336547795f9SHong Zhang         }
1337547795f9SHong Zhang       }
1338547795f9SHong Zhang     }
1339547795f9SHong Zhang   }
1340547795f9SHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
1341547795f9SHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
1342547795f9SHong Zhang   ierr = PetscLogFlops(2.0*a->nz*a->bs2);CHKERRQ(ierr);
1343547795f9SHong Zhang   PetscFunctionReturn(0);
1344547795f9SHong Zhang }
1345547795f9SHong Zhang 
1346547795f9SHong Zhang #undef __FUNCT__
13474a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_SeqBAIJ"
1348dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_SeqBAIJ(Mat A,Vec xx,Vec yy,Vec zz)
13492d61bbb3SSatish Balay 
13502d61bbb3SSatish Balay {
13512d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
1352dcf5cc72SBarry Smith   PetscScalar       *zb,*x,*z,*xb = 0,x1,x2,x3,x4,x5;
13533f1db9ecSBarry Smith   MatScalar         *v;
13546849ba73SBarry Smith   PetscErrorCode    ierr;
1355d0f46423SBarry Smith   PetscInt          mbs,i,*idx,*ii,rval,bs=A->rmap->bs,j,n,bs2=a->bs2,*ib,*ridx=PETSC_NULL;
13563447b6efSHong Zhang   Mat_CompressedRow cprow = a->compressedrow;
13573447b6efSHong Zhang   PetscTruth        usecprow=cprow.use;
13582d61bbb3SSatish Balay 
13592d61bbb3SSatish Balay   PetscFunctionBegin;
13606a65c766SHong Zhang   if (yy != zz) { ierr = VecCopy_Seq(yy,zz);CHKERRQ(ierr); }
13613447b6efSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
13623447b6efSHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
13632d61bbb3SSatish Balay 
13642d61bbb3SSatish Balay   idx = a->j;
13652d61bbb3SSatish Balay   v   = a->a;
13663447b6efSHong Zhang   if (usecprow){
13673447b6efSHong Zhang     mbs  = cprow.nrows;
13683447b6efSHong Zhang     ii   = cprow.i;
13697b2bb3b9SHong Zhang     ridx = cprow.rindex;
13703447b6efSHong Zhang   } else {
13713447b6efSHong Zhang     mbs=a->mbs;
13722d61bbb3SSatish Balay     ii = a->i;
1373f1af5d2fSBarry Smith     xb = x;
13743447b6efSHong Zhang   }
13752d61bbb3SSatish Balay 
13762d61bbb3SSatish Balay   switch (bs) {
13772d61bbb3SSatish Balay   case 1:
13782d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
13797b2bb3b9SHong Zhang       if (usecprow) xb = x + ridx[i];
1380f1af5d2fSBarry Smith       x1 = xb[0];
13813447b6efSHong Zhang       ib = idx + ii[0];
13823447b6efSHong Zhang       n  = ii[1] - ii[0]; ii++;
13832d61bbb3SSatish Balay       for (j=0; j<n; j++) {
13842d61bbb3SSatish Balay         rval    = ib[j];
1385f1af5d2fSBarry Smith         z[rval] += *v * x1;
1386f1af5d2fSBarry Smith         v++;
13872d61bbb3SSatish Balay       }
13883447b6efSHong Zhang       if (!usecprow) xb++;
13892d61bbb3SSatish Balay     }
13902d61bbb3SSatish Balay     break;
13912d61bbb3SSatish Balay   case 2:
13922d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
13937b2bb3b9SHong Zhang       if (usecprow) xb = x + 2*ridx[i];
1394f1af5d2fSBarry Smith       x1 = xb[0]; x2 = xb[1];
13953447b6efSHong Zhang       ib = idx + ii[0];
13963447b6efSHong Zhang       n  = ii[1] - ii[0]; ii++;
13972d61bbb3SSatish Balay       for (j=0; j<n; j++) {
13982d61bbb3SSatish Balay         rval      = ib[j]*2;
13992d61bbb3SSatish Balay         z[rval++] += v[0]*x1 + v[1]*x2;
14002d61bbb3SSatish Balay         z[rval++] += v[2]*x1 + v[3]*x2;
14012d61bbb3SSatish Balay         v  += 4;
14022d61bbb3SSatish Balay       }
14033447b6efSHong Zhang       if (!usecprow) xb += 2;
14042d61bbb3SSatish Balay     }
14052d61bbb3SSatish Balay     break;
14062d61bbb3SSatish Balay   case 3:
14072d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
14087b2bb3b9SHong Zhang       if (usecprow) xb = x + 3*ridx[i];
1409f1af5d2fSBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
14103447b6efSHong Zhang       ib = idx + ii[0];
14113447b6efSHong Zhang       n  = ii[1] - ii[0]; ii++;
14122d61bbb3SSatish Balay       for (j=0; j<n; j++) {
14132d61bbb3SSatish Balay         rval      = ib[j]*3;
14142d61bbb3SSatish Balay         z[rval++] += v[0]*x1 + v[1]*x2 + v[2]*x3;
14152d61bbb3SSatish Balay         z[rval++] += v[3]*x1 + v[4]*x2 + v[5]*x3;
14162d61bbb3SSatish Balay         z[rval++] += v[6]*x1 + v[7]*x2 + v[8]*x3;
14172d61bbb3SSatish Balay         v  += 9;
14182d61bbb3SSatish Balay       }
14193447b6efSHong Zhang       if (!usecprow) xb += 3;
14202d61bbb3SSatish Balay     }
14212d61bbb3SSatish Balay     break;
14222d61bbb3SSatish Balay   case 4:
14232d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
14247b2bb3b9SHong Zhang       if (usecprow) xb = x + 4*ridx[i];
1425f1af5d2fSBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
14263447b6efSHong Zhang       ib = idx + ii[0];
14273447b6efSHong Zhang       n  = ii[1] - ii[0]; ii++;
14282d61bbb3SSatish Balay       for (j=0; j<n; j++) {
14292d61bbb3SSatish Balay         rval      = ib[j]*4;
14302d61bbb3SSatish Balay         z[rval++] +=  v[0]*x1 +  v[1]*x2 +  v[2]*x3 +  v[3]*x4;
14312d61bbb3SSatish Balay         z[rval++] +=  v[4]*x1 +  v[5]*x2 +  v[6]*x3 +  v[7]*x4;
14322d61bbb3SSatish Balay         z[rval++] +=  v[8]*x1 +  v[9]*x2 + v[10]*x3 + v[11]*x4;
14332d61bbb3SSatish Balay         z[rval++] += v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4;
14342d61bbb3SSatish Balay         v  += 16;
14352d61bbb3SSatish Balay       }
14363447b6efSHong Zhang       if (!usecprow) xb += 4;
14372d61bbb3SSatish Balay     }
14382d61bbb3SSatish Balay     break;
14392d61bbb3SSatish Balay   case 5:
14402d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
14417b2bb3b9SHong Zhang       if (usecprow) xb = x + 5*ridx[i];
1442f1af5d2fSBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
14432d61bbb3SSatish Balay       x4 = xb[3]; x5 = xb[4];
14443447b6efSHong Zhang       ib = idx + ii[0];
14453447b6efSHong Zhang       n  = ii[1] - ii[0]; ii++;
14462d61bbb3SSatish Balay       for (j=0; j<n; j++) {
14472d61bbb3SSatish Balay         rval      = ib[j]*5;
14482d61bbb3SSatish Balay         z[rval++] +=  v[0]*x1 +  v[1]*x2 +  v[2]*x3 +  v[3]*x4 +  v[4]*x5;
14492d61bbb3SSatish Balay         z[rval++] +=  v[5]*x1 +  v[6]*x2 +  v[7]*x3 +  v[8]*x4 +  v[9]*x5;
14502d61bbb3SSatish Balay         z[rval++] += v[10]*x1 + v[11]*x2 + v[12]*x3 + v[13]*x4 + v[14]*x5;
14512d61bbb3SSatish Balay         z[rval++] += v[15]*x1 + v[16]*x2 + v[17]*x3 + v[18]*x4 + v[19]*x5;
14522d61bbb3SSatish Balay         z[rval++] += v[20]*x1 + v[21]*x2 + v[22]*x3 + v[23]*x4 + v[24]*x5;
14532d61bbb3SSatish Balay         v  += 25;
14542d61bbb3SSatish Balay       }
14553447b6efSHong Zhang       if (!usecprow) xb += 5;
14562d61bbb3SSatish Balay     }
14572d61bbb3SSatish Balay     break;
1458f1af5d2fSBarry Smith   default: {      /* block sizes larger then 5 by 5 are handled by BLAS */
1459690b6cddSBarry Smith       PetscInt     ncols,k;
14603447b6efSHong Zhang       PetscScalar  *work,*workt,*xtmp;
14613f1db9ecSBarry Smith 
14622d61bbb3SSatish Balay       if (!a->mult_work) {
1463d0f46423SBarry Smith         k = PetscMax(A->rmap->n,A->cmap->n);
146487828ca2SBarry Smith         ierr = PetscMalloc((k+1)*sizeof(PetscScalar),&a->mult_work);CHKERRQ(ierr);
14652d61bbb3SSatish Balay       }
14662d61bbb3SSatish Balay       work = a->mult_work;
14673447b6efSHong Zhang       xtmp = x;
14682d61bbb3SSatish Balay       for (i=0; i<mbs; i++) {
14692d61bbb3SSatish Balay         n     = ii[1] - ii[0]; ii++;
14702d61bbb3SSatish Balay         ncols = n*bs;
147187828ca2SBarry Smith         ierr  = PetscMemzero(work,ncols*sizeof(PetscScalar));CHKERRQ(ierr);
14723447b6efSHong Zhang         if (usecprow) {
14737b2bb3b9SHong Zhang           xtmp = x + bs*ridx[i];
14743447b6efSHong Zhang         }
14753447b6efSHong Zhang         Kernel_w_gets_w_plus_trans_Ar_times_v(bs,ncols,xtmp,v,work);
147671044d3cSBarry Smith         /* BLASgemv_("T",&bs,&ncols,&_DOne,v,&bs,xtmp,&_One,&_DOne,work,&_One); */
14772d61bbb3SSatish Balay         v += n*bs2;
14783447b6efSHong Zhang         if (!usecprow) xtmp += bs;
14792d61bbb3SSatish Balay         workt = work;
14802d61bbb3SSatish Balay         for (j=0; j<n; j++) {
14812d61bbb3SSatish Balay           zb = z + bs*(*idx++);
14822d61bbb3SSatish Balay           for (k=0; k<bs; k++) zb[k] += workt[k] ;
14832d61bbb3SSatish Balay           workt += bs;
14842d61bbb3SSatish Balay         }
14852d61bbb3SSatish Balay       }
14862d61bbb3SSatish Balay     }
14872d61bbb3SSatish Balay   }
14883447b6efSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
14893447b6efSHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
1490dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*a->nz*a->bs2);CHKERRQ(ierr);
14912d61bbb3SSatish Balay   PetscFunctionReturn(0);
14922d61bbb3SSatish Balay }
14932d61bbb3SSatish Balay 
14944a2ae208SSatish Balay #undef __FUNCT__
14954a2ae208SSatish Balay #define __FUNCT__ "MatScale_SeqBAIJ"
1496f4df32b1SMatthew Knepley PetscErrorCode MatScale_SeqBAIJ(Mat inA,PetscScalar alpha)
14972d61bbb3SSatish Balay {
14982d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)inA->data;
1499690b6cddSBarry Smith   PetscInt       totalnz = a->bs2*a->nz;
1500f4df32b1SMatthew Knepley   PetscScalar    oalpha = alpha;
1501efee365bSSatish Balay   PetscErrorCode ierr;
15020805154bSBarry Smith   PetscBLASInt   one = 1,tnz = PetscBLASIntCast(totalnz);
15032d61bbb3SSatish Balay 
15042d61bbb3SSatish Balay   PetscFunctionBegin;
1505f4df32b1SMatthew Knepley   BLASscal_(&tnz,&oalpha,a->a,&one);
1506efee365bSSatish Balay   ierr = PetscLogFlops(totalnz);CHKERRQ(ierr);
15072d61bbb3SSatish Balay   PetscFunctionReturn(0);
15082d61bbb3SSatish Balay }
15092d61bbb3SSatish Balay 
15104a2ae208SSatish Balay #undef __FUNCT__
15114a2ae208SSatish Balay #define __FUNCT__ "MatNorm_SeqBAIJ"
1512dfbe8321SBarry Smith PetscErrorCode MatNorm_SeqBAIJ(Mat A,NormType type,PetscReal *norm)
15132d61bbb3SSatish Balay {
15148a62d963SHong Zhang   PetscErrorCode ierr;
15152d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
15163f1db9ecSBarry Smith   MatScalar      *v = a->a;
1517329f5518SBarry Smith   PetscReal      sum = 0.0;
1518d0f46423SBarry Smith   PetscInt       i,j,k,bs=A->rmap->bs,nz=a->nz,bs2=a->bs2,k1;
15192d61bbb3SSatish Balay 
15202d61bbb3SSatish Balay   PetscFunctionBegin;
15212d61bbb3SSatish Balay   if (type == NORM_FROBENIUS) {
15222d61bbb3SSatish Balay     for (i=0; i< bs2*nz; i++) {
1523aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
1524329f5518SBarry Smith       sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
15252d61bbb3SSatish Balay #else
15262d61bbb3SSatish Balay       sum += (*v)*(*v); v++;
15272d61bbb3SSatish Balay #endif
15282d61bbb3SSatish Balay     }
15292d61bbb3SSatish Balay     *norm = sqrt(sum);
15308a62d963SHong Zhang   } else if (type == NORM_1) { /* maximum column sum */
15318a62d963SHong Zhang     PetscReal *tmp;
15328a62d963SHong Zhang     PetscInt  *bcol = a->j;
1533d0f46423SBarry Smith     ierr = PetscMalloc((A->cmap->n+1)*sizeof(PetscReal),&tmp);CHKERRQ(ierr);
1534d0f46423SBarry Smith     ierr = PetscMemzero(tmp,A->cmap->n*sizeof(PetscReal));CHKERRQ(ierr);
15358a62d963SHong Zhang     for (i=0; i<nz; i++){
15368a62d963SHong Zhang       for (j=0; j<bs; j++){
15378a62d963SHong Zhang         k1 = bs*(*bcol) + j; /* column index */
15388a62d963SHong Zhang         for (k=0; k<bs; k++){
15398a62d963SHong Zhang           tmp[k1] += PetscAbsScalar(*v); v++;
15408a62d963SHong Zhang         }
15418a62d963SHong Zhang       }
15428a62d963SHong Zhang       bcol++;
15438a62d963SHong Zhang     }
15448a62d963SHong Zhang     *norm = 0.0;
1545d0f46423SBarry Smith     for (j=0; j<A->cmap->n; j++) {
15468a62d963SHong Zhang       if (tmp[j] > *norm) *norm = tmp[j];
15478a62d963SHong Zhang     }
15488a62d963SHong Zhang     ierr = PetscFree(tmp);CHKERRQ(ierr);
1549596552b5SBarry Smith   } else if (type == NORM_INFINITY) { /* maximum row sum */
1550596552b5SBarry Smith     *norm = 0.0;
1551596552b5SBarry Smith     for (k=0; k<bs; k++) {
155274f84c7bSSatish Balay       for (j=0; j<a->mbs; j++) {
1553596552b5SBarry Smith         v = a->a + bs2*a->i[j] + k;
1554596552b5SBarry Smith         sum = 0.0;
1555596552b5SBarry Smith         for (i=0; i<a->i[j+1]-a->i[j]; i++) {
15560e90e235SBarry Smith           for (k1=0; k1<bs; k1++){
1557596552b5SBarry Smith             sum += PetscAbsScalar(*v);
1558596552b5SBarry Smith             v   += bs;
15592d61bbb3SSatish Balay           }
15600e90e235SBarry Smith         }
1561596552b5SBarry Smith         if (sum > *norm) *norm = sum;
1562596552b5SBarry Smith       }
1563596552b5SBarry Smith     }
1564596552b5SBarry Smith   } else {
156529bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"No support for this norm yet");
15662d61bbb3SSatish Balay   }
15672d61bbb3SSatish Balay   PetscFunctionReturn(0);
15682d61bbb3SSatish Balay }
15692d61bbb3SSatish Balay 
15702d61bbb3SSatish Balay 
15714a2ae208SSatish Balay #undef __FUNCT__
15724a2ae208SSatish Balay #define __FUNCT__ "MatEqual_SeqBAIJ"
1573dfbe8321SBarry Smith PetscErrorCode MatEqual_SeqBAIJ(Mat A,Mat B,PetscTruth* flg)
15742d61bbb3SSatish Balay {
15752d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ *)A->data,*b = (Mat_SeqBAIJ *)B->data;
1576dfbe8321SBarry Smith   PetscErrorCode ierr;
15772d61bbb3SSatish Balay 
15782d61bbb3SSatish Balay   PetscFunctionBegin;
15792d61bbb3SSatish Balay   /* If the  matrix/block dimensions are not equal, or no of nonzeros or shift */
1580d0f46423SBarry Smith   if ((A->rmap->N != B->rmap->N) || (A->cmap->n != B->cmap->n) || (A->rmap->bs != B->rmap->bs)|| (a->nz != b->nz)) {
1581273d9f13SBarry Smith     *flg = PETSC_FALSE;
1582273d9f13SBarry Smith     PetscFunctionReturn(0);
15832d61bbb3SSatish Balay   }
15842d61bbb3SSatish Balay 
15852d61bbb3SSatish Balay   /* if the a->i are the same */
1586690b6cddSBarry Smith   ierr = PetscMemcmp(a->i,b->i,(a->mbs+1)*sizeof(PetscInt),flg);CHKERRQ(ierr);
1587abc0a331SBarry Smith   if (!*flg) {
15880f5bd95cSBarry Smith     PetscFunctionReturn(0);
15892d61bbb3SSatish Balay   }
15902d61bbb3SSatish Balay 
15912d61bbb3SSatish Balay   /* if a->j are the same */
1592690b6cddSBarry Smith   ierr = PetscMemcmp(a->j,b->j,(a->nz)*sizeof(PetscInt),flg);CHKERRQ(ierr);
1593abc0a331SBarry Smith   if (!*flg) {
15940f5bd95cSBarry Smith     PetscFunctionReturn(0);
15952d61bbb3SSatish Balay   }
15962d61bbb3SSatish Balay   /* if a->a are the same */
1597d0f46423SBarry Smith   ierr = PetscMemcmp(a->a,b->a,(a->nz)*(A->rmap->bs)*(B->rmap->bs)*sizeof(PetscScalar),flg);CHKERRQ(ierr);
15982d61bbb3SSatish Balay   PetscFunctionReturn(0);
15992d61bbb3SSatish Balay 
16002d61bbb3SSatish Balay }
16012d61bbb3SSatish Balay 
16024a2ae208SSatish Balay #undef __FUNCT__
16034a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_SeqBAIJ"
1604dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_SeqBAIJ(Mat A,Vec v)
16052d61bbb3SSatish Balay {
16062d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
1607dfbe8321SBarry Smith   PetscErrorCode ierr;
1608690b6cddSBarry Smith   PetscInt       i,j,k,n,row,bs,*ai,*aj,ambs,bs2;
160987828ca2SBarry Smith   PetscScalar    *x,zero = 0.0;
16103f1db9ecSBarry Smith   MatScalar      *aa,*aa_j;
16112d61bbb3SSatish Balay 
16122d61bbb3SSatish Balay   PetscFunctionBegin;
161329bbc08cSBarry Smith   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1614d0f46423SBarry Smith   bs   = A->rmap->bs;
16152d61bbb3SSatish Balay   aa   = a->a;
16162d61bbb3SSatish Balay   ai   = a->i;
16172d61bbb3SSatish Balay   aj   = a->j;
16182d61bbb3SSatish Balay   ambs = a->mbs;
16192d61bbb3SSatish Balay   bs2  = a->bs2;
16202d61bbb3SSatish Balay 
16212dcb1b2aSMatthew Knepley   ierr = VecSet(v,zero);CHKERRQ(ierr);
16221ebc52fbSHong Zhang   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
1623e1311b90SBarry Smith   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
1624d0f46423SBarry Smith   if (n != A->rmap->N) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
16252d61bbb3SSatish Balay   for (i=0; i<ambs; i++) {
16262d61bbb3SSatish Balay     for (j=ai[i]; j<ai[i+1]; j++) {
16272d61bbb3SSatish Balay       if (aj[j] == i) {
16282d61bbb3SSatish Balay         row  = i*bs;
16292d61bbb3SSatish Balay         aa_j = aa+j*bs2;
16302d61bbb3SSatish Balay         for (k=0; k<bs2; k+=(bs+1),row++) x[row] = aa_j[k];
16312d61bbb3SSatish Balay         break;
16322d61bbb3SSatish Balay       }
16332d61bbb3SSatish Balay     }
16342d61bbb3SSatish Balay   }
16351ebc52fbSHong Zhang   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
16362d61bbb3SSatish Balay   PetscFunctionReturn(0);
16372d61bbb3SSatish Balay }
16382d61bbb3SSatish Balay 
16394a2ae208SSatish Balay #undef __FUNCT__
16404a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_SeqBAIJ"
1641dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_SeqBAIJ(Mat A,Vec ll,Vec rr)
16422d61bbb3SSatish Balay {
16432d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
164487828ca2SBarry Smith   PetscScalar    *l,*r,x,*li,*ri;
16453f1db9ecSBarry Smith   MatScalar      *aa,*v;
1646dfbe8321SBarry Smith   PetscErrorCode ierr;
1647690b6cddSBarry Smith   PetscInt       i,j,k,lm,rn,M,m,n,*ai,*aj,mbs,tmp,bs,bs2;
16482d61bbb3SSatish Balay 
16492d61bbb3SSatish Balay   PetscFunctionBegin;
16502d61bbb3SSatish Balay   ai  = a->i;
16512d61bbb3SSatish Balay   aj  = a->j;
16522d61bbb3SSatish Balay   aa  = a->a;
1653d0f46423SBarry Smith   m   = A->rmap->n;
1654d0f46423SBarry Smith   n   = A->cmap->n;
1655d0f46423SBarry Smith   bs  = A->rmap->bs;
16562d61bbb3SSatish Balay   mbs = a->mbs;
16572d61bbb3SSatish Balay   bs2 = a->bs2;
16582d61bbb3SSatish Balay   if (ll) {
16591ebc52fbSHong Zhang     ierr = VecGetArray(ll,&l);CHKERRQ(ierr);
1660e1311b90SBarry Smith     ierr = VecGetLocalSize(ll,&lm);CHKERRQ(ierr);
166129bbc08cSBarry Smith     if (lm != m) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vector wrong length");
16622d61bbb3SSatish Balay     for (i=0; i<mbs; i++) { /* for each block row */
16632d61bbb3SSatish Balay       M  = ai[i+1] - ai[i];
16642d61bbb3SSatish Balay       li = l + i*bs;
16652d61bbb3SSatish Balay       v  = aa + bs2*ai[i];
16662d61bbb3SSatish Balay       for (j=0; j<M; j++) { /* for each block */
16672d61bbb3SSatish Balay         for (k=0; k<bs2; k++) {
16682d61bbb3SSatish Balay           (*v++) *= li[k%bs];
16692d61bbb3SSatish Balay         }
16702d61bbb3SSatish Balay       }
16712d61bbb3SSatish Balay     }
16721ebc52fbSHong Zhang     ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr);
1673efee365bSSatish Balay     ierr = PetscLogFlops(a->nz);CHKERRQ(ierr);
16742d61bbb3SSatish Balay   }
16752d61bbb3SSatish Balay 
16762d61bbb3SSatish Balay   if (rr) {
16771ebc52fbSHong Zhang     ierr = VecGetArray(rr,&r);CHKERRQ(ierr);
1678e1311b90SBarry Smith     ierr = VecGetLocalSize(rr,&rn);CHKERRQ(ierr);
167929bbc08cSBarry Smith     if (rn != n) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vector wrong length");
16802d61bbb3SSatish Balay     for (i=0; i<mbs; i++) { /* for each block row */
16812d61bbb3SSatish Balay       M  = ai[i+1] - ai[i];
16822d61bbb3SSatish Balay       v  = aa + bs2*ai[i];
16832d61bbb3SSatish Balay       for (j=0; j<M; j++) { /* for each block */
16842d61bbb3SSatish Balay         ri = r + bs*aj[ai[i]+j];
16852d61bbb3SSatish Balay         for (k=0; k<bs; k++) {
16862d61bbb3SSatish Balay           x = ri[k];
16872d61bbb3SSatish Balay           for (tmp=0; tmp<bs; tmp++) (*v++) *= x;
16882d61bbb3SSatish Balay         }
16892d61bbb3SSatish Balay       }
16902d61bbb3SSatish Balay     }
16911ebc52fbSHong Zhang     ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr);
1692efee365bSSatish Balay     ierr = PetscLogFlops(a->nz);CHKERRQ(ierr);
16932d61bbb3SSatish Balay   }
16942d61bbb3SSatish Balay   PetscFunctionReturn(0);
16952d61bbb3SSatish Balay }
16962d61bbb3SSatish Balay 
16972d61bbb3SSatish Balay 
16984a2ae208SSatish Balay #undef __FUNCT__
16994a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_SeqBAIJ"
1700dfbe8321SBarry Smith PetscErrorCode MatGetInfo_SeqBAIJ(Mat A,MatInfoType flag,MatInfo *info)
17012d61bbb3SSatish Balay {
17022d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
17032d61bbb3SSatish Balay 
17042d61bbb3SSatish Balay   PetscFunctionBegin;
17052d61bbb3SSatish Balay   info->block_size     = a->bs2;
17062d61bbb3SSatish Balay   info->nz_allocated   = a->maxnz;
17072d61bbb3SSatish Balay   info->nz_used        = a->bs2*a->nz;
17082d61bbb3SSatish Balay   info->nz_unneeded    = (double)(info->nz_allocated - info->nz_used);
17092d61bbb3SSatish Balay   info->assemblies   = A->num_ass;
17102d61bbb3SSatish Balay   info->mallocs      = a->reallocs;
17117adad957SLisandro Dalcin   info->memory       = ((PetscObject)A)->mem;
17122d61bbb3SSatish Balay   if (A->factor) {
17132d61bbb3SSatish Balay     info->fill_ratio_given  = A->info.fill_ratio_given;
17142d61bbb3SSatish Balay     info->fill_ratio_needed = A->info.fill_ratio_needed;
17152d61bbb3SSatish Balay     info->factor_mallocs    = A->info.factor_mallocs;
17162d61bbb3SSatish Balay   } else {
17172d61bbb3SSatish Balay     info->fill_ratio_given  = 0;
17182d61bbb3SSatish Balay     info->fill_ratio_needed = 0;
17192d61bbb3SSatish Balay     info->factor_mallocs    = 0;
17202d61bbb3SSatish Balay   }
17212d61bbb3SSatish Balay   PetscFunctionReturn(0);
17222d61bbb3SSatish Balay }
17232d61bbb3SSatish Balay 
17242d61bbb3SSatish Balay 
17254a2ae208SSatish Balay #undef __FUNCT__
17264a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_SeqBAIJ"
1727dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_SeqBAIJ(Mat A)
17282d61bbb3SSatish Balay {
17292d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
1730dfbe8321SBarry Smith   PetscErrorCode ierr;
17312d61bbb3SSatish Balay 
17322d61bbb3SSatish Balay   PetscFunctionBegin;
1733549d3d68SSatish Balay   ierr = PetscMemzero(a->a,a->bs2*a->i[a->mbs]*sizeof(MatScalar));CHKERRQ(ierr);
17342d61bbb3SSatish Balay   PetscFunctionReturn(0);
17352d61bbb3SSatish Balay }
1736