1cac129eeSSatish Balay 2c6db04a5SJed Brown #include <../src/mat/impls/baij/seq/baij.h> 3af0996ceSBarry Smith #include <petsc/private/kernels/blockinvert.h> 4c6db04a5SJed Brown #include <petscbt.h> 5c6db04a5SJed Brown #include <petscblaslapack.h> 6cac129eeSSatish Balay 74a2ae208SSatish Balay #undef __FUNCT__ 84a2ae208SSatish Balay #define __FUNCT__ "MatIncreaseOverlap_SeqBAIJ" 9690b6cddSBarry Smith PetscErrorCode MatIncreaseOverlap_SeqBAIJ(Mat A,PetscInt is_max,IS is[],PetscInt ov) 10a3192f15SSatish Balay { 11a3192f15SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 126849ba73SBarry Smith PetscErrorCode ierr; 135d0c19d7SBarry Smith PetscInt row,i,j,k,l,m,n,*nidx,isz,val,ival; 145d0c19d7SBarry Smith const PetscInt *idx; 15690b6cddSBarry Smith PetscInt start,end,*ai,*aj,bs,*nidx2; 16f1af5d2fSBarry Smith PetscBT table; 17a3192f15SSatish Balay 183a40ed3dSBarry Smith PetscFunctionBegin; 19a3192f15SSatish Balay m = a->mbs; 20a3192f15SSatish Balay ai = a->i; 21a3192f15SSatish Balay aj = a->j; 22d0f46423SBarry Smith bs = A->rmap->bs; 23a3192f15SSatish Balay 24e32f2f54SBarry Smith if (ov < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap specified"); 25a3192f15SSatish Balay 2653b8de81SBarry Smith ierr = PetscBTCreate(m,&table);CHKERRQ(ierr); 27854ce69bSBarry Smith ierr = PetscMalloc1(m+1,&nidx);CHKERRQ(ierr); 28854ce69bSBarry Smith ierr = PetscMalloc1(A->rmap->N+1,&nidx2);CHKERRQ(ierr); 29a3192f15SSatish Balay 30a3192f15SSatish Balay for (i=0; i<is_max; i++) { 31a3192f15SSatish Balay /* Initialise the two local arrays */ 32a3192f15SSatish Balay isz = 0; 336831982aSBarry Smith ierr = PetscBTMemzero(m,table);CHKERRQ(ierr); 34a3192f15SSatish Balay 35a3192f15SSatish Balay /* Extract the indices, assume there can be duplicate entries */ 36a3192f15SSatish Balay ierr = ISGetIndices(is[i],&idx);CHKERRQ(ierr); 37b9b97703SBarry Smith ierr = ISGetLocalSize(is[i],&n);CHKERRQ(ierr); 38a3192f15SSatish Balay 39a3192f15SSatish Balay /* Enter these into the temp arrays i.e mark table[row], enter row into new index */ 40a3192f15SSatish Balay for (j=0; j<n; ++j) { 41218c64b6SSatish Balay ival = idx[j]/bs; /* convert the indices into block indices */ 42e32f2f54SBarry Smith if (ival>=m) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"index greater than mat-dim"); 4326fbe8dcSKarl Rupp if (!PetscBTLookupSet(table,ival)) nidx[isz++] = ival; 44a3192f15SSatish Balay } 45a3192f15SSatish Balay ierr = ISRestoreIndices(is[i],&idx);CHKERRQ(ierr); 466bf464f9SBarry Smith ierr = ISDestroy(&is[i]);CHKERRQ(ierr); 47a3192f15SSatish Balay 48a3192f15SSatish Balay k = 0; 49a3192f15SSatish Balay for (j=0; j<ov; j++) { /* for each overlap*/ 50a3192f15SSatish Balay n = isz; 51a3192f15SSatish Balay for (; k<n; k++) { /* do only those rows in nidx[k], which are not done yet */ 52a3192f15SSatish Balay row = nidx[k]; 53a3192f15SSatish Balay start = ai[row]; 54a3192f15SSatish Balay end = ai[row+1]; 55a3192f15SSatish Balay for (l = start; l<end; l++) { 56a3192f15SSatish Balay val = aj[l]; 5726fbe8dcSKarl Rupp if (!PetscBTLookupSet(table,val)) nidx[isz++] = val; 58a3192f15SSatish Balay } 59a3192f15SSatish Balay } 60a3192f15SSatish Balay } 61218c64b6SSatish Balay /* expand the Index Set */ 62218c64b6SSatish Balay for (j=0; j<isz; j++) { 6326fbe8dcSKarl Rupp for (k=0; k<bs; k++) nidx2[j*bs+k] = nidx[j]*bs+k; 64218c64b6SSatish Balay } 6570b3c8c7SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,isz*bs,nidx2,PETSC_COPY_VALUES,is+i);CHKERRQ(ierr); 66a3192f15SSatish Balay } 6794bacf5dSBarry Smith ierr = PetscBTDestroy(&table);CHKERRQ(ierr); 68606d414cSSatish Balay ierr = PetscFree(nidx);CHKERRQ(ierr); 69606d414cSSatish Balay ierr = PetscFree(nidx2);CHKERRQ(ierr); 703a40ed3dSBarry Smith PetscFunctionReturn(0); 71a3192f15SSatish Balay } 721c351548SSatish Balay 734a2ae208SSatish Balay #undef __FUNCT__ 744a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqBAIJ_Private" 754aa3045dSJed Brown PetscErrorCode MatGetSubMatrix_SeqBAIJ_Private(Mat A,IS isrow,IS iscol,MatReuse scall,Mat *B) 76736121d4SSatish Balay { 77736121d4SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*c; 786849ba73SBarry Smith PetscErrorCode ierr; 79690b6cddSBarry Smith PetscInt *smap,i,k,kstart,kend,oldcols = a->nbs,*lens; 80690b6cddSBarry Smith PetscInt row,mat_i,*mat_j,tcol,*mat_ilen; 815d0c19d7SBarry Smith const PetscInt *irow,*icol; 825d0c19d7SBarry Smith PetscInt nrows,ncols,*ssmap,bs=A->rmap->bs,bs2=a->bs2; 83690b6cddSBarry Smith PetscInt *aj = a->j,*ai = a->i; 843f1db9ecSBarry Smith MatScalar *mat_a; 85736121d4SSatish Balay Mat C; 866041f1b1SToby Isaac PetscBool flag; 87736121d4SSatish Balay 883a40ed3dSBarry Smith PetscFunctionBegin; 89736121d4SSatish Balay ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr); 90218c64b6SSatish Balay ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr); 91b9b97703SBarry Smith ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr); 92b9b97703SBarry Smith ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr); 93736121d4SSatish Balay 94854ce69bSBarry Smith ierr = PetscCalloc1(1+oldcols,&smap);CHKERRQ(ierr); 95736121d4SSatish Balay ssmap = smap; 96854ce69bSBarry Smith ierr = PetscMalloc1(1+nrows,&lens);CHKERRQ(ierr); 97736121d4SSatish Balay for (i=0; i<ncols; i++) smap[icol[i]] = i+1; 98736121d4SSatish Balay /* determine lens of each row */ 99736121d4SSatish Balay for (i=0; i<nrows; i++) { 100736121d4SSatish Balay kstart = ai[irow[i]]; 101736121d4SSatish Balay kend = kstart + a->ilen[irow[i]]; 102736121d4SSatish Balay lens[i] = 0; 103736121d4SSatish Balay for (k=kstart; k<kend; k++) { 10426fbe8dcSKarl Rupp if (ssmap[aj[k]]) lens[i]++; 105736121d4SSatish Balay } 106736121d4SSatish Balay } 107736121d4SSatish Balay /* Create and fill new matrix */ 108736121d4SSatish Balay if (scall == MAT_REUSE_MATRIX) { 109736121d4SSatish Balay c = (Mat_SeqBAIJ*)((*B)->data); 110736121d4SSatish Balay 111e32f2f54SBarry Smith if (c->mbs!=nrows || c->nbs!=ncols || (*B)->rmap->bs!=bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Submatrix wrong size"); 112690b6cddSBarry Smith ierr = PetscMemcmp(c->ilen,lens,c->mbs *sizeof(PetscInt),&flag);CHKERRQ(ierr); 113e7e72b3dSBarry Smith if (!flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong no of nonzeros"); 114690b6cddSBarry Smith ierr = PetscMemzero(c->ilen,c->mbs*sizeof(PetscInt));CHKERRQ(ierr); 115736121d4SSatish Balay C = *B; 1163a40ed3dSBarry Smith } else { 117ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)A),&C);CHKERRQ(ierr); 118f69a0ea3SMatthew Knepley ierr = MatSetSizes(C,nrows*bs,ncols*bs,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 1197adad957SLisandro Dalcin ierr = MatSetType(C,((PetscObject)A)->type_name);CHKERRQ(ierr); 120*367daffbSBarry Smith ierr = MatSeqBAIJSetPreallocation(C,bs,0,lens);CHKERRQ(ierr); 121736121d4SSatish Balay } 122736121d4SSatish Balay c = (Mat_SeqBAIJ*)(C->data); 123736121d4SSatish Balay for (i=0; i<nrows; i++) { 124736121d4SSatish Balay row = irow[i]; 125736121d4SSatish Balay kstart = ai[row]; 126736121d4SSatish Balay kend = kstart + a->ilen[row]; 127736121d4SSatish Balay mat_i = c->i[i]; 128736121d4SSatish Balay mat_j = c->j + mat_i; 129218c64b6SSatish Balay mat_a = c->a + mat_i*bs2; 130736121d4SSatish Balay mat_ilen = c->ilen + i; 131736121d4SSatish Balay for (k=kstart; k<kend; k++) { 132736121d4SSatish Balay if ((tcol=ssmap[a->j[k]])) { 133736121d4SSatish Balay *mat_j++ = tcol - 1; 134549d3d68SSatish Balay ierr = PetscMemcpy(mat_a,a->a+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 135549d3d68SSatish Balay mat_a += bs2; 136736121d4SSatish Balay (*mat_ilen)++; 137736121d4SSatish Balay } 138736121d4SSatish Balay } 139736121d4SSatish Balay } 140cdc6f3adSToby Isaac /* sort */ 141cdc6f3adSToby Isaac { 142cdc6f3adSToby Isaac MatScalar *work; 143cdc6f3adSToby Isaac ierr = PetscMalloc1(bs2,&work);CHKERRQ(ierr); 144cdc6f3adSToby Isaac for (i=0; i<nrows; i++) { 145cdc6f3adSToby Isaac PetscInt ilen; 146cdc6f3adSToby Isaac mat_i = c->i[i]; 147cdc6f3adSToby Isaac mat_j = c->j + mat_i; 148cdc6f3adSToby Isaac mat_a = c->a + mat_i*bs2; 149cdc6f3adSToby Isaac ilen = c->ilen[i]; 150cdc6f3adSToby Isaac ierr = PetscSortIntWithDataArray(ilen,mat_j,mat_a,bs2*sizeof(MatScalar),work);CHKERRQ(ierr); 151cdc6f3adSToby Isaac } 152cdc6f3adSToby Isaac ierr = PetscFree(work);CHKERRQ(ierr); 153cdc6f3adSToby Isaac } 154218c64b6SSatish Balay 155736121d4SSatish Balay /* Free work space */ 156736121d4SSatish Balay ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr); 157606d414cSSatish Balay ierr = PetscFree(smap);CHKERRQ(ierr); 158606d414cSSatish Balay ierr = PetscFree(lens);CHKERRQ(ierr); 1596d4a8577SBarry Smith ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1606d4a8577SBarry Smith ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 161736121d4SSatish Balay 162736121d4SSatish Balay ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr); 163736121d4SSatish Balay *B = C; 1643a40ed3dSBarry Smith PetscFunctionReturn(0); 165736121d4SSatish Balay } 166736121d4SSatish Balay 1674a2ae208SSatish Balay #undef __FUNCT__ 1684a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqBAIJ" 1694aa3045dSJed Brown PetscErrorCode MatGetSubMatrix_SeqBAIJ(Mat A,IS isrow,IS iscol,MatReuse scall,Mat *B) 170218c64b6SSatish Balay { 171218c64b6SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 172218c64b6SSatish Balay IS is1,is2; 1736849ba73SBarry Smith PetscErrorCode ierr; 174afebec48SHong Zhang PetscInt *vary,*iary,nrows,ncols,i,bs=A->rmap->bs,count,maxmnbs,j; 1755d0c19d7SBarry Smith const PetscInt *irow,*icol; 176218c64b6SSatish Balay 1773a40ed3dSBarry Smith PetscFunctionBegin; 178218c64b6SSatish Balay ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr); 179218c64b6SSatish Balay ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr); 180b9b97703SBarry Smith ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr); 181b9b97703SBarry Smith ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr); 182218c64b6SSatish Balay 183218c64b6SSatish Balay /* Verify if the indices corespond to each element in a block 184218c64b6SSatish Balay and form the IS with compressed IS */ 185f8ecb639SStefano Zampini maxmnbs = PetscMax(a->mbs,a->nbs); 186f8ecb639SStefano Zampini ierr = PetscMalloc2(maxmnbs,&vary,maxmnbs,&iary);CHKERRQ(ierr); 187fca92195SBarry Smith ierr = PetscMemzero(vary,a->mbs*sizeof(PetscInt));CHKERRQ(ierr); 188218c64b6SSatish Balay for (i=0; i<nrows; i++) vary[irow[i]/bs]++; 189218c64b6SSatish Balay for (i=0; i<a->mbs; i++) { 190e32f2f54SBarry Smith if (vary[i]!=0 && vary[i]!=bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Index set does not match blocks"); 1916041f1b1SToby Isaac } 1926041f1b1SToby Isaac count = 0; 1936041f1b1SToby Isaac for (i=0; i<nrows; i++) { 194afebec48SHong Zhang j = irow[i] / bs; 1956041f1b1SToby Isaac if ((vary[j]--)==bs) iary[count++] = j; 196218c64b6SSatish Balay } 19770b3c8c7SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,count,iary,PETSC_COPY_VALUES,&is1);CHKERRQ(ierr); 198218c64b6SSatish Balay 199f8ecb639SStefano Zampini ierr = PetscMemzero(vary,(a->nbs)*sizeof(PetscInt));CHKERRQ(ierr); 200218c64b6SSatish Balay for (i=0; i<ncols; i++) vary[icol[i]/bs]++; 201f8ecb639SStefano Zampini for (i=0; i<a->nbs; i++) { 202e32f2f54SBarry Smith if (vary[i]!=0 && vary[i]!=bs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Internal error in PETSc"); 203218c64b6SSatish Balay } 2046041f1b1SToby Isaac count = 0; 2056041f1b1SToby Isaac for (i=0; i<ncols; i++) { 206afebec48SHong Zhang j = icol[i] / bs; 2076041f1b1SToby Isaac if ((vary[j]--)==bs) iary[count++] = j; 2086041f1b1SToby Isaac } 20970b3c8c7SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,count,iary,PETSC_COPY_VALUES,&is2);CHKERRQ(ierr); 210218c64b6SSatish Balay ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr); 211218c64b6SSatish Balay ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr); 212fca92195SBarry Smith ierr = PetscFree2(vary,iary);CHKERRQ(ierr); 213218c64b6SSatish Balay 2144aa3045dSJed Brown ierr = MatGetSubMatrix_SeqBAIJ_Private(A,is1,is2,scall,B);CHKERRQ(ierr); 2156bf464f9SBarry Smith ierr = ISDestroy(&is1);CHKERRQ(ierr); 2166bf464f9SBarry Smith ierr = ISDestroy(&is2);CHKERRQ(ierr); 2173a40ed3dSBarry Smith PetscFunctionReturn(0); 218218c64b6SSatish Balay } 219218c64b6SSatish Balay 2204a2ae208SSatish Balay #undef __FUNCT__ 2214a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_SeqBAIJ" 222690b6cddSBarry Smith PetscErrorCode MatGetSubMatrices_SeqBAIJ(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[]) 223736121d4SSatish Balay { 2246849ba73SBarry Smith PetscErrorCode ierr; 225690b6cddSBarry Smith PetscInt i; 226736121d4SSatish Balay 2273a40ed3dSBarry Smith PetscFunctionBegin; 228736121d4SSatish Balay if (scall == MAT_INITIAL_MATRIX) { 229854ce69bSBarry Smith ierr = PetscMalloc1(n+1,B);CHKERRQ(ierr); 230736121d4SSatish Balay } 231736121d4SSatish Balay 232736121d4SSatish Balay for (i=0; i<n; i++) { 2334aa3045dSJed Brown ierr = MatGetSubMatrix_SeqBAIJ(A,irow[i],icol[i],scall,&(*B)[i]);CHKERRQ(ierr); 234736121d4SSatish Balay } 2353a40ed3dSBarry Smith PetscFunctionReturn(0); 236736121d4SSatish Balay } 237218c64b6SSatish Balay 238218c64b6SSatish Balay 2392d61bbb3SSatish Balay /* -------------------------------------------------------*/ 2402d61bbb3SSatish Balay /* Should check that shapes of vectors and matrices match */ 2412d61bbb3SSatish Balay /* -------------------------------------------------------*/ 2422d61bbb3SSatish Balay 2434a2ae208SSatish Balay #undef __FUNCT__ 2444a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_1" 245dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_1(Mat A,Vec xx,Vec zz) 2462d61bbb3SSatish Balay { 2472d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 248d9fead3dSBarry Smith PetscScalar *z,sum; 249d9fead3dSBarry Smith const PetscScalar *x; 250d9fead3dSBarry Smith const MatScalar *v; 2516849ba73SBarry Smith PetscErrorCode ierr; 2527c565772SBarry Smith PetscInt mbs,i,n; 2530298fd71SBarry Smith const PetscInt *idx,*ii,*ridx=NULL; 254ace3abfcSBarry Smith PetscBool usecprow=a->compressedrow.use; 2552d61bbb3SSatish Balay 2562d61bbb3SSatish Balay PetscFunctionBegin; 2573649974fSBarry Smith ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr); 2581ebc52fbSHong Zhang ierr = VecGetArray(zz,&z);CHKERRQ(ierr); 2592d61bbb3SSatish Balay 26026e093fcSHong Zhang if (usecprow) { 26126e093fcSHong Zhang mbs = a->compressedrow.nrows; 26226e093fcSHong Zhang ii = a->compressedrow.i; 2637b2bb3b9SHong Zhang ridx = a->compressedrow.rindex; 264a1c3900fSBarry Smith ierr = PetscMemzero(z,mbs*sizeof(PetscScalar));CHKERRQ(ierr); 26526e093fcSHong Zhang } else { 26626e093fcSHong Zhang mbs = a->mbs; 2672d61bbb3SSatish Balay ii = a->i; 26826e093fcSHong Zhang } 2692d61bbb3SSatish Balay 2702d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 271ee54c7eeSHong Zhang n = ii[1] - ii[0]; 272ee54c7eeSHong Zhang v = a->a + ii[0]; 273ee54c7eeSHong Zhang idx = a->j + ii[0]; 274ee54c7eeSHong Zhang ii++; 275444d8c10SJed Brown PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA); /* Indices for the next row (assumes same size as this one) */ 276444d8c10SJed Brown PetscPrefetchBlock(v+1*n,1*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */ 2772d61bbb3SSatish Balay sum = 0.0; 2782162cab8SBarry Smith PetscSparseDensePlusDot(sum,x,v,idx,n); 27926e093fcSHong Zhang if (usecprow) { 2807b2bb3b9SHong Zhang z[ridx[i]] = sum; 28126e093fcSHong Zhang } else { 2822d61bbb3SSatish Balay z[i] = sum; 2832d61bbb3SSatish Balay } 28426e093fcSHong Zhang } 2853649974fSBarry Smith ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr); 2861ebc52fbSHong Zhang ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr); 2877c565772SBarry Smith ierr = PetscLogFlops(2.0*a->nz - a->nonzerorowcnt);CHKERRQ(ierr); 2882d61bbb3SSatish Balay PetscFunctionReturn(0); 2892d61bbb3SSatish Balay } 2902d61bbb3SSatish Balay 2914a2ae208SSatish Balay #undef __FUNCT__ 2924a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_2" 293dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_2(Mat A,Vec xx,Vec zz) 2942d61bbb3SSatish Balay { 2952d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 296d9fead3dSBarry Smith PetscScalar *z = 0,sum1,sum2,*zarray; 297d9fead3dSBarry Smith const PetscScalar *x,*xb; 29887828ca2SBarry Smith PetscScalar x1,x2; 299d9fead3dSBarry Smith const MatScalar *v; 300dfbe8321SBarry Smith PetscErrorCode ierr; 3017c565772SBarry Smith PetscInt mbs,i,*idx,*ii,j,n,*ridx=NULL; 302ace3abfcSBarry Smith PetscBool usecprow=a->compressedrow.use; 3032d61bbb3SSatish Balay 3042d61bbb3SSatish Balay PetscFunctionBegin; 3053649974fSBarry Smith ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr); 30626e093fcSHong Zhang ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr); 3072d61bbb3SSatish Balay 3082d61bbb3SSatish Balay idx = a->j; 3092d61bbb3SSatish Balay v = a->a; 31026e093fcSHong Zhang if (usecprow) { 31126e093fcSHong Zhang mbs = a->compressedrow.nrows; 31226e093fcSHong Zhang ii = a->compressedrow.i; 3137b2bb3b9SHong Zhang ridx = a->compressedrow.rindex; 31426e093fcSHong Zhang } else { 31526e093fcSHong Zhang mbs = a->mbs; 3162d61bbb3SSatish Balay ii = a->i; 31726e093fcSHong Zhang z = zarray; 31826e093fcSHong Zhang } 3192d61bbb3SSatish Balay 3202d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 3212d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 3222d61bbb3SSatish Balay sum1 = 0.0; sum2 = 0.0; 323444d8c10SJed Brown PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA); /* Indices for the next row (assumes same size as this one) */ 324444d8c10SJed Brown PetscPrefetchBlock(v+4*n,4*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */ 3252d61bbb3SSatish Balay for (j=0; j<n; j++) { 3262d61bbb3SSatish Balay xb = x + 2*(*idx++); x1 = xb[0]; x2 = xb[1]; 3272d61bbb3SSatish Balay sum1 += v[0]*x1 + v[2]*x2; 3282d61bbb3SSatish Balay sum2 += v[1]*x1 + v[3]*x2; 3292d61bbb3SSatish Balay v += 4; 3302d61bbb3SSatish Balay } 3317b2bb3b9SHong Zhang if (usecprow) z = zarray + 2*ridx[i]; 3322d61bbb3SSatish Balay z[0] = sum1; z[1] = sum2; 33326e093fcSHong Zhang if (!usecprow) z += 2; 3342d61bbb3SSatish Balay } 3353649974fSBarry Smith ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr); 33626e093fcSHong Zhang ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr); 3377c565772SBarry Smith ierr = PetscLogFlops(8.0*a->nz - 2.0*a->nonzerorowcnt);CHKERRQ(ierr); 3382d61bbb3SSatish Balay PetscFunctionReturn(0); 3392d61bbb3SSatish Balay } 3402d61bbb3SSatish Balay 3414a2ae208SSatish Balay #undef __FUNCT__ 3424a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_3" 343dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_3(Mat A,Vec xx,Vec zz) 3442d61bbb3SSatish Balay { 3452d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 346d9fead3dSBarry Smith PetscScalar *z = 0,sum1,sum2,sum3,x1,x2,x3,*zarray; 347d9fead3dSBarry Smith const PetscScalar *x,*xb; 348d9fead3dSBarry Smith const MatScalar *v; 349dfbe8321SBarry Smith PetscErrorCode ierr; 3507c565772SBarry Smith PetscInt mbs,i,*idx,*ii,j,n,*ridx=NULL; 351ace3abfcSBarry Smith PetscBool usecprow=a->compressedrow.use; 35226e093fcSHong Zhang 353b6410449SSatish Balay #if defined(PETSC_HAVE_PRAGMA_DISJOINT) 354fee21e36SBarry Smith #pragma disjoint(*v,*z,*xb) 355fee21e36SBarry Smith #endif 356fee21e36SBarry Smith 3572d61bbb3SSatish Balay PetscFunctionBegin; 3583649974fSBarry Smith ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr); 35926e093fcSHong Zhang ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr); 3602d61bbb3SSatish Balay 3612d61bbb3SSatish Balay idx = a->j; 3622d61bbb3SSatish Balay v = a->a; 36326e093fcSHong Zhang if (usecprow) { 36426e093fcSHong Zhang mbs = a->compressedrow.nrows; 36526e093fcSHong Zhang ii = a->compressedrow.i; 3667b2bb3b9SHong Zhang ridx = a->compressedrow.rindex; 36726e093fcSHong Zhang } else { 36826e093fcSHong Zhang mbs = a->mbs; 3692d61bbb3SSatish Balay ii = a->i; 37026e093fcSHong Zhang z = zarray; 37126e093fcSHong Zhang } 3722d61bbb3SSatish Balay 3732d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 3742d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 3752d61bbb3SSatish Balay sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; 376444d8c10SJed Brown PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA); /* Indices for the next row (assumes same size as this one) */ 377444d8c10SJed Brown PetscPrefetchBlock(v+9*n,9*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */ 3782d61bbb3SSatish Balay for (j=0; j<n; j++) { 37926fbe8dcSKarl Rupp xb = x + 3*(*idx++); 38026fbe8dcSKarl Rupp x1 = xb[0]; 38126fbe8dcSKarl Rupp x2 = xb[1]; 38226fbe8dcSKarl Rupp x3 = xb[2]; 38326fbe8dcSKarl Rupp 3842d61bbb3SSatish Balay sum1 += v[0]*x1 + v[3]*x2 + v[6]*x3; 3852d61bbb3SSatish Balay sum2 += v[1]*x1 + v[4]*x2 + v[7]*x3; 3862d61bbb3SSatish Balay sum3 += v[2]*x1 + v[5]*x2 + v[8]*x3; 3872d61bbb3SSatish Balay v += 9; 3882d61bbb3SSatish Balay } 3897b2bb3b9SHong Zhang if (usecprow) z = zarray + 3*ridx[i]; 3902d61bbb3SSatish Balay z[0] = sum1; z[1] = sum2; z[2] = sum3; 39126e093fcSHong Zhang if (!usecprow) z += 3; 3922d61bbb3SSatish Balay } 3933649974fSBarry Smith ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr); 39426e093fcSHong Zhang ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr); 3957c565772SBarry Smith ierr = PetscLogFlops(18.0*a->nz - 3.0*a->nonzerorowcnt);CHKERRQ(ierr); 3962d61bbb3SSatish Balay PetscFunctionReturn(0); 3972d61bbb3SSatish Balay } 3982d61bbb3SSatish Balay 3994a2ae208SSatish Balay #undef __FUNCT__ 4004a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_4" 401dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_4(Mat A,Vec xx,Vec zz) 4022d61bbb3SSatish Balay { 4032d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 404d9fead3dSBarry Smith PetscScalar *z = 0,sum1,sum2,sum3,sum4,x1,x2,x3,x4,*zarray; 405d9fead3dSBarry Smith const PetscScalar *x,*xb; 406d9fead3dSBarry Smith const MatScalar *v; 407dfbe8321SBarry Smith PetscErrorCode ierr; 4087c565772SBarry Smith PetscInt mbs,i,*idx,*ii,j,n,*ridx=NULL; 409ace3abfcSBarry Smith PetscBool usecprow=a->compressedrow.use; 4102d61bbb3SSatish Balay 4112d61bbb3SSatish Balay PetscFunctionBegin; 4123649974fSBarry Smith ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr); 41326e093fcSHong Zhang ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr); 4142d61bbb3SSatish Balay 4152d61bbb3SSatish Balay idx = a->j; 4162d61bbb3SSatish Balay v = a->a; 41726e093fcSHong Zhang if (usecprow) { 41826e093fcSHong Zhang mbs = a->compressedrow.nrows; 41926e093fcSHong Zhang ii = a->compressedrow.i; 4207b2bb3b9SHong Zhang ridx = a->compressedrow.rindex; 42126e093fcSHong Zhang } else { 42226e093fcSHong Zhang mbs = a->mbs; 4232d61bbb3SSatish Balay ii = a->i; 42426e093fcSHong Zhang z = zarray; 42526e093fcSHong Zhang } 4262d61bbb3SSatish Balay 4272d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 42826fbe8dcSKarl Rupp n = ii[1] - ii[0]; 42926fbe8dcSKarl Rupp ii++; 43026fbe8dcSKarl Rupp sum1 = 0.0; 43126fbe8dcSKarl Rupp sum2 = 0.0; 43226fbe8dcSKarl Rupp sum3 = 0.0; 43326fbe8dcSKarl Rupp sum4 = 0.0; 43426fbe8dcSKarl Rupp 435444d8c10SJed Brown PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA); /* Indices for the next row (assumes same size as this one) */ 436444d8c10SJed Brown PetscPrefetchBlock(v+16*n,16*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */ 4372d61bbb3SSatish Balay for (j=0; j<n; j++) { 4382d61bbb3SSatish Balay xb = x + 4*(*idx++); 4392d61bbb3SSatish Balay x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; 4402d61bbb3SSatish Balay sum1 += v[0]*x1 + v[4]*x2 + v[8]*x3 + v[12]*x4; 4412d61bbb3SSatish Balay sum2 += v[1]*x1 + v[5]*x2 + v[9]*x3 + v[13]*x4; 4422d61bbb3SSatish Balay sum3 += v[2]*x1 + v[6]*x2 + v[10]*x3 + v[14]*x4; 4432d61bbb3SSatish Balay sum4 += v[3]*x1 + v[7]*x2 + v[11]*x3 + v[15]*x4; 4442d61bbb3SSatish Balay v += 16; 4452d61bbb3SSatish Balay } 4467b2bb3b9SHong Zhang if (usecprow) z = zarray + 4*ridx[i]; 4472d61bbb3SSatish Balay z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; 44826e093fcSHong Zhang if (!usecprow) z += 4; 4492d61bbb3SSatish Balay } 4503649974fSBarry Smith ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr); 45126e093fcSHong Zhang ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr); 4527c565772SBarry Smith ierr = PetscLogFlops(32.0*a->nz - 4.0*a->nonzerorowcnt);CHKERRQ(ierr); 4532d61bbb3SSatish Balay PetscFunctionReturn(0); 4542d61bbb3SSatish Balay } 4552d61bbb3SSatish Balay 4564a2ae208SSatish Balay #undef __FUNCT__ 4574a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_5" 458dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_5(Mat A,Vec xx,Vec zz) 4592d61bbb3SSatish Balay { 4602d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 461d9fead3dSBarry Smith PetscScalar sum1,sum2,sum3,sum4,sum5,x1,x2,x3,x4,x5,*z = 0,*zarray; 462d9fead3dSBarry Smith const PetscScalar *xb,*x; 463d9fead3dSBarry Smith const MatScalar *v; 464dfbe8321SBarry Smith PetscErrorCode ierr; 4650298fd71SBarry Smith const PetscInt *idx,*ii,*ridx=NULL; 4667c565772SBarry Smith PetscInt mbs,i,j,n; 467ace3abfcSBarry Smith PetscBool usecprow=a->compressedrow.use; 4682d61bbb3SSatish Balay 469433994e6SBarry Smith PetscFunctionBegin; 4703649974fSBarry Smith ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr); 47126e093fcSHong Zhang ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr); 4722d61bbb3SSatish Balay 4732d61bbb3SSatish Balay idx = a->j; 4742d61bbb3SSatish Balay v = a->a; 47526e093fcSHong Zhang if (usecprow) { 47626e093fcSHong Zhang mbs = a->compressedrow.nrows; 47726e093fcSHong Zhang ii = a->compressedrow.i; 4787b2bb3b9SHong Zhang ridx = a->compressedrow.rindex; 47926e093fcSHong Zhang } else { 48026e093fcSHong Zhang mbs = a->mbs; 4812d61bbb3SSatish Balay ii = a->i; 48226e093fcSHong Zhang z = zarray; 48326e093fcSHong Zhang } 4842d61bbb3SSatish Balay 4852d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 4862d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 4872d61bbb3SSatish Balay sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; 488444d8c10SJed Brown PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA); /* Indices for the next row (assumes same size as this one) */ 489444d8c10SJed Brown PetscPrefetchBlock(v+25*n,25*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */ 4902d61bbb3SSatish Balay for (j=0; j<n; j++) { 4912d61bbb3SSatish Balay xb = x + 5*(*idx++); 4922d61bbb3SSatish Balay x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; 4932d61bbb3SSatish Balay sum1 += v[0]*x1 + v[5]*x2 + v[10]*x3 + v[15]*x4 + v[20]*x5; 4942d61bbb3SSatish Balay sum2 += v[1]*x1 + v[6]*x2 + v[11]*x3 + v[16]*x4 + v[21]*x5; 4952d61bbb3SSatish Balay sum3 += v[2]*x1 + v[7]*x2 + v[12]*x3 + v[17]*x4 + v[22]*x5; 4962d61bbb3SSatish Balay sum4 += v[3]*x1 + v[8]*x2 + v[13]*x3 + v[18]*x4 + v[23]*x5; 4972d61bbb3SSatish Balay sum5 += v[4]*x1 + v[9]*x2 + v[14]*x3 + v[19]*x4 + v[24]*x5; 4982d61bbb3SSatish Balay v += 25; 4992d61bbb3SSatish Balay } 5007b2bb3b9SHong Zhang if (usecprow) z = zarray + 5*ridx[i]; 5012d61bbb3SSatish Balay z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; 50226e093fcSHong Zhang if (!usecprow) z += 5; 5032d61bbb3SSatish Balay } 5043649974fSBarry Smith ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr); 50526e093fcSHong Zhang ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr); 5067c565772SBarry Smith ierr = PetscLogFlops(50.0*a->nz - 5.0*a->nonzerorowcnt);CHKERRQ(ierr); 5072d61bbb3SSatish Balay PetscFunctionReturn(0); 5082d61bbb3SSatish Balay } 5092d61bbb3SSatish Balay 51015091d37SBarry Smith 511d6232bceSMichael Lange 5124a2ae208SSatish Balay #undef __FUNCT__ 5134a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_6" 514dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_6(Mat A,Vec xx,Vec zz) 51515091d37SBarry Smith { 51615091d37SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 517d9fead3dSBarry Smith PetscScalar *z = 0,sum1,sum2,sum3,sum4,sum5,sum6; 518d9fead3dSBarry Smith const PetscScalar *x,*xb; 51926e093fcSHong Zhang PetscScalar x1,x2,x3,x4,x5,x6,*zarray; 520d9fead3dSBarry Smith const MatScalar *v; 521dfbe8321SBarry Smith PetscErrorCode ierr; 5227c565772SBarry Smith PetscInt mbs,i,*idx,*ii,j,n,*ridx=NULL; 523ace3abfcSBarry Smith PetscBool usecprow=a->compressedrow.use; 52415091d37SBarry Smith 525433994e6SBarry Smith PetscFunctionBegin; 5263649974fSBarry Smith ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr); 52726e093fcSHong Zhang ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr); 52815091d37SBarry Smith 52915091d37SBarry Smith idx = a->j; 53015091d37SBarry Smith v = a->a; 53126e093fcSHong Zhang if (usecprow) { 53226e093fcSHong Zhang mbs = a->compressedrow.nrows; 53326e093fcSHong Zhang ii = a->compressedrow.i; 5347b2bb3b9SHong Zhang ridx = a->compressedrow.rindex; 53526e093fcSHong Zhang } else { 53626e093fcSHong Zhang mbs = a->mbs; 53715091d37SBarry Smith ii = a->i; 53826e093fcSHong Zhang z = zarray; 53926e093fcSHong Zhang } 54015091d37SBarry Smith 54115091d37SBarry Smith for (i=0; i<mbs; i++) { 54226fbe8dcSKarl Rupp n = ii[1] - ii[0]; 54326fbe8dcSKarl Rupp ii++; 54426fbe8dcSKarl Rupp sum1 = 0.0; 54526fbe8dcSKarl Rupp sum2 = 0.0; 54626fbe8dcSKarl Rupp sum3 = 0.0; 54726fbe8dcSKarl Rupp sum4 = 0.0; 54826fbe8dcSKarl Rupp sum5 = 0.0; 54926fbe8dcSKarl Rupp sum6 = 0.0; 55026fbe8dcSKarl Rupp 551444d8c10SJed Brown PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA); /* Indices for the next row (assumes same size as this one) */ 552444d8c10SJed Brown PetscPrefetchBlock(v+36*n,36*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */ 55315091d37SBarry Smith for (j=0; j<n; j++) { 55415091d37SBarry Smith xb = x + 6*(*idx++); 55515091d37SBarry Smith x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; 55615091d37SBarry Smith sum1 += v[0]*x1 + v[6]*x2 + v[12]*x3 + v[18]*x4 + v[24]*x5 + v[30]*x6; 55715091d37SBarry Smith sum2 += v[1]*x1 + v[7]*x2 + v[13]*x3 + v[19]*x4 + v[25]*x5 + v[31]*x6; 55815091d37SBarry Smith sum3 += v[2]*x1 + v[8]*x2 + v[14]*x3 + v[20]*x4 + v[26]*x5 + v[32]*x6; 55915091d37SBarry Smith sum4 += v[3]*x1 + v[9]*x2 + v[15]*x3 + v[21]*x4 + v[27]*x5 + v[33]*x6; 56015091d37SBarry Smith sum5 += v[4]*x1 + v[10]*x2 + v[16]*x3 + v[22]*x4 + v[28]*x5 + v[34]*x6; 56115091d37SBarry Smith sum6 += v[5]*x1 + v[11]*x2 + v[17]*x3 + v[23]*x4 + v[29]*x5 + v[35]*x6; 56215091d37SBarry Smith v += 36; 56315091d37SBarry Smith } 5647b2bb3b9SHong Zhang if (usecprow) z = zarray + 6*ridx[i]; 56515091d37SBarry Smith z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; 56626e093fcSHong Zhang if (!usecprow) z += 6; 56715091d37SBarry Smith } 56815091d37SBarry Smith 5693649974fSBarry Smith ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr); 57026e093fcSHong Zhang ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr); 5717c565772SBarry Smith ierr = PetscLogFlops(72.0*a->nz - 6.0*a->nonzerorowcnt);CHKERRQ(ierr); 57215091d37SBarry Smith PetscFunctionReturn(0); 57315091d37SBarry Smith } 5748ab949d8SShri Abhyankar 5754a2ae208SSatish Balay #undef __FUNCT__ 5764a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_7" 577dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_7(Mat A,Vec xx,Vec zz) 5782d61bbb3SSatish Balay { 5792d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 580d9fead3dSBarry Smith PetscScalar *z = 0,sum1,sum2,sum3,sum4,sum5,sum6,sum7; 581d9fead3dSBarry Smith const PetscScalar *x,*xb; 58226e093fcSHong Zhang PetscScalar x1,x2,x3,x4,x5,x6,x7,*zarray; 583d9fead3dSBarry Smith const MatScalar *v; 584dfbe8321SBarry Smith PetscErrorCode ierr; 5857c565772SBarry Smith PetscInt mbs,i,*idx,*ii,j,n,*ridx=NULL; 586ace3abfcSBarry Smith PetscBool usecprow=a->compressedrow.use; 5872d61bbb3SSatish Balay 588433994e6SBarry Smith PetscFunctionBegin; 5893649974fSBarry Smith ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr); 59026e093fcSHong Zhang ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr); 5912d61bbb3SSatish Balay 5922d61bbb3SSatish Balay idx = a->j; 5932d61bbb3SSatish Balay v = a->a; 59426e093fcSHong Zhang if (usecprow) { 59526e093fcSHong Zhang mbs = a->compressedrow.nrows; 59626e093fcSHong Zhang ii = a->compressedrow.i; 5977b2bb3b9SHong Zhang ridx = a->compressedrow.rindex; 59826e093fcSHong Zhang } else { 59926e093fcSHong Zhang mbs = a->mbs; 6002d61bbb3SSatish Balay ii = a->i; 60126e093fcSHong Zhang z = zarray; 60226e093fcSHong Zhang } 6032d61bbb3SSatish Balay 6042d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 60526fbe8dcSKarl Rupp n = ii[1] - ii[0]; 60626fbe8dcSKarl Rupp ii++; 60726fbe8dcSKarl Rupp sum1 = 0.0; 60826fbe8dcSKarl Rupp sum2 = 0.0; 60926fbe8dcSKarl Rupp sum3 = 0.0; 61026fbe8dcSKarl Rupp sum4 = 0.0; 61126fbe8dcSKarl Rupp sum5 = 0.0; 61226fbe8dcSKarl Rupp sum6 = 0.0; 61326fbe8dcSKarl Rupp sum7 = 0.0; 61426fbe8dcSKarl Rupp 615444d8c10SJed Brown PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA); /* Indices for the next row (assumes same size as this one) */ 616444d8c10SJed Brown PetscPrefetchBlock(v+49*n,49*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */ 6172d61bbb3SSatish Balay for (j=0; j<n; j++) { 6182d61bbb3SSatish Balay xb = x + 7*(*idx++); 6192d61bbb3SSatish Balay x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6]; 6202d61bbb3SSatish Balay sum1 += v[0]*x1 + v[7]*x2 + v[14]*x3 + v[21]*x4 + v[28]*x5 + v[35]*x6 + v[42]*x7; 6212d61bbb3SSatish Balay sum2 += v[1]*x1 + v[8]*x2 + v[15]*x3 + v[22]*x4 + v[29]*x5 + v[36]*x6 + v[43]*x7; 6222d61bbb3SSatish Balay sum3 += v[2]*x1 + v[9]*x2 + v[16]*x3 + v[23]*x4 + v[30]*x5 + v[37]*x6 + v[44]*x7; 6232d61bbb3SSatish Balay sum4 += v[3]*x1 + v[10]*x2 + v[17]*x3 + v[24]*x4 + v[31]*x5 + v[38]*x6 + v[45]*x7; 6242d61bbb3SSatish Balay sum5 += v[4]*x1 + v[11]*x2 + v[18]*x3 + v[25]*x4 + v[32]*x5 + v[39]*x6 + v[46]*x7; 6252d61bbb3SSatish Balay sum6 += v[5]*x1 + v[12]*x2 + v[19]*x3 + v[26]*x4 + v[33]*x5 + v[40]*x6 + v[47]*x7; 6262d61bbb3SSatish Balay sum7 += v[6]*x1 + v[13]*x2 + v[20]*x3 + v[27]*x4 + v[34]*x5 + v[41]*x6 + v[48]*x7; 6272d61bbb3SSatish Balay v += 49; 6282d61bbb3SSatish Balay } 6297b2bb3b9SHong Zhang if (usecprow) z = zarray + 7*ridx[i]; 6302d61bbb3SSatish Balay z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7; 63126e093fcSHong Zhang if (!usecprow) z += 7; 6322d61bbb3SSatish Balay } 6332d61bbb3SSatish Balay 6343649974fSBarry Smith ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr); 63526e093fcSHong Zhang ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr); 6367c565772SBarry Smith ierr = PetscLogFlops(98.0*a->nz - 7.0*a->nonzerorowcnt);CHKERRQ(ierr); 6372d61bbb3SSatish Balay PetscFunctionReturn(0); 6382d61bbb3SSatish Balay } 6392d61bbb3SSatish Balay 6408ab949d8SShri Abhyankar /* MatMult_SeqBAIJ_15 version 1: Columns in the block are accessed one at a time */ 641832cc040SShri Abhyankar /* Default MatMult for block size 15 */ 6428ab949d8SShri Abhyankar 6438ab949d8SShri Abhyankar #undef __FUNCT__ 6448ab949d8SShri Abhyankar #define __FUNCT__ "MatMult_SeqBAIJ_15_ver1" 6458ab949d8SShri Abhyankar PetscErrorCode MatMult_SeqBAIJ_15_ver1(Mat A,Vec xx,Vec zz) 6468ab949d8SShri Abhyankar { 6478ab949d8SShri Abhyankar Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 6488ab949d8SShri Abhyankar PetscScalar *z = 0,sum1,sum2,sum3,sum4,sum5,sum6,sum7,sum8,sum9,sum10,sum11,sum12,sum13,sum14,sum15; 6498ab949d8SShri Abhyankar const PetscScalar *x,*xb; 65053ef36baSBarry Smith PetscScalar *zarray,xv; 6518ab949d8SShri Abhyankar const MatScalar *v; 6528ab949d8SShri Abhyankar PetscErrorCode ierr; 6538ab949d8SShri Abhyankar const PetscInt *ii,*ij=a->j,*idx; 6547c565772SBarry Smith PetscInt mbs,i,j,k,n,*ridx=NULL; 655ace3abfcSBarry Smith PetscBool usecprow=a->compressedrow.use; 6568ab949d8SShri Abhyankar 6578ab949d8SShri Abhyankar PetscFunctionBegin; 6583649974fSBarry Smith ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr); 6598ab949d8SShri Abhyankar ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr); 6608ab949d8SShri Abhyankar 6618ab949d8SShri Abhyankar v = a->a; 6628ab949d8SShri Abhyankar if (usecprow) { 6638ab949d8SShri Abhyankar mbs = a->compressedrow.nrows; 6648ab949d8SShri Abhyankar ii = a->compressedrow.i; 6658ab949d8SShri Abhyankar ridx = a->compressedrow.rindex; 6668ab949d8SShri Abhyankar } else { 6678ab949d8SShri Abhyankar mbs = a->mbs; 6688ab949d8SShri Abhyankar ii = a->i; 6698ab949d8SShri Abhyankar z = zarray; 6708ab949d8SShri Abhyankar } 6718ab949d8SShri Abhyankar 6728ab949d8SShri Abhyankar for (i=0; i<mbs; i++) { 6738ab949d8SShri Abhyankar n = ii[i+1] - ii[i]; 6748ab949d8SShri Abhyankar idx = ij + ii[i]; 6758ab949d8SShri Abhyankar sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0; sum7 = 0.0; 6768ab949d8SShri 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; 6778ab949d8SShri Abhyankar 6788ab949d8SShri Abhyankar for (j=0; j<n; j++) { 6798ab949d8SShri Abhyankar xb = x + 15*(idx[j]); 6808ab949d8SShri Abhyankar 6818ab949d8SShri Abhyankar for (k=0; k<15; k++) { 68253ef36baSBarry Smith xv = xb[k]; 68353ef36baSBarry Smith sum1 += v[0]*xv; 68453ef36baSBarry Smith sum2 += v[1]*xv; 68553ef36baSBarry Smith sum3 += v[2]*xv; 68653ef36baSBarry Smith sum4 += v[3]*xv; 68753ef36baSBarry Smith sum5 += v[4]*xv; 68853ef36baSBarry Smith sum6 += v[5]*xv; 68953ef36baSBarry Smith sum7 += v[6]*xv; 69053ef36baSBarry Smith sum8 += v[7]*xv; 69153ef36baSBarry Smith sum9 += v[8]*xv; 69253ef36baSBarry Smith sum10 += v[9]*xv; 69353ef36baSBarry Smith sum11 += v[10]*xv; 69453ef36baSBarry Smith sum12 += v[11]*xv; 69553ef36baSBarry Smith sum13 += v[12]*xv; 69653ef36baSBarry Smith sum14 += v[13]*xv; 69753ef36baSBarry Smith sum15 += v[14]*xv; 6988ab949d8SShri Abhyankar v += 15; 6998ab949d8SShri Abhyankar } 7008ab949d8SShri Abhyankar } 7018ab949d8SShri Abhyankar if (usecprow) z = zarray + 15*ridx[i]; 7028ab949d8SShri Abhyankar z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7; 7038ab949d8SShri 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; 7048ab949d8SShri Abhyankar 7058ab949d8SShri Abhyankar if (!usecprow) z += 15; 7068ab949d8SShri Abhyankar } 7078ab949d8SShri Abhyankar 7083649974fSBarry Smith ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr); 7098ab949d8SShri Abhyankar ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr); 7107c565772SBarry Smith ierr = PetscLogFlops(450.0*a->nz - 15.0*a->nonzerorowcnt);CHKERRQ(ierr); 7118ab949d8SShri Abhyankar PetscFunctionReturn(0); 7128ab949d8SShri Abhyankar } 7138ab949d8SShri Abhyankar 7148ab949d8SShri Abhyankar /* MatMult_SeqBAIJ_15_ver2 : Columns in the block are accessed in sets of 4,4,4,3 */ 7158ab949d8SShri Abhyankar #undef __FUNCT__ 7168ab949d8SShri Abhyankar #define __FUNCT__ "MatMult_SeqBAIJ_15_ver2" 7178ab949d8SShri Abhyankar PetscErrorCode MatMult_SeqBAIJ_15_ver2(Mat A,Vec xx,Vec zz) 7188ab949d8SShri Abhyankar { 7198ab949d8SShri Abhyankar Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 7208ab949d8SShri Abhyankar PetscScalar *z = 0,sum1,sum2,sum3,sum4,sum5,sum6,sum7,sum8,sum9,sum10,sum11,sum12,sum13,sum14,sum15; 7218ab949d8SShri Abhyankar const PetscScalar *x,*xb; 7220b8f6341SShri Abhyankar PetscScalar x1,x2,x3,x4,*zarray; 7238ab949d8SShri Abhyankar const MatScalar *v; 7248ab949d8SShri Abhyankar PetscErrorCode ierr; 7258ab949d8SShri Abhyankar const PetscInt *ii,*ij=a->j,*idx; 7267c565772SBarry Smith PetscInt mbs,i,j,n,*ridx=NULL; 727ace3abfcSBarry Smith PetscBool usecprow=a->compressedrow.use; 7288ab949d8SShri Abhyankar 7298ab949d8SShri Abhyankar PetscFunctionBegin; 7303649974fSBarry Smith ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr); 7318ab949d8SShri Abhyankar ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr); 7328ab949d8SShri Abhyankar 7338ab949d8SShri Abhyankar v = a->a; 7348ab949d8SShri Abhyankar if (usecprow) { 7358ab949d8SShri Abhyankar mbs = a->compressedrow.nrows; 7368ab949d8SShri Abhyankar ii = a->compressedrow.i; 7378ab949d8SShri Abhyankar ridx = a->compressedrow.rindex; 7388ab949d8SShri Abhyankar } else { 7398ab949d8SShri Abhyankar mbs = a->mbs; 7408ab949d8SShri Abhyankar ii = a->i; 7418ab949d8SShri Abhyankar z = zarray; 7428ab949d8SShri Abhyankar } 7438ab949d8SShri Abhyankar 7448ab949d8SShri Abhyankar for (i=0; i<mbs; i++) { 7458ab949d8SShri Abhyankar n = ii[i+1] - ii[i]; 7468ab949d8SShri Abhyankar idx = ij + ii[i]; 7478ab949d8SShri Abhyankar sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0; sum7 = 0.0; 7488ab949d8SShri 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; 7498ab949d8SShri Abhyankar 7508ab949d8SShri Abhyankar for (j=0; j<n; j++) { 7518ab949d8SShri Abhyankar xb = x + 15*(idx[j]); 7520b8f6341SShri Abhyankar x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; 7538ab949d8SShri Abhyankar 7548ab949d8SShri Abhyankar sum1 += v[0]*x1 + v[15]*x2 + v[30]*x3 + v[45]*x4; 7558ab949d8SShri Abhyankar sum2 += v[1]*x1 + v[16]*x2 + v[31]*x3 + v[46]*x4; 7568ab949d8SShri Abhyankar sum3 += v[2]*x1 + v[17]*x2 + v[32]*x3 + v[47]*x4; 7578ab949d8SShri Abhyankar sum4 += v[3]*x1 + v[18]*x2 + v[33]*x3 + v[48]*x4; 7588ab949d8SShri Abhyankar sum5 += v[4]*x1 + v[19]*x2 + v[34]*x3 + v[49]*x4; 7598ab949d8SShri Abhyankar sum6 += v[5]*x1 + v[20]*x2 + v[35]*x3 + v[50]*x4; 7608ab949d8SShri Abhyankar sum7 += v[6]*x1 + v[21]*x2 + v[36]*x3 + v[51]*x4; 7618ab949d8SShri Abhyankar sum8 += v[7]*x1 + v[22]*x2 + v[37]*x3 + v[52]*x4; 7628ab949d8SShri Abhyankar sum9 += v[8]*x1 + v[23]*x2 + v[38]*x3 + v[53]*x4; 7638ab949d8SShri Abhyankar sum10 += v[9]*x1 + v[24]*x2 + v[39]*x3 + v[54]*x4; 7648ab949d8SShri Abhyankar sum11 += v[10]*x1 + v[25]*x2 + v[40]*x3 + v[55]*x4; 7658ab949d8SShri Abhyankar sum12 += v[11]*x1 + v[26]*x2 + v[41]*x3 + v[56]*x4; 7668ab949d8SShri Abhyankar sum13 += v[12]*x1 + v[27]*x2 + v[42]*x3 + v[57]*x4; 7678ab949d8SShri Abhyankar sum14 += v[13]*x1 + v[28]*x2 + v[43]*x3 + v[58]*x4; 7688ab949d8SShri Abhyankar sum15 += v[14]*x1 + v[29]*x2 + v[44]*x3 + v[59]*x4; 7698ab949d8SShri Abhyankar 7708ab949d8SShri Abhyankar v += 60; 7718ab949d8SShri Abhyankar 7720b8f6341SShri Abhyankar x1 = xb[4]; x2 = xb[5]; x3 = xb[6]; x4 = xb[7]; 7738ab949d8SShri Abhyankar 7748ab949d8SShri Abhyankar sum1 += v[0]*x1 + v[15]*x2 + v[30]*x3 + v[45]*x4; 7758ab949d8SShri Abhyankar sum2 += v[1]*x1 + v[16]*x2 + v[31]*x3 + v[46]*x4; 7768ab949d8SShri Abhyankar sum3 += v[2]*x1 + v[17]*x2 + v[32]*x3 + v[47]*x4; 7778ab949d8SShri Abhyankar sum4 += v[3]*x1 + v[18]*x2 + v[33]*x3 + v[48]*x4; 7788ab949d8SShri Abhyankar sum5 += v[4]*x1 + v[19]*x2 + v[34]*x3 + v[49]*x4; 7798ab949d8SShri Abhyankar sum6 += v[5]*x1 + v[20]*x2 + v[35]*x3 + v[50]*x4; 7808ab949d8SShri Abhyankar sum7 += v[6]*x1 + v[21]*x2 + v[36]*x3 + v[51]*x4; 7818ab949d8SShri Abhyankar sum8 += v[7]*x1 + v[22]*x2 + v[37]*x3 + v[52]*x4; 7828ab949d8SShri Abhyankar sum9 += v[8]*x1 + v[23]*x2 + v[38]*x3 + v[53]*x4; 7838ab949d8SShri Abhyankar sum10 += v[9]*x1 + v[24]*x2 + v[39]*x3 + v[54]*x4; 7848ab949d8SShri Abhyankar sum11 += v[10]*x1 + v[25]*x2 + v[40]*x3 + v[55]*x4; 7858ab949d8SShri Abhyankar sum12 += v[11]*x1 + v[26]*x2 + v[41]*x3 + v[56]*x4; 7868ab949d8SShri Abhyankar sum13 += v[12]*x1 + v[27]*x2 + v[42]*x3 + v[57]*x4; 7878ab949d8SShri Abhyankar sum14 += v[13]*x1 + v[28]*x2 + v[43]*x3 + v[58]*x4; 7888ab949d8SShri Abhyankar sum15 += v[14]*x1 + v[29]*x2 + v[44]*x3 + v[59]*x4; 7898ab949d8SShri Abhyankar v += 60; 7908ab949d8SShri Abhyankar 7910b8f6341SShri Abhyankar x1 = xb[8]; x2 = xb[9]; x3 = xb[10]; x4 = xb[11]; 7920b8f6341SShri Abhyankar sum1 += v[0]*x1 + v[15]*x2 + v[30]*x3 + v[45]*x4; 7930b8f6341SShri Abhyankar sum2 += v[1]*x1 + v[16]*x2 + v[31]*x3 + v[46]*x4; 7940b8f6341SShri Abhyankar sum3 += v[2]*x1 + v[17]*x2 + v[32]*x3 + v[47]*x4; 7950b8f6341SShri Abhyankar sum4 += v[3]*x1 + v[18]*x2 + v[33]*x3 + v[48]*x4; 7960b8f6341SShri Abhyankar sum5 += v[4]*x1 + v[19]*x2 + v[34]*x3 + v[49]*x4; 7970b8f6341SShri Abhyankar sum6 += v[5]*x1 + v[20]*x2 + v[35]*x3 + v[50]*x4; 7980b8f6341SShri Abhyankar sum7 += v[6]*x1 + v[21]*x2 + v[36]*x3 + v[51]*x4; 7990b8f6341SShri Abhyankar sum8 += v[7]*x1 + v[22]*x2 + v[37]*x3 + v[52]*x4; 8000b8f6341SShri Abhyankar sum9 += v[8]*x1 + v[23]*x2 + v[38]*x3 + v[53]*x4; 8010b8f6341SShri Abhyankar sum10 += v[9]*x1 + v[24]*x2 + v[39]*x3 + v[54]*x4; 8020b8f6341SShri Abhyankar sum11 += v[10]*x1 + v[25]*x2 + v[40]*x3 + v[55]*x4; 8030b8f6341SShri Abhyankar sum12 += v[11]*x1 + v[26]*x2 + v[41]*x3 + v[56]*x4; 8040b8f6341SShri Abhyankar sum13 += v[12]*x1 + v[27]*x2 + v[42]*x3 + v[57]*x4; 8050b8f6341SShri Abhyankar sum14 += v[13]*x1 + v[28]*x2 + v[43]*x3 + v[58]*x4; 8060b8f6341SShri Abhyankar sum15 += v[14]*x1 + v[29]*x2 + v[44]*x3 + v[59]*x4; 8070b8f6341SShri Abhyankar v += 60; 8080b8f6341SShri Abhyankar 8090b8f6341SShri Abhyankar x1 = xb[12]; x2 = xb[13]; x3 = xb[14]; 8108ab949d8SShri Abhyankar sum1 += v[0]*x1 + v[15]*x2 + v[30]*x3; 8118ab949d8SShri Abhyankar sum2 += v[1]*x1 + v[16]*x2 + v[31]*x3; 8128ab949d8SShri Abhyankar sum3 += v[2]*x1 + v[17]*x2 + v[32]*x3; 8138ab949d8SShri Abhyankar sum4 += v[3]*x1 + v[18]*x2 + v[33]*x3; 8148ab949d8SShri Abhyankar sum5 += v[4]*x1 + v[19]*x2 + v[34]*x3; 8158ab949d8SShri Abhyankar sum6 += v[5]*x1 + v[20]*x2 + v[35]*x3; 8168ab949d8SShri Abhyankar sum7 += v[6]*x1 + v[21]*x2 + v[36]*x3; 8178ab949d8SShri Abhyankar sum8 += v[7]*x1 + v[22]*x2 + v[37]*x3; 8188ab949d8SShri Abhyankar sum9 += v[8]*x1 + v[23]*x2 + v[38]*x3; 8198ab949d8SShri Abhyankar sum10 += v[9]*x1 + v[24]*x2 + v[39]*x3; 8208ab949d8SShri Abhyankar sum11 += v[10]*x1 + v[25]*x2 + v[40]*x3; 8218ab949d8SShri Abhyankar sum12 += v[11]*x1 + v[26]*x2 + v[41]*x3; 8228ab949d8SShri Abhyankar sum13 += v[12]*x1 + v[27]*x2 + v[42]*x3; 8238ab949d8SShri Abhyankar sum14 += v[13]*x1 + v[28]*x2 + v[43]*x3; 8248ab949d8SShri Abhyankar sum15 += v[14]*x1 + v[29]*x2 + v[44]*x3; 8258ab949d8SShri Abhyankar v += 45; 8268ab949d8SShri Abhyankar } 8278ab949d8SShri Abhyankar if (usecprow) z = zarray + 15*ridx[i]; 8288ab949d8SShri Abhyankar z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7; 8298ab949d8SShri 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; 8308ab949d8SShri Abhyankar 8318ab949d8SShri Abhyankar if (!usecprow) z += 15; 8328ab949d8SShri Abhyankar } 8338ab949d8SShri Abhyankar 8343649974fSBarry Smith ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr); 8358ab949d8SShri Abhyankar ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr); 8367c565772SBarry Smith ierr = PetscLogFlops(450.0*a->nz - 15.0*a->nonzerorowcnt);CHKERRQ(ierr); 8378ab949d8SShri Abhyankar PetscFunctionReturn(0); 8388ab949d8SShri Abhyankar } 8398ab949d8SShri Abhyankar 8408ab949d8SShri Abhyankar /* MatMult_SeqBAIJ_15_ver3 : Columns in the block are accessed in sets of 8,7 */ 8418ab949d8SShri Abhyankar #undef __FUNCT__ 8428ab949d8SShri Abhyankar #define __FUNCT__ "MatMult_SeqBAIJ_15_ver3" 8438ab949d8SShri Abhyankar PetscErrorCode MatMult_SeqBAIJ_15_ver3(Mat A,Vec xx,Vec zz) 8448ab949d8SShri Abhyankar { 8458ab949d8SShri Abhyankar Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 8468ab949d8SShri Abhyankar PetscScalar *z = 0,sum1,sum2,sum3,sum4,sum5,sum6,sum7,sum8,sum9,sum10,sum11,sum12,sum13,sum14,sum15; 8478ab949d8SShri Abhyankar const PetscScalar *x,*xb; 8480b8f6341SShri Abhyankar PetscScalar x1,x2,x3,x4,x5,x6,x7,x8,*zarray; 8498ab949d8SShri Abhyankar const MatScalar *v; 8508ab949d8SShri Abhyankar PetscErrorCode ierr; 8518ab949d8SShri Abhyankar const PetscInt *ii,*ij=a->j,*idx; 8527c565772SBarry Smith PetscInt mbs,i,j,n,*ridx=NULL; 853ace3abfcSBarry Smith PetscBool usecprow=a->compressedrow.use; 8548ab949d8SShri Abhyankar 8558ab949d8SShri Abhyankar PetscFunctionBegin; 8563649974fSBarry Smith ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr); 8578ab949d8SShri Abhyankar ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr); 8588ab949d8SShri Abhyankar 8598ab949d8SShri Abhyankar v = a->a; 8608ab949d8SShri Abhyankar if (usecprow) { 8618ab949d8SShri Abhyankar mbs = a->compressedrow.nrows; 8628ab949d8SShri Abhyankar ii = a->compressedrow.i; 8638ab949d8SShri Abhyankar ridx = a->compressedrow.rindex; 8648ab949d8SShri Abhyankar } else { 8658ab949d8SShri Abhyankar mbs = a->mbs; 8668ab949d8SShri Abhyankar ii = a->i; 8678ab949d8SShri Abhyankar z = zarray; 8688ab949d8SShri Abhyankar } 8698ab949d8SShri Abhyankar 8708ab949d8SShri Abhyankar for (i=0; i<mbs; i++) { 8718ab949d8SShri Abhyankar n = ii[i+1] - ii[i]; 8728ab949d8SShri Abhyankar idx = ij + ii[i]; 8738ab949d8SShri Abhyankar sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0; sum7 = 0.0; 8748ab949d8SShri 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; 8758ab949d8SShri Abhyankar 8768ab949d8SShri Abhyankar for (j=0; j<n; j++) { 8778ab949d8SShri Abhyankar xb = x + 15*(idx[j]); 8788ab949d8SShri Abhyankar x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6]; 8790b8f6341SShri Abhyankar x8 = xb[7]; 8808ab949d8SShri Abhyankar 8818ab949d8SShri 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; 8828ab949d8SShri 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; 8838ab949d8SShri 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; 8848ab949d8SShri 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; 8858ab949d8SShri 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; 8868ab949d8SShri 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; 8878ab949d8SShri 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; 8888ab949d8SShri 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; 8898ab949d8SShri 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; 8908ab949d8SShri 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; 8918ab949d8SShri 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; 8928ab949d8SShri 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; 8938ab949d8SShri 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; 8948ab949d8SShri 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; 8958ab949d8SShri 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; 8968ab949d8SShri Abhyankar v += 120; 8978ab949d8SShri Abhyankar 8980b8f6341SShri Abhyankar x1 = xb[8]; x2 = xb[9]; x3 = xb[10]; x4 = xb[11]; x5 = xb[12]; x6 = xb[13]; x7 = xb[14]; 8990b8f6341SShri Abhyankar 9008ab949d8SShri Abhyankar sum1 += v[0]*x1 + v[15]*x2 + v[30]*x3 + v[45]*x4 + v[60]*x5 + v[75]*x6 + v[90]*x7; 9018ab949d8SShri Abhyankar sum2 += v[1]*x1 + v[16]*x2 + v[31]*x3 + v[46]*x4 + v[61]*x5 + v[76]*x6 + v[91]*x7; 9028ab949d8SShri Abhyankar sum3 += v[2]*x1 + v[17]*x2 + v[32]*x3 + v[47]*x4 + v[62]*x5 + v[77]*x6 + v[92]*x7; 9038ab949d8SShri Abhyankar sum4 += v[3]*x1 + v[18]*x2 + v[33]*x3 + v[48]*x4 + v[63]*x5 + v[78]*x6 + v[93]*x7; 9048ab949d8SShri Abhyankar sum5 += v[4]*x1 + v[19]*x2 + v[34]*x3 + v[49]*x4 + v[64]*x5 + v[79]*x6 + v[94]*x7; 9058ab949d8SShri Abhyankar sum6 += v[5]*x1 + v[20]*x2 + v[35]*x3 + v[50]*x4 + v[65]*x5 + v[80]*x6 + v[95]*x7; 9068ab949d8SShri Abhyankar sum7 += v[6]*x1 + v[21]*x2 + v[36]*x3 + v[51]*x4 + v[66]*x5 + v[81]*x6 + v[96]*x7; 9078ab949d8SShri Abhyankar sum8 += v[7]*x1 + v[22]*x2 + v[37]*x3 + v[52]*x4 + v[67]*x5 + v[82]*x6 + v[97]*x7; 9088ab949d8SShri Abhyankar sum9 += v[8]*x1 + v[23]*x2 + v[38]*x3 + v[53]*x4 + v[68]*x5 + v[83]*x6 + v[98]*x7; 9098ab949d8SShri Abhyankar sum10 += v[9]*x1 + v[24]*x2 + v[39]*x3 + v[54]*x4 + v[69]*x5 + v[84]*x6 + v[99]*x7; 9108ab949d8SShri Abhyankar sum11 += v[10]*x1 + v[25]*x2 + v[40]*x3 + v[55]*x4 + v[70]*x5 + v[85]*x6 + v[100]*x7; 9118ab949d8SShri Abhyankar sum12 += v[11]*x1 + v[26]*x2 + v[41]*x3 + v[56]*x4 + v[71]*x5 + v[86]*x6 + v[101]*x7; 9128ab949d8SShri Abhyankar sum13 += v[12]*x1 + v[27]*x2 + v[42]*x3 + v[57]*x4 + v[72]*x5 + v[87]*x6 + v[102]*x7; 9138ab949d8SShri Abhyankar sum14 += v[13]*x1 + v[28]*x2 + v[43]*x3 + v[58]*x4 + v[73]*x5 + v[88]*x6 + v[103]*x7; 9148ab949d8SShri Abhyankar sum15 += v[14]*x1 + v[29]*x2 + v[44]*x3 + v[59]*x4 + v[74]*x5 + v[89]*x6 + v[104]*x7; 9158ab949d8SShri Abhyankar v += 105; 9168ab949d8SShri Abhyankar } 9178ab949d8SShri Abhyankar if (usecprow) z = zarray + 15*ridx[i]; 9188ab949d8SShri Abhyankar z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7; 9198ab949d8SShri 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; 9208ab949d8SShri Abhyankar 9218ab949d8SShri Abhyankar if (!usecprow) z += 15; 9228ab949d8SShri Abhyankar } 9238ab949d8SShri Abhyankar 9243649974fSBarry Smith ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr); 9258ab949d8SShri Abhyankar ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr); 9267c565772SBarry Smith ierr = PetscLogFlops(450.0*a->nz - 15.0*a->nonzerorowcnt);CHKERRQ(ierr); 9278ab949d8SShri Abhyankar PetscFunctionReturn(0); 9288ab949d8SShri Abhyankar } 9298ab949d8SShri Abhyankar 9308ab949d8SShri Abhyankar /* MatMult_SeqBAIJ_15_ver4 : All columns in the block are accessed at once */ 9318ab949d8SShri Abhyankar 9328ab949d8SShri Abhyankar #undef __FUNCT__ 9338ab949d8SShri Abhyankar #define __FUNCT__ "MatMult_SeqBAIJ_15_ver4" 9348ab949d8SShri Abhyankar PetscErrorCode MatMult_SeqBAIJ_15_ver4(Mat A,Vec xx,Vec zz) 9358ab949d8SShri Abhyankar { 9368ab949d8SShri Abhyankar Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 9378ab949d8SShri Abhyankar PetscScalar *z = 0,sum1,sum2,sum3,sum4,sum5,sum6,sum7,sum8,sum9,sum10,sum11,sum12,sum13,sum14,sum15; 9388ab949d8SShri Abhyankar const PetscScalar *x,*xb; 9398ab949d8SShri Abhyankar PetscScalar x1,x2,x3,x4,x5,x6,x7,x8,x9,x10,x11,x12,x13,x14,x15,*zarray; 9408ab949d8SShri Abhyankar const MatScalar *v; 9418ab949d8SShri Abhyankar PetscErrorCode ierr; 9428ab949d8SShri Abhyankar const PetscInt *ii,*ij=a->j,*idx; 9437c565772SBarry Smith PetscInt mbs,i,j,n,*ridx=NULL; 944ace3abfcSBarry Smith PetscBool usecprow=a->compressedrow.use; 9458ab949d8SShri Abhyankar 9468ab949d8SShri Abhyankar PetscFunctionBegin; 9473649974fSBarry Smith ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr); 9488ab949d8SShri Abhyankar ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr); 9498ab949d8SShri Abhyankar 9508ab949d8SShri Abhyankar v = a->a; 9518ab949d8SShri Abhyankar if (usecprow) { 9528ab949d8SShri Abhyankar mbs = a->compressedrow.nrows; 9538ab949d8SShri Abhyankar ii = a->compressedrow.i; 9548ab949d8SShri Abhyankar ridx = a->compressedrow.rindex; 9558ab949d8SShri Abhyankar } else { 9568ab949d8SShri Abhyankar mbs = a->mbs; 9578ab949d8SShri Abhyankar ii = a->i; 9588ab949d8SShri Abhyankar z = zarray; 9598ab949d8SShri Abhyankar } 9608ab949d8SShri Abhyankar 9618ab949d8SShri Abhyankar for (i=0; i<mbs; i++) { 9628ab949d8SShri Abhyankar n = ii[i+1] - ii[i]; 9638ab949d8SShri Abhyankar idx = ij + ii[i]; 9648ab949d8SShri Abhyankar sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0; sum7 = 0.0; 9658ab949d8SShri 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; 9668ab949d8SShri Abhyankar 9678ab949d8SShri Abhyankar for (j=0; j<n; j++) { 9688ab949d8SShri Abhyankar xb = x + 15*(idx[j]); 9698ab949d8SShri Abhyankar x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6]; 9708ab949d8SShri 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]; 9718ab949d8SShri Abhyankar 9728ab949d8SShri 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; 9738ab949d8SShri 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; 9748ab949d8SShri 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; 9758ab949d8SShri 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; 9768ab949d8SShri 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; 9778ab949d8SShri 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; 9788ab949d8SShri 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; 9798ab949d8SShri 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; 9808ab949d8SShri 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; 9818ab949d8SShri 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; 9828ab949d8SShri 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; 9838ab949d8SShri 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; 9848ab949d8SShri 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; 9858ab949d8SShri 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; 9868ab949d8SShri 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; 9878ab949d8SShri Abhyankar v += 225; 9888ab949d8SShri Abhyankar } 9898ab949d8SShri Abhyankar if (usecprow) z = zarray + 15*ridx[i]; 9908ab949d8SShri Abhyankar z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7; 9918ab949d8SShri 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; 9928ab949d8SShri Abhyankar 9938ab949d8SShri Abhyankar if (!usecprow) z += 15; 9948ab949d8SShri Abhyankar } 9958ab949d8SShri Abhyankar 9963649974fSBarry Smith ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr); 9978ab949d8SShri Abhyankar ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr); 9987c565772SBarry Smith ierr = PetscLogFlops(450.0*a->nz - 15.0*a->nonzerorowcnt);CHKERRQ(ierr); 9998ab949d8SShri Abhyankar PetscFunctionReturn(0); 10008ab949d8SShri Abhyankar } 10018ab949d8SShri Abhyankar 10028ab949d8SShri Abhyankar 10033f1db9ecSBarry Smith /* 10043f1db9ecSBarry Smith This will not work with MatScalar == float because it calls the BLAS 10053f1db9ecSBarry Smith */ 10064a2ae208SSatish Balay #undef __FUNCT__ 10074a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_N" 1008dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_N(Mat A,Vec xx,Vec zz) 10092d61bbb3SSatish Balay { 10102d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1011d9ca1df4SBarry Smith PetscScalar *z = 0,*work,*workt,*zarray; 1012d9ca1df4SBarry Smith const PetscScalar *x,*xb; 1013d9ca1df4SBarry Smith const MatScalar *v; 1014dfbe8321SBarry Smith PetscErrorCode ierr; 1015d9ca1df4SBarry Smith PetscInt mbs,i,bs=A->rmap->bs,j,n,bs2=a->bs2; 1016d9ca1df4SBarry Smith const PetscInt *idx,*ii,*ridx=NULL; 1017d9ca1df4SBarry Smith PetscInt ncols,k; 1018ace3abfcSBarry Smith PetscBool usecprow=a->compressedrow.use; 10192d61bbb3SSatish Balay 10202d61bbb3SSatish Balay PetscFunctionBegin; 1021d9ca1df4SBarry Smith ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr); 102226e093fcSHong Zhang ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr); 10232d61bbb3SSatish Balay 10242d61bbb3SSatish Balay idx = a->j; 10252d61bbb3SSatish Balay v = a->a; 102626e093fcSHong Zhang if (usecprow) { 102726e093fcSHong Zhang mbs = a->compressedrow.nrows; 102826e093fcSHong Zhang ii = a->compressedrow.i; 10297b2bb3b9SHong Zhang ridx = a->compressedrow.rindex; 103026e093fcSHong Zhang } else { 103126e093fcSHong Zhang mbs = a->mbs; 10322d61bbb3SSatish Balay ii = a->i; 103326e093fcSHong Zhang z = zarray; 103426e093fcSHong Zhang } 1035218c64b6SSatish Balay 10362d61bbb3SSatish Balay if (!a->mult_work) { 1037d0f46423SBarry Smith k = PetscMax(A->rmap->n,A->cmap->n); 1038854ce69bSBarry Smith ierr = PetscMalloc1(k+1,&a->mult_work);CHKERRQ(ierr); 10392d61bbb3SSatish Balay } 10402d61bbb3SSatish Balay work = a->mult_work; 10412d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 10422d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 10432d61bbb3SSatish Balay ncols = n*bs; 10442d61bbb3SSatish Balay workt = work; 10452d61bbb3SSatish Balay for (j=0; j<n; j++) { 10462d61bbb3SSatish Balay xb = x + bs*(*idx++); 10472d61bbb3SSatish Balay for (k=0; k<bs; k++) workt[k] = xb[k]; 10482d61bbb3SSatish Balay workt += bs; 10492d61bbb3SSatish Balay } 10507b2bb3b9SHong Zhang if (usecprow) z = zarray + bs*ridx[i]; 105196b95a6bSBarry Smith PetscKernel_w_gets_Ar_times_v(bs,ncols,work,v,z); 105271044d3cSBarry Smith /* BLASgemv_("N",&bs,&ncols,&_DOne,v,&bs,work,&_One,&_DZero,z,&_One); */ 10532d61bbb3SSatish Balay v += n*bs2; 105426e093fcSHong Zhang if (!usecprow) z += bs; 10552d61bbb3SSatish Balay } 1056d9ca1df4SBarry Smith ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr); 105726e093fcSHong Zhang ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr); 10587c565772SBarry Smith ierr = PetscLogFlops(2.0*a->nz*bs2 - bs*a->nonzerorowcnt);CHKERRQ(ierr); 10592d61bbb3SSatish Balay PetscFunctionReturn(0); 10602d61bbb3SSatish Balay } 10612d61bbb3SSatish Balay 10624a2ae208SSatish Balay #undef __FUNCT__ 10634a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_1" 1064dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_1(Mat A,Vec xx,Vec yy,Vec zz) 10652d61bbb3SSatish Balay { 10662d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1067122f12eaSBarry Smith const PetscScalar *x; 1068122f12eaSBarry Smith PetscScalar *y,*z,sum; 1069122f12eaSBarry Smith const MatScalar *v; 1070dfbe8321SBarry Smith PetscErrorCode ierr; 10717c565772SBarry Smith PetscInt mbs=a->mbs,i,n,*ridx=NULL; 1072122f12eaSBarry Smith const PetscInt *idx,*ii; 1073ace3abfcSBarry Smith PetscBool usecprow=a->compressedrow.use; 10742d61bbb3SSatish Balay 10752d61bbb3SSatish Balay PetscFunctionBegin; 10763649974fSBarry Smith ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr); 1077d9ca1df4SBarry Smith ierr = VecGetArrayPair(yy,zz,&y,&z);CHKERRQ(ierr); 10782d61bbb3SSatish Balay 10792d61bbb3SSatish Balay idx = a->j; 10802d61bbb3SSatish Balay v = a->a; 108126e093fcSHong Zhang if (usecprow) { 10824eb6d288SHong Zhang if (zz != yy) { 1083122f12eaSBarry Smith ierr = PetscMemcpy(z,y,mbs*sizeof(PetscScalar));CHKERRQ(ierr); 10844eb6d288SHong Zhang } 108526e093fcSHong Zhang mbs = a->compressedrow.nrows; 108626e093fcSHong Zhang ii = a->compressedrow.i; 10877b2bb3b9SHong Zhang ridx = a->compressedrow.rindex; 108826e093fcSHong Zhang } else { 10892d61bbb3SSatish Balay ii = a->i; 109026e093fcSHong Zhang } 10912d61bbb3SSatish Balay 10922d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 1093122f12eaSBarry Smith n = ii[1] - ii[0]; 1094122f12eaSBarry Smith ii++; 109526e093fcSHong Zhang if (!usecprow) { 1096122f12eaSBarry Smith sum = y[i]; 1097122f12eaSBarry Smith } else { 1098122f12eaSBarry Smith sum = y[ridx[i]]; 1099122f12eaSBarry Smith } 1100444d8c10SJed Brown PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA); /* Indices for the next row (assumes same size as this one) */ 1101444d8c10SJed Brown PetscPrefetchBlock(v+n,n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */ 1102122f12eaSBarry Smith PetscSparseDensePlusDot(sum,x,v,idx,n); 1103122f12eaSBarry Smith v += n; 1104122f12eaSBarry Smith idx += n; 1105122f12eaSBarry Smith if (usecprow) { 1106122f12eaSBarry Smith z[ridx[i]] = sum; 1107122f12eaSBarry Smith } else { 1108122f12eaSBarry Smith z[i] = sum; 110926e093fcSHong Zhang } 11102d61bbb3SSatish Balay } 11113649974fSBarry Smith ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr); 1112d9ca1df4SBarry Smith ierr = VecRestoreArrayPair(yy,zz,&y,&z);CHKERRQ(ierr); 11137c565772SBarry Smith ierr = PetscLogFlops(2.0*a->nz);CHKERRQ(ierr); 11142d61bbb3SSatish Balay PetscFunctionReturn(0); 11152d61bbb3SSatish Balay } 11162d61bbb3SSatish Balay 11174a2ae208SSatish Balay #undef __FUNCT__ 11184a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_2" 1119dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_2(Mat A,Vec xx,Vec yy,Vec zz) 11202d61bbb3SSatish Balay { 11212d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1122d9ca1df4SBarry Smith PetscScalar *y = 0,*z = 0,sum1,sum2; 1123d9ca1df4SBarry Smith const PetscScalar *x,*xb; 112426e093fcSHong Zhang PetscScalar x1,x2,*yarray,*zarray; 1125d9ca1df4SBarry Smith const MatScalar *v; 1126dfbe8321SBarry Smith PetscErrorCode ierr; 1127d9ca1df4SBarry Smith PetscInt mbs = a->mbs,i,n,j; 1128d9ca1df4SBarry Smith const PetscInt *idx,*ii,*ridx = NULL; 1129ace3abfcSBarry Smith PetscBool usecprow = a->compressedrow.use; 11302d61bbb3SSatish Balay 11312d61bbb3SSatish Balay PetscFunctionBegin; 1132d9ca1df4SBarry Smith ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr); 1133d9ca1df4SBarry Smith ierr = VecGetArrayPair(yy,zz,&yarray,&zarray);CHKERRQ(ierr); 11342d61bbb3SSatish Balay 11352d61bbb3SSatish Balay idx = a->j; 11362d61bbb3SSatish Balay v = a->a; 113726e093fcSHong Zhang if (usecprow) { 11384eb6d288SHong Zhang if (zz != yy) { 11394eb6d288SHong Zhang ierr = PetscMemcpy(zarray,yarray,2*mbs*sizeof(PetscScalar));CHKERRQ(ierr); 11404eb6d288SHong Zhang } 114126e093fcSHong Zhang mbs = a->compressedrow.nrows; 114226e093fcSHong Zhang ii = a->compressedrow.i; 11437b2bb3b9SHong Zhang ridx = a->compressedrow.rindex; 11444eb6d288SHong Zhang if (zz != yy) { 11454eb6d288SHong Zhang ierr = PetscMemcpy(zarray,yarray,a->mbs*sizeof(PetscScalar));CHKERRQ(ierr); 11464eb6d288SHong Zhang } 114726e093fcSHong Zhang } else { 11482d61bbb3SSatish Balay ii = a->i; 114926e093fcSHong Zhang y = yarray; 115026e093fcSHong Zhang z = zarray; 115126e093fcSHong Zhang } 11522d61bbb3SSatish Balay 11532d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 11542d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 115526e093fcSHong Zhang if (usecprow) { 11567b2bb3b9SHong Zhang z = zarray + 2*ridx[i]; 11577b2bb3b9SHong Zhang y = yarray + 2*ridx[i]; 115826e093fcSHong Zhang } 11592d61bbb3SSatish Balay sum1 = y[0]; sum2 = y[1]; 1160444d8c10SJed Brown PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA); /* Indices for the next row (assumes same size as this one) */ 1161444d8c10SJed Brown PetscPrefetchBlock(v+4*n,4*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */ 11622d61bbb3SSatish Balay for (j=0; j<n; j++) { 116326fbe8dcSKarl Rupp xb = x + 2*(*idx++); 116426fbe8dcSKarl Rupp x1 = xb[0]; 116526fbe8dcSKarl Rupp x2 = xb[1]; 116626fbe8dcSKarl Rupp 11672d61bbb3SSatish Balay sum1 += v[0]*x1 + v[2]*x2; 11682d61bbb3SSatish Balay sum2 += v[1]*x1 + v[3]*x2; 11692d61bbb3SSatish Balay v += 4; 11702d61bbb3SSatish Balay } 11712d61bbb3SSatish Balay z[0] = sum1; z[1] = sum2; 117226e093fcSHong Zhang if (!usecprow) { 11732d61bbb3SSatish Balay z += 2; y += 2; 11742d61bbb3SSatish Balay } 117526e093fcSHong Zhang } 1176d9ca1df4SBarry Smith ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr); 1177d9ca1df4SBarry Smith ierr = VecRestoreArrayPair(yy,zz,&yarray,&zarray);CHKERRQ(ierr); 1178dc0b31edSSatish Balay ierr = PetscLogFlops(4.0*a->nz);CHKERRQ(ierr); 11792d61bbb3SSatish Balay PetscFunctionReturn(0); 11802d61bbb3SSatish Balay } 11812d61bbb3SSatish Balay 11824a2ae208SSatish Balay #undef __FUNCT__ 11834a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_3" 1184dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_3(Mat A,Vec xx,Vec yy,Vec zz) 11852d61bbb3SSatish Balay { 11862d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1187d9ca1df4SBarry Smith PetscScalar *y = 0,*z = 0,sum1,sum2,sum3,x1,x2,x3,*yarray,*zarray; 1188d9ca1df4SBarry Smith const PetscScalar *x,*xb; 1189d9ca1df4SBarry Smith const MatScalar *v; 1190dfbe8321SBarry Smith PetscErrorCode ierr; 1191d9ca1df4SBarry Smith PetscInt mbs = a->mbs,i,j,n; 1192d9ca1df4SBarry Smith const PetscInt *idx,*ii,*ridx = NULL; 1193ace3abfcSBarry Smith PetscBool usecprow = a->compressedrow.use; 11942d61bbb3SSatish Balay 11952d61bbb3SSatish Balay PetscFunctionBegin; 1196d9ca1df4SBarry Smith ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr); 1197d9ca1df4SBarry Smith ierr = VecGetArrayPair(yy,zz,&yarray,&zarray);CHKERRQ(ierr); 11982d61bbb3SSatish Balay 11992d61bbb3SSatish Balay idx = a->j; 12002d61bbb3SSatish Balay v = a->a; 120126e093fcSHong Zhang if (usecprow) { 12024eb6d288SHong Zhang if (zz != yy) { 12034eb6d288SHong Zhang ierr = PetscMemcpy(zarray,yarray,3*mbs*sizeof(PetscScalar));CHKERRQ(ierr); 12044eb6d288SHong Zhang } 120526e093fcSHong Zhang mbs = a->compressedrow.nrows; 120626e093fcSHong Zhang ii = a->compressedrow.i; 12077b2bb3b9SHong Zhang ridx = a->compressedrow.rindex; 120826e093fcSHong Zhang } else { 12092d61bbb3SSatish Balay ii = a->i; 121026e093fcSHong Zhang y = yarray; 121126e093fcSHong Zhang z = zarray; 121226e093fcSHong Zhang } 12132d61bbb3SSatish Balay 12142d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 12152d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 121626e093fcSHong Zhang if (usecprow) { 12177b2bb3b9SHong Zhang z = zarray + 3*ridx[i]; 12187b2bb3b9SHong Zhang y = yarray + 3*ridx[i]; 121926e093fcSHong Zhang } 12202d61bbb3SSatish Balay sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; 1221444d8c10SJed Brown PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA); /* Indices for the next row (assumes same size as this one) */ 1222444d8c10SJed Brown PetscPrefetchBlock(v+9*n,9*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */ 12232d61bbb3SSatish Balay for (j=0; j<n; j++) { 12242d61bbb3SSatish Balay xb = x + 3*(*idx++); x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; 12252d61bbb3SSatish Balay sum1 += v[0]*x1 + v[3]*x2 + v[6]*x3; 12262d61bbb3SSatish Balay sum2 += v[1]*x1 + v[4]*x2 + v[7]*x3; 12272d61bbb3SSatish Balay sum3 += v[2]*x1 + v[5]*x2 + v[8]*x3; 12282d61bbb3SSatish Balay v += 9; 12292d61bbb3SSatish Balay } 12302d61bbb3SSatish Balay z[0] = sum1; z[1] = sum2; z[2] = sum3; 123126e093fcSHong Zhang if (!usecprow) { 12322d61bbb3SSatish Balay z += 3; y += 3; 12332d61bbb3SSatish Balay } 123426e093fcSHong Zhang } 1235d9ca1df4SBarry Smith ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr); 1236d9ca1df4SBarry Smith ierr = VecRestoreArrayPair(yy,zz,&yarray,&zarray);CHKERRQ(ierr); 1237dc0b31edSSatish Balay ierr = PetscLogFlops(18.0*a->nz);CHKERRQ(ierr); 12382d61bbb3SSatish Balay PetscFunctionReturn(0); 12392d61bbb3SSatish Balay } 12402d61bbb3SSatish Balay 12414a2ae208SSatish Balay #undef __FUNCT__ 12424a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_4" 1243dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_4(Mat A,Vec xx,Vec yy,Vec zz) 12442d61bbb3SSatish Balay { 12452d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1246d9ca1df4SBarry Smith PetscScalar *y = 0,*z = 0,sum1,sum2,sum3,sum4,x1,x2,x3,x4,*yarray,*zarray; 1247d9ca1df4SBarry Smith const PetscScalar *x,*xb; 1248d9ca1df4SBarry Smith const MatScalar *v; 1249dfbe8321SBarry Smith PetscErrorCode ierr; 1250d9ca1df4SBarry Smith PetscInt mbs = a->mbs,i,j,n; 1251d9ca1df4SBarry Smith const PetscInt *idx,*ii,*ridx=NULL; 1252ace3abfcSBarry Smith PetscBool usecprow=a->compressedrow.use; 12532d61bbb3SSatish Balay 12542d61bbb3SSatish Balay PetscFunctionBegin; 1255d9ca1df4SBarry Smith ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr); 1256d9ca1df4SBarry Smith ierr = VecGetArrayPair(yy,zz,&yarray,&zarray);CHKERRQ(ierr); 12572d61bbb3SSatish Balay 12582d61bbb3SSatish Balay idx = a->j; 12592d61bbb3SSatish Balay v = a->a; 126026e093fcSHong Zhang if (usecprow) { 12614eb6d288SHong Zhang if (zz != yy) { 12624eb6d288SHong Zhang ierr = PetscMemcpy(zarray,yarray,4*mbs*sizeof(PetscScalar));CHKERRQ(ierr); 12634eb6d288SHong Zhang } 126426e093fcSHong Zhang mbs = a->compressedrow.nrows; 126526e093fcSHong Zhang ii = a->compressedrow.i; 12667b2bb3b9SHong Zhang ridx = a->compressedrow.rindex; 126726e093fcSHong Zhang } else { 12682d61bbb3SSatish Balay ii = a->i; 126926e093fcSHong Zhang y = yarray; 127026e093fcSHong Zhang z = zarray; 127126e093fcSHong Zhang } 12722d61bbb3SSatish Balay 12732d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 12742d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 127526e093fcSHong Zhang if (usecprow) { 12767b2bb3b9SHong Zhang z = zarray + 4*ridx[i]; 12777b2bb3b9SHong Zhang y = yarray + 4*ridx[i]; 127826e093fcSHong Zhang } 12792d61bbb3SSatish Balay sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; 1280444d8c10SJed Brown PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA); /* Indices for the next row (assumes same size as this one) */ 1281444d8c10SJed Brown PetscPrefetchBlock(v+16*n,16*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */ 12822d61bbb3SSatish Balay for (j=0; j<n; j++) { 12832d61bbb3SSatish Balay xb = x + 4*(*idx++); 12842d61bbb3SSatish Balay x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; 12852d61bbb3SSatish Balay sum1 += v[0]*x1 + v[4]*x2 + v[8]*x3 + v[12]*x4; 12862d61bbb3SSatish Balay sum2 += v[1]*x1 + v[5]*x2 + v[9]*x3 + v[13]*x4; 12872d61bbb3SSatish Balay sum3 += v[2]*x1 + v[6]*x2 + v[10]*x3 + v[14]*x4; 12882d61bbb3SSatish Balay sum4 += v[3]*x1 + v[7]*x2 + v[11]*x3 + v[15]*x4; 12892d61bbb3SSatish Balay v += 16; 12902d61bbb3SSatish Balay } 12912d61bbb3SSatish Balay z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; 129226e093fcSHong Zhang if (!usecprow) { 12932d61bbb3SSatish Balay z += 4; y += 4; 12942d61bbb3SSatish Balay } 129526e093fcSHong Zhang } 1296d9ca1df4SBarry Smith ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr); 1297d9ca1df4SBarry Smith ierr = VecRestoreArrayPair(yy,zz,&yarray,&zarray);CHKERRQ(ierr); 1298dc0b31edSSatish Balay ierr = PetscLogFlops(32.0*a->nz);CHKERRQ(ierr); 12992d61bbb3SSatish Balay PetscFunctionReturn(0); 13002d61bbb3SSatish Balay } 13012d61bbb3SSatish Balay 13024a2ae208SSatish Balay #undef __FUNCT__ 13034a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_5" 1304dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_5(Mat A,Vec xx,Vec yy,Vec zz) 13052d61bbb3SSatish Balay { 13062d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1307d9ca1df4SBarry Smith PetscScalar *y = 0,*z = 0,sum1,sum2,sum3,sum4,sum5,x1,x2,x3,x4,x5; 1308d9ca1df4SBarry Smith const PetscScalar *x,*xb; 130926e093fcSHong Zhang PetscScalar *yarray,*zarray; 1310d9ca1df4SBarry Smith const MatScalar *v; 1311dfbe8321SBarry Smith PetscErrorCode ierr; 1312d9ca1df4SBarry Smith PetscInt mbs = a->mbs,i,j,n; 1313d9ca1df4SBarry Smith const PetscInt *idx,*ii,*ridx = NULL; 1314ace3abfcSBarry Smith PetscBool usecprow=a->compressedrow.use; 13152d61bbb3SSatish Balay 13162d61bbb3SSatish Balay PetscFunctionBegin; 1317d9ca1df4SBarry Smith ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr); 1318d9ca1df4SBarry Smith ierr = VecGetArrayPair(yy,zz,&yarray,&zarray);CHKERRQ(ierr); 13192d61bbb3SSatish Balay 13202d61bbb3SSatish Balay idx = a->j; 13212d61bbb3SSatish Balay v = a->a; 132226e093fcSHong Zhang if (usecprow) { 13234eb6d288SHong Zhang if (zz != yy) { 13244eb6d288SHong Zhang ierr = PetscMemcpy(zarray,yarray,5*mbs*sizeof(PetscScalar));CHKERRQ(ierr); 13254eb6d288SHong Zhang } 132626e093fcSHong Zhang mbs = a->compressedrow.nrows; 132726e093fcSHong Zhang ii = a->compressedrow.i; 13287b2bb3b9SHong Zhang ridx = a->compressedrow.rindex; 132926e093fcSHong Zhang } else { 13302d61bbb3SSatish Balay ii = a->i; 133126e093fcSHong Zhang y = yarray; 133226e093fcSHong Zhang z = zarray; 133326e093fcSHong Zhang } 13342d61bbb3SSatish Balay 13352d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 13362d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 133726e093fcSHong Zhang if (usecprow) { 13387b2bb3b9SHong Zhang z = zarray + 5*ridx[i]; 13397b2bb3b9SHong Zhang y = yarray + 5*ridx[i]; 134026e093fcSHong Zhang } 13412d61bbb3SSatish Balay sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4]; 1342444d8c10SJed Brown PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA); /* Indices for the next row (assumes same size as this one) */ 1343444d8c10SJed Brown PetscPrefetchBlock(v+25*n,25*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */ 13442d61bbb3SSatish Balay for (j=0; j<n; j++) { 13452d61bbb3SSatish Balay xb = x + 5*(*idx++); 13462d61bbb3SSatish Balay x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; 13472d61bbb3SSatish Balay sum1 += v[0]*x1 + v[5]*x2 + v[10]*x3 + v[15]*x4 + v[20]*x5; 13482d61bbb3SSatish Balay sum2 += v[1]*x1 + v[6]*x2 + v[11]*x3 + v[16]*x4 + v[21]*x5; 13492d61bbb3SSatish Balay sum3 += v[2]*x1 + v[7]*x2 + v[12]*x3 + v[17]*x4 + v[22]*x5; 13502d61bbb3SSatish Balay sum4 += v[3]*x1 + v[8]*x2 + v[13]*x3 + v[18]*x4 + v[23]*x5; 13512d61bbb3SSatish Balay sum5 += v[4]*x1 + v[9]*x2 + v[14]*x3 + v[19]*x4 + v[24]*x5; 13522d61bbb3SSatish Balay v += 25; 13532d61bbb3SSatish Balay } 13542d61bbb3SSatish Balay z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; 135526e093fcSHong Zhang if (!usecprow) { 13562d61bbb3SSatish Balay z += 5; y += 5; 13572d61bbb3SSatish Balay } 135826e093fcSHong Zhang } 1359d9ca1df4SBarry Smith ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr); 1360d9ca1df4SBarry Smith ierr = VecRestoreArrayPair(yy,zz,&yarray,&zarray);CHKERRQ(ierr); 1361dc0b31edSSatish Balay ierr = PetscLogFlops(50.0*a->nz);CHKERRQ(ierr); 13622d61bbb3SSatish Balay PetscFunctionReturn(0); 13632d61bbb3SSatish Balay } 13644a2ae208SSatish Balay #undef __FUNCT__ 13654a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_6" 1366dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_6(Mat A,Vec xx,Vec yy,Vec zz) 136715091d37SBarry Smith { 136815091d37SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1369d9ca1df4SBarry Smith PetscScalar *y = 0,*z = 0,sum1,sum2,sum3,sum4,sum5,sum6; 1370d9ca1df4SBarry Smith const PetscScalar *x,*xb; 137126e093fcSHong Zhang PetscScalar x1,x2,x3,x4,x5,x6,*yarray,*zarray; 1372d9ca1df4SBarry Smith const MatScalar *v; 1373dfbe8321SBarry Smith PetscErrorCode ierr; 1374d9ca1df4SBarry Smith PetscInt mbs = a->mbs,i,j,n; 1375d9ca1df4SBarry Smith const PetscInt *idx,*ii,*ridx=NULL; 1376ace3abfcSBarry Smith PetscBool usecprow=a->compressedrow.use; 137715091d37SBarry Smith 137815091d37SBarry Smith PetscFunctionBegin; 1379d9ca1df4SBarry Smith ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr); 1380d9ca1df4SBarry Smith ierr = VecGetArrayPair(yy,zz,&yarray,&zarray);CHKERRQ(ierr); 138115091d37SBarry Smith 138215091d37SBarry Smith idx = a->j; 138315091d37SBarry Smith v = a->a; 138426e093fcSHong Zhang if (usecprow) { 13854eb6d288SHong Zhang if (zz != yy) { 13864eb6d288SHong Zhang ierr = PetscMemcpy(zarray,yarray,6*mbs*sizeof(PetscScalar));CHKERRQ(ierr); 13874eb6d288SHong Zhang } 138826e093fcSHong Zhang mbs = a->compressedrow.nrows; 138926e093fcSHong Zhang ii = a->compressedrow.i; 13907b2bb3b9SHong Zhang ridx = a->compressedrow.rindex; 139126e093fcSHong Zhang } else { 139215091d37SBarry Smith ii = a->i; 139326e093fcSHong Zhang y = yarray; 139426e093fcSHong Zhang z = zarray; 139526e093fcSHong Zhang } 139615091d37SBarry Smith 139715091d37SBarry Smith for (i=0; i<mbs; i++) { 139815091d37SBarry Smith n = ii[1] - ii[0]; ii++; 139926e093fcSHong Zhang if (usecprow) { 14007b2bb3b9SHong Zhang z = zarray + 6*ridx[i]; 14017b2bb3b9SHong Zhang y = yarray + 6*ridx[i]; 140226e093fcSHong Zhang } 140315091d37SBarry Smith sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4]; sum6 = y[5]; 1404444d8c10SJed Brown PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA); /* Indices for the next row (assumes same size as this one) */ 1405444d8c10SJed Brown PetscPrefetchBlock(v+36*n,36*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */ 140615091d37SBarry Smith for (j=0; j<n; j++) { 14073b95cb0eSSatish Balay xb = x + 6*(*idx++); 140815091d37SBarry Smith x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; 140915091d37SBarry Smith sum1 += v[0]*x1 + v[6]*x2 + v[12]*x3 + v[18]*x4 + v[24]*x5 + v[30]*x6; 141015091d37SBarry Smith sum2 += v[1]*x1 + v[7]*x2 + v[13]*x3 + v[19]*x4 + v[25]*x5 + v[31]*x6; 141115091d37SBarry Smith sum3 += v[2]*x1 + v[8]*x2 + v[14]*x3 + v[20]*x4 + v[26]*x5 + v[32]*x6; 141215091d37SBarry Smith sum4 += v[3]*x1 + v[9]*x2 + v[15]*x3 + v[21]*x4 + v[27]*x5 + v[33]*x6; 141315091d37SBarry Smith sum5 += v[4]*x1 + v[10]*x2 + v[16]*x3 + v[22]*x4 + v[28]*x5 + v[34]*x6; 141415091d37SBarry Smith sum6 += v[5]*x1 + v[11]*x2 + v[17]*x3 + v[23]*x4 + v[29]*x5 + v[35]*x6; 141515091d37SBarry Smith v += 36; 141615091d37SBarry Smith } 141715091d37SBarry Smith z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; 141826e093fcSHong Zhang if (!usecprow) { 141915091d37SBarry Smith z += 6; y += 6; 142015091d37SBarry Smith } 142126e093fcSHong Zhang } 1422d9ca1df4SBarry Smith ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr); 1423d9ca1df4SBarry Smith ierr = VecRestoreArrayPair(yy,zz,&yarray,&zarray);CHKERRQ(ierr); 1424dc0b31edSSatish Balay ierr = PetscLogFlops(72.0*a->nz);CHKERRQ(ierr); 142515091d37SBarry Smith PetscFunctionReturn(0); 142615091d37SBarry Smith } 14272d61bbb3SSatish Balay 14284a2ae208SSatish Balay #undef __FUNCT__ 14294a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_7" 1430dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_7(Mat A,Vec xx,Vec yy,Vec zz) 14312d61bbb3SSatish Balay { 14322d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1433d9ca1df4SBarry Smith PetscScalar *y = 0,*z = 0,sum1,sum2,sum3,sum4,sum5,sum6,sum7; 1434d9ca1df4SBarry Smith const PetscScalar *x,*xb; 143526e093fcSHong Zhang PetscScalar x1,x2,x3,x4,x5,x6,x7,*yarray,*zarray; 1436d9ca1df4SBarry Smith const MatScalar *v; 1437dfbe8321SBarry Smith PetscErrorCode ierr; 1438d9ca1df4SBarry Smith PetscInt mbs = a->mbs,i,j,n; 1439d9ca1df4SBarry Smith const PetscInt *idx,*ii,*ridx = NULL; 1440ace3abfcSBarry Smith PetscBool usecprow=a->compressedrow.use; 14412d61bbb3SSatish Balay 14422d61bbb3SSatish Balay PetscFunctionBegin; 1443d9ca1df4SBarry Smith ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr); 1444d9ca1df4SBarry Smith ierr = VecGetArrayPair(yy,zz,&yarray,&zarray);CHKERRQ(ierr); 14452d61bbb3SSatish Balay 14462d61bbb3SSatish Balay idx = a->j; 14472d61bbb3SSatish Balay v = a->a; 144826e093fcSHong Zhang if (usecprow) { 14494eb6d288SHong Zhang if (zz != yy) { 14504eb6d288SHong Zhang ierr = PetscMemcpy(zarray,yarray,7*mbs*sizeof(PetscScalar));CHKERRQ(ierr); 14514eb6d288SHong Zhang } 145226e093fcSHong Zhang mbs = a->compressedrow.nrows; 145326e093fcSHong Zhang ii = a->compressedrow.i; 14547b2bb3b9SHong Zhang ridx = a->compressedrow.rindex; 145526e093fcSHong Zhang } else { 14562d61bbb3SSatish Balay ii = a->i; 145726e093fcSHong Zhang y = yarray; 145826e093fcSHong Zhang z = zarray; 145926e093fcSHong Zhang } 14602d61bbb3SSatish Balay 14612d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 14622d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 146326e093fcSHong Zhang if (usecprow) { 14647b2bb3b9SHong Zhang z = zarray + 7*ridx[i]; 14657b2bb3b9SHong Zhang y = yarray + 7*ridx[i]; 146626e093fcSHong Zhang } 14672d61bbb3SSatish Balay sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4]; sum6 = y[5]; sum7 = y[6]; 1468444d8c10SJed Brown PetscPrefetchBlock(idx+n,n,0,PETSC_PREFETCH_HINT_NTA); /* Indices for the next row (assumes same size as this one) */ 1469444d8c10SJed Brown PetscPrefetchBlock(v+49*n,49*n,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */ 14702d61bbb3SSatish Balay for (j=0; j<n; j++) { 14712d61bbb3SSatish Balay xb = x + 7*(*idx++); 14722d61bbb3SSatish Balay x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6]; 14732d61bbb3SSatish Balay sum1 += v[0]*x1 + v[7]*x2 + v[14]*x3 + v[21]*x4 + v[28]*x5 + v[35]*x6 + v[42]*x7; 14742d61bbb3SSatish Balay sum2 += v[1]*x1 + v[8]*x2 + v[15]*x3 + v[22]*x4 + v[29]*x5 + v[36]*x6 + v[43]*x7; 14752d61bbb3SSatish Balay sum3 += v[2]*x1 + v[9]*x2 + v[16]*x3 + v[23]*x4 + v[30]*x5 + v[37]*x6 + v[44]*x7; 14762d61bbb3SSatish Balay sum4 += v[3]*x1 + v[10]*x2 + v[17]*x3 + v[24]*x4 + v[31]*x5 + v[38]*x6 + v[45]*x7; 14772d61bbb3SSatish Balay sum5 += v[4]*x1 + v[11]*x2 + v[18]*x3 + v[25]*x4 + v[32]*x5 + v[39]*x6 + v[46]*x7; 14782d61bbb3SSatish Balay sum6 += v[5]*x1 + v[12]*x2 + v[19]*x3 + v[26]*x4 + v[33]*x5 + v[40]*x6 + v[47]*x7; 14792d61bbb3SSatish Balay sum7 += v[6]*x1 + v[13]*x2 + v[20]*x3 + v[27]*x4 + v[34]*x5 + v[41]*x6 + v[48]*x7; 14802d61bbb3SSatish Balay v += 49; 14812d61bbb3SSatish Balay } 14822d61bbb3SSatish Balay z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7; 148326e093fcSHong Zhang if (!usecprow) { 14842d61bbb3SSatish Balay z += 7; y += 7; 14852d61bbb3SSatish Balay } 148626e093fcSHong Zhang } 1487d9ca1df4SBarry Smith ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr); 1488d9ca1df4SBarry Smith ierr = VecRestoreArrayPair(yy,zz,&yarray,&zarray);CHKERRQ(ierr); 1489dc0b31edSSatish Balay ierr = PetscLogFlops(98.0*a->nz);CHKERRQ(ierr); 14902d61bbb3SSatish Balay PetscFunctionReturn(0); 14912d61bbb3SSatish Balay } 1492218c64b6SSatish Balay 14934a2ae208SSatish Balay #undef __FUNCT__ 14944a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_N" 1495dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_N(Mat A,Vec xx,Vec yy,Vec zz) 14962d61bbb3SSatish Balay { 14972d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1498d9ca1df4SBarry Smith PetscScalar *z = 0,*work,*workt,*zarray; 1499d9ca1df4SBarry Smith const PetscScalar *x,*xb; 1500d9ca1df4SBarry Smith const MatScalar *v; 15016849ba73SBarry Smith PetscErrorCode ierr; 1502d9ca1df4SBarry Smith PetscInt mbs,i,bs=A->rmap->bs,j,n,bs2=a->bs2; 1503d9ca1df4SBarry Smith PetscInt ncols,k; 1504d9ca1df4SBarry Smith const PetscInt *ridx = NULL,*idx,*ii; 1505ace3abfcSBarry Smith PetscBool usecprow = a->compressedrow.use; 1506218c64b6SSatish Balay 15072d61bbb3SSatish Balay PetscFunctionBegin; 1508343551c4SBarry Smith ierr = VecCopy(yy,zz);CHKERRQ(ierr); 1509d9ca1df4SBarry Smith ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr); 151026e093fcSHong Zhang ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr); 15112d61bbb3SSatish Balay 15122d61bbb3SSatish Balay idx = a->j; 15132d61bbb3SSatish Balay v = a->a; 151426e093fcSHong Zhang if (usecprow) { 151526e093fcSHong Zhang mbs = a->compressedrow.nrows; 151626e093fcSHong Zhang ii = a->compressedrow.i; 15177b2bb3b9SHong Zhang ridx = a->compressedrow.rindex; 151826e093fcSHong Zhang } else { 151926e093fcSHong Zhang mbs = a->mbs; 15202d61bbb3SSatish Balay ii = a->i; 152126e093fcSHong Zhang z = zarray; 152226e093fcSHong Zhang } 15232d61bbb3SSatish Balay 15242d61bbb3SSatish Balay if (!a->mult_work) { 1525d0f46423SBarry Smith k = PetscMax(A->rmap->n,A->cmap->n); 1526854ce69bSBarry Smith ierr = PetscMalloc1(k+1,&a->mult_work);CHKERRQ(ierr); 15272d61bbb3SSatish Balay } 15282d61bbb3SSatish Balay work = a->mult_work; 15292d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 15302d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 15312d61bbb3SSatish Balay ncols = n*bs; 15322d61bbb3SSatish Balay workt = work; 15332d61bbb3SSatish Balay for (j=0; j<n; j++) { 15342d61bbb3SSatish Balay xb = x + bs*(*idx++); 15352d61bbb3SSatish Balay for (k=0; k<bs; k++) workt[k] = xb[k]; 15362d61bbb3SSatish Balay workt += bs; 15372d61bbb3SSatish Balay } 15387b2bb3b9SHong Zhang if (usecprow) z = zarray + bs*ridx[i]; 153996b95a6bSBarry Smith PetscKernel_w_gets_w_plus_Ar_times_v(bs,ncols,work,v,z); 154071044d3cSBarry Smith /* BLASgemv_("N",&bs,&ncols,&_DOne,v,&bs,work,&_One,&_DOne,z,&_One); */ 15412d61bbb3SSatish Balay v += n*bs2; 154226fbe8dcSKarl Rupp if (!usecprow) z += bs; 154326e093fcSHong Zhang } 1544d9ca1df4SBarry Smith ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr); 154526e093fcSHong Zhang ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr); 1546dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*a->nz*bs2);CHKERRQ(ierr); 15472d61bbb3SSatish Balay PetscFunctionReturn(0); 15482d61bbb3SSatish Balay } 15492d61bbb3SSatish Balay 15504a2ae208SSatish Balay #undef __FUNCT__ 1551547795f9SHong Zhang #define __FUNCT__ "MatMultHermitianTranspose_SeqBAIJ" 1552547795f9SHong Zhang PetscErrorCode MatMultHermitianTranspose_SeqBAIJ(Mat A,Vec xx,Vec zz) 1553547795f9SHong Zhang { 1554547795f9SHong Zhang PetscScalar zero = 0.0; 1555547795f9SHong Zhang PetscErrorCode ierr; 1556547795f9SHong Zhang 1557547795f9SHong Zhang PetscFunctionBegin; 1558547795f9SHong Zhang ierr = VecSet(zz,zero);CHKERRQ(ierr); 1559547795f9SHong Zhang ierr = MatMultHermitianTransposeAdd_SeqBAIJ(A,xx,zz,zz);CHKERRQ(ierr); 1560547795f9SHong Zhang PetscFunctionReturn(0); 1561547795f9SHong Zhang } 1562547795f9SHong Zhang 1563547795f9SHong Zhang #undef __FUNCT__ 15644a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_SeqBAIJ" 1565dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_SeqBAIJ(Mat A,Vec xx,Vec zz) 15662d61bbb3SSatish Balay { 15673447b6efSHong Zhang PetscScalar zero = 0.0; 15686849ba73SBarry Smith PetscErrorCode ierr; 15692d61bbb3SSatish Balay 15702d61bbb3SSatish Balay PetscFunctionBegin; 15712dcb1b2aSMatthew Knepley ierr = VecSet(zz,zero);CHKERRQ(ierr); 15723447b6efSHong Zhang ierr = MatMultTransposeAdd_SeqBAIJ(A,xx,zz,zz);CHKERRQ(ierr); 15732d61bbb3SSatish Balay PetscFunctionReturn(0); 15742d61bbb3SSatish Balay } 15752d61bbb3SSatish Balay 15764a2ae208SSatish Balay #undef __FUNCT__ 1577547795f9SHong Zhang #define __FUNCT__ "MatMultHermitianTransposeAdd_SeqBAIJ" 1578547795f9SHong Zhang PetscErrorCode MatMultHermitianTransposeAdd_SeqBAIJ(Mat A,Vec xx,Vec yy,Vec zz) 1579547795f9SHong Zhang { 1580547795f9SHong Zhang Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1581b8c08b77SHong Zhang PetscScalar *z,x1,x2,x3,x4,x5; 1582d9ca1df4SBarry Smith const PetscScalar *x,*xb = NULL; 1583d9ca1df4SBarry Smith const MatScalar *v; 1584547795f9SHong Zhang PetscErrorCode ierr; 1585b8c08b77SHong Zhang PetscInt mbs,i,rval,bs=A->rmap->bs,j,n; 1586d9ca1df4SBarry Smith const PetscInt *idx,*ii,*ib,*ridx = NULL; 1587547795f9SHong Zhang Mat_CompressedRow cprow = a->compressedrow; 1588ace3abfcSBarry Smith PetscBool usecprow = cprow.use; 1589547795f9SHong Zhang 1590547795f9SHong Zhang PetscFunctionBegin; 1591343551c4SBarry Smith if (yy != zz) { ierr = VecCopy(yy,zz);CHKERRQ(ierr); } 1592d9ca1df4SBarry Smith ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr); 1593547795f9SHong Zhang ierr = VecGetArray(zz,&z);CHKERRQ(ierr); 1594547795f9SHong Zhang 1595547795f9SHong Zhang idx = a->j; 1596547795f9SHong Zhang v = a->a; 1597547795f9SHong Zhang if (usecprow) { 1598547795f9SHong Zhang mbs = cprow.nrows; 1599547795f9SHong Zhang ii = cprow.i; 1600547795f9SHong Zhang ridx = cprow.rindex; 1601547795f9SHong Zhang } else { 1602547795f9SHong Zhang mbs=a->mbs; 1603547795f9SHong Zhang ii = a->i; 1604547795f9SHong Zhang xb = x; 1605547795f9SHong Zhang } 1606547795f9SHong Zhang 1607547795f9SHong Zhang switch (bs) { 1608547795f9SHong Zhang case 1: 1609547795f9SHong Zhang for (i=0; i<mbs; i++) { 1610547795f9SHong Zhang if (usecprow) xb = x + ridx[i]; 1611547795f9SHong Zhang x1 = xb[0]; 1612547795f9SHong Zhang ib = idx + ii[0]; 1613547795f9SHong Zhang n = ii[1] - ii[0]; ii++; 1614547795f9SHong Zhang for (j=0; j<n; j++) { 1615547795f9SHong Zhang rval = ib[j]; 1616547795f9SHong Zhang z[rval] += PetscConj(*v) * x1; 1617547795f9SHong Zhang v++; 1618547795f9SHong Zhang } 1619547795f9SHong Zhang if (!usecprow) xb++; 1620547795f9SHong Zhang } 1621547795f9SHong Zhang break; 1622547795f9SHong Zhang case 2: 1623547795f9SHong Zhang for (i=0; i<mbs; i++) { 1624547795f9SHong Zhang if (usecprow) xb = x + 2*ridx[i]; 1625547795f9SHong Zhang x1 = xb[0]; x2 = xb[1]; 1626547795f9SHong Zhang ib = idx + ii[0]; 1627547795f9SHong Zhang n = ii[1] - ii[0]; ii++; 1628547795f9SHong Zhang for (j=0; j<n; j++) { 1629547795f9SHong Zhang rval = ib[j]*2; 1630547795f9SHong Zhang z[rval++] += PetscConj(v[0])*x1 + PetscConj(v[1])*x2; 1631547795f9SHong Zhang z[rval++] += PetscConj(v[2])*x1 + PetscConj(v[3])*x2; 1632547795f9SHong Zhang v += 4; 1633547795f9SHong Zhang } 1634547795f9SHong Zhang if (!usecprow) xb += 2; 1635547795f9SHong Zhang } 1636547795f9SHong Zhang break; 1637547795f9SHong Zhang case 3: 1638547795f9SHong Zhang for (i=0; i<mbs; i++) { 1639547795f9SHong Zhang if (usecprow) xb = x + 3*ridx[i]; 1640547795f9SHong Zhang x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; 1641547795f9SHong Zhang ib = idx + ii[0]; 1642547795f9SHong Zhang n = ii[1] - ii[0]; ii++; 1643547795f9SHong Zhang for (j=0; j<n; j++) { 1644547795f9SHong Zhang rval = ib[j]*3; 1645547795f9SHong Zhang z[rval++] += PetscConj(v[0])*x1 + PetscConj(v[1])*x2 + PetscConj(v[2])*x3; 1646547795f9SHong Zhang z[rval++] += PetscConj(v[3])*x1 + PetscConj(v[4])*x2 + PetscConj(v[5])*x3; 1647547795f9SHong Zhang z[rval++] += PetscConj(v[6])*x1 + PetscConj(v[7])*x2 + PetscConj(v[8])*x3; 1648547795f9SHong Zhang v += 9; 1649547795f9SHong Zhang } 1650547795f9SHong Zhang if (!usecprow) xb += 3; 1651547795f9SHong Zhang } 1652547795f9SHong Zhang break; 1653547795f9SHong Zhang case 4: 1654547795f9SHong Zhang for (i=0; i<mbs; i++) { 1655547795f9SHong Zhang if (usecprow) xb = x + 4*ridx[i]; 1656547795f9SHong Zhang x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; 1657547795f9SHong Zhang ib = idx + ii[0]; 1658547795f9SHong Zhang n = ii[1] - ii[0]; ii++; 1659547795f9SHong Zhang for (j=0; j<n; j++) { 1660547795f9SHong Zhang rval = ib[j]*4; 1661547795f9SHong Zhang z[rval++] += PetscConj(v[0])*x1 + PetscConj(v[1])*x2 + PetscConj(v[2])*x3 + PetscConj(v[3])*x4; 1662547795f9SHong Zhang z[rval++] += PetscConj(v[4])*x1 + PetscConj(v[5])*x2 + PetscConj(v[6])*x3 + PetscConj(v[7])*x4; 1663547795f9SHong Zhang z[rval++] += PetscConj(v[8])*x1 + PetscConj(v[9])*x2 + PetscConj(v[10])*x3 + PetscConj(v[11])*x4; 1664547795f9SHong Zhang z[rval++] += PetscConj(v[12])*x1 + PetscConj(v[13])*x2 + PetscConj(v[14])*x3 + PetscConj(v[15])*x4; 1665547795f9SHong Zhang v += 16; 1666547795f9SHong Zhang } 1667547795f9SHong Zhang if (!usecprow) xb += 4; 1668547795f9SHong Zhang } 1669547795f9SHong Zhang break; 1670547795f9SHong Zhang case 5: 1671547795f9SHong Zhang for (i=0; i<mbs; i++) { 1672547795f9SHong Zhang if (usecprow) xb = x + 5*ridx[i]; 1673547795f9SHong Zhang x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; 1674547795f9SHong Zhang x4 = xb[3]; x5 = xb[4]; 1675547795f9SHong Zhang ib = idx + ii[0]; 1676547795f9SHong Zhang n = ii[1] - ii[0]; ii++; 1677547795f9SHong Zhang for (j=0; j<n; j++) { 1678547795f9SHong Zhang rval = ib[j]*5; 1679547795f9SHong Zhang z[rval++] += PetscConj(v[0])*x1 + PetscConj(v[1])*x2 + PetscConj(v[2])*x3 + PetscConj(v[3])*x4 + PetscConj(v[4])*x5; 1680547795f9SHong Zhang z[rval++] += PetscConj(v[5])*x1 + PetscConj(v[6])*x2 + PetscConj(v[7])*x3 + PetscConj(v[8])*x4 + PetscConj(v[9])*x5; 1681547795f9SHong Zhang z[rval++] += PetscConj(v[10])*x1 + PetscConj(v[11])*x2 + PetscConj(v[12])*x3 + PetscConj(v[13])*x4 + PetscConj(v[14])*x5; 1682547795f9SHong Zhang z[rval++] += PetscConj(v[15])*x1 + PetscConj(v[16])*x2 + PetscConj(v[17])*x3 + PetscConj(v[18])*x4 + PetscConj(v[19])*x5; 1683547795f9SHong Zhang z[rval++] += PetscConj(v[20])*x1 + PetscConj(v[21])*x2 + PetscConj(v[22])*x3 + PetscConj(v[23])*x4 + PetscConj(v[24])*x5; 1684547795f9SHong Zhang v += 25; 1685547795f9SHong Zhang } 1686547795f9SHong Zhang if (!usecprow) xb += 5; 1687547795f9SHong Zhang } 1688547795f9SHong Zhang break; 1689968ae2c8SSatish Balay default: /* block sizes larger than 5 by 5 are handled by BLAS */ 1690968ae2c8SSatish Balay SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"block size larger than 5 is not supported yet"); 1691968ae2c8SSatish Balay #if 0 1692968ae2c8SSatish Balay { 1693b8c08b77SHong Zhang PetscInt ncols,k,bs2=a->bs2; 1694b8c08b77SHong Zhang PetscScalar *work,*workt,zb; 1695d9ca1df4SBarry Smith const PetscScalar *xtmp; 1696547795f9SHong Zhang if (!a->mult_work) { 1697547795f9SHong Zhang k = PetscMax(A->rmap->n,A->cmap->n); 1698854ce69bSBarry Smith ierr = PetscMalloc1(k+1,&a->mult_work);CHKERRQ(ierr); 1699547795f9SHong Zhang } 1700547795f9SHong Zhang work = a->mult_work; 1701547795f9SHong Zhang xtmp = x; 1702547795f9SHong Zhang for (i=0; i<mbs; i++) { 1703547795f9SHong Zhang n = ii[1] - ii[0]; ii++; 1704547795f9SHong Zhang ncols = n*bs; 1705547795f9SHong Zhang ierr = PetscMemzero(work,ncols*sizeof(PetscScalar));CHKERRQ(ierr); 170626fbe8dcSKarl Rupp if (usecprow) xtmp = x + bs*ridx[i]; 170796b95a6bSBarry Smith PetscKernel_w_gets_w_plus_trans_Ar_times_v(bs,ncols,xtmp,v,work); 1708547795f9SHong Zhang /* BLASgemv_("T",&bs,&ncols,&_DOne,v,&bs,xtmp,&_One,&_DOne,work,&_One); */ 1709547795f9SHong Zhang v += n*bs2; 1710547795f9SHong Zhang if (!usecprow) xtmp += bs; 1711547795f9SHong Zhang workt = work; 1712547795f9SHong Zhang for (j=0; j<n; j++) { 1713547795f9SHong Zhang zb = z + bs*(*idx++); 1714547795f9SHong Zhang for (k=0; k<bs; k++) zb[k] += workt[k] ; 1715547795f9SHong Zhang workt += bs; 1716547795f9SHong Zhang } 1717547795f9SHong Zhang } 1718547795f9SHong Zhang } 1719968ae2c8SSatish Balay #endif 1720547795f9SHong Zhang } 1721d9ca1df4SBarry Smith ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr); 1722547795f9SHong Zhang ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr); 1723547795f9SHong Zhang ierr = PetscLogFlops(2.0*a->nz*a->bs2);CHKERRQ(ierr); 1724547795f9SHong Zhang PetscFunctionReturn(0); 1725547795f9SHong Zhang } 1726547795f9SHong Zhang 1727547795f9SHong Zhang #undef __FUNCT__ 17284a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_SeqBAIJ" 1729dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_SeqBAIJ(Mat A,Vec xx,Vec yy,Vec zz) 17302d61bbb3SSatish Balay { 17312d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1732d9ca1df4SBarry Smith PetscScalar *zb,*z,x1,x2,x3,x4,x5; 1733d9ca1df4SBarry Smith const PetscScalar *x,*xb = 0; 1734d9ca1df4SBarry Smith const MatScalar *v; 17356849ba73SBarry Smith PetscErrorCode ierr; 1736d9ca1df4SBarry Smith PetscInt mbs,i,rval,bs=A->rmap->bs,j,n,bs2=a->bs2; 1737d9ca1df4SBarry Smith const PetscInt *idx,*ii,*ib,*ridx = NULL; 17383447b6efSHong Zhang Mat_CompressedRow cprow = a->compressedrow; 1739ace3abfcSBarry Smith PetscBool usecprow=cprow.use; 17402d61bbb3SSatish Balay 17412d61bbb3SSatish Balay PetscFunctionBegin; 1742343551c4SBarry Smith if (yy != zz) { ierr = VecCopy(yy,zz);CHKERRQ(ierr); } 1743d9ca1df4SBarry Smith ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr); 17443447b6efSHong Zhang ierr = VecGetArray(zz,&z);CHKERRQ(ierr); 17452d61bbb3SSatish Balay 17462d61bbb3SSatish Balay idx = a->j; 17472d61bbb3SSatish Balay v = a->a; 17483447b6efSHong Zhang if (usecprow) { 17493447b6efSHong Zhang mbs = cprow.nrows; 17503447b6efSHong Zhang ii = cprow.i; 17517b2bb3b9SHong Zhang ridx = cprow.rindex; 17523447b6efSHong Zhang } else { 17533447b6efSHong Zhang mbs=a->mbs; 17542d61bbb3SSatish Balay ii = a->i; 1755f1af5d2fSBarry Smith xb = x; 17563447b6efSHong Zhang } 17572d61bbb3SSatish Balay 17582d61bbb3SSatish Balay switch (bs) { 17592d61bbb3SSatish Balay case 1: 17602d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 17617b2bb3b9SHong Zhang if (usecprow) xb = x + ridx[i]; 1762f1af5d2fSBarry Smith x1 = xb[0]; 17633447b6efSHong Zhang ib = idx + ii[0]; 17643447b6efSHong Zhang n = ii[1] - ii[0]; ii++; 17652d61bbb3SSatish Balay for (j=0; j<n; j++) { 17662d61bbb3SSatish Balay rval = ib[j]; 1767f1af5d2fSBarry Smith z[rval] += *v * x1; 1768f1af5d2fSBarry Smith v++; 17692d61bbb3SSatish Balay } 17703447b6efSHong Zhang if (!usecprow) xb++; 17712d61bbb3SSatish Balay } 17722d61bbb3SSatish Balay break; 17732d61bbb3SSatish Balay case 2: 17742d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 17757b2bb3b9SHong Zhang if (usecprow) xb = x + 2*ridx[i]; 1776f1af5d2fSBarry Smith x1 = xb[0]; x2 = xb[1]; 17773447b6efSHong Zhang ib = idx + ii[0]; 17783447b6efSHong Zhang n = ii[1] - ii[0]; ii++; 17792d61bbb3SSatish Balay for (j=0; j<n; j++) { 17802d61bbb3SSatish Balay rval = ib[j]*2; 17812d61bbb3SSatish Balay z[rval++] += v[0]*x1 + v[1]*x2; 17822d61bbb3SSatish Balay z[rval++] += v[2]*x1 + v[3]*x2; 17832d61bbb3SSatish Balay v += 4; 17842d61bbb3SSatish Balay } 17853447b6efSHong Zhang if (!usecprow) xb += 2; 17862d61bbb3SSatish Balay } 17872d61bbb3SSatish Balay break; 17882d61bbb3SSatish Balay case 3: 17892d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 17907b2bb3b9SHong Zhang if (usecprow) xb = x + 3*ridx[i]; 1791f1af5d2fSBarry Smith x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; 17923447b6efSHong Zhang ib = idx + ii[0]; 17933447b6efSHong Zhang n = ii[1] - ii[0]; ii++; 17942d61bbb3SSatish Balay for (j=0; j<n; j++) { 17952d61bbb3SSatish Balay rval = ib[j]*3; 17962d61bbb3SSatish Balay z[rval++] += v[0]*x1 + v[1]*x2 + v[2]*x3; 17972d61bbb3SSatish Balay z[rval++] += v[3]*x1 + v[4]*x2 + v[5]*x3; 17982d61bbb3SSatish Balay z[rval++] += v[6]*x1 + v[7]*x2 + v[8]*x3; 17992d61bbb3SSatish Balay v += 9; 18002d61bbb3SSatish Balay } 18013447b6efSHong Zhang if (!usecprow) xb += 3; 18022d61bbb3SSatish Balay } 18032d61bbb3SSatish Balay break; 18042d61bbb3SSatish Balay case 4: 18052d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 18067b2bb3b9SHong Zhang if (usecprow) xb = x + 4*ridx[i]; 1807f1af5d2fSBarry Smith x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; 18083447b6efSHong Zhang ib = idx + ii[0]; 18093447b6efSHong Zhang n = ii[1] - ii[0]; ii++; 18102d61bbb3SSatish Balay for (j=0; j<n; j++) { 18112d61bbb3SSatish Balay rval = ib[j]*4; 18122d61bbb3SSatish Balay z[rval++] += v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4; 18132d61bbb3SSatish Balay z[rval++] += v[4]*x1 + v[5]*x2 + v[6]*x3 + v[7]*x4; 18142d61bbb3SSatish Balay z[rval++] += v[8]*x1 + v[9]*x2 + v[10]*x3 + v[11]*x4; 18152d61bbb3SSatish Balay z[rval++] += v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4; 18162d61bbb3SSatish Balay v += 16; 18172d61bbb3SSatish Balay } 18183447b6efSHong Zhang if (!usecprow) xb += 4; 18192d61bbb3SSatish Balay } 18202d61bbb3SSatish Balay break; 18212d61bbb3SSatish Balay case 5: 18222d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 18237b2bb3b9SHong Zhang if (usecprow) xb = x + 5*ridx[i]; 1824f1af5d2fSBarry Smith x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; 18252d61bbb3SSatish Balay x4 = xb[3]; x5 = xb[4]; 18263447b6efSHong Zhang ib = idx + ii[0]; 18273447b6efSHong Zhang n = ii[1] - ii[0]; ii++; 18282d61bbb3SSatish Balay for (j=0; j<n; j++) { 18292d61bbb3SSatish Balay rval = ib[j]*5; 18302d61bbb3SSatish Balay z[rval++] += v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4 + v[4]*x5; 18312d61bbb3SSatish Balay z[rval++] += v[5]*x1 + v[6]*x2 + v[7]*x3 + v[8]*x4 + v[9]*x5; 18322d61bbb3SSatish Balay z[rval++] += v[10]*x1 + v[11]*x2 + v[12]*x3 + v[13]*x4 + v[14]*x5; 18332d61bbb3SSatish Balay z[rval++] += v[15]*x1 + v[16]*x2 + v[17]*x3 + v[18]*x4 + v[19]*x5; 18342d61bbb3SSatish Balay z[rval++] += v[20]*x1 + v[21]*x2 + v[22]*x3 + v[23]*x4 + v[24]*x5; 18352d61bbb3SSatish Balay v += 25; 18362d61bbb3SSatish Balay } 18373447b6efSHong Zhang if (!usecprow) xb += 5; 18382d61bbb3SSatish Balay } 18392d61bbb3SSatish Balay break; 1840f1af5d2fSBarry Smith default: { /* block sizes larger then 5 by 5 are handled by BLAS */ 1841690b6cddSBarry Smith PetscInt ncols,k; 1842d9ca1df4SBarry Smith PetscScalar *work,*workt; 1843d9ca1df4SBarry Smith const PetscScalar *xtmp; 18442d61bbb3SSatish Balay if (!a->mult_work) { 1845d0f46423SBarry Smith k = PetscMax(A->rmap->n,A->cmap->n); 1846854ce69bSBarry Smith ierr = PetscMalloc1(k+1,&a->mult_work);CHKERRQ(ierr); 18472d61bbb3SSatish Balay } 18482d61bbb3SSatish Balay work = a->mult_work; 18493447b6efSHong Zhang xtmp = x; 18502d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 18512d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 18522d61bbb3SSatish Balay ncols = n*bs; 185387828ca2SBarry Smith ierr = PetscMemzero(work,ncols*sizeof(PetscScalar));CHKERRQ(ierr); 185426fbe8dcSKarl Rupp if (usecprow) xtmp = x + bs*ridx[i]; 185596b95a6bSBarry Smith PetscKernel_w_gets_w_plus_trans_Ar_times_v(bs,ncols,xtmp,v,work); 185671044d3cSBarry Smith /* BLASgemv_("T",&bs,&ncols,&_DOne,v,&bs,xtmp,&_One,&_DOne,work,&_One); */ 18572d61bbb3SSatish Balay v += n*bs2; 18583447b6efSHong Zhang if (!usecprow) xtmp += bs; 18592d61bbb3SSatish Balay workt = work; 18602d61bbb3SSatish Balay for (j=0; j<n; j++) { 18612d61bbb3SSatish Balay zb = z + bs*(*idx++); 18622d61bbb3SSatish Balay for (k=0; k<bs; k++) zb[k] += workt[k]; 18632d61bbb3SSatish Balay workt += bs; 18642d61bbb3SSatish Balay } 18652d61bbb3SSatish Balay } 18662d61bbb3SSatish Balay } 18672d61bbb3SSatish Balay } 1868d9ca1df4SBarry Smith ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr); 18693447b6efSHong Zhang ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr); 1870dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*a->nz*a->bs2);CHKERRQ(ierr); 18712d61bbb3SSatish Balay PetscFunctionReturn(0); 18722d61bbb3SSatish Balay } 18732d61bbb3SSatish Balay 18744a2ae208SSatish Balay #undef __FUNCT__ 18754a2ae208SSatish Balay #define __FUNCT__ "MatScale_SeqBAIJ" 1876f4df32b1SMatthew Knepley PetscErrorCode MatScale_SeqBAIJ(Mat inA,PetscScalar alpha) 18772d61bbb3SSatish Balay { 18782d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)inA->data; 1879690b6cddSBarry Smith PetscInt totalnz = a->bs2*a->nz; 1880f4df32b1SMatthew Knepley PetscScalar oalpha = alpha; 1881efee365bSSatish Balay PetscErrorCode ierr; 1882c5df96a5SBarry Smith PetscBLASInt one = 1,tnz; 18832d61bbb3SSatish Balay 18842d61bbb3SSatish Balay PetscFunctionBegin; 1885c5df96a5SBarry Smith ierr = PetscBLASIntCast(totalnz,&tnz);CHKERRQ(ierr); 18868b83055fSJed Brown PetscStackCallBLAS("BLASscal",BLASscal_(&tnz,&oalpha,a->a,&one)); 1887efee365bSSatish Balay ierr = PetscLogFlops(totalnz);CHKERRQ(ierr); 18882d61bbb3SSatish Balay PetscFunctionReturn(0); 18892d61bbb3SSatish Balay } 18902d61bbb3SSatish Balay 18914a2ae208SSatish Balay #undef __FUNCT__ 18924a2ae208SSatish Balay #define __FUNCT__ "MatNorm_SeqBAIJ" 1893dfbe8321SBarry Smith PetscErrorCode MatNorm_SeqBAIJ(Mat A,NormType type,PetscReal *norm) 18942d61bbb3SSatish Balay { 18958a62d963SHong Zhang PetscErrorCode ierr; 18962d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 18973f1db9ecSBarry Smith MatScalar *v = a->a; 1898329f5518SBarry Smith PetscReal sum = 0.0; 1899d0f46423SBarry Smith PetscInt i,j,k,bs=A->rmap->bs,nz=a->nz,bs2=a->bs2,k1; 19002d61bbb3SSatish Balay 19012d61bbb3SSatish Balay PetscFunctionBegin; 19022d61bbb3SSatish Balay if (type == NORM_FROBENIUS) { 19032d61bbb3SSatish Balay for (i=0; i< bs2*nz; i++) { 1904329f5518SBarry Smith sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 19052d61bbb3SSatish Balay } 19068f1a2a5eSBarry Smith *norm = PetscSqrtReal(sum); 190751f70360SJed Brown ierr = PetscLogFlops(2*bs2*nz);CHKERRQ(ierr); 19088a62d963SHong Zhang } else if (type == NORM_1) { /* maximum column sum */ 19098a62d963SHong Zhang PetscReal *tmp; 19108a62d963SHong Zhang PetscInt *bcol = a->j; 1911854ce69bSBarry Smith ierr = PetscCalloc1(A->cmap->n+1,&tmp);CHKERRQ(ierr); 19128a62d963SHong Zhang for (i=0; i<nz; i++) { 19138a62d963SHong Zhang for (j=0; j<bs; j++) { 19148a62d963SHong Zhang k1 = bs*(*bcol) + j; /* column index */ 19158a62d963SHong Zhang for (k=0; k<bs; k++) { 19168a62d963SHong Zhang tmp[k1] += PetscAbsScalar(*v); v++; 19178a62d963SHong Zhang } 19188a62d963SHong Zhang } 19198a62d963SHong Zhang bcol++; 19208a62d963SHong Zhang } 19218a62d963SHong Zhang *norm = 0.0; 1922d0f46423SBarry Smith for (j=0; j<A->cmap->n; j++) { 19238a62d963SHong Zhang if (tmp[j] > *norm) *norm = tmp[j]; 19248a62d963SHong Zhang } 19258a62d963SHong Zhang ierr = PetscFree(tmp);CHKERRQ(ierr); 192651f70360SJed Brown ierr = PetscLogFlops(PetscMax(bs2*nz-1,0));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 } 194251f70360SJed Brown ierr = PetscLogFlops(PetscMax(bs2*nz-1,0));CHKERRQ(ierr); 1943e7e72b3dSBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for this norm yet"); 19442d61bbb3SSatish Balay PetscFunctionReturn(0); 19452d61bbb3SSatish Balay } 19462d61bbb3SSatish Balay 19472d61bbb3SSatish Balay 19484a2ae208SSatish Balay #undef __FUNCT__ 19494a2ae208SSatish Balay #define __FUNCT__ "MatEqual_SeqBAIJ" 1950ace3abfcSBarry Smith PetscErrorCode MatEqual_SeqBAIJ(Mat A,Mat B,PetscBool * flg) 19512d61bbb3SSatish Balay { 19522d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)B->data; 1953dfbe8321SBarry Smith PetscErrorCode ierr; 19542d61bbb3SSatish Balay 19552d61bbb3SSatish Balay PetscFunctionBegin; 19562d61bbb3SSatish Balay /* If the matrix/block dimensions are not equal, or no of nonzeros or shift */ 1957d0f46423SBarry Smith if ((A->rmap->N != B->rmap->N) || (A->cmap->n != B->cmap->n) || (A->rmap->bs != B->rmap->bs)|| (a->nz != b->nz)) { 1958273d9f13SBarry Smith *flg = PETSC_FALSE; 1959273d9f13SBarry Smith PetscFunctionReturn(0); 19602d61bbb3SSatish Balay } 19612d61bbb3SSatish Balay 19622d61bbb3SSatish Balay /* if the a->i are the same */ 1963690b6cddSBarry Smith ierr = PetscMemcmp(a->i,b->i,(a->mbs+1)*sizeof(PetscInt),flg);CHKERRQ(ierr); 196426fbe8dcSKarl Rupp if (!*flg) PetscFunctionReturn(0); 19652d61bbb3SSatish Balay 19662d61bbb3SSatish Balay /* if a->j are the same */ 1967690b6cddSBarry Smith ierr = PetscMemcmp(a->j,b->j,(a->nz)*sizeof(PetscInt),flg);CHKERRQ(ierr); 196826fbe8dcSKarl Rupp if (!*flg) PetscFunctionReturn(0); 196926fbe8dcSKarl Rupp 19702d61bbb3SSatish Balay /* if a->a are the same */ 1971d0f46423SBarry Smith ierr = PetscMemcmp(a->a,b->a,(a->nz)*(A->rmap->bs)*(B->rmap->bs)*sizeof(PetscScalar),flg);CHKERRQ(ierr); 19722d61bbb3SSatish Balay PetscFunctionReturn(0); 19732d61bbb3SSatish Balay 19742d61bbb3SSatish Balay } 19752d61bbb3SSatish Balay 19764a2ae208SSatish Balay #undef __FUNCT__ 19774a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_SeqBAIJ" 1978dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_SeqBAIJ(Mat A,Vec v) 19792d61bbb3SSatish Balay { 19802d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1981dfbe8321SBarry Smith PetscErrorCode ierr; 1982690b6cddSBarry Smith PetscInt i,j,k,n,row,bs,*ai,*aj,ambs,bs2; 198387828ca2SBarry Smith PetscScalar *x,zero = 0.0; 19843f1db9ecSBarry Smith MatScalar *aa,*aa_j; 19852d61bbb3SSatish Balay 19862d61bbb3SSatish Balay PetscFunctionBegin; 1987e32f2f54SBarry Smith if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 1988d0f46423SBarry Smith bs = A->rmap->bs; 19892d61bbb3SSatish Balay aa = a->a; 19902d61bbb3SSatish Balay ai = a->i; 19912d61bbb3SSatish Balay aj = a->j; 19922d61bbb3SSatish Balay ambs = a->mbs; 19932d61bbb3SSatish Balay bs2 = a->bs2; 19942d61bbb3SSatish Balay 19952dcb1b2aSMatthew Knepley ierr = VecSet(v,zero);CHKERRQ(ierr); 19961ebc52fbSHong Zhang ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1997e1311b90SBarry Smith ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr); 1998e32f2f54SBarry Smith if (n != A->rmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector"); 19992d61bbb3SSatish Balay for (i=0; i<ambs; i++) { 20002d61bbb3SSatish Balay for (j=ai[i]; j<ai[i+1]; j++) { 20012d61bbb3SSatish Balay if (aj[j] == i) { 20022d61bbb3SSatish Balay row = i*bs; 20032d61bbb3SSatish Balay aa_j = aa+j*bs2; 20042d61bbb3SSatish Balay for (k=0; k<bs2; k+=(bs+1),row++) x[row] = aa_j[k]; 20052d61bbb3SSatish Balay break; 20062d61bbb3SSatish Balay } 20072d61bbb3SSatish Balay } 20082d61bbb3SSatish Balay } 20091ebc52fbSHong Zhang ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 20102d61bbb3SSatish Balay PetscFunctionReturn(0); 20112d61bbb3SSatish Balay } 20122d61bbb3SSatish Balay 20134a2ae208SSatish Balay #undef __FUNCT__ 20144a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_SeqBAIJ" 2015dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_SeqBAIJ(Mat A,Vec ll,Vec rr) 20162d61bbb3SSatish Balay { 20172d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 201853ef36baSBarry Smith const PetscScalar *l,*r,*li,*ri; 201953ef36baSBarry Smith PetscScalar x; 20203f1db9ecSBarry Smith MatScalar *aa, *v; 2021dfbe8321SBarry Smith PetscErrorCode ierr; 202253ef36baSBarry Smith PetscInt i,j,k,lm,rn,M,m,n,mbs,tmp,bs,bs2,iai; 202353ef36baSBarry Smith const PetscInt *ai,*aj; 20242d61bbb3SSatish Balay 20252d61bbb3SSatish Balay PetscFunctionBegin; 20262d61bbb3SSatish Balay ai = a->i; 20272d61bbb3SSatish Balay aj = a->j; 20282d61bbb3SSatish Balay aa = a->a; 2029d0f46423SBarry Smith m = A->rmap->n; 2030d0f46423SBarry Smith n = A->cmap->n; 2031d0f46423SBarry Smith bs = A->rmap->bs; 20322d61bbb3SSatish Balay mbs = a->mbs; 20332d61bbb3SSatish Balay bs2 = a->bs2; 20342d61bbb3SSatish Balay if (ll) { 203553ef36baSBarry Smith ierr = VecGetArrayRead(ll,&l);CHKERRQ(ierr); 2036e1311b90SBarry Smith ierr = VecGetLocalSize(ll,&lm);CHKERRQ(ierr); 2037e32f2f54SBarry Smith if (lm != m) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Left scaling vector wrong length"); 20382d61bbb3SSatish Balay for (i=0; i<mbs; i++) { /* for each block row */ 20392d61bbb3SSatish Balay M = ai[i+1] - ai[i]; 20402d61bbb3SSatish Balay li = l + i*bs; 20412d61bbb3SSatish Balay v = aa + bs2*ai[i]; 20422d61bbb3SSatish Balay for (j=0; j<M; j++) { /* for each block */ 20432d61bbb3SSatish Balay for (k=0; k<bs2; k++) { 20442d61bbb3SSatish Balay (*v++) *= li[k%bs]; 20452d61bbb3SSatish Balay } 20462d61bbb3SSatish Balay } 20472d61bbb3SSatish Balay } 204853ef36baSBarry Smith ierr = VecRestoreArrayRead(ll,&l);CHKERRQ(ierr); 2049efee365bSSatish Balay ierr = PetscLogFlops(a->nz);CHKERRQ(ierr); 20502d61bbb3SSatish Balay } 20512d61bbb3SSatish Balay 20522d61bbb3SSatish Balay if (rr) { 205353ef36baSBarry Smith ierr = VecGetArrayRead(rr,&r);CHKERRQ(ierr); 2054e1311b90SBarry Smith ierr = VecGetLocalSize(rr,&rn);CHKERRQ(ierr); 2055e32f2f54SBarry Smith if (rn != n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Right scaling vector wrong length"); 20562d61bbb3SSatish Balay for (i=0; i<mbs; i++) { /* for each block row */ 205753ef36baSBarry Smith iai = ai[i]; 205853ef36baSBarry Smith M = ai[i+1] - iai; 205953ef36baSBarry Smith v = aa + bs2*iai; 20602d61bbb3SSatish Balay for (j=0; j<M; j++) { /* for each block */ 206153ef36baSBarry Smith ri = r + bs*aj[iai+j]; 20622d61bbb3SSatish Balay for (k=0; k<bs; k++) { 20632d61bbb3SSatish Balay x = ri[k]; 206453ef36baSBarry Smith for (tmp=0; tmp<bs; tmp++) v[tmp] *= x; 206553ef36baSBarry Smith v += bs; 20662d61bbb3SSatish Balay } 20672d61bbb3SSatish Balay } 20682d61bbb3SSatish Balay } 206953ef36baSBarry Smith ierr = VecRestoreArrayRead(rr,&r);CHKERRQ(ierr); 2070efee365bSSatish Balay ierr = PetscLogFlops(a->nz);CHKERRQ(ierr); 20712d61bbb3SSatish Balay } 20722d61bbb3SSatish Balay PetscFunctionReturn(0); 20732d61bbb3SSatish Balay } 20742d61bbb3SSatish Balay 20752d61bbb3SSatish Balay 20764a2ae208SSatish Balay #undef __FUNCT__ 20774a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_SeqBAIJ" 2078dfbe8321SBarry Smith PetscErrorCode MatGetInfo_SeqBAIJ(Mat A,MatInfoType flag,MatInfo *info) 20792d61bbb3SSatish Balay { 20802d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 20812d61bbb3SSatish Balay 20822d61bbb3SSatish Balay PetscFunctionBegin; 20832d61bbb3SSatish Balay info->block_size = a->bs2; 2084ceed8ce5SJed Brown info->nz_allocated = a->bs2*a->maxnz; 20852d61bbb3SSatish Balay info->nz_used = a->bs2*a->nz; 20862d61bbb3SSatish Balay info->nz_unneeded = (double)(info->nz_allocated - info->nz_used); 20872d61bbb3SSatish Balay info->assemblies = A->num_ass; 20888e58a170SBarry Smith info->mallocs = A->info.mallocs; 20897adad957SLisandro Dalcin info->memory = ((PetscObject)A)->mem; 2090d5f3da31SBarry Smith if (A->factortype) { 20912d61bbb3SSatish Balay info->fill_ratio_given = A->info.fill_ratio_given; 20922d61bbb3SSatish Balay info->fill_ratio_needed = A->info.fill_ratio_needed; 20932d61bbb3SSatish Balay info->factor_mallocs = A->info.factor_mallocs; 20942d61bbb3SSatish Balay } else { 20952d61bbb3SSatish Balay info->fill_ratio_given = 0; 20962d61bbb3SSatish Balay info->fill_ratio_needed = 0; 20972d61bbb3SSatish Balay info->factor_mallocs = 0; 20982d61bbb3SSatish Balay } 20992d61bbb3SSatish Balay PetscFunctionReturn(0); 21002d61bbb3SSatish Balay } 21012d61bbb3SSatish Balay 21024a2ae208SSatish Balay #undef __FUNCT__ 21034a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_SeqBAIJ" 2104dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_SeqBAIJ(Mat A) 21052d61bbb3SSatish Balay { 21062d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 2107dfbe8321SBarry Smith PetscErrorCode ierr; 21082d61bbb3SSatish Balay 21092d61bbb3SSatish Balay PetscFunctionBegin; 2110549d3d68SSatish Balay ierr = PetscMemzero(a->a,a->bs2*a->i[a->mbs]*sizeof(MatScalar));CHKERRQ(ierr); 21112d61bbb3SSatish Balay PetscFunctionReturn(0); 21122d61bbb3SSatish Balay } 2113