xref: /petsc/src/mat/impls/baij/seq/baij2.c (revision ef16e67128e275b9fd6dc747c1e3ef0ca8beb92f)
1cac129eeSSatish Balay 
2c6db04a5SJed Brown #include <../src/mat/impls/baij/seq/baij.h>
3af0996ceSBarry 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);
27854ce69bSBarry Smith   ierr = PetscMalloc1(m+1,&nidx);CHKERRQ(ierr);
28854ce69bSBarry Smith   ierr = PetscMalloc1(A->rmap->N+1,&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;
866041f1b1SToby Isaac   PetscBool      flag;
87736121d4SSatish Balay 
883a40ed3dSBarry Smith   PetscFunctionBegin;
89736121d4SSatish Balay   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
90218c64b6SSatish Balay   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
91b9b97703SBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
92b9b97703SBarry Smith   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
93736121d4SSatish Balay 
94854ce69bSBarry Smith   ierr  = PetscCalloc1(1+oldcols,&smap);CHKERRQ(ierr);
95736121d4SSatish Balay   ssmap = smap;
96854ce69bSBarry Smith   ierr  = PetscMalloc1(1+nrows,&lens);CHKERRQ(ierr);
97736121d4SSatish Balay   for (i=0; i<ncols; i++) smap[icol[i]] = i+1;
98736121d4SSatish Balay   /* determine lens of each row */
99736121d4SSatish Balay   for (i=0; i<nrows; i++) {
100736121d4SSatish Balay     kstart  = ai[irow[i]];
101736121d4SSatish Balay     kend    = kstart + a->ilen[irow[i]];
102736121d4SSatish Balay     lens[i] = 0;
103736121d4SSatish Balay     for (k=kstart; k<kend; k++) {
10426fbe8dcSKarl Rupp       if (ssmap[aj[k]]) lens[i]++;
105736121d4SSatish Balay     }
106736121d4SSatish Balay   }
107736121d4SSatish Balay   /* Create and fill new matrix */
108736121d4SSatish Balay   if (scall == MAT_REUSE_MATRIX) {
109736121d4SSatish Balay     c = (Mat_SeqBAIJ*)((*B)->data);
110736121d4SSatish Balay 
111e32f2f54SBarry Smith     if (c->mbs!=nrows || c->nbs!=ncols || (*B)->rmap->bs!=bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Submatrix wrong size");
112690b6cddSBarry Smith     ierr = PetscMemcmp(c->ilen,lens,c->mbs *sizeof(PetscInt),&flag);CHKERRQ(ierr);
113e7e72b3dSBarry Smith     if (!flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong no of nonzeros");
114690b6cddSBarry Smith     ierr = PetscMemzero(c->ilen,c->mbs*sizeof(PetscInt));CHKERRQ(ierr);
115736121d4SSatish Balay     C    = *B;
1163a40ed3dSBarry Smith   } else {
117ce94432eSBarry Smith     ierr = MatCreate(PetscObjectComm((PetscObject)A),&C);CHKERRQ(ierr);
118f69a0ea3SMatthew Knepley     ierr = MatSetSizes(C,nrows*bs,ncols*bs,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1197adad957SLisandro Dalcin     ierr = MatSetType(C,((PetscObject)A)->type_name);CHKERRQ(ierr);
120e50a8114SHong Zhang     ierr = MatSeqBAIJSetPreallocation_SeqBAIJ(C,bs,0,lens);CHKERRQ(ierr);
121736121d4SSatish Balay   }
122736121d4SSatish Balay   c = (Mat_SeqBAIJ*)(C->data);
123736121d4SSatish Balay   for (i=0; i<nrows; i++) {
124736121d4SSatish Balay     row      = irow[i];
125736121d4SSatish Balay     kstart   = ai[row];
126736121d4SSatish Balay     kend     = kstart + a->ilen[row];
127736121d4SSatish Balay     mat_i    = c->i[i];
128736121d4SSatish Balay     mat_j    = c->j + mat_i;
129218c64b6SSatish Balay     mat_a    = c->a + mat_i*bs2;
130736121d4SSatish Balay     mat_ilen = c->ilen + i;
131736121d4SSatish Balay     for (k=kstart; k<kend; k++) {
132736121d4SSatish Balay       if ((tcol=ssmap[a->j[k]])) {
133736121d4SSatish Balay         *mat_j++ = tcol - 1;
134549d3d68SSatish Balay         ierr     = PetscMemcpy(mat_a,a->a+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
135549d3d68SSatish Balay         mat_a   += bs2;
136736121d4SSatish Balay         (*mat_ilen)++;
137736121d4SSatish Balay       }
138736121d4SSatish Balay     }
139736121d4SSatish Balay   }
140cdc6f3adSToby Isaac   /* sort */
141cdc6f3adSToby Isaac   {
142cdc6f3adSToby Isaac     MatScalar *work;
143cdc6f3adSToby Isaac     ierr = PetscMalloc1(bs2,&work);CHKERRQ(ierr);
144cdc6f3adSToby Isaac     for (i=0; i<nrows; i++) {
145cdc6f3adSToby Isaac       PetscInt ilen;
146cdc6f3adSToby Isaac       mat_i = c->i[i];
147cdc6f3adSToby Isaac       mat_j = c->j + mat_i;
148cdc6f3adSToby Isaac       mat_a = c->a + mat_i*bs2;
149cdc6f3adSToby Isaac       ilen  = c->ilen[i];
150cdc6f3adSToby Isaac       ierr  = PetscSortIntWithDataArray(ilen,mat_j,mat_a,bs2*sizeof(MatScalar),work);CHKERRQ(ierr);
151cdc6f3adSToby Isaac     }
152cdc6f3adSToby Isaac     ierr = PetscFree(work);CHKERRQ(ierr);
153cdc6f3adSToby Isaac   }
154218c64b6SSatish Balay 
155736121d4SSatish Balay   /* Free work space */
156736121d4SSatish Balay   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
157606d414cSSatish Balay   ierr = PetscFree(smap);CHKERRQ(ierr);
158606d414cSSatish Balay   ierr = PetscFree(lens);CHKERRQ(ierr);
1596d4a8577SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1606d4a8577SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
161736121d4SSatish Balay 
162736121d4SSatish Balay   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
163736121d4SSatish Balay   *B   = C;
1643a40ed3dSBarry Smith   PetscFunctionReturn(0);
165736121d4SSatish Balay }
166736121d4SSatish Balay 
1674a2ae208SSatish Balay #undef __FUNCT__
1684a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqBAIJ"
1694aa3045dSJed Brown PetscErrorCode MatGetSubMatrix_SeqBAIJ(Mat A,IS isrow,IS iscol,MatReuse scall,Mat *B)
170218c64b6SSatish Balay {
171218c64b6SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
172218c64b6SSatish Balay   IS             is1,is2;
1736849ba73SBarry Smith   PetscErrorCode ierr;
174afebec48SHong Zhang   PetscInt       *vary,*iary,nrows,ncols,i,bs=A->rmap->bs,count,maxmnbs,j;
1755d0c19d7SBarry Smith   const PetscInt *irow,*icol;
176218c64b6SSatish Balay 
1773a40ed3dSBarry Smith   PetscFunctionBegin;
178218c64b6SSatish Balay   ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr);
179218c64b6SSatish Balay   ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr);
180b9b97703SBarry Smith   ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr);
181b9b97703SBarry Smith   ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr);
182218c64b6SSatish Balay 
183218c64b6SSatish Balay   /* Verify if the indices corespond to each element in a block
184218c64b6SSatish Balay    and form the IS with compressed IS */
185f8ecb639SStefano Zampini   maxmnbs = PetscMax(a->mbs,a->nbs);
186f8ecb639SStefano Zampini   ierr = PetscMalloc2(maxmnbs,&vary,maxmnbs,&iary);CHKERRQ(ierr);
187fca92195SBarry Smith   ierr = PetscMemzero(vary,a->mbs*sizeof(PetscInt));CHKERRQ(ierr);
188218c64b6SSatish Balay   for (i=0; i<nrows; i++) vary[irow[i]/bs]++;
189218c64b6SSatish Balay   for (i=0; i<a->mbs; i++) {
190e32f2f54SBarry Smith     if (vary[i]!=0 && vary[i]!=bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Index set does not match blocks");
1916041f1b1SToby Isaac   }
1926041f1b1SToby Isaac   count = 0;
1936041f1b1SToby Isaac   for (i=0; i<nrows; i++) {
194afebec48SHong Zhang     j = irow[i] / bs;
1956041f1b1SToby Isaac     if ((vary[j]--)==bs) iary[count++] = j;
196218c64b6SSatish Balay   }
19770b3c8c7SBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,count,iary,PETSC_COPY_VALUES,&is1);CHKERRQ(ierr);
198218c64b6SSatish Balay 
199f8ecb639SStefano Zampini   ierr = PetscMemzero(vary,(a->nbs)*sizeof(PetscInt));CHKERRQ(ierr);
200218c64b6SSatish Balay   for (i=0; i<ncols; i++) vary[icol[i]/bs]++;
201f8ecb639SStefano Zampini   for (i=0; i<a->nbs; i++) {
202e32f2f54SBarry Smith     if (vary[i]!=0 && vary[i]!=bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Internal error in PETSc");
203218c64b6SSatish Balay   }
2046041f1b1SToby Isaac   count = 0;
2056041f1b1SToby Isaac   for (i=0; i<ncols; i++) {
206afebec48SHong Zhang     j = icol[i] / bs;
2076041f1b1SToby Isaac     if ((vary[j]--)==bs) iary[count++] = j;
2086041f1b1SToby Isaac   }
20970b3c8c7SBarry Smith   ierr = ISCreateGeneral(PETSC_COMM_SELF,count,iary,PETSC_COPY_VALUES,&is2);CHKERRQ(ierr);
210218c64b6SSatish Balay   ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr);
211218c64b6SSatish Balay   ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr);
212fca92195SBarry Smith   ierr = PetscFree2(vary,iary);CHKERRQ(ierr);
213218c64b6SSatish Balay 
2144aa3045dSJed Brown   ierr = MatGetSubMatrix_SeqBAIJ_Private(A,is1,is2,scall,B);CHKERRQ(ierr);
2156bf464f9SBarry Smith   ierr = ISDestroy(&is1);CHKERRQ(ierr);
2166bf464f9SBarry Smith   ierr = ISDestroy(&is2);CHKERRQ(ierr);
2173a40ed3dSBarry Smith   PetscFunctionReturn(0);
218218c64b6SSatish Balay }
219218c64b6SSatish Balay 
2204a2ae208SSatish Balay #undef __FUNCT__
2214a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_SeqBAIJ"
222690b6cddSBarry Smith PetscErrorCode MatGetSubMatrices_SeqBAIJ(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[])
223736121d4SSatish Balay {
2246849ba73SBarry Smith   PetscErrorCode ierr;
225690b6cddSBarry Smith   PetscInt       i;
226736121d4SSatish Balay 
2273a40ed3dSBarry Smith   PetscFunctionBegin;
228736121d4SSatish Balay   if (scall == MAT_INITIAL_MATRIX) {
229854ce69bSBarry Smith     ierr = PetscMalloc1(n+1,B);CHKERRQ(ierr);
230736121d4SSatish Balay   }
231736121d4SSatish Balay 
232736121d4SSatish Balay   for (i=0; i<n; i++) {
2334aa3045dSJed Brown     ierr = MatGetSubMatrix_SeqBAIJ(A,irow[i],icol[i],scall,&(*B)[i]);CHKERRQ(ierr);
234736121d4SSatish Balay   }
2353a40ed3dSBarry Smith   PetscFunctionReturn(0);
236736121d4SSatish Balay }
237218c64b6SSatish Balay 
238218c64b6SSatish Balay 
2392d61bbb3SSatish Balay /* -------------------------------------------------------*/
2402d61bbb3SSatish Balay /* Should check that shapes of vectors and matrices match */
2412d61bbb3SSatish Balay /* -------------------------------------------------------*/
2422d61bbb3SSatish Balay 
2434a2ae208SSatish Balay #undef __FUNCT__
2444a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_1"
245dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_1(Mat A,Vec xx,Vec zz)
2462d61bbb3SSatish Balay {
2472d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
248d9fead3dSBarry Smith   PetscScalar       *z,sum;
249d9fead3dSBarry Smith   const PetscScalar *x;
250d9fead3dSBarry Smith   const MatScalar   *v;
2516849ba73SBarry Smith   PetscErrorCode    ierr;
2527c565772SBarry Smith   PetscInt          mbs,i,n;
2530298fd71SBarry Smith   const PetscInt    *idx,*ii,*ridx=NULL;
254ace3abfcSBarry Smith   PetscBool         usecprow=a->compressedrow.use;
2552d61bbb3SSatish Balay 
2562d61bbb3SSatish Balay   PetscFunctionBegin;
2573649974fSBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
2581ebc52fbSHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
2592d61bbb3SSatish Balay 
26026e093fcSHong Zhang   if (usecprow) {
26126e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
26226e093fcSHong Zhang     ii   = a->compressedrow.i;
2637b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
264*ef16e671SStefano Zampini     ierr = PetscMemzero(z,a->mbs*sizeof(PetscScalar));CHKERRQ(ierr);
26526e093fcSHong Zhang   } else {
26626e093fcSHong Zhang     mbs = a->mbs;
2672d61bbb3SSatish Balay     ii  = a->i;
26826e093fcSHong Zhang   }
2692d61bbb3SSatish Balay 
2702d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
271ee54c7eeSHong Zhang     n   = ii[1] - ii[0];
272ee54c7eeSHong Zhang     v   = a->a + ii[0];
273ee54c7eeSHong Zhang     idx = a->j + ii[0];
274ee54c7eeSHong Zhang     ii++;
275444d8c10SJed Brown     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);   /* Indices for the next row (assumes same size as this one) */
276444d8c10SJed Brown     PetscPrefetchBlock(v+1*n,1*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
2772d61bbb3SSatish Balay     sum = 0.0;
2782162cab8SBarry Smith     PetscSparseDensePlusDot(sum,x,v,idx,n);
27926e093fcSHong Zhang     if (usecprow) {
2807b2bb3b9SHong Zhang       z[ridx[i]] = sum;
28126e093fcSHong Zhang     } else {
2822d61bbb3SSatish Balay       z[i]       = sum;
2832d61bbb3SSatish Balay     }
28426e093fcSHong Zhang   }
2853649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
2861ebc52fbSHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
2877c565772SBarry Smith   ierr = PetscLogFlops(2.0*a->nz - a->nonzerorowcnt);CHKERRQ(ierr);
2882d61bbb3SSatish Balay   PetscFunctionReturn(0);
2892d61bbb3SSatish Balay }
2902d61bbb3SSatish Balay 
2914a2ae208SSatish Balay #undef __FUNCT__
2924a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_2"
293dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_2(Mat A,Vec xx,Vec zz)
2942d61bbb3SSatish Balay {
2952d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
296d9fead3dSBarry Smith   PetscScalar       *z = 0,sum1,sum2,*zarray;
297d9fead3dSBarry Smith   const PetscScalar *x,*xb;
29887828ca2SBarry Smith   PetscScalar       x1,x2;
299d9fead3dSBarry Smith   const MatScalar   *v;
300dfbe8321SBarry Smith   PetscErrorCode    ierr;
3017c565772SBarry Smith   PetscInt          mbs,i,*idx,*ii,j,n,*ridx=NULL;
302ace3abfcSBarry Smith   PetscBool         usecprow=a->compressedrow.use;
3032d61bbb3SSatish Balay 
3042d61bbb3SSatish Balay   PetscFunctionBegin;
3053649974fSBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
30626e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
3072d61bbb3SSatish Balay 
3082d61bbb3SSatish Balay   idx = a->j;
3092d61bbb3SSatish Balay   v   = a->a;
31026e093fcSHong Zhang   if (usecprow) {
31126e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
31226e093fcSHong Zhang     ii   = a->compressedrow.i;
3137b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
314*ef16e671SStefano Zampini     ierr = PetscMemzero(zarray,2*a->mbs*sizeof(PetscScalar));CHKERRQ(ierr);
31526e093fcSHong Zhang   } else {
31626e093fcSHong Zhang     mbs = a->mbs;
3172d61bbb3SSatish Balay     ii  = a->i;
31826e093fcSHong Zhang     z   = zarray;
31926e093fcSHong Zhang   }
3202d61bbb3SSatish Balay 
3212d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
3222d61bbb3SSatish Balay     n           = ii[1] - ii[0]; ii++;
3232d61bbb3SSatish Balay     sum1        = 0.0; sum2 = 0.0;
324444d8c10SJed Brown     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);   /* Indices for the next row (assumes same size as this one) */
325444d8c10SJed Brown     PetscPrefetchBlock(v+4*n,4*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
3262d61bbb3SSatish Balay     for (j=0; j<n; j++) {
3272d61bbb3SSatish Balay       xb    = x + 2*(*idx++); x1 = xb[0]; x2 = xb[1];
3282d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[2]*x2;
3292d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[3]*x2;
3302d61bbb3SSatish Balay       v    += 4;
3312d61bbb3SSatish Balay     }
3327b2bb3b9SHong Zhang     if (usecprow) z = zarray + 2*ridx[i];
3332d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2;
33426e093fcSHong Zhang     if (!usecprow) z += 2;
3352d61bbb3SSatish Balay   }
3363649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
33726e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
3387c565772SBarry Smith   ierr = PetscLogFlops(8.0*a->nz - 2.0*a->nonzerorowcnt);CHKERRQ(ierr);
3392d61bbb3SSatish Balay   PetscFunctionReturn(0);
3402d61bbb3SSatish Balay }
3412d61bbb3SSatish Balay 
3424a2ae208SSatish Balay #undef __FUNCT__
3434a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_3"
344dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_3(Mat A,Vec xx,Vec zz)
3452d61bbb3SSatish Balay {
3462d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
347d9fead3dSBarry Smith   PetscScalar       *z = 0,sum1,sum2,sum3,x1,x2,x3,*zarray;
348d9fead3dSBarry Smith   const PetscScalar *x,*xb;
349d9fead3dSBarry Smith   const MatScalar   *v;
350dfbe8321SBarry Smith   PetscErrorCode    ierr;
3517c565772SBarry Smith   PetscInt          mbs,i,*idx,*ii,j,n,*ridx=NULL;
352ace3abfcSBarry Smith   PetscBool         usecprow=a->compressedrow.use;
35326e093fcSHong Zhang 
354b6410449SSatish Balay #if defined(PETSC_HAVE_PRAGMA_DISJOINT)
355fee21e36SBarry Smith #pragma disjoint(*v,*z,*xb)
356fee21e36SBarry Smith #endif
357fee21e36SBarry Smith 
3582d61bbb3SSatish Balay   PetscFunctionBegin;
3593649974fSBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
36026e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
3612d61bbb3SSatish Balay 
3622d61bbb3SSatish Balay   idx = a->j;
3632d61bbb3SSatish Balay   v   = a->a;
36426e093fcSHong Zhang   if (usecprow) {
36526e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
36626e093fcSHong Zhang     ii   = a->compressedrow.i;
3677b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
368*ef16e671SStefano Zampini     ierr = PetscMemzero(zarray,3*a->mbs*sizeof(PetscScalar));CHKERRQ(ierr);
36926e093fcSHong Zhang   } else {
37026e093fcSHong Zhang     mbs = a->mbs;
3712d61bbb3SSatish Balay     ii  = a->i;
37226e093fcSHong Zhang     z   = zarray;
37326e093fcSHong Zhang   }
3742d61bbb3SSatish Balay 
3752d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
3762d61bbb3SSatish Balay     n           = ii[1] - ii[0]; ii++;
3772d61bbb3SSatish Balay     sum1        = 0.0; sum2 = 0.0; sum3 = 0.0;
378444d8c10SJed Brown     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);   /* Indices for the next row (assumes same size as this one) */
379444d8c10SJed Brown     PetscPrefetchBlock(v+9*n,9*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
3802d61bbb3SSatish Balay     for (j=0; j<n; j++) {
38126fbe8dcSKarl Rupp       xb = x + 3*(*idx++);
38226fbe8dcSKarl Rupp       x1 = xb[0];
38326fbe8dcSKarl Rupp       x2 = xb[1];
38426fbe8dcSKarl Rupp       x3 = xb[2];
38526fbe8dcSKarl Rupp 
3862d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[3]*x2 + v[6]*x3;
3872d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[4]*x2 + v[7]*x3;
3882d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[5]*x2 + v[8]*x3;
3892d61bbb3SSatish Balay       v    += 9;
3902d61bbb3SSatish Balay     }
3917b2bb3b9SHong Zhang     if (usecprow) z = zarray + 3*ridx[i];
3922d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3;
39326e093fcSHong Zhang     if (!usecprow) z += 3;
3942d61bbb3SSatish Balay   }
3953649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
39626e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
3977c565772SBarry Smith   ierr = PetscLogFlops(18.0*a->nz - 3.0*a->nonzerorowcnt);CHKERRQ(ierr);
3982d61bbb3SSatish Balay   PetscFunctionReturn(0);
3992d61bbb3SSatish Balay }
4002d61bbb3SSatish Balay 
4014a2ae208SSatish Balay #undef __FUNCT__
4024a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_4"
403dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_4(Mat A,Vec xx,Vec zz)
4042d61bbb3SSatish Balay {
4052d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
406d9fead3dSBarry Smith   PetscScalar       *z = 0,sum1,sum2,sum3,sum4,x1,x2,x3,x4,*zarray;
407d9fead3dSBarry Smith   const PetscScalar *x,*xb;
408d9fead3dSBarry Smith   const MatScalar   *v;
409dfbe8321SBarry Smith   PetscErrorCode    ierr;
4107c565772SBarry Smith   PetscInt          mbs,i,*idx,*ii,j,n,*ridx=NULL;
411ace3abfcSBarry Smith   PetscBool         usecprow=a->compressedrow.use;
4122d61bbb3SSatish Balay 
4132d61bbb3SSatish Balay   PetscFunctionBegin;
4143649974fSBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
41526e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
4162d61bbb3SSatish Balay 
4172d61bbb3SSatish Balay   idx = a->j;
4182d61bbb3SSatish Balay   v   = a->a;
41926e093fcSHong Zhang   if (usecprow) {
42026e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
42126e093fcSHong Zhang     ii   = a->compressedrow.i;
4227b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
423*ef16e671SStefano Zampini     ierr = PetscMemzero(zarray,4*a->mbs*sizeof(PetscScalar));CHKERRQ(ierr);
42426e093fcSHong Zhang   } else {
42526e093fcSHong Zhang     mbs = a->mbs;
4262d61bbb3SSatish Balay     ii  = a->i;
42726e093fcSHong Zhang     z   = zarray;
42826e093fcSHong Zhang   }
4292d61bbb3SSatish Balay 
4302d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
43126fbe8dcSKarl Rupp     n = ii[1] - ii[0];
43226fbe8dcSKarl Rupp     ii++;
43326fbe8dcSKarl Rupp     sum1 = 0.0;
43426fbe8dcSKarl Rupp     sum2 = 0.0;
43526fbe8dcSKarl Rupp     sum3 = 0.0;
43626fbe8dcSKarl Rupp     sum4 = 0.0;
43726fbe8dcSKarl Rupp 
438444d8c10SJed Brown     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);     /* Indices for the next row (assumes same size as this one) */
439444d8c10SJed Brown     PetscPrefetchBlock(v+16*n,16*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
4402d61bbb3SSatish Balay     for (j=0; j<n; j++) {
4412d61bbb3SSatish Balay       xb    = x + 4*(*idx++);
4422d61bbb3SSatish Balay       x1    = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
4432d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[4]*x2 + v[8]*x3   + v[12]*x4;
4442d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[5]*x2 + v[9]*x3   + v[13]*x4;
4452d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[6]*x2 + v[10]*x3  + v[14]*x4;
4462d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[7]*x2 + v[11]*x3  + v[15]*x4;
4472d61bbb3SSatish Balay       v    += 16;
4482d61bbb3SSatish Balay     }
4497b2bb3b9SHong Zhang     if (usecprow) z = zarray + 4*ridx[i];
4502d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4;
45126e093fcSHong Zhang     if (!usecprow) z += 4;
4522d61bbb3SSatish Balay   }
4533649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
45426e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
4557c565772SBarry Smith   ierr = PetscLogFlops(32.0*a->nz - 4.0*a->nonzerorowcnt);CHKERRQ(ierr);
4562d61bbb3SSatish Balay   PetscFunctionReturn(0);
4572d61bbb3SSatish Balay }
4582d61bbb3SSatish Balay 
4594a2ae208SSatish Balay #undef __FUNCT__
4604a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_5"
461dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_5(Mat A,Vec xx,Vec zz)
4622d61bbb3SSatish Balay {
4632d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
464d9fead3dSBarry Smith   PetscScalar       sum1,sum2,sum3,sum4,sum5,x1,x2,x3,x4,x5,*z = 0,*zarray;
465d9fead3dSBarry Smith   const PetscScalar *xb,*x;
466d9fead3dSBarry Smith   const MatScalar   *v;
467dfbe8321SBarry Smith   PetscErrorCode    ierr;
4680298fd71SBarry Smith   const PetscInt    *idx,*ii,*ridx=NULL;
4697c565772SBarry Smith   PetscInt          mbs,i,j,n;
470ace3abfcSBarry Smith   PetscBool         usecprow=a->compressedrow.use;
4712d61bbb3SSatish Balay 
472433994e6SBarry Smith   PetscFunctionBegin;
4733649974fSBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
47426e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
4752d61bbb3SSatish Balay 
4762d61bbb3SSatish Balay   idx = a->j;
4772d61bbb3SSatish Balay   v   = a->a;
47826e093fcSHong Zhang   if (usecprow) {
47926e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
48026e093fcSHong Zhang     ii   = a->compressedrow.i;
4817b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
482*ef16e671SStefano Zampini     ierr = PetscMemzero(zarray,5*a->mbs*sizeof(PetscScalar));CHKERRQ(ierr);
48326e093fcSHong Zhang   } else {
48426e093fcSHong Zhang     mbs = a->mbs;
4852d61bbb3SSatish Balay     ii  = a->i;
48626e093fcSHong Zhang     z   = zarray;
48726e093fcSHong Zhang   }
4882d61bbb3SSatish Balay 
4892d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
4902d61bbb3SSatish Balay     n           = ii[1] - ii[0]; ii++;
4912d61bbb3SSatish Balay     sum1        = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0;
492444d8c10SJed Brown     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);     /* Indices for the next row (assumes same size as this one) */
493444d8c10SJed Brown     PetscPrefetchBlock(v+25*n,25*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
4942d61bbb3SSatish Balay     for (j=0; j<n; j++) {
4952d61bbb3SSatish Balay       xb    = x + 5*(*idx++);
4962d61bbb3SSatish Balay       x1    = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4];
4972d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[5]*x2 + v[10]*x3  + v[15]*x4 + v[20]*x5;
4982d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[6]*x2 + v[11]*x3  + v[16]*x4 + v[21]*x5;
4992d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[7]*x2 + v[12]*x3  + v[17]*x4 + v[22]*x5;
5002d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[8]*x2 + v[13]*x3  + v[18]*x4 + v[23]*x5;
5012d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[9]*x2 + v[14]*x3  + v[19]*x4 + v[24]*x5;
5022d61bbb3SSatish Balay       v    += 25;
5032d61bbb3SSatish Balay     }
5047b2bb3b9SHong Zhang     if (usecprow) z = zarray + 5*ridx[i];
5052d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5;
50626e093fcSHong Zhang     if (!usecprow) z += 5;
5072d61bbb3SSatish Balay   }
5083649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
50926e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
5107c565772SBarry Smith   ierr = PetscLogFlops(50.0*a->nz - 5.0*a->nonzerorowcnt);CHKERRQ(ierr);
5112d61bbb3SSatish Balay   PetscFunctionReturn(0);
5122d61bbb3SSatish Balay }
5132d61bbb3SSatish Balay 
51415091d37SBarry Smith 
515d6232bceSMichael Lange 
5164a2ae208SSatish Balay #undef __FUNCT__
5174a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_6"
518dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_6(Mat A,Vec xx,Vec zz)
51915091d37SBarry Smith {
52015091d37SBarry Smith   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
521d9fead3dSBarry Smith   PetscScalar       *z = 0,sum1,sum2,sum3,sum4,sum5,sum6;
522d9fead3dSBarry Smith   const PetscScalar *x,*xb;
52326e093fcSHong Zhang   PetscScalar       x1,x2,x3,x4,x5,x6,*zarray;
524d9fead3dSBarry Smith   const MatScalar   *v;
525dfbe8321SBarry Smith   PetscErrorCode    ierr;
5267c565772SBarry Smith   PetscInt          mbs,i,*idx,*ii,j,n,*ridx=NULL;
527ace3abfcSBarry Smith   PetscBool         usecprow=a->compressedrow.use;
52815091d37SBarry Smith 
529433994e6SBarry Smith   PetscFunctionBegin;
5303649974fSBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
53126e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
53215091d37SBarry Smith 
53315091d37SBarry Smith   idx = a->j;
53415091d37SBarry Smith   v   = a->a;
53526e093fcSHong Zhang   if (usecprow) {
53626e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
53726e093fcSHong Zhang     ii   = a->compressedrow.i;
5387b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
539*ef16e671SStefano Zampini     ierr = PetscMemzero(zarray,6*a->mbs*sizeof(PetscScalar));CHKERRQ(ierr);
54026e093fcSHong Zhang   } else {
54126e093fcSHong Zhang     mbs = a->mbs;
54215091d37SBarry Smith     ii  = a->i;
54326e093fcSHong Zhang     z   = zarray;
54426e093fcSHong Zhang   }
54515091d37SBarry Smith 
54615091d37SBarry Smith   for (i=0; i<mbs; i++) {
54726fbe8dcSKarl Rupp     n  = ii[1] - ii[0];
54826fbe8dcSKarl Rupp     ii++;
54926fbe8dcSKarl Rupp     sum1 = 0.0;
55026fbe8dcSKarl Rupp     sum2 = 0.0;
55126fbe8dcSKarl Rupp     sum3 = 0.0;
55226fbe8dcSKarl Rupp     sum4 = 0.0;
55326fbe8dcSKarl Rupp     sum5 = 0.0;
55426fbe8dcSKarl Rupp     sum6 = 0.0;
55526fbe8dcSKarl Rupp 
556444d8c10SJed Brown     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);     /* Indices for the next row (assumes same size as this one) */
557444d8c10SJed Brown     PetscPrefetchBlock(v+36*n,36*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
55815091d37SBarry Smith     for (j=0; j<n; j++) {
55915091d37SBarry Smith       xb    = x + 6*(*idx++);
56015091d37SBarry Smith       x1    = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5];
56115091d37SBarry Smith       sum1 += v[0]*x1 + v[6]*x2  + v[12]*x3  + v[18]*x4 + v[24]*x5 + v[30]*x6;
56215091d37SBarry Smith       sum2 += v[1]*x1 + v[7]*x2  + v[13]*x3  + v[19]*x4 + v[25]*x5 + v[31]*x6;
56315091d37SBarry Smith       sum3 += v[2]*x1 + v[8]*x2  + v[14]*x3  + v[20]*x4 + v[26]*x5 + v[32]*x6;
56415091d37SBarry Smith       sum4 += v[3]*x1 + v[9]*x2  + v[15]*x3  + v[21]*x4 + v[27]*x5 + v[33]*x6;
56515091d37SBarry Smith       sum5 += v[4]*x1 + v[10]*x2 + v[16]*x3  + v[22]*x4 + v[28]*x5 + v[34]*x6;
56615091d37SBarry Smith       sum6 += v[5]*x1 + v[11]*x2 + v[17]*x3  + v[23]*x4 + v[29]*x5 + v[35]*x6;
56715091d37SBarry Smith       v    += 36;
56815091d37SBarry Smith     }
5697b2bb3b9SHong Zhang     if (usecprow) z = zarray + 6*ridx[i];
57015091d37SBarry Smith     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6;
57126e093fcSHong Zhang     if (!usecprow) z += 6;
57215091d37SBarry Smith   }
57315091d37SBarry Smith 
5743649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
57526e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
5767c565772SBarry Smith   ierr = PetscLogFlops(72.0*a->nz - 6.0*a->nonzerorowcnt);CHKERRQ(ierr);
57715091d37SBarry Smith   PetscFunctionReturn(0);
57815091d37SBarry Smith }
5798ab949d8SShri Abhyankar 
5804a2ae208SSatish Balay #undef __FUNCT__
5814a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_7"
582dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_7(Mat A,Vec xx,Vec zz)
5832d61bbb3SSatish Balay {
5842d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
585d9fead3dSBarry Smith   PetscScalar       *z = 0,sum1,sum2,sum3,sum4,sum5,sum6,sum7;
586d9fead3dSBarry Smith   const PetscScalar *x,*xb;
58726e093fcSHong Zhang   PetscScalar       x1,x2,x3,x4,x5,x6,x7,*zarray;
588d9fead3dSBarry Smith   const MatScalar   *v;
589dfbe8321SBarry Smith   PetscErrorCode    ierr;
5907c565772SBarry Smith   PetscInt          mbs,i,*idx,*ii,j,n,*ridx=NULL;
591ace3abfcSBarry Smith   PetscBool         usecprow=a->compressedrow.use;
5922d61bbb3SSatish Balay 
593433994e6SBarry Smith   PetscFunctionBegin;
5943649974fSBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
59526e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
5962d61bbb3SSatish Balay 
5972d61bbb3SSatish Balay   idx = a->j;
5982d61bbb3SSatish Balay   v   = a->a;
59926e093fcSHong Zhang   if (usecprow) {
60026e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
60126e093fcSHong Zhang     ii   = a->compressedrow.i;
6027b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
603*ef16e671SStefano Zampini     ierr = PetscMemzero(zarray,7*a->mbs*sizeof(PetscScalar));CHKERRQ(ierr);
60426e093fcSHong Zhang   } else {
60526e093fcSHong Zhang     mbs = a->mbs;
6062d61bbb3SSatish Balay     ii  = a->i;
60726e093fcSHong Zhang     z   = zarray;
60826e093fcSHong Zhang   }
6092d61bbb3SSatish Balay 
6102d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
61126fbe8dcSKarl Rupp     n  = ii[1] - ii[0];
61226fbe8dcSKarl Rupp     ii++;
61326fbe8dcSKarl Rupp     sum1 = 0.0;
61426fbe8dcSKarl Rupp     sum2 = 0.0;
61526fbe8dcSKarl Rupp     sum3 = 0.0;
61626fbe8dcSKarl Rupp     sum4 = 0.0;
61726fbe8dcSKarl Rupp     sum5 = 0.0;
61826fbe8dcSKarl Rupp     sum6 = 0.0;
61926fbe8dcSKarl Rupp     sum7 = 0.0;
62026fbe8dcSKarl Rupp 
621444d8c10SJed Brown     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);     /* Indices for the next row (assumes same size as this one) */
622444d8c10SJed Brown     PetscPrefetchBlock(v+49*n,49*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
6232d61bbb3SSatish Balay     for (j=0; j<n; j++) {
6242d61bbb3SSatish Balay       xb    = x + 7*(*idx++);
6252d61bbb3SSatish Balay       x1    = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6];
6262d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[7]*x2  + v[14]*x3  + v[21]*x4 + v[28]*x5 + v[35]*x6 + v[42]*x7;
6272d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[8]*x2  + v[15]*x3  + v[22]*x4 + v[29]*x5 + v[36]*x6 + v[43]*x7;
6282d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[9]*x2  + v[16]*x3  + v[23]*x4 + v[30]*x5 + v[37]*x6 + v[44]*x7;
6292d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[10]*x2 + v[17]*x3  + v[24]*x4 + v[31]*x5 + v[38]*x6 + v[45]*x7;
6302d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[11]*x2 + v[18]*x3  + v[25]*x4 + v[32]*x5 + v[39]*x6 + v[46]*x7;
6312d61bbb3SSatish Balay       sum6 += v[5]*x1 + v[12]*x2 + v[19]*x3  + v[26]*x4 + v[33]*x5 + v[40]*x6 + v[47]*x7;
6322d61bbb3SSatish Balay       sum7 += v[6]*x1 + v[13]*x2 + v[20]*x3  + v[27]*x4 + v[34]*x5 + v[41]*x6 + v[48]*x7;
6332d61bbb3SSatish Balay       v    += 49;
6342d61bbb3SSatish Balay     }
6357b2bb3b9SHong Zhang     if (usecprow) z = zarray + 7*ridx[i];
6362d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
63726e093fcSHong Zhang     if (!usecprow) z += 7;
6382d61bbb3SSatish Balay   }
6392d61bbb3SSatish Balay 
6403649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
64126e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
6427c565772SBarry Smith   ierr = PetscLogFlops(98.0*a->nz - 7.0*a->nonzerorowcnt);CHKERRQ(ierr);
6432d61bbb3SSatish Balay   PetscFunctionReturn(0);
6442d61bbb3SSatish Balay }
6452d61bbb3SSatish Balay 
6468ab949d8SShri Abhyankar /* MatMult_SeqBAIJ_15 version 1: Columns in the block are accessed one at a time */
647832cc040SShri Abhyankar /* Default MatMult for block size 15 */
6488ab949d8SShri Abhyankar 
6498ab949d8SShri Abhyankar #undef __FUNCT__
6508ab949d8SShri Abhyankar #define __FUNCT__ "MatMult_SeqBAIJ_15_ver1"
6518ab949d8SShri Abhyankar PetscErrorCode MatMult_SeqBAIJ_15_ver1(Mat A,Vec xx,Vec zz)
6528ab949d8SShri Abhyankar {
6538ab949d8SShri Abhyankar   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
6548ab949d8SShri Abhyankar   PetscScalar       *z = 0,sum1,sum2,sum3,sum4,sum5,sum6,sum7,sum8,sum9,sum10,sum11,sum12,sum13,sum14,sum15;
6558ab949d8SShri Abhyankar   const PetscScalar *x,*xb;
65653ef36baSBarry Smith   PetscScalar       *zarray,xv;
6578ab949d8SShri Abhyankar   const MatScalar   *v;
6588ab949d8SShri Abhyankar   PetscErrorCode    ierr;
6598ab949d8SShri Abhyankar   const PetscInt    *ii,*ij=a->j,*idx;
6607c565772SBarry Smith   PetscInt          mbs,i,j,k,n,*ridx=NULL;
661ace3abfcSBarry Smith   PetscBool         usecprow=a->compressedrow.use;
6628ab949d8SShri Abhyankar 
6638ab949d8SShri Abhyankar   PetscFunctionBegin;
6643649974fSBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
6658ab949d8SShri Abhyankar   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
6668ab949d8SShri Abhyankar 
6678ab949d8SShri Abhyankar   v = a->a;
6688ab949d8SShri Abhyankar   if (usecprow) {
6698ab949d8SShri Abhyankar     mbs  = a->compressedrow.nrows;
6708ab949d8SShri Abhyankar     ii   = a->compressedrow.i;
6718ab949d8SShri Abhyankar     ridx = a->compressedrow.rindex;
672*ef16e671SStefano Zampini     ierr = PetscMemzero(zarray,15*a->mbs*sizeof(PetscScalar));CHKERRQ(ierr);
6738ab949d8SShri Abhyankar   } else {
6748ab949d8SShri Abhyankar     mbs = a->mbs;
6758ab949d8SShri Abhyankar     ii  = a->i;
6768ab949d8SShri Abhyankar     z   = zarray;
6778ab949d8SShri Abhyankar   }
6788ab949d8SShri Abhyankar 
6798ab949d8SShri Abhyankar   for (i=0; i<mbs; i++) {
6808ab949d8SShri Abhyankar     n    = ii[i+1] - ii[i];
6818ab949d8SShri Abhyankar     idx  = ij + ii[i];
6828ab949d8SShri Abhyankar     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0; sum7 = 0.0;
6838ab949d8SShri 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;
6848ab949d8SShri Abhyankar 
6858ab949d8SShri Abhyankar     for (j=0; j<n; j++) {
6868ab949d8SShri Abhyankar       xb = x + 15*(idx[j]);
6878ab949d8SShri Abhyankar 
6888ab949d8SShri Abhyankar       for (k=0; k<15; k++) {
68953ef36baSBarry Smith         xv     =  xb[k];
69053ef36baSBarry Smith         sum1  += v[0]*xv;
69153ef36baSBarry Smith         sum2  += v[1]*xv;
69253ef36baSBarry Smith         sum3  += v[2]*xv;
69353ef36baSBarry Smith         sum4  += v[3]*xv;
69453ef36baSBarry Smith         sum5  += v[4]*xv;
69553ef36baSBarry Smith         sum6  += v[5]*xv;
69653ef36baSBarry Smith         sum7  += v[6]*xv;
69753ef36baSBarry Smith         sum8  += v[7]*xv;
69853ef36baSBarry Smith         sum9  += v[8]*xv;
69953ef36baSBarry Smith         sum10 += v[9]*xv;
70053ef36baSBarry Smith         sum11 += v[10]*xv;
70153ef36baSBarry Smith         sum12 += v[11]*xv;
70253ef36baSBarry Smith         sum13 += v[12]*xv;
70353ef36baSBarry Smith         sum14 += v[13]*xv;
70453ef36baSBarry Smith         sum15 += v[14]*xv;
7058ab949d8SShri Abhyankar         v     += 15;
7068ab949d8SShri Abhyankar       }
7078ab949d8SShri Abhyankar     }
7088ab949d8SShri Abhyankar     if (usecprow) z = zarray + 15*ridx[i];
7098ab949d8SShri Abhyankar     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
7108ab949d8SShri 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;
7118ab949d8SShri Abhyankar 
7128ab949d8SShri Abhyankar     if (!usecprow) z += 15;
7138ab949d8SShri Abhyankar   }
7148ab949d8SShri Abhyankar 
7153649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
7168ab949d8SShri Abhyankar   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
7177c565772SBarry Smith   ierr = PetscLogFlops(450.0*a->nz - 15.0*a->nonzerorowcnt);CHKERRQ(ierr);
7188ab949d8SShri Abhyankar   PetscFunctionReturn(0);
7198ab949d8SShri Abhyankar }
7208ab949d8SShri Abhyankar 
7218ab949d8SShri Abhyankar /* MatMult_SeqBAIJ_15_ver2 : Columns in the block are accessed in sets of 4,4,4,3 */
7228ab949d8SShri Abhyankar #undef __FUNCT__
7238ab949d8SShri Abhyankar #define __FUNCT__ "MatMult_SeqBAIJ_15_ver2"
7248ab949d8SShri Abhyankar PetscErrorCode MatMult_SeqBAIJ_15_ver2(Mat A,Vec xx,Vec zz)
7258ab949d8SShri Abhyankar {
7268ab949d8SShri Abhyankar   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
7278ab949d8SShri Abhyankar   PetscScalar       *z = 0,sum1,sum2,sum3,sum4,sum5,sum6,sum7,sum8,sum9,sum10,sum11,sum12,sum13,sum14,sum15;
7288ab949d8SShri Abhyankar   const PetscScalar *x,*xb;
7290b8f6341SShri Abhyankar   PetscScalar       x1,x2,x3,x4,*zarray;
7308ab949d8SShri Abhyankar   const MatScalar   *v;
7318ab949d8SShri Abhyankar   PetscErrorCode    ierr;
7328ab949d8SShri Abhyankar   const PetscInt    *ii,*ij=a->j,*idx;
7337c565772SBarry Smith   PetscInt          mbs,i,j,n,*ridx=NULL;
734ace3abfcSBarry Smith   PetscBool         usecprow=a->compressedrow.use;
7358ab949d8SShri Abhyankar 
7368ab949d8SShri Abhyankar   PetscFunctionBegin;
7373649974fSBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
7388ab949d8SShri Abhyankar   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
7398ab949d8SShri Abhyankar 
7408ab949d8SShri Abhyankar   v = a->a;
7418ab949d8SShri Abhyankar   if (usecprow) {
7428ab949d8SShri Abhyankar     mbs  = a->compressedrow.nrows;
7438ab949d8SShri Abhyankar     ii   = a->compressedrow.i;
7448ab949d8SShri Abhyankar     ridx = a->compressedrow.rindex;
745*ef16e671SStefano Zampini     ierr = PetscMemzero(zarray,15*a->mbs*sizeof(PetscScalar));CHKERRQ(ierr);
7468ab949d8SShri Abhyankar   } else {
7478ab949d8SShri Abhyankar     mbs = a->mbs;
7488ab949d8SShri Abhyankar     ii  = a->i;
7498ab949d8SShri Abhyankar     z   = zarray;
7508ab949d8SShri Abhyankar   }
7518ab949d8SShri Abhyankar 
7528ab949d8SShri Abhyankar   for (i=0; i<mbs; i++) {
7538ab949d8SShri Abhyankar     n    = ii[i+1] - ii[i];
7548ab949d8SShri Abhyankar     idx  = ij + ii[i];
7558ab949d8SShri Abhyankar     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0; sum7 = 0.0;
7568ab949d8SShri 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;
7578ab949d8SShri Abhyankar 
7588ab949d8SShri Abhyankar     for (j=0; j<n; j++) {
7598ab949d8SShri Abhyankar       xb = x + 15*(idx[j]);
7600b8f6341SShri Abhyankar       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
7618ab949d8SShri Abhyankar 
7628ab949d8SShri Abhyankar       sum1  += v[0]*x1 + v[15]*x2 + v[30]*x3   + v[45]*x4;
7638ab949d8SShri Abhyankar       sum2  += v[1]*x1 + v[16]*x2 + v[31]*x3   + v[46]*x4;
7648ab949d8SShri Abhyankar       sum3  += v[2]*x1 + v[17]*x2 + v[32]*x3  + v[47]*x4;
7658ab949d8SShri Abhyankar       sum4  += v[3]*x1 + v[18]*x2 + v[33]*x3  + v[48]*x4;
7668ab949d8SShri Abhyankar       sum5  += v[4]*x1 + v[19]*x2 + v[34]*x3   + v[49]*x4;
7678ab949d8SShri Abhyankar       sum6  += v[5]*x1 + v[20]*x2 + v[35]*x3   + v[50]*x4;
7688ab949d8SShri Abhyankar       sum7  += v[6]*x1 + v[21]*x2 + v[36]*x3  + v[51]*x4;
7698ab949d8SShri Abhyankar       sum8  += v[7]*x1 + v[22]*x2 + v[37]*x3  + v[52]*x4;
7708ab949d8SShri Abhyankar       sum9  += v[8]*x1 + v[23]*x2 + v[38]*x3   + v[53]*x4;
7718ab949d8SShri Abhyankar       sum10 += v[9]*x1 + v[24]*x2 + v[39]*x3   + v[54]*x4;
7728ab949d8SShri Abhyankar       sum11 += v[10]*x1 + v[25]*x2 + v[40]*x3  + v[55]*x4;
7738ab949d8SShri Abhyankar       sum12 += v[11]*x1 + v[26]*x2 + v[41]*x3  + v[56]*x4;
7748ab949d8SShri Abhyankar       sum13 += v[12]*x1 + v[27]*x2 + v[42]*x3   + v[57]*x4;
7758ab949d8SShri Abhyankar       sum14 += v[13]*x1 + v[28]*x2 + v[43]*x3   + v[58]*x4;
7768ab949d8SShri Abhyankar       sum15 += v[14]*x1 + v[29]*x2 + v[44]*x3  + v[59]*x4;
7778ab949d8SShri Abhyankar 
7788ab949d8SShri Abhyankar       v += 60;
7798ab949d8SShri Abhyankar 
7800b8f6341SShri Abhyankar       x1 = xb[4]; x2 = xb[5]; x3 = xb[6]; x4 = xb[7];
7818ab949d8SShri Abhyankar 
7828ab949d8SShri Abhyankar       sum1  += v[0]*x1 + v[15]*x2 + v[30]*x3   + v[45]*x4;
7838ab949d8SShri Abhyankar       sum2  += v[1]*x1 + v[16]*x2 + v[31]*x3   + v[46]*x4;
7848ab949d8SShri Abhyankar       sum3  += v[2]*x1 + v[17]*x2 + v[32]*x3  + v[47]*x4;
7858ab949d8SShri Abhyankar       sum4  += v[3]*x1 + v[18]*x2 + v[33]*x3  + v[48]*x4;
7868ab949d8SShri Abhyankar       sum5  += v[4]*x1 + v[19]*x2 + v[34]*x3   + v[49]*x4;
7878ab949d8SShri Abhyankar       sum6  += v[5]*x1 + v[20]*x2 + v[35]*x3   + v[50]*x4;
7888ab949d8SShri Abhyankar       sum7  += v[6]*x1 + v[21]*x2 + v[36]*x3  + v[51]*x4;
7898ab949d8SShri Abhyankar       sum8  += v[7]*x1 + v[22]*x2 + v[37]*x3  + v[52]*x4;
7908ab949d8SShri Abhyankar       sum9  += v[8]*x1 + v[23]*x2 + v[38]*x3   + v[53]*x4;
7918ab949d8SShri Abhyankar       sum10 += v[9]*x1 + v[24]*x2 + v[39]*x3   + v[54]*x4;
7928ab949d8SShri Abhyankar       sum11 += v[10]*x1 + v[25]*x2 + v[40]*x3  + v[55]*x4;
7938ab949d8SShri Abhyankar       sum12 += v[11]*x1 + v[26]*x2 + v[41]*x3  + v[56]*x4;
7948ab949d8SShri Abhyankar       sum13 += v[12]*x1 + v[27]*x2 + v[42]*x3   + v[57]*x4;
7958ab949d8SShri Abhyankar       sum14 += v[13]*x1 + v[28]*x2 + v[43]*x3   + v[58]*x4;
7968ab949d8SShri Abhyankar       sum15 += v[14]*x1 + v[29]*x2 + v[44]*x3  + v[59]*x4;
7978ab949d8SShri Abhyankar       v     += 60;
7988ab949d8SShri Abhyankar 
7990b8f6341SShri Abhyankar       x1     = xb[8]; x2 = xb[9]; x3 = xb[10]; x4 = xb[11];
8000b8f6341SShri Abhyankar       sum1  += v[0]*x1 + v[15]*x2 + v[30]*x3   + v[45]*x4;
8010b8f6341SShri Abhyankar       sum2  += v[1]*x1 + v[16]*x2 + v[31]*x3   + v[46]*x4;
8020b8f6341SShri Abhyankar       sum3  += v[2]*x1 + v[17]*x2 + v[32]*x3  + v[47]*x4;
8030b8f6341SShri Abhyankar       sum4  += v[3]*x1 + v[18]*x2 + v[33]*x3  + v[48]*x4;
8040b8f6341SShri Abhyankar       sum5  += v[4]*x1 + v[19]*x2 + v[34]*x3   + v[49]*x4;
8050b8f6341SShri Abhyankar       sum6  += v[5]*x1 + v[20]*x2 + v[35]*x3   + v[50]*x4;
8060b8f6341SShri Abhyankar       sum7  += v[6]*x1 + v[21]*x2 + v[36]*x3  + v[51]*x4;
8070b8f6341SShri Abhyankar       sum8  += v[7]*x1 + v[22]*x2 + v[37]*x3  + v[52]*x4;
8080b8f6341SShri Abhyankar       sum9  += v[8]*x1 + v[23]*x2 + v[38]*x3   + v[53]*x4;
8090b8f6341SShri Abhyankar       sum10 += v[9]*x1 + v[24]*x2 + v[39]*x3   + v[54]*x4;
8100b8f6341SShri Abhyankar       sum11 += v[10]*x1 + v[25]*x2 + v[40]*x3  + v[55]*x4;
8110b8f6341SShri Abhyankar       sum12 += v[11]*x1 + v[26]*x2 + v[41]*x3  + v[56]*x4;
8120b8f6341SShri Abhyankar       sum13 += v[12]*x1 + v[27]*x2 + v[42]*x3   + v[57]*x4;
8130b8f6341SShri Abhyankar       sum14 += v[13]*x1 + v[28]*x2 + v[43]*x3   + v[58]*x4;
8140b8f6341SShri Abhyankar       sum15 += v[14]*x1 + v[29]*x2 + v[44]*x3  + v[59]*x4;
8150b8f6341SShri Abhyankar       v     += 60;
8160b8f6341SShri Abhyankar 
8170b8f6341SShri Abhyankar       x1     = xb[12]; x2 = xb[13]; x3 = xb[14];
8188ab949d8SShri Abhyankar       sum1  += v[0]*x1 + v[15]*x2 + v[30]*x3;
8198ab949d8SShri Abhyankar       sum2  += v[1]*x1 + v[16]*x2 + v[31]*x3;
8208ab949d8SShri Abhyankar       sum3  += v[2]*x1 + v[17]*x2 + v[32]*x3;
8218ab949d8SShri Abhyankar       sum4  += v[3]*x1 + v[18]*x2 + v[33]*x3;
8228ab949d8SShri Abhyankar       sum5  += v[4]*x1 + v[19]*x2 + v[34]*x3;
8238ab949d8SShri Abhyankar       sum6  += v[5]*x1 + v[20]*x2 + v[35]*x3;
8248ab949d8SShri Abhyankar       sum7  += v[6]*x1 + v[21]*x2 + v[36]*x3;
8258ab949d8SShri Abhyankar       sum8  += v[7]*x1 + v[22]*x2 + v[37]*x3;
8268ab949d8SShri Abhyankar       sum9  += v[8]*x1 + v[23]*x2 + v[38]*x3;
8278ab949d8SShri Abhyankar       sum10 += v[9]*x1 + v[24]*x2 + v[39]*x3;
8288ab949d8SShri Abhyankar       sum11 += v[10]*x1 + v[25]*x2 + v[40]*x3;
8298ab949d8SShri Abhyankar       sum12 += v[11]*x1 + v[26]*x2 + v[41]*x3;
8308ab949d8SShri Abhyankar       sum13 += v[12]*x1 + v[27]*x2 + v[42]*x3;
8318ab949d8SShri Abhyankar       sum14 += v[13]*x1 + v[28]*x2 + v[43]*x3;
8328ab949d8SShri Abhyankar       sum15 += v[14]*x1 + v[29]*x2 + v[44]*x3;
8338ab949d8SShri Abhyankar       v     += 45;
8348ab949d8SShri Abhyankar     }
8358ab949d8SShri Abhyankar     if (usecprow) z = zarray + 15*ridx[i];
8368ab949d8SShri Abhyankar     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
8378ab949d8SShri 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;
8388ab949d8SShri Abhyankar 
8398ab949d8SShri Abhyankar     if (!usecprow) z += 15;
8408ab949d8SShri Abhyankar   }
8418ab949d8SShri Abhyankar 
8423649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
8438ab949d8SShri Abhyankar   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
8447c565772SBarry Smith   ierr = PetscLogFlops(450.0*a->nz - 15.0*a->nonzerorowcnt);CHKERRQ(ierr);
8458ab949d8SShri Abhyankar   PetscFunctionReturn(0);
8468ab949d8SShri Abhyankar }
8478ab949d8SShri Abhyankar 
8488ab949d8SShri Abhyankar /* MatMult_SeqBAIJ_15_ver3 : Columns in the block are accessed in sets of 8,7 */
8498ab949d8SShri Abhyankar #undef __FUNCT__
8508ab949d8SShri Abhyankar #define __FUNCT__ "MatMult_SeqBAIJ_15_ver3"
8518ab949d8SShri Abhyankar PetscErrorCode MatMult_SeqBAIJ_15_ver3(Mat A,Vec xx,Vec zz)
8528ab949d8SShri Abhyankar {
8538ab949d8SShri Abhyankar   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
8548ab949d8SShri Abhyankar   PetscScalar       *z = 0,sum1,sum2,sum3,sum4,sum5,sum6,sum7,sum8,sum9,sum10,sum11,sum12,sum13,sum14,sum15;
8558ab949d8SShri Abhyankar   const PetscScalar *x,*xb;
8560b8f6341SShri Abhyankar   PetscScalar       x1,x2,x3,x4,x5,x6,x7,x8,*zarray;
8578ab949d8SShri Abhyankar   const MatScalar   *v;
8588ab949d8SShri Abhyankar   PetscErrorCode    ierr;
8598ab949d8SShri Abhyankar   const PetscInt    *ii,*ij=a->j,*idx;
8607c565772SBarry Smith   PetscInt          mbs,i,j,n,*ridx=NULL;
861ace3abfcSBarry Smith   PetscBool         usecprow=a->compressedrow.use;
8628ab949d8SShri Abhyankar 
8638ab949d8SShri Abhyankar   PetscFunctionBegin;
8643649974fSBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
8658ab949d8SShri Abhyankar   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
8668ab949d8SShri Abhyankar 
8678ab949d8SShri Abhyankar   v = a->a;
8688ab949d8SShri Abhyankar   if (usecprow) {
8698ab949d8SShri Abhyankar     mbs  = a->compressedrow.nrows;
8708ab949d8SShri Abhyankar     ii   = a->compressedrow.i;
8718ab949d8SShri Abhyankar     ridx = a->compressedrow.rindex;
872*ef16e671SStefano Zampini     ierr = PetscMemzero(zarray,15*a->mbs*sizeof(PetscScalar));CHKERRQ(ierr);
8738ab949d8SShri Abhyankar   } else {
8748ab949d8SShri Abhyankar     mbs = a->mbs;
8758ab949d8SShri Abhyankar     ii  = a->i;
8768ab949d8SShri Abhyankar     z   = zarray;
8778ab949d8SShri Abhyankar   }
8788ab949d8SShri Abhyankar 
8798ab949d8SShri Abhyankar   for (i=0; i<mbs; i++) {
8808ab949d8SShri Abhyankar     n    = ii[i+1] - ii[i];
8818ab949d8SShri Abhyankar     idx  = ij + ii[i];
8828ab949d8SShri Abhyankar     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0; sum7 = 0.0;
8838ab949d8SShri 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;
8848ab949d8SShri Abhyankar 
8858ab949d8SShri Abhyankar     for (j=0; j<n; j++) {
8868ab949d8SShri Abhyankar       xb = x + 15*(idx[j]);
8878ab949d8SShri Abhyankar       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6];
8880b8f6341SShri Abhyankar       x8 = xb[7];
8898ab949d8SShri Abhyankar 
8908ab949d8SShri 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;
8918ab949d8SShri 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;
8928ab949d8SShri 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;
8938ab949d8SShri 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;
8948ab949d8SShri 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;
8958ab949d8SShri 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;
8968ab949d8SShri 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;
8978ab949d8SShri 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;
8988ab949d8SShri 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;
8998ab949d8SShri 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;
9008ab949d8SShri 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;
9018ab949d8SShri 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;
9028ab949d8SShri 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;
9038ab949d8SShri 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;
9048ab949d8SShri 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;
9058ab949d8SShri Abhyankar       v     += 120;
9068ab949d8SShri Abhyankar 
9070b8f6341SShri Abhyankar       x1 = xb[8]; x2 = xb[9]; x3 = xb[10]; x4 = xb[11]; x5 = xb[12]; x6 = xb[13]; x7 = xb[14];
9080b8f6341SShri Abhyankar 
9098ab949d8SShri Abhyankar       sum1  += v[0]*x1 + v[15]*x2  + v[30]*x3  + v[45]*x4 + v[60]*x5 + v[75]*x6 + v[90]*x7;
9108ab949d8SShri Abhyankar       sum2  += v[1]*x1 + v[16]*x2  + v[31]*x3  + v[46]*x4 + v[61]*x5 + v[76]*x6 + v[91]*x7;
9118ab949d8SShri Abhyankar       sum3  += v[2]*x1 + v[17]*x2  + v[32]*x3  + v[47]*x4 + v[62]*x5 + v[77]*x6 + v[92]*x7;
9128ab949d8SShri Abhyankar       sum4  += v[3]*x1 + v[18]*x2 + v[33]*x3  + v[48]*x4 + v[63]*x5 + v[78]*x6 + v[93]*x7;
9138ab949d8SShri Abhyankar       sum5  += v[4]*x1 + v[19]*x2 + v[34]*x3  + v[49]*x4 + v[64]*x5 + v[79]*x6 + v[94]*x7;
9148ab949d8SShri Abhyankar       sum6  += v[5]*x1 + v[20]*x2 + v[35]*x3  + v[50]*x4 + v[65]*x5 + v[80]*x6 + v[95]*x7;
9158ab949d8SShri Abhyankar       sum7  += v[6]*x1 + v[21]*x2 + v[36]*x3  + v[51]*x4 + v[66]*x5 + v[81]*x6 + v[96]*x7;
9168ab949d8SShri Abhyankar       sum8  += v[7]*x1 + v[22]*x2  + v[37]*x3  + v[52]*x4 + v[67]*x5 + v[82]*x6 + v[97]*x7;
9178ab949d8SShri Abhyankar       sum9  += v[8]*x1 + v[23]*x2  + v[38]*x3  + v[53]*x4 + v[68]*x5 + v[83]*x6 + v[98]*x7;
9188ab949d8SShri Abhyankar       sum10 += v[9]*x1 + v[24]*x2  + v[39]*x3  + v[54]*x4 + v[69]*x5 + v[84]*x6 + v[99]*x7;
9198ab949d8SShri Abhyankar       sum11 += v[10]*x1 + v[25]*x2 + v[40]*x3  + v[55]*x4 + v[70]*x5 + v[85]*x6 + v[100]*x7;
9208ab949d8SShri Abhyankar       sum12 += v[11]*x1 + v[26]*x2 + v[41]*x3  + v[56]*x4 + v[71]*x5 + v[86]*x6 + v[101]*x7;
9218ab949d8SShri Abhyankar       sum13 += v[12]*x1 + v[27]*x2 + v[42]*x3  + v[57]*x4 + v[72]*x5 + v[87]*x6 + v[102]*x7;
9228ab949d8SShri Abhyankar       sum14 += v[13]*x1 + v[28]*x2 + v[43]*x3  + v[58]*x4 + v[73]*x5 + v[88]*x6 + v[103]*x7;
9238ab949d8SShri Abhyankar       sum15 += v[14]*x1 + v[29]*x2 + v[44]*x3  + v[59]*x4 + v[74]*x5 + v[89]*x6 + v[104]*x7;
9248ab949d8SShri Abhyankar       v     += 105;
9258ab949d8SShri Abhyankar     }
9268ab949d8SShri Abhyankar     if (usecprow) z = zarray + 15*ridx[i];
9278ab949d8SShri Abhyankar     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
9288ab949d8SShri 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;
9298ab949d8SShri Abhyankar 
9308ab949d8SShri Abhyankar     if (!usecprow) z += 15;
9318ab949d8SShri Abhyankar   }
9328ab949d8SShri Abhyankar 
9333649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
9348ab949d8SShri Abhyankar   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
9357c565772SBarry Smith   ierr = PetscLogFlops(450.0*a->nz - 15.0*a->nonzerorowcnt);CHKERRQ(ierr);
9368ab949d8SShri Abhyankar   PetscFunctionReturn(0);
9378ab949d8SShri Abhyankar }
9388ab949d8SShri Abhyankar 
9398ab949d8SShri Abhyankar /* MatMult_SeqBAIJ_15_ver4 : All columns in the block are accessed at once */
9408ab949d8SShri Abhyankar 
9418ab949d8SShri Abhyankar #undef __FUNCT__
9428ab949d8SShri Abhyankar #define __FUNCT__ "MatMult_SeqBAIJ_15_ver4"
9438ab949d8SShri Abhyankar PetscErrorCode MatMult_SeqBAIJ_15_ver4(Mat A,Vec xx,Vec zz)
9448ab949d8SShri Abhyankar {
9458ab949d8SShri Abhyankar   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
9468ab949d8SShri Abhyankar   PetscScalar       *z = 0,sum1,sum2,sum3,sum4,sum5,sum6,sum7,sum8,sum9,sum10,sum11,sum12,sum13,sum14,sum15;
9478ab949d8SShri Abhyankar   const PetscScalar *x,*xb;
9488ab949d8SShri Abhyankar   PetscScalar       x1,x2,x3,x4,x5,x6,x7,x8,x9,x10,x11,x12,x13,x14,x15,*zarray;
9498ab949d8SShri Abhyankar   const MatScalar   *v;
9508ab949d8SShri Abhyankar   PetscErrorCode    ierr;
9518ab949d8SShri Abhyankar   const PetscInt    *ii,*ij=a->j,*idx;
9527c565772SBarry Smith   PetscInt          mbs,i,j,n,*ridx=NULL;
953ace3abfcSBarry Smith   PetscBool         usecprow=a->compressedrow.use;
9548ab949d8SShri Abhyankar 
9558ab949d8SShri Abhyankar   PetscFunctionBegin;
9563649974fSBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
9578ab949d8SShri Abhyankar   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
9588ab949d8SShri Abhyankar 
9598ab949d8SShri Abhyankar   v = a->a;
9608ab949d8SShri Abhyankar   if (usecprow) {
9618ab949d8SShri Abhyankar     mbs  = a->compressedrow.nrows;
9628ab949d8SShri Abhyankar     ii   = a->compressedrow.i;
9638ab949d8SShri Abhyankar     ridx = a->compressedrow.rindex;
964*ef16e671SStefano Zampini     ierr = PetscMemzero(zarray,15*a->mbs*sizeof(PetscScalar));CHKERRQ(ierr);
9658ab949d8SShri Abhyankar   } else {
9668ab949d8SShri Abhyankar     mbs = a->mbs;
9678ab949d8SShri Abhyankar     ii  = a->i;
9688ab949d8SShri Abhyankar     z   = zarray;
9698ab949d8SShri Abhyankar   }
9708ab949d8SShri Abhyankar 
9718ab949d8SShri Abhyankar   for (i=0; i<mbs; i++) {
9728ab949d8SShri Abhyankar     n    = ii[i+1] - ii[i];
9738ab949d8SShri Abhyankar     idx  = ij + ii[i];
9748ab949d8SShri Abhyankar     sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0; sum7 = 0.0;
9758ab949d8SShri 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;
9768ab949d8SShri Abhyankar 
9778ab949d8SShri Abhyankar     for (j=0; j<n; j++) {
9788ab949d8SShri Abhyankar       xb = x + 15*(idx[j]);
9798ab949d8SShri Abhyankar       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6];
9808ab949d8SShri 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];
9818ab949d8SShri Abhyankar 
9828ab949d8SShri 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;
9838ab949d8SShri 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;
9848ab949d8SShri 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;
9858ab949d8SShri 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;
9868ab949d8SShri 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;
9878ab949d8SShri 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;
9888ab949d8SShri 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;
9898ab949d8SShri 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;
9908ab949d8SShri 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;
9918ab949d8SShri 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;
9928ab949d8SShri 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;
9938ab949d8SShri 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;
9948ab949d8SShri 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;
9958ab949d8SShri 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;
9968ab949d8SShri 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;
9978ab949d8SShri Abhyankar       v     += 225;
9988ab949d8SShri Abhyankar     }
9998ab949d8SShri Abhyankar     if (usecprow) z = zarray + 15*ridx[i];
10008ab949d8SShri Abhyankar     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
10018ab949d8SShri 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;
10028ab949d8SShri Abhyankar 
10038ab949d8SShri Abhyankar     if (!usecprow) z += 15;
10048ab949d8SShri Abhyankar   }
10058ab949d8SShri Abhyankar 
10063649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
10078ab949d8SShri Abhyankar   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
10087c565772SBarry Smith   ierr = PetscLogFlops(450.0*a->nz - 15.0*a->nonzerorowcnt);CHKERRQ(ierr);
10098ab949d8SShri Abhyankar   PetscFunctionReturn(0);
10108ab949d8SShri Abhyankar }
10118ab949d8SShri Abhyankar 
10128ab949d8SShri Abhyankar 
10133f1db9ecSBarry Smith /*
10143f1db9ecSBarry Smith     This will not work with MatScalar == float because it calls the BLAS
10153f1db9ecSBarry Smith */
10164a2ae208SSatish Balay #undef __FUNCT__
10174a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_N"
1018dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_N(Mat A,Vec xx,Vec zz)
10192d61bbb3SSatish Balay {
10202d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
1021d9ca1df4SBarry Smith   PetscScalar       *z = 0,*work,*workt,*zarray;
1022d9ca1df4SBarry Smith   const PetscScalar *x,*xb;
1023d9ca1df4SBarry Smith   const MatScalar   *v;
1024dfbe8321SBarry Smith   PetscErrorCode    ierr;
1025d9ca1df4SBarry Smith   PetscInt          mbs,i,bs=A->rmap->bs,j,n,bs2=a->bs2;
1026d9ca1df4SBarry Smith   const PetscInt    *idx,*ii,*ridx=NULL;
1027d9ca1df4SBarry Smith   PetscInt          ncols,k;
1028ace3abfcSBarry Smith   PetscBool         usecprow=a->compressedrow.use;
10292d61bbb3SSatish Balay 
10302d61bbb3SSatish Balay   PetscFunctionBegin;
1031d9ca1df4SBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
103226e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
10332d61bbb3SSatish Balay 
10342d61bbb3SSatish Balay   idx = a->j;
10352d61bbb3SSatish Balay   v   = a->a;
103626e093fcSHong Zhang   if (usecprow) {
103726e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
103826e093fcSHong Zhang     ii   = a->compressedrow.i;
10397b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
1040*ef16e671SStefano Zampini     ierr = PetscMemzero(zarray,bs*a->mbs*sizeof(PetscScalar));CHKERRQ(ierr);
104126e093fcSHong Zhang   } else {
104226e093fcSHong Zhang     mbs = a->mbs;
10432d61bbb3SSatish Balay     ii  = a->i;
104426e093fcSHong Zhang     z   = zarray;
104526e093fcSHong Zhang   }
1046218c64b6SSatish Balay 
10472d61bbb3SSatish Balay   if (!a->mult_work) {
1048d0f46423SBarry Smith     k    = PetscMax(A->rmap->n,A->cmap->n);
1049854ce69bSBarry Smith     ierr = PetscMalloc1(k+1,&a->mult_work);CHKERRQ(ierr);
10502d61bbb3SSatish Balay   }
10512d61bbb3SSatish Balay   work = a->mult_work;
10522d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
10532d61bbb3SSatish Balay     n           = ii[1] - ii[0]; ii++;
10542d61bbb3SSatish Balay     ncols       = n*bs;
10552d61bbb3SSatish Balay     workt       = work;
10562d61bbb3SSatish Balay     for (j=0; j<n; j++) {
10572d61bbb3SSatish Balay       xb = x + bs*(*idx++);
10582d61bbb3SSatish Balay       for (k=0; k<bs; k++) workt[k] = xb[k];
10592d61bbb3SSatish Balay       workt += bs;
10602d61bbb3SSatish Balay     }
10617b2bb3b9SHong Zhang     if (usecprow) z = zarray + bs*ridx[i];
106296b95a6bSBarry Smith     PetscKernel_w_gets_Ar_times_v(bs,ncols,work,v,z);
106371044d3cSBarry Smith     /* BLASgemv_("N",&bs,&ncols,&_DOne,v,&bs,work,&_One,&_DZero,z,&_One); */
10642d61bbb3SSatish Balay     v += n*bs2;
106526e093fcSHong Zhang     if (!usecprow) z += bs;
10662d61bbb3SSatish Balay   }
1067d9ca1df4SBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
106826e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
10697c565772SBarry Smith   ierr = PetscLogFlops(2.0*a->nz*bs2 - bs*a->nonzerorowcnt);CHKERRQ(ierr);
10702d61bbb3SSatish Balay   PetscFunctionReturn(0);
10712d61bbb3SSatish Balay }
10722d61bbb3SSatish Balay 
10734a2ae208SSatish Balay #undef __FUNCT__
10744a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_1"
1075dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_1(Mat A,Vec xx,Vec yy,Vec zz)
10762d61bbb3SSatish Balay {
10772d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
1078122f12eaSBarry Smith   const PetscScalar *x;
1079122f12eaSBarry Smith   PetscScalar       *y,*z,sum;
1080122f12eaSBarry Smith   const MatScalar   *v;
1081dfbe8321SBarry Smith   PetscErrorCode    ierr;
10827c565772SBarry Smith   PetscInt          mbs=a->mbs,i,n,*ridx=NULL;
1083122f12eaSBarry Smith   const PetscInt    *idx,*ii;
1084ace3abfcSBarry Smith   PetscBool         usecprow=a->compressedrow.use;
10852d61bbb3SSatish Balay 
10862d61bbb3SSatish Balay   PetscFunctionBegin;
10873649974fSBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
1088d9ca1df4SBarry Smith   ierr = VecGetArrayPair(yy,zz,&y,&z);CHKERRQ(ierr);
10892d61bbb3SSatish Balay 
10902d61bbb3SSatish Balay   idx = a->j;
10912d61bbb3SSatish Balay   v   = a->a;
109226e093fcSHong Zhang   if (usecprow) {
10934eb6d288SHong Zhang     if (zz != yy) {
1094122f12eaSBarry Smith       ierr = PetscMemcpy(z,y,mbs*sizeof(PetscScalar));CHKERRQ(ierr);
10954eb6d288SHong Zhang     }
109626e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
109726e093fcSHong Zhang     ii   = a->compressedrow.i;
10987b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
109926e093fcSHong Zhang   } else {
11002d61bbb3SSatish Balay     ii = a->i;
110126e093fcSHong Zhang   }
11022d61bbb3SSatish Balay 
11032d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
1104122f12eaSBarry Smith     n = ii[1] - ii[0];
1105122f12eaSBarry Smith     ii++;
110626e093fcSHong Zhang     if (!usecprow) {
1107122f12eaSBarry Smith       sum         = y[i];
1108122f12eaSBarry Smith     } else {
1109122f12eaSBarry Smith       sum = y[ridx[i]];
1110122f12eaSBarry Smith     }
1111444d8c10SJed Brown     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA); /* Indices for the next row (assumes same size as this one) */
1112444d8c10SJed Brown     PetscPrefetchBlock(v+n,n,0,PETSC_PREFETCH_HINT_NTA);   /* Entries for the next row */
1113122f12eaSBarry Smith     PetscSparseDensePlusDot(sum,x,v,idx,n);
1114122f12eaSBarry Smith     v   += n;
1115122f12eaSBarry Smith     idx += n;
1116122f12eaSBarry Smith     if (usecprow) {
1117122f12eaSBarry Smith       z[ridx[i]] = sum;
1118122f12eaSBarry Smith     } else {
1119122f12eaSBarry Smith       z[i] = sum;
112026e093fcSHong Zhang     }
11212d61bbb3SSatish Balay   }
11223649974fSBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
1123d9ca1df4SBarry Smith   ierr = VecRestoreArrayPair(yy,zz,&y,&z);CHKERRQ(ierr);
11247c565772SBarry Smith   ierr = PetscLogFlops(2.0*a->nz);CHKERRQ(ierr);
11252d61bbb3SSatish Balay   PetscFunctionReturn(0);
11262d61bbb3SSatish Balay }
11272d61bbb3SSatish Balay 
11284a2ae208SSatish Balay #undef __FUNCT__
11294a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_2"
1130dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_2(Mat A,Vec xx,Vec yy,Vec zz)
11312d61bbb3SSatish Balay {
11322d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
1133d9ca1df4SBarry Smith   PetscScalar       *y = 0,*z = 0,sum1,sum2;
1134d9ca1df4SBarry Smith   const PetscScalar *x,*xb;
113526e093fcSHong Zhang   PetscScalar       x1,x2,*yarray,*zarray;
1136d9ca1df4SBarry Smith   const MatScalar   *v;
1137dfbe8321SBarry Smith   PetscErrorCode    ierr;
1138d9ca1df4SBarry Smith   PetscInt          mbs = a->mbs,i,n,j;
1139d9ca1df4SBarry Smith   const PetscInt    *idx,*ii,*ridx = NULL;
1140ace3abfcSBarry Smith   PetscBool         usecprow = a->compressedrow.use;
11412d61bbb3SSatish Balay 
11422d61bbb3SSatish Balay   PetscFunctionBegin;
1143d9ca1df4SBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
1144d9ca1df4SBarry Smith   ierr = VecGetArrayPair(yy,zz,&yarray,&zarray);CHKERRQ(ierr);
11452d61bbb3SSatish Balay 
11462d61bbb3SSatish Balay   idx = a->j;
11472d61bbb3SSatish Balay   v   = a->a;
114826e093fcSHong Zhang   if (usecprow) {
11494eb6d288SHong Zhang     if (zz != yy) {
11504eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,2*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
11514eb6d288SHong Zhang     }
115226e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
115326e093fcSHong Zhang     ii   = a->compressedrow.i;
11547b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
115526e093fcSHong Zhang   } else {
11562d61bbb3SSatish Balay     ii = a->i;
115726e093fcSHong Zhang     y  = yarray;
115826e093fcSHong Zhang     z  = zarray;
115926e093fcSHong Zhang   }
11602d61bbb3SSatish Balay 
11612d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
11622d61bbb3SSatish Balay     n = ii[1] - ii[0]; ii++;
116326e093fcSHong Zhang     if (usecprow) {
11647b2bb3b9SHong Zhang       z = zarray + 2*ridx[i];
11657b2bb3b9SHong Zhang       y = yarray + 2*ridx[i];
116626e093fcSHong Zhang     }
11672d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1];
1168444d8c10SJed Brown     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);   /* Indices for the next row (assumes same size as this one) */
1169444d8c10SJed Brown     PetscPrefetchBlock(v+4*n,4*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
11702d61bbb3SSatish Balay     for (j=0; j<n; j++) {
117126fbe8dcSKarl Rupp       xb = x + 2*(*idx++);
117226fbe8dcSKarl Rupp       x1 = xb[0];
117326fbe8dcSKarl Rupp       x2 = xb[1];
117426fbe8dcSKarl Rupp 
11752d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[2]*x2;
11762d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[3]*x2;
11772d61bbb3SSatish Balay       v    += 4;
11782d61bbb3SSatish Balay     }
11792d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2;
118026e093fcSHong Zhang     if (!usecprow) {
11812d61bbb3SSatish Balay       z += 2; y += 2;
11822d61bbb3SSatish Balay     }
118326e093fcSHong Zhang   }
1184d9ca1df4SBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
1185d9ca1df4SBarry Smith   ierr = VecRestoreArrayPair(yy,zz,&yarray,&zarray);CHKERRQ(ierr);
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;
1195d9ca1df4SBarry Smith   PetscScalar       *y = 0,*z = 0,sum1,sum2,sum3,x1,x2,x3,*yarray,*zarray;
1196d9ca1df4SBarry Smith   const PetscScalar *x,*xb;
1197d9ca1df4SBarry Smith   const MatScalar   *v;
1198dfbe8321SBarry Smith   PetscErrorCode    ierr;
1199d9ca1df4SBarry Smith   PetscInt          mbs = a->mbs,i,j,n;
1200d9ca1df4SBarry Smith   const PetscInt    *idx,*ii,*ridx = NULL;
1201ace3abfcSBarry Smith   PetscBool         usecprow = a->compressedrow.use;
12022d61bbb3SSatish Balay 
12032d61bbb3SSatish Balay   PetscFunctionBegin;
1204d9ca1df4SBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
1205d9ca1df4SBarry Smith   ierr = VecGetArrayPair(yy,zz,&yarray,&zarray);CHKERRQ(ierr);
12062d61bbb3SSatish Balay 
12072d61bbb3SSatish Balay   idx = a->j;
12082d61bbb3SSatish Balay   v   = a->a;
120926e093fcSHong Zhang   if (usecprow) {
12104eb6d288SHong Zhang     if (zz != yy) {
12114eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,3*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
12124eb6d288SHong Zhang     }
121326e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
121426e093fcSHong Zhang     ii   = a->compressedrow.i;
12157b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
121626e093fcSHong Zhang   } else {
12172d61bbb3SSatish Balay     ii = a->i;
121826e093fcSHong Zhang     y  = yarray;
121926e093fcSHong Zhang     z  = zarray;
122026e093fcSHong Zhang   }
12212d61bbb3SSatish Balay 
12222d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
12232d61bbb3SSatish Balay     n = ii[1] - ii[0]; ii++;
122426e093fcSHong Zhang     if (usecprow) {
12257b2bb3b9SHong Zhang       z = zarray + 3*ridx[i];
12267b2bb3b9SHong Zhang       y = yarray + 3*ridx[i];
122726e093fcSHong Zhang     }
12282d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2];
1229444d8c10SJed Brown     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);   /* Indices for the next row (assumes same size as this one) */
1230444d8c10SJed Brown     PetscPrefetchBlock(v+9*n,9*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
12312d61bbb3SSatish Balay     for (j=0; j<n; j++) {
12322d61bbb3SSatish Balay       xb    = x + 3*(*idx++); x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
12332d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[3]*x2 + v[6]*x3;
12342d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[4]*x2 + v[7]*x3;
12352d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[5]*x2 + v[8]*x3;
12362d61bbb3SSatish Balay       v    += 9;
12372d61bbb3SSatish Balay     }
12382d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3;
123926e093fcSHong Zhang     if (!usecprow) {
12402d61bbb3SSatish Balay       z += 3; y += 3;
12412d61bbb3SSatish Balay     }
124226e093fcSHong Zhang   }
1243d9ca1df4SBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
1244d9ca1df4SBarry Smith   ierr = VecRestoreArrayPair(yy,zz,&yarray,&zarray);CHKERRQ(ierr);
1245dc0b31edSSatish Balay   ierr = PetscLogFlops(18.0*a->nz);CHKERRQ(ierr);
12462d61bbb3SSatish Balay   PetscFunctionReturn(0);
12472d61bbb3SSatish Balay }
12482d61bbb3SSatish Balay 
12494a2ae208SSatish Balay #undef __FUNCT__
12504a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_4"
1251dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_4(Mat A,Vec xx,Vec yy,Vec zz)
12522d61bbb3SSatish Balay {
12532d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
1254d9ca1df4SBarry Smith   PetscScalar       *y = 0,*z = 0,sum1,sum2,sum3,sum4,x1,x2,x3,x4,*yarray,*zarray;
1255d9ca1df4SBarry Smith   const PetscScalar *x,*xb;
1256d9ca1df4SBarry Smith   const MatScalar   *v;
1257dfbe8321SBarry Smith   PetscErrorCode    ierr;
1258d9ca1df4SBarry Smith   PetscInt          mbs = a->mbs,i,j,n;
1259d9ca1df4SBarry Smith   const PetscInt    *idx,*ii,*ridx=NULL;
1260ace3abfcSBarry Smith   PetscBool         usecprow=a->compressedrow.use;
12612d61bbb3SSatish Balay 
12622d61bbb3SSatish Balay   PetscFunctionBegin;
1263d9ca1df4SBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
1264d9ca1df4SBarry Smith   ierr = VecGetArrayPair(yy,zz,&yarray,&zarray);CHKERRQ(ierr);
12652d61bbb3SSatish Balay 
12662d61bbb3SSatish Balay   idx = a->j;
12672d61bbb3SSatish Balay   v   = a->a;
126826e093fcSHong Zhang   if (usecprow) {
12694eb6d288SHong Zhang     if (zz != yy) {
12704eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,4*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
12714eb6d288SHong Zhang     }
127226e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
127326e093fcSHong Zhang     ii   = a->compressedrow.i;
12747b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
127526e093fcSHong Zhang   } else {
12762d61bbb3SSatish Balay     ii = a->i;
127726e093fcSHong Zhang     y  = yarray;
127826e093fcSHong Zhang     z  = zarray;
127926e093fcSHong Zhang   }
12802d61bbb3SSatish Balay 
12812d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
12822d61bbb3SSatish Balay     n = ii[1] - ii[0]; ii++;
128326e093fcSHong Zhang     if (usecprow) {
12847b2bb3b9SHong Zhang       z = zarray + 4*ridx[i];
12857b2bb3b9SHong Zhang       y = yarray + 4*ridx[i];
128626e093fcSHong Zhang     }
12872d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3];
1288444d8c10SJed Brown     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);     /* Indices for the next row (assumes same size as this one) */
1289444d8c10SJed Brown     PetscPrefetchBlock(v+16*n,16*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
12902d61bbb3SSatish Balay     for (j=0; j<n; j++) {
12912d61bbb3SSatish Balay       xb    = x + 4*(*idx++);
12922d61bbb3SSatish Balay       x1    = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
12932d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[4]*x2 + v[8]*x3   + v[12]*x4;
12942d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[5]*x2 + v[9]*x3   + v[13]*x4;
12952d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[6]*x2 + v[10]*x3  + v[14]*x4;
12962d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[7]*x2 + v[11]*x3  + v[15]*x4;
12972d61bbb3SSatish Balay       v    += 16;
12982d61bbb3SSatish Balay     }
12992d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4;
130026e093fcSHong Zhang     if (!usecprow) {
13012d61bbb3SSatish Balay       z += 4; y += 4;
13022d61bbb3SSatish Balay     }
130326e093fcSHong Zhang   }
1304d9ca1df4SBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
1305d9ca1df4SBarry Smith   ierr = VecRestoreArrayPair(yy,zz,&yarray,&zarray);CHKERRQ(ierr);
1306dc0b31edSSatish Balay   ierr = PetscLogFlops(32.0*a->nz);CHKERRQ(ierr);
13072d61bbb3SSatish Balay   PetscFunctionReturn(0);
13082d61bbb3SSatish Balay }
13092d61bbb3SSatish Balay 
13104a2ae208SSatish Balay #undef __FUNCT__
13114a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_5"
1312dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_5(Mat A,Vec xx,Vec yy,Vec zz)
13132d61bbb3SSatish Balay {
13142d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
1315d9ca1df4SBarry Smith   PetscScalar       *y = 0,*z = 0,sum1,sum2,sum3,sum4,sum5,x1,x2,x3,x4,x5;
1316d9ca1df4SBarry Smith   const PetscScalar *x,*xb;
131726e093fcSHong Zhang   PetscScalar       *yarray,*zarray;
1318d9ca1df4SBarry Smith   const MatScalar   *v;
1319dfbe8321SBarry Smith   PetscErrorCode    ierr;
1320d9ca1df4SBarry Smith   PetscInt          mbs = a->mbs,i,j,n;
1321d9ca1df4SBarry Smith   const PetscInt    *idx,*ii,*ridx = NULL;
1322ace3abfcSBarry Smith   PetscBool         usecprow=a->compressedrow.use;
13232d61bbb3SSatish Balay 
13242d61bbb3SSatish Balay   PetscFunctionBegin;
1325d9ca1df4SBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
1326d9ca1df4SBarry Smith   ierr = VecGetArrayPair(yy,zz,&yarray,&zarray);CHKERRQ(ierr);
13272d61bbb3SSatish Balay 
13282d61bbb3SSatish Balay   idx = a->j;
13292d61bbb3SSatish Balay   v   = a->a;
133026e093fcSHong Zhang   if (usecprow) {
13314eb6d288SHong Zhang     if (zz != yy) {
13324eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,5*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
13334eb6d288SHong Zhang     }
133426e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
133526e093fcSHong Zhang     ii   = a->compressedrow.i;
13367b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
133726e093fcSHong Zhang   } else {
13382d61bbb3SSatish Balay     ii = a->i;
133926e093fcSHong Zhang     y  = yarray;
134026e093fcSHong Zhang     z  = zarray;
134126e093fcSHong Zhang   }
13422d61bbb3SSatish Balay 
13432d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
13442d61bbb3SSatish Balay     n = ii[1] - ii[0]; ii++;
134526e093fcSHong Zhang     if (usecprow) {
13467b2bb3b9SHong Zhang       z = zarray + 5*ridx[i];
13477b2bb3b9SHong Zhang       y = yarray + 5*ridx[i];
134826e093fcSHong Zhang     }
13492d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4];
1350444d8c10SJed Brown     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);     /* Indices for the next row (assumes same size as this one) */
1351444d8c10SJed Brown     PetscPrefetchBlock(v+25*n,25*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
13522d61bbb3SSatish Balay     for (j=0; j<n; j++) {
13532d61bbb3SSatish Balay       xb    = x + 5*(*idx++);
13542d61bbb3SSatish Balay       x1    = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4];
13552d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[5]*x2 + v[10]*x3  + v[15]*x4 + v[20]*x5;
13562d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[6]*x2 + v[11]*x3  + v[16]*x4 + v[21]*x5;
13572d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[7]*x2 + v[12]*x3  + v[17]*x4 + v[22]*x5;
13582d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[8]*x2 + v[13]*x3  + v[18]*x4 + v[23]*x5;
13592d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[9]*x2 + v[14]*x3  + v[19]*x4 + v[24]*x5;
13602d61bbb3SSatish Balay       v    += 25;
13612d61bbb3SSatish Balay     }
13622d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5;
136326e093fcSHong Zhang     if (!usecprow) {
13642d61bbb3SSatish Balay       z += 5; y += 5;
13652d61bbb3SSatish Balay     }
136626e093fcSHong Zhang   }
1367d9ca1df4SBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
1368d9ca1df4SBarry Smith   ierr = VecRestoreArrayPair(yy,zz,&yarray,&zarray);CHKERRQ(ierr);
1369dc0b31edSSatish Balay   ierr = PetscLogFlops(50.0*a->nz);CHKERRQ(ierr);
13702d61bbb3SSatish Balay   PetscFunctionReturn(0);
13712d61bbb3SSatish Balay }
13724a2ae208SSatish Balay #undef __FUNCT__
13734a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_6"
1374dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_6(Mat A,Vec xx,Vec yy,Vec zz)
137515091d37SBarry Smith {
137615091d37SBarry Smith   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
1377d9ca1df4SBarry Smith   PetscScalar       *y = 0,*z = 0,sum1,sum2,sum3,sum4,sum5,sum6;
1378d9ca1df4SBarry Smith   const PetscScalar *x,*xb;
137926e093fcSHong Zhang   PetscScalar       x1,x2,x3,x4,x5,x6,*yarray,*zarray;
1380d9ca1df4SBarry Smith   const MatScalar   *v;
1381dfbe8321SBarry Smith   PetscErrorCode    ierr;
1382d9ca1df4SBarry Smith   PetscInt          mbs = a->mbs,i,j,n;
1383d9ca1df4SBarry Smith   const PetscInt    *idx,*ii,*ridx=NULL;
1384ace3abfcSBarry Smith   PetscBool         usecprow=a->compressedrow.use;
138515091d37SBarry Smith 
138615091d37SBarry Smith   PetscFunctionBegin;
1387d9ca1df4SBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
1388d9ca1df4SBarry Smith   ierr = VecGetArrayPair(yy,zz,&yarray,&zarray);CHKERRQ(ierr);
138915091d37SBarry Smith 
139015091d37SBarry Smith   idx = a->j;
139115091d37SBarry Smith   v   = a->a;
139226e093fcSHong Zhang   if (usecprow) {
13934eb6d288SHong Zhang     if (zz != yy) {
13944eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,6*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
13954eb6d288SHong Zhang     }
139626e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
139726e093fcSHong Zhang     ii   = a->compressedrow.i;
13987b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
139926e093fcSHong Zhang   } else {
140015091d37SBarry Smith     ii = a->i;
140126e093fcSHong Zhang     y  = yarray;
140226e093fcSHong Zhang     z  = zarray;
140326e093fcSHong Zhang   }
140415091d37SBarry Smith 
140515091d37SBarry Smith   for (i=0; i<mbs; i++) {
140615091d37SBarry Smith     n = ii[1] - ii[0]; ii++;
140726e093fcSHong Zhang     if (usecprow) {
14087b2bb3b9SHong Zhang       z = zarray + 6*ridx[i];
14097b2bb3b9SHong Zhang       y = yarray + 6*ridx[i];
141026e093fcSHong Zhang     }
141115091d37SBarry Smith     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4]; sum6 = y[5];
1412444d8c10SJed Brown     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);     /* Indices for the next row (assumes same size as this one) */
1413444d8c10SJed Brown     PetscPrefetchBlock(v+36*n,36*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
141415091d37SBarry Smith     for (j=0; j<n; j++) {
14153b95cb0eSSatish Balay       xb    = x + 6*(*idx++);
141615091d37SBarry Smith       x1    = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5];
141715091d37SBarry Smith       sum1 += v[0]*x1 + v[6]*x2  + v[12]*x3  + v[18]*x4 + v[24]*x5 + v[30]*x6;
141815091d37SBarry Smith       sum2 += v[1]*x1 + v[7]*x2  + v[13]*x3  + v[19]*x4 + v[25]*x5 + v[31]*x6;
141915091d37SBarry Smith       sum3 += v[2]*x1 + v[8]*x2  + v[14]*x3  + v[20]*x4 + v[26]*x5 + v[32]*x6;
142015091d37SBarry Smith       sum4 += v[3]*x1 + v[9]*x2  + v[15]*x3  + v[21]*x4 + v[27]*x5 + v[33]*x6;
142115091d37SBarry Smith       sum5 += v[4]*x1 + v[10]*x2 + v[16]*x3  + v[22]*x4 + v[28]*x5 + v[34]*x6;
142215091d37SBarry Smith       sum6 += v[5]*x1 + v[11]*x2 + v[17]*x3  + v[23]*x4 + v[29]*x5 + v[35]*x6;
142315091d37SBarry Smith       v    += 36;
142415091d37SBarry Smith     }
142515091d37SBarry Smith     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6;
142626e093fcSHong Zhang     if (!usecprow) {
142715091d37SBarry Smith       z += 6; y += 6;
142815091d37SBarry Smith     }
142926e093fcSHong Zhang   }
1430d9ca1df4SBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
1431d9ca1df4SBarry Smith   ierr = VecRestoreArrayPair(yy,zz,&yarray,&zarray);CHKERRQ(ierr);
1432dc0b31edSSatish Balay   ierr = PetscLogFlops(72.0*a->nz);CHKERRQ(ierr);
143315091d37SBarry Smith   PetscFunctionReturn(0);
143415091d37SBarry Smith }
14352d61bbb3SSatish Balay 
14364a2ae208SSatish Balay #undef __FUNCT__
14374a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_7"
1438dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_7(Mat A,Vec xx,Vec yy,Vec zz)
14392d61bbb3SSatish Balay {
14402d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
1441d9ca1df4SBarry Smith   PetscScalar       *y = 0,*z = 0,sum1,sum2,sum3,sum4,sum5,sum6,sum7;
1442d9ca1df4SBarry Smith   const PetscScalar *x,*xb;
144326e093fcSHong Zhang   PetscScalar       x1,x2,x3,x4,x5,x6,x7,*yarray,*zarray;
1444d9ca1df4SBarry Smith   const MatScalar   *v;
1445dfbe8321SBarry Smith   PetscErrorCode    ierr;
1446d9ca1df4SBarry Smith   PetscInt          mbs = a->mbs,i,j,n;
1447d9ca1df4SBarry Smith   const PetscInt    *idx,*ii,*ridx = NULL;
1448ace3abfcSBarry Smith   PetscBool         usecprow=a->compressedrow.use;
14492d61bbb3SSatish Balay 
14502d61bbb3SSatish Balay   PetscFunctionBegin;
1451d9ca1df4SBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
1452d9ca1df4SBarry Smith   ierr = VecGetArrayPair(yy,zz,&yarray,&zarray);CHKERRQ(ierr);
14532d61bbb3SSatish Balay 
14542d61bbb3SSatish Balay   idx = a->j;
14552d61bbb3SSatish Balay   v   = a->a;
145626e093fcSHong Zhang   if (usecprow) {
14574eb6d288SHong Zhang     if (zz != yy) {
14584eb6d288SHong Zhang       ierr = PetscMemcpy(zarray,yarray,7*mbs*sizeof(PetscScalar));CHKERRQ(ierr);
14594eb6d288SHong Zhang     }
146026e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
146126e093fcSHong Zhang     ii   = a->compressedrow.i;
14627b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
146326e093fcSHong Zhang   } else {
14642d61bbb3SSatish Balay     ii = a->i;
146526e093fcSHong Zhang     y  = yarray;
146626e093fcSHong Zhang     z  = zarray;
146726e093fcSHong Zhang   }
14682d61bbb3SSatish Balay 
14692d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
14702d61bbb3SSatish Balay     n = ii[1] - ii[0]; ii++;
147126e093fcSHong Zhang     if (usecprow) {
14727b2bb3b9SHong Zhang       z = zarray + 7*ridx[i];
14737b2bb3b9SHong Zhang       y = yarray + 7*ridx[i];
147426e093fcSHong Zhang     }
14752d61bbb3SSatish Balay     sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4]; sum6 = y[5]; sum7 = y[6];
1476444d8c10SJed Brown     PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA);     /* Indices for the next row (assumes same size as this one) */
1477444d8c10SJed Brown     PetscPrefetchBlock(v+49*n,49*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */
14782d61bbb3SSatish Balay     for (j=0; j<n; j++) {
14792d61bbb3SSatish Balay       xb    = x + 7*(*idx++);
14802d61bbb3SSatish Balay       x1    = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6];
14812d61bbb3SSatish Balay       sum1 += v[0]*x1 + v[7]*x2  + v[14]*x3  + v[21]*x4 + v[28]*x5 + v[35]*x6 + v[42]*x7;
14822d61bbb3SSatish Balay       sum2 += v[1]*x1 + v[8]*x2  + v[15]*x3  + v[22]*x4 + v[29]*x5 + v[36]*x6 + v[43]*x7;
14832d61bbb3SSatish Balay       sum3 += v[2]*x1 + v[9]*x2  + v[16]*x3  + v[23]*x4 + v[30]*x5 + v[37]*x6 + v[44]*x7;
14842d61bbb3SSatish Balay       sum4 += v[3]*x1 + v[10]*x2 + v[17]*x3  + v[24]*x4 + v[31]*x5 + v[38]*x6 + v[45]*x7;
14852d61bbb3SSatish Balay       sum5 += v[4]*x1 + v[11]*x2 + v[18]*x3  + v[25]*x4 + v[32]*x5 + v[39]*x6 + v[46]*x7;
14862d61bbb3SSatish Balay       sum6 += v[5]*x1 + v[12]*x2 + v[19]*x3  + v[26]*x4 + v[33]*x5 + v[40]*x6 + v[47]*x7;
14872d61bbb3SSatish Balay       sum7 += v[6]*x1 + v[13]*x2 + v[20]*x3  + v[27]*x4 + v[34]*x5 + v[41]*x6 + v[48]*x7;
14882d61bbb3SSatish Balay       v    += 49;
14892d61bbb3SSatish Balay     }
14902d61bbb3SSatish Balay     z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7;
149126e093fcSHong Zhang     if (!usecprow) {
14922d61bbb3SSatish Balay       z += 7; y += 7;
14932d61bbb3SSatish Balay     }
149426e093fcSHong Zhang   }
1495d9ca1df4SBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
1496d9ca1df4SBarry Smith   ierr = VecRestoreArrayPair(yy,zz,&yarray,&zarray);CHKERRQ(ierr);
1497dc0b31edSSatish Balay   ierr = PetscLogFlops(98.0*a->nz);CHKERRQ(ierr);
14982d61bbb3SSatish Balay   PetscFunctionReturn(0);
14992d61bbb3SSatish Balay }
1500218c64b6SSatish Balay 
15014a2ae208SSatish Balay #undef __FUNCT__
15024a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_N"
1503dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_N(Mat A,Vec xx,Vec yy,Vec zz)
15042d61bbb3SSatish Balay {
15052d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
1506d9ca1df4SBarry Smith   PetscScalar       *z = 0,*work,*workt,*zarray;
1507d9ca1df4SBarry Smith   const PetscScalar *x,*xb;
1508d9ca1df4SBarry Smith   const MatScalar   *v;
15096849ba73SBarry Smith   PetscErrorCode    ierr;
1510d9ca1df4SBarry Smith   PetscInt          mbs,i,bs=A->rmap->bs,j,n,bs2=a->bs2;
1511d9ca1df4SBarry Smith   PetscInt          ncols,k;
1512d9ca1df4SBarry Smith   const PetscInt    *ridx = NULL,*idx,*ii;
1513ace3abfcSBarry Smith   PetscBool         usecprow = a->compressedrow.use;
1514218c64b6SSatish Balay 
15152d61bbb3SSatish Balay   PetscFunctionBegin;
1516343551c4SBarry Smith   ierr = VecCopy(yy,zz);CHKERRQ(ierr);
1517d9ca1df4SBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
151826e093fcSHong Zhang   ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr);
15192d61bbb3SSatish Balay 
15202d61bbb3SSatish Balay   idx = a->j;
15212d61bbb3SSatish Balay   v   = a->a;
152226e093fcSHong Zhang   if (usecprow) {
152326e093fcSHong Zhang     mbs  = a->compressedrow.nrows;
152426e093fcSHong Zhang     ii   = a->compressedrow.i;
15257b2bb3b9SHong Zhang     ridx = a->compressedrow.rindex;
152626e093fcSHong Zhang   } else {
152726e093fcSHong Zhang     mbs = a->mbs;
15282d61bbb3SSatish Balay     ii  = a->i;
152926e093fcSHong Zhang     z   = zarray;
153026e093fcSHong Zhang   }
15312d61bbb3SSatish Balay 
15322d61bbb3SSatish Balay   if (!a->mult_work) {
1533d0f46423SBarry Smith     k    = PetscMax(A->rmap->n,A->cmap->n);
1534854ce69bSBarry Smith     ierr = PetscMalloc1(k+1,&a->mult_work);CHKERRQ(ierr);
15352d61bbb3SSatish Balay   }
15362d61bbb3SSatish Balay   work = a->mult_work;
15372d61bbb3SSatish Balay   for (i=0; i<mbs; i++) {
15382d61bbb3SSatish Balay     n     = ii[1] - ii[0]; ii++;
15392d61bbb3SSatish Balay     ncols = n*bs;
15402d61bbb3SSatish Balay     workt = work;
15412d61bbb3SSatish Balay     for (j=0; j<n; j++) {
15422d61bbb3SSatish Balay       xb = x + bs*(*idx++);
15432d61bbb3SSatish Balay       for (k=0; k<bs; k++) workt[k] = xb[k];
15442d61bbb3SSatish Balay       workt += bs;
15452d61bbb3SSatish Balay     }
15467b2bb3b9SHong Zhang     if (usecprow) z = zarray + bs*ridx[i];
154796b95a6bSBarry Smith     PetscKernel_w_gets_w_plus_Ar_times_v(bs,ncols,work,v,z);
154871044d3cSBarry Smith     /* BLASgemv_("N",&bs,&ncols,&_DOne,v,&bs,work,&_One,&_DOne,z,&_One); */
15492d61bbb3SSatish Balay     v += n*bs2;
155026fbe8dcSKarl Rupp     if (!usecprow) z += bs;
155126e093fcSHong Zhang   }
1552d9ca1df4SBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
155326e093fcSHong Zhang   ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr);
1554dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*a->nz*bs2);CHKERRQ(ierr);
15552d61bbb3SSatish Balay   PetscFunctionReturn(0);
15562d61bbb3SSatish Balay }
15572d61bbb3SSatish Balay 
15584a2ae208SSatish Balay #undef __FUNCT__
1559547795f9SHong Zhang #define __FUNCT__ "MatMultHermitianTranspose_SeqBAIJ"
1560547795f9SHong Zhang PetscErrorCode MatMultHermitianTranspose_SeqBAIJ(Mat A,Vec xx,Vec zz)
1561547795f9SHong Zhang {
1562547795f9SHong Zhang   PetscScalar    zero = 0.0;
1563547795f9SHong Zhang   PetscErrorCode ierr;
1564547795f9SHong Zhang 
1565547795f9SHong Zhang   PetscFunctionBegin;
1566547795f9SHong Zhang   ierr = VecSet(zz,zero);CHKERRQ(ierr);
1567547795f9SHong Zhang   ierr = MatMultHermitianTransposeAdd_SeqBAIJ(A,xx,zz,zz);CHKERRQ(ierr);
1568547795f9SHong Zhang   PetscFunctionReturn(0);
1569547795f9SHong Zhang }
1570547795f9SHong Zhang 
1571547795f9SHong Zhang #undef __FUNCT__
15724a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_SeqBAIJ"
1573dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_SeqBAIJ(Mat A,Vec xx,Vec zz)
15742d61bbb3SSatish Balay {
15753447b6efSHong Zhang   PetscScalar    zero = 0.0;
15766849ba73SBarry Smith   PetscErrorCode ierr;
15772d61bbb3SSatish Balay 
15782d61bbb3SSatish Balay   PetscFunctionBegin;
15792dcb1b2aSMatthew Knepley   ierr = VecSet(zz,zero);CHKERRQ(ierr);
15803447b6efSHong Zhang   ierr = MatMultTransposeAdd_SeqBAIJ(A,xx,zz,zz);CHKERRQ(ierr);
15812d61bbb3SSatish Balay   PetscFunctionReturn(0);
15822d61bbb3SSatish Balay }
15832d61bbb3SSatish Balay 
15844a2ae208SSatish Balay #undef __FUNCT__
1585547795f9SHong Zhang #define __FUNCT__ "MatMultHermitianTransposeAdd_SeqBAIJ"
1586547795f9SHong Zhang PetscErrorCode MatMultHermitianTransposeAdd_SeqBAIJ(Mat A,Vec xx,Vec yy,Vec zz)
1587547795f9SHong Zhang {
1588547795f9SHong Zhang   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
1589b8c08b77SHong Zhang   PetscScalar       *z,x1,x2,x3,x4,x5;
1590d9ca1df4SBarry Smith   const PetscScalar *x,*xb = NULL;
1591d9ca1df4SBarry Smith   const MatScalar   *v;
1592547795f9SHong Zhang   PetscErrorCode    ierr;
1593b8c08b77SHong Zhang   PetscInt          mbs,i,rval,bs=A->rmap->bs,j,n;
1594d9ca1df4SBarry Smith   const PetscInt    *idx,*ii,*ib,*ridx = NULL;
1595547795f9SHong Zhang   Mat_CompressedRow cprow = a->compressedrow;
1596ace3abfcSBarry Smith   PetscBool         usecprow = cprow.use;
1597547795f9SHong Zhang 
1598547795f9SHong Zhang   PetscFunctionBegin;
1599343551c4SBarry Smith   if (yy != zz) { ierr = VecCopy(yy,zz);CHKERRQ(ierr); }
1600d9ca1df4SBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
1601547795f9SHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
1602547795f9SHong Zhang 
1603547795f9SHong Zhang   idx = a->j;
1604547795f9SHong Zhang   v   = a->a;
1605547795f9SHong Zhang   if (usecprow) {
1606547795f9SHong Zhang     mbs  = cprow.nrows;
1607547795f9SHong Zhang     ii   = cprow.i;
1608547795f9SHong Zhang     ridx = cprow.rindex;
1609547795f9SHong Zhang   } else {
1610547795f9SHong Zhang     mbs=a->mbs;
1611547795f9SHong Zhang     ii = a->i;
1612547795f9SHong Zhang     xb = x;
1613547795f9SHong Zhang   }
1614547795f9SHong Zhang 
1615547795f9SHong Zhang   switch (bs) {
1616547795f9SHong Zhang   case 1:
1617547795f9SHong Zhang     for (i=0; i<mbs; i++) {
1618547795f9SHong Zhang       if (usecprow) xb = x + ridx[i];
1619547795f9SHong Zhang       x1 = xb[0];
1620547795f9SHong Zhang       ib = idx + ii[0];
1621547795f9SHong Zhang       n  = ii[1] - ii[0]; ii++;
1622547795f9SHong Zhang       for (j=0; j<n; j++) {
1623547795f9SHong Zhang         rval     = ib[j];
1624547795f9SHong Zhang         z[rval] += PetscConj(*v) * x1;
1625547795f9SHong Zhang         v++;
1626547795f9SHong Zhang       }
1627547795f9SHong Zhang       if (!usecprow) xb++;
1628547795f9SHong Zhang     }
1629547795f9SHong Zhang     break;
1630547795f9SHong Zhang   case 2:
1631547795f9SHong Zhang     for (i=0; i<mbs; i++) {
1632547795f9SHong Zhang       if (usecprow) xb = x + 2*ridx[i];
1633547795f9SHong Zhang       x1 = xb[0]; x2 = xb[1];
1634547795f9SHong Zhang       ib = idx + ii[0];
1635547795f9SHong Zhang       n  = ii[1] - ii[0]; ii++;
1636547795f9SHong Zhang       for (j=0; j<n; j++) {
1637547795f9SHong Zhang         rval       = ib[j]*2;
1638547795f9SHong Zhang         z[rval++] += PetscConj(v[0])*x1 + PetscConj(v[1])*x2;
1639547795f9SHong Zhang         z[rval++] += PetscConj(v[2])*x1 + PetscConj(v[3])*x2;
1640547795f9SHong Zhang         v         += 4;
1641547795f9SHong Zhang       }
1642547795f9SHong Zhang       if (!usecprow) xb += 2;
1643547795f9SHong Zhang     }
1644547795f9SHong Zhang     break;
1645547795f9SHong Zhang   case 3:
1646547795f9SHong Zhang     for (i=0; i<mbs; i++) {
1647547795f9SHong Zhang       if (usecprow) xb = x + 3*ridx[i];
1648547795f9SHong Zhang       x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
1649547795f9SHong Zhang       ib = idx + ii[0];
1650547795f9SHong Zhang       n  = ii[1] - ii[0]; ii++;
1651547795f9SHong Zhang       for (j=0; j<n; j++) {
1652547795f9SHong Zhang         rval       = ib[j]*3;
1653547795f9SHong Zhang         z[rval++] += PetscConj(v[0])*x1 + PetscConj(v[1])*x2 + PetscConj(v[2])*x3;
1654547795f9SHong Zhang         z[rval++] += PetscConj(v[3])*x1 + PetscConj(v[4])*x2 + PetscConj(v[5])*x3;
1655547795f9SHong Zhang         z[rval++] += PetscConj(v[6])*x1 + PetscConj(v[7])*x2 + PetscConj(v[8])*x3;
1656547795f9SHong Zhang         v         += 9;
1657547795f9SHong Zhang       }
1658547795f9SHong Zhang       if (!usecprow) xb += 3;
1659547795f9SHong Zhang     }
1660547795f9SHong Zhang     break;
1661547795f9SHong Zhang   case 4:
1662547795f9SHong Zhang     for (i=0; i<mbs; i++) {
1663547795f9SHong Zhang       if (usecprow) xb = x + 4*ridx[i];
1664547795f9SHong Zhang       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
1665547795f9SHong Zhang       ib = idx + ii[0];
1666547795f9SHong Zhang       n  = ii[1] - ii[0]; ii++;
1667547795f9SHong Zhang       for (j=0; j<n; j++) {
1668547795f9SHong Zhang         rval       = ib[j]*4;
1669547795f9SHong Zhang         z[rval++] +=  PetscConj(v[0])*x1 + PetscConj(v[1])*x2  + PetscConj(v[2])*x3  + PetscConj(v[3])*x4;
1670547795f9SHong Zhang         z[rval++] +=  PetscConj(v[4])*x1 + PetscConj(v[5])*x2  + PetscConj(v[6])*x3  + PetscConj(v[7])*x4;
1671547795f9SHong Zhang         z[rval++] +=  PetscConj(v[8])*x1 + PetscConj(v[9])*x2  + PetscConj(v[10])*x3 + PetscConj(v[11])*x4;
1672547795f9SHong Zhang         z[rval++] += PetscConj(v[12])*x1 + PetscConj(v[13])*x2 + PetscConj(v[14])*x3 + PetscConj(v[15])*x4;
1673547795f9SHong Zhang         v         += 16;
1674547795f9SHong Zhang       }
1675547795f9SHong Zhang       if (!usecprow) xb += 4;
1676547795f9SHong Zhang     }
1677547795f9SHong Zhang     break;
1678547795f9SHong Zhang   case 5:
1679547795f9SHong Zhang     for (i=0; i<mbs; i++) {
1680547795f9SHong Zhang       if (usecprow) xb = x + 5*ridx[i];
1681547795f9SHong Zhang       x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
1682547795f9SHong Zhang       x4 = xb[3]; x5 = xb[4];
1683547795f9SHong Zhang       ib = idx + ii[0];
1684547795f9SHong Zhang       n  = ii[1] - ii[0]; ii++;
1685547795f9SHong Zhang       for (j=0; j<n; j++) {
1686547795f9SHong Zhang         rval       = ib[j]*5;
1687547795f9SHong Zhang         z[rval++] +=  PetscConj(v[0])*x1 +  PetscConj(v[1])*x2 +  PetscConj(v[2])*x3 +  PetscConj(v[3])*x4 +  PetscConj(v[4])*x5;
1688547795f9SHong Zhang         z[rval++] +=  PetscConj(v[5])*x1 +  PetscConj(v[6])*x2 +  PetscConj(v[7])*x3 +  PetscConj(v[8])*x4 +  PetscConj(v[9])*x5;
1689547795f9SHong Zhang         z[rval++] += PetscConj(v[10])*x1 + PetscConj(v[11])*x2 + PetscConj(v[12])*x3 + PetscConj(v[13])*x4 + PetscConj(v[14])*x5;
1690547795f9SHong Zhang         z[rval++] += PetscConj(v[15])*x1 + PetscConj(v[16])*x2 + PetscConj(v[17])*x3 + PetscConj(v[18])*x4 + PetscConj(v[19])*x5;
1691547795f9SHong Zhang         z[rval++] += PetscConj(v[20])*x1 + PetscConj(v[21])*x2 + PetscConj(v[22])*x3 + PetscConj(v[23])*x4 + PetscConj(v[24])*x5;
1692547795f9SHong Zhang         v         += 25;
1693547795f9SHong Zhang       }
1694547795f9SHong Zhang       if (!usecprow) xb += 5;
1695547795f9SHong Zhang     }
1696547795f9SHong Zhang     break;
1697968ae2c8SSatish Balay   default: /* block sizes larger than 5 by 5 are handled by BLAS */
1698968ae2c8SSatish Balay     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"block size larger than 5 is not supported yet");
1699968ae2c8SSatish Balay #if 0
1700968ae2c8SSatish Balay     {
1701b8c08b77SHong Zhang       PetscInt          ncols,k,bs2=a->bs2;
1702b8c08b77SHong Zhang       PetscScalar       *work,*workt,zb;
1703d9ca1df4SBarry Smith       const PetscScalar *xtmp;
1704547795f9SHong Zhang       if (!a->mult_work) {
1705547795f9SHong Zhang         k    = PetscMax(A->rmap->n,A->cmap->n);
1706854ce69bSBarry Smith         ierr = PetscMalloc1(k+1,&a->mult_work);CHKERRQ(ierr);
1707547795f9SHong Zhang       }
1708547795f9SHong Zhang       work = a->mult_work;
1709547795f9SHong Zhang       xtmp = x;
1710547795f9SHong Zhang       for (i=0; i<mbs; i++) {
1711547795f9SHong Zhang         n     = ii[1] - ii[0]; ii++;
1712547795f9SHong Zhang         ncols = n*bs;
1713547795f9SHong Zhang         ierr  = PetscMemzero(work,ncols*sizeof(PetscScalar));CHKERRQ(ierr);
171426fbe8dcSKarl Rupp         if (usecprow) xtmp = x + bs*ridx[i];
171596b95a6bSBarry Smith         PetscKernel_w_gets_w_plus_trans_Ar_times_v(bs,ncols,xtmp,v,work);
1716547795f9SHong Zhang         /* BLASgemv_("T",&bs,&ncols,&_DOne,v,&bs,xtmp,&_One,&_DOne,work,&_One); */
1717547795f9SHong Zhang         v += n*bs2;
1718547795f9SHong Zhang         if (!usecprow) xtmp += bs;
1719547795f9SHong Zhang         workt = work;
1720547795f9SHong Zhang         for (j=0; j<n; j++) {
1721547795f9SHong Zhang           zb = z + bs*(*idx++);
1722547795f9SHong Zhang           for (k=0; k<bs; k++) zb[k] += workt[k] ;
1723547795f9SHong Zhang           workt += bs;
1724547795f9SHong Zhang         }
1725547795f9SHong Zhang       }
1726547795f9SHong Zhang     }
1727968ae2c8SSatish Balay #endif
1728547795f9SHong Zhang   }
1729d9ca1df4SBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
1730547795f9SHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
1731547795f9SHong Zhang   ierr = PetscLogFlops(2.0*a->nz*a->bs2);CHKERRQ(ierr);
1732547795f9SHong Zhang   PetscFunctionReturn(0);
1733547795f9SHong Zhang }
1734547795f9SHong Zhang 
1735547795f9SHong Zhang #undef __FUNCT__
17364a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_SeqBAIJ"
1737dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_SeqBAIJ(Mat A,Vec xx,Vec yy,Vec zz)
17382d61bbb3SSatish Balay {
17392d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
1740d9ca1df4SBarry Smith   PetscScalar       *zb,*z,x1,x2,x3,x4,x5;
1741d9ca1df4SBarry Smith   const PetscScalar *x,*xb = 0;
1742d9ca1df4SBarry Smith   const MatScalar   *v;
17436849ba73SBarry Smith   PetscErrorCode    ierr;
1744d9ca1df4SBarry Smith   PetscInt          mbs,i,rval,bs=A->rmap->bs,j,n,bs2=a->bs2;
1745d9ca1df4SBarry Smith   const PetscInt    *idx,*ii,*ib,*ridx = NULL;
17463447b6efSHong Zhang   Mat_CompressedRow cprow   = a->compressedrow;
1747ace3abfcSBarry Smith   PetscBool         usecprow=cprow.use;
17482d61bbb3SSatish Balay 
17492d61bbb3SSatish Balay   PetscFunctionBegin;
1750343551c4SBarry Smith   if (yy != zz) { ierr = VecCopy(yy,zz);CHKERRQ(ierr); }
1751d9ca1df4SBarry Smith   ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr);
17523447b6efSHong Zhang   ierr = VecGetArray(zz,&z);CHKERRQ(ierr);
17532d61bbb3SSatish Balay 
17542d61bbb3SSatish Balay   idx = a->j;
17552d61bbb3SSatish Balay   v   = a->a;
17563447b6efSHong Zhang   if (usecprow) {
17573447b6efSHong Zhang     mbs  = cprow.nrows;
17583447b6efSHong Zhang     ii   = cprow.i;
17597b2bb3b9SHong Zhang     ridx = cprow.rindex;
17603447b6efSHong Zhang   } else {
17613447b6efSHong Zhang     mbs=a->mbs;
17622d61bbb3SSatish Balay     ii = a->i;
1763f1af5d2fSBarry Smith     xb = x;
17643447b6efSHong Zhang   }
17652d61bbb3SSatish Balay 
17662d61bbb3SSatish Balay   switch (bs) {
17672d61bbb3SSatish Balay   case 1:
17682d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
17697b2bb3b9SHong Zhang       if (usecprow) xb = x + ridx[i];
1770f1af5d2fSBarry Smith       x1 = xb[0];
17713447b6efSHong Zhang       ib = idx + ii[0];
17723447b6efSHong Zhang       n  = ii[1] - ii[0]; ii++;
17732d61bbb3SSatish Balay       for (j=0; j<n; j++) {
17742d61bbb3SSatish Balay         rval     = ib[j];
1775f1af5d2fSBarry Smith         z[rval] += *v * x1;
1776f1af5d2fSBarry Smith         v++;
17772d61bbb3SSatish Balay       }
17783447b6efSHong Zhang       if (!usecprow) xb++;
17792d61bbb3SSatish Balay     }
17802d61bbb3SSatish Balay     break;
17812d61bbb3SSatish Balay   case 2:
17822d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
17837b2bb3b9SHong Zhang       if (usecprow) xb = x + 2*ridx[i];
1784f1af5d2fSBarry Smith       x1 = xb[0]; x2 = xb[1];
17853447b6efSHong Zhang       ib = idx + ii[0];
17863447b6efSHong Zhang       n  = ii[1] - ii[0]; ii++;
17872d61bbb3SSatish Balay       for (j=0; j<n; j++) {
17882d61bbb3SSatish Balay         rval       = ib[j]*2;
17892d61bbb3SSatish Balay         z[rval++] += v[0]*x1 + v[1]*x2;
17902d61bbb3SSatish Balay         z[rval++] += v[2]*x1 + v[3]*x2;
17912d61bbb3SSatish Balay         v         += 4;
17922d61bbb3SSatish Balay       }
17933447b6efSHong Zhang       if (!usecprow) xb += 2;
17942d61bbb3SSatish Balay     }
17952d61bbb3SSatish Balay     break;
17962d61bbb3SSatish Balay   case 3:
17972d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
17987b2bb3b9SHong Zhang       if (usecprow) xb = x + 3*ridx[i];
1799f1af5d2fSBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
18003447b6efSHong Zhang       ib = idx + ii[0];
18013447b6efSHong Zhang       n  = ii[1] - ii[0]; ii++;
18022d61bbb3SSatish Balay       for (j=0; j<n; j++) {
18032d61bbb3SSatish Balay         rval       = ib[j]*3;
18042d61bbb3SSatish Balay         z[rval++] += v[0]*x1 + v[1]*x2 + v[2]*x3;
18052d61bbb3SSatish Balay         z[rval++] += v[3]*x1 + v[4]*x2 + v[5]*x3;
18062d61bbb3SSatish Balay         z[rval++] += v[6]*x1 + v[7]*x2 + v[8]*x3;
18072d61bbb3SSatish Balay         v         += 9;
18082d61bbb3SSatish Balay       }
18093447b6efSHong Zhang       if (!usecprow) xb += 3;
18102d61bbb3SSatish Balay     }
18112d61bbb3SSatish Balay     break;
18122d61bbb3SSatish Balay   case 4:
18132d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
18147b2bb3b9SHong Zhang       if (usecprow) xb = x + 4*ridx[i];
1815f1af5d2fSBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3];
18163447b6efSHong Zhang       ib = idx + ii[0];
18173447b6efSHong Zhang       n  = ii[1] - ii[0]; ii++;
18182d61bbb3SSatish Balay       for (j=0; j<n; j++) {
18192d61bbb3SSatish Balay         rval       = ib[j]*4;
18202d61bbb3SSatish Balay         z[rval++] +=  v[0]*x1 +  v[1]*x2 +  v[2]*x3 +  v[3]*x4;
18212d61bbb3SSatish Balay         z[rval++] +=  v[4]*x1 +  v[5]*x2 +  v[6]*x3 +  v[7]*x4;
18222d61bbb3SSatish Balay         z[rval++] +=  v[8]*x1 +  v[9]*x2 + v[10]*x3 + v[11]*x4;
18232d61bbb3SSatish Balay         z[rval++] += v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4;
18242d61bbb3SSatish Balay         v         += 16;
18252d61bbb3SSatish Balay       }
18263447b6efSHong Zhang       if (!usecprow) xb += 4;
18272d61bbb3SSatish Balay     }
18282d61bbb3SSatish Balay     break;
18292d61bbb3SSatish Balay   case 5:
18302d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
18317b2bb3b9SHong Zhang       if (usecprow) xb = x + 5*ridx[i];
1832f1af5d2fSBarry Smith       x1 = xb[0]; x2 = xb[1]; x3 = xb[2];
18332d61bbb3SSatish Balay       x4 = xb[3]; x5 = xb[4];
18343447b6efSHong Zhang       ib = idx + ii[0];
18353447b6efSHong Zhang       n  = ii[1] - ii[0]; ii++;
18362d61bbb3SSatish Balay       for (j=0; j<n; j++) {
18372d61bbb3SSatish Balay         rval       = ib[j]*5;
18382d61bbb3SSatish Balay         z[rval++] +=  v[0]*x1 +  v[1]*x2 +  v[2]*x3 +  v[3]*x4 +  v[4]*x5;
18392d61bbb3SSatish Balay         z[rval++] +=  v[5]*x1 +  v[6]*x2 +  v[7]*x3 +  v[8]*x4 +  v[9]*x5;
18402d61bbb3SSatish Balay         z[rval++] += v[10]*x1 + v[11]*x2 + v[12]*x3 + v[13]*x4 + v[14]*x5;
18412d61bbb3SSatish Balay         z[rval++] += v[15]*x1 + v[16]*x2 + v[17]*x3 + v[18]*x4 + v[19]*x5;
18422d61bbb3SSatish Balay         z[rval++] += v[20]*x1 + v[21]*x2 + v[22]*x3 + v[23]*x4 + v[24]*x5;
18432d61bbb3SSatish Balay         v         += 25;
18442d61bbb3SSatish Balay       }
18453447b6efSHong Zhang       if (!usecprow) xb += 5;
18462d61bbb3SSatish Balay     }
18472d61bbb3SSatish Balay     break;
1848f1af5d2fSBarry Smith   default: {      /* block sizes larger then 5 by 5 are handled by BLAS */
1849690b6cddSBarry Smith     PetscInt          ncols,k;
1850d9ca1df4SBarry Smith     PetscScalar       *work,*workt;
1851d9ca1df4SBarry Smith     const PetscScalar *xtmp;
18522d61bbb3SSatish Balay     if (!a->mult_work) {
1853d0f46423SBarry Smith       k    = PetscMax(A->rmap->n,A->cmap->n);
1854854ce69bSBarry Smith       ierr = PetscMalloc1(k+1,&a->mult_work);CHKERRQ(ierr);
18552d61bbb3SSatish Balay     }
18562d61bbb3SSatish Balay     work = a->mult_work;
18573447b6efSHong Zhang     xtmp = x;
18582d61bbb3SSatish Balay     for (i=0; i<mbs; i++) {
18592d61bbb3SSatish Balay       n     = ii[1] - ii[0]; ii++;
18602d61bbb3SSatish Balay       ncols = n*bs;
186187828ca2SBarry Smith       ierr  = PetscMemzero(work,ncols*sizeof(PetscScalar));CHKERRQ(ierr);
186226fbe8dcSKarl Rupp       if (usecprow) xtmp = x + bs*ridx[i];
186396b95a6bSBarry Smith       PetscKernel_w_gets_w_plus_trans_Ar_times_v(bs,ncols,xtmp,v,work);
186471044d3cSBarry Smith       /* BLASgemv_("T",&bs,&ncols,&_DOne,v,&bs,xtmp,&_One,&_DOne,work,&_One); */
18652d61bbb3SSatish Balay       v += n*bs2;
18663447b6efSHong Zhang       if (!usecprow) xtmp += bs;
18672d61bbb3SSatish Balay       workt = work;
18682d61bbb3SSatish Balay       for (j=0; j<n; j++) {
18692d61bbb3SSatish Balay         zb = z + bs*(*idx++);
18702d61bbb3SSatish Balay         for (k=0; k<bs; k++) zb[k] += workt[k];
18712d61bbb3SSatish Balay         workt += bs;
18722d61bbb3SSatish Balay       }
18732d61bbb3SSatish Balay     }
18742d61bbb3SSatish Balay     }
18752d61bbb3SSatish Balay   }
1876d9ca1df4SBarry Smith   ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr);
18773447b6efSHong Zhang   ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr);
1878dc0b31edSSatish Balay   ierr = PetscLogFlops(2.0*a->nz*a->bs2);CHKERRQ(ierr);
18792d61bbb3SSatish Balay   PetscFunctionReturn(0);
18802d61bbb3SSatish Balay }
18812d61bbb3SSatish Balay 
18824a2ae208SSatish Balay #undef __FUNCT__
18834a2ae208SSatish Balay #define __FUNCT__ "MatScale_SeqBAIJ"
1884f4df32b1SMatthew Knepley PetscErrorCode MatScale_SeqBAIJ(Mat inA,PetscScalar alpha)
18852d61bbb3SSatish Balay {
18862d61bbb3SSatish Balay   Mat_SeqBAIJ    *a      = (Mat_SeqBAIJ*)inA->data;
1887690b6cddSBarry Smith   PetscInt       totalnz = a->bs2*a->nz;
1888f4df32b1SMatthew Knepley   PetscScalar    oalpha  = alpha;
1889efee365bSSatish Balay   PetscErrorCode ierr;
1890c5df96a5SBarry Smith   PetscBLASInt   one = 1,tnz;
18912d61bbb3SSatish Balay 
18922d61bbb3SSatish Balay   PetscFunctionBegin;
1893c5df96a5SBarry Smith   ierr = PetscBLASIntCast(totalnz,&tnz);CHKERRQ(ierr);
18948b83055fSJed Brown   PetscStackCallBLAS("BLASscal",BLASscal_(&tnz,&oalpha,a->a,&one));
1895efee365bSSatish Balay   ierr = PetscLogFlops(totalnz);CHKERRQ(ierr);
18962d61bbb3SSatish Balay   PetscFunctionReturn(0);
18972d61bbb3SSatish Balay }
18982d61bbb3SSatish Balay 
18994a2ae208SSatish Balay #undef __FUNCT__
19004a2ae208SSatish Balay #define __FUNCT__ "MatNorm_SeqBAIJ"
1901dfbe8321SBarry Smith PetscErrorCode MatNorm_SeqBAIJ(Mat A,NormType type,PetscReal *norm)
19022d61bbb3SSatish Balay {
19038a62d963SHong Zhang   PetscErrorCode ierr;
19042d61bbb3SSatish Balay   Mat_SeqBAIJ    *a  = (Mat_SeqBAIJ*)A->data;
19053f1db9ecSBarry Smith   MatScalar      *v  = a->a;
1906329f5518SBarry Smith   PetscReal      sum = 0.0;
1907d0f46423SBarry Smith   PetscInt       i,j,k,bs=A->rmap->bs,nz=a->nz,bs2=a->bs2,k1;
19082d61bbb3SSatish Balay 
19092d61bbb3SSatish Balay   PetscFunctionBegin;
19102d61bbb3SSatish Balay   if (type == NORM_FROBENIUS) {
19112d61bbb3SSatish Balay     for (i=0; i< bs2*nz; i++) {
1912329f5518SBarry Smith       sum += PetscRealPart(PetscConj(*v)*(*v)); v++;
19132d61bbb3SSatish Balay     }
19148f1a2a5eSBarry Smith     *norm = PetscSqrtReal(sum);
191551f70360SJed Brown     ierr = PetscLogFlops(2*bs2*nz);CHKERRQ(ierr);
19168a62d963SHong Zhang   } else if (type == NORM_1) { /* maximum column sum */
19178a62d963SHong Zhang     PetscReal *tmp;
19188a62d963SHong Zhang     PetscInt  *bcol = a->j;
1919854ce69bSBarry Smith     ierr = PetscCalloc1(A->cmap->n+1,&tmp);CHKERRQ(ierr);
19208a62d963SHong Zhang     for (i=0; i<nz; i++) {
19218a62d963SHong Zhang       for (j=0; j<bs; j++) {
19228a62d963SHong Zhang         k1 = bs*(*bcol) + j; /* column index */
19238a62d963SHong Zhang         for (k=0; k<bs; k++) {
19248a62d963SHong Zhang           tmp[k1] += PetscAbsScalar(*v); v++;
19258a62d963SHong Zhang         }
19268a62d963SHong Zhang       }
19278a62d963SHong Zhang       bcol++;
19288a62d963SHong Zhang     }
19298a62d963SHong Zhang     *norm = 0.0;
1930d0f46423SBarry Smith     for (j=0; j<A->cmap->n; j++) {
19318a62d963SHong Zhang       if (tmp[j] > *norm) *norm = tmp[j];
19328a62d963SHong Zhang     }
19338a62d963SHong Zhang     ierr = PetscFree(tmp);CHKERRQ(ierr);
193451f70360SJed Brown     ierr = PetscLogFlops(PetscMax(bs2*nz-1,0));CHKERRQ(ierr);
1935596552b5SBarry Smith   } else if (type == NORM_INFINITY) { /* maximum row sum */
1936596552b5SBarry Smith     *norm = 0.0;
1937596552b5SBarry Smith     for (k=0; k<bs; k++) {
193874f84c7bSSatish Balay       for (j=0; j<a->mbs; j++) {
1939596552b5SBarry Smith         v   = a->a + bs2*a->i[j] + k;
1940596552b5SBarry Smith         sum = 0.0;
1941596552b5SBarry Smith         for (i=0; i<a->i[j+1]-a->i[j]; i++) {
19420e90e235SBarry Smith           for (k1=0; k1<bs; k1++) {
1943596552b5SBarry Smith             sum += PetscAbsScalar(*v);
1944596552b5SBarry Smith             v   += bs;
19452d61bbb3SSatish Balay           }
19460e90e235SBarry Smith         }
1947596552b5SBarry Smith         if (sum > *norm) *norm = sum;
1948596552b5SBarry Smith       }
1949596552b5SBarry Smith     }
195051f70360SJed Brown     ierr = PetscLogFlops(PetscMax(bs2*nz-1,0));CHKERRQ(ierr);
1951e7e72b3dSBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for this norm yet");
19522d61bbb3SSatish Balay   PetscFunctionReturn(0);
19532d61bbb3SSatish Balay }
19542d61bbb3SSatish Balay 
19552d61bbb3SSatish Balay 
19564a2ae208SSatish Balay #undef __FUNCT__
19574a2ae208SSatish Balay #define __FUNCT__ "MatEqual_SeqBAIJ"
1958ace3abfcSBarry Smith PetscErrorCode MatEqual_SeqBAIJ(Mat A,Mat B,PetscBool * flg)
19592d61bbb3SSatish Balay {
19602d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)B->data;
1961dfbe8321SBarry Smith   PetscErrorCode ierr;
19622d61bbb3SSatish Balay 
19632d61bbb3SSatish Balay   PetscFunctionBegin;
19642d61bbb3SSatish Balay   /* If the  matrix/block dimensions are not equal, or no of nonzeros or shift */
1965d0f46423SBarry Smith   if ((A->rmap->N != B->rmap->N) || (A->cmap->n != B->cmap->n) || (A->rmap->bs != B->rmap->bs)|| (a->nz != b->nz)) {
1966273d9f13SBarry Smith     *flg = PETSC_FALSE;
1967273d9f13SBarry Smith     PetscFunctionReturn(0);
19682d61bbb3SSatish Balay   }
19692d61bbb3SSatish Balay 
19702d61bbb3SSatish Balay   /* if the a->i are the same */
1971690b6cddSBarry Smith   ierr = PetscMemcmp(a->i,b->i,(a->mbs+1)*sizeof(PetscInt),flg);CHKERRQ(ierr);
197226fbe8dcSKarl Rupp   if (!*flg) PetscFunctionReturn(0);
19732d61bbb3SSatish Balay 
19742d61bbb3SSatish Balay   /* if a->j are the same */
1975690b6cddSBarry Smith   ierr = PetscMemcmp(a->j,b->j,(a->nz)*sizeof(PetscInt),flg);CHKERRQ(ierr);
197626fbe8dcSKarl Rupp   if (!*flg) PetscFunctionReturn(0);
197726fbe8dcSKarl Rupp 
19782d61bbb3SSatish Balay   /* if a->a are the same */
1979d0f46423SBarry Smith   ierr = PetscMemcmp(a->a,b->a,(a->nz)*(A->rmap->bs)*(B->rmap->bs)*sizeof(PetscScalar),flg);CHKERRQ(ierr);
19802d61bbb3SSatish Balay   PetscFunctionReturn(0);
19812d61bbb3SSatish Balay 
19822d61bbb3SSatish Balay }
19832d61bbb3SSatish Balay 
19844a2ae208SSatish Balay #undef __FUNCT__
19854a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_SeqBAIJ"
1986dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_SeqBAIJ(Mat A,Vec v)
19872d61bbb3SSatish Balay {
19882d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
1989dfbe8321SBarry Smith   PetscErrorCode ierr;
1990690b6cddSBarry Smith   PetscInt       i,j,k,n,row,bs,*ai,*aj,ambs,bs2;
199187828ca2SBarry Smith   PetscScalar    *x,zero = 0.0;
19923f1db9ecSBarry Smith   MatScalar      *aa,*aa_j;
19932d61bbb3SSatish Balay 
19942d61bbb3SSatish Balay   PetscFunctionBegin;
1995e32f2f54SBarry Smith   if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
1996d0f46423SBarry Smith   bs   = A->rmap->bs;
19972d61bbb3SSatish Balay   aa   = a->a;
19982d61bbb3SSatish Balay   ai   = a->i;
19992d61bbb3SSatish Balay   aj   = a->j;
20002d61bbb3SSatish Balay   ambs = a->mbs;
20012d61bbb3SSatish Balay   bs2  = a->bs2;
20022d61bbb3SSatish Balay 
20032dcb1b2aSMatthew Knepley   ierr = VecSet(v,zero);CHKERRQ(ierr);
20041ebc52fbSHong Zhang   ierr = VecGetArray(v,&x);CHKERRQ(ierr);
2005e1311b90SBarry Smith   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
2006e32f2f54SBarry Smith   if (n != A->rmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector");
20072d61bbb3SSatish Balay   for (i=0; i<ambs; i++) {
20082d61bbb3SSatish Balay     for (j=ai[i]; j<ai[i+1]; j++) {
20092d61bbb3SSatish Balay       if (aj[j] == i) {
20102d61bbb3SSatish Balay         row  = i*bs;
20112d61bbb3SSatish Balay         aa_j = aa+j*bs2;
20122d61bbb3SSatish Balay         for (k=0; k<bs2; k+=(bs+1),row++) x[row] = aa_j[k];
20132d61bbb3SSatish Balay         break;
20142d61bbb3SSatish Balay       }
20152d61bbb3SSatish Balay     }
20162d61bbb3SSatish Balay   }
20171ebc52fbSHong Zhang   ierr = VecRestoreArray(v,&x);CHKERRQ(ierr);
20182d61bbb3SSatish Balay   PetscFunctionReturn(0);
20192d61bbb3SSatish Balay }
20202d61bbb3SSatish Balay 
20214a2ae208SSatish Balay #undef __FUNCT__
20224a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_SeqBAIJ"
2023dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_SeqBAIJ(Mat A,Vec ll,Vec rr)
20242d61bbb3SSatish Balay {
20252d61bbb3SSatish Balay   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
202653ef36baSBarry Smith   const PetscScalar *l,*r,*li,*ri;
202753ef36baSBarry Smith   PetscScalar       x;
20283f1db9ecSBarry Smith   MatScalar         *aa, *v;
2029dfbe8321SBarry Smith   PetscErrorCode    ierr;
203053ef36baSBarry Smith   PetscInt          i,j,k,lm,rn,M,m,n,mbs,tmp,bs,bs2,iai;
203153ef36baSBarry Smith   const PetscInt    *ai,*aj;
20322d61bbb3SSatish Balay 
20332d61bbb3SSatish Balay   PetscFunctionBegin;
20342d61bbb3SSatish Balay   ai  = a->i;
20352d61bbb3SSatish Balay   aj  = a->j;
20362d61bbb3SSatish Balay   aa  = a->a;
2037d0f46423SBarry Smith   m   = A->rmap->n;
2038d0f46423SBarry Smith   n   = A->cmap->n;
2039d0f46423SBarry Smith   bs  = A->rmap->bs;
20402d61bbb3SSatish Balay   mbs = a->mbs;
20412d61bbb3SSatish Balay   bs2 = a->bs2;
20422d61bbb3SSatish Balay   if (ll) {
204353ef36baSBarry Smith     ierr = VecGetArrayRead(ll,&l);CHKERRQ(ierr);
2044e1311b90SBarry Smith     ierr = VecGetLocalSize(ll,&lm);CHKERRQ(ierr);
2045e32f2f54SBarry Smith     if (lm != m) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Left scaling vector wrong length");
20462d61bbb3SSatish Balay     for (i=0; i<mbs; i++) { /* for each block row */
20472d61bbb3SSatish Balay       M  = ai[i+1] - ai[i];
20482d61bbb3SSatish Balay       li = l + i*bs;
20492d61bbb3SSatish Balay       v  = aa + bs2*ai[i];
20502d61bbb3SSatish Balay       for (j=0; j<M; j++) { /* for each block */
20512d61bbb3SSatish Balay         for (k=0; k<bs2; k++) {
20522d61bbb3SSatish Balay           (*v++) *= li[k%bs];
20532d61bbb3SSatish Balay         }
20542d61bbb3SSatish Balay       }
20552d61bbb3SSatish Balay     }
205653ef36baSBarry Smith     ierr = VecRestoreArrayRead(ll,&l);CHKERRQ(ierr);
2057efee365bSSatish Balay     ierr = PetscLogFlops(a->nz);CHKERRQ(ierr);
20582d61bbb3SSatish Balay   }
20592d61bbb3SSatish Balay 
20602d61bbb3SSatish Balay   if (rr) {
206153ef36baSBarry Smith     ierr = VecGetArrayRead(rr,&r);CHKERRQ(ierr);
2062e1311b90SBarry Smith     ierr = VecGetLocalSize(rr,&rn);CHKERRQ(ierr);
2063e32f2f54SBarry Smith     if (rn != n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Right scaling vector wrong length");
20642d61bbb3SSatish Balay     for (i=0; i<mbs; i++) { /* for each block row */
206553ef36baSBarry Smith       iai = ai[i];
206653ef36baSBarry Smith       M   = ai[i+1] - iai;
206753ef36baSBarry Smith       v   = aa + bs2*iai;
20682d61bbb3SSatish Balay       for (j=0; j<M; j++) { /* for each block */
206953ef36baSBarry Smith         ri = r + bs*aj[iai+j];
20702d61bbb3SSatish Balay         for (k=0; k<bs; k++) {
20712d61bbb3SSatish Balay           x = ri[k];
207253ef36baSBarry Smith           for (tmp=0; tmp<bs; tmp++) v[tmp] *= x;
207353ef36baSBarry Smith           v += bs;
20742d61bbb3SSatish Balay         }
20752d61bbb3SSatish Balay       }
20762d61bbb3SSatish Balay     }
207753ef36baSBarry Smith     ierr = VecRestoreArrayRead(rr,&r);CHKERRQ(ierr);
2078efee365bSSatish Balay     ierr = PetscLogFlops(a->nz);CHKERRQ(ierr);
20792d61bbb3SSatish Balay   }
20802d61bbb3SSatish Balay   PetscFunctionReturn(0);
20812d61bbb3SSatish Balay }
20822d61bbb3SSatish Balay 
20832d61bbb3SSatish Balay 
20844a2ae208SSatish Balay #undef __FUNCT__
20854a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_SeqBAIJ"
2086dfbe8321SBarry Smith PetscErrorCode MatGetInfo_SeqBAIJ(Mat A,MatInfoType flag,MatInfo *info)
20872d61bbb3SSatish Balay {
20882d61bbb3SSatish Balay   Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
20892d61bbb3SSatish Balay 
20902d61bbb3SSatish Balay   PetscFunctionBegin;
20912d61bbb3SSatish Balay   info->block_size   = a->bs2;
2092ceed8ce5SJed Brown   info->nz_allocated = a->bs2*a->maxnz;
20932d61bbb3SSatish Balay   info->nz_used      = a->bs2*a->nz;
20942d61bbb3SSatish Balay   info->nz_unneeded  = (double)(info->nz_allocated - info->nz_used);
20952d61bbb3SSatish Balay   info->assemblies   = A->num_ass;
20968e58a170SBarry Smith   info->mallocs      = A->info.mallocs;
20977adad957SLisandro Dalcin   info->memory       = ((PetscObject)A)->mem;
2098d5f3da31SBarry Smith   if (A->factortype) {
20992d61bbb3SSatish Balay     info->fill_ratio_given  = A->info.fill_ratio_given;
21002d61bbb3SSatish Balay     info->fill_ratio_needed = A->info.fill_ratio_needed;
21012d61bbb3SSatish Balay     info->factor_mallocs    = A->info.factor_mallocs;
21022d61bbb3SSatish Balay   } else {
21032d61bbb3SSatish Balay     info->fill_ratio_given  = 0;
21042d61bbb3SSatish Balay     info->fill_ratio_needed = 0;
21052d61bbb3SSatish Balay     info->factor_mallocs    = 0;
21062d61bbb3SSatish Balay   }
21072d61bbb3SSatish Balay   PetscFunctionReturn(0);
21082d61bbb3SSatish Balay }
21092d61bbb3SSatish Balay 
21104a2ae208SSatish Balay #undef __FUNCT__
21114a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_SeqBAIJ"
2112dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_SeqBAIJ(Mat A)
21132d61bbb3SSatish Balay {
21142d61bbb3SSatish Balay   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data;
2115dfbe8321SBarry Smith   PetscErrorCode ierr;
21162d61bbb3SSatish Balay 
21172d61bbb3SSatish Balay   PetscFunctionBegin;
2118549d3d68SSatish Balay   ierr = PetscMemzero(a->a,a->bs2*a->i[a->mbs]*sizeof(MatScalar));CHKERRQ(ierr);
21192d61bbb3SSatish Balay   PetscFunctionReturn(0);
21202d61bbb3SSatish Balay }
2121