xref: /petsc/src/mat/impls/baij/seq/baij2.c (revision 8b83055f287c8deb1a16325ae3faee24d0648b7a)
1cac129eeSSatish Balay 
2c6db04a5SJed Brown #include <../src/mat/impls/baij/seq/baij.h>
306873bf2SBarry Smith #include <petsc-private/kernels/blockinvert.h>
4c6db04a5SJed Brown #include <petscbt.h>
5c6db04a5SJed Brown #include <petscblaslapack.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 
24e32f2f54SBarry Smith   if (ov < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap specified");
25a3192f15SSatish Balay 
2653b8de81SBarry 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 */
42e32f2f54SBarry Smith       if (ival>=m) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"index greater than mat-dim");
4326fbe8dcSKarl Rupp       if (!PetscBTLookupSet(table,ival)) nidx[isz++] = ival;
44a3192f15SSatish Balay     }
45a3192f15SSatish Balay     ierr = ISRestoreIndices(is[i],&idx);CHKERRQ(ierr);
466bf464f9SBarry Smith     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];
5726fbe8dcSKarl Rupp           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++) {
6326fbe8dcSKarl Rupp       for (k=0; k<bs; k++) nidx2[j*bs+k] = nidx[j]*bs+k;
64218c64b6SSatish Balay     }
6570b3c8c7SBarry Smith     ierr = ISCreateGeneral(PETSC_COMM_SELF,isz*bs,nidx2,PETSC_COPY_VALUES,is+i);CHKERRQ(ierr);
66a3192f15SSatish Balay   }
6794bacf5dSBarry Smith   ierr = PetscBTDestroy(&table);CHKERRQ(ierr);
68606d414cSSatish Balay   ierr = PetscFree(nidx);CHKERRQ(ierr);
69606d414cSSatish Balay   ierr = PetscFree(nidx2);CHKERRQ(ierr);
703a40ed3dSBarry Smith   PetscFunctionReturn(0);
71a3192f15SSatish Balay }
721c351548SSatish Balay 
734a2ae208SSatish Balay #undef __FUNCT__
744a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqBAIJ_Private"
754aa3045dSJed Brown PetscErrorCode MatGetSubMatrix_SeqBAIJ_Private(Mat A,IS isrow,IS iscol,MatReuse scall,Mat *B)
76736121d4SSatish Balay {
77736121d4SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*c;
786849ba73SBarry Smith   PetscErrorCode ierr;
79690b6cddSBarry Smith   PetscInt       *smap,i,k,kstart,kend,oldcols = a->nbs,*lens;
80690b6cddSBarry Smith   PetscInt       row,mat_i,*mat_j,tcol,*mat_ilen;
815d0c19d7SBarry Smith   const PetscInt *irow,*icol;
825d0c19d7SBarry Smith   PetscInt       nrows,ncols,*ssmap,bs=A->rmap->bs,bs2=a->bs2;
83690b6cddSBarry Smith   PetscInt       *aj = a->j,*ai = a->i;
843f1db9ecSBarry Smith   MatScalar      *mat_a;
85736121d4SSatish Balay   Mat            C;
86ace3abfcSBarry Smith   PetscBool      flag,sorted;
87736121d4SSatish Balay 
883a40ed3dSBarry Smith   PetscFunctionBegin;
8914ca34e6SBarry Smith   ierr = ISSorted(iscol,&sorted);CHKERRQ(ierr);
90e32f2f54SBarry Smith   if (!sorted) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"IS is not sorted");
91736121d4SSatish Balay 
92736121d4SSatish Balay   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
93218c64b6SSatish Balay   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
94b9b97703SBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
95b9b97703SBarry Smith   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
96736121d4SSatish Balay 
97690b6cddSBarry Smith   ierr  = PetscMalloc((1+oldcols)*sizeof(PetscInt),&smap);CHKERRQ(ierr);
98736121d4SSatish Balay   ssmap = smap;
99690b6cddSBarry Smith   ierr  = PetscMalloc((1+nrows)*sizeof(PetscInt),&lens);CHKERRQ(ierr);
100690b6cddSBarry Smith   ierr  = PetscMemzero(smap,oldcols*sizeof(PetscInt));CHKERRQ(ierr);
101736121d4SSatish Balay   for (i=0; i<ncols; i++) smap[icol[i]] = i+1;
102736121d4SSatish Balay   /* determine lens of each row */
103736121d4SSatish Balay   for (i=0; i<nrows; i++) {
104736121d4SSatish Balay     kstart  = ai[irow[i]];
105736121d4SSatish Balay     kend    = kstart + a->ilen[irow[i]];
106736121d4SSatish Balay     lens[i] = 0;
107736121d4SSatish Balay     for (k=kstart; k<kend; k++) {
10826fbe8dcSKarl Rupp       if (ssmap[aj[k]]) lens[i]++;
109736121d4SSatish Balay     }
110736121d4SSatish Balay   }
111736121d4SSatish Balay   /* Create and fill new matrix */
112736121d4SSatish Balay   if (scall == MAT_REUSE_MATRIX) {
113736121d4SSatish Balay     c = (Mat_SeqBAIJ*)((*B)->data);
114736121d4SSatish Balay 
115e32f2f54SBarry Smith     if (c->mbs!=nrows || c->nbs!=ncols || (*B)->rmap->bs!=bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Submatrix wrong size");
116690b6cddSBarry Smith     ierr = PetscMemcmp(c->ilen,lens,c->mbs *sizeof(PetscInt),&flag);CHKERRQ(ierr);
117e7e72b3dSBarry Smith     if (!flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong no of nonzeros");
118690b6cddSBarry Smith     ierr = PetscMemzero(c->ilen,c->mbs*sizeof(PetscInt));CHKERRQ(ierr);
119736121d4SSatish Balay     C    = *B;
1203a40ed3dSBarry Smith   } else {
121ce94432eSBarry Smith     ierr = MatCreate(PetscObjectComm((PetscObject)A),&C);CHKERRQ(ierr);
122f69a0ea3SMatthew Knepley     ierr = MatSetSizes(C,nrows*bs,ncols*bs,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1237adad957SLisandro Dalcin     ierr = MatSetType(C,((PetscObject)A)->type_name);CHKERRQ(ierr);
124ab93d7beSBarry Smith     ierr = MatSeqBAIJSetPreallocation_SeqBAIJ(C,bs,0,lens);CHKERRQ(ierr);
125736121d4SSatish Balay   }
126736121d4SSatish Balay   c = (Mat_SeqBAIJ*)(C->data);
127736121d4SSatish Balay   for (i=0; i<nrows; i++) {
128736121d4SSatish Balay     row      = irow[i];
129736121d4SSatish Balay     kstart   = ai[row];
130736121d4SSatish Balay     kend     = kstart + a->ilen[row];
131736121d4SSatish Balay     mat_i    = c->i[i];
132736121d4SSatish Balay     mat_j    = c->j + mat_i;
133218c64b6SSatish Balay     mat_a    = c->a + mat_i*bs2;
134736121d4SSatish Balay     mat_ilen = c->ilen + i;
135736121d4SSatish Balay     for (k=kstart; k<kend; k++) {
136736121d4SSatish Balay       if ((tcol=ssmap[a->j[k]])) {
137736121d4SSatish Balay         *mat_j++ = tcol - 1;
138549d3d68SSatish Balay         ierr     = PetscMemcpy(mat_a,a->a+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
139549d3d68SSatish Balay         mat_a   += bs2;
140736121d4SSatish Balay         (*mat_ilen)++;
141736121d4SSatish Balay       }
142736121d4SSatish Balay     }
143736121d4SSatish Balay   }
144218c64b6SSatish Balay 
145736121d4SSatish Balay   /* Free work space */
146736121d4SSatish Balay   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
147606d414cSSatish Balay   ierr = PetscFree(smap);CHKERRQ(ierr);
148606d414cSSatish Balay   ierr = PetscFree(lens);CHKERRQ(ierr);
1496d4a8577SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1506d4a8577SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
151736121d4SSatish Balay 
152736121d4SSatish Balay   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
153736121d4SSatish Balay   *B   = C;
1543a40ed3dSBarry Smith   PetscFunctionReturn(0);
155736121d4SSatish Balay }
156736121d4SSatish Balay 
1574a2ae208SSatish Balay #undef __FUNCT__
1584a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqBAIJ"
1594aa3045dSJed Brown PetscErrorCode MatGetSubMatrix_SeqBAIJ(Mat A,IS isrow,IS iscol,MatReuse scall,Mat *B)
160218c64b6SSatish Balay {
161218c64b6SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
162218c64b6SSatish Balay   IS             is1,is2;
1636849ba73SBarry Smith   PetscErrorCode ierr;
1645d0c19d7SBarry Smith   PetscInt       *vary,*iary,nrows,ncols,i,bs=A->rmap->bs,count;
1655d0c19d7SBarry Smith   const PetscInt *irow,*icol;
166218c64b6SSatish Balay 
1673a40ed3dSBarry Smith   PetscFunctionBegin;
168218c64b6SSatish Balay   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
169218c64b6SSatish Balay   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
170b9b97703SBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
171b9b97703SBarry Smith   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
172218c64b6SSatish Balay 
173218c64b6SSatish Balay   /* Verify if the indices corespond to each element in a block
174218c64b6SSatish Balay    and form the IS with compressed IS */
175fca92195SBarry Smith   ierr = PetscMalloc2(a->mbs,PetscInt,&vary,a->mbs,PetscInt,&iary);CHKERRQ(ierr);
176fca92195SBarry Smith   ierr = PetscMemzero(vary,a->mbs*sizeof(PetscInt));CHKERRQ(ierr);
177218c64b6SSatish Balay   for (i=0; i<nrows; i++) vary[irow[i]/bs]++;
178218c64b6SSatish Balay   count = 0;
179218c64b6SSatish Balay   for (i=0; i<a->mbs; i++) {
180e32f2f54SBarry Smith     if (vary[i]!=0 && vary[i]!=bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Index set does not match blocks");
181218c64b6SSatish Balay     if (vary[i]==bs) iary[count++] = i;
182218c64b6SSatish Balay   }
18370b3c8c7SBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,count,iary,PETSC_COPY_VALUES,&is1);CHKERRQ(ierr);
184218c64b6SSatish Balay 
185690b6cddSBarry Smith   ierr = PetscMemzero(vary,(a->mbs)*sizeof(PetscInt));CHKERRQ(ierr);
186218c64b6SSatish Balay   for (i=0; i<ncols; i++) vary[icol[i]/bs]++;
187218c64b6SSatish Balay   count = 0;
188218c64b6SSatish Balay   for (i=0; i<a->mbs; i++) {
189e32f2f54SBarry Smith     if (vary[i]!=0 && vary[i]!=bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Internal error in PETSc");
190218c64b6SSatish Balay     if (vary[i]==bs) iary[count++] = i;
191218c64b6SSatish Balay   }
19270b3c8c7SBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,count,iary,PETSC_COPY_VALUES,&is2);CHKERRQ(ierr);
193218c64b6SSatish Balay   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
194218c64b6SSatish Balay   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
195fca92195SBarry Smith   ierr = PetscFree2(vary,iary);CHKERRQ(ierr);
196218c64b6SSatish Balay 
1974aa3045dSJed Brown   ierr = MatGetSubMatrix_SeqBAIJ_Private(A,is1,is2,scall,B);CHKERRQ(ierr);
1986bf464f9SBarry Smith   ierr = ISDestroy(&is1);CHKERRQ(ierr);
1996bf464f9SBarry Smith   ierr = ISDestroy(&is2);CHKERRQ(ierr);
2003a40ed3dSBarry Smith   PetscFunctionReturn(0);
201218c64b6SSatish Balay }
202218c64b6SSatish Balay 
2034a2ae208SSatish Balay #undef __FUNCT__
2044a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_SeqBAIJ"
205690b6cddSBarry Smith PetscErrorCode MatGetSubMatrices_SeqBAIJ(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
206736121d4SSatish Balay {
2076849ba73SBarry Smith   PetscErrorCode ierr;
208690b6cddSBarry Smith   PetscInt       i;
209736121d4SSatish Balay 
2103a40ed3dSBarry Smith   PetscFunctionBegin;
211736121d4SSatish Balay   if (scall == MAT_INITIAL_MATRIX) {
21282502324SSatish Balay     ierr = PetscMalloc((n+1)*sizeof(Mat),B);CHKERRQ(ierr);
213736121d4SSatish Balay   }
214736121d4SSatish Balay 
215736121d4SSatish Balay   for (i=0; i<n; i++) {
2164aa3045dSJed Brown     ierr = MatGetSubMatrix_SeqBAIJ(A,irow[i],icol[i],scall,&(*B)[i]);CHKERRQ(ierr);
217736121d4SSatish Balay   }
2183a40ed3dSBarry Smith   PetscFunctionReturn(0);
219736121d4SSatish Balay }
220218c64b6SSatish Balay 
221218c64b6SSatish Balay 
2222d61bbb3SSatish Balay /* -------------------------------------------------------*/
2232d61bbb3SSatish Balay /* Should check that shapes of vectors and matrices match */
2242d61bbb3SSatish Balay /* -------------------------------------------------------*/
2252d61bbb3SSatish Balay 
2264a2ae208SSatish Balay #undef __FUNCT__
2274a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_1"
228dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_1(Mat A,Vec xx,Vec zz)
2292d61bbb3SSatish Balay {
2302d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
231d9fead3dSBarry Smith   PetscScalar       *z,sum;
232d9fead3dSBarry Smith   const PetscScalar *x;
233d9fead3dSBarry Smith   const MatScalar   *v;
2346849ba73SBarry Smith   PetscErrorCode    ierr;
2352162cab8SBarry Smith   PetscInt          mbs,i,n,nonzerorow=0;
2360298fd71SBarry Smith   const PetscInt    *idx,*ii,*ridx=NULL;
237ace3abfcSBarry Smith   PetscBool         usecprow=a->compressedrow.use;
2382d61bbb3SSatish Balay 
2392d61bbb3SSatish Balay   PetscFunctionBegin;
2403649974fSBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
2411ebc52fbSHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
2422d61bbb3SSatish Balay 
24326e093fcSHong Zhang   if (usecprow) {
24426e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
24526e093fcSHong Zhang     ii   = a->compressedrow.i;
2467b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
247a1c3900fSBarry Smith     ierr = PetscMemzero(z,mbs*sizeof(PetscScalar));CHKERRQ(ierr);
24826e093fcSHong Zhang   } else {
24926e093fcSHong Zhang     mbs = a->mbs;
2502d61bbb3SSatish Balay     ii  = a->i;
25126e093fcSHong Zhang   }
2522d61bbb3SSatish Balay 
2532d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
254ee54c7eeSHong Zhang     n   = ii[1] - ii[0];
255ee54c7eeSHong Zhang     v   = a->a + ii[0];
256ee54c7eeSHong Zhang     idx = a->j + ii[0];
257ee54c7eeSHong Zhang     ii++;
258444d8c10SJed Brown     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);   /* Indices for the next row (assumes same size as this one) */
259444d8c10SJed Brown     PetscPrefetchBlock(v+1*n,1*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
2602d61bbb3SSatish Balay     sum = 0.0;
2612162cab8SBarry Smith     PetscSparseDensePlusDot(sum,x,v,idx,n);
26226e093fcSHong Zhang     if (usecprow) {
2637b2bb3b9SHong Zhang       z[ridx[i]] = sum;
26426e093fcSHong Zhang     } else {
265122f12eaSBarry Smith       nonzerorow += (n>0);
2662d61bbb3SSatish Balay       z[i]        = sum;
2672d61bbb3SSatish Balay     }
26826e093fcSHong Zhang   }
2693649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
2701ebc52fbSHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
271dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*a->nz - nonzerorow);CHKERRQ(ierr);
2722d61bbb3SSatish Balay   PetscFunctionReturn(0);
2732d61bbb3SSatish Balay }
2742d61bbb3SSatish Balay 
2754a2ae208SSatish Balay #undef __FUNCT__
2764a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_2"
277dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_2(Mat A,Vec xx,Vec zz)
2782d61bbb3SSatish Balay {
2792d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
280d9fead3dSBarry Smith   PetscScalar       *z = 0,sum1,sum2,*zarray;
281d9fead3dSBarry Smith   const PetscScalar *x,*xb;
28287828ca2SBarry Smith   PetscScalar       x1,x2;
283d9fead3dSBarry Smith   const MatScalar   *v;
284dfbe8321SBarry Smith   PetscErrorCode    ierr;
2850298fd71SBarry Smith   PetscInt          mbs,i,*idx,*ii,j,n,*ridx=NULL,nonzerorow=0;
286ace3abfcSBarry Smith   PetscBool         usecprow=a->compressedrow.use;
2872d61bbb3SSatish Balay 
2882d61bbb3SSatish Balay   PetscFunctionBegin;
2893649974fSBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
29026e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
2912d61bbb3SSatish Balay 
2922d61bbb3SSatish Balay   idx = a->j;
2932d61bbb3SSatish Balay   v   = a->a;
29426e093fcSHong Zhang   if (usecprow) {
29526e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
29626e093fcSHong Zhang     ii   = a->compressedrow.i;
2977b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
29826e093fcSHong Zhang   } else {
29926e093fcSHong Zhang     mbs = a->mbs;
3002d61bbb3SSatish Balay     ii  = a->i;
30126e093fcSHong Zhang     z   = zarray;
30226e093fcSHong Zhang   }
3032d61bbb3SSatish Balay 
3042d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
3052d61bbb3SSatish Balay     n           = ii[1] - ii[0]; ii++;
3062d61bbb3SSatish Balay     sum1        = 0.0; sum2 = 0.0;
30798c9bda7SSatish Balay     nonzerorow += (n>0);
308444d8c10SJed Brown     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);   /* Indices for the next row (assumes same size as this one) */
309444d8c10SJed Brown     PetscPrefetchBlock(v+4*n,4*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
3102d61bbb3SSatish Balay     for (j=0; j<n; j++) {
3112d61bbb3SSatish Balay       xb    = x + 2*(*idx++); x1 = xb[0]; x2 = xb[1];
3122d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[2]*x2;
3132d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[3]*x2;
3142d61bbb3SSatish Balay       v    += 4;
3152d61bbb3SSatish Balay     }
3167b2bb3b9SHong Zhang     if (usecprow) z = zarray + 2*ridx[i];
3172d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2;
31826e093fcSHong Zhang     if (!usecprow) z += 2;
3192d61bbb3SSatish Balay   }
3203649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
32126e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
322dc0b31edSSatish Balay   ierr = PetscLogFlops(8.0*a->nz - 2.0*nonzerorow);CHKERRQ(ierr);
3232d61bbb3SSatish Balay   PetscFunctionReturn(0);
3242d61bbb3SSatish Balay }
3252d61bbb3SSatish Balay 
3264a2ae208SSatish Balay #undef __FUNCT__
3274a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_3"
328dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_3(Mat A,Vec xx,Vec zz)
3292d61bbb3SSatish Balay {
3302d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
331d9fead3dSBarry Smith   PetscScalar       *z = 0,sum1,sum2,sum3,x1,x2,x3,*zarray;
332d9fead3dSBarry Smith   const PetscScalar *x,*xb;
333d9fead3dSBarry Smith   const MatScalar   *v;
334dfbe8321SBarry Smith   PetscErrorCode    ierr;
3350298fd71SBarry Smith   PetscInt          mbs,i,*idx,*ii,j,n,*ridx=NULL,nonzerorow=0;
336ace3abfcSBarry Smith   PetscBool         usecprow=a->compressedrow.use;
33726e093fcSHong Zhang 
3382d61bbb3SSatish Balay 
339b6410449SSatish Balay #if defined(PETSC_HAVE_PRAGMA_DISJOINT)
340fee21e36SBarry Smith #pragma disjoint(*v,*z,*xb)
341fee21e36SBarry Smith #endif
342fee21e36SBarry Smith 
3432d61bbb3SSatish Balay   PetscFunctionBegin;
3443649974fSBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
34526e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
3462d61bbb3SSatish Balay 
3472d61bbb3SSatish Balay   idx = a->j;
3482d61bbb3SSatish Balay   v   = a->a;
34926e093fcSHong Zhang   if (usecprow) {
35026e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
35126e093fcSHong Zhang     ii   = a->compressedrow.i;
3527b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
35326e093fcSHong Zhang   } else {
35426e093fcSHong Zhang     mbs = a->mbs;
3552d61bbb3SSatish Balay     ii  = a->i;
35626e093fcSHong Zhang     z   = zarray;
35726e093fcSHong Zhang   }
3582d61bbb3SSatish Balay 
3592d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
3602d61bbb3SSatish Balay     n           = ii[1] - ii[0]; ii++;
3612d61bbb3SSatish Balay     sum1        = 0.0; sum2 = 0.0; sum3 = 0.0;
36298c9bda7SSatish Balay     nonzerorow += (n>0);
363444d8c10SJed Brown     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);   /* Indices for the next row (assumes same size as this one) */
364444d8c10SJed Brown     PetscPrefetchBlock(v+9*n,9*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
3652d61bbb3SSatish Balay     for (j=0; j<n; j++) {
36626fbe8dcSKarl Rupp       xb = x + 3*(*idx++);
36726fbe8dcSKarl Rupp       x1 = xb[0];
36826fbe8dcSKarl Rupp       x2 = xb[1];
36926fbe8dcSKarl Rupp       x3 = xb[2];
37026fbe8dcSKarl Rupp 
3712d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[3]*x2 + v[6]*x3;
3722d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[4]*x2 + v[7]*x3;
3732d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[5]*x2 + v[8]*x3;
3742d61bbb3SSatish Balay       v    += 9;
3752d61bbb3SSatish Balay     }
3767b2bb3b9SHong Zhang     if (usecprow) z = zarray + 3*ridx[i];
3772d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3;
37826e093fcSHong Zhang     if (!usecprow) z += 3;
3792d61bbb3SSatish Balay   }
3803649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
38126e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
382dc0b31edSSatish Balay   ierr = PetscLogFlops(18.0*a->nz - 3.0*nonzerorow);CHKERRQ(ierr);
3832d61bbb3SSatish Balay   PetscFunctionReturn(0);
3842d61bbb3SSatish Balay }
3852d61bbb3SSatish Balay 
3864a2ae208SSatish Balay #undef __FUNCT__
3874a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_4"
388dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_4(Mat A,Vec xx,Vec zz)
3892d61bbb3SSatish Balay {
3902d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
391d9fead3dSBarry Smith   PetscScalar       *z = 0,sum1,sum2,sum3,sum4,x1,x2,x3,x4,*zarray;
392d9fead3dSBarry Smith   const PetscScalar *x,*xb;
393d9fead3dSBarry Smith   const MatScalar   *v;
394dfbe8321SBarry Smith   PetscErrorCode    ierr;
3950298fd71SBarry Smith   PetscInt          mbs,i,*idx,*ii,j,n,*ridx=NULL,nonzerorow=0;
396ace3abfcSBarry Smith   PetscBool         usecprow=a->compressedrow.use;
3972d61bbb3SSatish Balay 
3982d61bbb3SSatish Balay   PetscFunctionBegin;
3993649974fSBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
40026e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
4012d61bbb3SSatish Balay 
4022d61bbb3SSatish Balay   idx = a->j;
4032d61bbb3SSatish Balay   v   = a->a;
40426e093fcSHong Zhang   if (usecprow) {
40526e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
40626e093fcSHong Zhang     ii   = a->compressedrow.i;
4077b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
40826e093fcSHong Zhang   } else {
40926e093fcSHong Zhang     mbs = a->mbs;
4102d61bbb3SSatish Balay     ii  = a->i;
41126e093fcSHong Zhang     z   = zarray;
41226e093fcSHong Zhang   }
4132d61bbb3SSatish Balay 
4142d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
41526fbe8dcSKarl Rupp     n = ii[1] - ii[0];
41626fbe8dcSKarl Rupp     ii++;
41726fbe8dcSKarl Rupp     sum1 = 0.0;
41826fbe8dcSKarl Rupp     sum2 = 0.0;
41926fbe8dcSKarl Rupp     sum3 = 0.0;
42026fbe8dcSKarl Rupp     sum4 = 0.0;
42126fbe8dcSKarl Rupp 
42298c9bda7SSatish Balay     nonzerorow += (n>0);
423444d8c10SJed Brown     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);     /* Indices for the next row (assumes same size as this one) */
424444d8c10SJed Brown     PetscPrefetchBlock(v+16*n,16*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
4252d61bbb3SSatish Balay     for (j=0; j<n; j++) {
4262d61bbb3SSatish Balay       xb    = x + 4*(*idx++);
4272d61bbb3SSatish Balay       x1    = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
4282d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[4]*x2 + v[8]*x3   + v[12]*x4;
4292d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[5]*x2 + v[9]*x3   + v[13]*x4;
4302d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[6]*x2 + v[10]*x3  + v[14]*x4;
4312d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[7]*x2 + v[11]*x3  + v[15]*x4;
4322d61bbb3SSatish Balay       v    += 16;
4332d61bbb3SSatish Balay     }
4347b2bb3b9SHong Zhang     if (usecprow) z = zarray + 4*ridx[i];
4352d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4;
43626e093fcSHong Zhang     if (!usecprow) z += 4;
4372d61bbb3SSatish Balay   }
4383649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
43926e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
440dc0b31edSSatish Balay   ierr = PetscLogFlops(32.0*a->nz - 4.0*nonzerorow);CHKERRQ(ierr);
4412d61bbb3SSatish Balay   PetscFunctionReturn(0);
4422d61bbb3SSatish Balay }
4432d61bbb3SSatish Balay 
4444a2ae208SSatish Balay #undef __FUNCT__
4454a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_5"
446dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_5(Mat A,Vec xx,Vec zz)
4472d61bbb3SSatish Balay {
4482d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
449d9fead3dSBarry Smith   PetscScalar       sum1,sum2,sum3,sum4,sum5,x1,x2,x3,x4,x5,*z = 0,*zarray;
450d9fead3dSBarry Smith   const PetscScalar *xb,*x;
451d9fead3dSBarry Smith   const MatScalar   *v;
452dfbe8321SBarry Smith   PetscErrorCode    ierr;
4530298fd71SBarry Smith   const PetscInt    *idx,*ii,*ridx=NULL;
454766f9fbaSBarry Smith   PetscInt          mbs,i,j,n,nonzerorow=0;
455ace3abfcSBarry Smith   PetscBool         usecprow=a->compressedrow.use;
4562d61bbb3SSatish Balay 
457433994e6SBarry Smith   PetscFunctionBegin;
4583649974fSBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
45926e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
4602d61bbb3SSatish Balay 
4612d61bbb3SSatish Balay   idx = a->j;
4622d61bbb3SSatish Balay   v   = a->a;
46326e093fcSHong Zhang   if (usecprow) {
46426e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
46526e093fcSHong Zhang     ii   = a->compressedrow.i;
4667b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
46726e093fcSHong Zhang   } else {
46826e093fcSHong Zhang     mbs = a->mbs;
4692d61bbb3SSatish Balay     ii  = a->i;
47026e093fcSHong Zhang     z   = zarray;
47126e093fcSHong Zhang   }
4722d61bbb3SSatish Balay 
4732d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
4742d61bbb3SSatish Balay     n           = ii[1] - ii[0]; ii++;
4752d61bbb3SSatish Balay     sum1        = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0;
47698c9bda7SSatish Balay     nonzerorow += (n>0);
477444d8c10SJed Brown     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);     /* Indices for the next row (assumes same size as this one) */
478444d8c10SJed Brown     PetscPrefetchBlock(v+25*n,25*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
4792d61bbb3SSatish Balay     for (j=0; j<n; j++) {
4802d61bbb3SSatish Balay       xb    = x + 5*(*idx++);
4812d61bbb3SSatish Balay       x1    = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4];
4822d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[5]*x2 + v[10]*x3  + v[15]*x4 + v[20]*x5;
4832d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[6]*x2 + v[11]*x3  + v[16]*x4 + v[21]*x5;
4842d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[7]*x2 + v[12]*x3  + v[17]*x4 + v[22]*x5;
4852d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[8]*x2 + v[13]*x3  + v[18]*x4 + v[23]*x5;
4862d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[9]*x2 + v[14]*x3  + v[19]*x4 + v[24]*x5;
4872d61bbb3SSatish Balay       v    += 25;
4882d61bbb3SSatish Balay     }
4897b2bb3b9SHong Zhang     if (usecprow) z = zarray + 5*ridx[i];
4902d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5;
49126e093fcSHong Zhang     if (!usecprow) z += 5;
4922d61bbb3SSatish Balay   }
4933649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
49426e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
495dc0b31edSSatish Balay   ierr = PetscLogFlops(50.0*a->nz - 5.0*nonzerorow);CHKERRQ(ierr);
4962d61bbb3SSatish Balay   PetscFunctionReturn(0);
4972d61bbb3SSatish Balay }
4982d61bbb3SSatish Balay 
49915091d37SBarry Smith 
5004a2ae208SSatish Balay #undef __FUNCT__
5014a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_6"
502dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_6(Mat A,Vec xx,Vec zz)
50315091d37SBarry Smith {
50415091d37SBarry Smith   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
505d9fead3dSBarry Smith   PetscScalar       *z = 0,sum1,sum2,sum3,sum4,sum5,sum6;
506d9fead3dSBarry Smith   const PetscScalar *x,*xb;
50726e093fcSHong Zhang   PetscScalar       x1,x2,x3,x4,x5,x6,*zarray;
508d9fead3dSBarry Smith   const MatScalar   *v;
509dfbe8321SBarry Smith   PetscErrorCode    ierr;
5100298fd71SBarry Smith   PetscInt          mbs,i,*idx,*ii,j,n,*ridx=NULL,nonzerorow=0;
511ace3abfcSBarry Smith   PetscBool         usecprow=a->compressedrow.use;
51215091d37SBarry Smith 
513433994e6SBarry Smith   PetscFunctionBegin;
5143649974fSBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
51526e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
51615091d37SBarry Smith 
51715091d37SBarry Smith   idx = a->j;
51815091d37SBarry Smith   v   = a->a;
51926e093fcSHong Zhang   if (usecprow) {
52026e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
52126e093fcSHong Zhang     ii   = a->compressedrow.i;
5227b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
52326e093fcSHong Zhang   } else {
52426e093fcSHong Zhang     mbs = a->mbs;
52515091d37SBarry Smith     ii  = a->i;
52626e093fcSHong Zhang     z   = zarray;
52726e093fcSHong Zhang   }
52815091d37SBarry Smith 
52915091d37SBarry Smith   for (i=0; i<mbs; i++) {
53026fbe8dcSKarl Rupp     n  = ii[1] - ii[0];
53126fbe8dcSKarl Rupp     ii++;
53226fbe8dcSKarl Rupp     sum1 = 0.0;
53326fbe8dcSKarl Rupp     sum2 = 0.0;
53426fbe8dcSKarl Rupp     sum3 = 0.0;
53526fbe8dcSKarl Rupp     sum4 = 0.0;
53626fbe8dcSKarl Rupp     sum5 = 0.0;
53726fbe8dcSKarl Rupp     sum6 = 0.0;
53826fbe8dcSKarl Rupp 
53998c9bda7SSatish Balay     nonzerorow += (n>0);
540444d8c10SJed Brown     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);     /* Indices for the next row (assumes same size as this one) */
541444d8c10SJed Brown     PetscPrefetchBlock(v+36*n,36*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
54215091d37SBarry Smith     for (j=0; j<n; j++) {
54315091d37SBarry Smith       xb    = x + 6*(*idx++);
54415091d37SBarry Smith       x1    = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5];
54515091d37SBarry Smith       sum1 += v[0]*x1 + v[6]*x2  + v[12]*x3  + v[18]*x4 + v[24]*x5 + v[30]*x6;
54615091d37SBarry Smith       sum2 += v[1]*x1 + v[7]*x2  + v[13]*x3  + v[19]*x4 + v[25]*x5 + v[31]*x6;
54715091d37SBarry Smith       sum3 += v[2]*x1 + v[8]*x2  + v[14]*x3  + v[20]*x4 + v[26]*x5 + v[32]*x6;
54815091d37SBarry Smith       sum4 += v[3]*x1 + v[9]*x2  + v[15]*x3  + v[21]*x4 + v[27]*x5 + v[33]*x6;
54915091d37SBarry Smith       sum5 += v[4]*x1 + v[10]*x2 + v[16]*x3  + v[22]*x4 + v[28]*x5 + v[34]*x6;
55015091d37SBarry Smith       sum6 += v[5]*x1 + v[11]*x2 + v[17]*x3  + v[23]*x4 + v[29]*x5 + v[35]*x6;
55115091d37SBarry Smith       v    += 36;
55215091d37SBarry Smith     }
5537b2bb3b9SHong Zhang     if (usecprow) z = zarray + 6*ridx[i];
55415091d37SBarry Smith     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6;
55526e093fcSHong Zhang     if (!usecprow) z += 6;
55615091d37SBarry Smith   }
55715091d37SBarry Smith 
5583649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
55926e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
560dc0b31edSSatish Balay   ierr = PetscLogFlops(72.0*a->nz - 6.0*nonzerorow);CHKERRQ(ierr);
56115091d37SBarry Smith   PetscFunctionReturn(0);
56215091d37SBarry Smith }
5638ab949d8SShri Abhyankar 
5644a2ae208SSatish Balay #undef __FUNCT__
5654a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_7"
566dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_7(Mat A,Vec xx,Vec zz)
5672d61bbb3SSatish Balay {
5682d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
569d9fead3dSBarry Smith   PetscScalar       *z = 0,sum1,sum2,sum3,sum4,sum5,sum6,sum7;
570d9fead3dSBarry Smith   const PetscScalar *x,*xb;
57126e093fcSHong Zhang   PetscScalar       x1,x2,x3,x4,x5,x6,x7,*zarray;
572d9fead3dSBarry Smith   const MatScalar   *v;
573dfbe8321SBarry Smith   PetscErrorCode    ierr;
5740298fd71SBarry Smith   PetscInt          mbs,i,*idx,*ii,j,n,*ridx=NULL,nonzerorow=0;
575ace3abfcSBarry Smith   PetscBool         usecprow=a->compressedrow.use;
5762d61bbb3SSatish Balay 
577433994e6SBarry Smith   PetscFunctionBegin;
5783649974fSBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
57926e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
5802d61bbb3SSatish Balay 
5812d61bbb3SSatish Balay   idx = a->j;
5822d61bbb3SSatish Balay   v   = a->a;
58326e093fcSHong Zhang   if (usecprow) {
58426e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
58526e093fcSHong Zhang     ii   = a->compressedrow.i;
5867b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
58726e093fcSHong Zhang   } else {
58826e093fcSHong Zhang     mbs = a->mbs;
5892d61bbb3SSatish Balay     ii  = a->i;
59026e093fcSHong Zhang     z   = zarray;
59126e093fcSHong Zhang   }
5922d61bbb3SSatish Balay 
5932d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
59426fbe8dcSKarl Rupp     n  = ii[1] - ii[0];
59526fbe8dcSKarl Rupp     ii++;
59626fbe8dcSKarl Rupp     sum1 = 0.0;
59726fbe8dcSKarl Rupp     sum2 = 0.0;
59826fbe8dcSKarl Rupp     sum3 = 0.0;
59926fbe8dcSKarl Rupp     sum4 = 0.0;
60026fbe8dcSKarl Rupp     sum5 = 0.0;
60126fbe8dcSKarl Rupp     sum6 = 0.0;
60226fbe8dcSKarl Rupp     sum7 = 0.0;
60326fbe8dcSKarl Rupp 
60498c9bda7SSatish Balay     nonzerorow += (n>0);
605444d8c10SJed Brown     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);     /* Indices for the next row (assumes same size as this one) */
606444d8c10SJed Brown     PetscPrefetchBlock(v+49*n,49*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
6072d61bbb3SSatish Balay     for (j=0; j<n; j++) {
6082d61bbb3SSatish Balay       xb    = x + 7*(*idx++);
6092d61bbb3SSatish Balay       x1    = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6];
6102d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[7]*x2  + v[14]*x3  + v[21]*x4 + v[28]*x5 + v[35]*x6 + v[42]*x7;
6112d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[8]*x2  + v[15]*x3  + v[22]*x4 + v[29]*x5 + v[36]*x6 + v[43]*x7;
6122d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[9]*x2  + v[16]*x3  + v[23]*x4 + v[30]*x5 + v[37]*x6 + v[44]*x7;
6132d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[10]*x2 + v[17]*x3  + v[24]*x4 + v[31]*x5 + v[38]*x6 + v[45]*x7;
6142d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[11]*x2 + v[18]*x3  + v[25]*x4 + v[32]*x5 + v[39]*x6 + v[46]*x7;
6152d61bbb3SSatish Balay       sum6 += v[5]*x1 + v[12]*x2 + v[19]*x3  + v[26]*x4 + v[33]*x5 + v[40]*x6 + v[47]*x7;
6162d61bbb3SSatish Balay       sum7 += v[6]*x1 + v[13]*x2 + v[20]*x3  + v[27]*x4 + v[34]*x5 + v[41]*x6 + v[48]*x7;
6172d61bbb3SSatish Balay       v    += 49;
6182d61bbb3SSatish Balay     }
6197b2bb3b9SHong Zhang     if (usecprow) z = zarray + 7*ridx[i];
6202d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
62126e093fcSHong Zhang     if (!usecprow) z += 7;
6222d61bbb3SSatish Balay   }
6232d61bbb3SSatish Balay 
6243649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
62526e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
626dc0b31edSSatish Balay   ierr = PetscLogFlops(98.0*a->nz - 7.0*nonzerorow);CHKERRQ(ierr);
6272d61bbb3SSatish Balay   PetscFunctionReturn(0);
6282d61bbb3SSatish Balay }
6292d61bbb3SSatish Balay 
6308ab949d8SShri Abhyankar /* MatMult_SeqBAIJ_15 version 1: Columns in the block are accessed one at a time */
631832cc040SShri Abhyankar /* Default MatMult for block size 15 */
6328ab949d8SShri Abhyankar 
6338ab949d8SShri Abhyankar #undef __FUNCT__
6348ab949d8SShri Abhyankar #define __FUNCT__ "MatMult_SeqBAIJ_15_ver1"
6358ab949d8SShri Abhyankar PetscErrorCode MatMult_SeqBAIJ_15_ver1(Mat A,Vec xx,Vec zz)
6368ab949d8SShri Abhyankar {
6378ab949d8SShri Abhyankar   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
6388ab949d8SShri Abhyankar   PetscScalar       *z = 0,sum1,sum2,sum3,sum4,sum5,sum6,sum7,sum8,sum9,sum10,sum11,sum12,sum13,sum14,sum15;
6398ab949d8SShri Abhyankar   const PetscScalar *x,*xb;
64053ef36baSBarry Smith   PetscScalar       *zarray,xv;
6418ab949d8SShri Abhyankar   const MatScalar   *v;
6428ab949d8SShri Abhyankar   PetscErrorCode    ierr;
6438ab949d8SShri Abhyankar   const PetscInt    *ii,*ij=a->j,*idx;
6440298fd71SBarry Smith   PetscInt          mbs,i,j,k,n,*ridx=NULL,nonzerorow=0;
645ace3abfcSBarry Smith   PetscBool         usecprow=a->compressedrow.use;
6468ab949d8SShri Abhyankar 
6478ab949d8SShri Abhyankar   PetscFunctionBegin;
6483649974fSBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
6498ab949d8SShri Abhyankar   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
6508ab949d8SShri Abhyankar 
6518ab949d8SShri Abhyankar   v = a->a;
6528ab949d8SShri Abhyankar   if (usecprow) {
6538ab949d8SShri Abhyankar     mbs  = a->compressedrow.nrows;
6548ab949d8SShri Abhyankar     ii   = a->compressedrow.i;
6558ab949d8SShri Abhyankar     ridx = a->compressedrow.rindex;
6568ab949d8SShri Abhyankar   } else {
6578ab949d8SShri Abhyankar     mbs = a->mbs;
6588ab949d8SShri Abhyankar     ii  = a->i;
6598ab949d8SShri Abhyankar     z   = zarray;
6608ab949d8SShri Abhyankar   }
6618ab949d8SShri Abhyankar 
6628ab949d8SShri Abhyankar   for (i=0; i<mbs; i++) {
6638ab949d8SShri Abhyankar     n    = ii[i+1] - ii[i];
6648ab949d8SShri Abhyankar     idx  = ij + ii[i];
6658ab949d8SShri Abhyankar     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0; sum7 = 0.0;
6668ab949d8SShri Abhyankar     sum8 = 0.0; sum9 = 0.0; sum10 = 0.0; sum11 = 0.0; sum12 = 0.0; sum13 = 0.0; sum14 = 0.0;sum15 = 0.0;
6678ab949d8SShri Abhyankar 
6688ab949d8SShri Abhyankar     nonzerorow += (n>0);
6698ab949d8SShri Abhyankar     for (j=0; j<n; j++) {
6708ab949d8SShri Abhyankar       xb = x + 15*(idx[j]);
6718ab949d8SShri Abhyankar 
6728ab949d8SShri Abhyankar       for (k=0; k<15; k++) {
67353ef36baSBarry Smith         xv     =  xb[k];
67453ef36baSBarry Smith         sum1  += v[0]*xv;
67553ef36baSBarry Smith         sum2  += v[1]*xv;
67653ef36baSBarry Smith         sum3  += v[2]*xv;
67753ef36baSBarry Smith         sum4  += v[3]*xv;
67853ef36baSBarry Smith         sum5  += v[4]*xv;
67953ef36baSBarry Smith         sum6  += v[5]*xv;
68053ef36baSBarry Smith         sum7  += v[6]*xv;
68153ef36baSBarry Smith         sum8  += v[7]*xv;
68253ef36baSBarry Smith         sum9  += v[8]*xv;
68353ef36baSBarry Smith         sum10 += v[9]*xv;
68453ef36baSBarry Smith         sum11 += v[10]*xv;
68553ef36baSBarry Smith         sum12 += v[11]*xv;
68653ef36baSBarry Smith         sum13 += v[12]*xv;
68753ef36baSBarry Smith         sum14 += v[13]*xv;
68853ef36baSBarry Smith         sum15 += v[14]*xv;
6898ab949d8SShri Abhyankar         v     += 15;
6908ab949d8SShri Abhyankar       }
6918ab949d8SShri Abhyankar     }
6928ab949d8SShri Abhyankar     if (usecprow) z = zarray + 15*ridx[i];
6938ab949d8SShri Abhyankar     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
6948ab949d8SShri Abhyankar     z[7] = sum8; z[8] = sum9; z[9] = sum10; z[10] = sum11; z[11] = sum12; z[12] = sum13; z[13] = sum14;z[14] = sum15;
6958ab949d8SShri Abhyankar 
6968ab949d8SShri Abhyankar     if (!usecprow) z += 15;
6978ab949d8SShri Abhyankar   }
6988ab949d8SShri Abhyankar 
6993649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
7008ab949d8SShri Abhyankar   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
7018ab949d8SShri Abhyankar   ierr = PetscLogFlops(450.0*a->nz - 15.0*nonzerorow);CHKERRQ(ierr);
7028ab949d8SShri Abhyankar   PetscFunctionReturn(0);
7038ab949d8SShri Abhyankar }
7048ab949d8SShri Abhyankar 
7058ab949d8SShri Abhyankar /* MatMult_SeqBAIJ_15_ver2 : Columns in the block are accessed in sets of 4,4,4,3 */
7068ab949d8SShri Abhyankar #undef __FUNCT__
7078ab949d8SShri Abhyankar #define __FUNCT__ "MatMult_SeqBAIJ_15_ver2"
7088ab949d8SShri Abhyankar PetscErrorCode MatMult_SeqBAIJ_15_ver2(Mat A,Vec xx,Vec zz)
7098ab949d8SShri Abhyankar {
7108ab949d8SShri Abhyankar   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
7118ab949d8SShri Abhyankar   PetscScalar       *z = 0,sum1,sum2,sum3,sum4,sum5,sum6,sum7,sum8,sum9,sum10,sum11,sum12,sum13,sum14,sum15;
7128ab949d8SShri Abhyankar   const PetscScalar *x,*xb;
7130b8f6341SShri Abhyankar   PetscScalar       x1,x2,x3,x4,*zarray;
7148ab949d8SShri Abhyankar   const MatScalar   *v;
7158ab949d8SShri Abhyankar   PetscErrorCode    ierr;
7168ab949d8SShri Abhyankar   const PetscInt    *ii,*ij=a->j,*idx;
7170298fd71SBarry Smith   PetscInt          mbs,i,j,n,*ridx=NULL,nonzerorow=0;
718ace3abfcSBarry Smith   PetscBool         usecprow=a->compressedrow.use;
7198ab949d8SShri Abhyankar 
7208ab949d8SShri Abhyankar   PetscFunctionBegin;
7213649974fSBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
7228ab949d8SShri Abhyankar   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
7238ab949d8SShri Abhyankar 
7248ab949d8SShri Abhyankar   v = a->a;
7258ab949d8SShri Abhyankar   if (usecprow) {
7268ab949d8SShri Abhyankar     mbs  = a->compressedrow.nrows;
7278ab949d8SShri Abhyankar     ii   = a->compressedrow.i;
7288ab949d8SShri Abhyankar     ridx = a->compressedrow.rindex;
7298ab949d8SShri Abhyankar   } else {
7308ab949d8SShri Abhyankar     mbs = a->mbs;
7318ab949d8SShri Abhyankar     ii  = a->i;
7328ab949d8SShri Abhyankar     z   = zarray;
7338ab949d8SShri Abhyankar   }
7348ab949d8SShri Abhyankar 
7358ab949d8SShri Abhyankar   for (i=0; i<mbs; i++) {
7368ab949d8SShri Abhyankar     n    = ii[i+1] - ii[i];
7378ab949d8SShri Abhyankar     idx  = ij + ii[i];
7388ab949d8SShri Abhyankar     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0; sum7 = 0.0;
7398ab949d8SShri Abhyankar     sum8 = 0.0; sum9 = 0.0; sum10 = 0.0; sum11 = 0.0; sum12 = 0.0; sum13 = 0.0; sum14 = 0.0;sum15 = 0.0;
7408ab949d8SShri Abhyankar 
7418ab949d8SShri Abhyankar     nonzerorow += (n>0);
7428ab949d8SShri Abhyankar     for (j=0; j<n; j++) {
7438ab949d8SShri Abhyankar       xb = x + 15*(idx[j]);
7440b8f6341SShri Abhyankar       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
7458ab949d8SShri Abhyankar 
7468ab949d8SShri Abhyankar       sum1  += v[0]*x1 + v[15]*x2 + v[30]*x3   + v[45]*x4;
7478ab949d8SShri Abhyankar       sum2  += v[1]*x1 + v[16]*x2 + v[31]*x3   + v[46]*x4;
7488ab949d8SShri Abhyankar       sum3  += v[2]*x1 + v[17]*x2 + v[32]*x3  + v[47]*x4;
7498ab949d8SShri Abhyankar       sum4  += v[3]*x1 + v[18]*x2 + v[33]*x3  + v[48]*x4;
7508ab949d8SShri Abhyankar       sum5  += v[4]*x1 + v[19]*x2 + v[34]*x3   + v[49]*x4;
7518ab949d8SShri Abhyankar       sum6  += v[5]*x1 + v[20]*x2 + v[35]*x3   + v[50]*x4;
7528ab949d8SShri Abhyankar       sum7  += v[6]*x1 + v[21]*x2 + v[36]*x3  + v[51]*x4;
7538ab949d8SShri Abhyankar       sum8  += v[7]*x1 + v[22]*x2 + v[37]*x3  + v[52]*x4;
7548ab949d8SShri Abhyankar       sum9  += v[8]*x1 + v[23]*x2 + v[38]*x3   + v[53]*x4;
7558ab949d8SShri Abhyankar       sum10 += v[9]*x1 + v[24]*x2 + v[39]*x3   + v[54]*x4;
7568ab949d8SShri Abhyankar       sum11 += v[10]*x1 + v[25]*x2 + v[40]*x3  + v[55]*x4;
7578ab949d8SShri Abhyankar       sum12 += v[11]*x1 + v[26]*x2 + v[41]*x3  + v[56]*x4;
7588ab949d8SShri Abhyankar       sum13 += v[12]*x1 + v[27]*x2 + v[42]*x3   + v[57]*x4;
7598ab949d8SShri Abhyankar       sum14 += v[13]*x1 + v[28]*x2 + v[43]*x3   + v[58]*x4;
7608ab949d8SShri Abhyankar       sum15 += v[14]*x1 + v[29]*x2 + v[44]*x3  + v[59]*x4;
7618ab949d8SShri Abhyankar 
7628ab949d8SShri Abhyankar       v += 60;
7638ab949d8SShri Abhyankar 
7640b8f6341SShri Abhyankar       x1 = xb[4]; x2 = xb[5]; x3 = xb[6]; x4 = xb[7];
7658ab949d8SShri Abhyankar 
7668ab949d8SShri Abhyankar       sum1  += v[0]*x1 + v[15]*x2 + v[30]*x3   + v[45]*x4;
7678ab949d8SShri Abhyankar       sum2  += v[1]*x1 + v[16]*x2 + v[31]*x3   + v[46]*x4;
7688ab949d8SShri Abhyankar       sum3  += v[2]*x1 + v[17]*x2 + v[32]*x3  + v[47]*x4;
7698ab949d8SShri Abhyankar       sum4  += v[3]*x1 + v[18]*x2 + v[33]*x3  + v[48]*x4;
7708ab949d8SShri Abhyankar       sum5  += v[4]*x1 + v[19]*x2 + v[34]*x3   + v[49]*x4;
7718ab949d8SShri Abhyankar       sum6  += v[5]*x1 + v[20]*x2 + v[35]*x3   + v[50]*x4;
7728ab949d8SShri Abhyankar       sum7  += v[6]*x1 + v[21]*x2 + v[36]*x3  + v[51]*x4;
7738ab949d8SShri Abhyankar       sum8  += v[7]*x1 + v[22]*x2 + v[37]*x3  + v[52]*x4;
7748ab949d8SShri Abhyankar       sum9  += v[8]*x1 + v[23]*x2 + v[38]*x3   + v[53]*x4;
7758ab949d8SShri Abhyankar       sum10 += v[9]*x1 + v[24]*x2 + v[39]*x3   + v[54]*x4;
7768ab949d8SShri Abhyankar       sum11 += v[10]*x1 + v[25]*x2 + v[40]*x3  + v[55]*x4;
7778ab949d8SShri Abhyankar       sum12 += v[11]*x1 + v[26]*x2 + v[41]*x3  + v[56]*x4;
7788ab949d8SShri Abhyankar       sum13 += v[12]*x1 + v[27]*x2 + v[42]*x3   + v[57]*x4;
7798ab949d8SShri Abhyankar       sum14 += v[13]*x1 + v[28]*x2 + v[43]*x3   + v[58]*x4;
7808ab949d8SShri Abhyankar       sum15 += v[14]*x1 + v[29]*x2 + v[44]*x3  + v[59]*x4;
7818ab949d8SShri Abhyankar       v     += 60;
7828ab949d8SShri Abhyankar 
7830b8f6341SShri Abhyankar       x1     = xb[8]; x2 = xb[9]; x3 = xb[10]; x4 = xb[11];
7840b8f6341SShri Abhyankar       sum1  += v[0]*x1 + v[15]*x2 + v[30]*x3   + v[45]*x4;
7850b8f6341SShri Abhyankar       sum2  += v[1]*x1 + v[16]*x2 + v[31]*x3   + v[46]*x4;
7860b8f6341SShri Abhyankar       sum3  += v[2]*x1 + v[17]*x2 + v[32]*x3  + v[47]*x4;
7870b8f6341SShri Abhyankar       sum4  += v[3]*x1 + v[18]*x2 + v[33]*x3  + v[48]*x4;
7880b8f6341SShri Abhyankar       sum5  += v[4]*x1 + v[19]*x2 + v[34]*x3   + v[49]*x4;
7890b8f6341SShri Abhyankar       sum6  += v[5]*x1 + v[20]*x2 + v[35]*x3   + v[50]*x4;
7900b8f6341SShri Abhyankar       sum7  += v[6]*x1 + v[21]*x2 + v[36]*x3  + v[51]*x4;
7910b8f6341SShri Abhyankar       sum8  += v[7]*x1 + v[22]*x2 + v[37]*x3  + v[52]*x4;
7920b8f6341SShri Abhyankar       sum9  += v[8]*x1 + v[23]*x2 + v[38]*x3   + v[53]*x4;
7930b8f6341SShri Abhyankar       sum10 += v[9]*x1 + v[24]*x2 + v[39]*x3   + v[54]*x4;
7940b8f6341SShri Abhyankar       sum11 += v[10]*x1 + v[25]*x2 + v[40]*x3  + v[55]*x4;
7950b8f6341SShri Abhyankar       sum12 += v[11]*x1 + v[26]*x2 + v[41]*x3  + v[56]*x4;
7960b8f6341SShri Abhyankar       sum13 += v[12]*x1 + v[27]*x2 + v[42]*x3   + v[57]*x4;
7970b8f6341SShri Abhyankar       sum14 += v[13]*x1 + v[28]*x2 + v[43]*x3   + v[58]*x4;
7980b8f6341SShri Abhyankar       sum15 += v[14]*x1 + v[29]*x2 + v[44]*x3  + v[59]*x4;
7990b8f6341SShri Abhyankar       v     += 60;
8000b8f6341SShri Abhyankar 
8010b8f6341SShri Abhyankar       x1     = xb[12]; x2 = xb[13]; x3 = xb[14];
8028ab949d8SShri Abhyankar       sum1  += v[0]*x1 + v[15]*x2 + v[30]*x3;
8038ab949d8SShri Abhyankar       sum2  += v[1]*x1 + v[16]*x2 + v[31]*x3;
8048ab949d8SShri Abhyankar       sum3  += v[2]*x1 + v[17]*x2 + v[32]*x3;
8058ab949d8SShri Abhyankar       sum4  += v[3]*x1 + v[18]*x2 + v[33]*x3;
8068ab949d8SShri Abhyankar       sum5  += v[4]*x1 + v[19]*x2 + v[34]*x3;
8078ab949d8SShri Abhyankar       sum6  += v[5]*x1 + v[20]*x2 + v[35]*x3;
8088ab949d8SShri Abhyankar       sum7  += v[6]*x1 + v[21]*x2 + v[36]*x3;
8098ab949d8SShri Abhyankar       sum8  += v[7]*x1 + v[22]*x2 + v[37]*x3;
8108ab949d8SShri Abhyankar       sum9  += v[8]*x1 + v[23]*x2 + v[38]*x3;
8118ab949d8SShri Abhyankar       sum10 += v[9]*x1 + v[24]*x2 + v[39]*x3;
8128ab949d8SShri Abhyankar       sum11 += v[10]*x1 + v[25]*x2 + v[40]*x3;
8138ab949d8SShri Abhyankar       sum12 += v[11]*x1 + v[26]*x2 + v[41]*x3;
8148ab949d8SShri Abhyankar       sum13 += v[12]*x1 + v[27]*x2 + v[42]*x3;
8158ab949d8SShri Abhyankar       sum14 += v[13]*x1 + v[28]*x2 + v[43]*x3;
8168ab949d8SShri Abhyankar       sum15 += v[14]*x1 + v[29]*x2 + v[44]*x3;
8178ab949d8SShri Abhyankar       v     += 45;
8188ab949d8SShri Abhyankar     }
8198ab949d8SShri Abhyankar     if (usecprow) z = zarray + 15*ridx[i];
8208ab949d8SShri Abhyankar     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
8218ab949d8SShri Abhyankar     z[7] = sum8; z[8] = sum9; z[9] = sum10; z[10] = sum11; z[11] = sum12; z[12] = sum13; z[13] = sum14;z[14] = sum15;
8228ab949d8SShri Abhyankar 
8238ab949d8SShri Abhyankar     if (!usecprow) z += 15;
8248ab949d8SShri Abhyankar   }
8258ab949d8SShri Abhyankar 
8263649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
8278ab949d8SShri Abhyankar   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
8288ab949d8SShri Abhyankar   ierr = PetscLogFlops(450.0*a->nz - 15.0*nonzerorow);CHKERRQ(ierr);
8298ab949d8SShri Abhyankar   PetscFunctionReturn(0);
8308ab949d8SShri Abhyankar }
8318ab949d8SShri Abhyankar 
8328ab949d8SShri Abhyankar /* MatMult_SeqBAIJ_15_ver3 : Columns in the block are accessed in sets of 8,7 */
8338ab949d8SShri Abhyankar #undef __FUNCT__
8348ab949d8SShri Abhyankar #define __FUNCT__ "MatMult_SeqBAIJ_15_ver3"
8358ab949d8SShri Abhyankar PetscErrorCode MatMult_SeqBAIJ_15_ver3(Mat A,Vec xx,Vec zz)
8368ab949d8SShri Abhyankar {
8378ab949d8SShri Abhyankar   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
8388ab949d8SShri Abhyankar   PetscScalar       *z = 0,sum1,sum2,sum3,sum4,sum5,sum6,sum7,sum8,sum9,sum10,sum11,sum12,sum13,sum14,sum15;
8398ab949d8SShri Abhyankar   const PetscScalar *x,*xb;
8400b8f6341SShri Abhyankar   PetscScalar       x1,x2,x3,x4,x5,x6,x7,x8,*zarray;
8418ab949d8SShri Abhyankar   const MatScalar   *v;
8428ab949d8SShri Abhyankar   PetscErrorCode    ierr;
8438ab949d8SShri Abhyankar   const PetscInt    *ii,*ij=a->j,*idx;
8440298fd71SBarry Smith   PetscInt          mbs,i,j,n,*ridx=NULL,nonzerorow=0;
845ace3abfcSBarry Smith   PetscBool         usecprow=a->compressedrow.use;
8468ab949d8SShri Abhyankar 
8478ab949d8SShri Abhyankar   PetscFunctionBegin;
8483649974fSBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
8498ab949d8SShri Abhyankar   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
8508ab949d8SShri Abhyankar 
8518ab949d8SShri Abhyankar   v = a->a;
8528ab949d8SShri Abhyankar   if (usecprow) {
8538ab949d8SShri Abhyankar     mbs  = a->compressedrow.nrows;
8548ab949d8SShri Abhyankar     ii   = a->compressedrow.i;
8558ab949d8SShri Abhyankar     ridx = a->compressedrow.rindex;
8568ab949d8SShri Abhyankar   } else {
8578ab949d8SShri Abhyankar     mbs = a->mbs;
8588ab949d8SShri Abhyankar     ii  = a->i;
8598ab949d8SShri Abhyankar     z   = zarray;
8608ab949d8SShri Abhyankar   }
8618ab949d8SShri Abhyankar 
8628ab949d8SShri Abhyankar   for (i=0; i<mbs; i++) {
8638ab949d8SShri Abhyankar     n    = ii[i+1] - ii[i];
8648ab949d8SShri Abhyankar     idx  = ij + ii[i];
8658ab949d8SShri Abhyankar     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0; sum7 = 0.0;
8668ab949d8SShri Abhyankar     sum8 = 0.0; sum9 = 0.0; sum10 = 0.0; sum11 = 0.0; sum12 = 0.0; sum13 = 0.0; sum14 = 0.0;sum15 = 0.0;
8678ab949d8SShri Abhyankar 
8688ab949d8SShri Abhyankar     nonzerorow += (n>0);
8698ab949d8SShri Abhyankar     for (j=0; j<n; j++) {
8708ab949d8SShri Abhyankar       xb = x + 15*(idx[j]);
8718ab949d8SShri Abhyankar       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6];
8720b8f6341SShri Abhyankar       x8 = xb[7];
8738ab949d8SShri Abhyankar 
8748ab949d8SShri Abhyankar       sum1  += v[0]*x1 + v[15]*x2  + v[30]*x3  + v[45]*x4 + v[60]*x5 + v[75]*x6 + v[90]*x7 + v[105]*x8;
8758ab949d8SShri Abhyankar       sum2  += v[1]*x1 + v[16]*x2  + v[31]*x3  + v[46]*x4 + v[61]*x5 + v[76]*x6 + v[91]*x7 + v[106]*x8;
8768ab949d8SShri Abhyankar       sum3  += v[2]*x1 + v[17]*x2  + v[32]*x3  + v[47]*x4 + v[62]*x5 + v[77]*x6 + v[92]*x7 + v[107]*x8;
8778ab949d8SShri Abhyankar       sum4  += v[3]*x1 + v[18]*x2 + v[33]*x3  + v[48]*x4 + v[63]*x5 + v[78]*x6 + v[93]*x7 + v[108]*x8;
8788ab949d8SShri Abhyankar       sum5  += v[4]*x1 + v[19]*x2 + v[34]*x3  + v[49]*x4 + v[64]*x5 + v[79]*x6 + v[94]*x7 + v[109]*x8;
8798ab949d8SShri Abhyankar       sum6  += v[5]*x1 + v[20]*x2 + v[35]*x3  + v[50]*x4 + v[65]*x5 + v[80]*x6 + v[95]*x7 + v[110]*x8;
8808ab949d8SShri Abhyankar       sum7  += v[6]*x1 + v[21]*x2 + v[36]*x3  + v[51]*x4 + v[66]*x5 + v[81]*x6 + v[96]*x7 + v[111]*x8;
8818ab949d8SShri Abhyankar       sum8  += v[7]*x1 + v[22]*x2  + v[37]*x3  + v[52]*x4 + v[67]*x5 + v[82]*x6 + v[97]*x7 + v[112]*x8;
8828ab949d8SShri Abhyankar       sum9  += v[8]*x1 + v[23]*x2  + v[38]*x3  + v[53]*x4 + v[68]*x5 + v[83]*x6 + v[98]*x7 + v[113]*x8;
8838ab949d8SShri Abhyankar       sum10 += v[9]*x1 + v[24]*x2  + v[39]*x3  + v[54]*x4 + v[69]*x5 + v[84]*x6 + v[99]*x7 + v[114]*x8;
8848ab949d8SShri Abhyankar       sum11 += v[10]*x1 + v[25]*x2 + v[40]*x3  + v[55]*x4 + v[70]*x5 + v[85]*x6 + v[100]*x7 + v[115]*x8;
8858ab949d8SShri Abhyankar       sum12 += v[11]*x1 + v[26]*x2 + v[41]*x3  + v[56]*x4 + v[71]*x5 + v[86]*x6 + v[101]*x7 + v[116]*x8;
8868ab949d8SShri Abhyankar       sum13 += v[12]*x1 + v[27]*x2 + v[42]*x3  + v[57]*x4 + v[72]*x5 + v[87]*x6 + v[102]*x7 + v[117]*x8;
8878ab949d8SShri Abhyankar       sum14 += v[13]*x1 + v[28]*x2 + v[43]*x3  + v[58]*x4 + v[73]*x5 + v[88]*x6 + v[103]*x7 + v[118]*x8;
8888ab949d8SShri Abhyankar       sum15 += v[14]*x1 + v[29]*x2 + v[44]*x3  + v[59]*x4 + v[74]*x5 + v[89]*x6 + v[104]*x7 + v[119]*x8;
8898ab949d8SShri Abhyankar       v     += 120;
8908ab949d8SShri Abhyankar 
8910b8f6341SShri Abhyankar       x1 = xb[8]; x2 = xb[9]; x3 = xb[10]; x4 = xb[11]; x5 = xb[12]; x6 = xb[13]; x7 = xb[14];
8920b8f6341SShri Abhyankar 
8938ab949d8SShri Abhyankar       sum1  += v[0]*x1 + v[15]*x2  + v[30]*x3  + v[45]*x4 + v[60]*x5 + v[75]*x6 + v[90]*x7;
8948ab949d8SShri Abhyankar       sum2  += v[1]*x1 + v[16]*x2  + v[31]*x3  + v[46]*x4 + v[61]*x5 + v[76]*x6 + v[91]*x7;
8958ab949d8SShri Abhyankar       sum3  += v[2]*x1 + v[17]*x2  + v[32]*x3  + v[47]*x4 + v[62]*x5 + v[77]*x6 + v[92]*x7;
8968ab949d8SShri Abhyankar       sum4  += v[3]*x1 + v[18]*x2 + v[33]*x3  + v[48]*x4 + v[63]*x5 + v[78]*x6 + v[93]*x7;
8978ab949d8SShri Abhyankar       sum5  += v[4]*x1 + v[19]*x2 + v[34]*x3  + v[49]*x4 + v[64]*x5 + v[79]*x6 + v[94]*x7;
8988ab949d8SShri Abhyankar       sum6  += v[5]*x1 + v[20]*x2 + v[35]*x3  + v[50]*x4 + v[65]*x5 + v[80]*x6 + v[95]*x7;
8998ab949d8SShri Abhyankar       sum7  += v[6]*x1 + v[21]*x2 + v[36]*x3  + v[51]*x4 + v[66]*x5 + v[81]*x6 + v[96]*x7;
9008ab949d8SShri Abhyankar       sum8  += v[7]*x1 + v[22]*x2  + v[37]*x3  + v[52]*x4 + v[67]*x5 + v[82]*x6 + v[97]*x7;
9018ab949d8SShri Abhyankar       sum9  += v[8]*x1 + v[23]*x2  + v[38]*x3  + v[53]*x4 + v[68]*x5 + v[83]*x6 + v[98]*x7;
9028ab949d8SShri Abhyankar       sum10 += v[9]*x1 + v[24]*x2  + v[39]*x3  + v[54]*x4 + v[69]*x5 + v[84]*x6 + v[99]*x7;
9038ab949d8SShri Abhyankar       sum11 += v[10]*x1 + v[25]*x2 + v[40]*x3  + v[55]*x4 + v[70]*x5 + v[85]*x6 + v[100]*x7;
9048ab949d8SShri Abhyankar       sum12 += v[11]*x1 + v[26]*x2 + v[41]*x3  + v[56]*x4 + v[71]*x5 + v[86]*x6 + v[101]*x7;
9058ab949d8SShri Abhyankar       sum13 += v[12]*x1 + v[27]*x2 + v[42]*x3  + v[57]*x4 + v[72]*x5 + v[87]*x6 + v[102]*x7;
9068ab949d8SShri Abhyankar       sum14 += v[13]*x1 + v[28]*x2 + v[43]*x3  + v[58]*x4 + v[73]*x5 + v[88]*x6 + v[103]*x7;
9078ab949d8SShri Abhyankar       sum15 += v[14]*x1 + v[29]*x2 + v[44]*x3  + v[59]*x4 + v[74]*x5 + v[89]*x6 + v[104]*x7;
9088ab949d8SShri Abhyankar       v     += 105;
9098ab949d8SShri Abhyankar     }
9108ab949d8SShri Abhyankar     if (usecprow) z = zarray + 15*ridx[i];
9118ab949d8SShri Abhyankar     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
9128ab949d8SShri Abhyankar     z[7] = sum8; z[8] = sum9; z[9] = sum10; z[10] = sum11; z[11] = sum12; z[12] = sum13; z[13] = sum14;z[14] = sum15;
9138ab949d8SShri Abhyankar 
9148ab949d8SShri Abhyankar     if (!usecprow) z += 15;
9158ab949d8SShri Abhyankar   }
9168ab949d8SShri Abhyankar 
9173649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
9188ab949d8SShri Abhyankar   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
9198ab949d8SShri Abhyankar   ierr = PetscLogFlops(450.0*a->nz - 15.0*nonzerorow);CHKERRQ(ierr);
9208ab949d8SShri Abhyankar   PetscFunctionReturn(0);
9218ab949d8SShri Abhyankar }
9228ab949d8SShri Abhyankar 
9238ab949d8SShri Abhyankar /* MatMult_SeqBAIJ_15_ver4 : All columns in the block are accessed at once */
9248ab949d8SShri Abhyankar 
9258ab949d8SShri Abhyankar #undef __FUNCT__
9268ab949d8SShri Abhyankar #define __FUNCT__ "MatMult_SeqBAIJ_15_ver4"
9278ab949d8SShri Abhyankar PetscErrorCode MatMult_SeqBAIJ_15_ver4(Mat A,Vec xx,Vec zz)
9288ab949d8SShri Abhyankar {
9298ab949d8SShri Abhyankar   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
9308ab949d8SShri Abhyankar   PetscScalar       *z = 0,sum1,sum2,sum3,sum4,sum5,sum6,sum7,sum8,sum9,sum10,sum11,sum12,sum13,sum14,sum15;
9318ab949d8SShri Abhyankar   const PetscScalar *x,*xb;
9328ab949d8SShri Abhyankar   PetscScalar       x1,x2,x3,x4,x5,x6,x7,x8,x9,x10,x11,x12,x13,x14,x15,*zarray;
9338ab949d8SShri Abhyankar   const MatScalar   *v;
9348ab949d8SShri Abhyankar   PetscErrorCode    ierr;
9358ab949d8SShri Abhyankar   const PetscInt    *ii,*ij=a->j,*idx;
9360298fd71SBarry Smith   PetscInt          mbs,i,j,n,*ridx=NULL,nonzerorow=0;
937ace3abfcSBarry Smith   PetscBool         usecprow=a->compressedrow.use;
9388ab949d8SShri Abhyankar 
9398ab949d8SShri Abhyankar   PetscFunctionBegin;
9403649974fSBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
9418ab949d8SShri Abhyankar   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
9428ab949d8SShri Abhyankar 
9438ab949d8SShri Abhyankar   v = a->a;
9448ab949d8SShri Abhyankar   if (usecprow) {
9458ab949d8SShri Abhyankar     mbs  = a->compressedrow.nrows;
9468ab949d8SShri Abhyankar     ii   = a->compressedrow.i;
9478ab949d8SShri Abhyankar     ridx = a->compressedrow.rindex;
9488ab949d8SShri Abhyankar   } else {
9498ab949d8SShri Abhyankar     mbs = a->mbs;
9508ab949d8SShri Abhyankar     ii  = a->i;
9518ab949d8SShri Abhyankar     z   = zarray;
9528ab949d8SShri Abhyankar   }
9538ab949d8SShri Abhyankar 
9548ab949d8SShri Abhyankar   for (i=0; i<mbs; i++) {
9558ab949d8SShri Abhyankar     n    = ii[i+1] - ii[i];
9568ab949d8SShri Abhyankar     idx  = ij + ii[i];
9578ab949d8SShri Abhyankar     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0; sum7 = 0.0;
9588ab949d8SShri Abhyankar     sum8 = 0.0; sum9 = 0.0; sum10 = 0.0; sum11 = 0.0; sum12 = 0.0; sum13 = 0.0; sum14 = 0.0;sum15 = 0.0;
9598ab949d8SShri Abhyankar 
9608ab949d8SShri Abhyankar     nonzerorow += (n>0);
9618ab949d8SShri Abhyankar     for (j=0; j<n; j++) {
9628ab949d8SShri Abhyankar       xb = x + 15*(idx[j]);
9638ab949d8SShri Abhyankar       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6];
9648ab949d8SShri Abhyankar       x8 = xb[7]; x9 = xb[8]; x10 = xb[9]; x11 = xb[10]; x12 = xb[11]; x13 = xb[12]; x14 = xb[13];x15 = xb[14];
9658ab949d8SShri Abhyankar 
9668ab949d8SShri Abhyankar       sum1  +=  v[0]*x1  + v[15]*x2 + v[30]*x3 + v[45]*x4 + v[60]*x5 + v[75]*x6 + v[90]*x7  + v[105]*x8 + v[120]*x9 + v[135]*x10 + v[150]*x11 + v[165]*x12 + v[180]*x13 + v[195]*x14 + v[210]*x15;
9678ab949d8SShri Abhyankar       sum2  +=  v[1]*x1  + v[16]*x2 + v[31]*x3 + v[46]*x4 + v[61]*x5 + v[76]*x6 + v[91]*x7  + v[106]*x8 + v[121]*x9 + v[136]*x10 + v[151]*x11 + v[166]*x12 + v[181]*x13 + v[196]*x14 + v[211]*x15;
9688ab949d8SShri Abhyankar       sum3  +=  v[2]*x1  + v[17]*x2 + v[32]*x3 + v[47]*x4 + v[62]*x5 + v[77]*x6 + v[92]*x7  + v[107]*x8 + v[122]*x9 + v[137]*x10 + v[152]*x11 + v[167]*x12 + v[182]*x13 + v[197]*x14 + v[212]*x15;
9698ab949d8SShri Abhyankar       sum4  +=  v[3]*x1  + v[18]*x2 + v[33]*x3 + v[48]*x4 + v[63]*x5 + v[78]*x6 + v[93]*x7  + v[108]*x8 + v[123]*x9 + v[138]*x10 + v[153]*x11 + v[168]*x12 + v[183]*x13 + v[198]*x14 + v[213]*x15;
9708ab949d8SShri Abhyankar       sum5  += v[4]*x1  + v[19]*x2 + v[34]*x3 + v[49]*x4 + v[64]*x5 + v[79]*x6 + v[94]*x7  + v[109]*x8 + v[124]*x9 + v[139]*x10 + v[154]*x11 + v[169]*x12 + v[184]*x13 + v[199]*x14 + v[214]*x15;
9718ab949d8SShri Abhyankar       sum6  += v[5]*x1  + v[20]*x2 + v[35]*x3 + v[50]*x4 + v[65]*x5 + v[80]*x6 + v[95]*x7  + v[110]*x8 + v[125]*x9 + v[140]*x10 + v[155]*x11 + v[170]*x12 + v[185]*x13 + v[200]*x14 + v[215]*x15;
9728ab949d8SShri Abhyankar       sum7  += v[6]*x1  + v[21]*x2 + v[36]*x3 + v[51]*x4 + v[66]*x5 + v[81]*x6 + v[96]*x7  + v[111]*x8 + v[126]*x9 + v[141]*x10 + v[156]*x11 + v[171]*x12 + v[186]*x13 + v[201]*x14 + v[216]*x15;
9738ab949d8SShri Abhyankar       sum8  += v[7]*x1  + v[22]*x2 + v[37]*x3 + v[52]*x4 + v[67]*x5 + v[82]*x6 + v[97]*x7  + v[112]*x8 + v[127]*x9 + v[142]*x10 + v[157]*x11 + v[172]*x12 + v[187]*x13 + v[202]*x14 + v[217]*x15;
9748ab949d8SShri Abhyankar       sum9  += v[8]*x1  + v[23]*x2 + v[38]*x3 + v[53]*x4 + v[68]*x5 + v[83]*x6 + v[98]*x7  + v[113]*x8 + v[128]*x9 + v[143]*x10 + v[158]*x11 + v[173]*x12 + v[188]*x13 + v[203]*x14 + v[218]*x15;
9758ab949d8SShri Abhyankar       sum10 += v[9]*x1  + v[24]*x2 + v[39]*x3 + v[54]*x4 + v[69]*x5 + v[84]*x6 + v[99]*x7  + v[114]*x8 + v[129]*x9 + v[144]*x10 + v[159]*x11 + v[174]*x12 + v[189]*x13 + v[204]*x14 + v[219]*x15;
9768ab949d8SShri Abhyankar       sum11 += v[10]*x1 + v[25]*x2 + v[40]*x3 + v[55]*x4 + v[70]*x5 + v[85]*x6 + v[100]*x7 + v[115]*x8 + v[130]*x9 + v[145]*x10 + v[160]*x11 + v[175]*x12 + v[190]*x13 + v[205]*x14 + v[220]*x15;
9778ab949d8SShri Abhyankar       sum12 += v[11]*x1 + v[26]*x2 + v[41]*x3 + v[56]*x4 + v[71]*x5 + v[86]*x6 + v[101]*x7 + v[116]*x8 + v[131]*x9 + v[146]*x10 + v[161]*x11 + v[176]*x12 + v[191]*x13 + v[206]*x14 + v[221]*x15;
9788ab949d8SShri Abhyankar       sum13 += v[12]*x1 + v[27]*x2 + v[42]*x3 + v[57]*x4 + v[72]*x5 + v[87]*x6 + v[102]*x7 + v[117]*x8 + v[132]*x9 + v[147]*x10 + v[162]*x11 + v[177]*x12 + v[192]*x13 + v[207]*x14 + v[222]*x15;
9798ab949d8SShri Abhyankar       sum14 += v[13]*x1 + v[28]*x2 + v[43]*x3 + v[58]*x4 + v[73]*x5 + v[88]*x6 + v[103]*x7 + v[118]*x8 + v[133]*x9 + v[148]*x10 + v[163]*x11 + v[178]*x12 + v[193]*x13 + v[208]*x14 + v[223]*x15;
9808ab949d8SShri Abhyankar       sum15 += v[14]*x1 + v[29]*x2 + v[44]*x3 + v[59]*x4 + v[74]*x5 + v[89]*x6 + v[104]*x7 + v[119]*x8 + v[134]*x9 + v[149]*x10 + v[164]*x11 + v[179]*x12 + v[194]*x13 + v[209]*x14 + v[224]*x15;
9818ab949d8SShri Abhyankar       v     += 225;
9828ab949d8SShri Abhyankar     }
9838ab949d8SShri Abhyankar     if (usecprow) z = zarray + 15*ridx[i];
9848ab949d8SShri Abhyankar     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
9858ab949d8SShri Abhyankar     z[7] = sum8; z[8] = sum9; z[9] = sum10; z[10] = sum11; z[11] = sum12; z[12] = sum13; z[13] = sum14;z[14] = sum15;
9868ab949d8SShri Abhyankar 
9878ab949d8SShri Abhyankar     if (!usecprow) z += 15;
9888ab949d8SShri Abhyankar   }
9898ab949d8SShri Abhyankar 
9903649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
9918ab949d8SShri Abhyankar   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
9928ab949d8SShri Abhyankar   ierr = PetscLogFlops(450.0*a->nz - 15.0*nonzerorow);CHKERRQ(ierr);
9938ab949d8SShri Abhyankar   PetscFunctionReturn(0);
9948ab949d8SShri Abhyankar }
9958ab949d8SShri Abhyankar 
9968ab949d8SShri Abhyankar 
9973f1db9ecSBarry Smith /*
9983f1db9ecSBarry Smith     This will not work with MatScalar == float because it calls the BLAS
9993f1db9ecSBarry Smith */
10004a2ae208SSatish Balay #undef __FUNCT__
10014a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_N"
1002dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_N(Mat A,Vec xx,Vec zz)
10032d61bbb3SSatish Balay {
10042d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
1005dcf5cc72SBarry Smith   PetscScalar    *x,*z = 0,*xb,*work,*workt,*zarray;
10063f1db9ecSBarry Smith   MatScalar      *v;
1007dfbe8321SBarry Smith   PetscErrorCode ierr;
1008a2ea699eSBarry Smith   PetscInt       mbs,i,*idx,*ii,bs=A->rmap->bs,j,n,bs2=a->bs2;
10090298fd71SBarry Smith   PetscInt       ncols,k,*ridx=NULL,nonzerorow=0;
1010ace3abfcSBarry Smith   PetscBool      usecprow=a->compressedrow.use;
10112d61bbb3SSatish Balay 
10122d61bbb3SSatish Balay   PetscFunctionBegin;
10131ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
101426e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
10152d61bbb3SSatish Balay 
10162d61bbb3SSatish Balay   idx = a->j;
10172d61bbb3SSatish Balay   v   = a->a;
101826e093fcSHong Zhang   if (usecprow) {
101926e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
102026e093fcSHong Zhang     ii   = a->compressedrow.i;
10217b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
102226e093fcSHong Zhang   } else {
102326e093fcSHong Zhang     mbs = a->mbs;
10242d61bbb3SSatish Balay     ii  = a->i;
102526e093fcSHong Zhang     z   = zarray;
102626e093fcSHong Zhang   }
1027218c64b6SSatish Balay 
10282d61bbb3SSatish Balay   if (!a->mult_work) {
1029d0f46423SBarry Smith     k    = PetscMax(A->rmap->n,A->cmap->n);
103087828ca2SBarry Smith     ierr = PetscMalloc((k+1)*sizeof(PetscScalar),&a->mult_work);CHKERRQ(ierr);
10312d61bbb3SSatish Balay   }
10322d61bbb3SSatish Balay   work = a->mult_work;
10332d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
10342d61bbb3SSatish Balay     n           = ii[1] - ii[0]; ii++;
10352d61bbb3SSatish Balay     ncols       = n*bs;
10362d61bbb3SSatish Balay     workt       = work;
103798c9bda7SSatish Balay     nonzerorow += (n>0);
10382d61bbb3SSatish Balay     for (j=0; j<n; j++) {
10392d61bbb3SSatish Balay       xb = x + bs*(*idx++);
10402d61bbb3SSatish Balay       for (k=0; k<bs; k++) workt[k] = xb[k];
10412d61bbb3SSatish Balay       workt += bs;
10422d61bbb3SSatish Balay     }
10437b2bb3b9SHong Zhang     if (usecprow) z = zarray + bs*ridx[i];
104496b95a6bSBarry Smith     PetscKernel_w_gets_Ar_times_v(bs,ncols,work,v,z);
104571044d3cSBarry Smith     /* BLASgemv_("N",&bs,&ncols,&_DOne,v,&bs,work,&_One,&_DZero,z,&_One); */
10462d61bbb3SSatish Balay     v += n*bs2;
104726e093fcSHong Zhang     if (!usecprow) z += bs;
10482d61bbb3SSatish Balay   }
10491ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
105026e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
1051dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*a->nz*bs2 - bs*nonzerorow);CHKERRQ(ierr);
10522d61bbb3SSatish Balay   PetscFunctionReturn(0);
10532d61bbb3SSatish Balay }
10542d61bbb3SSatish Balay 
10554a2ae208SSatish Balay #undef __FUNCT__
10564a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_1"
1057dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_1(Mat A,Vec xx,Vec yy,Vec zz)
10582d61bbb3SSatish Balay {
10592d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
1060122f12eaSBarry Smith   const PetscScalar *x;
1061122f12eaSBarry Smith   PetscScalar       *y,*z,sum;
1062122f12eaSBarry Smith   const MatScalar   *v;
1063dfbe8321SBarry Smith   PetscErrorCode    ierr;
10640298fd71SBarry Smith   PetscInt          mbs=a->mbs,i,n,*ridx=NULL,nonzerorow=0;
1065122f12eaSBarry Smith   const PetscInt    *idx,*ii;
1066ace3abfcSBarry Smith   PetscBool         usecprow=a->compressedrow.use;
10672d61bbb3SSatish Balay 
10682d61bbb3SSatish Balay   PetscFunctionBegin;
10693649974fSBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
1070122f12eaSBarry Smith   ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
10712e8a6d31SBarry Smith   if (zz != yy) {
1072122f12eaSBarry Smith     ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
10732e8a6d31SBarry Smith   } else {
1074122f12eaSBarry Smith     z = y;
10752e8a6d31SBarry Smith   }
10762d61bbb3SSatish Balay 
10772d61bbb3SSatish Balay   idx = a->j;
10782d61bbb3SSatish Balay   v   = a->a;
107926e093fcSHong Zhang   if (usecprow) {
10804eb6d288SHong Zhang     if (zz != yy) {
1081122f12eaSBarry Smith       ierr = PetscMemcpy(z,y,mbs*sizeof(PetscScalar));CHKERRQ(ierr);
10824eb6d288SHong Zhang     }
108326e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
108426e093fcSHong Zhang     ii   = a->compressedrow.i;
10857b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
108626e093fcSHong Zhang   } else {
10872d61bbb3SSatish Balay     ii = a->i;
108826e093fcSHong Zhang   }
10892d61bbb3SSatish Balay 
10902d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
1091122f12eaSBarry Smith     n = ii[1] - ii[0];
1092122f12eaSBarry Smith     ii++;
109326e093fcSHong Zhang     if (!usecprow) {
1094122f12eaSBarry Smith       nonzerorow += (n>0);
1095122f12eaSBarry Smith       sum         = y[i];
1096122f12eaSBarry Smith     } else {
1097122f12eaSBarry Smith       sum = y[ridx[i]];
1098122f12eaSBarry Smith     }
1099444d8c10SJed Brown     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA); /* Indices for the next row (assumes same size as this one) */
1100444d8c10SJed Brown     PetscPrefetchBlock(v+n,n,0,PETSC_PREFETCH_HINT_NTA);   /* Entries for the next row */
1101122f12eaSBarry Smith     PetscSparseDensePlusDot(sum,x,v,idx,n);
1102122f12eaSBarry Smith     v   += n;
1103122f12eaSBarry Smith     idx += n;
1104122f12eaSBarry Smith     if (usecprow) {
1105122f12eaSBarry Smith       z[ridx[i]] = sum;
1106122f12eaSBarry Smith     } else {
1107122f12eaSBarry Smith       z[i] = sum;
110826e093fcSHong Zhang     }
11092d61bbb3SSatish Balay   }
11103649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
1111122f12eaSBarry Smith   ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
11122e8a6d31SBarry Smith   if (zz != yy) {
1113122f12eaSBarry Smith     ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
11142e8a6d31SBarry Smith   }
1115122f12eaSBarry Smith   ierr = PetscLogFlops(2.0*a->nz - nonzerorow);CHKERRQ(ierr);
11162d61bbb3SSatish Balay   PetscFunctionReturn(0);
11172d61bbb3SSatish Balay }
11182d61bbb3SSatish Balay 
11194a2ae208SSatish Balay #undef __FUNCT__
11204a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_2"
1121dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_2(Mat A,Vec xx,Vec yy,Vec zz)
11222d61bbb3SSatish Balay {
11232d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
1124dcf5cc72SBarry Smith   PetscScalar    *x,*y = 0,*z = 0,*xb,sum1,sum2;
112526e093fcSHong Zhang   PetscScalar    x1,x2,*yarray,*zarray;
11263f1db9ecSBarry Smith   MatScalar      *v;
1127dfbe8321SBarry Smith   PetscErrorCode ierr;
11280298fd71SBarry Smith   PetscInt       mbs     =a->mbs,i,*idx,*ii,j,n,*ridx=NULL;
1129ace3abfcSBarry Smith   PetscBool      usecprow=a->compressedrow.use;
11302d61bbb3SSatish Balay 
11312d61bbb3SSatish Balay   PetscFunctionBegin;
11321ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
113326e093fcSHong Zhang   ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr);
11342e8a6d31SBarry Smith   if (zz != yy) {
113526e093fcSHong Zhang     ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
11362e8a6d31SBarry Smith   } else {
113726e093fcSHong Zhang     zarray = yarray;
11382e8a6d31SBarry Smith   }
11392d61bbb3SSatish Balay 
11402d61bbb3SSatish Balay   idx = a->j;
11412d61bbb3SSatish Balay   v   = a->a;
114226e093fcSHong Zhang   if (usecprow) {
11434eb6d288SHong Zhang     if (zz != yy) {
11444eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,2*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
11454eb6d288SHong Zhang     }
114626e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
114726e093fcSHong Zhang     ii   = a->compressedrow.i;
11487b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
11494eb6d288SHong Zhang     if (zz != yy) {
11504eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,a->mbs*sizeof(PetscScalar));CHKERRQ(ierr);
11514eb6d288SHong Zhang     }
115226e093fcSHong Zhang   } else {
11532d61bbb3SSatish Balay     ii = a->i;
115426e093fcSHong Zhang     y  = yarray;
115526e093fcSHong Zhang     z  = zarray;
115626e093fcSHong Zhang   }
11572d61bbb3SSatish Balay 
11582d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
11592d61bbb3SSatish Balay     n = ii[1] - ii[0]; ii++;
116026e093fcSHong Zhang     if (usecprow) {
11617b2bb3b9SHong Zhang       z = zarray + 2*ridx[i];
11627b2bb3b9SHong Zhang       y = yarray + 2*ridx[i];
116326e093fcSHong Zhang     }
11642d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1];
1165444d8c10SJed Brown     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);   /* Indices for the next row (assumes same size as this one) */
1166444d8c10SJed Brown     PetscPrefetchBlock(v+4*n,4*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
11672d61bbb3SSatish Balay     for (j=0; j<n; j++) {
116826fbe8dcSKarl Rupp       xb = x + 2*(*idx++);
116926fbe8dcSKarl Rupp       x1 = xb[0];
117026fbe8dcSKarl Rupp       x2 = xb[1];
117126fbe8dcSKarl Rupp 
11722d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[2]*x2;
11732d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[3]*x2;
11742d61bbb3SSatish Balay       v    += 4;
11752d61bbb3SSatish Balay     }
11762d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2;
117726e093fcSHong Zhang     if (!usecprow) {
11782d61bbb3SSatish Balay       z += 2; y += 2;
11792d61bbb3SSatish Balay     }
118026e093fcSHong Zhang   }
11811ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
118226e093fcSHong Zhang   ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr);
11832e8a6d31SBarry Smith   if (zz != yy) {
118426e093fcSHong Zhang     ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
11852e8a6d31SBarry Smith   }
1186dc0b31edSSatish Balay   ierr = PetscLogFlops(4.0*a->nz);CHKERRQ(ierr);
11872d61bbb3SSatish Balay   PetscFunctionReturn(0);
11882d61bbb3SSatish Balay }
11892d61bbb3SSatish Balay 
11904a2ae208SSatish Balay #undef __FUNCT__
11914a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_3"
1192dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_3(Mat A,Vec xx,Vec yy,Vec zz)
11932d61bbb3SSatish Balay {
11942d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
1195dcf5cc72SBarry Smith   PetscScalar    *x,*y = 0,*z = 0,*xb,sum1,sum2,sum3,x1,x2,x3,*yarray,*zarray;
11963f1db9ecSBarry Smith   MatScalar      *v;
1197dfbe8321SBarry Smith   PetscErrorCode ierr;
11980298fd71SBarry Smith   PetscInt       mbs     =a->mbs,i,*idx,*ii,j,n,*ridx=NULL;
1199ace3abfcSBarry Smith   PetscBool      usecprow=a->compressedrow.use;
12002d61bbb3SSatish Balay 
12012d61bbb3SSatish Balay   PetscFunctionBegin;
12021ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
120326e093fcSHong Zhang   ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr);
12042e8a6d31SBarry Smith   if (zz != yy) {
120526e093fcSHong Zhang     ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
12062e8a6d31SBarry Smith   } else {
120726e093fcSHong Zhang     zarray = yarray;
12082e8a6d31SBarry Smith   }
12092d61bbb3SSatish Balay 
12102d61bbb3SSatish Balay   idx = a->j;
12112d61bbb3SSatish Balay   v   = a->a;
121226e093fcSHong Zhang   if (usecprow) {
12134eb6d288SHong Zhang     if (zz != yy) {
12144eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,3*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
12154eb6d288SHong Zhang     }
121626e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
121726e093fcSHong Zhang     ii   = a->compressedrow.i;
12187b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
121926e093fcSHong Zhang   } else {
12202d61bbb3SSatish Balay     ii = a->i;
122126e093fcSHong Zhang     y  = yarray;
122226e093fcSHong Zhang     z  = zarray;
122326e093fcSHong Zhang   }
12242d61bbb3SSatish Balay 
12252d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
12262d61bbb3SSatish Balay     n = ii[1] - ii[0]; ii++;
122726e093fcSHong Zhang     if (usecprow) {
12287b2bb3b9SHong Zhang       z = zarray + 3*ridx[i];
12297b2bb3b9SHong Zhang       y = yarray + 3*ridx[i];
123026e093fcSHong Zhang     }
12312d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2];
1232444d8c10SJed Brown     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);   /* Indices for the next row (assumes same size as this one) */
1233444d8c10SJed Brown     PetscPrefetchBlock(v+9*n,9*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
12342d61bbb3SSatish Balay     for (j=0; j<n; j++) {
12352d61bbb3SSatish Balay       xb    = x + 3*(*idx++); x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
12362d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[3]*x2 + v[6]*x3;
12372d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[4]*x2 + v[7]*x3;
12382d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[5]*x2 + v[8]*x3;
12392d61bbb3SSatish Balay       v    += 9;
12402d61bbb3SSatish Balay     }
12412d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3;
124226e093fcSHong Zhang     if (!usecprow) {
12432d61bbb3SSatish Balay       z += 3; y += 3;
12442d61bbb3SSatish Balay     }
124526e093fcSHong Zhang   }
12461ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
124726e093fcSHong Zhang   ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr);
12482e8a6d31SBarry Smith   if (zz != yy) {
124926e093fcSHong Zhang     ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
12502e8a6d31SBarry Smith   }
1251dc0b31edSSatish Balay   ierr = PetscLogFlops(18.0*a->nz);CHKERRQ(ierr);
12522d61bbb3SSatish Balay   PetscFunctionReturn(0);
12532d61bbb3SSatish Balay }
12542d61bbb3SSatish Balay 
12554a2ae208SSatish Balay #undef __FUNCT__
12564a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_4"
1257dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_4(Mat A,Vec xx,Vec yy,Vec zz)
12582d61bbb3SSatish Balay {
12592d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
1260dcf5cc72SBarry Smith   PetscScalar    *x,*y = 0,*z = 0,*xb,sum1,sum2,sum3,sum4,x1,x2,x3,x4,*yarray,*zarray;
12613f1db9ecSBarry Smith   MatScalar      *v;
1262dfbe8321SBarry Smith   PetscErrorCode ierr;
12630298fd71SBarry Smith   PetscInt       mbs     =a->mbs,i,*idx,*ii,j,n,*ridx=NULL;
1264ace3abfcSBarry Smith   PetscBool      usecprow=a->compressedrow.use;
12652d61bbb3SSatish Balay 
12662d61bbb3SSatish Balay   PetscFunctionBegin;
12671ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
126826e093fcSHong Zhang   ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr);
12692e8a6d31SBarry Smith   if (zz != yy) {
127026e093fcSHong Zhang     ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
12712e8a6d31SBarry Smith   } else {
127226e093fcSHong Zhang     zarray = yarray;
12732e8a6d31SBarry Smith   }
12742d61bbb3SSatish Balay 
12752d61bbb3SSatish Balay   idx = a->j;
12762d61bbb3SSatish Balay   v   = a->a;
127726e093fcSHong Zhang   if (usecprow) {
12784eb6d288SHong Zhang     if (zz != yy) {
12794eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,4*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
12804eb6d288SHong Zhang     }
128126e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
128226e093fcSHong Zhang     ii   = a->compressedrow.i;
12837b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
128426e093fcSHong Zhang   } else {
12852d61bbb3SSatish Balay     ii = a->i;
128626e093fcSHong Zhang     y  = yarray;
128726e093fcSHong Zhang     z  = zarray;
128826e093fcSHong Zhang   }
12892d61bbb3SSatish Balay 
12902d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
12912d61bbb3SSatish Balay     n = ii[1] - ii[0]; ii++;
129226e093fcSHong Zhang     if (usecprow) {
12937b2bb3b9SHong Zhang       z = zarray + 4*ridx[i];
12947b2bb3b9SHong Zhang       y = yarray + 4*ridx[i];
129526e093fcSHong Zhang     }
12962d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3];
1297444d8c10SJed Brown     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);     /* Indices for the next row (assumes same size as this one) */
1298444d8c10SJed Brown     PetscPrefetchBlock(v+16*n,16*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
12992d61bbb3SSatish Balay     for (j=0; j<n; j++) {
13002d61bbb3SSatish Balay       xb    = x + 4*(*idx++);
13012d61bbb3SSatish Balay       x1    = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
13022d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[4]*x2 + v[8]*x3   + v[12]*x4;
13032d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[5]*x2 + v[9]*x3   + v[13]*x4;
13042d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[6]*x2 + v[10]*x3  + v[14]*x4;
13052d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[7]*x2 + v[11]*x3  + v[15]*x4;
13062d61bbb3SSatish Balay       v    += 16;
13072d61bbb3SSatish Balay     }
13082d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4;
130926e093fcSHong Zhang     if (!usecprow) {
13102d61bbb3SSatish Balay       z += 4; y += 4;
13112d61bbb3SSatish Balay     }
131226e093fcSHong Zhang   }
13131ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
131426e093fcSHong Zhang   ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr);
13152e8a6d31SBarry Smith   if (zz != yy) {
131626e093fcSHong Zhang     ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
13172e8a6d31SBarry Smith   }
1318dc0b31edSSatish Balay   ierr = PetscLogFlops(32.0*a->nz);CHKERRQ(ierr);
13192d61bbb3SSatish Balay   PetscFunctionReturn(0);
13202d61bbb3SSatish Balay }
13212d61bbb3SSatish Balay 
13224a2ae208SSatish Balay #undef __FUNCT__
13234a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_5"
1324dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_5(Mat A,Vec xx,Vec yy,Vec zz)
13252d61bbb3SSatish Balay {
13262d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
1327dcf5cc72SBarry Smith   PetscScalar    *x,*y = 0,*z = 0,*xb,sum1,sum2,sum3,sum4,sum5,x1,x2,x3,x4,x5;
132826e093fcSHong Zhang   PetscScalar    *yarray,*zarray;
13293f1db9ecSBarry Smith   MatScalar      *v;
1330dfbe8321SBarry Smith   PetscErrorCode ierr;
13310298fd71SBarry Smith   PetscInt       mbs     =a->mbs,i,*idx,*ii,j,n,*ridx=NULL;
1332ace3abfcSBarry Smith   PetscBool      usecprow=a->compressedrow.use;
13332d61bbb3SSatish Balay 
13342d61bbb3SSatish Balay   PetscFunctionBegin;
13351ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
133626e093fcSHong Zhang   ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr);
13372e8a6d31SBarry Smith   if (zz != yy) {
133826e093fcSHong Zhang     ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
13392e8a6d31SBarry Smith   } else {
134026e093fcSHong Zhang     zarray = yarray;
13412e8a6d31SBarry Smith   }
13422d61bbb3SSatish Balay 
13432d61bbb3SSatish Balay   idx = a->j;
13442d61bbb3SSatish Balay   v   = a->a;
134526e093fcSHong Zhang   if (usecprow) {
13464eb6d288SHong Zhang     if (zz != yy) {
13474eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,5*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
13484eb6d288SHong Zhang     }
134926e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
135026e093fcSHong Zhang     ii   = a->compressedrow.i;
13517b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
135226e093fcSHong Zhang   } else {
13532d61bbb3SSatish Balay     ii = a->i;
135426e093fcSHong Zhang     y  = yarray;
135526e093fcSHong Zhang     z  = zarray;
135626e093fcSHong Zhang   }
13572d61bbb3SSatish Balay 
13582d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
13592d61bbb3SSatish Balay     n = ii[1] - ii[0]; ii++;
136026e093fcSHong Zhang     if (usecprow) {
13617b2bb3b9SHong Zhang       z = zarray + 5*ridx[i];
13627b2bb3b9SHong Zhang       y = yarray + 5*ridx[i];
136326e093fcSHong Zhang     }
13642d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4];
1365444d8c10SJed Brown     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);     /* Indices for the next row (assumes same size as this one) */
1366444d8c10SJed Brown     PetscPrefetchBlock(v+25*n,25*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
13672d61bbb3SSatish Balay     for (j=0; j<n; j++) {
13682d61bbb3SSatish Balay       xb    = x + 5*(*idx++);
13692d61bbb3SSatish Balay       x1    = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4];
13702d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[5]*x2 + v[10]*x3  + v[15]*x4 + v[20]*x5;
13712d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[6]*x2 + v[11]*x3  + v[16]*x4 + v[21]*x5;
13722d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[7]*x2 + v[12]*x3  + v[17]*x4 + v[22]*x5;
13732d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[8]*x2 + v[13]*x3  + v[18]*x4 + v[23]*x5;
13742d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[9]*x2 + v[14]*x3  + v[19]*x4 + v[24]*x5;
13752d61bbb3SSatish Balay       v    += 25;
13762d61bbb3SSatish Balay     }
13772d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5;
137826e093fcSHong Zhang     if (!usecprow) {
13792d61bbb3SSatish Balay       z += 5; y += 5;
13802d61bbb3SSatish Balay     }
138126e093fcSHong Zhang   }
13821ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
138326e093fcSHong Zhang   ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr);
13842e8a6d31SBarry Smith   if (zz != yy) {
138526e093fcSHong Zhang     ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
13862e8a6d31SBarry Smith   }
1387dc0b31edSSatish Balay   ierr = PetscLogFlops(50.0*a->nz);CHKERRQ(ierr);
13882d61bbb3SSatish Balay   PetscFunctionReturn(0);
13892d61bbb3SSatish Balay }
13904a2ae208SSatish Balay #undef __FUNCT__
13914a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_6"
1392dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_6(Mat A,Vec xx,Vec yy,Vec zz)
139315091d37SBarry Smith {
139415091d37SBarry Smith   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
1395dcf5cc72SBarry Smith   PetscScalar    *x,*y = 0,*z = 0,*xb,sum1,sum2,sum3,sum4,sum5,sum6;
139626e093fcSHong Zhang   PetscScalar    x1,x2,x3,x4,x5,x6,*yarray,*zarray;
139715091d37SBarry Smith   MatScalar      *v;
1398dfbe8321SBarry Smith   PetscErrorCode ierr;
13990298fd71SBarry Smith   PetscInt       mbs     =a->mbs,i,*idx,*ii,j,n,*ridx=NULL;
1400ace3abfcSBarry Smith   PetscBool      usecprow=a->compressedrow.use;
140115091d37SBarry Smith 
140215091d37SBarry Smith   PetscFunctionBegin;
14031ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
140426e093fcSHong Zhang   ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr);
140515091d37SBarry Smith   if (zz != yy) {
140626e093fcSHong Zhang     ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
140715091d37SBarry Smith   } else {
140826e093fcSHong Zhang     zarray = yarray;
140915091d37SBarry Smith   }
141015091d37SBarry Smith 
141115091d37SBarry Smith   idx = a->j;
141215091d37SBarry Smith   v   = a->a;
141326e093fcSHong Zhang   if (usecprow) {
14144eb6d288SHong Zhang     if (zz != yy) {
14154eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,6*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
14164eb6d288SHong Zhang     }
141726e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
141826e093fcSHong Zhang     ii   = a->compressedrow.i;
14197b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
142026e093fcSHong Zhang   } else {
142115091d37SBarry Smith     ii = a->i;
142226e093fcSHong Zhang     y  = yarray;
142326e093fcSHong Zhang     z  = zarray;
142426e093fcSHong Zhang   }
142515091d37SBarry Smith 
142615091d37SBarry Smith   for (i=0; i<mbs; i++) {
142715091d37SBarry Smith     n = ii[1] - ii[0]; ii++;
142826e093fcSHong Zhang     if (usecprow) {
14297b2bb3b9SHong Zhang       z = zarray + 6*ridx[i];
14307b2bb3b9SHong Zhang       y = yarray + 6*ridx[i];
143126e093fcSHong Zhang     }
143215091d37SBarry Smith     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4]; sum6 = y[5];
1433444d8c10SJed Brown     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);     /* Indices for the next row (assumes same size as this one) */
1434444d8c10SJed Brown     PetscPrefetchBlock(v+36*n,36*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
143515091d37SBarry Smith     for (j=0; j<n; j++) {
14363b95cb0eSSatish Balay       xb    = x + 6*(*idx++);
143715091d37SBarry Smith       x1    = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5];
143815091d37SBarry Smith       sum1 += v[0]*x1 + v[6]*x2  + v[12]*x3  + v[18]*x4 + v[24]*x5 + v[30]*x6;
143915091d37SBarry Smith       sum2 += v[1]*x1 + v[7]*x2  + v[13]*x3  + v[19]*x4 + v[25]*x5 + v[31]*x6;
144015091d37SBarry Smith       sum3 += v[2]*x1 + v[8]*x2  + v[14]*x3  + v[20]*x4 + v[26]*x5 + v[32]*x6;
144115091d37SBarry Smith       sum4 += v[3]*x1 + v[9]*x2  + v[15]*x3  + v[21]*x4 + v[27]*x5 + v[33]*x6;
144215091d37SBarry Smith       sum5 += v[4]*x1 + v[10]*x2 + v[16]*x3  + v[22]*x4 + v[28]*x5 + v[34]*x6;
144315091d37SBarry Smith       sum6 += v[5]*x1 + v[11]*x2 + v[17]*x3  + v[23]*x4 + v[29]*x5 + v[35]*x6;
144415091d37SBarry Smith       v    += 36;
144515091d37SBarry Smith     }
144615091d37SBarry Smith     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6;
144726e093fcSHong Zhang     if (!usecprow) {
144815091d37SBarry Smith       z += 6; y += 6;
144915091d37SBarry Smith     }
145026e093fcSHong Zhang   }
14511ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
145226e093fcSHong Zhang   ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr);
145315091d37SBarry Smith   if (zz != yy) {
145426e093fcSHong Zhang     ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
145515091d37SBarry Smith   }
1456dc0b31edSSatish Balay   ierr = PetscLogFlops(72.0*a->nz);CHKERRQ(ierr);
145715091d37SBarry Smith   PetscFunctionReturn(0);
145815091d37SBarry Smith }
14592d61bbb3SSatish Balay 
14604a2ae208SSatish Balay #undef __FUNCT__
14614a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_7"
1462dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_7(Mat A,Vec xx,Vec yy,Vec zz)
14632d61bbb3SSatish Balay {
14642d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
1465dcf5cc72SBarry Smith   PetscScalar    *x,*y = 0,*z = 0,*xb,sum1,sum2,sum3,sum4,sum5,sum6,sum7;
146626e093fcSHong Zhang   PetscScalar    x1,x2,x3,x4,x5,x6,x7,*yarray,*zarray;
14673f1db9ecSBarry Smith   MatScalar      *v;
1468dfbe8321SBarry Smith   PetscErrorCode ierr;
14690298fd71SBarry Smith   PetscInt       mbs     =a->mbs,i,*idx,*ii,j,n,*ridx=NULL;
1470ace3abfcSBarry Smith   PetscBool      usecprow=a->compressedrow.use;
14712d61bbb3SSatish Balay 
14722d61bbb3SSatish Balay   PetscFunctionBegin;
14731ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
147426e093fcSHong Zhang   ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr);
14752e8a6d31SBarry Smith   if (zz != yy) {
147626e093fcSHong Zhang     ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
14772e8a6d31SBarry Smith   } else {
147826e093fcSHong Zhang     zarray = yarray;
14792e8a6d31SBarry Smith   }
14802d61bbb3SSatish Balay 
14812d61bbb3SSatish Balay   idx = a->j;
14822d61bbb3SSatish Balay   v   = a->a;
148326e093fcSHong Zhang   if (usecprow) {
14844eb6d288SHong Zhang     if (zz != yy) {
14854eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,7*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
14864eb6d288SHong Zhang     }
148726e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
148826e093fcSHong Zhang     ii   = a->compressedrow.i;
14897b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
149026e093fcSHong Zhang   } else {
14912d61bbb3SSatish Balay     ii = a->i;
149226e093fcSHong Zhang     y  = yarray;
149326e093fcSHong Zhang     z  = zarray;
149426e093fcSHong Zhang   }
14952d61bbb3SSatish Balay 
14962d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
14972d61bbb3SSatish Balay     n = ii[1] - ii[0]; ii++;
149826e093fcSHong Zhang     if (usecprow) {
14997b2bb3b9SHong Zhang       z = zarray + 7*ridx[i];
15007b2bb3b9SHong Zhang       y = yarray + 7*ridx[i];
150126e093fcSHong Zhang     }
15022d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4]; sum6 = y[5]; sum7 = y[6];
1503444d8c10SJed Brown     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);     /* Indices for the next row (assumes same size as this one) */
1504444d8c10SJed Brown     PetscPrefetchBlock(v+49*n,49*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
15052d61bbb3SSatish Balay     for (j=0; j<n; j++) {
15062d61bbb3SSatish Balay       xb    = x + 7*(*idx++);
15072d61bbb3SSatish Balay       x1    = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6];
15082d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[7]*x2  + v[14]*x3  + v[21]*x4 + v[28]*x5 + v[35]*x6 + v[42]*x7;
15092d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[8]*x2  + v[15]*x3  + v[22]*x4 + v[29]*x5 + v[36]*x6 + v[43]*x7;
15102d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[9]*x2  + v[16]*x3  + v[23]*x4 + v[30]*x5 + v[37]*x6 + v[44]*x7;
15112d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[10]*x2 + v[17]*x3  + v[24]*x4 + v[31]*x5 + v[38]*x6 + v[45]*x7;
15122d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[11]*x2 + v[18]*x3  + v[25]*x4 + v[32]*x5 + v[39]*x6 + v[46]*x7;
15132d61bbb3SSatish Balay       sum6 += v[5]*x1 + v[12]*x2 + v[19]*x3  + v[26]*x4 + v[33]*x5 + v[40]*x6 + v[47]*x7;
15142d61bbb3SSatish Balay       sum7 += v[6]*x1 + v[13]*x2 + v[20]*x3  + v[27]*x4 + v[34]*x5 + v[41]*x6 + v[48]*x7;
15152d61bbb3SSatish Balay       v    += 49;
15162d61bbb3SSatish Balay     }
15172d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
151826e093fcSHong Zhang     if (!usecprow) {
15192d61bbb3SSatish Balay       z += 7; y += 7;
15202d61bbb3SSatish Balay     }
152126e093fcSHong Zhang   }
15221ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
152326e093fcSHong Zhang   ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr);
15242e8a6d31SBarry Smith   if (zz != yy) {
152526e093fcSHong Zhang     ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
15262e8a6d31SBarry Smith   }
1527dc0b31edSSatish Balay   ierr = PetscLogFlops(98.0*a->nz);CHKERRQ(ierr);
15282d61bbb3SSatish Balay   PetscFunctionReturn(0);
15292d61bbb3SSatish Balay }
1530218c64b6SSatish Balay 
15314a2ae208SSatish Balay #undef __FUNCT__
15324a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_N"
1533dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_N(Mat A,Vec xx,Vec yy,Vec zz)
15342d61bbb3SSatish Balay {
15352d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
1536bba805e6SBarry Smith   PetscScalar    *x,*z = 0,*xb,*work,*workt,*zarray;
15373f1db9ecSBarry Smith   MatScalar      *v;
15386849ba73SBarry Smith   PetscErrorCode ierr;
1539d0f46423SBarry Smith   PetscInt       mbs,i,*idx,*ii,bs=A->rmap->bs,j,n,bs2=a->bs2;
15400298fd71SBarry Smith   PetscInt       ncols,k,*ridx=NULL;
1541ace3abfcSBarry Smith   PetscBool      usecprow=a->compressedrow.use;
1542218c64b6SSatish Balay 
15432d61bbb3SSatish Balay   PetscFunctionBegin;
1544343551c4SBarry Smith   ierr = VecCopy(yy,zz);CHKERRQ(ierr);
15451ebc52fbSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
154626e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
15472d61bbb3SSatish Balay 
15482d61bbb3SSatish Balay   idx = a->j;
15492d61bbb3SSatish Balay   v   = a->a;
155026e093fcSHong Zhang   if (usecprow) {
155126e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
155226e093fcSHong Zhang     ii   = a->compressedrow.i;
15537b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
155426e093fcSHong Zhang   } else {
155526e093fcSHong Zhang     mbs = a->mbs;
15562d61bbb3SSatish Balay     ii  = a->i;
155726e093fcSHong Zhang     z   = zarray;
155826e093fcSHong Zhang   }
15592d61bbb3SSatish Balay 
15602d61bbb3SSatish Balay   if (!a->mult_work) {
1561d0f46423SBarry Smith     k    = PetscMax(A->rmap->n,A->cmap->n);
156287828ca2SBarry Smith     ierr = PetscMalloc((k+1)*sizeof(PetscScalar),&a->mult_work);CHKERRQ(ierr);
15632d61bbb3SSatish Balay   }
15642d61bbb3SSatish Balay   work = a->mult_work;
15652d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
15662d61bbb3SSatish Balay     n     = ii[1] - ii[0]; ii++;
15672d61bbb3SSatish Balay     ncols = n*bs;
15682d61bbb3SSatish Balay     workt = work;
15692d61bbb3SSatish Balay     for (j=0; j<n; j++) {
15702d61bbb3SSatish Balay       xb = x + bs*(*idx++);
15712d61bbb3SSatish Balay       for (k=0; k<bs; k++) workt[k] = xb[k];
15722d61bbb3SSatish Balay       workt += bs;
15732d61bbb3SSatish Balay     }
15747b2bb3b9SHong Zhang     if (usecprow) z = zarray + bs*ridx[i];
157596b95a6bSBarry Smith     PetscKernel_w_gets_w_plus_Ar_times_v(bs,ncols,work,v,z);
157671044d3cSBarry Smith     /* BLASgemv_("N",&bs,&ncols,&_DOne,v,&bs,work,&_One,&_DOne,z,&_One); */
15772d61bbb3SSatish Balay     v += n*bs2;
157826fbe8dcSKarl Rupp     if (!usecprow) z += bs;
157926e093fcSHong Zhang   }
15801ebc52fbSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
158126e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
1582dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*a->nz*bs2);CHKERRQ(ierr);
15832d61bbb3SSatish Balay   PetscFunctionReturn(0);
15842d61bbb3SSatish Balay }
15852d61bbb3SSatish Balay 
15864a2ae208SSatish Balay #undef __FUNCT__
1587547795f9SHong Zhang #define __FUNCT__ "MatMultHermitianTranspose_SeqBAIJ"
1588547795f9SHong Zhang PetscErrorCode MatMultHermitianTranspose_SeqBAIJ(Mat A,Vec xx,Vec zz)
1589547795f9SHong Zhang {
1590547795f9SHong Zhang   PetscScalar    zero = 0.0;
1591547795f9SHong Zhang   PetscErrorCode ierr;
1592547795f9SHong Zhang 
1593547795f9SHong Zhang   PetscFunctionBegin;
1594547795f9SHong Zhang   ierr = VecSet(zz,zero);CHKERRQ(ierr);
1595547795f9SHong Zhang   ierr = MatMultHermitianTransposeAdd_SeqBAIJ(A,xx,zz,zz);CHKERRQ(ierr);
1596547795f9SHong Zhang   PetscFunctionReturn(0);
1597547795f9SHong Zhang }
1598547795f9SHong Zhang 
1599547795f9SHong Zhang #undef __FUNCT__
16004a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_SeqBAIJ"
1601dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_SeqBAIJ(Mat A,Vec xx,Vec zz)
16022d61bbb3SSatish Balay {
16033447b6efSHong Zhang   PetscScalar    zero = 0.0;
16046849ba73SBarry Smith   PetscErrorCode ierr;
16052d61bbb3SSatish Balay 
16062d61bbb3SSatish Balay   PetscFunctionBegin;
16072dcb1b2aSMatthew Knepley   ierr = VecSet(zz,zero);CHKERRQ(ierr);
16083447b6efSHong Zhang   ierr = MatMultTransposeAdd_SeqBAIJ(A,xx,zz,zz);CHKERRQ(ierr);
16092d61bbb3SSatish Balay   PetscFunctionReturn(0);
16102d61bbb3SSatish Balay }
16112d61bbb3SSatish Balay 
16124a2ae208SSatish Balay #undef __FUNCT__
1613547795f9SHong Zhang #define __FUNCT__ "MatMultHermitianTransposeAdd_SeqBAIJ"
1614547795f9SHong Zhang PetscErrorCode MatMultHermitianTransposeAdd_SeqBAIJ(Mat A,Vec xx,Vec yy,Vec zz)
1615547795f9SHong Zhang {
1616547795f9SHong Zhang   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
1617547795f9SHong Zhang   PetscScalar       *zb,*x,*z,*xb = 0,x1,x2,x3,x4,x5;
1618547795f9SHong Zhang   MatScalar         *v;
1619547795f9SHong Zhang   PetscErrorCode    ierr;
16200298fd71SBarry Smith   PetscInt          mbs,i,*idx,*ii,rval,bs=A->rmap->bs,j,n,bs2=a->bs2,*ib,*ridx=NULL;
1621547795f9SHong Zhang   Mat_CompressedRow cprow   = a->compressedrow;
1622ace3abfcSBarry Smith   PetscBool         usecprow=cprow.use;
1623547795f9SHong Zhang 
1624547795f9SHong Zhang   PetscFunctionBegin;
1625343551c4SBarry Smith   if (yy != zz) { ierr = VecCopy(yy,zz);CHKERRQ(ierr); }
1626547795f9SHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
1627547795f9SHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
1628547795f9SHong Zhang 
1629547795f9SHong Zhang   idx = a->j;
1630547795f9SHong Zhang   v   = a->a;
1631547795f9SHong Zhang   if (usecprow) {
1632547795f9SHong Zhang     mbs  = cprow.nrows;
1633547795f9SHong Zhang     ii   = cprow.i;
1634547795f9SHong Zhang     ridx = cprow.rindex;
1635547795f9SHong Zhang   } else {
1636547795f9SHong Zhang     mbs=a->mbs;
1637547795f9SHong Zhang     ii = a->i;
1638547795f9SHong Zhang     xb = x;
1639547795f9SHong Zhang   }
1640547795f9SHong Zhang 
1641547795f9SHong Zhang   switch (bs) {
1642547795f9SHong Zhang   case 1:
1643547795f9SHong Zhang     for (i=0; i<mbs; i++) {
1644547795f9SHong Zhang       if (usecprow) xb = x + ridx[i];
1645547795f9SHong Zhang       x1 = xb[0];
1646547795f9SHong Zhang       ib = idx + ii[0];
1647547795f9SHong Zhang       n  = ii[1] - ii[0]; ii++;
1648547795f9SHong Zhang       for (j=0; j<n; j++) {
1649547795f9SHong Zhang         rval     = ib[j];
1650547795f9SHong Zhang         z[rval] += PetscConj(*v) * x1;
1651547795f9SHong Zhang         v++;
1652547795f9SHong Zhang       }
1653547795f9SHong Zhang       if (!usecprow) xb++;
1654547795f9SHong Zhang     }
1655547795f9SHong Zhang     break;
1656547795f9SHong Zhang   case 2:
1657547795f9SHong Zhang     for (i=0; i<mbs; i++) {
1658547795f9SHong Zhang       if (usecprow) xb = x + 2*ridx[i];
1659547795f9SHong Zhang       x1 = xb[0]; x2 = xb[1];
1660547795f9SHong Zhang       ib = idx + ii[0];
1661547795f9SHong Zhang       n  = ii[1] - ii[0]; ii++;
1662547795f9SHong Zhang       for (j=0; j<n; j++) {
1663547795f9SHong Zhang         rval       = ib[j]*2;
1664547795f9SHong Zhang         z[rval++] += PetscConj(v[0])*x1 + PetscConj(v[1])*x2;
1665547795f9SHong Zhang         z[rval++] += PetscConj(v[2])*x1 + PetscConj(v[3])*x2;
1666547795f9SHong Zhang         v         += 4;
1667547795f9SHong Zhang       }
1668547795f9SHong Zhang       if (!usecprow) xb += 2;
1669547795f9SHong Zhang     }
1670547795f9SHong Zhang     break;
1671547795f9SHong Zhang   case 3:
1672547795f9SHong Zhang     for (i=0; i<mbs; i++) {
1673547795f9SHong Zhang       if (usecprow) xb = x + 3*ridx[i];
1674547795f9SHong Zhang       x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
1675547795f9SHong Zhang       ib = idx + ii[0];
1676547795f9SHong Zhang       n  = ii[1] - ii[0]; ii++;
1677547795f9SHong Zhang       for (j=0; j<n; j++) {
1678547795f9SHong Zhang         rval       = ib[j]*3;
1679547795f9SHong Zhang         z[rval++] += PetscConj(v[0])*x1 + PetscConj(v[1])*x2 + PetscConj(v[2])*x3;
1680547795f9SHong Zhang         z[rval++] += PetscConj(v[3])*x1 + PetscConj(v[4])*x2 + PetscConj(v[5])*x3;
1681547795f9SHong Zhang         z[rval++] += PetscConj(v[6])*x1 + PetscConj(v[7])*x2 + PetscConj(v[8])*x3;
1682547795f9SHong Zhang         v         += 9;
1683547795f9SHong Zhang       }
1684547795f9SHong Zhang       if (!usecprow) xb += 3;
1685547795f9SHong Zhang     }
1686547795f9SHong Zhang     break;
1687547795f9SHong Zhang   case 4:
1688547795f9SHong Zhang     for (i=0; i<mbs; i++) {
1689547795f9SHong Zhang       if (usecprow) xb = x + 4*ridx[i];
1690547795f9SHong Zhang       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
1691547795f9SHong Zhang       ib = idx + ii[0];
1692547795f9SHong Zhang       n  = ii[1] - ii[0]; ii++;
1693547795f9SHong Zhang       for (j=0; j<n; j++) {
1694547795f9SHong Zhang         rval       = ib[j]*4;
1695547795f9SHong Zhang         z[rval++] +=  PetscConj(v[0])*x1 + PetscConj(v[1])*x2  + PetscConj(v[2])*x3  + PetscConj(v[3])*x4;
1696547795f9SHong Zhang         z[rval++] +=  PetscConj(v[4])*x1 + PetscConj(v[5])*x2  + PetscConj(v[6])*x3  + PetscConj(v[7])*x4;
1697547795f9SHong Zhang         z[rval++] +=  PetscConj(v[8])*x1 + PetscConj(v[9])*x2  + PetscConj(v[10])*x3 + PetscConj(v[11])*x4;
1698547795f9SHong Zhang         z[rval++] += PetscConj(v[12])*x1 + PetscConj(v[13])*x2 + PetscConj(v[14])*x3 + PetscConj(v[15])*x4;
1699547795f9SHong Zhang         v         += 16;
1700547795f9SHong Zhang       }
1701547795f9SHong Zhang       if (!usecprow) xb += 4;
1702547795f9SHong Zhang     }
1703547795f9SHong Zhang     break;
1704547795f9SHong Zhang   case 5:
1705547795f9SHong Zhang     for (i=0; i<mbs; i++) {
1706547795f9SHong Zhang       if (usecprow) xb = x + 5*ridx[i];
1707547795f9SHong Zhang       x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
1708547795f9SHong Zhang       x4 = xb[3]; x5 = xb[4];
1709547795f9SHong Zhang       ib = idx + ii[0];
1710547795f9SHong Zhang       n  = ii[1] - ii[0]; ii++;
1711547795f9SHong Zhang       for (j=0; j<n; j++) {
1712547795f9SHong Zhang         rval       = ib[j]*5;
1713547795f9SHong Zhang         z[rval++] +=  PetscConj(v[0])*x1 +  PetscConj(v[1])*x2 +  PetscConj(v[2])*x3 +  PetscConj(v[3])*x4 +  PetscConj(v[4])*x5;
1714547795f9SHong Zhang         z[rval++] +=  PetscConj(v[5])*x1 +  PetscConj(v[6])*x2 +  PetscConj(v[7])*x3 +  PetscConj(v[8])*x4 +  PetscConj(v[9])*x5;
1715547795f9SHong Zhang         z[rval++] += PetscConj(v[10])*x1 + PetscConj(v[11])*x2 + PetscConj(v[12])*x3 + PetscConj(v[13])*x4 + PetscConj(v[14])*x5;
1716547795f9SHong Zhang         z[rval++] += PetscConj(v[15])*x1 + PetscConj(v[16])*x2 + PetscConj(v[17])*x3 + PetscConj(v[18])*x4 + PetscConj(v[19])*x5;
1717547795f9SHong Zhang         z[rval++] += PetscConj(v[20])*x1 + PetscConj(v[21])*x2 + PetscConj(v[22])*x3 + PetscConj(v[23])*x4 + PetscConj(v[24])*x5;
1718547795f9SHong Zhang         v         += 25;
1719547795f9SHong Zhang       }
1720547795f9SHong Zhang       if (!usecprow) xb += 5;
1721547795f9SHong Zhang     }
1722547795f9SHong Zhang     break;
1723547795f9SHong Zhang   default: {      /* block sizes larger than 5 by 5 are handled by BLAS */
1724547795f9SHong Zhang     PetscInt    ncols,k;
1725547795f9SHong Zhang     PetscScalar *work,*workt,*xtmp;
1726547795f9SHong Zhang 
1727e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"block size larger than 5 is not supported yet");
1728547795f9SHong Zhang     if (!a->mult_work) {
1729547795f9SHong Zhang       k    = PetscMax(A->rmap->n,A->cmap->n);
1730547795f9SHong Zhang       ierr = PetscMalloc((k+1)*sizeof(PetscScalar),&a->mult_work);CHKERRQ(ierr);
1731547795f9SHong Zhang     }
1732547795f9SHong Zhang     work = a->mult_work;
1733547795f9SHong Zhang     xtmp = x;
1734547795f9SHong Zhang     for (i=0; i<mbs; i++) {
1735547795f9SHong Zhang       n     = ii[1] - ii[0]; ii++;
1736547795f9SHong Zhang       ncols = n*bs;
1737547795f9SHong Zhang       ierr  = PetscMemzero(work,ncols*sizeof(PetscScalar));CHKERRQ(ierr);
173826fbe8dcSKarl Rupp       if (usecprow) xtmp = x + bs*ridx[i];
173996b95a6bSBarry Smith       PetscKernel_w_gets_w_plus_trans_Ar_times_v(bs,ncols,xtmp,v,work);
1740547795f9SHong Zhang       /* BLASgemv_("T",&bs,&ncols,&_DOne,v,&bs,xtmp,&_One,&_DOne,work,&_One); */
1741547795f9SHong Zhang       v += n*bs2;
1742547795f9SHong Zhang       if (!usecprow) xtmp += bs;
1743547795f9SHong Zhang       workt = work;
1744547795f9SHong Zhang       for (j=0; j<n; j++) {
1745547795f9SHong Zhang         zb = z + bs*(*idx++);
1746547795f9SHong Zhang         for (k=0; k<bs; k++) zb[k] += workt[k] ;
1747547795f9SHong Zhang         workt += bs;
1748547795f9SHong Zhang       }
1749547795f9SHong Zhang     }
1750547795f9SHong Zhang     }
1751547795f9SHong Zhang   }
1752547795f9SHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
1753547795f9SHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
1754547795f9SHong Zhang   ierr = PetscLogFlops(2.0*a->nz*a->bs2);CHKERRQ(ierr);
1755547795f9SHong Zhang   PetscFunctionReturn(0);
1756547795f9SHong Zhang }
1757547795f9SHong Zhang 
1758547795f9SHong Zhang #undef __FUNCT__
17594a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_SeqBAIJ"
1760dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_SeqBAIJ(Mat A,Vec xx,Vec yy,Vec zz)
17612d61bbb3SSatish Balay {
17622d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
1763dcf5cc72SBarry Smith   PetscScalar       *zb,*x,*z,*xb = 0,x1,x2,x3,x4,x5;
17643f1db9ecSBarry Smith   MatScalar         *v;
17656849ba73SBarry Smith   PetscErrorCode    ierr;
17660298fd71SBarry Smith   PetscInt          mbs,i,*idx,*ii,rval,bs=A->rmap->bs,j,n,bs2=a->bs2,*ib,*ridx=NULL;
17673447b6efSHong Zhang   Mat_CompressedRow cprow   = a->compressedrow;
1768ace3abfcSBarry Smith   PetscBool         usecprow=cprow.use;
17692d61bbb3SSatish Balay 
17702d61bbb3SSatish Balay   PetscFunctionBegin;
1771343551c4SBarry Smith   if (yy != zz) { ierr = VecCopy(yy,zz);CHKERRQ(ierr); }
17723447b6efSHong Zhang   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
17733447b6efSHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
17742d61bbb3SSatish Balay 
17752d61bbb3SSatish Balay   idx = a->j;
17762d61bbb3SSatish Balay   v   = a->a;
17773447b6efSHong Zhang   if (usecprow) {
17783447b6efSHong Zhang     mbs  = cprow.nrows;
17793447b6efSHong Zhang     ii   = cprow.i;
17807b2bb3b9SHong Zhang     ridx = cprow.rindex;
17813447b6efSHong Zhang   } else {
17823447b6efSHong Zhang     mbs=a->mbs;
17832d61bbb3SSatish Balay     ii = a->i;
1784f1af5d2fSBarry Smith     xb = x;
17853447b6efSHong Zhang   }
17862d61bbb3SSatish Balay 
17872d61bbb3SSatish Balay   switch (bs) {
17882d61bbb3SSatish Balay   case 1:
17892d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
17907b2bb3b9SHong Zhang       if (usecprow) xb = x + ridx[i];
1791f1af5d2fSBarry Smith       x1 = xb[0];
17923447b6efSHong Zhang       ib = idx + ii[0];
17933447b6efSHong Zhang       n  = ii[1] - ii[0]; ii++;
17942d61bbb3SSatish Balay       for (j=0; j<n; j++) {
17952d61bbb3SSatish Balay         rval     = ib[j];
1796f1af5d2fSBarry Smith         z[rval] += *v * x1;
1797f1af5d2fSBarry Smith         v++;
17982d61bbb3SSatish Balay       }
17993447b6efSHong Zhang       if (!usecprow) xb++;
18002d61bbb3SSatish Balay     }
18012d61bbb3SSatish Balay     break;
18022d61bbb3SSatish Balay   case 2:
18032d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
18047b2bb3b9SHong Zhang       if (usecprow) xb = x + 2*ridx[i];
1805f1af5d2fSBarry Smith       x1 = xb[0]; x2 = xb[1];
18063447b6efSHong Zhang       ib = idx + ii[0];
18073447b6efSHong Zhang       n  = ii[1] - ii[0]; ii++;
18082d61bbb3SSatish Balay       for (j=0; j<n; j++) {
18092d61bbb3SSatish Balay         rval       = ib[j]*2;
18102d61bbb3SSatish Balay         z[rval++] += v[0]*x1 + v[1]*x2;
18112d61bbb3SSatish Balay         z[rval++] += v[2]*x1 + v[3]*x2;
18122d61bbb3SSatish Balay         v         += 4;
18132d61bbb3SSatish Balay       }
18143447b6efSHong Zhang       if (!usecprow) xb += 2;
18152d61bbb3SSatish Balay     }
18162d61bbb3SSatish Balay     break;
18172d61bbb3SSatish Balay   case 3:
18182d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
18197b2bb3b9SHong Zhang       if (usecprow) xb = x + 3*ridx[i];
1820f1af5d2fSBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
18213447b6efSHong Zhang       ib = idx + ii[0];
18223447b6efSHong Zhang       n  = ii[1] - ii[0]; ii++;
18232d61bbb3SSatish Balay       for (j=0; j<n; j++) {
18242d61bbb3SSatish Balay         rval       = ib[j]*3;
18252d61bbb3SSatish Balay         z[rval++] += v[0]*x1 + v[1]*x2 + v[2]*x3;
18262d61bbb3SSatish Balay         z[rval++] += v[3]*x1 + v[4]*x2 + v[5]*x3;
18272d61bbb3SSatish Balay         z[rval++] += v[6]*x1 + v[7]*x2 + v[8]*x3;
18282d61bbb3SSatish Balay         v         += 9;
18292d61bbb3SSatish Balay       }
18303447b6efSHong Zhang       if (!usecprow) xb += 3;
18312d61bbb3SSatish Balay     }
18322d61bbb3SSatish Balay     break;
18332d61bbb3SSatish Balay   case 4:
18342d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
18357b2bb3b9SHong Zhang       if (usecprow) xb = x + 4*ridx[i];
1836f1af5d2fSBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
18373447b6efSHong Zhang       ib = idx + ii[0];
18383447b6efSHong Zhang       n  = ii[1] - ii[0]; ii++;
18392d61bbb3SSatish Balay       for (j=0; j<n; j++) {
18402d61bbb3SSatish Balay         rval       = ib[j]*4;
18412d61bbb3SSatish Balay         z[rval++] +=  v[0]*x1 +  v[1]*x2 +  v[2]*x3 +  v[3]*x4;
18422d61bbb3SSatish Balay         z[rval++] +=  v[4]*x1 +  v[5]*x2 +  v[6]*x3 +  v[7]*x4;
18432d61bbb3SSatish Balay         z[rval++] +=  v[8]*x1 +  v[9]*x2 + v[10]*x3 + v[11]*x4;
18442d61bbb3SSatish Balay         z[rval++] += v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4;
18452d61bbb3SSatish Balay         v         += 16;
18462d61bbb3SSatish Balay       }
18473447b6efSHong Zhang       if (!usecprow) xb += 4;
18482d61bbb3SSatish Balay     }
18492d61bbb3SSatish Balay     break;
18502d61bbb3SSatish Balay   case 5:
18512d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
18527b2bb3b9SHong Zhang       if (usecprow) xb = x + 5*ridx[i];
1853f1af5d2fSBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
18542d61bbb3SSatish Balay       x4 = xb[3]; x5 = xb[4];
18553447b6efSHong Zhang       ib = idx + ii[0];
18563447b6efSHong Zhang       n  = ii[1] - ii[0]; ii++;
18572d61bbb3SSatish Balay       for (j=0; j<n; j++) {
18582d61bbb3SSatish Balay         rval       = ib[j]*5;
18592d61bbb3SSatish Balay         z[rval++] +=  v[0]*x1 +  v[1]*x2 +  v[2]*x3 +  v[3]*x4 +  v[4]*x5;
18602d61bbb3SSatish Balay         z[rval++] +=  v[5]*x1 +  v[6]*x2 +  v[7]*x3 +  v[8]*x4 +  v[9]*x5;
18612d61bbb3SSatish Balay         z[rval++] += v[10]*x1 + v[11]*x2 + v[12]*x3 + v[13]*x4 + v[14]*x5;
18622d61bbb3SSatish Balay         z[rval++] += v[15]*x1 + v[16]*x2 + v[17]*x3 + v[18]*x4 + v[19]*x5;
18632d61bbb3SSatish Balay         z[rval++] += v[20]*x1 + v[21]*x2 + v[22]*x3 + v[23]*x4 + v[24]*x5;
18642d61bbb3SSatish Balay         v         += 25;
18652d61bbb3SSatish Balay       }
18663447b6efSHong Zhang       if (!usecprow) xb += 5;
18672d61bbb3SSatish Balay     }
18682d61bbb3SSatish Balay     break;
1869f1af5d2fSBarry Smith   default: {      /* block sizes larger then 5 by 5 are handled by BLAS */
1870690b6cddSBarry Smith     PetscInt    ncols,k;
18713447b6efSHong Zhang     PetscScalar *work,*workt,*xtmp;
18723f1db9ecSBarry Smith 
18732d61bbb3SSatish Balay     if (!a->mult_work) {
1874d0f46423SBarry Smith       k    = PetscMax(A->rmap->n,A->cmap->n);
187587828ca2SBarry Smith       ierr = PetscMalloc((k+1)*sizeof(PetscScalar),&a->mult_work);CHKERRQ(ierr);
18762d61bbb3SSatish Balay     }
18772d61bbb3SSatish Balay     work = a->mult_work;
18783447b6efSHong Zhang     xtmp = x;
18792d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
18802d61bbb3SSatish Balay       n     = ii[1] - ii[0]; ii++;
18812d61bbb3SSatish Balay       ncols = n*bs;
188287828ca2SBarry Smith       ierr  = PetscMemzero(work,ncols*sizeof(PetscScalar));CHKERRQ(ierr);
188326fbe8dcSKarl Rupp       if (usecprow) xtmp = x + bs*ridx[i];
188496b95a6bSBarry Smith       PetscKernel_w_gets_w_plus_trans_Ar_times_v(bs,ncols,xtmp,v,work);
188571044d3cSBarry Smith       /* BLASgemv_("T",&bs,&ncols,&_DOne,v,&bs,xtmp,&_One,&_DOne,work,&_One); */
18862d61bbb3SSatish Balay       v += n*bs2;
18873447b6efSHong Zhang       if (!usecprow) xtmp += bs;
18882d61bbb3SSatish Balay       workt = work;
18892d61bbb3SSatish Balay       for (j=0; j<n; j++) {
18902d61bbb3SSatish Balay         zb = z + bs*(*idx++);
18912d61bbb3SSatish Balay         for (k=0; k<bs; k++) zb[k] += workt[k];
18922d61bbb3SSatish Balay         workt += bs;
18932d61bbb3SSatish Balay       }
18942d61bbb3SSatish Balay     }
18952d61bbb3SSatish Balay     }
18962d61bbb3SSatish Balay   }
18973447b6efSHong Zhang   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
18983447b6efSHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
1899dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*a->nz*a->bs2);CHKERRQ(ierr);
19002d61bbb3SSatish Balay   PetscFunctionReturn(0);
19012d61bbb3SSatish Balay }
19022d61bbb3SSatish Balay 
19034a2ae208SSatish Balay #undef __FUNCT__
19044a2ae208SSatish Balay #define __FUNCT__ "MatScale_SeqBAIJ"
1905f4df32b1SMatthew Knepley PetscErrorCode MatScale_SeqBAIJ(Mat inA,PetscScalar alpha)
19062d61bbb3SSatish Balay {
19072d61bbb3SSatish Balay   Mat_SeqBAIJ    *a      = (Mat_SeqBAIJ*)inA->data;
1908690b6cddSBarry Smith   PetscInt       totalnz = a->bs2*a->nz;
1909f4df32b1SMatthew Knepley   PetscScalar    oalpha  = alpha;
1910efee365bSSatish Balay   PetscErrorCode ierr;
1911c5df96a5SBarry Smith   PetscBLASInt   one = 1,tnz;
19122d61bbb3SSatish Balay 
19132d61bbb3SSatish Balay   PetscFunctionBegin;
1914c5df96a5SBarry Smith   ierr = PetscBLASIntCast(totalnz,&tnz);CHKERRQ(ierr);
1915*8b83055fSJed Brown   PetscStackCallBLAS("BLASscal",BLASscal_(&tnz,&oalpha,a->a,&one));
1916efee365bSSatish Balay   ierr = PetscLogFlops(totalnz);CHKERRQ(ierr);
19172d61bbb3SSatish Balay   PetscFunctionReturn(0);
19182d61bbb3SSatish Balay }
19192d61bbb3SSatish Balay 
19204a2ae208SSatish Balay #undef __FUNCT__
19214a2ae208SSatish Balay #define __FUNCT__ "MatNorm_SeqBAIJ"
1922dfbe8321SBarry Smith PetscErrorCode MatNorm_SeqBAIJ(Mat A,NormType type,PetscReal *norm)
19232d61bbb3SSatish Balay {
19248a62d963SHong Zhang   PetscErrorCode ierr;
19252d61bbb3SSatish Balay   Mat_SeqBAIJ    *a  = (Mat_SeqBAIJ*)A->data;
19263f1db9ecSBarry Smith   MatScalar      *v  = a->a;
1927329f5518SBarry Smith   PetscReal      sum = 0.0;
1928d0f46423SBarry Smith   PetscInt       i,j,k,bs=A->rmap->bs,nz=a->nz,bs2=a->bs2,k1;
19292d61bbb3SSatish Balay 
19302d61bbb3SSatish Balay   PetscFunctionBegin;
19312d61bbb3SSatish Balay   if (type == NORM_FROBENIUS) {
19322d61bbb3SSatish Balay     for (i=0; i< bs2*nz; i++) {
1933329f5518SBarry Smith       sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
19342d61bbb3SSatish Balay     }
19358f1a2a5eSBarry Smith     *norm = PetscSqrtReal(sum);
19368a62d963SHong Zhang   } else if (type == NORM_1) { /* maximum column sum */
19378a62d963SHong Zhang     PetscReal *tmp;
19388a62d963SHong Zhang     PetscInt  *bcol = a->j;
1939d0f46423SBarry Smith     ierr = PetscMalloc((A->cmap->n+1)*sizeof(PetscReal),&tmp);CHKERRQ(ierr);
1940d0f46423SBarry Smith     ierr = PetscMemzero(tmp,A->cmap->n*sizeof(PetscReal));CHKERRQ(ierr);
19418a62d963SHong Zhang     for (i=0; i<nz; i++) {
19428a62d963SHong Zhang       for (j=0; j<bs; j++) {
19438a62d963SHong Zhang         k1 = bs*(*bcol) + j; /* column index */
19448a62d963SHong Zhang         for (k=0; k<bs; k++) {
19458a62d963SHong Zhang           tmp[k1] += PetscAbsScalar(*v); v++;
19468a62d963SHong Zhang         }
19478a62d963SHong Zhang       }
19488a62d963SHong Zhang       bcol++;
19498a62d963SHong Zhang     }
19508a62d963SHong Zhang     *norm = 0.0;
1951d0f46423SBarry Smith     for (j=0; j<A->cmap->n; j++) {
19528a62d963SHong Zhang       if (tmp[j] > *norm) *norm = tmp[j];
19538a62d963SHong Zhang     }
19548a62d963SHong Zhang     ierr = PetscFree(tmp);CHKERRQ(ierr);
1955596552b5SBarry Smith   } else if (type == NORM_INFINITY) { /* maximum row sum */
1956596552b5SBarry Smith     *norm = 0.0;
1957596552b5SBarry Smith     for (k=0; k<bs; k++) {
195874f84c7bSSatish Balay       for (j=0; j<a->mbs; j++) {
1959596552b5SBarry Smith         v   = a->a + bs2*a->i[j] + k;
1960596552b5SBarry Smith         sum = 0.0;
1961596552b5SBarry Smith         for (i=0; i<a->i[j+1]-a->i[j]; i++) {
19620e90e235SBarry Smith           for (k1=0; k1<bs; k1++) {
1963596552b5SBarry Smith             sum += PetscAbsScalar(*v);
1964596552b5SBarry Smith             v   += bs;
19652d61bbb3SSatish Balay           }
19660e90e235SBarry Smith         }
1967596552b5SBarry Smith         if (sum > *norm) *norm = sum;
1968596552b5SBarry Smith       }
1969596552b5SBarry Smith     }
1970e7e72b3dSBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for this norm yet");
19712d61bbb3SSatish Balay   PetscFunctionReturn(0);
19722d61bbb3SSatish Balay }
19732d61bbb3SSatish Balay 
19742d61bbb3SSatish Balay 
19754a2ae208SSatish Balay #undef __FUNCT__
19764a2ae208SSatish Balay #define __FUNCT__ "MatEqual_SeqBAIJ"
1977ace3abfcSBarry Smith PetscErrorCode MatEqual_SeqBAIJ(Mat A,Mat B,PetscBool * flg)
19782d61bbb3SSatish Balay {
19792d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)B->data;
1980dfbe8321SBarry Smith   PetscErrorCode ierr;
19812d61bbb3SSatish Balay 
19822d61bbb3SSatish Balay   PetscFunctionBegin;
19832d61bbb3SSatish Balay   /* If the  matrix/block dimensions are not equal, or no of nonzeros or shift */
1984d0f46423SBarry Smith   if ((A->rmap->N != B->rmap->N) || (A->cmap->n != B->cmap->n) || (A->rmap->bs != B->rmap->bs)|| (a->nz != b->nz)) {
1985273d9f13SBarry Smith     *flg = PETSC_FALSE;
1986273d9f13SBarry Smith     PetscFunctionReturn(0);
19872d61bbb3SSatish Balay   }
19882d61bbb3SSatish Balay 
19892d61bbb3SSatish Balay   /* if the a->i are the same */
1990690b6cddSBarry Smith   ierr = PetscMemcmp(a->i,b->i,(a->mbs+1)*sizeof(PetscInt),flg);CHKERRQ(ierr);
199126fbe8dcSKarl Rupp   if (!*flg) PetscFunctionReturn(0);
19922d61bbb3SSatish Balay 
19932d61bbb3SSatish Balay   /* if a->j are the same */
1994690b6cddSBarry Smith   ierr = PetscMemcmp(a->j,b->j,(a->nz)*sizeof(PetscInt),flg);CHKERRQ(ierr);
199526fbe8dcSKarl Rupp   if (!*flg) PetscFunctionReturn(0);
199626fbe8dcSKarl Rupp 
19972d61bbb3SSatish Balay   /* if a->a are the same */
1998d0f46423SBarry Smith   ierr = PetscMemcmp(a->a,b->a,(a->nz)*(A->rmap->bs)*(B->rmap->bs)*sizeof(PetscScalar),flg);CHKERRQ(ierr);
19992d61bbb3SSatish Balay   PetscFunctionReturn(0);
20002d61bbb3SSatish Balay 
20012d61bbb3SSatish Balay }
20022d61bbb3SSatish Balay 
20034a2ae208SSatish Balay #undef __FUNCT__
20044a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_SeqBAIJ"
2005dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_SeqBAIJ(Mat A,Vec v)
20062d61bbb3SSatish Balay {
20072d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
2008dfbe8321SBarry Smith   PetscErrorCode ierr;
2009690b6cddSBarry Smith   PetscInt       i,j,k,n,row,bs,*ai,*aj,ambs,bs2;
201087828ca2SBarry Smith   PetscScalar    *x,zero = 0.0;
20113f1db9ecSBarry Smith   MatScalar      *aa,*aa_j;
20122d61bbb3SSatish Balay 
20132d61bbb3SSatish Balay   PetscFunctionBegin;
2014e32f2f54SBarry Smith   if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
2015d0f46423SBarry Smith   bs   = A->rmap->bs;
20162d61bbb3SSatish Balay   aa   = a->a;
20172d61bbb3SSatish Balay   ai   = a->i;
20182d61bbb3SSatish Balay   aj   = a->j;
20192d61bbb3SSatish Balay   ambs = a->mbs;
20202d61bbb3SSatish Balay   bs2  = a->bs2;
20212d61bbb3SSatish Balay 
20222dcb1b2aSMatthew Knepley   ierr = VecSet(v,zero);CHKERRQ(ierr);
20231ebc52fbSHong Zhang   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
2024e1311b90SBarry Smith   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
2025e32f2f54SBarry Smith   if (n != A->rmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
20262d61bbb3SSatish Balay   for (i=0; i<ambs; i++) {
20272d61bbb3SSatish Balay     for (j=ai[i]; j<ai[i+1]; j++) {
20282d61bbb3SSatish Balay       if (aj[j] == i) {
20292d61bbb3SSatish Balay         row  = i*bs;
20302d61bbb3SSatish Balay         aa_j = aa+j*bs2;
20312d61bbb3SSatish Balay         for (k=0; k<bs2; k+=(bs+1),row++) x[row] = aa_j[k];
20322d61bbb3SSatish Balay         break;
20332d61bbb3SSatish Balay       }
20342d61bbb3SSatish Balay     }
20352d61bbb3SSatish Balay   }
20361ebc52fbSHong Zhang   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
20372d61bbb3SSatish Balay   PetscFunctionReturn(0);
20382d61bbb3SSatish Balay }
20392d61bbb3SSatish Balay 
20404a2ae208SSatish Balay #undef __FUNCT__
20414a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_SeqBAIJ"
2042dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_SeqBAIJ(Mat A,Vec ll,Vec rr)
20432d61bbb3SSatish Balay {
20442d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
204553ef36baSBarry Smith   const PetscScalar *l,*r,*li,*ri;
204653ef36baSBarry Smith   PetscScalar       x;
20473f1db9ecSBarry Smith   MatScalar         *aa, *v;
2048dfbe8321SBarry Smith   PetscErrorCode    ierr;
204953ef36baSBarry Smith   PetscInt          i,j,k,lm,rn,M,m,n,mbs,tmp,bs,bs2,iai;
205053ef36baSBarry Smith   const PetscInt    *ai,*aj;
20512d61bbb3SSatish Balay 
20522d61bbb3SSatish Balay   PetscFunctionBegin;
20532d61bbb3SSatish Balay   ai  = a->i;
20542d61bbb3SSatish Balay   aj  = a->j;
20552d61bbb3SSatish Balay   aa  = a->a;
2056d0f46423SBarry Smith   m   = A->rmap->n;
2057d0f46423SBarry Smith   n   = A->cmap->n;
2058d0f46423SBarry Smith   bs  = A->rmap->bs;
20592d61bbb3SSatish Balay   mbs = a->mbs;
20602d61bbb3SSatish Balay   bs2 = a->bs2;
20612d61bbb3SSatish Balay   if (ll) {
206253ef36baSBarry Smith     ierr = VecGetArrayRead(ll,&l);CHKERRQ(ierr);
2063e1311b90SBarry Smith     ierr = VecGetLocalSize(ll,&lm);CHKERRQ(ierr);
2064e32f2f54SBarry Smith     if (lm != m) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Left scaling vector wrong length");
20652d61bbb3SSatish Balay     for (i=0; i<mbs; i++) { /* for each block row */
20662d61bbb3SSatish Balay       M  = ai[i+1] - ai[i];
20672d61bbb3SSatish Balay       li = l + i*bs;
20682d61bbb3SSatish Balay       v  = aa + bs2*ai[i];
20692d61bbb3SSatish Balay       for (j=0; j<M; j++) { /* for each block */
20702d61bbb3SSatish Balay         for (k=0; k<bs2; k++) {
20712d61bbb3SSatish Balay           (*v++) *= li[k%bs];
20722d61bbb3SSatish Balay         }
20732d61bbb3SSatish Balay       }
20742d61bbb3SSatish Balay     }
207553ef36baSBarry Smith     ierr = VecRestoreArrayRead(ll,&l);CHKERRQ(ierr);
2076efee365bSSatish Balay     ierr = PetscLogFlops(a->nz);CHKERRQ(ierr);
20772d61bbb3SSatish Balay   }
20782d61bbb3SSatish Balay 
20792d61bbb3SSatish Balay   if (rr) {
208053ef36baSBarry Smith     ierr = VecGetArrayRead(rr,&r);CHKERRQ(ierr);
2081e1311b90SBarry Smith     ierr = VecGetLocalSize(rr,&rn);CHKERRQ(ierr);
2082e32f2f54SBarry Smith     if (rn != n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Right scaling vector wrong length");
20832d61bbb3SSatish Balay     for (i=0; i<mbs; i++) { /* for each block row */
208453ef36baSBarry Smith       iai = ai[i];
208553ef36baSBarry Smith       M   = ai[i+1] - iai;
208653ef36baSBarry Smith       v   = aa + bs2*iai;
20872d61bbb3SSatish Balay       for (j=0; j<M; j++) { /* for each block */
208853ef36baSBarry Smith         ri = r + bs*aj[iai+j];
20892d61bbb3SSatish Balay         for (k=0; k<bs; k++) {
20902d61bbb3SSatish Balay           x = ri[k];
209153ef36baSBarry Smith           for (tmp=0; tmp<bs; tmp++) v[tmp] *= x;
209253ef36baSBarry Smith           v += bs;
20932d61bbb3SSatish Balay         }
20942d61bbb3SSatish Balay       }
20952d61bbb3SSatish Balay     }
209653ef36baSBarry Smith     ierr = VecRestoreArrayRead(rr,&r);CHKERRQ(ierr);
2097efee365bSSatish Balay     ierr = PetscLogFlops(a->nz);CHKERRQ(ierr);
20982d61bbb3SSatish Balay   }
20992d61bbb3SSatish Balay   PetscFunctionReturn(0);
21002d61bbb3SSatish Balay }
21012d61bbb3SSatish Balay 
21022d61bbb3SSatish Balay 
21034a2ae208SSatish Balay #undef __FUNCT__
21044a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_SeqBAIJ"
2105dfbe8321SBarry Smith PetscErrorCode MatGetInfo_SeqBAIJ(Mat A,MatInfoType flag,MatInfo *info)
21062d61bbb3SSatish Balay {
21072d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
21082d61bbb3SSatish Balay 
21092d61bbb3SSatish Balay   PetscFunctionBegin;
21102d61bbb3SSatish Balay   info->block_size   = a->bs2;
2111ceed8ce5SJed Brown   info->nz_allocated = a->bs2*a->maxnz;
21122d61bbb3SSatish Balay   info->nz_used      = a->bs2*a->nz;
21132d61bbb3SSatish Balay   info->nz_unneeded  = (double)(info->nz_allocated - info->nz_used);
21142d61bbb3SSatish Balay   info->assemblies   = A->num_ass;
21158e58a170SBarry Smith   info->mallocs      = A->info.mallocs;
21167adad957SLisandro Dalcin   info->memory       = ((PetscObject)A)->mem;
2117d5f3da31SBarry Smith   if (A->factortype) {
21182d61bbb3SSatish Balay     info->fill_ratio_given  = A->info.fill_ratio_given;
21192d61bbb3SSatish Balay     info->fill_ratio_needed = A->info.fill_ratio_needed;
21202d61bbb3SSatish Balay     info->factor_mallocs    = A->info.factor_mallocs;
21212d61bbb3SSatish Balay   } else {
21222d61bbb3SSatish Balay     info->fill_ratio_given  = 0;
21232d61bbb3SSatish Balay     info->fill_ratio_needed = 0;
21242d61bbb3SSatish Balay     info->factor_mallocs    = 0;
21252d61bbb3SSatish Balay   }
21262d61bbb3SSatish Balay   PetscFunctionReturn(0);
21272d61bbb3SSatish Balay }
21282d61bbb3SSatish Balay 
21292d61bbb3SSatish Balay 
21304a2ae208SSatish Balay #undef __FUNCT__
21314a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_SeqBAIJ"
2132dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_SeqBAIJ(Mat A)
21332d61bbb3SSatish Balay {
21342d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
2135dfbe8321SBarry Smith   PetscErrorCode ierr;
21362d61bbb3SSatish Balay 
21372d61bbb3SSatish Balay   PetscFunctionBegin;
2138549d3d68SSatish Balay   ierr = PetscMemzero(a->a,a->bs2*a->i[a->mbs]*sizeof(MatScalar));CHKERRQ(ierr);
21392d61bbb3SSatish Balay   PetscFunctionReturn(0);
21402d61bbb3SSatish Balay }
2141