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 90736121d4SSatish Balay ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr); 91218c64b6SSatish Balay ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr); 92b9b97703SBarry Smith ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr); 93b9b97703SBarry Smith ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr); 94736121d4SSatish Balay 95854ce69bSBarry Smith ierr = PetscCalloc1(1+oldcols,&smap);CHKERRQ(ierr); 96736121d4SSatish Balay ssmap = smap; 97854ce69bSBarry Smith ierr = PetscMalloc1(1+nrows,&lens);CHKERRQ(ierr); 98736121d4SSatish Balay for (i=0; i<ncols; i++) smap[icol[i]] = i+1; 99736121d4SSatish Balay /* determine lens of each row */ 100736121d4SSatish Balay for (i=0; i<nrows; i++) { 101736121d4SSatish Balay kstart = ai[irow[i]]; 102736121d4SSatish Balay kend = kstart + a->ilen[irow[i]]; 103736121d4SSatish Balay lens[i] = 0; 104736121d4SSatish Balay for (k=kstart; k<kend; k++) { 10526fbe8dcSKarl Rupp if (ssmap[aj[k]]) lens[i]++; 106736121d4SSatish Balay } 107736121d4SSatish Balay } 108736121d4SSatish Balay /* Create and fill new matrix */ 109736121d4SSatish Balay if (scall == MAT_REUSE_MATRIX) { 110736121d4SSatish Balay c = (Mat_SeqBAIJ*)((*B)->data); 111736121d4SSatish Balay 112e32f2f54SBarry Smith if (c->mbs!=nrows || c->nbs!=ncols || (*B)->rmap->bs!=bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Submatrix wrong size"); 113690b6cddSBarry Smith ierr = PetscMemcmp(c->ilen,lens,c->mbs *sizeof(PetscInt),&flag);CHKERRQ(ierr); 114e7e72b3dSBarry Smith if (!flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong no of nonzeros"); 115690b6cddSBarry Smith ierr = PetscMemzero(c->ilen,c->mbs*sizeof(PetscInt));CHKERRQ(ierr); 116736121d4SSatish Balay C = *B; 1173a40ed3dSBarry Smith } else { 118ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)A),&C);CHKERRQ(ierr); 119f69a0ea3SMatthew Knepley ierr = MatSetSizes(C,nrows*bs,ncols*bs,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 1207adad957SLisandro Dalcin ierr = MatSetType(C,((PetscObject)A)->type_name);CHKERRQ(ierr); 121ab93d7beSBarry Smith ierr = MatSeqBAIJSetPreallocation_SeqBAIJ(C,bs,0,lens);CHKERRQ(ierr); 122736121d4SSatish Balay } 123736121d4SSatish Balay c = (Mat_SeqBAIJ*)(C->data); 124736121d4SSatish Balay for (i=0; i<nrows; i++) { 125736121d4SSatish Balay row = irow[i]; 126736121d4SSatish Balay kstart = ai[row]; 127736121d4SSatish Balay kend = kstart + a->ilen[row]; 128736121d4SSatish Balay mat_i = c->i[i]; 129736121d4SSatish Balay mat_j = c->j + mat_i; 130218c64b6SSatish Balay mat_a = c->a + mat_i*bs2; 131736121d4SSatish Balay mat_ilen = c->ilen + i; 132736121d4SSatish Balay for (k=kstart; k<kend; k++) { 133736121d4SSatish Balay if ((tcol=ssmap[a->j[k]])) { 134736121d4SSatish Balay *mat_j++ = tcol - 1; 135549d3d68SSatish Balay ierr = PetscMemcpy(mat_a,a->a+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 136549d3d68SSatish Balay mat_a += bs2; 137736121d4SSatish Balay (*mat_ilen)++; 138736121d4SSatish Balay } 139736121d4SSatish Balay } 140736121d4SSatish Balay } 141cdc6f3adSToby Isaac /* sort */ 142cdc6f3adSToby Isaac { 143cdc6f3adSToby Isaac MatScalar *work; 144cdc6f3adSToby Isaac 145cdc6f3adSToby Isaac ierr = PetscMalloc1(bs2,&work);CHKERRQ(ierr); 146cdc6f3adSToby Isaac for (i=0; i<nrows; i++) { 147cdc6f3adSToby Isaac PetscInt ilen; 148cdc6f3adSToby Isaac mat_i = c->i[i]; 149cdc6f3adSToby Isaac mat_j = c->j + mat_i; 150cdc6f3adSToby Isaac mat_a = c->a + mat_i*bs2; 151cdc6f3adSToby Isaac ilen = c->ilen[i]; 152cdc6f3adSToby Isaac ierr = PetscSortIntWithDataArray(ilen,mat_j,mat_a,bs2*sizeof(MatScalar),work);CHKERRQ(ierr); 153cdc6f3adSToby Isaac } 154cdc6f3adSToby Isaac ierr = PetscFree(work);CHKERRQ(ierr); 155cdc6f3adSToby Isaac } 156218c64b6SSatish Balay 157736121d4SSatish Balay /* Free work space */ 158736121d4SSatish Balay ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr); 159606d414cSSatish Balay ierr = PetscFree(smap);CHKERRQ(ierr); 160606d414cSSatish Balay ierr = PetscFree(lens);CHKERRQ(ierr); 1616d4a8577SBarry Smith ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1626d4a8577SBarry Smith ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 163736121d4SSatish Balay 164736121d4SSatish Balay ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr); 165736121d4SSatish Balay *B = C; 1663a40ed3dSBarry Smith PetscFunctionReturn(0); 167736121d4SSatish Balay } 168736121d4SSatish Balay 1694a2ae208SSatish Balay #undef __FUNCT__ 1704a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqBAIJ" 1714aa3045dSJed Brown PetscErrorCode MatGetSubMatrix_SeqBAIJ(Mat A,IS isrow,IS iscol,MatReuse scall,Mat *B) 172218c64b6SSatish Balay { 173218c64b6SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 174218c64b6SSatish Balay IS is1,is2; 1756849ba73SBarry Smith PetscErrorCode ierr; 176*afebec48SHong Zhang PetscInt *vary,*iary,nrows,ncols,i,bs=A->rmap->bs,count,maxmnbs,j; 1775d0c19d7SBarry Smith const PetscInt *irow,*icol; 178218c64b6SSatish Balay 1793a40ed3dSBarry Smith PetscFunctionBegin; 180218c64b6SSatish Balay ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr); 181218c64b6SSatish Balay ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr); 182b9b97703SBarry Smith ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr); 183b9b97703SBarry Smith ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr); 184218c64b6SSatish Balay 185218c64b6SSatish Balay /* Verify if the indices corespond to each element in a block 186218c64b6SSatish Balay and form the IS with compressed IS */ 187f8ecb639SStefano Zampini maxmnbs = PetscMax(a->mbs,a->nbs); 188f8ecb639SStefano Zampini ierr = PetscMalloc2(maxmnbs,&vary,maxmnbs,&iary);CHKERRQ(ierr); 189fca92195SBarry Smith ierr = PetscMemzero(vary,a->mbs*sizeof(PetscInt));CHKERRQ(ierr); 190218c64b6SSatish Balay for (i=0; i<nrows; i++) vary[irow[i]/bs]++; 191218c64b6SSatish Balay for (i=0; i<a->mbs; i++) { 192e32f2f54SBarry Smith if (vary[i]!=0 && vary[i]!=bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Index set does not match blocks"); 1936041f1b1SToby Isaac } 1946041f1b1SToby Isaac count = 0; 1956041f1b1SToby Isaac for (i=0; i<nrows; i++) { 196*afebec48SHong Zhang j = irow[i] / bs; 1976041f1b1SToby Isaac if ((vary[j]--)==bs) iary[count++] = j; 198218c64b6SSatish Balay } 19970b3c8c7SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,count,iary,PETSC_COPY_VALUES,&is1);CHKERRQ(ierr); 200218c64b6SSatish Balay 201f8ecb639SStefano Zampini ierr = PetscMemzero(vary,(a->nbs)*sizeof(PetscInt));CHKERRQ(ierr); 202218c64b6SSatish Balay for (i=0; i<ncols; i++) vary[icol[i]/bs]++; 203f8ecb639SStefano Zampini for (i=0; i<a->nbs; i++) { 204e32f2f54SBarry Smith if (vary[i]!=0 && vary[i]!=bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Internal error in PETSc"); 205218c64b6SSatish Balay } 2066041f1b1SToby Isaac count = 0; 2076041f1b1SToby Isaac for (i=0; i<ncols; i++) { 208*afebec48SHong Zhang j = icol[i] / bs; 2096041f1b1SToby Isaac if ((vary[j]--)==bs) iary[count++] = j; 2106041f1b1SToby Isaac } 21170b3c8c7SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,count,iary,PETSC_COPY_VALUES,&is2);CHKERRQ(ierr); 212218c64b6SSatish Balay ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr); 213218c64b6SSatish Balay ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr); 214fca92195SBarry Smith ierr = PetscFree2(vary,iary);CHKERRQ(ierr); 215218c64b6SSatish Balay 2164aa3045dSJed Brown ierr = MatGetSubMatrix_SeqBAIJ_Private(A,is1,is2,scall,B);CHKERRQ(ierr); 2176bf464f9SBarry Smith ierr = ISDestroy(&is1);CHKERRQ(ierr); 2186bf464f9SBarry Smith ierr = ISDestroy(&is2);CHKERRQ(ierr); 2193a40ed3dSBarry Smith PetscFunctionReturn(0); 220218c64b6SSatish Balay } 221218c64b6SSatish Balay 2224a2ae208SSatish Balay #undef __FUNCT__ 2234a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_SeqBAIJ" 224690b6cddSBarry Smith PetscErrorCode MatGetSubMatrices_SeqBAIJ(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[]) 225736121d4SSatish Balay { 2266849ba73SBarry Smith PetscErrorCode ierr; 227690b6cddSBarry Smith PetscInt i; 228736121d4SSatish Balay 2293a40ed3dSBarry Smith PetscFunctionBegin; 230736121d4SSatish Balay if (scall == MAT_INITIAL_MATRIX) { 231854ce69bSBarry Smith ierr = PetscMalloc1(n+1,B);CHKERRQ(ierr); 232736121d4SSatish Balay } 233736121d4SSatish Balay 234736121d4SSatish Balay for (i=0; i<n; i++) { 2354aa3045dSJed Brown ierr = MatGetSubMatrix_SeqBAIJ(A,irow[i],icol[i],scall,&(*B)[i]);CHKERRQ(ierr); 236736121d4SSatish Balay } 2373a40ed3dSBarry Smith PetscFunctionReturn(0); 238736121d4SSatish Balay } 239218c64b6SSatish Balay 240218c64b6SSatish Balay 2412d61bbb3SSatish Balay /* -------------------------------------------------------*/ 2422d61bbb3SSatish Balay /* Should check that shapes of vectors and matrices match */ 2432d61bbb3SSatish Balay /* -------------------------------------------------------*/ 2442d61bbb3SSatish Balay 2454a2ae208SSatish Balay #undef __FUNCT__ 2464a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_1" 247dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_1(Mat A,Vec xx,Vec zz) 2482d61bbb3SSatish Balay { 2492d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 250d9fead3dSBarry Smith PetscScalar *z,sum; 251d9fead3dSBarry Smith const PetscScalar *x; 252d9fead3dSBarry Smith const MatScalar *v; 2536849ba73SBarry Smith PetscErrorCode ierr; 2547c565772SBarry Smith PetscInt mbs,i,n; 2550298fd71SBarry Smith const PetscInt *idx,*ii,*ridx=NULL; 256ace3abfcSBarry Smith PetscBool usecprow=a->compressedrow.use; 2572d61bbb3SSatish Balay 2582d61bbb3SSatish Balay PetscFunctionBegin; 2593649974fSBarry Smith ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr); 2601ebc52fbSHong Zhang ierr = VecGetArray(zz,&z);CHKERRQ(ierr); 2612d61bbb3SSatish Balay 26226e093fcSHong Zhang if (usecprow) { 26326e093fcSHong Zhang mbs = a->compressedrow.nrows; 26426e093fcSHong Zhang ii = a->compressedrow.i; 2657b2bb3b9SHong Zhang ridx = a->compressedrow.rindex; 266a1c3900fSBarry Smith ierr = PetscMemzero(z,mbs*sizeof(PetscScalar));CHKERRQ(ierr); 26726e093fcSHong Zhang } else { 26826e093fcSHong Zhang mbs = a->mbs; 2692d61bbb3SSatish Balay ii = a->i; 27026e093fcSHong Zhang } 2712d61bbb3SSatish Balay 2722d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 273ee54c7eeSHong Zhang n = ii[1] - ii[0]; 274ee54c7eeSHong Zhang v = a->a + ii[0]; 275ee54c7eeSHong Zhang idx = a->j + ii[0]; 276ee54c7eeSHong Zhang ii++; 277444d8c10SJed Brown PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA); /* Indices for the next row (assumes same size as this one) */ 278444d8c10SJed Brown PetscPrefetchBlock(v+1*n,1*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */ 2792d61bbb3SSatish Balay sum = 0.0; 2802162cab8SBarry Smith PetscSparseDensePlusDot(sum,x,v,idx,n); 28126e093fcSHong Zhang if (usecprow) { 2827b2bb3b9SHong Zhang z[ridx[i]] = sum; 28326e093fcSHong Zhang } else { 2842d61bbb3SSatish Balay z[i] = sum; 2852d61bbb3SSatish Balay } 28626e093fcSHong Zhang } 2873649974fSBarry Smith ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr); 2881ebc52fbSHong Zhang ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr); 2897c565772SBarry Smith ierr = PetscLogFlops(2.0*a->nz - a->nonzerorowcnt);CHKERRQ(ierr); 2902d61bbb3SSatish Balay PetscFunctionReturn(0); 2912d61bbb3SSatish Balay } 2922d61bbb3SSatish Balay 2934a2ae208SSatish Balay #undef __FUNCT__ 2944a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_2" 295dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_2(Mat A,Vec xx,Vec zz) 2962d61bbb3SSatish Balay { 2972d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 298d9fead3dSBarry Smith PetscScalar *z = 0,sum1,sum2,*zarray; 299d9fead3dSBarry Smith const PetscScalar *x,*xb; 30087828ca2SBarry Smith PetscScalar x1,x2; 301d9fead3dSBarry Smith const MatScalar *v; 302dfbe8321SBarry Smith PetscErrorCode ierr; 3037c565772SBarry Smith PetscInt mbs,i,*idx,*ii,j,n,*ridx=NULL; 304ace3abfcSBarry Smith PetscBool usecprow=a->compressedrow.use; 3052d61bbb3SSatish Balay 3062d61bbb3SSatish Balay PetscFunctionBegin; 3073649974fSBarry Smith ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr); 30826e093fcSHong Zhang ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr); 3092d61bbb3SSatish Balay 3102d61bbb3SSatish Balay idx = a->j; 3112d61bbb3SSatish Balay v = a->a; 31226e093fcSHong Zhang if (usecprow) { 31326e093fcSHong Zhang mbs = a->compressedrow.nrows; 31426e093fcSHong Zhang ii = a->compressedrow.i; 3157b2bb3b9SHong Zhang ridx = a->compressedrow.rindex; 31626e093fcSHong Zhang } else { 31726e093fcSHong Zhang mbs = a->mbs; 3182d61bbb3SSatish Balay ii = a->i; 31926e093fcSHong Zhang z = zarray; 32026e093fcSHong Zhang } 3212d61bbb3SSatish Balay 3222d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 3232d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 3242d61bbb3SSatish Balay sum1 = 0.0; sum2 = 0.0; 325444d8c10SJed Brown PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA); /* Indices for the next row (assumes same size as this one) */ 326444d8c10SJed Brown PetscPrefetchBlock(v+4*n,4*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */ 3272d61bbb3SSatish Balay for (j=0; j<n; j++) { 3282d61bbb3SSatish Balay xb = x + 2*(*idx++); x1 = xb[0]; x2 = xb[1]; 3292d61bbb3SSatish Balay sum1 += v[0]*x1 + v[2]*x2; 3302d61bbb3SSatish Balay sum2 += v[1]*x1 + v[3]*x2; 3312d61bbb3SSatish Balay v += 4; 3322d61bbb3SSatish Balay } 3337b2bb3b9SHong Zhang if (usecprow) z = zarray + 2*ridx[i]; 3342d61bbb3SSatish Balay z[0] = sum1; z[1] = sum2; 33526e093fcSHong Zhang if (!usecprow) z += 2; 3362d61bbb3SSatish Balay } 3373649974fSBarry Smith ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr); 33826e093fcSHong Zhang ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr); 3397c565772SBarry Smith ierr = PetscLogFlops(8.0*a->nz - 2.0*a->nonzerorowcnt);CHKERRQ(ierr); 3402d61bbb3SSatish Balay PetscFunctionReturn(0); 3412d61bbb3SSatish Balay } 3422d61bbb3SSatish Balay 3434a2ae208SSatish Balay #undef __FUNCT__ 3444a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_3" 345dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_3(Mat A,Vec xx,Vec zz) 3462d61bbb3SSatish Balay { 3472d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 348d9fead3dSBarry Smith PetscScalar *z = 0,sum1,sum2,sum3,x1,x2,x3,*zarray; 349d9fead3dSBarry Smith const PetscScalar *x,*xb; 350d9fead3dSBarry Smith const MatScalar *v; 351dfbe8321SBarry Smith PetscErrorCode ierr; 3527c565772SBarry Smith PetscInt mbs,i,*idx,*ii,j,n,*ridx=NULL; 353ace3abfcSBarry Smith PetscBool usecprow=a->compressedrow.use; 35426e093fcSHong Zhang 355b6410449SSatish Balay #if defined(PETSC_HAVE_PRAGMA_DISJOINT) 356fee21e36SBarry Smith #pragma disjoint(*v,*z,*xb) 357fee21e36SBarry Smith #endif 358fee21e36SBarry Smith 3592d61bbb3SSatish Balay PetscFunctionBegin; 3603649974fSBarry Smith ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr); 36126e093fcSHong Zhang ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr); 3622d61bbb3SSatish Balay 3632d61bbb3SSatish Balay idx = a->j; 3642d61bbb3SSatish Balay v = a->a; 36526e093fcSHong Zhang if (usecprow) { 36626e093fcSHong Zhang mbs = a->compressedrow.nrows; 36726e093fcSHong Zhang ii = a->compressedrow.i; 3687b2bb3b9SHong Zhang ridx = a->compressedrow.rindex; 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; 42326e093fcSHong Zhang } else { 42426e093fcSHong Zhang mbs = a->mbs; 4252d61bbb3SSatish Balay ii = a->i; 42626e093fcSHong Zhang z = zarray; 42726e093fcSHong Zhang } 4282d61bbb3SSatish Balay 4292d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 43026fbe8dcSKarl Rupp n = ii[1] - ii[0]; 43126fbe8dcSKarl Rupp ii++; 43226fbe8dcSKarl Rupp sum1 = 0.0; 43326fbe8dcSKarl Rupp sum2 = 0.0; 43426fbe8dcSKarl Rupp sum3 = 0.0; 43526fbe8dcSKarl Rupp sum4 = 0.0; 43626fbe8dcSKarl Rupp 437444d8c10SJed Brown PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA); /* Indices for the next row (assumes same size as this one) */ 438444d8c10SJed Brown PetscPrefetchBlock(v+16*n,16*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */ 4392d61bbb3SSatish Balay for (j=0; j<n; j++) { 4402d61bbb3SSatish Balay xb = x + 4*(*idx++); 4412d61bbb3SSatish Balay x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; 4422d61bbb3SSatish Balay sum1 += v[0]*x1 + v[4]*x2 + v[8]*x3 + v[12]*x4; 4432d61bbb3SSatish Balay sum2 += v[1]*x1 + v[5]*x2 + v[9]*x3 + v[13]*x4; 4442d61bbb3SSatish Balay sum3 += v[2]*x1 + v[6]*x2 + v[10]*x3 + v[14]*x4; 4452d61bbb3SSatish Balay sum4 += v[3]*x1 + v[7]*x2 + v[11]*x3 + v[15]*x4; 4462d61bbb3SSatish Balay v += 16; 4472d61bbb3SSatish Balay } 4487b2bb3b9SHong Zhang if (usecprow) z = zarray + 4*ridx[i]; 4492d61bbb3SSatish Balay z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; 45026e093fcSHong Zhang if (!usecprow) z += 4; 4512d61bbb3SSatish Balay } 4523649974fSBarry Smith ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr); 45326e093fcSHong Zhang ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr); 4547c565772SBarry Smith ierr = PetscLogFlops(32.0*a->nz - 4.0*a->nonzerorowcnt);CHKERRQ(ierr); 4552d61bbb3SSatish Balay PetscFunctionReturn(0); 4562d61bbb3SSatish Balay } 4572d61bbb3SSatish Balay 4584a2ae208SSatish Balay #undef __FUNCT__ 4594a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_5" 460dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_5(Mat A,Vec xx,Vec zz) 4612d61bbb3SSatish Balay { 4622d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 463d9fead3dSBarry Smith PetscScalar sum1,sum2,sum3,sum4,sum5,x1,x2,x3,x4,x5,*z = 0,*zarray; 464d9fead3dSBarry Smith const PetscScalar *xb,*x; 465d9fead3dSBarry Smith const MatScalar *v; 466dfbe8321SBarry Smith PetscErrorCode ierr; 4670298fd71SBarry Smith const PetscInt *idx,*ii,*ridx=NULL; 4687c565772SBarry Smith PetscInt mbs,i,j,n; 469ace3abfcSBarry Smith PetscBool usecprow=a->compressedrow.use; 4702d61bbb3SSatish Balay 471433994e6SBarry Smith PetscFunctionBegin; 4723649974fSBarry Smith ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr); 47326e093fcSHong Zhang ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr); 4742d61bbb3SSatish Balay 4752d61bbb3SSatish Balay idx = a->j; 4762d61bbb3SSatish Balay v = a->a; 47726e093fcSHong Zhang if (usecprow) { 47826e093fcSHong Zhang mbs = a->compressedrow.nrows; 47926e093fcSHong Zhang ii = a->compressedrow.i; 4807b2bb3b9SHong Zhang ridx = a->compressedrow.rindex; 48126e093fcSHong Zhang } else { 48226e093fcSHong Zhang mbs = a->mbs; 4832d61bbb3SSatish Balay ii = a->i; 48426e093fcSHong Zhang z = zarray; 48526e093fcSHong Zhang } 4862d61bbb3SSatish Balay 4872d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 4882d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 4892d61bbb3SSatish Balay sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; 490444d8c10SJed Brown PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA); /* Indices for the next row (assumes same size as this one) */ 491444d8c10SJed Brown PetscPrefetchBlock(v+25*n,25*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */ 4922d61bbb3SSatish Balay for (j=0; j<n; j++) { 4932d61bbb3SSatish Balay xb = x + 5*(*idx++); 4942d61bbb3SSatish Balay x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; 4952d61bbb3SSatish Balay sum1 += v[0]*x1 + v[5]*x2 + v[10]*x3 + v[15]*x4 + v[20]*x5; 4962d61bbb3SSatish Balay sum2 += v[1]*x1 + v[6]*x2 + v[11]*x3 + v[16]*x4 + v[21]*x5; 4972d61bbb3SSatish Balay sum3 += v[2]*x1 + v[7]*x2 + v[12]*x3 + v[17]*x4 + v[22]*x5; 4982d61bbb3SSatish Balay sum4 += v[3]*x1 + v[8]*x2 + v[13]*x3 + v[18]*x4 + v[23]*x5; 4992d61bbb3SSatish Balay sum5 += v[4]*x1 + v[9]*x2 + v[14]*x3 + v[19]*x4 + v[24]*x5; 5002d61bbb3SSatish Balay v += 25; 5012d61bbb3SSatish Balay } 5027b2bb3b9SHong Zhang if (usecprow) z = zarray + 5*ridx[i]; 5032d61bbb3SSatish Balay z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; 50426e093fcSHong Zhang if (!usecprow) z += 5; 5052d61bbb3SSatish Balay } 5063649974fSBarry Smith ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr); 50726e093fcSHong Zhang ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr); 5087c565772SBarry Smith ierr = PetscLogFlops(50.0*a->nz - 5.0*a->nonzerorowcnt);CHKERRQ(ierr); 5092d61bbb3SSatish Balay PetscFunctionReturn(0); 5102d61bbb3SSatish Balay } 5112d61bbb3SSatish Balay 51215091d37SBarry Smith 513d6232bceSMichael Lange 5144a2ae208SSatish Balay #undef __FUNCT__ 5154a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_6" 516dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_6(Mat A,Vec xx,Vec zz) 51715091d37SBarry Smith { 51815091d37SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 519d9fead3dSBarry Smith PetscScalar *z = 0,sum1,sum2,sum3,sum4,sum5,sum6; 520d9fead3dSBarry Smith const PetscScalar *x,*xb; 52126e093fcSHong Zhang PetscScalar x1,x2,x3,x4,x5,x6,*zarray; 522d9fead3dSBarry Smith const MatScalar *v; 523dfbe8321SBarry Smith PetscErrorCode ierr; 5247c565772SBarry Smith PetscInt mbs,i,*idx,*ii,j,n,*ridx=NULL; 525ace3abfcSBarry Smith PetscBool usecprow=a->compressedrow.use; 52615091d37SBarry Smith 527433994e6SBarry Smith PetscFunctionBegin; 5283649974fSBarry Smith ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr); 52926e093fcSHong Zhang ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr); 53015091d37SBarry Smith 53115091d37SBarry Smith idx = a->j; 53215091d37SBarry Smith v = a->a; 53326e093fcSHong Zhang if (usecprow) { 53426e093fcSHong Zhang mbs = a->compressedrow.nrows; 53526e093fcSHong Zhang ii = a->compressedrow.i; 5367b2bb3b9SHong Zhang ridx = a->compressedrow.rindex; 53726e093fcSHong Zhang } else { 53826e093fcSHong Zhang mbs = a->mbs; 53915091d37SBarry Smith ii = a->i; 54026e093fcSHong Zhang z = zarray; 54126e093fcSHong Zhang } 54215091d37SBarry Smith 54315091d37SBarry Smith for (i=0; i<mbs; i++) { 54426fbe8dcSKarl Rupp n = ii[1] - ii[0]; 54526fbe8dcSKarl Rupp ii++; 54626fbe8dcSKarl Rupp sum1 = 0.0; 54726fbe8dcSKarl Rupp sum2 = 0.0; 54826fbe8dcSKarl Rupp sum3 = 0.0; 54926fbe8dcSKarl Rupp sum4 = 0.0; 55026fbe8dcSKarl Rupp sum5 = 0.0; 55126fbe8dcSKarl Rupp sum6 = 0.0; 55226fbe8dcSKarl Rupp 553444d8c10SJed Brown PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA); /* Indices for the next row (assumes same size as this one) */ 554444d8c10SJed Brown PetscPrefetchBlock(v+36*n,36*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */ 55515091d37SBarry Smith for (j=0; j<n; j++) { 55615091d37SBarry Smith xb = x + 6*(*idx++); 55715091d37SBarry Smith x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; 55815091d37SBarry Smith sum1 += v[0]*x1 + v[6]*x2 + v[12]*x3 + v[18]*x4 + v[24]*x5 + v[30]*x6; 55915091d37SBarry Smith sum2 += v[1]*x1 + v[7]*x2 + v[13]*x3 + v[19]*x4 + v[25]*x5 + v[31]*x6; 56015091d37SBarry Smith sum3 += v[2]*x1 + v[8]*x2 + v[14]*x3 + v[20]*x4 + v[26]*x5 + v[32]*x6; 56115091d37SBarry Smith sum4 += v[3]*x1 + v[9]*x2 + v[15]*x3 + v[21]*x4 + v[27]*x5 + v[33]*x6; 56215091d37SBarry Smith sum5 += v[4]*x1 + v[10]*x2 + v[16]*x3 + v[22]*x4 + v[28]*x5 + v[34]*x6; 56315091d37SBarry Smith sum6 += v[5]*x1 + v[11]*x2 + v[17]*x3 + v[23]*x4 + v[29]*x5 + v[35]*x6; 56415091d37SBarry Smith v += 36; 56515091d37SBarry Smith } 5667b2bb3b9SHong Zhang if (usecprow) z = zarray + 6*ridx[i]; 56715091d37SBarry Smith z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; 56826e093fcSHong Zhang if (!usecprow) z += 6; 56915091d37SBarry Smith } 57015091d37SBarry Smith 5713649974fSBarry Smith ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr); 57226e093fcSHong Zhang ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr); 5737c565772SBarry Smith ierr = PetscLogFlops(72.0*a->nz - 6.0*a->nonzerorowcnt);CHKERRQ(ierr); 57415091d37SBarry Smith PetscFunctionReturn(0); 57515091d37SBarry Smith } 5768ab949d8SShri Abhyankar 5774a2ae208SSatish Balay #undef __FUNCT__ 5784a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_7" 579dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_7(Mat A,Vec xx,Vec zz) 5802d61bbb3SSatish Balay { 5812d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 582d9fead3dSBarry Smith PetscScalar *z = 0,sum1,sum2,sum3,sum4,sum5,sum6,sum7; 583d9fead3dSBarry Smith const PetscScalar *x,*xb; 58426e093fcSHong Zhang PetscScalar x1,x2,x3,x4,x5,x6,x7,*zarray; 585d9fead3dSBarry Smith const MatScalar *v; 586dfbe8321SBarry Smith PetscErrorCode ierr; 5877c565772SBarry Smith PetscInt mbs,i,*idx,*ii,j,n,*ridx=NULL; 588ace3abfcSBarry Smith PetscBool usecprow=a->compressedrow.use; 5892d61bbb3SSatish Balay 590433994e6SBarry Smith PetscFunctionBegin; 5913649974fSBarry Smith ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr); 59226e093fcSHong Zhang ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr); 5932d61bbb3SSatish Balay 5942d61bbb3SSatish Balay idx = a->j; 5952d61bbb3SSatish Balay v = a->a; 59626e093fcSHong Zhang if (usecprow) { 59726e093fcSHong Zhang mbs = a->compressedrow.nrows; 59826e093fcSHong Zhang ii = a->compressedrow.i; 5997b2bb3b9SHong Zhang ridx = a->compressedrow.rindex; 60026e093fcSHong Zhang } else { 60126e093fcSHong Zhang mbs = a->mbs; 6022d61bbb3SSatish Balay ii = a->i; 60326e093fcSHong Zhang z = zarray; 60426e093fcSHong Zhang } 6052d61bbb3SSatish Balay 6062d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 60726fbe8dcSKarl Rupp n = ii[1] - ii[0]; 60826fbe8dcSKarl Rupp ii++; 60926fbe8dcSKarl Rupp sum1 = 0.0; 61026fbe8dcSKarl Rupp sum2 = 0.0; 61126fbe8dcSKarl Rupp sum3 = 0.0; 61226fbe8dcSKarl Rupp sum4 = 0.0; 61326fbe8dcSKarl Rupp sum5 = 0.0; 61426fbe8dcSKarl Rupp sum6 = 0.0; 61526fbe8dcSKarl Rupp sum7 = 0.0; 61626fbe8dcSKarl Rupp 617444d8c10SJed Brown PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA); /* Indices for the next row (assumes same size as this one) */ 618444d8c10SJed Brown PetscPrefetchBlock(v+49*n,49*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */ 6192d61bbb3SSatish Balay for (j=0; j<n; j++) { 6202d61bbb3SSatish Balay xb = x + 7*(*idx++); 6212d61bbb3SSatish Balay x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6]; 6222d61bbb3SSatish Balay sum1 += v[0]*x1 + v[7]*x2 + v[14]*x3 + v[21]*x4 + v[28]*x5 + v[35]*x6 + v[42]*x7; 6232d61bbb3SSatish Balay sum2 += v[1]*x1 + v[8]*x2 + v[15]*x3 + v[22]*x4 + v[29]*x5 + v[36]*x6 + v[43]*x7; 6242d61bbb3SSatish Balay sum3 += v[2]*x1 + v[9]*x2 + v[16]*x3 + v[23]*x4 + v[30]*x5 + v[37]*x6 + v[44]*x7; 6252d61bbb3SSatish Balay sum4 += v[3]*x1 + v[10]*x2 + v[17]*x3 + v[24]*x4 + v[31]*x5 + v[38]*x6 + v[45]*x7; 6262d61bbb3SSatish Balay sum5 += v[4]*x1 + v[11]*x2 + v[18]*x3 + v[25]*x4 + v[32]*x5 + v[39]*x6 + v[46]*x7; 6272d61bbb3SSatish Balay sum6 += v[5]*x1 + v[12]*x2 + v[19]*x3 + v[26]*x4 + v[33]*x5 + v[40]*x6 + v[47]*x7; 6282d61bbb3SSatish Balay sum7 += v[6]*x1 + v[13]*x2 + v[20]*x3 + v[27]*x4 + v[34]*x5 + v[41]*x6 + v[48]*x7; 6292d61bbb3SSatish Balay v += 49; 6302d61bbb3SSatish Balay } 6317b2bb3b9SHong Zhang if (usecprow) z = zarray + 7*ridx[i]; 6322d61bbb3SSatish Balay z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7; 63326e093fcSHong Zhang if (!usecprow) z += 7; 6342d61bbb3SSatish Balay } 6352d61bbb3SSatish Balay 6363649974fSBarry Smith ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr); 63726e093fcSHong Zhang ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr); 6387c565772SBarry Smith ierr = PetscLogFlops(98.0*a->nz - 7.0*a->nonzerorowcnt);CHKERRQ(ierr); 6392d61bbb3SSatish Balay PetscFunctionReturn(0); 6402d61bbb3SSatish Balay } 6412d61bbb3SSatish Balay 6428ab949d8SShri Abhyankar /* MatMult_SeqBAIJ_15 version 1: Columns in the block are accessed one at a time */ 643832cc040SShri Abhyankar /* Default MatMult for block size 15 */ 6448ab949d8SShri Abhyankar 6458ab949d8SShri Abhyankar #undef __FUNCT__ 6468ab949d8SShri Abhyankar #define __FUNCT__ "MatMult_SeqBAIJ_15_ver1" 6478ab949d8SShri Abhyankar PetscErrorCode MatMult_SeqBAIJ_15_ver1(Mat A,Vec xx,Vec zz) 6488ab949d8SShri Abhyankar { 6498ab949d8SShri Abhyankar Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 6508ab949d8SShri Abhyankar PetscScalar *z = 0,sum1,sum2,sum3,sum4,sum5,sum6,sum7,sum8,sum9,sum10,sum11,sum12,sum13,sum14,sum15; 6518ab949d8SShri Abhyankar const PetscScalar *x,*xb; 65253ef36baSBarry Smith PetscScalar *zarray,xv; 6538ab949d8SShri Abhyankar const MatScalar *v; 6548ab949d8SShri Abhyankar PetscErrorCode ierr; 6558ab949d8SShri Abhyankar const PetscInt *ii,*ij=a->j,*idx; 6567c565772SBarry Smith PetscInt mbs,i,j,k,n,*ridx=NULL; 657ace3abfcSBarry Smith PetscBool usecprow=a->compressedrow.use; 6588ab949d8SShri Abhyankar 6598ab949d8SShri Abhyankar PetscFunctionBegin; 6603649974fSBarry Smith ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr); 6618ab949d8SShri Abhyankar ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr); 6628ab949d8SShri Abhyankar 6638ab949d8SShri Abhyankar v = a->a; 6648ab949d8SShri Abhyankar if (usecprow) { 6658ab949d8SShri Abhyankar mbs = a->compressedrow.nrows; 6668ab949d8SShri Abhyankar ii = a->compressedrow.i; 6678ab949d8SShri Abhyankar ridx = a->compressedrow.rindex; 6688ab949d8SShri Abhyankar } else { 6698ab949d8SShri Abhyankar mbs = a->mbs; 6708ab949d8SShri Abhyankar ii = a->i; 6718ab949d8SShri Abhyankar z = zarray; 6728ab949d8SShri Abhyankar } 6738ab949d8SShri Abhyankar 6748ab949d8SShri Abhyankar for (i=0; i<mbs; i++) { 6758ab949d8SShri Abhyankar n = ii[i+1] - ii[i]; 6768ab949d8SShri Abhyankar idx = ij + ii[i]; 6778ab949d8SShri Abhyankar sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0; sum7 = 0.0; 6788ab949d8SShri 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; 6798ab949d8SShri Abhyankar 6808ab949d8SShri Abhyankar for (j=0; j<n; j++) { 6818ab949d8SShri Abhyankar xb = x + 15*(idx[j]); 6828ab949d8SShri Abhyankar 6838ab949d8SShri Abhyankar for (k=0; k<15; k++) { 68453ef36baSBarry Smith xv = xb[k]; 68553ef36baSBarry Smith sum1 += v[0]*xv; 68653ef36baSBarry Smith sum2 += v[1]*xv; 68753ef36baSBarry Smith sum3 += v[2]*xv; 68853ef36baSBarry Smith sum4 += v[3]*xv; 68953ef36baSBarry Smith sum5 += v[4]*xv; 69053ef36baSBarry Smith sum6 += v[5]*xv; 69153ef36baSBarry Smith sum7 += v[6]*xv; 69253ef36baSBarry Smith sum8 += v[7]*xv; 69353ef36baSBarry Smith sum9 += v[8]*xv; 69453ef36baSBarry Smith sum10 += v[9]*xv; 69553ef36baSBarry Smith sum11 += v[10]*xv; 69653ef36baSBarry Smith sum12 += v[11]*xv; 69753ef36baSBarry Smith sum13 += v[12]*xv; 69853ef36baSBarry Smith sum14 += v[13]*xv; 69953ef36baSBarry Smith sum15 += v[14]*xv; 7008ab949d8SShri Abhyankar v += 15; 7018ab949d8SShri Abhyankar } 7028ab949d8SShri Abhyankar } 7038ab949d8SShri Abhyankar if (usecprow) z = zarray + 15*ridx[i]; 7048ab949d8SShri Abhyankar z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7; 7058ab949d8SShri 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; 7068ab949d8SShri Abhyankar 7078ab949d8SShri Abhyankar if (!usecprow) z += 15; 7088ab949d8SShri Abhyankar } 7098ab949d8SShri Abhyankar 7103649974fSBarry Smith ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr); 7118ab949d8SShri Abhyankar ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr); 7127c565772SBarry Smith ierr = PetscLogFlops(450.0*a->nz - 15.0*a->nonzerorowcnt);CHKERRQ(ierr); 7138ab949d8SShri Abhyankar PetscFunctionReturn(0); 7148ab949d8SShri Abhyankar } 7158ab949d8SShri Abhyankar 7168ab949d8SShri Abhyankar /* MatMult_SeqBAIJ_15_ver2 : Columns in the block are accessed in sets of 4,4,4,3 */ 7178ab949d8SShri Abhyankar #undef __FUNCT__ 7188ab949d8SShri Abhyankar #define __FUNCT__ "MatMult_SeqBAIJ_15_ver2" 7198ab949d8SShri Abhyankar PetscErrorCode MatMult_SeqBAIJ_15_ver2(Mat A,Vec xx,Vec zz) 7208ab949d8SShri Abhyankar { 7218ab949d8SShri Abhyankar Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 7228ab949d8SShri Abhyankar PetscScalar *z = 0,sum1,sum2,sum3,sum4,sum5,sum6,sum7,sum8,sum9,sum10,sum11,sum12,sum13,sum14,sum15; 7238ab949d8SShri Abhyankar const PetscScalar *x,*xb; 7240b8f6341SShri Abhyankar PetscScalar x1,x2,x3,x4,*zarray; 7258ab949d8SShri Abhyankar const MatScalar *v; 7268ab949d8SShri Abhyankar PetscErrorCode ierr; 7278ab949d8SShri Abhyankar const PetscInt *ii,*ij=a->j,*idx; 7287c565772SBarry Smith PetscInt mbs,i,j,n,*ridx=NULL; 729ace3abfcSBarry Smith PetscBool usecprow=a->compressedrow.use; 7308ab949d8SShri Abhyankar 7318ab949d8SShri Abhyankar PetscFunctionBegin; 7323649974fSBarry Smith ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr); 7338ab949d8SShri Abhyankar ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr); 7348ab949d8SShri Abhyankar 7358ab949d8SShri Abhyankar v = a->a; 7368ab949d8SShri Abhyankar if (usecprow) { 7378ab949d8SShri Abhyankar mbs = a->compressedrow.nrows; 7388ab949d8SShri Abhyankar ii = a->compressedrow.i; 7398ab949d8SShri Abhyankar ridx = a->compressedrow.rindex; 7408ab949d8SShri Abhyankar } else { 7418ab949d8SShri Abhyankar mbs = a->mbs; 7428ab949d8SShri Abhyankar ii = a->i; 7438ab949d8SShri Abhyankar z = zarray; 7448ab949d8SShri Abhyankar } 7458ab949d8SShri Abhyankar 7468ab949d8SShri Abhyankar for (i=0; i<mbs; i++) { 7478ab949d8SShri Abhyankar n = ii[i+1] - ii[i]; 7488ab949d8SShri Abhyankar idx = ij + ii[i]; 7498ab949d8SShri Abhyankar sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0; sum7 = 0.0; 7508ab949d8SShri 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; 7518ab949d8SShri Abhyankar 7528ab949d8SShri Abhyankar for (j=0; j<n; j++) { 7538ab949d8SShri Abhyankar xb = x + 15*(idx[j]); 7540b8f6341SShri Abhyankar x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; 7558ab949d8SShri Abhyankar 7568ab949d8SShri Abhyankar sum1 += v[0]*x1 + v[15]*x2 + v[30]*x3 + v[45]*x4; 7578ab949d8SShri Abhyankar sum2 += v[1]*x1 + v[16]*x2 + v[31]*x3 + v[46]*x4; 7588ab949d8SShri Abhyankar sum3 += v[2]*x1 + v[17]*x2 + v[32]*x3 + v[47]*x4; 7598ab949d8SShri Abhyankar sum4 += v[3]*x1 + v[18]*x2 + v[33]*x3 + v[48]*x4; 7608ab949d8SShri Abhyankar sum5 += v[4]*x1 + v[19]*x2 + v[34]*x3 + v[49]*x4; 7618ab949d8SShri Abhyankar sum6 += v[5]*x1 + v[20]*x2 + v[35]*x3 + v[50]*x4; 7628ab949d8SShri Abhyankar sum7 += v[6]*x1 + v[21]*x2 + v[36]*x3 + v[51]*x4; 7638ab949d8SShri Abhyankar sum8 += v[7]*x1 + v[22]*x2 + v[37]*x3 + v[52]*x4; 7648ab949d8SShri Abhyankar sum9 += v[8]*x1 + v[23]*x2 + v[38]*x3 + v[53]*x4; 7658ab949d8SShri Abhyankar sum10 += v[9]*x1 + v[24]*x2 + v[39]*x3 + v[54]*x4; 7668ab949d8SShri Abhyankar sum11 += v[10]*x1 + v[25]*x2 + v[40]*x3 + v[55]*x4; 7678ab949d8SShri Abhyankar sum12 += v[11]*x1 + v[26]*x2 + v[41]*x3 + v[56]*x4; 7688ab949d8SShri Abhyankar sum13 += v[12]*x1 + v[27]*x2 + v[42]*x3 + v[57]*x4; 7698ab949d8SShri Abhyankar sum14 += v[13]*x1 + v[28]*x2 + v[43]*x3 + v[58]*x4; 7708ab949d8SShri Abhyankar sum15 += v[14]*x1 + v[29]*x2 + v[44]*x3 + v[59]*x4; 7718ab949d8SShri Abhyankar 7728ab949d8SShri Abhyankar v += 60; 7738ab949d8SShri Abhyankar 7740b8f6341SShri Abhyankar x1 = xb[4]; x2 = xb[5]; x3 = xb[6]; x4 = xb[7]; 7758ab949d8SShri Abhyankar 7768ab949d8SShri Abhyankar sum1 += v[0]*x1 + v[15]*x2 + v[30]*x3 + v[45]*x4; 7778ab949d8SShri Abhyankar sum2 += v[1]*x1 + v[16]*x2 + v[31]*x3 + v[46]*x4; 7788ab949d8SShri Abhyankar sum3 += v[2]*x1 + v[17]*x2 + v[32]*x3 + v[47]*x4; 7798ab949d8SShri Abhyankar sum4 += v[3]*x1 + v[18]*x2 + v[33]*x3 + v[48]*x4; 7808ab949d8SShri Abhyankar sum5 += v[4]*x1 + v[19]*x2 + v[34]*x3 + v[49]*x4; 7818ab949d8SShri Abhyankar sum6 += v[5]*x1 + v[20]*x2 + v[35]*x3 + v[50]*x4; 7828ab949d8SShri Abhyankar sum7 += v[6]*x1 + v[21]*x2 + v[36]*x3 + v[51]*x4; 7838ab949d8SShri Abhyankar sum8 += v[7]*x1 + v[22]*x2 + v[37]*x3 + v[52]*x4; 7848ab949d8SShri Abhyankar sum9 += v[8]*x1 + v[23]*x2 + v[38]*x3 + v[53]*x4; 7858ab949d8SShri Abhyankar sum10 += v[9]*x1 + v[24]*x2 + v[39]*x3 + v[54]*x4; 7868ab949d8SShri Abhyankar sum11 += v[10]*x1 + v[25]*x2 + v[40]*x3 + v[55]*x4; 7878ab949d8SShri Abhyankar sum12 += v[11]*x1 + v[26]*x2 + v[41]*x3 + v[56]*x4; 7888ab949d8SShri Abhyankar sum13 += v[12]*x1 + v[27]*x2 + v[42]*x3 + v[57]*x4; 7898ab949d8SShri Abhyankar sum14 += v[13]*x1 + v[28]*x2 + v[43]*x3 + v[58]*x4; 7908ab949d8SShri Abhyankar sum15 += v[14]*x1 + v[29]*x2 + v[44]*x3 + v[59]*x4; 7918ab949d8SShri Abhyankar v += 60; 7928ab949d8SShri Abhyankar 7930b8f6341SShri Abhyankar x1 = xb[8]; x2 = xb[9]; x3 = xb[10]; x4 = xb[11]; 7940b8f6341SShri Abhyankar sum1 += v[0]*x1 + v[15]*x2 + v[30]*x3 + v[45]*x4; 7950b8f6341SShri Abhyankar sum2 += v[1]*x1 + v[16]*x2 + v[31]*x3 + v[46]*x4; 7960b8f6341SShri Abhyankar sum3 += v[2]*x1 + v[17]*x2 + v[32]*x3 + v[47]*x4; 7970b8f6341SShri Abhyankar sum4 += v[3]*x1 + v[18]*x2 + v[33]*x3 + v[48]*x4; 7980b8f6341SShri Abhyankar sum5 += v[4]*x1 + v[19]*x2 + v[34]*x3 + v[49]*x4; 7990b8f6341SShri Abhyankar sum6 += v[5]*x1 + v[20]*x2 + v[35]*x3 + v[50]*x4; 8000b8f6341SShri Abhyankar sum7 += v[6]*x1 + v[21]*x2 + v[36]*x3 + v[51]*x4; 8010b8f6341SShri Abhyankar sum8 += v[7]*x1 + v[22]*x2 + v[37]*x3 + v[52]*x4; 8020b8f6341SShri Abhyankar sum9 += v[8]*x1 + v[23]*x2 + v[38]*x3 + v[53]*x4; 8030b8f6341SShri Abhyankar sum10 += v[9]*x1 + v[24]*x2 + v[39]*x3 + v[54]*x4; 8040b8f6341SShri Abhyankar sum11 += v[10]*x1 + v[25]*x2 + v[40]*x3 + v[55]*x4; 8050b8f6341SShri Abhyankar sum12 += v[11]*x1 + v[26]*x2 + v[41]*x3 + v[56]*x4; 8060b8f6341SShri Abhyankar sum13 += v[12]*x1 + v[27]*x2 + v[42]*x3 + v[57]*x4; 8070b8f6341SShri Abhyankar sum14 += v[13]*x1 + v[28]*x2 + v[43]*x3 + v[58]*x4; 8080b8f6341SShri Abhyankar sum15 += v[14]*x1 + v[29]*x2 + v[44]*x3 + v[59]*x4; 8090b8f6341SShri Abhyankar v += 60; 8100b8f6341SShri Abhyankar 8110b8f6341SShri Abhyankar x1 = xb[12]; x2 = xb[13]; x3 = xb[14]; 8128ab949d8SShri Abhyankar sum1 += v[0]*x1 + v[15]*x2 + v[30]*x3; 8138ab949d8SShri Abhyankar sum2 += v[1]*x1 + v[16]*x2 + v[31]*x3; 8148ab949d8SShri Abhyankar sum3 += v[2]*x1 + v[17]*x2 + v[32]*x3; 8158ab949d8SShri Abhyankar sum4 += v[3]*x1 + v[18]*x2 + v[33]*x3; 8168ab949d8SShri Abhyankar sum5 += v[4]*x1 + v[19]*x2 + v[34]*x3; 8178ab949d8SShri Abhyankar sum6 += v[5]*x1 + v[20]*x2 + v[35]*x3; 8188ab949d8SShri Abhyankar sum7 += v[6]*x1 + v[21]*x2 + v[36]*x3; 8198ab949d8SShri Abhyankar sum8 += v[7]*x1 + v[22]*x2 + v[37]*x3; 8208ab949d8SShri Abhyankar sum9 += v[8]*x1 + v[23]*x2 + v[38]*x3; 8218ab949d8SShri Abhyankar sum10 += v[9]*x1 + v[24]*x2 + v[39]*x3; 8228ab949d8SShri Abhyankar sum11 += v[10]*x1 + v[25]*x2 + v[40]*x3; 8238ab949d8SShri Abhyankar sum12 += v[11]*x1 + v[26]*x2 + v[41]*x3; 8248ab949d8SShri Abhyankar sum13 += v[12]*x1 + v[27]*x2 + v[42]*x3; 8258ab949d8SShri Abhyankar sum14 += v[13]*x1 + v[28]*x2 + v[43]*x3; 8268ab949d8SShri Abhyankar sum15 += v[14]*x1 + v[29]*x2 + v[44]*x3; 8278ab949d8SShri Abhyankar v += 45; 8288ab949d8SShri Abhyankar } 8298ab949d8SShri Abhyankar if (usecprow) z = zarray + 15*ridx[i]; 8308ab949d8SShri Abhyankar z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7; 8318ab949d8SShri 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; 8328ab949d8SShri Abhyankar 8338ab949d8SShri Abhyankar if (!usecprow) z += 15; 8348ab949d8SShri Abhyankar } 8358ab949d8SShri Abhyankar 8363649974fSBarry Smith ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr); 8378ab949d8SShri Abhyankar ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr); 8387c565772SBarry Smith ierr = PetscLogFlops(450.0*a->nz - 15.0*a->nonzerorowcnt);CHKERRQ(ierr); 8398ab949d8SShri Abhyankar PetscFunctionReturn(0); 8408ab949d8SShri Abhyankar } 8418ab949d8SShri Abhyankar 8428ab949d8SShri Abhyankar /* MatMult_SeqBAIJ_15_ver3 : Columns in the block are accessed in sets of 8,7 */ 8438ab949d8SShri Abhyankar #undef __FUNCT__ 8448ab949d8SShri Abhyankar #define __FUNCT__ "MatMult_SeqBAIJ_15_ver3" 8458ab949d8SShri Abhyankar PetscErrorCode MatMult_SeqBAIJ_15_ver3(Mat A,Vec xx,Vec zz) 8468ab949d8SShri Abhyankar { 8478ab949d8SShri Abhyankar Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 8488ab949d8SShri Abhyankar PetscScalar *z = 0,sum1,sum2,sum3,sum4,sum5,sum6,sum7,sum8,sum9,sum10,sum11,sum12,sum13,sum14,sum15; 8498ab949d8SShri Abhyankar const PetscScalar *x,*xb; 8500b8f6341SShri Abhyankar PetscScalar x1,x2,x3,x4,x5,x6,x7,x8,*zarray; 8518ab949d8SShri Abhyankar const MatScalar *v; 8528ab949d8SShri Abhyankar PetscErrorCode ierr; 8538ab949d8SShri Abhyankar const PetscInt *ii,*ij=a->j,*idx; 8547c565772SBarry Smith PetscInt mbs,i,j,n,*ridx=NULL; 855ace3abfcSBarry Smith PetscBool usecprow=a->compressedrow.use; 8568ab949d8SShri Abhyankar 8578ab949d8SShri Abhyankar PetscFunctionBegin; 8583649974fSBarry Smith ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr); 8598ab949d8SShri Abhyankar ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr); 8608ab949d8SShri Abhyankar 8618ab949d8SShri Abhyankar v = a->a; 8628ab949d8SShri Abhyankar if (usecprow) { 8638ab949d8SShri Abhyankar mbs = a->compressedrow.nrows; 8648ab949d8SShri Abhyankar ii = a->compressedrow.i; 8658ab949d8SShri Abhyankar ridx = a->compressedrow.rindex; 8668ab949d8SShri Abhyankar } else { 8678ab949d8SShri Abhyankar mbs = a->mbs; 8688ab949d8SShri Abhyankar ii = a->i; 8698ab949d8SShri Abhyankar z = zarray; 8708ab949d8SShri Abhyankar } 8718ab949d8SShri Abhyankar 8728ab949d8SShri Abhyankar for (i=0; i<mbs; i++) { 8738ab949d8SShri Abhyankar n = ii[i+1] - ii[i]; 8748ab949d8SShri Abhyankar idx = ij + ii[i]; 8758ab949d8SShri Abhyankar sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0; sum7 = 0.0; 8768ab949d8SShri 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; 8778ab949d8SShri Abhyankar 8788ab949d8SShri Abhyankar for (j=0; j<n; j++) { 8798ab949d8SShri Abhyankar xb = x + 15*(idx[j]); 8808ab949d8SShri Abhyankar x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6]; 8810b8f6341SShri Abhyankar x8 = xb[7]; 8828ab949d8SShri Abhyankar 8838ab949d8SShri 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; 8848ab949d8SShri 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; 8858ab949d8SShri 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; 8868ab949d8SShri 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; 8878ab949d8SShri 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; 8888ab949d8SShri 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; 8898ab949d8SShri 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; 8908ab949d8SShri 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; 8918ab949d8SShri 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; 8928ab949d8SShri 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; 8938ab949d8SShri 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; 8948ab949d8SShri 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; 8958ab949d8SShri 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; 8968ab949d8SShri 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; 8978ab949d8SShri 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; 8988ab949d8SShri Abhyankar v += 120; 8998ab949d8SShri Abhyankar 9000b8f6341SShri Abhyankar x1 = xb[8]; x2 = xb[9]; x3 = xb[10]; x4 = xb[11]; x5 = xb[12]; x6 = xb[13]; x7 = xb[14]; 9010b8f6341SShri Abhyankar 9028ab949d8SShri Abhyankar sum1 += v[0]*x1 + v[15]*x2 + v[30]*x3 + v[45]*x4 + v[60]*x5 + v[75]*x6 + v[90]*x7; 9038ab949d8SShri Abhyankar sum2 += v[1]*x1 + v[16]*x2 + v[31]*x3 + v[46]*x4 + v[61]*x5 + v[76]*x6 + v[91]*x7; 9048ab949d8SShri Abhyankar sum3 += v[2]*x1 + v[17]*x2 + v[32]*x3 + v[47]*x4 + v[62]*x5 + v[77]*x6 + v[92]*x7; 9058ab949d8SShri Abhyankar sum4 += v[3]*x1 + v[18]*x2 + v[33]*x3 + v[48]*x4 + v[63]*x5 + v[78]*x6 + v[93]*x7; 9068ab949d8SShri Abhyankar sum5 += v[4]*x1 + v[19]*x2 + v[34]*x3 + v[49]*x4 + v[64]*x5 + v[79]*x6 + v[94]*x7; 9078ab949d8SShri Abhyankar sum6 += v[5]*x1 + v[20]*x2 + v[35]*x3 + v[50]*x4 + v[65]*x5 + v[80]*x6 + v[95]*x7; 9088ab949d8SShri Abhyankar sum7 += v[6]*x1 + v[21]*x2 + v[36]*x3 + v[51]*x4 + v[66]*x5 + v[81]*x6 + v[96]*x7; 9098ab949d8SShri Abhyankar sum8 += v[7]*x1 + v[22]*x2 + v[37]*x3 + v[52]*x4 + v[67]*x5 + v[82]*x6 + v[97]*x7; 9108ab949d8SShri Abhyankar sum9 += v[8]*x1 + v[23]*x2 + v[38]*x3 + v[53]*x4 + v[68]*x5 + v[83]*x6 + v[98]*x7; 9118ab949d8SShri Abhyankar sum10 += v[9]*x1 + v[24]*x2 + v[39]*x3 + v[54]*x4 + v[69]*x5 + v[84]*x6 + v[99]*x7; 9128ab949d8SShri Abhyankar sum11 += v[10]*x1 + v[25]*x2 + v[40]*x3 + v[55]*x4 + v[70]*x5 + v[85]*x6 + v[100]*x7; 9138ab949d8SShri Abhyankar sum12 += v[11]*x1 + v[26]*x2 + v[41]*x3 + v[56]*x4 + v[71]*x5 + v[86]*x6 + v[101]*x7; 9148ab949d8SShri Abhyankar sum13 += v[12]*x1 + v[27]*x2 + v[42]*x3 + v[57]*x4 + v[72]*x5 + v[87]*x6 + v[102]*x7; 9158ab949d8SShri Abhyankar sum14 += v[13]*x1 + v[28]*x2 + v[43]*x3 + v[58]*x4 + v[73]*x5 + v[88]*x6 + v[103]*x7; 9168ab949d8SShri Abhyankar sum15 += v[14]*x1 + v[29]*x2 + v[44]*x3 + v[59]*x4 + v[74]*x5 + v[89]*x6 + v[104]*x7; 9178ab949d8SShri Abhyankar v += 105; 9188ab949d8SShri Abhyankar } 9198ab949d8SShri Abhyankar if (usecprow) z = zarray + 15*ridx[i]; 9208ab949d8SShri Abhyankar z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7; 9218ab949d8SShri 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; 9228ab949d8SShri Abhyankar 9238ab949d8SShri Abhyankar if (!usecprow) z += 15; 9248ab949d8SShri Abhyankar } 9258ab949d8SShri Abhyankar 9263649974fSBarry Smith ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr); 9278ab949d8SShri Abhyankar ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr); 9287c565772SBarry Smith ierr = PetscLogFlops(450.0*a->nz - 15.0*a->nonzerorowcnt);CHKERRQ(ierr); 9298ab949d8SShri Abhyankar PetscFunctionReturn(0); 9308ab949d8SShri Abhyankar } 9318ab949d8SShri Abhyankar 9328ab949d8SShri Abhyankar /* MatMult_SeqBAIJ_15_ver4 : All columns in the block are accessed at once */ 9338ab949d8SShri Abhyankar 9348ab949d8SShri Abhyankar #undef __FUNCT__ 9358ab949d8SShri Abhyankar #define __FUNCT__ "MatMult_SeqBAIJ_15_ver4" 9368ab949d8SShri Abhyankar PetscErrorCode MatMult_SeqBAIJ_15_ver4(Mat A,Vec xx,Vec zz) 9378ab949d8SShri Abhyankar { 9388ab949d8SShri Abhyankar Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 9398ab949d8SShri Abhyankar PetscScalar *z = 0,sum1,sum2,sum3,sum4,sum5,sum6,sum7,sum8,sum9,sum10,sum11,sum12,sum13,sum14,sum15; 9408ab949d8SShri Abhyankar const PetscScalar *x,*xb; 9418ab949d8SShri Abhyankar PetscScalar x1,x2,x3,x4,x5,x6,x7,x8,x9,x10,x11,x12,x13,x14,x15,*zarray; 9428ab949d8SShri Abhyankar const MatScalar *v; 9438ab949d8SShri Abhyankar PetscErrorCode ierr; 9448ab949d8SShri Abhyankar const PetscInt *ii,*ij=a->j,*idx; 9457c565772SBarry Smith PetscInt mbs,i,j,n,*ridx=NULL; 946ace3abfcSBarry Smith PetscBool usecprow=a->compressedrow.use; 9478ab949d8SShri Abhyankar 9488ab949d8SShri Abhyankar PetscFunctionBegin; 9493649974fSBarry Smith ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr); 9508ab949d8SShri Abhyankar ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr); 9518ab949d8SShri Abhyankar 9528ab949d8SShri Abhyankar v = a->a; 9538ab949d8SShri Abhyankar if (usecprow) { 9548ab949d8SShri Abhyankar mbs = a->compressedrow.nrows; 9558ab949d8SShri Abhyankar ii = a->compressedrow.i; 9568ab949d8SShri Abhyankar ridx = a->compressedrow.rindex; 9578ab949d8SShri Abhyankar } else { 9588ab949d8SShri Abhyankar mbs = a->mbs; 9598ab949d8SShri Abhyankar ii = a->i; 9608ab949d8SShri Abhyankar z = zarray; 9618ab949d8SShri Abhyankar } 9628ab949d8SShri Abhyankar 9638ab949d8SShri Abhyankar for (i=0; i<mbs; i++) { 9648ab949d8SShri Abhyankar n = ii[i+1] - ii[i]; 9658ab949d8SShri Abhyankar idx = ij + ii[i]; 9668ab949d8SShri Abhyankar sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0; sum7 = 0.0; 9678ab949d8SShri 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; 9688ab949d8SShri Abhyankar 9698ab949d8SShri Abhyankar for (j=0; j<n; j++) { 9708ab949d8SShri Abhyankar xb = x + 15*(idx[j]); 9718ab949d8SShri Abhyankar x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6]; 9728ab949d8SShri 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]; 9738ab949d8SShri Abhyankar 9748ab949d8SShri 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; 9758ab949d8SShri 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; 9768ab949d8SShri 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; 9778ab949d8SShri 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; 9788ab949d8SShri 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; 9798ab949d8SShri 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; 9808ab949d8SShri 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; 9818ab949d8SShri 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; 9828ab949d8SShri 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; 9838ab949d8SShri 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; 9848ab949d8SShri 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; 9858ab949d8SShri 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; 9868ab949d8SShri 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; 9878ab949d8SShri 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; 9888ab949d8SShri 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; 9898ab949d8SShri Abhyankar v += 225; 9908ab949d8SShri Abhyankar } 9918ab949d8SShri Abhyankar if (usecprow) z = zarray + 15*ridx[i]; 9928ab949d8SShri Abhyankar z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7; 9938ab949d8SShri 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; 9948ab949d8SShri Abhyankar 9958ab949d8SShri Abhyankar if (!usecprow) z += 15; 9968ab949d8SShri Abhyankar } 9978ab949d8SShri Abhyankar 9983649974fSBarry Smith ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr); 9998ab949d8SShri Abhyankar ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr); 10007c565772SBarry Smith ierr = PetscLogFlops(450.0*a->nz - 15.0*a->nonzerorowcnt);CHKERRQ(ierr); 10018ab949d8SShri Abhyankar PetscFunctionReturn(0); 10028ab949d8SShri Abhyankar } 10038ab949d8SShri Abhyankar 10048ab949d8SShri Abhyankar 10053f1db9ecSBarry Smith /* 10063f1db9ecSBarry Smith This will not work with MatScalar == float because it calls the BLAS 10073f1db9ecSBarry Smith */ 10084a2ae208SSatish Balay #undef __FUNCT__ 10094a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_N" 1010dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_N(Mat A,Vec xx,Vec zz) 10112d61bbb3SSatish Balay { 10122d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1013d9ca1df4SBarry Smith PetscScalar *z = 0,*work,*workt,*zarray; 1014d9ca1df4SBarry Smith const PetscScalar *x,*xb; 1015d9ca1df4SBarry Smith const MatScalar *v; 1016dfbe8321SBarry Smith PetscErrorCode ierr; 1017d9ca1df4SBarry Smith PetscInt mbs,i,bs=A->rmap->bs,j,n,bs2=a->bs2; 1018d9ca1df4SBarry Smith const PetscInt *idx,*ii,*ridx=NULL; 1019d9ca1df4SBarry Smith PetscInt ncols,k; 1020ace3abfcSBarry Smith PetscBool usecprow=a->compressedrow.use; 10212d61bbb3SSatish Balay 10222d61bbb3SSatish Balay PetscFunctionBegin; 1023d9ca1df4SBarry Smith ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr); 102426e093fcSHong Zhang ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr); 10252d61bbb3SSatish Balay 10262d61bbb3SSatish Balay idx = a->j; 10272d61bbb3SSatish Balay v = a->a; 102826e093fcSHong Zhang if (usecprow) { 102926e093fcSHong Zhang mbs = a->compressedrow.nrows; 103026e093fcSHong Zhang ii = a->compressedrow.i; 10317b2bb3b9SHong Zhang ridx = a->compressedrow.rindex; 103226e093fcSHong Zhang } else { 103326e093fcSHong Zhang mbs = a->mbs; 10342d61bbb3SSatish Balay ii = a->i; 103526e093fcSHong Zhang z = zarray; 103626e093fcSHong Zhang } 1037218c64b6SSatish Balay 10382d61bbb3SSatish Balay if (!a->mult_work) { 1039d0f46423SBarry Smith k = PetscMax(A->rmap->n,A->cmap->n); 1040854ce69bSBarry Smith ierr = PetscMalloc1(k+1,&a->mult_work);CHKERRQ(ierr); 10412d61bbb3SSatish Balay } 10422d61bbb3SSatish Balay work = a->mult_work; 10432d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 10442d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 10452d61bbb3SSatish Balay ncols = n*bs; 10462d61bbb3SSatish Balay workt = work; 10472d61bbb3SSatish Balay for (j=0; j<n; j++) { 10482d61bbb3SSatish Balay xb = x + bs*(*idx++); 10492d61bbb3SSatish Balay for (k=0; k<bs; k++) workt[k] = xb[k]; 10502d61bbb3SSatish Balay workt += bs; 10512d61bbb3SSatish Balay } 10527b2bb3b9SHong Zhang if (usecprow) z = zarray + bs*ridx[i]; 105396b95a6bSBarry Smith PetscKernel_w_gets_Ar_times_v(bs,ncols,work,v,z); 105471044d3cSBarry Smith /* BLASgemv_("N",&bs,&ncols,&_DOne,v,&bs,work,&_One,&_DZero,z,&_One); */ 10552d61bbb3SSatish Balay v += n*bs2; 105626e093fcSHong Zhang if (!usecprow) z += bs; 10572d61bbb3SSatish Balay } 1058d9ca1df4SBarry Smith ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr); 105926e093fcSHong Zhang ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr); 10607c565772SBarry Smith ierr = PetscLogFlops(2.0*a->nz*bs2 - bs*a->nonzerorowcnt);CHKERRQ(ierr); 10612d61bbb3SSatish Balay PetscFunctionReturn(0); 10622d61bbb3SSatish Balay } 10632d61bbb3SSatish Balay 10644a2ae208SSatish Balay #undef __FUNCT__ 10654a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_1" 1066dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_1(Mat A,Vec xx,Vec yy,Vec zz) 10672d61bbb3SSatish Balay { 10682d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1069122f12eaSBarry Smith const PetscScalar *x; 1070122f12eaSBarry Smith PetscScalar *y,*z,sum; 1071122f12eaSBarry Smith const MatScalar *v; 1072dfbe8321SBarry Smith PetscErrorCode ierr; 10737c565772SBarry Smith PetscInt mbs=a->mbs,i,n,*ridx=NULL; 1074122f12eaSBarry Smith const PetscInt *idx,*ii; 1075ace3abfcSBarry Smith PetscBool usecprow=a->compressedrow.use; 10762d61bbb3SSatish Balay 10772d61bbb3SSatish Balay PetscFunctionBegin; 10783649974fSBarry Smith ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr); 1079d9ca1df4SBarry Smith ierr = VecGetArrayPair(yy,zz,&y,&z);CHKERRQ(ierr); 10802d61bbb3SSatish Balay 10812d61bbb3SSatish Balay idx = a->j; 10822d61bbb3SSatish Balay v = a->a; 108326e093fcSHong Zhang if (usecprow) { 10844eb6d288SHong Zhang if (zz != yy) { 1085122f12eaSBarry Smith ierr = PetscMemcpy(z,y,mbs*sizeof(PetscScalar));CHKERRQ(ierr); 10864eb6d288SHong Zhang } 108726e093fcSHong Zhang mbs = a->compressedrow.nrows; 108826e093fcSHong Zhang ii = a->compressedrow.i; 10897b2bb3b9SHong Zhang ridx = a->compressedrow.rindex; 109026e093fcSHong Zhang } else { 10912d61bbb3SSatish Balay ii = a->i; 109226e093fcSHong Zhang } 10932d61bbb3SSatish Balay 10942d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 1095122f12eaSBarry Smith n = ii[1] - ii[0]; 1096122f12eaSBarry Smith ii++; 109726e093fcSHong Zhang if (!usecprow) { 1098122f12eaSBarry Smith sum = y[i]; 1099122f12eaSBarry Smith } else { 1100122f12eaSBarry Smith sum = y[ridx[i]]; 1101122f12eaSBarry Smith } 1102444d8c10SJed Brown PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA); /* Indices for the next row (assumes same size as this one) */ 1103444d8c10SJed Brown PetscPrefetchBlock(v+n,n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */ 1104122f12eaSBarry Smith PetscSparseDensePlusDot(sum,x,v,idx,n); 1105122f12eaSBarry Smith v += n; 1106122f12eaSBarry Smith idx += n; 1107122f12eaSBarry Smith if (usecprow) { 1108122f12eaSBarry Smith z[ridx[i]] = sum; 1109122f12eaSBarry Smith } else { 1110122f12eaSBarry Smith z[i] = sum; 111126e093fcSHong Zhang } 11122d61bbb3SSatish Balay } 11133649974fSBarry Smith ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr); 1114d9ca1df4SBarry Smith ierr = VecRestoreArrayPair(yy,zz,&y,&z);CHKERRQ(ierr); 11157c565772SBarry Smith ierr = PetscLogFlops(2.0*a->nz);CHKERRQ(ierr); 11162d61bbb3SSatish Balay PetscFunctionReturn(0); 11172d61bbb3SSatish Balay } 11182d61bbb3SSatish Balay 11194a2ae208SSatish Balay #undef __FUNCT__ 11204a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_2" 1121dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_2(Mat A,Vec xx,Vec yy,Vec zz) 11222d61bbb3SSatish Balay { 11232d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1124d9ca1df4SBarry Smith PetscScalar *y = 0,*z = 0,sum1,sum2; 1125d9ca1df4SBarry Smith const PetscScalar *x,*xb; 112626e093fcSHong Zhang PetscScalar x1,x2,*yarray,*zarray; 1127d9ca1df4SBarry Smith const MatScalar *v; 1128dfbe8321SBarry Smith PetscErrorCode ierr; 1129d9ca1df4SBarry Smith PetscInt mbs = a->mbs,i,n,j; 1130d9ca1df4SBarry Smith const PetscInt *idx,*ii,*ridx = NULL; 1131ace3abfcSBarry Smith PetscBool usecprow = a->compressedrow.use; 11322d61bbb3SSatish Balay 11332d61bbb3SSatish Balay PetscFunctionBegin; 1134d9ca1df4SBarry Smith ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr); 1135d9ca1df4SBarry Smith ierr = VecGetArrayPair(yy,zz,&yarray,&zarray);CHKERRQ(ierr); 11362d61bbb3SSatish Balay 11372d61bbb3SSatish Balay idx = a->j; 11382d61bbb3SSatish Balay v = a->a; 113926e093fcSHong Zhang if (usecprow) { 11404eb6d288SHong Zhang if (zz != yy) { 11414eb6d288SHong Zhang ierr = PetscMemcpy(zarray,yarray,2*mbs*sizeof(PetscScalar));CHKERRQ(ierr); 11424eb6d288SHong Zhang } 114326e093fcSHong Zhang mbs = a->compressedrow.nrows; 114426e093fcSHong Zhang ii = a->compressedrow.i; 11457b2bb3b9SHong Zhang ridx = a->compressedrow.rindex; 11464eb6d288SHong Zhang if (zz != yy) { 11474eb6d288SHong Zhang ierr = PetscMemcpy(zarray,yarray,a->mbs*sizeof(PetscScalar));CHKERRQ(ierr); 11484eb6d288SHong Zhang } 114926e093fcSHong Zhang } else { 11502d61bbb3SSatish Balay ii = a->i; 115126e093fcSHong Zhang y = yarray; 115226e093fcSHong Zhang z = zarray; 115326e093fcSHong Zhang } 11542d61bbb3SSatish Balay 11552d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 11562d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 115726e093fcSHong Zhang if (usecprow) { 11587b2bb3b9SHong Zhang z = zarray + 2*ridx[i]; 11597b2bb3b9SHong Zhang y = yarray + 2*ridx[i]; 116026e093fcSHong Zhang } 11612d61bbb3SSatish Balay sum1 = y[0]; sum2 = y[1]; 1162444d8c10SJed Brown PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA); /* Indices for the next row (assumes same size as this one) */ 1163444d8c10SJed Brown PetscPrefetchBlock(v+4*n,4*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */ 11642d61bbb3SSatish Balay for (j=0; j<n; j++) { 116526fbe8dcSKarl Rupp xb = x + 2*(*idx++); 116626fbe8dcSKarl Rupp x1 = xb[0]; 116726fbe8dcSKarl Rupp x2 = xb[1]; 116826fbe8dcSKarl Rupp 11692d61bbb3SSatish Balay sum1 += v[0]*x1 + v[2]*x2; 11702d61bbb3SSatish Balay sum2 += v[1]*x1 + v[3]*x2; 11712d61bbb3SSatish Balay v += 4; 11722d61bbb3SSatish Balay } 11732d61bbb3SSatish Balay z[0] = sum1; z[1] = sum2; 117426e093fcSHong Zhang if (!usecprow) { 11752d61bbb3SSatish Balay z += 2; y += 2; 11762d61bbb3SSatish Balay } 117726e093fcSHong Zhang } 1178d9ca1df4SBarry Smith ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr); 1179d9ca1df4SBarry Smith ierr = VecRestoreArrayPair(yy,zz,&yarray,&zarray);CHKERRQ(ierr); 1180dc0b31edSSatish Balay ierr = PetscLogFlops(4.0*a->nz);CHKERRQ(ierr); 11812d61bbb3SSatish Balay PetscFunctionReturn(0); 11822d61bbb3SSatish Balay } 11832d61bbb3SSatish Balay 11844a2ae208SSatish Balay #undef __FUNCT__ 11854a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_3" 1186dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_3(Mat A,Vec xx,Vec yy,Vec zz) 11872d61bbb3SSatish Balay { 11882d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1189d9ca1df4SBarry Smith PetscScalar *y = 0,*z = 0,sum1,sum2,sum3,x1,x2,x3,*yarray,*zarray; 1190d9ca1df4SBarry Smith const PetscScalar *x,*xb; 1191d9ca1df4SBarry Smith const MatScalar *v; 1192dfbe8321SBarry Smith PetscErrorCode ierr; 1193d9ca1df4SBarry Smith PetscInt mbs = a->mbs,i,j,n; 1194d9ca1df4SBarry Smith const PetscInt *idx,*ii,*ridx = NULL; 1195ace3abfcSBarry Smith PetscBool usecprow = a->compressedrow.use; 11962d61bbb3SSatish Balay 11972d61bbb3SSatish Balay PetscFunctionBegin; 1198d9ca1df4SBarry Smith ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr); 1199d9ca1df4SBarry Smith ierr = VecGetArrayPair(yy,zz,&yarray,&zarray);CHKERRQ(ierr); 12002d61bbb3SSatish Balay 12012d61bbb3SSatish Balay idx = a->j; 12022d61bbb3SSatish Balay v = a->a; 120326e093fcSHong Zhang if (usecprow) { 12044eb6d288SHong Zhang if (zz != yy) { 12054eb6d288SHong Zhang ierr = PetscMemcpy(zarray,yarray,3*mbs*sizeof(PetscScalar));CHKERRQ(ierr); 12064eb6d288SHong Zhang } 120726e093fcSHong Zhang mbs = a->compressedrow.nrows; 120826e093fcSHong Zhang ii = a->compressedrow.i; 12097b2bb3b9SHong Zhang ridx = a->compressedrow.rindex; 121026e093fcSHong Zhang } else { 12112d61bbb3SSatish Balay ii = a->i; 121226e093fcSHong Zhang y = yarray; 121326e093fcSHong Zhang z = zarray; 121426e093fcSHong Zhang } 12152d61bbb3SSatish Balay 12162d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 12172d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 121826e093fcSHong Zhang if (usecprow) { 12197b2bb3b9SHong Zhang z = zarray + 3*ridx[i]; 12207b2bb3b9SHong Zhang y = yarray + 3*ridx[i]; 122126e093fcSHong Zhang } 12222d61bbb3SSatish Balay sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; 1223444d8c10SJed Brown PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA); /* Indices for the next row (assumes same size as this one) */ 1224444d8c10SJed Brown PetscPrefetchBlock(v+9*n,9*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */ 12252d61bbb3SSatish Balay for (j=0; j<n; j++) { 12262d61bbb3SSatish Balay xb = x + 3*(*idx++); x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; 12272d61bbb3SSatish Balay sum1 += v[0]*x1 + v[3]*x2 + v[6]*x3; 12282d61bbb3SSatish Balay sum2 += v[1]*x1 + v[4]*x2 + v[7]*x3; 12292d61bbb3SSatish Balay sum3 += v[2]*x1 + v[5]*x2 + v[8]*x3; 12302d61bbb3SSatish Balay v += 9; 12312d61bbb3SSatish Balay } 12322d61bbb3SSatish Balay z[0] = sum1; z[1] = sum2; z[2] = sum3; 123326e093fcSHong Zhang if (!usecprow) { 12342d61bbb3SSatish Balay z += 3; y += 3; 12352d61bbb3SSatish Balay } 123626e093fcSHong Zhang } 1237d9ca1df4SBarry Smith ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr); 1238d9ca1df4SBarry Smith ierr = VecRestoreArrayPair(yy,zz,&yarray,&zarray);CHKERRQ(ierr); 1239dc0b31edSSatish Balay ierr = PetscLogFlops(18.0*a->nz);CHKERRQ(ierr); 12402d61bbb3SSatish Balay PetscFunctionReturn(0); 12412d61bbb3SSatish Balay } 12422d61bbb3SSatish Balay 12434a2ae208SSatish Balay #undef __FUNCT__ 12444a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_4" 1245dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_4(Mat A,Vec xx,Vec yy,Vec zz) 12462d61bbb3SSatish Balay { 12472d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1248d9ca1df4SBarry Smith PetscScalar *y = 0,*z = 0,sum1,sum2,sum3,sum4,x1,x2,x3,x4,*yarray,*zarray; 1249d9ca1df4SBarry Smith const PetscScalar *x,*xb; 1250d9ca1df4SBarry Smith const MatScalar *v; 1251dfbe8321SBarry Smith PetscErrorCode ierr; 1252d9ca1df4SBarry Smith PetscInt mbs = a->mbs,i,j,n; 1253d9ca1df4SBarry Smith const PetscInt *idx,*ii,*ridx=NULL; 1254ace3abfcSBarry Smith PetscBool usecprow=a->compressedrow.use; 12552d61bbb3SSatish Balay 12562d61bbb3SSatish Balay PetscFunctionBegin; 1257d9ca1df4SBarry Smith ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr); 1258d9ca1df4SBarry Smith ierr = VecGetArrayPair(yy,zz,&yarray,&zarray);CHKERRQ(ierr); 12592d61bbb3SSatish Balay 12602d61bbb3SSatish Balay idx = a->j; 12612d61bbb3SSatish Balay v = a->a; 126226e093fcSHong Zhang if (usecprow) { 12634eb6d288SHong Zhang if (zz != yy) { 12644eb6d288SHong Zhang ierr = PetscMemcpy(zarray,yarray,4*mbs*sizeof(PetscScalar));CHKERRQ(ierr); 12654eb6d288SHong Zhang } 126626e093fcSHong Zhang mbs = a->compressedrow.nrows; 126726e093fcSHong Zhang ii = a->compressedrow.i; 12687b2bb3b9SHong Zhang ridx = a->compressedrow.rindex; 126926e093fcSHong Zhang } else { 12702d61bbb3SSatish Balay ii = a->i; 127126e093fcSHong Zhang y = yarray; 127226e093fcSHong Zhang z = zarray; 127326e093fcSHong Zhang } 12742d61bbb3SSatish Balay 12752d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 12762d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 127726e093fcSHong Zhang if (usecprow) { 12787b2bb3b9SHong Zhang z = zarray + 4*ridx[i]; 12797b2bb3b9SHong Zhang y = yarray + 4*ridx[i]; 128026e093fcSHong Zhang } 12812d61bbb3SSatish Balay sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; 1282444d8c10SJed Brown PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA); /* Indices for the next row (assumes same size as this one) */ 1283444d8c10SJed Brown PetscPrefetchBlock(v+16*n,16*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */ 12842d61bbb3SSatish Balay for (j=0; j<n; j++) { 12852d61bbb3SSatish Balay xb = x + 4*(*idx++); 12862d61bbb3SSatish Balay x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; 12872d61bbb3SSatish Balay sum1 += v[0]*x1 + v[4]*x2 + v[8]*x3 + v[12]*x4; 12882d61bbb3SSatish Balay sum2 += v[1]*x1 + v[5]*x2 + v[9]*x3 + v[13]*x4; 12892d61bbb3SSatish Balay sum3 += v[2]*x1 + v[6]*x2 + v[10]*x3 + v[14]*x4; 12902d61bbb3SSatish Balay sum4 += v[3]*x1 + v[7]*x2 + v[11]*x3 + v[15]*x4; 12912d61bbb3SSatish Balay v += 16; 12922d61bbb3SSatish Balay } 12932d61bbb3SSatish Balay z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; 129426e093fcSHong Zhang if (!usecprow) { 12952d61bbb3SSatish Balay z += 4; y += 4; 12962d61bbb3SSatish Balay } 129726e093fcSHong Zhang } 1298d9ca1df4SBarry Smith ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr); 1299d9ca1df4SBarry Smith ierr = VecRestoreArrayPair(yy,zz,&yarray,&zarray);CHKERRQ(ierr); 1300dc0b31edSSatish Balay ierr = PetscLogFlops(32.0*a->nz);CHKERRQ(ierr); 13012d61bbb3SSatish Balay PetscFunctionReturn(0); 13022d61bbb3SSatish Balay } 13032d61bbb3SSatish Balay 13044a2ae208SSatish Balay #undef __FUNCT__ 13054a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_5" 1306dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_5(Mat A,Vec xx,Vec yy,Vec zz) 13072d61bbb3SSatish Balay { 13082d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1309d9ca1df4SBarry Smith PetscScalar *y = 0,*z = 0,sum1,sum2,sum3,sum4,sum5,x1,x2,x3,x4,x5; 1310d9ca1df4SBarry Smith const PetscScalar *x,*xb; 131126e093fcSHong Zhang PetscScalar *yarray,*zarray; 1312d9ca1df4SBarry Smith const MatScalar *v; 1313dfbe8321SBarry Smith PetscErrorCode ierr; 1314d9ca1df4SBarry Smith PetscInt mbs = a->mbs,i,j,n; 1315d9ca1df4SBarry Smith const PetscInt *idx,*ii,*ridx = NULL; 1316ace3abfcSBarry Smith PetscBool usecprow=a->compressedrow.use; 13172d61bbb3SSatish Balay 13182d61bbb3SSatish Balay PetscFunctionBegin; 1319d9ca1df4SBarry Smith ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr); 1320d9ca1df4SBarry Smith ierr = VecGetArrayPair(yy,zz,&yarray,&zarray);CHKERRQ(ierr); 13212d61bbb3SSatish Balay 13222d61bbb3SSatish Balay idx = a->j; 13232d61bbb3SSatish Balay v = a->a; 132426e093fcSHong Zhang if (usecprow) { 13254eb6d288SHong Zhang if (zz != yy) { 13264eb6d288SHong Zhang ierr = PetscMemcpy(zarray,yarray,5*mbs*sizeof(PetscScalar));CHKERRQ(ierr); 13274eb6d288SHong Zhang } 132826e093fcSHong Zhang mbs = a->compressedrow.nrows; 132926e093fcSHong Zhang ii = a->compressedrow.i; 13307b2bb3b9SHong Zhang ridx = a->compressedrow.rindex; 133126e093fcSHong Zhang } else { 13322d61bbb3SSatish Balay ii = a->i; 133326e093fcSHong Zhang y = yarray; 133426e093fcSHong Zhang z = zarray; 133526e093fcSHong Zhang } 13362d61bbb3SSatish Balay 13372d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 13382d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 133926e093fcSHong Zhang if (usecprow) { 13407b2bb3b9SHong Zhang z = zarray + 5*ridx[i]; 13417b2bb3b9SHong Zhang y = yarray + 5*ridx[i]; 134226e093fcSHong Zhang } 13432d61bbb3SSatish Balay sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4]; 1344444d8c10SJed Brown PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA); /* Indices for the next row (assumes same size as this one) */ 1345444d8c10SJed Brown PetscPrefetchBlock(v+25*n,25*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */ 13462d61bbb3SSatish Balay for (j=0; j<n; j++) { 13472d61bbb3SSatish Balay xb = x + 5*(*idx++); 13482d61bbb3SSatish Balay x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; 13492d61bbb3SSatish Balay sum1 += v[0]*x1 + v[5]*x2 + v[10]*x3 + v[15]*x4 + v[20]*x5; 13502d61bbb3SSatish Balay sum2 += v[1]*x1 + v[6]*x2 + v[11]*x3 + v[16]*x4 + v[21]*x5; 13512d61bbb3SSatish Balay sum3 += v[2]*x1 + v[7]*x2 + v[12]*x3 + v[17]*x4 + v[22]*x5; 13522d61bbb3SSatish Balay sum4 += v[3]*x1 + v[8]*x2 + v[13]*x3 + v[18]*x4 + v[23]*x5; 13532d61bbb3SSatish Balay sum5 += v[4]*x1 + v[9]*x2 + v[14]*x3 + v[19]*x4 + v[24]*x5; 13542d61bbb3SSatish Balay v += 25; 13552d61bbb3SSatish Balay } 13562d61bbb3SSatish Balay z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; 135726e093fcSHong Zhang if (!usecprow) { 13582d61bbb3SSatish Balay z += 5; y += 5; 13592d61bbb3SSatish Balay } 136026e093fcSHong Zhang } 1361d9ca1df4SBarry Smith ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr); 1362d9ca1df4SBarry Smith ierr = VecRestoreArrayPair(yy,zz,&yarray,&zarray);CHKERRQ(ierr); 1363dc0b31edSSatish Balay ierr = PetscLogFlops(50.0*a->nz);CHKERRQ(ierr); 13642d61bbb3SSatish Balay PetscFunctionReturn(0); 13652d61bbb3SSatish Balay } 13664a2ae208SSatish Balay #undef __FUNCT__ 13674a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_6" 1368dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_6(Mat A,Vec xx,Vec yy,Vec zz) 136915091d37SBarry Smith { 137015091d37SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1371d9ca1df4SBarry Smith PetscScalar *y = 0,*z = 0,sum1,sum2,sum3,sum4,sum5,sum6; 1372d9ca1df4SBarry Smith const PetscScalar *x,*xb; 137326e093fcSHong Zhang PetscScalar x1,x2,x3,x4,x5,x6,*yarray,*zarray; 1374d9ca1df4SBarry Smith const MatScalar *v; 1375dfbe8321SBarry Smith PetscErrorCode ierr; 1376d9ca1df4SBarry Smith PetscInt mbs = a->mbs,i,j,n; 1377d9ca1df4SBarry Smith const PetscInt *idx,*ii,*ridx=NULL; 1378ace3abfcSBarry Smith PetscBool usecprow=a->compressedrow.use; 137915091d37SBarry Smith 138015091d37SBarry Smith PetscFunctionBegin; 1381d9ca1df4SBarry Smith ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr); 1382d9ca1df4SBarry Smith ierr = VecGetArrayPair(yy,zz,&yarray,&zarray);CHKERRQ(ierr); 138315091d37SBarry Smith 138415091d37SBarry Smith idx = a->j; 138515091d37SBarry Smith v = a->a; 138626e093fcSHong Zhang if (usecprow) { 13874eb6d288SHong Zhang if (zz != yy) { 13884eb6d288SHong Zhang ierr = PetscMemcpy(zarray,yarray,6*mbs*sizeof(PetscScalar));CHKERRQ(ierr); 13894eb6d288SHong Zhang } 139026e093fcSHong Zhang mbs = a->compressedrow.nrows; 139126e093fcSHong Zhang ii = a->compressedrow.i; 13927b2bb3b9SHong Zhang ridx = a->compressedrow.rindex; 139326e093fcSHong Zhang } else { 139415091d37SBarry Smith ii = a->i; 139526e093fcSHong Zhang y = yarray; 139626e093fcSHong Zhang z = zarray; 139726e093fcSHong Zhang } 139815091d37SBarry Smith 139915091d37SBarry Smith for (i=0; i<mbs; i++) { 140015091d37SBarry Smith n = ii[1] - ii[0]; ii++; 140126e093fcSHong Zhang if (usecprow) { 14027b2bb3b9SHong Zhang z = zarray + 6*ridx[i]; 14037b2bb3b9SHong Zhang y = yarray + 6*ridx[i]; 140426e093fcSHong Zhang } 140515091d37SBarry Smith sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4]; sum6 = y[5]; 1406444d8c10SJed Brown PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA); /* Indices for the next row (assumes same size as this one) */ 1407444d8c10SJed Brown PetscPrefetchBlock(v+36*n,36*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */ 140815091d37SBarry Smith for (j=0; j<n; j++) { 14093b95cb0eSSatish Balay xb = x + 6*(*idx++); 141015091d37SBarry Smith x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; 141115091d37SBarry Smith sum1 += v[0]*x1 + v[6]*x2 + v[12]*x3 + v[18]*x4 + v[24]*x5 + v[30]*x6; 141215091d37SBarry Smith sum2 += v[1]*x1 + v[7]*x2 + v[13]*x3 + v[19]*x4 + v[25]*x5 + v[31]*x6; 141315091d37SBarry Smith sum3 += v[2]*x1 + v[8]*x2 + v[14]*x3 + v[20]*x4 + v[26]*x5 + v[32]*x6; 141415091d37SBarry Smith sum4 += v[3]*x1 + v[9]*x2 + v[15]*x3 + v[21]*x4 + v[27]*x5 + v[33]*x6; 141515091d37SBarry Smith sum5 += v[4]*x1 + v[10]*x2 + v[16]*x3 + v[22]*x4 + v[28]*x5 + v[34]*x6; 141615091d37SBarry Smith sum6 += v[5]*x1 + v[11]*x2 + v[17]*x3 + v[23]*x4 + v[29]*x5 + v[35]*x6; 141715091d37SBarry Smith v += 36; 141815091d37SBarry Smith } 141915091d37SBarry Smith z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; 142026e093fcSHong Zhang if (!usecprow) { 142115091d37SBarry Smith z += 6; y += 6; 142215091d37SBarry Smith } 142326e093fcSHong Zhang } 1424d9ca1df4SBarry Smith ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr); 1425d9ca1df4SBarry Smith ierr = VecRestoreArrayPair(yy,zz,&yarray,&zarray);CHKERRQ(ierr); 1426dc0b31edSSatish Balay ierr = PetscLogFlops(72.0*a->nz);CHKERRQ(ierr); 142715091d37SBarry Smith PetscFunctionReturn(0); 142815091d37SBarry Smith } 14292d61bbb3SSatish Balay 14304a2ae208SSatish Balay #undef __FUNCT__ 14314a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_7" 1432dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_7(Mat A,Vec xx,Vec yy,Vec zz) 14332d61bbb3SSatish Balay { 14342d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1435d9ca1df4SBarry Smith PetscScalar *y = 0,*z = 0,sum1,sum2,sum3,sum4,sum5,sum6,sum7; 1436d9ca1df4SBarry Smith const PetscScalar *x,*xb; 143726e093fcSHong Zhang PetscScalar x1,x2,x3,x4,x5,x6,x7,*yarray,*zarray; 1438d9ca1df4SBarry Smith const MatScalar *v; 1439dfbe8321SBarry Smith PetscErrorCode ierr; 1440d9ca1df4SBarry Smith PetscInt mbs = a->mbs,i,j,n; 1441d9ca1df4SBarry Smith const PetscInt *idx,*ii,*ridx = NULL; 1442ace3abfcSBarry Smith PetscBool usecprow=a->compressedrow.use; 14432d61bbb3SSatish Balay 14442d61bbb3SSatish Balay PetscFunctionBegin; 1445d9ca1df4SBarry Smith ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr); 1446d9ca1df4SBarry Smith ierr = VecGetArrayPair(yy,zz,&yarray,&zarray);CHKERRQ(ierr); 14472d61bbb3SSatish Balay 14482d61bbb3SSatish Balay idx = a->j; 14492d61bbb3SSatish Balay v = a->a; 145026e093fcSHong Zhang if (usecprow) { 14514eb6d288SHong Zhang if (zz != yy) { 14524eb6d288SHong Zhang ierr = PetscMemcpy(zarray,yarray,7*mbs*sizeof(PetscScalar));CHKERRQ(ierr); 14534eb6d288SHong Zhang } 145426e093fcSHong Zhang mbs = a->compressedrow.nrows; 145526e093fcSHong Zhang ii = a->compressedrow.i; 14567b2bb3b9SHong Zhang ridx = a->compressedrow.rindex; 145726e093fcSHong Zhang } else { 14582d61bbb3SSatish Balay ii = a->i; 145926e093fcSHong Zhang y = yarray; 146026e093fcSHong Zhang z = zarray; 146126e093fcSHong Zhang } 14622d61bbb3SSatish Balay 14632d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 14642d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 146526e093fcSHong Zhang if (usecprow) { 14667b2bb3b9SHong Zhang z = zarray + 7*ridx[i]; 14677b2bb3b9SHong Zhang y = yarray + 7*ridx[i]; 146826e093fcSHong Zhang } 14692d61bbb3SSatish Balay sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4]; sum6 = y[5]; sum7 = y[6]; 1470444d8c10SJed Brown PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA); /* Indices for the next row (assumes same size as this one) */ 1471444d8c10SJed Brown PetscPrefetchBlock(v+49*n,49*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */ 14722d61bbb3SSatish Balay for (j=0; j<n; j++) { 14732d61bbb3SSatish Balay xb = x + 7*(*idx++); 14742d61bbb3SSatish Balay x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6]; 14752d61bbb3SSatish Balay sum1 += v[0]*x1 + v[7]*x2 + v[14]*x3 + v[21]*x4 + v[28]*x5 + v[35]*x6 + v[42]*x7; 14762d61bbb3SSatish Balay sum2 += v[1]*x1 + v[8]*x2 + v[15]*x3 + v[22]*x4 + v[29]*x5 + v[36]*x6 + v[43]*x7; 14772d61bbb3SSatish Balay sum3 += v[2]*x1 + v[9]*x2 + v[16]*x3 + v[23]*x4 + v[30]*x5 + v[37]*x6 + v[44]*x7; 14782d61bbb3SSatish Balay sum4 += v[3]*x1 + v[10]*x2 + v[17]*x3 + v[24]*x4 + v[31]*x5 + v[38]*x6 + v[45]*x7; 14792d61bbb3SSatish Balay sum5 += v[4]*x1 + v[11]*x2 + v[18]*x3 + v[25]*x4 + v[32]*x5 + v[39]*x6 + v[46]*x7; 14802d61bbb3SSatish Balay sum6 += v[5]*x1 + v[12]*x2 + v[19]*x3 + v[26]*x4 + v[33]*x5 + v[40]*x6 + v[47]*x7; 14812d61bbb3SSatish Balay sum7 += v[6]*x1 + v[13]*x2 + v[20]*x3 + v[27]*x4 + v[34]*x5 + v[41]*x6 + v[48]*x7; 14822d61bbb3SSatish Balay v += 49; 14832d61bbb3SSatish Balay } 14842d61bbb3SSatish Balay z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7; 148526e093fcSHong Zhang if (!usecprow) { 14862d61bbb3SSatish Balay z += 7; y += 7; 14872d61bbb3SSatish Balay } 148826e093fcSHong Zhang } 1489d9ca1df4SBarry Smith ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr); 1490d9ca1df4SBarry Smith ierr = VecRestoreArrayPair(yy,zz,&yarray,&zarray);CHKERRQ(ierr); 1491dc0b31edSSatish Balay ierr = PetscLogFlops(98.0*a->nz);CHKERRQ(ierr); 14922d61bbb3SSatish Balay PetscFunctionReturn(0); 14932d61bbb3SSatish Balay } 1494218c64b6SSatish Balay 14954a2ae208SSatish Balay #undef __FUNCT__ 14964a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_N" 1497dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_N(Mat A,Vec xx,Vec yy,Vec zz) 14982d61bbb3SSatish Balay { 14992d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1500d9ca1df4SBarry Smith PetscScalar *z = 0,*work,*workt,*zarray; 1501d9ca1df4SBarry Smith const PetscScalar *x,*xb; 1502d9ca1df4SBarry Smith const MatScalar *v; 15036849ba73SBarry Smith PetscErrorCode ierr; 1504d9ca1df4SBarry Smith PetscInt mbs,i,bs=A->rmap->bs,j,n,bs2=a->bs2; 1505d9ca1df4SBarry Smith PetscInt ncols,k; 1506d9ca1df4SBarry Smith const PetscInt *ridx = NULL,*idx,*ii; 1507ace3abfcSBarry Smith PetscBool usecprow = a->compressedrow.use; 1508218c64b6SSatish Balay 15092d61bbb3SSatish Balay PetscFunctionBegin; 1510343551c4SBarry Smith ierr = VecCopy(yy,zz);CHKERRQ(ierr); 1511d9ca1df4SBarry Smith ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr); 151226e093fcSHong Zhang ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr); 15132d61bbb3SSatish Balay 15142d61bbb3SSatish Balay idx = a->j; 15152d61bbb3SSatish Balay v = a->a; 151626e093fcSHong Zhang if (usecprow) { 151726e093fcSHong Zhang mbs = a->compressedrow.nrows; 151826e093fcSHong Zhang ii = a->compressedrow.i; 15197b2bb3b9SHong Zhang ridx = a->compressedrow.rindex; 152026e093fcSHong Zhang } else { 152126e093fcSHong Zhang mbs = a->mbs; 15222d61bbb3SSatish Balay ii = a->i; 152326e093fcSHong Zhang z = zarray; 152426e093fcSHong Zhang } 15252d61bbb3SSatish Balay 15262d61bbb3SSatish Balay if (!a->mult_work) { 1527d0f46423SBarry Smith k = PetscMax(A->rmap->n,A->cmap->n); 1528854ce69bSBarry Smith ierr = PetscMalloc1(k+1,&a->mult_work);CHKERRQ(ierr); 15292d61bbb3SSatish Balay } 15302d61bbb3SSatish Balay work = a->mult_work; 15312d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 15322d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 15332d61bbb3SSatish Balay ncols = n*bs; 15342d61bbb3SSatish Balay workt = work; 15352d61bbb3SSatish Balay for (j=0; j<n; j++) { 15362d61bbb3SSatish Balay xb = x + bs*(*idx++); 15372d61bbb3SSatish Balay for (k=0; k<bs; k++) workt[k] = xb[k]; 15382d61bbb3SSatish Balay workt += bs; 15392d61bbb3SSatish Balay } 15407b2bb3b9SHong Zhang if (usecprow) z = zarray + bs*ridx[i]; 154196b95a6bSBarry Smith PetscKernel_w_gets_w_plus_Ar_times_v(bs,ncols,work,v,z); 154271044d3cSBarry Smith /* BLASgemv_("N",&bs,&ncols,&_DOne,v,&bs,work,&_One,&_DOne,z,&_One); */ 15432d61bbb3SSatish Balay v += n*bs2; 154426fbe8dcSKarl Rupp if (!usecprow) z += bs; 154526e093fcSHong Zhang } 1546d9ca1df4SBarry Smith ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr); 154726e093fcSHong Zhang ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr); 1548dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*a->nz*bs2);CHKERRQ(ierr); 15492d61bbb3SSatish Balay PetscFunctionReturn(0); 15502d61bbb3SSatish Balay } 15512d61bbb3SSatish Balay 15524a2ae208SSatish Balay #undef __FUNCT__ 1553547795f9SHong Zhang #define __FUNCT__ "MatMultHermitianTranspose_SeqBAIJ" 1554547795f9SHong Zhang PetscErrorCode MatMultHermitianTranspose_SeqBAIJ(Mat A,Vec xx,Vec zz) 1555547795f9SHong Zhang { 1556547795f9SHong Zhang PetscScalar zero = 0.0; 1557547795f9SHong Zhang PetscErrorCode ierr; 1558547795f9SHong Zhang 1559547795f9SHong Zhang PetscFunctionBegin; 1560547795f9SHong Zhang ierr = VecSet(zz,zero);CHKERRQ(ierr); 1561547795f9SHong Zhang ierr = MatMultHermitianTransposeAdd_SeqBAIJ(A,xx,zz,zz);CHKERRQ(ierr); 1562547795f9SHong Zhang PetscFunctionReturn(0); 1563547795f9SHong Zhang } 1564547795f9SHong Zhang 1565547795f9SHong Zhang #undef __FUNCT__ 15664a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_SeqBAIJ" 1567dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_SeqBAIJ(Mat A,Vec xx,Vec zz) 15682d61bbb3SSatish Balay { 15693447b6efSHong Zhang PetscScalar zero = 0.0; 15706849ba73SBarry Smith PetscErrorCode ierr; 15712d61bbb3SSatish Balay 15722d61bbb3SSatish Balay PetscFunctionBegin; 15732dcb1b2aSMatthew Knepley ierr = VecSet(zz,zero);CHKERRQ(ierr); 15743447b6efSHong Zhang ierr = MatMultTransposeAdd_SeqBAIJ(A,xx,zz,zz);CHKERRQ(ierr); 15752d61bbb3SSatish Balay PetscFunctionReturn(0); 15762d61bbb3SSatish Balay } 15772d61bbb3SSatish Balay 15784a2ae208SSatish Balay #undef __FUNCT__ 1579547795f9SHong Zhang #define __FUNCT__ "MatMultHermitianTransposeAdd_SeqBAIJ" 1580547795f9SHong Zhang PetscErrorCode MatMultHermitianTransposeAdd_SeqBAIJ(Mat A,Vec xx,Vec yy,Vec zz) 1581547795f9SHong Zhang { 1582547795f9SHong Zhang Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1583b8c08b77SHong Zhang PetscScalar *z,x1,x2,x3,x4,x5; 1584d9ca1df4SBarry Smith const PetscScalar *x,*xb = NULL; 1585d9ca1df4SBarry Smith const MatScalar *v; 1586547795f9SHong Zhang PetscErrorCode ierr; 1587b8c08b77SHong Zhang PetscInt mbs,i,rval,bs=A->rmap->bs,j,n; 1588d9ca1df4SBarry Smith const PetscInt *idx,*ii,*ib,*ridx = NULL; 1589547795f9SHong Zhang Mat_CompressedRow cprow = a->compressedrow; 1590ace3abfcSBarry Smith PetscBool usecprow = cprow.use; 1591547795f9SHong Zhang 1592547795f9SHong Zhang PetscFunctionBegin; 1593343551c4SBarry Smith if (yy != zz) { ierr = VecCopy(yy,zz);CHKERRQ(ierr); } 1594d9ca1df4SBarry Smith ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr); 1595547795f9SHong Zhang ierr = VecGetArray(zz,&z);CHKERRQ(ierr); 1596547795f9SHong Zhang 1597547795f9SHong Zhang idx = a->j; 1598547795f9SHong Zhang v = a->a; 1599547795f9SHong Zhang if (usecprow) { 1600547795f9SHong Zhang mbs = cprow.nrows; 1601547795f9SHong Zhang ii = cprow.i; 1602547795f9SHong Zhang ridx = cprow.rindex; 1603547795f9SHong Zhang } else { 1604547795f9SHong Zhang mbs=a->mbs; 1605547795f9SHong Zhang ii = a->i; 1606547795f9SHong Zhang xb = x; 1607547795f9SHong Zhang } 1608547795f9SHong Zhang 1609547795f9SHong Zhang switch (bs) { 1610547795f9SHong Zhang case 1: 1611547795f9SHong Zhang for (i=0; i<mbs; i++) { 1612547795f9SHong Zhang if (usecprow) xb = x + ridx[i]; 1613547795f9SHong Zhang x1 = xb[0]; 1614547795f9SHong Zhang ib = idx + ii[0]; 1615547795f9SHong Zhang n = ii[1] - ii[0]; ii++; 1616547795f9SHong Zhang for (j=0; j<n; j++) { 1617547795f9SHong Zhang rval = ib[j]; 1618547795f9SHong Zhang z[rval] += PetscConj(*v) * x1; 1619547795f9SHong Zhang v++; 1620547795f9SHong Zhang } 1621547795f9SHong Zhang if (!usecprow) xb++; 1622547795f9SHong Zhang } 1623547795f9SHong Zhang break; 1624547795f9SHong Zhang case 2: 1625547795f9SHong Zhang for (i=0; i<mbs; i++) { 1626547795f9SHong Zhang if (usecprow) xb = x + 2*ridx[i]; 1627547795f9SHong Zhang x1 = xb[0]; x2 = xb[1]; 1628547795f9SHong Zhang ib = idx + ii[0]; 1629547795f9SHong Zhang n = ii[1] - ii[0]; ii++; 1630547795f9SHong Zhang for (j=0; j<n; j++) { 1631547795f9SHong Zhang rval = ib[j]*2; 1632547795f9SHong Zhang z[rval++] += PetscConj(v[0])*x1 + PetscConj(v[1])*x2; 1633547795f9SHong Zhang z[rval++] += PetscConj(v[2])*x1 + PetscConj(v[3])*x2; 1634547795f9SHong Zhang v += 4; 1635547795f9SHong Zhang } 1636547795f9SHong Zhang if (!usecprow) xb += 2; 1637547795f9SHong Zhang } 1638547795f9SHong Zhang break; 1639547795f9SHong Zhang case 3: 1640547795f9SHong Zhang for (i=0; i<mbs; i++) { 1641547795f9SHong Zhang if (usecprow) xb = x + 3*ridx[i]; 1642547795f9SHong Zhang x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; 1643547795f9SHong Zhang ib = idx + ii[0]; 1644547795f9SHong Zhang n = ii[1] - ii[0]; ii++; 1645547795f9SHong Zhang for (j=0; j<n; j++) { 1646547795f9SHong Zhang rval = ib[j]*3; 1647547795f9SHong Zhang z[rval++] += PetscConj(v[0])*x1 + PetscConj(v[1])*x2 + PetscConj(v[2])*x3; 1648547795f9SHong Zhang z[rval++] += PetscConj(v[3])*x1 + PetscConj(v[4])*x2 + PetscConj(v[5])*x3; 1649547795f9SHong Zhang z[rval++] += PetscConj(v[6])*x1 + PetscConj(v[7])*x2 + PetscConj(v[8])*x3; 1650547795f9SHong Zhang v += 9; 1651547795f9SHong Zhang } 1652547795f9SHong Zhang if (!usecprow) xb += 3; 1653547795f9SHong Zhang } 1654547795f9SHong Zhang break; 1655547795f9SHong Zhang case 4: 1656547795f9SHong Zhang for (i=0; i<mbs; i++) { 1657547795f9SHong Zhang if (usecprow) xb = x + 4*ridx[i]; 1658547795f9SHong Zhang x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; 1659547795f9SHong Zhang ib = idx + ii[0]; 1660547795f9SHong Zhang n = ii[1] - ii[0]; ii++; 1661547795f9SHong Zhang for (j=0; j<n; j++) { 1662547795f9SHong Zhang rval = ib[j]*4; 1663547795f9SHong Zhang z[rval++] += PetscConj(v[0])*x1 + PetscConj(v[1])*x2 + PetscConj(v[2])*x3 + PetscConj(v[3])*x4; 1664547795f9SHong Zhang z[rval++] += PetscConj(v[4])*x1 + PetscConj(v[5])*x2 + PetscConj(v[6])*x3 + PetscConj(v[7])*x4; 1665547795f9SHong Zhang z[rval++] += PetscConj(v[8])*x1 + PetscConj(v[9])*x2 + PetscConj(v[10])*x3 + PetscConj(v[11])*x4; 1666547795f9SHong Zhang z[rval++] += PetscConj(v[12])*x1 + PetscConj(v[13])*x2 + PetscConj(v[14])*x3 + PetscConj(v[15])*x4; 1667547795f9SHong Zhang v += 16; 1668547795f9SHong Zhang } 1669547795f9SHong Zhang if (!usecprow) xb += 4; 1670547795f9SHong Zhang } 1671547795f9SHong Zhang break; 1672547795f9SHong Zhang case 5: 1673547795f9SHong Zhang for (i=0; i<mbs; i++) { 1674547795f9SHong Zhang if (usecprow) xb = x + 5*ridx[i]; 1675547795f9SHong Zhang x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; 1676547795f9SHong Zhang x4 = xb[3]; x5 = xb[4]; 1677547795f9SHong Zhang ib = idx + ii[0]; 1678547795f9SHong Zhang n = ii[1] - ii[0]; ii++; 1679547795f9SHong Zhang for (j=0; j<n; j++) { 1680547795f9SHong Zhang rval = ib[j]*5; 1681547795f9SHong Zhang z[rval++] += PetscConj(v[0])*x1 + PetscConj(v[1])*x2 + PetscConj(v[2])*x3 + PetscConj(v[3])*x4 + PetscConj(v[4])*x5; 1682547795f9SHong Zhang z[rval++] += PetscConj(v[5])*x1 + PetscConj(v[6])*x2 + PetscConj(v[7])*x3 + PetscConj(v[8])*x4 + PetscConj(v[9])*x5; 1683547795f9SHong Zhang z[rval++] += PetscConj(v[10])*x1 + PetscConj(v[11])*x2 + PetscConj(v[12])*x3 + PetscConj(v[13])*x4 + PetscConj(v[14])*x5; 1684547795f9SHong Zhang z[rval++] += PetscConj(v[15])*x1 + PetscConj(v[16])*x2 + PetscConj(v[17])*x3 + PetscConj(v[18])*x4 + PetscConj(v[19])*x5; 1685547795f9SHong Zhang z[rval++] += PetscConj(v[20])*x1 + PetscConj(v[21])*x2 + PetscConj(v[22])*x3 + PetscConj(v[23])*x4 + PetscConj(v[24])*x5; 1686547795f9SHong Zhang v += 25; 1687547795f9SHong Zhang } 1688547795f9SHong Zhang if (!usecprow) xb += 5; 1689547795f9SHong Zhang } 1690547795f9SHong Zhang break; 1691968ae2c8SSatish Balay default: /* block sizes larger than 5 by 5 are handled by BLAS */ 1692968ae2c8SSatish Balay SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"block size larger than 5 is not supported yet"); 1693968ae2c8SSatish Balay #if 0 1694968ae2c8SSatish Balay { 1695b8c08b77SHong Zhang PetscInt ncols,k,bs2=a->bs2; 1696b8c08b77SHong Zhang PetscScalar *work,*workt,zb; 1697d9ca1df4SBarry Smith const PetscScalar *xtmp; 1698547795f9SHong Zhang if (!a->mult_work) { 1699547795f9SHong Zhang k = PetscMax(A->rmap->n,A->cmap->n); 1700854ce69bSBarry Smith ierr = PetscMalloc1(k+1,&a->mult_work);CHKERRQ(ierr); 1701547795f9SHong Zhang } 1702547795f9SHong Zhang work = a->mult_work; 1703547795f9SHong Zhang xtmp = x; 1704547795f9SHong Zhang for (i=0; i<mbs; i++) { 1705547795f9SHong Zhang n = ii[1] - ii[0]; ii++; 1706547795f9SHong Zhang ncols = n*bs; 1707547795f9SHong Zhang ierr = PetscMemzero(work,ncols*sizeof(PetscScalar));CHKERRQ(ierr); 170826fbe8dcSKarl Rupp if (usecprow) xtmp = x + bs*ridx[i]; 170996b95a6bSBarry Smith PetscKernel_w_gets_w_plus_trans_Ar_times_v(bs,ncols,xtmp,v,work); 1710547795f9SHong Zhang /* BLASgemv_("T",&bs,&ncols,&_DOne,v,&bs,xtmp,&_One,&_DOne,work,&_One); */ 1711547795f9SHong Zhang v += n*bs2; 1712547795f9SHong Zhang if (!usecprow) xtmp += bs; 1713547795f9SHong Zhang workt = work; 1714547795f9SHong Zhang for (j=0; j<n; j++) { 1715547795f9SHong Zhang zb = z + bs*(*idx++); 1716547795f9SHong Zhang for (k=0; k<bs; k++) zb[k] += workt[k] ; 1717547795f9SHong Zhang workt += bs; 1718547795f9SHong Zhang } 1719547795f9SHong Zhang } 1720547795f9SHong Zhang } 1721968ae2c8SSatish Balay #endif 1722547795f9SHong Zhang } 1723d9ca1df4SBarry Smith ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr); 1724547795f9SHong Zhang ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr); 1725547795f9SHong Zhang ierr = PetscLogFlops(2.0*a->nz*a->bs2);CHKERRQ(ierr); 1726547795f9SHong Zhang PetscFunctionReturn(0); 1727547795f9SHong Zhang } 1728547795f9SHong Zhang 1729547795f9SHong Zhang #undef __FUNCT__ 17304a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_SeqBAIJ" 1731dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_SeqBAIJ(Mat A,Vec xx,Vec yy,Vec zz) 17322d61bbb3SSatish Balay { 17332d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1734d9ca1df4SBarry Smith PetscScalar *zb,*z,x1,x2,x3,x4,x5; 1735d9ca1df4SBarry Smith const PetscScalar *x,*xb = 0; 1736d9ca1df4SBarry Smith const MatScalar *v; 17376849ba73SBarry Smith PetscErrorCode ierr; 1738d9ca1df4SBarry Smith PetscInt mbs,i,rval,bs=A->rmap->bs,j,n,bs2=a->bs2; 1739d9ca1df4SBarry Smith const PetscInt *idx,*ii,*ib,*ridx = NULL; 17403447b6efSHong Zhang Mat_CompressedRow cprow = a->compressedrow; 1741ace3abfcSBarry Smith PetscBool usecprow=cprow.use; 17422d61bbb3SSatish Balay 17432d61bbb3SSatish Balay PetscFunctionBegin; 1744343551c4SBarry Smith if (yy != zz) { ierr = VecCopy(yy,zz);CHKERRQ(ierr); } 1745d9ca1df4SBarry Smith ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr); 17463447b6efSHong Zhang ierr = VecGetArray(zz,&z);CHKERRQ(ierr); 17472d61bbb3SSatish Balay 17482d61bbb3SSatish Balay idx = a->j; 17492d61bbb3SSatish Balay v = a->a; 17503447b6efSHong Zhang if (usecprow) { 17513447b6efSHong Zhang mbs = cprow.nrows; 17523447b6efSHong Zhang ii = cprow.i; 17537b2bb3b9SHong Zhang ridx = cprow.rindex; 17543447b6efSHong Zhang } else { 17553447b6efSHong Zhang mbs=a->mbs; 17562d61bbb3SSatish Balay ii = a->i; 1757f1af5d2fSBarry Smith xb = x; 17583447b6efSHong Zhang } 17592d61bbb3SSatish Balay 17602d61bbb3SSatish Balay switch (bs) { 17612d61bbb3SSatish Balay case 1: 17622d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 17637b2bb3b9SHong Zhang if (usecprow) xb = x + ridx[i]; 1764f1af5d2fSBarry Smith x1 = xb[0]; 17653447b6efSHong Zhang ib = idx + ii[0]; 17663447b6efSHong Zhang n = ii[1] - ii[0]; ii++; 17672d61bbb3SSatish Balay for (j=0; j<n; j++) { 17682d61bbb3SSatish Balay rval = ib[j]; 1769f1af5d2fSBarry Smith z[rval] += *v * x1; 1770f1af5d2fSBarry Smith v++; 17712d61bbb3SSatish Balay } 17723447b6efSHong Zhang if (!usecprow) xb++; 17732d61bbb3SSatish Balay } 17742d61bbb3SSatish Balay break; 17752d61bbb3SSatish Balay case 2: 17762d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 17777b2bb3b9SHong Zhang if (usecprow) xb = x + 2*ridx[i]; 1778f1af5d2fSBarry Smith x1 = xb[0]; x2 = xb[1]; 17793447b6efSHong Zhang ib = idx + ii[0]; 17803447b6efSHong Zhang n = ii[1] - ii[0]; ii++; 17812d61bbb3SSatish Balay for (j=0; j<n; j++) { 17822d61bbb3SSatish Balay rval = ib[j]*2; 17832d61bbb3SSatish Balay z[rval++] += v[0]*x1 + v[1]*x2; 17842d61bbb3SSatish Balay z[rval++] += v[2]*x1 + v[3]*x2; 17852d61bbb3SSatish Balay v += 4; 17862d61bbb3SSatish Balay } 17873447b6efSHong Zhang if (!usecprow) xb += 2; 17882d61bbb3SSatish Balay } 17892d61bbb3SSatish Balay break; 17902d61bbb3SSatish Balay case 3: 17912d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 17927b2bb3b9SHong Zhang if (usecprow) xb = x + 3*ridx[i]; 1793f1af5d2fSBarry Smith x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; 17943447b6efSHong Zhang ib = idx + ii[0]; 17953447b6efSHong Zhang n = ii[1] - ii[0]; ii++; 17962d61bbb3SSatish Balay for (j=0; j<n; j++) { 17972d61bbb3SSatish Balay rval = ib[j]*3; 17982d61bbb3SSatish Balay z[rval++] += v[0]*x1 + v[1]*x2 + v[2]*x3; 17992d61bbb3SSatish Balay z[rval++] += v[3]*x1 + v[4]*x2 + v[5]*x3; 18002d61bbb3SSatish Balay z[rval++] += v[6]*x1 + v[7]*x2 + v[8]*x3; 18012d61bbb3SSatish Balay v += 9; 18022d61bbb3SSatish Balay } 18033447b6efSHong Zhang if (!usecprow) xb += 3; 18042d61bbb3SSatish Balay } 18052d61bbb3SSatish Balay break; 18062d61bbb3SSatish Balay case 4: 18072d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 18087b2bb3b9SHong Zhang if (usecprow) xb = x + 4*ridx[i]; 1809f1af5d2fSBarry Smith x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; 18103447b6efSHong Zhang ib = idx + ii[0]; 18113447b6efSHong Zhang n = ii[1] - ii[0]; ii++; 18122d61bbb3SSatish Balay for (j=0; j<n; j++) { 18132d61bbb3SSatish Balay rval = ib[j]*4; 18142d61bbb3SSatish Balay z[rval++] += v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4; 18152d61bbb3SSatish Balay z[rval++] += v[4]*x1 + v[5]*x2 + v[6]*x3 + v[7]*x4; 18162d61bbb3SSatish Balay z[rval++] += v[8]*x1 + v[9]*x2 + v[10]*x3 + v[11]*x4; 18172d61bbb3SSatish Balay z[rval++] += v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4; 18182d61bbb3SSatish Balay v += 16; 18192d61bbb3SSatish Balay } 18203447b6efSHong Zhang if (!usecprow) xb += 4; 18212d61bbb3SSatish Balay } 18222d61bbb3SSatish Balay break; 18232d61bbb3SSatish Balay case 5: 18242d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 18257b2bb3b9SHong Zhang if (usecprow) xb = x + 5*ridx[i]; 1826f1af5d2fSBarry Smith x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; 18272d61bbb3SSatish Balay x4 = xb[3]; x5 = xb[4]; 18283447b6efSHong Zhang ib = idx + ii[0]; 18293447b6efSHong Zhang n = ii[1] - ii[0]; ii++; 18302d61bbb3SSatish Balay for (j=0; j<n; j++) { 18312d61bbb3SSatish Balay rval = ib[j]*5; 18322d61bbb3SSatish Balay z[rval++] += v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4 + v[4]*x5; 18332d61bbb3SSatish Balay z[rval++] += v[5]*x1 + v[6]*x2 + v[7]*x3 + v[8]*x4 + v[9]*x5; 18342d61bbb3SSatish Balay z[rval++] += v[10]*x1 + v[11]*x2 + v[12]*x3 + v[13]*x4 + v[14]*x5; 18352d61bbb3SSatish Balay z[rval++] += v[15]*x1 + v[16]*x2 + v[17]*x3 + v[18]*x4 + v[19]*x5; 18362d61bbb3SSatish Balay z[rval++] += v[20]*x1 + v[21]*x2 + v[22]*x3 + v[23]*x4 + v[24]*x5; 18372d61bbb3SSatish Balay v += 25; 18382d61bbb3SSatish Balay } 18393447b6efSHong Zhang if (!usecprow) xb += 5; 18402d61bbb3SSatish Balay } 18412d61bbb3SSatish Balay break; 1842f1af5d2fSBarry Smith default: { /* block sizes larger then 5 by 5 are handled by BLAS */ 1843690b6cddSBarry Smith PetscInt ncols,k; 1844d9ca1df4SBarry Smith PetscScalar *work,*workt; 1845d9ca1df4SBarry Smith const PetscScalar *xtmp; 18462d61bbb3SSatish Balay if (!a->mult_work) { 1847d0f46423SBarry Smith k = PetscMax(A->rmap->n,A->cmap->n); 1848854ce69bSBarry Smith ierr = PetscMalloc1(k+1,&a->mult_work);CHKERRQ(ierr); 18492d61bbb3SSatish Balay } 18502d61bbb3SSatish Balay work = a->mult_work; 18513447b6efSHong Zhang xtmp = x; 18522d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 18532d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 18542d61bbb3SSatish Balay ncols = n*bs; 185587828ca2SBarry Smith ierr = PetscMemzero(work,ncols*sizeof(PetscScalar));CHKERRQ(ierr); 185626fbe8dcSKarl Rupp if (usecprow) xtmp = x + bs*ridx[i]; 185796b95a6bSBarry Smith PetscKernel_w_gets_w_plus_trans_Ar_times_v(bs,ncols,xtmp,v,work); 185871044d3cSBarry Smith /* BLASgemv_("T",&bs,&ncols,&_DOne,v,&bs,xtmp,&_One,&_DOne,work,&_One); */ 18592d61bbb3SSatish Balay v += n*bs2; 18603447b6efSHong Zhang if (!usecprow) xtmp += bs; 18612d61bbb3SSatish Balay workt = work; 18622d61bbb3SSatish Balay for (j=0; j<n; j++) { 18632d61bbb3SSatish Balay zb = z + bs*(*idx++); 18642d61bbb3SSatish Balay for (k=0; k<bs; k++) zb[k] += workt[k]; 18652d61bbb3SSatish Balay workt += bs; 18662d61bbb3SSatish Balay } 18672d61bbb3SSatish Balay } 18682d61bbb3SSatish Balay } 18692d61bbb3SSatish Balay } 1870d9ca1df4SBarry Smith ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr); 18713447b6efSHong Zhang ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr); 1872dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*a->nz*a->bs2);CHKERRQ(ierr); 18732d61bbb3SSatish Balay PetscFunctionReturn(0); 18742d61bbb3SSatish Balay } 18752d61bbb3SSatish Balay 18764a2ae208SSatish Balay #undef __FUNCT__ 18774a2ae208SSatish Balay #define __FUNCT__ "MatScale_SeqBAIJ" 1878f4df32b1SMatthew Knepley PetscErrorCode MatScale_SeqBAIJ(Mat inA,PetscScalar alpha) 18792d61bbb3SSatish Balay { 18802d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)inA->data; 1881690b6cddSBarry Smith PetscInt totalnz = a->bs2*a->nz; 1882f4df32b1SMatthew Knepley PetscScalar oalpha = alpha; 1883efee365bSSatish Balay PetscErrorCode ierr; 1884c5df96a5SBarry Smith PetscBLASInt one = 1,tnz; 18852d61bbb3SSatish Balay 18862d61bbb3SSatish Balay PetscFunctionBegin; 1887c5df96a5SBarry Smith ierr = PetscBLASIntCast(totalnz,&tnz);CHKERRQ(ierr); 18888b83055fSJed Brown PetscStackCallBLAS("BLASscal",BLASscal_(&tnz,&oalpha,a->a,&one)); 1889efee365bSSatish Balay ierr = PetscLogFlops(totalnz);CHKERRQ(ierr); 18902d61bbb3SSatish Balay PetscFunctionReturn(0); 18912d61bbb3SSatish Balay } 18922d61bbb3SSatish Balay 18934a2ae208SSatish Balay #undef __FUNCT__ 18944a2ae208SSatish Balay #define __FUNCT__ "MatNorm_SeqBAIJ" 1895dfbe8321SBarry Smith PetscErrorCode MatNorm_SeqBAIJ(Mat A,NormType type,PetscReal *norm) 18962d61bbb3SSatish Balay { 18978a62d963SHong Zhang PetscErrorCode ierr; 18982d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 18993f1db9ecSBarry Smith MatScalar *v = a->a; 1900329f5518SBarry Smith PetscReal sum = 0.0; 1901d0f46423SBarry Smith PetscInt i,j,k,bs=A->rmap->bs,nz=a->nz,bs2=a->bs2,k1; 19022d61bbb3SSatish Balay 19032d61bbb3SSatish Balay PetscFunctionBegin; 19042d61bbb3SSatish Balay if (type == NORM_FROBENIUS) { 19052d61bbb3SSatish Balay for (i=0; i< bs2*nz; i++) { 1906329f5518SBarry Smith sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 19072d61bbb3SSatish Balay } 19088f1a2a5eSBarry Smith *norm = PetscSqrtReal(sum); 19098a62d963SHong Zhang } else if (type == NORM_1) { /* maximum column sum */ 19108a62d963SHong Zhang PetscReal *tmp; 19118a62d963SHong Zhang PetscInt *bcol = a->j; 1912854ce69bSBarry Smith ierr = PetscCalloc1(A->cmap->n+1,&tmp);CHKERRQ(ierr); 19138a62d963SHong Zhang for (i=0; i<nz; i++) { 19148a62d963SHong Zhang for (j=0; j<bs; j++) { 19158a62d963SHong Zhang k1 = bs*(*bcol) + j; /* column index */ 19168a62d963SHong Zhang for (k=0; k<bs; k++) { 19178a62d963SHong Zhang tmp[k1] += PetscAbsScalar(*v); v++; 19188a62d963SHong Zhang } 19198a62d963SHong Zhang } 19208a62d963SHong Zhang bcol++; 19218a62d963SHong Zhang } 19228a62d963SHong Zhang *norm = 0.0; 1923d0f46423SBarry Smith for (j=0; j<A->cmap->n; j++) { 19248a62d963SHong Zhang if (tmp[j] > *norm) *norm = tmp[j]; 19258a62d963SHong Zhang } 19268a62d963SHong Zhang ierr = PetscFree(tmp);CHKERRQ(ierr); 1927596552b5SBarry Smith } else if (type == NORM_INFINITY) { /* maximum row sum */ 1928596552b5SBarry Smith *norm = 0.0; 1929596552b5SBarry Smith for (k=0; k<bs; k++) { 193074f84c7bSSatish Balay for (j=0; j<a->mbs; j++) { 1931596552b5SBarry Smith v = a->a + bs2*a->i[j] + k; 1932596552b5SBarry Smith sum = 0.0; 1933596552b5SBarry Smith for (i=0; i<a->i[j+1]-a->i[j]; i++) { 19340e90e235SBarry Smith for (k1=0; k1<bs; k1++) { 1935596552b5SBarry Smith sum += PetscAbsScalar(*v); 1936596552b5SBarry Smith v += bs; 19372d61bbb3SSatish Balay } 19380e90e235SBarry Smith } 1939596552b5SBarry Smith if (sum > *norm) *norm = sum; 1940596552b5SBarry Smith } 1941596552b5SBarry Smith } 1942e7e72b3dSBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for this norm yet"); 19432d61bbb3SSatish Balay PetscFunctionReturn(0); 19442d61bbb3SSatish Balay } 19452d61bbb3SSatish Balay 19462d61bbb3SSatish Balay 19474a2ae208SSatish Balay #undef __FUNCT__ 19484a2ae208SSatish Balay #define __FUNCT__ "MatEqual_SeqBAIJ" 1949ace3abfcSBarry Smith PetscErrorCode MatEqual_SeqBAIJ(Mat A,Mat B,PetscBool * flg) 19502d61bbb3SSatish Balay { 19512d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)B->data; 1952dfbe8321SBarry Smith PetscErrorCode ierr; 19532d61bbb3SSatish Balay 19542d61bbb3SSatish Balay PetscFunctionBegin; 19552d61bbb3SSatish Balay /* If the matrix/block dimensions are not equal, or no of nonzeros or shift */ 1956d0f46423SBarry Smith if ((A->rmap->N != B->rmap->N) || (A->cmap->n != B->cmap->n) || (A->rmap->bs != B->rmap->bs)|| (a->nz != b->nz)) { 1957273d9f13SBarry Smith *flg = PETSC_FALSE; 1958273d9f13SBarry Smith PetscFunctionReturn(0); 19592d61bbb3SSatish Balay } 19602d61bbb3SSatish Balay 19612d61bbb3SSatish Balay /* if the a->i are the same */ 1962690b6cddSBarry Smith ierr = PetscMemcmp(a->i,b->i,(a->mbs+1)*sizeof(PetscInt),flg);CHKERRQ(ierr); 196326fbe8dcSKarl Rupp if (!*flg) PetscFunctionReturn(0); 19642d61bbb3SSatish Balay 19652d61bbb3SSatish Balay /* if a->j are the same */ 1966690b6cddSBarry Smith ierr = PetscMemcmp(a->j,b->j,(a->nz)*sizeof(PetscInt),flg);CHKERRQ(ierr); 196726fbe8dcSKarl Rupp if (!*flg) PetscFunctionReturn(0); 196826fbe8dcSKarl Rupp 19692d61bbb3SSatish Balay /* if a->a are the same */ 1970d0f46423SBarry Smith ierr = PetscMemcmp(a->a,b->a,(a->nz)*(A->rmap->bs)*(B->rmap->bs)*sizeof(PetscScalar),flg);CHKERRQ(ierr); 19712d61bbb3SSatish Balay PetscFunctionReturn(0); 19722d61bbb3SSatish Balay 19732d61bbb3SSatish Balay } 19742d61bbb3SSatish Balay 19754a2ae208SSatish Balay #undef __FUNCT__ 19764a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_SeqBAIJ" 1977dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_SeqBAIJ(Mat A,Vec v) 19782d61bbb3SSatish Balay { 19792d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1980dfbe8321SBarry Smith PetscErrorCode ierr; 1981690b6cddSBarry Smith PetscInt i,j,k,n,row,bs,*ai,*aj,ambs,bs2; 198287828ca2SBarry Smith PetscScalar *x,zero = 0.0; 19833f1db9ecSBarry Smith MatScalar *aa,*aa_j; 19842d61bbb3SSatish Balay 19852d61bbb3SSatish Balay PetscFunctionBegin; 1986e32f2f54SBarry Smith if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 1987d0f46423SBarry Smith bs = A->rmap->bs; 19882d61bbb3SSatish Balay aa = a->a; 19892d61bbb3SSatish Balay ai = a->i; 19902d61bbb3SSatish Balay aj = a->j; 19912d61bbb3SSatish Balay ambs = a->mbs; 19922d61bbb3SSatish Balay bs2 = a->bs2; 19932d61bbb3SSatish Balay 19942dcb1b2aSMatthew Knepley ierr = VecSet(v,zero);CHKERRQ(ierr); 19951ebc52fbSHong Zhang ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1996e1311b90SBarry Smith ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr); 1997e32f2f54SBarry Smith if (n != A->rmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector"); 19982d61bbb3SSatish Balay for (i=0; i<ambs; i++) { 19992d61bbb3SSatish Balay for (j=ai[i]; j<ai[i+1]; j++) { 20002d61bbb3SSatish Balay if (aj[j] == i) { 20012d61bbb3SSatish Balay row = i*bs; 20022d61bbb3SSatish Balay aa_j = aa+j*bs2; 20032d61bbb3SSatish Balay for (k=0; k<bs2; k+=(bs+1),row++) x[row] = aa_j[k]; 20042d61bbb3SSatish Balay break; 20052d61bbb3SSatish Balay } 20062d61bbb3SSatish Balay } 20072d61bbb3SSatish Balay } 20081ebc52fbSHong Zhang ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 20092d61bbb3SSatish Balay PetscFunctionReturn(0); 20102d61bbb3SSatish Balay } 20112d61bbb3SSatish Balay 20124a2ae208SSatish Balay #undef __FUNCT__ 20134a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_SeqBAIJ" 2014dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_SeqBAIJ(Mat A,Vec ll,Vec rr) 20152d61bbb3SSatish Balay { 20162d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 201753ef36baSBarry Smith const PetscScalar *l,*r,*li,*ri; 201853ef36baSBarry Smith PetscScalar x; 20193f1db9ecSBarry Smith MatScalar *aa, *v; 2020dfbe8321SBarry Smith PetscErrorCode ierr; 202153ef36baSBarry Smith PetscInt i,j,k,lm,rn,M,m,n,mbs,tmp,bs,bs2,iai; 202253ef36baSBarry Smith const PetscInt *ai,*aj; 20232d61bbb3SSatish Balay 20242d61bbb3SSatish Balay PetscFunctionBegin; 20252d61bbb3SSatish Balay ai = a->i; 20262d61bbb3SSatish Balay aj = a->j; 20272d61bbb3SSatish Balay aa = a->a; 2028d0f46423SBarry Smith m = A->rmap->n; 2029d0f46423SBarry Smith n = A->cmap->n; 2030d0f46423SBarry Smith bs = A->rmap->bs; 20312d61bbb3SSatish Balay mbs = a->mbs; 20322d61bbb3SSatish Balay bs2 = a->bs2; 20332d61bbb3SSatish Balay if (ll) { 203453ef36baSBarry Smith ierr = VecGetArrayRead(ll,&l);CHKERRQ(ierr); 2035e1311b90SBarry Smith ierr = VecGetLocalSize(ll,&lm);CHKERRQ(ierr); 2036e32f2f54SBarry Smith if (lm != m) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Left scaling vector wrong length"); 20372d61bbb3SSatish Balay for (i=0; i<mbs; i++) { /* for each block row */ 20382d61bbb3SSatish Balay M = ai[i+1] - ai[i]; 20392d61bbb3SSatish Balay li = l + i*bs; 20402d61bbb3SSatish Balay v = aa + bs2*ai[i]; 20412d61bbb3SSatish Balay for (j=0; j<M; j++) { /* for each block */ 20422d61bbb3SSatish Balay for (k=0; k<bs2; k++) { 20432d61bbb3SSatish Balay (*v++) *= li[k%bs]; 20442d61bbb3SSatish Balay } 20452d61bbb3SSatish Balay } 20462d61bbb3SSatish Balay } 204753ef36baSBarry Smith ierr = VecRestoreArrayRead(ll,&l);CHKERRQ(ierr); 2048efee365bSSatish Balay ierr = PetscLogFlops(a->nz);CHKERRQ(ierr); 20492d61bbb3SSatish Balay } 20502d61bbb3SSatish Balay 20512d61bbb3SSatish Balay if (rr) { 205253ef36baSBarry Smith ierr = VecGetArrayRead(rr,&r);CHKERRQ(ierr); 2053e1311b90SBarry Smith ierr = VecGetLocalSize(rr,&rn);CHKERRQ(ierr); 2054e32f2f54SBarry Smith if (rn != n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Right scaling vector wrong length"); 20552d61bbb3SSatish Balay for (i=0; i<mbs; i++) { /* for each block row */ 205653ef36baSBarry Smith iai = ai[i]; 205753ef36baSBarry Smith M = ai[i+1] - iai; 205853ef36baSBarry Smith v = aa + bs2*iai; 20592d61bbb3SSatish Balay for (j=0; j<M; j++) { /* for each block */ 206053ef36baSBarry Smith ri = r + bs*aj[iai+j]; 20612d61bbb3SSatish Balay for (k=0; k<bs; k++) { 20622d61bbb3SSatish Balay x = ri[k]; 206353ef36baSBarry Smith for (tmp=0; tmp<bs; tmp++) v[tmp] *= x; 206453ef36baSBarry Smith v += bs; 20652d61bbb3SSatish Balay } 20662d61bbb3SSatish Balay } 20672d61bbb3SSatish Balay } 206853ef36baSBarry Smith ierr = VecRestoreArrayRead(rr,&r);CHKERRQ(ierr); 2069efee365bSSatish Balay ierr = PetscLogFlops(a->nz);CHKERRQ(ierr); 20702d61bbb3SSatish Balay } 20712d61bbb3SSatish Balay PetscFunctionReturn(0); 20722d61bbb3SSatish Balay } 20732d61bbb3SSatish Balay 20742d61bbb3SSatish Balay 20754a2ae208SSatish Balay #undef __FUNCT__ 20764a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_SeqBAIJ" 2077dfbe8321SBarry Smith PetscErrorCode MatGetInfo_SeqBAIJ(Mat A,MatInfoType flag,MatInfo *info) 20782d61bbb3SSatish Balay { 20792d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 20802d61bbb3SSatish Balay 20812d61bbb3SSatish Balay PetscFunctionBegin; 20822d61bbb3SSatish Balay info->block_size = a->bs2; 2083ceed8ce5SJed Brown info->nz_allocated = a->bs2*a->maxnz; 20842d61bbb3SSatish Balay info->nz_used = a->bs2*a->nz; 20852d61bbb3SSatish Balay info->nz_unneeded = (double)(info->nz_allocated - info->nz_used); 20862d61bbb3SSatish Balay info->assemblies = A->num_ass; 20878e58a170SBarry Smith info->mallocs = A->info.mallocs; 20887adad957SLisandro Dalcin info->memory = ((PetscObject)A)->mem; 2089d5f3da31SBarry Smith if (A->factortype) { 20902d61bbb3SSatish Balay info->fill_ratio_given = A->info.fill_ratio_given; 20912d61bbb3SSatish Balay info->fill_ratio_needed = A->info.fill_ratio_needed; 20922d61bbb3SSatish Balay info->factor_mallocs = A->info.factor_mallocs; 20932d61bbb3SSatish Balay } else { 20942d61bbb3SSatish Balay info->fill_ratio_given = 0; 20952d61bbb3SSatish Balay info->fill_ratio_needed = 0; 20962d61bbb3SSatish Balay info->factor_mallocs = 0; 20972d61bbb3SSatish Balay } 20982d61bbb3SSatish Balay PetscFunctionReturn(0); 20992d61bbb3SSatish Balay } 21002d61bbb3SSatish Balay 21014a2ae208SSatish Balay #undef __FUNCT__ 21024a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_SeqBAIJ" 2103dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_SeqBAIJ(Mat A) 21042d61bbb3SSatish Balay { 21052d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 2106dfbe8321SBarry Smith PetscErrorCode ierr; 21072d61bbb3SSatish Balay 21082d61bbb3SSatish Balay PetscFunctionBegin; 2109549d3d68SSatish Balay ierr = PetscMemzero(a->a,a->bs2*a->i[a->mbs]*sizeof(MatScalar));CHKERRQ(ierr); 21102d61bbb3SSatish Balay PetscFunctionReturn(0); 21112d61bbb3SSatish Balay } 2112