1be1d678aSKris Buschelman #define PETSCMAT_DLL 2cac129eeSSatish Balay 37c4f633dSBarry Smith #include "../src/mat/impls/baij/seq/baij.h" 4c60f0209SBarry Smith #include "../src/mat/blockinvert.h" 50a835dfdSSatish Balay #include "petscbt.h" 6cac129eeSSatish Balay 74a2ae208SSatish Balay #undef __FUNCT__ 84a2ae208SSatish Balay #define __FUNCT__ "MatIncreaseOverlap_SeqBAIJ" 9690b6cddSBarry Smith PetscErrorCode MatIncreaseOverlap_SeqBAIJ(Mat A,PetscInt is_max,IS is[],PetscInt ov) 10a3192f15SSatish Balay { 11a3192f15SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 126849ba73SBarry Smith PetscErrorCode ierr; 135d0c19d7SBarry Smith PetscInt row,i,j,k,l,m,n,*nidx,isz,val,ival; 145d0c19d7SBarry Smith const PetscInt *idx; 15690b6cddSBarry Smith PetscInt start,end,*ai,*aj,bs,*nidx2; 16f1af5d2fSBarry Smith PetscBT table; 17a3192f15SSatish Balay 183a40ed3dSBarry Smith PetscFunctionBegin; 19a3192f15SSatish Balay m = a->mbs; 20a3192f15SSatish Balay ai = a->i; 21a3192f15SSatish Balay aj = a->j; 22d0f46423SBarry Smith bs = A->rmap->bs; 23a3192f15SSatish Balay 2429bbc08cSBarry Smith if (ov < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap specified"); 25a3192f15SSatish Balay 266831982aSBarry Smith ierr = PetscBTCreate(m,table);CHKERRQ(ierr); 27690b6cddSBarry Smith ierr = PetscMalloc((m+1)*sizeof(PetscInt),&nidx);CHKERRQ(ierr); 28d0f46423SBarry Smith ierr = PetscMalloc((A->rmap->N+1)*sizeof(PetscInt),&nidx2);CHKERRQ(ierr); 29a3192f15SSatish Balay 30a3192f15SSatish Balay for (i=0; i<is_max; i++) { 31a3192f15SSatish Balay /* Initialise the two local arrays */ 32a3192f15SSatish Balay isz = 0; 336831982aSBarry Smith ierr = PetscBTMemzero(m,table);CHKERRQ(ierr); 34a3192f15SSatish Balay 35a3192f15SSatish Balay /* Extract the indices, assume there can be duplicate entries */ 36a3192f15SSatish Balay ierr = ISGetIndices(is[i],&idx);CHKERRQ(ierr); 37b9b97703SBarry Smith ierr = ISGetLocalSize(is[i],&n);CHKERRQ(ierr); 38a3192f15SSatish Balay 39a3192f15SSatish Balay /* Enter these into the temp arrays i.e mark table[row], enter row into new index */ 40a3192f15SSatish Balay for (j=0; j<n ; ++j){ 41218c64b6SSatish Balay ival = idx[j]/bs; /* convert the indices into block indices */ 425eee224dSHong Zhang if (ival>=m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"index greater than mat-dim"); 43f1af5d2fSBarry Smith if(!PetscBTLookupSet(table,ival)) { nidx[isz++] = ival;} 44a3192f15SSatish Balay } 45a3192f15SSatish Balay ierr = ISRestoreIndices(is[i],&idx);CHKERRQ(ierr); 46a3192f15SSatish Balay ierr = ISDestroy(is[i]);CHKERRQ(ierr); 47a3192f15SSatish Balay 48a3192f15SSatish Balay k = 0; 49a3192f15SSatish Balay for (j=0; j<ov; j++){ /* for each overlap*/ 50a3192f15SSatish Balay n = isz; 51a3192f15SSatish Balay for (; k<n ; k++){ /* do only those rows in nidx[k], which are not done yet */ 52a3192f15SSatish Balay row = nidx[k]; 53a3192f15SSatish Balay start = ai[row]; 54a3192f15SSatish Balay end = ai[row+1]; 55a3192f15SSatish Balay for (l = start; l<end ; l++){ 56a3192f15SSatish Balay val = aj[l]; 57f1af5d2fSBarry Smith if (!PetscBTLookupSet(table,val)) {nidx[isz++] = val;} 58a3192f15SSatish Balay } 59a3192f15SSatish Balay } 60a3192f15SSatish Balay } 61218c64b6SSatish Balay /* expand the Index Set */ 62218c64b6SSatish Balay for (j=0; j<isz; j++) { 63218c64b6SSatish Balay for (k=0; k<bs; k++) 64218c64b6SSatish Balay nidx2[j*bs+k] = nidx[j]*bs+k; 65218c64b6SSatish Balay } 66329f5518SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,isz*bs,nidx2,is+i);CHKERRQ(ierr); 67a3192f15SSatish Balay } 686831982aSBarry Smith ierr = PetscBTDestroy(table);CHKERRQ(ierr); 69606d414cSSatish Balay ierr = PetscFree(nidx);CHKERRQ(ierr); 70606d414cSSatish Balay ierr = PetscFree(nidx2);CHKERRQ(ierr); 713a40ed3dSBarry Smith PetscFunctionReturn(0); 72a3192f15SSatish Balay } 731c351548SSatish Balay 744a2ae208SSatish Balay #undef __FUNCT__ 754a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqBAIJ_Private" 764aa3045dSJed Brown PetscErrorCode MatGetSubMatrix_SeqBAIJ_Private(Mat A,IS isrow,IS iscol,MatReuse scall,Mat *B) 77736121d4SSatish Balay { 78736121d4SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*c; 796849ba73SBarry Smith PetscErrorCode ierr; 80690b6cddSBarry Smith PetscInt *smap,i,k,kstart,kend,oldcols = a->nbs,*lens; 81690b6cddSBarry Smith PetscInt row,mat_i,*mat_j,tcol,*mat_ilen; 825d0c19d7SBarry Smith const PetscInt *irow,*icol; 835d0c19d7SBarry Smith PetscInt nrows,ncols,*ssmap,bs=A->rmap->bs,bs2=a->bs2; 84690b6cddSBarry Smith PetscInt *aj = a->j,*ai = a->i; 853f1db9ecSBarry Smith MatScalar *mat_a; 86736121d4SSatish Balay Mat C; 8714ca34e6SBarry Smith PetscTruth flag,sorted; 88736121d4SSatish Balay 893a40ed3dSBarry Smith PetscFunctionBegin; 9014ca34e6SBarry Smith ierr = ISSorted(iscol,&sorted);CHKERRQ(ierr); 9114ca34e6SBarry Smith if (!sorted) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"IS is not sorted"); 92736121d4SSatish Balay 93736121d4SSatish Balay ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr); 94218c64b6SSatish Balay ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr); 95b9b97703SBarry Smith ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr); 96b9b97703SBarry Smith ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr); 97736121d4SSatish Balay 98690b6cddSBarry Smith ierr = PetscMalloc((1+oldcols)*sizeof(PetscInt),&smap);CHKERRQ(ierr); 99736121d4SSatish Balay ssmap = smap; 100690b6cddSBarry Smith ierr = PetscMalloc((1+nrows)*sizeof(PetscInt),&lens);CHKERRQ(ierr); 101690b6cddSBarry Smith ierr = PetscMemzero(smap,oldcols*sizeof(PetscInt));CHKERRQ(ierr); 102736121d4SSatish Balay for (i=0; i<ncols; i++) smap[icol[i]] = i+1; 103736121d4SSatish Balay /* determine lens of each row */ 104736121d4SSatish Balay for (i=0; i<nrows; i++) { 105736121d4SSatish Balay kstart = ai[irow[i]]; 106736121d4SSatish Balay kend = kstart + a->ilen[irow[i]]; 107736121d4SSatish Balay lens[i] = 0; 108736121d4SSatish Balay for (k=kstart; k<kend; k++) { 109736121d4SSatish Balay if (ssmap[aj[k]]) { 110736121d4SSatish Balay lens[i]++; 111736121d4SSatish Balay } 112736121d4SSatish Balay } 113736121d4SSatish Balay } 114736121d4SSatish Balay /* Create and fill new matrix */ 115736121d4SSatish Balay if (scall == MAT_REUSE_MATRIX) { 116736121d4SSatish Balay c = (Mat_SeqBAIJ *)((*B)->data); 117736121d4SSatish Balay 118d0f46423SBarry Smith if (c->mbs!=nrows || c->nbs!=ncols || (*B)->rmap->bs!=bs) SETERRQ(PETSC_ERR_ARG_SIZ,"Submatrix wrong size"); 119690b6cddSBarry Smith ierr = PetscMemcmp(c->ilen,lens,c->mbs *sizeof(PetscInt),&flag);CHKERRQ(ierr); 120abc0a331SBarry Smith if (!flag) { 12129bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong no of nonzeros"); 122736121d4SSatish Balay } 123690b6cddSBarry Smith ierr = PetscMemzero(c->ilen,c->mbs*sizeof(PetscInt));CHKERRQ(ierr); 124736121d4SSatish Balay C = *B; 1253a40ed3dSBarry Smith } else { 1267adad957SLisandro Dalcin ierr = MatCreate(((PetscObject)A)->comm,&C);CHKERRQ(ierr); 127f69a0ea3SMatthew Knepley ierr = MatSetSizes(C,nrows*bs,ncols*bs,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 1287adad957SLisandro Dalcin ierr = MatSetType(C,((PetscObject)A)->type_name);CHKERRQ(ierr); 129ab93d7beSBarry Smith ierr = MatSeqBAIJSetPreallocation_SeqBAIJ(C,bs,0,lens);CHKERRQ(ierr); 130736121d4SSatish Balay } 131736121d4SSatish Balay c = (Mat_SeqBAIJ *)(C->data); 132736121d4SSatish Balay for (i=0; i<nrows; i++) { 133736121d4SSatish Balay row = irow[i]; 134736121d4SSatish Balay kstart = ai[row]; 135736121d4SSatish Balay kend = kstart + a->ilen[row]; 136736121d4SSatish Balay mat_i = c->i[i]; 137736121d4SSatish Balay mat_j = c->j + mat_i; 138218c64b6SSatish Balay mat_a = c->a + mat_i*bs2; 139736121d4SSatish Balay mat_ilen = c->ilen + i; 140736121d4SSatish Balay for (k=kstart; k<kend; k++) { 141736121d4SSatish Balay if ((tcol=ssmap[a->j[k]])) { 142736121d4SSatish Balay *mat_j++ = tcol - 1; 143549d3d68SSatish Balay ierr = PetscMemcpy(mat_a,a->a+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 144549d3d68SSatish Balay mat_a += bs2; 145736121d4SSatish Balay (*mat_ilen)++; 146736121d4SSatish Balay } 147736121d4SSatish Balay } 148736121d4SSatish Balay } 149218c64b6SSatish Balay 150736121d4SSatish Balay /* Free work space */ 151736121d4SSatish Balay ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr); 152606d414cSSatish Balay ierr = PetscFree(smap);CHKERRQ(ierr); 153606d414cSSatish Balay ierr = PetscFree(lens);CHKERRQ(ierr); 1546d4a8577SBarry Smith ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1556d4a8577SBarry Smith ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 156736121d4SSatish Balay 157736121d4SSatish Balay ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr); 158736121d4SSatish Balay *B = C; 1593a40ed3dSBarry Smith PetscFunctionReturn(0); 160736121d4SSatish Balay } 161736121d4SSatish Balay 1624a2ae208SSatish Balay #undef __FUNCT__ 1634a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqBAIJ" 1644aa3045dSJed Brown PetscErrorCode MatGetSubMatrix_SeqBAIJ(Mat A,IS isrow,IS iscol,MatReuse scall,Mat *B) 165218c64b6SSatish Balay { 166218c64b6SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 167218c64b6SSatish Balay IS is1,is2; 1686849ba73SBarry Smith PetscErrorCode ierr; 1695d0c19d7SBarry Smith PetscInt *vary,*iary,nrows,ncols,i,bs=A->rmap->bs,count; 1705d0c19d7SBarry Smith const PetscInt *irow,*icol; 171218c64b6SSatish Balay 1723a40ed3dSBarry Smith PetscFunctionBegin; 173218c64b6SSatish Balay ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr); 174218c64b6SSatish Balay ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr); 175b9b97703SBarry Smith ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr); 176b9b97703SBarry Smith ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr); 177218c64b6SSatish Balay 178218c64b6SSatish Balay /* Verify if the indices corespond to each element in a block 179218c64b6SSatish Balay and form the IS with compressed IS */ 180fca92195SBarry Smith ierr = PetscMalloc2(a->mbs,PetscInt,&vary,a->mbs,PetscInt,&iary);CHKERRQ(ierr); 181fca92195SBarry Smith ierr = PetscMemzero(vary,a->mbs*sizeof(PetscInt));CHKERRQ(ierr); 182218c64b6SSatish Balay for (i=0; i<nrows; i++) vary[irow[i]/bs]++; 183218c64b6SSatish Balay count = 0; 184218c64b6SSatish Balay for (i=0; i<a->mbs; i++) { 185634064b4SBarry Smith if (vary[i]!=0 && vary[i]!=bs) SETERRQ(PETSC_ERR_ARG_SIZ,"Index set does not match blocks"); 186218c64b6SSatish Balay if (vary[i]==bs) iary[count++] = i; 187218c64b6SSatish Balay } 188029af93fSBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,count,iary,&is1);CHKERRQ(ierr); 189218c64b6SSatish Balay 190690b6cddSBarry Smith ierr = PetscMemzero(vary,(a->mbs)*sizeof(PetscInt));CHKERRQ(ierr); 191218c64b6SSatish Balay for (i=0; i<ncols; i++) vary[icol[i]/bs]++; 192218c64b6SSatish Balay count = 0; 193218c64b6SSatish Balay for (i=0; i<a->mbs; i++) { 194634064b4SBarry Smith if (vary[i]!=0 && vary[i]!=bs) SETERRQ(PETSC_ERR_PLIB,"Internal error in PETSc"); 195218c64b6SSatish Balay if (vary[i]==bs) iary[count++] = i; 196218c64b6SSatish Balay } 197029af93fSBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,count,iary,&is2);CHKERRQ(ierr); 198218c64b6SSatish Balay ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr); 199218c64b6SSatish Balay ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr); 200fca92195SBarry Smith ierr = PetscFree2(vary,iary);CHKERRQ(ierr); 201218c64b6SSatish Balay 2024aa3045dSJed Brown ierr = MatGetSubMatrix_SeqBAIJ_Private(A,is1,is2,scall,B);CHKERRQ(ierr); 203218c64b6SSatish Balay ISDestroy(is1); 204218c64b6SSatish Balay ISDestroy(is2); 2053a40ed3dSBarry Smith PetscFunctionReturn(0); 206218c64b6SSatish Balay } 207218c64b6SSatish Balay 2084a2ae208SSatish Balay #undef __FUNCT__ 2094a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_SeqBAIJ" 210690b6cddSBarry Smith PetscErrorCode MatGetSubMatrices_SeqBAIJ(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[]) 211736121d4SSatish Balay { 2126849ba73SBarry Smith PetscErrorCode ierr; 213690b6cddSBarry Smith PetscInt i; 214736121d4SSatish Balay 2153a40ed3dSBarry Smith PetscFunctionBegin; 216736121d4SSatish Balay if (scall == MAT_INITIAL_MATRIX) { 21782502324SSatish Balay ierr = PetscMalloc((n+1)*sizeof(Mat),B);CHKERRQ(ierr); 218736121d4SSatish Balay } 219736121d4SSatish Balay 220736121d4SSatish Balay for (i=0; i<n; i++) { 2214aa3045dSJed Brown ierr = MatGetSubMatrix_SeqBAIJ(A,irow[i],icol[i],scall,&(*B)[i]);CHKERRQ(ierr); 222736121d4SSatish Balay } 2233a40ed3dSBarry Smith PetscFunctionReturn(0); 224736121d4SSatish Balay } 225218c64b6SSatish Balay 226218c64b6SSatish Balay 2272d61bbb3SSatish Balay /* -------------------------------------------------------*/ 2282d61bbb3SSatish Balay /* Should check that shapes of vectors and matrices match */ 2292d61bbb3SSatish Balay /* -------------------------------------------------------*/ 230d9eff348SSatish Balay #include "petscblaslapack.h" 2312d61bbb3SSatish Balay 2324a2ae208SSatish Balay #undef __FUNCT__ 2334a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_1" 234dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_1(Mat A,Vec xx,Vec zz) 2352d61bbb3SSatish Balay { 2362d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 237d9fead3dSBarry Smith PetscScalar *z,sum; 238d9fead3dSBarry Smith const PetscScalar *x; 239d9fead3dSBarry Smith const MatScalar *v; 2406849ba73SBarry Smith PetscErrorCode ierr; 2412162cab8SBarry Smith PetscInt mbs,i,n,nonzerorow=0; 2422162cab8SBarry Smith const PetscInt *idx,*ii,*ridx=PETSC_NULL; 24326e093fcSHong Zhang PetscTruth usecprow=a->compressedrow.use; 2442d61bbb3SSatish Balay 2452d61bbb3SSatish Balay PetscFunctionBegin; 246d9fead3dSBarry Smith ierr = VecGetArray(xx,(PetscScalar**)&x);CHKERRQ(ierr); 2471ebc52fbSHong Zhang ierr = VecGetArray(zz,&z);CHKERRQ(ierr); 2482d61bbb3SSatish Balay 24926e093fcSHong Zhang if (usecprow){ 25026e093fcSHong Zhang mbs = a->compressedrow.nrows; 25126e093fcSHong Zhang ii = a->compressedrow.i; 2527b2bb3b9SHong Zhang ridx = a->compressedrow.rindex; 253a1c3900fSBarry Smith ierr = PetscMemzero(z,mbs*sizeof(PetscScalar));CHKERRQ(ierr); 25426e093fcSHong Zhang } else { 25526e093fcSHong Zhang mbs = a->mbs; 2562d61bbb3SSatish Balay ii = a->i; 25726e093fcSHong Zhang } 2582d61bbb3SSatish Balay 2592d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 260ee54c7eeSHong Zhang n = ii[1] - ii[0]; 261ee54c7eeSHong Zhang v = a->a + ii[0]; 262ee54c7eeSHong Zhang idx = a->j + ii[0]; 263ee54c7eeSHong Zhang ii++; 2642d61bbb3SSatish Balay sum = 0.0; 2652162cab8SBarry Smith PetscSparseDensePlusDot(sum,x,v,idx,n); 26626e093fcSHong Zhang if (usecprow){ 2677b2bb3b9SHong Zhang z[ridx[i]] = sum; 26826e093fcSHong Zhang } else { 269122f12eaSBarry Smith nonzerorow += (n>0); 2702d61bbb3SSatish Balay z[i] = sum; 2712d61bbb3SSatish Balay } 27226e093fcSHong Zhang } 273d9fead3dSBarry Smith ierr = VecRestoreArray(xx,(PetscScalar**)&x);CHKERRQ(ierr); 2741ebc52fbSHong Zhang ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr); 275dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*a->nz - nonzerorow);CHKERRQ(ierr); 2762d61bbb3SSatish Balay PetscFunctionReturn(0); 2772d61bbb3SSatish Balay } 2782d61bbb3SSatish Balay 2794a2ae208SSatish Balay #undef __FUNCT__ 2804a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_2" 281dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_2(Mat A,Vec xx,Vec zz) 2822d61bbb3SSatish Balay { 2832d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 284d9fead3dSBarry Smith PetscScalar *z = 0,sum1,sum2,*zarray; 285d9fead3dSBarry Smith const PetscScalar *x,*xb; 28687828ca2SBarry Smith PetscScalar x1,x2; 287d9fead3dSBarry Smith const MatScalar *v; 288dfbe8321SBarry Smith PetscErrorCode ierr; 28998c9bda7SSatish Balay PetscInt mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL,nonzerorow=0; 29026e093fcSHong Zhang PetscTruth usecprow=a->compressedrow.use; 2912d61bbb3SSatish Balay 2922d61bbb3SSatish Balay PetscFunctionBegin; 293d9fead3dSBarry Smith ierr = VecGetArray(xx,(PetscScalar**)&x);CHKERRQ(ierr); 29426e093fcSHong Zhang ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr); 2952d61bbb3SSatish Balay 2962d61bbb3SSatish Balay idx = a->j; 2972d61bbb3SSatish Balay v = a->a; 29826e093fcSHong Zhang if (usecprow){ 29926e093fcSHong Zhang mbs = a->compressedrow.nrows; 30026e093fcSHong Zhang ii = a->compressedrow.i; 3017b2bb3b9SHong Zhang ridx = a->compressedrow.rindex; 30226e093fcSHong Zhang } else { 30326e093fcSHong Zhang mbs = a->mbs; 3042d61bbb3SSatish Balay ii = a->i; 30526e093fcSHong Zhang z = zarray; 30626e093fcSHong Zhang } 3072d61bbb3SSatish Balay 3082d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 3092d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 3102d61bbb3SSatish Balay sum1 = 0.0; sum2 = 0.0; 31198c9bda7SSatish Balay nonzerorow += (n>0); 3122d61bbb3SSatish Balay for (j=0; j<n; j++) { 3132d61bbb3SSatish Balay xb = x + 2*(*idx++); x1 = xb[0]; x2 = xb[1]; 3142d61bbb3SSatish Balay sum1 += v[0]*x1 + v[2]*x2; 3152d61bbb3SSatish Balay sum2 += v[1]*x1 + v[3]*x2; 3162d61bbb3SSatish Balay v += 4; 3172d61bbb3SSatish Balay } 3187b2bb3b9SHong Zhang if (usecprow) z = zarray + 2*ridx[i]; 3192d61bbb3SSatish Balay z[0] = sum1; z[1] = sum2; 32026e093fcSHong Zhang if (!usecprow) z += 2; 3212d61bbb3SSatish Balay } 322d9fead3dSBarry Smith ierr = VecRestoreArray(xx,(PetscScalar**)&x);CHKERRQ(ierr); 32326e093fcSHong Zhang ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr); 324dc0b31edSSatish Balay ierr = PetscLogFlops(8.0*a->nz - 2.0*nonzerorow);CHKERRQ(ierr); 3252d61bbb3SSatish Balay PetscFunctionReturn(0); 3262d61bbb3SSatish Balay } 3272d61bbb3SSatish Balay 3284a2ae208SSatish Balay #undef __FUNCT__ 3294a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_3" 330dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_3(Mat A,Vec xx,Vec zz) 3312d61bbb3SSatish Balay { 3322d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 333d9fead3dSBarry Smith PetscScalar *z = 0,sum1,sum2,sum3,x1,x2,x3,*zarray; 334d9fead3dSBarry Smith const PetscScalar *x,*xb; 335d9fead3dSBarry Smith const MatScalar *v; 336dfbe8321SBarry Smith PetscErrorCode ierr; 337028cd4eaSSatish Balay PetscInt mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL,nonzerorow=0; 338028cd4eaSSatish Balay PetscTruth usecprow=a->compressedrow.use; 33926e093fcSHong Zhang 3402d61bbb3SSatish Balay 341b6410449SSatish Balay #if defined(PETSC_HAVE_PRAGMA_DISJOINT) 342fee21e36SBarry Smith #pragma disjoint(*v,*z,*xb) 343fee21e36SBarry Smith #endif 344fee21e36SBarry Smith 3452d61bbb3SSatish Balay PetscFunctionBegin; 346d9fead3dSBarry Smith ierr = VecGetArray(xx,(PetscScalar**)&x);CHKERRQ(ierr); 34726e093fcSHong Zhang ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr); 3482d61bbb3SSatish Balay 3492d61bbb3SSatish Balay idx = a->j; 3502d61bbb3SSatish Balay v = a->a; 35126e093fcSHong Zhang if (usecprow){ 35226e093fcSHong Zhang mbs = a->compressedrow.nrows; 35326e093fcSHong Zhang ii = a->compressedrow.i; 3547b2bb3b9SHong Zhang ridx = a->compressedrow.rindex; 35526e093fcSHong Zhang } else { 35626e093fcSHong Zhang mbs = a->mbs; 3572d61bbb3SSatish Balay ii = a->i; 35826e093fcSHong Zhang z = zarray; 35926e093fcSHong Zhang } 3602d61bbb3SSatish Balay 3612d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 3622d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 3632d61bbb3SSatish Balay sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; 36498c9bda7SSatish Balay nonzerorow += (n>0); 3652d61bbb3SSatish Balay for (j=0; j<n; j++) { 3662d61bbb3SSatish Balay xb = x + 3*(*idx++); x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; 3672d61bbb3SSatish Balay sum1 += v[0]*x1 + v[3]*x2 + v[6]*x3; 3682d61bbb3SSatish Balay sum2 += v[1]*x1 + v[4]*x2 + v[7]*x3; 3692d61bbb3SSatish Balay sum3 += v[2]*x1 + v[5]*x2 + v[8]*x3; 3702d61bbb3SSatish Balay v += 9; 3712d61bbb3SSatish Balay } 3727b2bb3b9SHong Zhang if (usecprow) z = zarray + 3*ridx[i]; 3732d61bbb3SSatish Balay z[0] = sum1; z[1] = sum2; z[2] = sum3; 37426e093fcSHong Zhang if (!usecprow) z += 3; 3752d61bbb3SSatish Balay } 376d9fead3dSBarry Smith ierr = VecRestoreArray(xx,(PetscScalar**)&x);CHKERRQ(ierr); 37726e093fcSHong Zhang ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr); 378dc0b31edSSatish Balay ierr = PetscLogFlops(18.0*a->nz - 3.0*nonzerorow);CHKERRQ(ierr); 3792d61bbb3SSatish Balay PetscFunctionReturn(0); 3802d61bbb3SSatish Balay } 3812d61bbb3SSatish Balay 3824a2ae208SSatish Balay #undef __FUNCT__ 3834a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_4" 384dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_4(Mat A,Vec xx,Vec zz) 3852d61bbb3SSatish Balay { 3862d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 387d9fead3dSBarry Smith PetscScalar *z = 0,sum1,sum2,sum3,sum4,x1,x2,x3,x4,*zarray; 388d9fead3dSBarry Smith const PetscScalar *x,*xb; 389d9fead3dSBarry Smith const MatScalar *v; 390dfbe8321SBarry Smith PetscErrorCode ierr; 39198c9bda7SSatish Balay PetscInt mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL,nonzerorow=0; 39226e093fcSHong Zhang PetscTruth usecprow=a->compressedrow.use; 3932d61bbb3SSatish Balay 3942d61bbb3SSatish Balay PetscFunctionBegin; 395d9fead3dSBarry Smith ierr = VecGetArray(xx,(PetscScalar**)&x);CHKERRQ(ierr); 39626e093fcSHong Zhang ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr); 3972d61bbb3SSatish Balay 3982d61bbb3SSatish Balay idx = a->j; 3992d61bbb3SSatish Balay v = a->a; 40026e093fcSHong Zhang if (usecprow){ 40126e093fcSHong Zhang mbs = a->compressedrow.nrows; 40226e093fcSHong Zhang ii = a->compressedrow.i; 4037b2bb3b9SHong Zhang ridx = a->compressedrow.rindex; 40426e093fcSHong Zhang } else { 40526e093fcSHong Zhang mbs = a->mbs; 4062d61bbb3SSatish Balay ii = a->i; 40726e093fcSHong Zhang z = zarray; 40826e093fcSHong Zhang } 4092d61bbb3SSatish Balay 4102d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 4112d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 4122d61bbb3SSatish Balay sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; 41398c9bda7SSatish Balay nonzerorow += (n>0); 4142d61bbb3SSatish Balay for (j=0; j<n; j++) { 4152d61bbb3SSatish Balay xb = x + 4*(*idx++); 4162d61bbb3SSatish Balay x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; 4172d61bbb3SSatish Balay sum1 += v[0]*x1 + v[4]*x2 + v[8]*x3 + v[12]*x4; 4182d61bbb3SSatish Balay sum2 += v[1]*x1 + v[5]*x2 + v[9]*x3 + v[13]*x4; 4192d61bbb3SSatish Balay sum3 += v[2]*x1 + v[6]*x2 + v[10]*x3 + v[14]*x4; 4202d61bbb3SSatish Balay sum4 += v[3]*x1 + v[7]*x2 + v[11]*x3 + v[15]*x4; 4212d61bbb3SSatish Balay v += 16; 4222d61bbb3SSatish Balay } 4237b2bb3b9SHong Zhang if (usecprow) z = zarray + 4*ridx[i]; 4242d61bbb3SSatish Balay z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; 42526e093fcSHong Zhang if (!usecprow) z += 4; 4262d61bbb3SSatish Balay } 427d9fead3dSBarry Smith ierr = VecRestoreArray(xx,(PetscScalar**)&x);CHKERRQ(ierr); 42826e093fcSHong Zhang ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr); 429dc0b31edSSatish Balay ierr = PetscLogFlops(32.0*a->nz - 4.0*nonzerorow);CHKERRQ(ierr); 4302d61bbb3SSatish Balay PetscFunctionReturn(0); 4312d61bbb3SSatish Balay } 4322d61bbb3SSatish Balay 4334a2ae208SSatish Balay #undef __FUNCT__ 4344a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_5" 435dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_5(Mat A,Vec xx,Vec zz) 4362d61bbb3SSatish Balay { 4372d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 438d9fead3dSBarry Smith PetscScalar sum1,sum2,sum3,sum4,sum5,x1,x2,x3,x4,x5,*z = 0,*zarray; 439d9fead3dSBarry Smith const PetscScalar *xb,*x; 440d9fead3dSBarry Smith const MatScalar *v; 441dfbe8321SBarry Smith PetscErrorCode ierr; 442766f9fbaSBarry Smith const PetscInt *idx,*ii,*ridx=PETSC_NULL; 443766f9fbaSBarry Smith PetscInt mbs,i,j,n,nonzerorow=0; 44426e093fcSHong Zhang PetscTruth usecprow=a->compressedrow.use; 4452d61bbb3SSatish Balay 446433994e6SBarry Smith PetscFunctionBegin; 447d9fead3dSBarry Smith ierr = VecGetArray(xx,(PetscScalar**)&x);CHKERRQ(ierr); 44826e093fcSHong Zhang ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr); 4492d61bbb3SSatish Balay 4502d61bbb3SSatish Balay idx = a->j; 4512d61bbb3SSatish Balay v = a->a; 45226e093fcSHong Zhang if (usecprow){ 45326e093fcSHong Zhang mbs = a->compressedrow.nrows; 45426e093fcSHong Zhang ii = a->compressedrow.i; 4557b2bb3b9SHong Zhang ridx = a->compressedrow.rindex; 45626e093fcSHong Zhang } else { 45726e093fcSHong Zhang mbs = a->mbs; 4582d61bbb3SSatish Balay ii = a->i; 45926e093fcSHong Zhang z = zarray; 46026e093fcSHong Zhang } 4612d61bbb3SSatish Balay 4622d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 4632d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 4642d61bbb3SSatish Balay sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; 46598c9bda7SSatish Balay nonzerorow += (n>0); 4662d61bbb3SSatish Balay for (j=0; j<n; j++) { 4672d61bbb3SSatish Balay xb = x + 5*(*idx++); 4682d61bbb3SSatish Balay x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; 4692d61bbb3SSatish Balay sum1 += v[0]*x1 + v[5]*x2 + v[10]*x3 + v[15]*x4 + v[20]*x5; 4702d61bbb3SSatish Balay sum2 += v[1]*x1 + v[6]*x2 + v[11]*x3 + v[16]*x4 + v[21]*x5; 4712d61bbb3SSatish Balay sum3 += v[2]*x1 + v[7]*x2 + v[12]*x3 + v[17]*x4 + v[22]*x5; 4722d61bbb3SSatish Balay sum4 += v[3]*x1 + v[8]*x2 + v[13]*x3 + v[18]*x4 + v[23]*x5; 4732d61bbb3SSatish Balay sum5 += v[4]*x1 + v[9]*x2 + v[14]*x3 + v[19]*x4 + v[24]*x5; 4742d61bbb3SSatish Balay v += 25; 4752d61bbb3SSatish Balay } 4767b2bb3b9SHong Zhang if (usecprow) z = zarray + 5*ridx[i]; 4772d61bbb3SSatish Balay z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; 47826e093fcSHong Zhang if (!usecprow) z += 5; 4792d61bbb3SSatish Balay } 480d9fead3dSBarry Smith ierr = VecRestoreArray(xx,(PetscScalar**)&x);CHKERRQ(ierr); 48126e093fcSHong Zhang ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr); 482dc0b31edSSatish Balay ierr = PetscLogFlops(50.0*a->nz - 5.0*nonzerorow);CHKERRQ(ierr); 4832d61bbb3SSatish Balay PetscFunctionReturn(0); 4842d61bbb3SSatish Balay } 4852d61bbb3SSatish Balay 48615091d37SBarry Smith 4874a2ae208SSatish Balay #undef __FUNCT__ 4884a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_6" 489dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_6(Mat A,Vec xx,Vec zz) 49015091d37SBarry Smith { 49115091d37SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 492d9fead3dSBarry Smith PetscScalar *z = 0,sum1,sum2,sum3,sum4,sum5,sum6; 493d9fead3dSBarry Smith const PetscScalar *x,*xb; 49426e093fcSHong Zhang PetscScalar x1,x2,x3,x4,x5,x6,*zarray; 495d9fead3dSBarry Smith const MatScalar *v; 496dfbe8321SBarry Smith PetscErrorCode ierr; 49798c9bda7SSatish Balay PetscInt mbs=a->mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL,nonzerorow=0; 49826e093fcSHong Zhang PetscTruth usecprow=a->compressedrow.use; 49915091d37SBarry Smith 500433994e6SBarry Smith PetscFunctionBegin; 501d9fead3dSBarry Smith ierr = VecGetArray(xx,(PetscScalar**)&x);CHKERRQ(ierr); 50226e093fcSHong Zhang ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr); 50315091d37SBarry Smith 50415091d37SBarry Smith idx = a->j; 50515091d37SBarry Smith v = a->a; 50626e093fcSHong Zhang if (usecprow){ 50726e093fcSHong Zhang mbs = a->compressedrow.nrows; 50826e093fcSHong Zhang ii = a->compressedrow.i; 5097b2bb3b9SHong Zhang ridx = a->compressedrow.rindex; 51026e093fcSHong Zhang } else { 51126e093fcSHong Zhang mbs = a->mbs; 51215091d37SBarry Smith ii = a->i; 51326e093fcSHong Zhang z = zarray; 51426e093fcSHong Zhang } 51515091d37SBarry Smith 51615091d37SBarry Smith for (i=0; i<mbs; i++) { 51715091d37SBarry Smith n = ii[1] - ii[0]; ii++; 51815091d37SBarry Smith sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0; 51998c9bda7SSatish Balay nonzerorow += (n>0); 52015091d37SBarry Smith for (j=0; j<n; j++) { 52115091d37SBarry Smith xb = x + 6*(*idx++); 52215091d37SBarry Smith x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; 52315091d37SBarry Smith sum1 += v[0]*x1 + v[6]*x2 + v[12]*x3 + v[18]*x4 + v[24]*x5 + v[30]*x6; 52415091d37SBarry Smith sum2 += v[1]*x1 + v[7]*x2 + v[13]*x3 + v[19]*x4 + v[25]*x5 + v[31]*x6; 52515091d37SBarry Smith sum3 += v[2]*x1 + v[8]*x2 + v[14]*x3 + v[20]*x4 + v[26]*x5 + v[32]*x6; 52615091d37SBarry Smith sum4 += v[3]*x1 + v[9]*x2 + v[15]*x3 + v[21]*x4 + v[27]*x5 + v[33]*x6; 52715091d37SBarry Smith sum5 += v[4]*x1 + v[10]*x2 + v[16]*x3 + v[22]*x4 + v[28]*x5 + v[34]*x6; 52815091d37SBarry Smith sum6 += v[5]*x1 + v[11]*x2 + v[17]*x3 + v[23]*x4 + v[29]*x5 + v[35]*x6; 52915091d37SBarry Smith v += 36; 53015091d37SBarry Smith } 5317b2bb3b9SHong Zhang if (usecprow) z = zarray + 6*ridx[i]; 53215091d37SBarry Smith z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; 53326e093fcSHong Zhang if (!usecprow) z += 6; 53415091d37SBarry Smith } 53515091d37SBarry Smith 536d9fead3dSBarry Smith ierr = VecRestoreArray(xx,(PetscScalar**)&x);CHKERRQ(ierr); 53726e093fcSHong Zhang ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr); 538dc0b31edSSatish Balay ierr = PetscLogFlops(72.0*a->nz - 6.0*nonzerorow);CHKERRQ(ierr); 53915091d37SBarry Smith PetscFunctionReturn(0); 54015091d37SBarry Smith } 5418ab949d8SShri Abhyankar 5424a2ae208SSatish Balay #undef __FUNCT__ 5434a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_7" 544dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_7(Mat A,Vec xx,Vec zz) 5452d61bbb3SSatish Balay { 5462d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 547d9fead3dSBarry Smith PetscScalar *z = 0,sum1,sum2,sum3,sum4,sum5,sum6,sum7; 548d9fead3dSBarry Smith const PetscScalar *x,*xb; 54926e093fcSHong Zhang PetscScalar x1,x2,x3,x4,x5,x6,x7,*zarray; 550d9fead3dSBarry Smith const MatScalar *v; 551dfbe8321SBarry Smith PetscErrorCode ierr; 55298c9bda7SSatish Balay PetscInt mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL,nonzerorow=0; 55326e093fcSHong Zhang PetscTruth usecprow=a->compressedrow.use; 5542d61bbb3SSatish Balay 555433994e6SBarry Smith PetscFunctionBegin; 556d9fead3dSBarry Smith ierr = VecGetArray(xx,(PetscScalar**)&x);CHKERRQ(ierr); 55726e093fcSHong Zhang ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr); 5582d61bbb3SSatish Balay 5592d61bbb3SSatish Balay idx = a->j; 5602d61bbb3SSatish Balay v = a->a; 56126e093fcSHong Zhang if (usecprow){ 56226e093fcSHong Zhang mbs = a->compressedrow.nrows; 56326e093fcSHong Zhang ii = a->compressedrow.i; 5647b2bb3b9SHong Zhang ridx = a->compressedrow.rindex; 56526e093fcSHong Zhang } else { 56626e093fcSHong Zhang mbs = a->mbs; 5672d61bbb3SSatish Balay ii = a->i; 56826e093fcSHong Zhang z = zarray; 56926e093fcSHong Zhang } 5702d61bbb3SSatish Balay 5712d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 5722d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 5732d61bbb3SSatish Balay sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0; sum7 = 0.0; 57498c9bda7SSatish Balay nonzerorow += (n>0); 5752d61bbb3SSatish Balay for (j=0; j<n; j++) { 5762d61bbb3SSatish Balay xb = x + 7*(*idx++); 5772d61bbb3SSatish Balay x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6]; 5782d61bbb3SSatish Balay sum1 += v[0]*x1 + v[7]*x2 + v[14]*x3 + v[21]*x4 + v[28]*x5 + v[35]*x6 + v[42]*x7; 5792d61bbb3SSatish Balay sum2 += v[1]*x1 + v[8]*x2 + v[15]*x3 + v[22]*x4 + v[29]*x5 + v[36]*x6 + v[43]*x7; 5802d61bbb3SSatish Balay sum3 += v[2]*x1 + v[9]*x2 + v[16]*x3 + v[23]*x4 + v[30]*x5 + v[37]*x6 + v[44]*x7; 5812d61bbb3SSatish Balay sum4 += v[3]*x1 + v[10]*x2 + v[17]*x3 + v[24]*x4 + v[31]*x5 + v[38]*x6 + v[45]*x7; 5822d61bbb3SSatish Balay sum5 += v[4]*x1 + v[11]*x2 + v[18]*x3 + v[25]*x4 + v[32]*x5 + v[39]*x6 + v[46]*x7; 5832d61bbb3SSatish Balay sum6 += v[5]*x1 + v[12]*x2 + v[19]*x3 + v[26]*x4 + v[33]*x5 + v[40]*x6 + v[47]*x7; 5842d61bbb3SSatish Balay sum7 += v[6]*x1 + v[13]*x2 + v[20]*x3 + v[27]*x4 + v[34]*x5 + v[41]*x6 + v[48]*x7; 5852d61bbb3SSatish Balay v += 49; 5862d61bbb3SSatish Balay } 5877b2bb3b9SHong Zhang if (usecprow) z = zarray + 7*ridx[i]; 5882d61bbb3SSatish Balay z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7; 58926e093fcSHong Zhang if (!usecprow) z += 7; 5902d61bbb3SSatish Balay } 5912d61bbb3SSatish Balay 592d9fead3dSBarry Smith ierr = VecRestoreArray(xx,(PetscScalar**)&x);CHKERRQ(ierr); 59326e093fcSHong Zhang ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr); 594dc0b31edSSatish Balay ierr = PetscLogFlops(98.0*a->nz - 7.0*nonzerorow);CHKERRQ(ierr); 5952d61bbb3SSatish Balay PetscFunctionReturn(0); 5962d61bbb3SSatish Balay } 5972d61bbb3SSatish Balay 5988ab949d8SShri Abhyankar /* MatMult_SeqBAIJ_15 version 1: Columns in the block are accessed one at a time */ 5998ab949d8SShri Abhyankar 6008ab949d8SShri Abhyankar #undef __FUNCT__ 6018ab949d8SShri Abhyankar #define __FUNCT__ "MatMult_SeqBAIJ_15_ver1" 6028ab949d8SShri Abhyankar PetscErrorCode MatMult_SeqBAIJ_15_ver1(Mat A,Vec xx,Vec zz) 6038ab949d8SShri Abhyankar { 6048ab949d8SShri Abhyankar Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 6058ab949d8SShri Abhyankar PetscScalar *z = 0,sum1,sum2,sum3,sum4,sum5,sum6,sum7,sum8,sum9,sum10,sum11,sum12,sum13,sum14,sum15; 6068ab949d8SShri Abhyankar const PetscScalar *x,*xb; 6078ab949d8SShri Abhyankar PetscScalar x1,x2,x3,x4,x5,x6,x7,x8,x9,x10,x11,x12,x13,x14,x15,*zarray; 6088ab949d8SShri Abhyankar const MatScalar *v; 6098ab949d8SShri Abhyankar PetscErrorCode ierr; 6108ab949d8SShri Abhyankar const PetscInt *ii,*ij=a->j,*idx; 6118ab949d8SShri Abhyankar PetscInt mbs,i,j,k,n,*ridx=PETSC_NULL,nonzerorow=0; 6128ab949d8SShri Abhyankar PetscTruth usecprow=a->compressedrow.use; 6138ab949d8SShri Abhyankar 6148ab949d8SShri Abhyankar PetscFunctionBegin; 6158ab949d8SShri Abhyankar ierr = VecGetArray(xx,(PetscScalar**)&x);CHKERRQ(ierr); 6168ab949d8SShri Abhyankar ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr); 6178ab949d8SShri Abhyankar 6188ab949d8SShri Abhyankar v = a->a; 6198ab949d8SShri Abhyankar if (usecprow){ 6208ab949d8SShri Abhyankar mbs = a->compressedrow.nrows; 6218ab949d8SShri Abhyankar ii = a->compressedrow.i; 6228ab949d8SShri Abhyankar ridx = a->compressedrow.rindex; 6238ab949d8SShri Abhyankar } else { 6248ab949d8SShri Abhyankar mbs = a->mbs; 6258ab949d8SShri Abhyankar ii = a->i; 6268ab949d8SShri Abhyankar z = zarray; 6278ab949d8SShri Abhyankar } 6288ab949d8SShri Abhyankar 6298ab949d8SShri Abhyankar for (i=0; i<mbs; i++) { 6308ab949d8SShri Abhyankar n = ii[i+1] - ii[i]; 6318ab949d8SShri Abhyankar idx = ij + ii[i]; 6328ab949d8SShri Abhyankar sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0; sum7 = 0.0; 6338ab949d8SShri 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; 6348ab949d8SShri Abhyankar 6358ab949d8SShri Abhyankar nonzerorow += (n>0); 6368ab949d8SShri Abhyankar for (j=0; j<n; j++) { 6378ab949d8SShri Abhyankar xb = x + 15*(idx[j]); 6388ab949d8SShri Abhyankar x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6]; 6398ab949d8SShri 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]; 6408ab949d8SShri Abhyankar 6418ab949d8SShri Abhyankar for(k=0;k<15;k++){ 6428ab949d8SShri Abhyankar sum1 += v[0]*xb[k]; 6438ab949d8SShri Abhyankar sum2 += v[1]*xb[k]; 6448ab949d8SShri Abhyankar sum3 += v[2]*xb[k]; 6458ab949d8SShri Abhyankar sum4 += v[3]*xb[k]; 6468ab949d8SShri Abhyankar sum5 += v[4]*xb[k]; 6478ab949d8SShri Abhyankar sum6 += v[5]*xb[k]; 6488ab949d8SShri Abhyankar sum7 += v[6]*xb[k]; 6498ab949d8SShri Abhyankar sum8 += v[7]*xb[k]; 6508ab949d8SShri Abhyankar sum9 += v[8]*xb[k]; 6518ab949d8SShri Abhyankar sum10 += v[9]*xb[k]; 6528ab949d8SShri Abhyankar sum11 += v[10]*xb[k]; 6538ab949d8SShri Abhyankar sum12 += v[11]*xb[k]; 6548ab949d8SShri Abhyankar sum13 += v[12]*xb[k]; 6558ab949d8SShri Abhyankar sum14 += v[13]*xb[k]; 6568ab949d8SShri Abhyankar sum15 += v[14]*xb[k]; 6578ab949d8SShri Abhyankar v += 15; 6588ab949d8SShri Abhyankar } 6598ab949d8SShri Abhyankar } 6608ab949d8SShri Abhyankar if (usecprow) z = zarray + 15*ridx[i]; 6618ab949d8SShri Abhyankar z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7; 6628ab949d8SShri 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; 6638ab949d8SShri Abhyankar 6648ab949d8SShri Abhyankar if (!usecprow) z += 15; 6658ab949d8SShri Abhyankar } 6668ab949d8SShri Abhyankar 6678ab949d8SShri Abhyankar ierr = VecRestoreArray(xx,(PetscScalar**)&x);CHKERRQ(ierr); 6688ab949d8SShri Abhyankar ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr); 6698ab949d8SShri Abhyankar ierr = PetscLogFlops(450.0*a->nz - 15.0*nonzerorow);CHKERRQ(ierr); 6708ab949d8SShri Abhyankar PetscFunctionReturn(0); 6718ab949d8SShri Abhyankar } 6728ab949d8SShri Abhyankar 6738ab949d8SShri Abhyankar /* MatMult_SeqBAIJ_15_ver2 : Columns in the block are accessed in sets of 4,4,4,3 */ 6748ab949d8SShri Abhyankar #undef __FUNCT__ 6758ab949d8SShri Abhyankar #define __FUNCT__ "MatMult_SeqBAIJ_15_ver2" 6768ab949d8SShri Abhyankar PetscErrorCode MatMult_SeqBAIJ_15_ver2(Mat A,Vec xx,Vec zz) 6778ab949d8SShri Abhyankar { 6788ab949d8SShri Abhyankar Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 6798ab949d8SShri Abhyankar PetscScalar *z = 0,sum1,sum2,sum3,sum4,sum5,sum6,sum7,sum8,sum9,sum10,sum11,sum12,sum13,sum14,sum15; 6808ab949d8SShri Abhyankar const PetscScalar *x,*xb; 681*0b8f6341SShri Abhyankar PetscScalar x1,x2,x3,x4,*zarray; 6828ab949d8SShri Abhyankar const MatScalar *v; 6838ab949d8SShri Abhyankar PetscErrorCode ierr; 6848ab949d8SShri Abhyankar const PetscInt *ii,*ij=a->j,*idx; 6858ab949d8SShri Abhyankar PetscInt mbs,i,j,n,*ridx=PETSC_NULL,nonzerorow=0; 6868ab949d8SShri Abhyankar PetscTruth usecprow=a->compressedrow.use; 6878ab949d8SShri Abhyankar 6888ab949d8SShri Abhyankar PetscFunctionBegin; 6898ab949d8SShri Abhyankar ierr = VecGetArray(xx,(PetscScalar**)&x);CHKERRQ(ierr); 6908ab949d8SShri Abhyankar ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr); 6918ab949d8SShri Abhyankar 6928ab949d8SShri Abhyankar v = a->a; 6938ab949d8SShri Abhyankar if (usecprow){ 6948ab949d8SShri Abhyankar mbs = a->compressedrow.nrows; 6958ab949d8SShri Abhyankar ii = a->compressedrow.i; 6968ab949d8SShri Abhyankar ridx = a->compressedrow.rindex; 6978ab949d8SShri Abhyankar } else { 6988ab949d8SShri Abhyankar mbs = a->mbs; 6998ab949d8SShri Abhyankar ii = a->i; 7008ab949d8SShri Abhyankar z = zarray; 7018ab949d8SShri Abhyankar } 7028ab949d8SShri Abhyankar 7038ab949d8SShri Abhyankar for (i=0; i<mbs; i++) { 7048ab949d8SShri Abhyankar n = ii[i+1] - ii[i]; 7058ab949d8SShri Abhyankar idx = ij + ii[i]; 7068ab949d8SShri Abhyankar sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0; sum7 = 0.0; 7078ab949d8SShri 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; 7088ab949d8SShri Abhyankar 7098ab949d8SShri Abhyankar nonzerorow += (n>0); 7108ab949d8SShri Abhyankar for (j=0; j<n; j++) { 7118ab949d8SShri Abhyankar xb = x + 15*(idx[j]); 712*0b8f6341SShri Abhyankar x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; 7138ab949d8SShri Abhyankar 7148ab949d8SShri Abhyankar sum1 += v[0]*x1 + v[15]*x2 + v[30]*x3 + v[45]*x4; 7158ab949d8SShri Abhyankar sum2 += v[1]*x1 + v[16]*x2 + v[31]*x3 + v[46]*x4; 7168ab949d8SShri Abhyankar sum3 += v[2]*x1 + v[17]*x2 + v[32]*x3 + v[47]*x4; 7178ab949d8SShri Abhyankar sum4 += v[3]*x1 + v[18]*x2 + v[33]*x3 + v[48]*x4; 7188ab949d8SShri Abhyankar sum5 += v[4]*x1 + v[19]*x2 + v[34]*x3 + v[49]*x4; 7198ab949d8SShri Abhyankar sum6 += v[5]*x1 + v[20]*x2 + v[35]*x3 + v[50]*x4; 7208ab949d8SShri Abhyankar sum7 += v[6]*x1 + v[21]*x2 + v[36]*x3 + v[51]*x4; 7218ab949d8SShri Abhyankar sum8 += v[7]*x1 + v[22]*x2 + v[37]*x3 + v[52]*x4; 7228ab949d8SShri Abhyankar sum9 += v[8]*x1 + v[23]*x2 + v[38]*x3 + v[53]*x4; 7238ab949d8SShri Abhyankar sum10 += v[9]*x1 + v[24]*x2 + v[39]*x3 + v[54]*x4; 7248ab949d8SShri Abhyankar sum11 += v[10]*x1 + v[25]*x2 + v[40]*x3 + v[55]*x4; 7258ab949d8SShri Abhyankar sum12 += v[11]*x1 + v[26]*x2 + v[41]*x3 + v[56]*x4; 7268ab949d8SShri Abhyankar sum13 += v[12]*x1 + v[27]*x2 + v[42]*x3 + v[57]*x4; 7278ab949d8SShri Abhyankar sum14 += v[13]*x1 + v[28]*x2 + v[43]*x3 + v[58]*x4; 7288ab949d8SShri Abhyankar sum15 += v[14]*x1 + v[29]*x2 + v[44]*x3 + v[59]*x4; 7298ab949d8SShri Abhyankar 7308ab949d8SShri Abhyankar v += 60; 7318ab949d8SShri Abhyankar 732*0b8f6341SShri Abhyankar x1 = xb[4]; x2 = xb[5]; x3 = xb[6]; x4 = xb[7]; 7338ab949d8SShri Abhyankar 7348ab949d8SShri Abhyankar sum1 += v[0]*x1 + v[15]*x2 + v[30]*x3 + v[45]*x4; 7358ab949d8SShri Abhyankar sum2 += v[1]*x1 + v[16]*x2 + v[31]*x3 + v[46]*x4; 7368ab949d8SShri Abhyankar sum3 += v[2]*x1 + v[17]*x2 + v[32]*x3 + v[47]*x4; 7378ab949d8SShri Abhyankar sum4 += v[3]*x1 + v[18]*x2 + v[33]*x3 + v[48]*x4; 7388ab949d8SShri Abhyankar sum5 += v[4]*x1 + v[19]*x2 + v[34]*x3 + v[49]*x4; 7398ab949d8SShri Abhyankar sum6 += v[5]*x1 + v[20]*x2 + v[35]*x3 + v[50]*x4; 7408ab949d8SShri Abhyankar sum7 += v[6]*x1 + v[21]*x2 + v[36]*x3 + v[51]*x4; 7418ab949d8SShri Abhyankar sum8 += v[7]*x1 + v[22]*x2 + v[37]*x3 + v[52]*x4; 7428ab949d8SShri Abhyankar sum9 += v[8]*x1 + v[23]*x2 + v[38]*x3 + v[53]*x4; 7438ab949d8SShri Abhyankar sum10 += v[9]*x1 + v[24]*x2 + v[39]*x3 + v[54]*x4; 7448ab949d8SShri Abhyankar sum11 += v[10]*x1 + v[25]*x2 + v[40]*x3 + v[55]*x4; 7458ab949d8SShri Abhyankar sum12 += v[11]*x1 + v[26]*x2 + v[41]*x3 + v[56]*x4; 7468ab949d8SShri Abhyankar sum13 += v[12]*x1 + v[27]*x2 + v[42]*x3 + v[57]*x4; 7478ab949d8SShri Abhyankar sum14 += v[13]*x1 + v[28]*x2 + v[43]*x3 + v[58]*x4; 7488ab949d8SShri Abhyankar sum15 += v[14]*x1 + v[29]*x2 + v[44]*x3 + v[59]*x4; 7498ab949d8SShri Abhyankar v += 60; 7508ab949d8SShri Abhyankar 751*0b8f6341SShri Abhyankar x1 = xb[8]; x2 = xb[9]; x3 = xb[10]; x4 = xb[11]; 752*0b8f6341SShri Abhyankar sum1 += v[0]*x1 + v[15]*x2 + v[30]*x3 + v[45]*x4; 753*0b8f6341SShri Abhyankar sum2 += v[1]*x1 + v[16]*x2 + v[31]*x3 + v[46]*x4; 754*0b8f6341SShri Abhyankar sum3 += v[2]*x1 + v[17]*x2 + v[32]*x3 + v[47]*x4; 755*0b8f6341SShri Abhyankar sum4 += v[3]*x1 + v[18]*x2 + v[33]*x3 + v[48]*x4; 756*0b8f6341SShri Abhyankar sum5 += v[4]*x1 + v[19]*x2 + v[34]*x3 + v[49]*x4; 757*0b8f6341SShri Abhyankar sum6 += v[5]*x1 + v[20]*x2 + v[35]*x3 + v[50]*x4; 758*0b8f6341SShri Abhyankar sum7 += v[6]*x1 + v[21]*x2 + v[36]*x3 + v[51]*x4; 759*0b8f6341SShri Abhyankar sum8 += v[7]*x1 + v[22]*x2 + v[37]*x3 + v[52]*x4; 760*0b8f6341SShri Abhyankar sum9 += v[8]*x1 + v[23]*x2 + v[38]*x3 + v[53]*x4; 761*0b8f6341SShri Abhyankar sum10 += v[9]*x1 + v[24]*x2 + v[39]*x3 + v[54]*x4; 762*0b8f6341SShri Abhyankar sum11 += v[10]*x1 + v[25]*x2 + v[40]*x3 + v[55]*x4; 763*0b8f6341SShri Abhyankar sum12 += v[11]*x1 + v[26]*x2 + v[41]*x3 + v[56]*x4; 764*0b8f6341SShri Abhyankar sum13 += v[12]*x1 + v[27]*x2 + v[42]*x3 + v[57]*x4; 765*0b8f6341SShri Abhyankar sum14 += v[13]*x1 + v[28]*x2 + v[43]*x3 + v[58]*x4; 766*0b8f6341SShri Abhyankar sum15 += v[14]*x1 + v[29]*x2 + v[44]*x3 + v[59]*x4; 767*0b8f6341SShri Abhyankar v += 60; 768*0b8f6341SShri Abhyankar 769*0b8f6341SShri Abhyankar x1 = xb[12]; x2 = xb[13]; x3 = xb[14]; 7708ab949d8SShri Abhyankar sum1 += v[0]*x1 + v[15]*x2 + v[30]*x3; 7718ab949d8SShri Abhyankar sum2 += v[1]*x1 + v[16]*x2 + v[31]*x3; 7728ab949d8SShri Abhyankar sum3 += v[2]*x1 + v[17]*x2 + v[32]*x3; 7738ab949d8SShri Abhyankar sum4 += v[3]*x1 + v[18]*x2 + v[33]*x3; 7748ab949d8SShri Abhyankar sum5 += v[4]*x1 + v[19]*x2 + v[34]*x3; 7758ab949d8SShri Abhyankar sum6 += v[5]*x1 + v[20]*x2 + v[35]*x3; 7768ab949d8SShri Abhyankar sum7 += v[6]*x1 + v[21]*x2 + v[36]*x3; 7778ab949d8SShri Abhyankar sum8 += v[7]*x1 + v[22]*x2 + v[37]*x3; 7788ab949d8SShri Abhyankar sum9 += v[8]*x1 + v[23]*x2 + v[38]*x3; 7798ab949d8SShri Abhyankar sum10 += v[9]*x1 + v[24]*x2 + v[39]*x3; 7808ab949d8SShri Abhyankar sum11 += v[10]*x1 + v[25]*x2 + v[40]*x3; 7818ab949d8SShri Abhyankar sum12 += v[11]*x1 + v[26]*x2 + v[41]*x3; 7828ab949d8SShri Abhyankar sum13 += v[12]*x1 + v[27]*x2 + v[42]*x3; 7838ab949d8SShri Abhyankar sum14 += v[13]*x1 + v[28]*x2 + v[43]*x3; 7848ab949d8SShri Abhyankar sum15 += v[14]*x1 + v[29]*x2 + v[44]*x3; 7858ab949d8SShri Abhyankar v += 45; 7868ab949d8SShri Abhyankar } 7878ab949d8SShri Abhyankar if (usecprow) z = zarray + 15*ridx[i]; 7888ab949d8SShri Abhyankar z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7; 7898ab949d8SShri 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; 7908ab949d8SShri Abhyankar 7918ab949d8SShri Abhyankar if (!usecprow) z += 15; 7928ab949d8SShri Abhyankar } 7938ab949d8SShri Abhyankar 7948ab949d8SShri Abhyankar ierr = VecRestoreArray(xx,(PetscScalar**)&x);CHKERRQ(ierr); 7958ab949d8SShri Abhyankar ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr); 7968ab949d8SShri Abhyankar ierr = PetscLogFlops(450.0*a->nz - 15.0*nonzerorow);CHKERRQ(ierr); 7978ab949d8SShri Abhyankar PetscFunctionReturn(0); 7988ab949d8SShri Abhyankar } 7998ab949d8SShri Abhyankar 8008ab949d8SShri Abhyankar /* MatMult_SeqBAIJ_15_ver3 : Columns in the block are accessed in sets of 8,7 */ 8018ab949d8SShri Abhyankar #undef __FUNCT__ 8028ab949d8SShri Abhyankar #define __FUNCT__ "MatMult_SeqBAIJ_15_ver3" 8038ab949d8SShri Abhyankar PetscErrorCode MatMult_SeqBAIJ_15_ver3(Mat A,Vec xx,Vec zz) 8048ab949d8SShri Abhyankar { 8058ab949d8SShri Abhyankar Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 8068ab949d8SShri Abhyankar PetscScalar *z = 0,sum1,sum2,sum3,sum4,sum5,sum6,sum7,sum8,sum9,sum10,sum11,sum12,sum13,sum14,sum15; 8078ab949d8SShri Abhyankar const PetscScalar *x,*xb; 808*0b8f6341SShri Abhyankar PetscScalar x1,x2,x3,x4,x5,x6,x7,x8,*zarray; 8098ab949d8SShri Abhyankar const MatScalar *v; 8108ab949d8SShri Abhyankar PetscErrorCode ierr; 8118ab949d8SShri Abhyankar const PetscInt *ii,*ij=a->j,*idx; 8128ab949d8SShri Abhyankar PetscInt mbs,i,j,n,*ridx=PETSC_NULL,nonzerorow=0; 8138ab949d8SShri Abhyankar PetscTruth usecprow=a->compressedrow.use; 8148ab949d8SShri Abhyankar 8158ab949d8SShri Abhyankar PetscFunctionBegin; 8168ab949d8SShri Abhyankar ierr = VecGetArray(xx,(PetscScalar**)&x);CHKERRQ(ierr); 8178ab949d8SShri Abhyankar ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr); 8188ab949d8SShri Abhyankar 8198ab949d8SShri Abhyankar v = a->a; 8208ab949d8SShri Abhyankar if (usecprow){ 8218ab949d8SShri Abhyankar mbs = a->compressedrow.nrows; 8228ab949d8SShri Abhyankar ii = a->compressedrow.i; 8238ab949d8SShri Abhyankar ridx = a->compressedrow.rindex; 8248ab949d8SShri Abhyankar } else { 8258ab949d8SShri Abhyankar mbs = a->mbs; 8268ab949d8SShri Abhyankar ii = a->i; 8278ab949d8SShri Abhyankar z = zarray; 8288ab949d8SShri Abhyankar } 8298ab949d8SShri Abhyankar 8308ab949d8SShri Abhyankar for (i=0; i<mbs; i++) { 8318ab949d8SShri Abhyankar n = ii[i+1] - ii[i]; 8328ab949d8SShri Abhyankar idx = ij + ii[i]; 8338ab949d8SShri Abhyankar sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0; sum7 = 0.0; 8348ab949d8SShri 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; 8358ab949d8SShri Abhyankar 8368ab949d8SShri Abhyankar nonzerorow += (n>0); 8378ab949d8SShri Abhyankar for (j=0; j<n; j++) { 8388ab949d8SShri Abhyankar xb = x + 15*(idx[j]); 8398ab949d8SShri Abhyankar x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6]; 840*0b8f6341SShri Abhyankar x8 = xb[7]; 8418ab949d8SShri Abhyankar 8428ab949d8SShri 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; 8438ab949d8SShri 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; 8448ab949d8SShri 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; 8458ab949d8SShri 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; 8468ab949d8SShri 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; 8478ab949d8SShri 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; 8488ab949d8SShri 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; 8498ab949d8SShri 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; 8508ab949d8SShri 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; 8518ab949d8SShri 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; 8528ab949d8SShri 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; 8538ab949d8SShri 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; 8548ab949d8SShri 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; 8558ab949d8SShri 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; 8568ab949d8SShri 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; 8578ab949d8SShri Abhyankar v += 120; 8588ab949d8SShri Abhyankar 859*0b8f6341SShri Abhyankar x1 = xb[8]; x2 = xb[9]; x3 = xb[10]; x4 = xb[11]; x5 = xb[12]; x6 = xb[13]; x7 = xb[14]; 860*0b8f6341SShri Abhyankar 8618ab949d8SShri Abhyankar sum1 += v[0]*x1 + v[15]*x2 + v[30]*x3 + v[45]*x4 + v[60]*x5 + v[75]*x6 + v[90]*x7; 8628ab949d8SShri Abhyankar sum2 += v[1]*x1 + v[16]*x2 + v[31]*x3 + v[46]*x4 + v[61]*x5 + v[76]*x6 + v[91]*x7; 8638ab949d8SShri Abhyankar sum3 += v[2]*x1 + v[17]*x2 + v[32]*x3 + v[47]*x4 + v[62]*x5 + v[77]*x6 + v[92]*x7; 8648ab949d8SShri Abhyankar sum4 += v[3]*x1 + v[18]*x2 + v[33]*x3 + v[48]*x4 + v[63]*x5 + v[78]*x6 + v[93]*x7; 8658ab949d8SShri Abhyankar sum5 += v[4]*x1 + v[19]*x2 + v[34]*x3 + v[49]*x4 + v[64]*x5 + v[79]*x6 + v[94]*x7; 8668ab949d8SShri Abhyankar sum6 += v[5]*x1 + v[20]*x2 + v[35]*x3 + v[50]*x4 + v[65]*x5 + v[80]*x6 + v[95]*x7; 8678ab949d8SShri Abhyankar sum7 += v[6]*x1 + v[21]*x2 + v[36]*x3 + v[51]*x4 + v[66]*x5 + v[81]*x6 + v[96]*x7; 8688ab949d8SShri Abhyankar sum8 += v[7]*x1 + v[22]*x2 + v[37]*x3 + v[52]*x4 + v[67]*x5 + v[82]*x6 + v[97]*x7; 8698ab949d8SShri Abhyankar sum9 += v[8]*x1 + v[23]*x2 + v[38]*x3 + v[53]*x4 + v[68]*x5 + v[83]*x6 + v[98]*x7; 8708ab949d8SShri Abhyankar sum10 += v[9]*x1 + v[24]*x2 + v[39]*x3 + v[54]*x4 + v[69]*x5 + v[84]*x6 + v[99]*x7; 8718ab949d8SShri Abhyankar sum11 += v[10]*x1 + v[25]*x2 + v[40]*x3 + v[55]*x4 + v[70]*x5 + v[85]*x6 + v[100]*x7; 8728ab949d8SShri Abhyankar sum12 += v[11]*x1 + v[26]*x2 + v[41]*x3 + v[56]*x4 + v[71]*x5 + v[86]*x6 + v[101]*x7; 8738ab949d8SShri Abhyankar sum13 += v[12]*x1 + v[27]*x2 + v[42]*x3 + v[57]*x4 + v[72]*x5 + v[87]*x6 + v[102]*x7; 8748ab949d8SShri Abhyankar sum14 += v[13]*x1 + v[28]*x2 + v[43]*x3 + v[58]*x4 + v[73]*x5 + v[88]*x6 + v[103]*x7; 8758ab949d8SShri Abhyankar sum15 += v[14]*x1 + v[29]*x2 + v[44]*x3 + v[59]*x4 + v[74]*x5 + v[89]*x6 + v[104]*x7; 8768ab949d8SShri Abhyankar v += 105; 8778ab949d8SShri Abhyankar } 8788ab949d8SShri Abhyankar if (usecprow) z = zarray + 15*ridx[i]; 8798ab949d8SShri Abhyankar z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7; 8808ab949d8SShri 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; 8818ab949d8SShri Abhyankar 8828ab949d8SShri Abhyankar if (!usecprow) z += 15; 8838ab949d8SShri Abhyankar } 8848ab949d8SShri Abhyankar 8858ab949d8SShri Abhyankar ierr = VecRestoreArray(xx,(PetscScalar**)&x);CHKERRQ(ierr); 8868ab949d8SShri Abhyankar ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr); 8878ab949d8SShri Abhyankar ierr = PetscLogFlops(450.0*a->nz - 15.0*nonzerorow);CHKERRQ(ierr); 8888ab949d8SShri Abhyankar PetscFunctionReturn(0); 8898ab949d8SShri Abhyankar } 8908ab949d8SShri Abhyankar 8918ab949d8SShri Abhyankar /* MatMult_SeqBAIJ_15_ver4 : All columns in the block are accessed at once */ 8928ab949d8SShri Abhyankar 8938ab949d8SShri Abhyankar #undef __FUNCT__ 8948ab949d8SShri Abhyankar #define __FUNCT__ "MatMult_SeqBAIJ_15_ver4" 8958ab949d8SShri Abhyankar PetscErrorCode MatMult_SeqBAIJ_15_ver4(Mat A,Vec xx,Vec zz) 8968ab949d8SShri Abhyankar { 8978ab949d8SShri Abhyankar Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 8988ab949d8SShri Abhyankar PetscScalar *z = 0,sum1,sum2,sum3,sum4,sum5,sum6,sum7,sum8,sum9,sum10,sum11,sum12,sum13,sum14,sum15; 8998ab949d8SShri Abhyankar const PetscScalar *x,*xb; 9008ab949d8SShri Abhyankar PetscScalar x1,x2,x3,x4,x5,x6,x7,x8,x9,x10,x11,x12,x13,x14,x15,*zarray; 9018ab949d8SShri Abhyankar const MatScalar *v; 9028ab949d8SShri Abhyankar PetscErrorCode ierr; 9038ab949d8SShri Abhyankar const PetscInt *ii,*ij=a->j,*idx; 9048ab949d8SShri Abhyankar PetscInt mbs,i,j,n,*ridx=PETSC_NULL,nonzerorow=0; 9058ab949d8SShri Abhyankar PetscTruth usecprow=a->compressedrow.use; 9068ab949d8SShri Abhyankar 9078ab949d8SShri Abhyankar PetscFunctionBegin; 9088ab949d8SShri Abhyankar ierr = VecGetArray(xx,(PetscScalar**)&x);CHKERRQ(ierr); 9098ab949d8SShri Abhyankar ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr); 9108ab949d8SShri Abhyankar 9118ab949d8SShri Abhyankar v = a->a; 9128ab949d8SShri Abhyankar if (usecprow){ 9138ab949d8SShri Abhyankar mbs = a->compressedrow.nrows; 9148ab949d8SShri Abhyankar ii = a->compressedrow.i; 9158ab949d8SShri Abhyankar ridx = a->compressedrow.rindex; 9168ab949d8SShri Abhyankar } else { 9178ab949d8SShri Abhyankar mbs = a->mbs; 9188ab949d8SShri Abhyankar ii = a->i; 9198ab949d8SShri Abhyankar z = zarray; 9208ab949d8SShri Abhyankar } 9218ab949d8SShri Abhyankar 9228ab949d8SShri Abhyankar for (i=0; i<mbs; i++) { 9238ab949d8SShri Abhyankar n = ii[i+1] - ii[i]; 9248ab949d8SShri Abhyankar idx = ij + ii[i]; 9258ab949d8SShri Abhyankar sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0; sum7 = 0.0; 9268ab949d8SShri 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; 9278ab949d8SShri Abhyankar 9288ab949d8SShri Abhyankar nonzerorow += (n>0); 9298ab949d8SShri Abhyankar for (j=0; j<n; j++) { 9308ab949d8SShri Abhyankar xb = x + 15*(idx[j]); 9318ab949d8SShri Abhyankar x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6]; 9328ab949d8SShri 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]; 9338ab949d8SShri Abhyankar 9348ab949d8SShri 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; 9358ab949d8SShri 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; 9368ab949d8SShri 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; 9378ab949d8SShri 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; 9388ab949d8SShri 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; 9398ab949d8SShri 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; 9408ab949d8SShri 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; 9418ab949d8SShri 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; 9428ab949d8SShri 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; 9438ab949d8SShri 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; 9448ab949d8SShri 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; 9458ab949d8SShri 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; 9468ab949d8SShri 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; 9478ab949d8SShri 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; 9488ab949d8SShri 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; 9498ab949d8SShri Abhyankar v += 225; 9508ab949d8SShri Abhyankar } 9518ab949d8SShri Abhyankar if (usecprow) z = zarray + 15*ridx[i]; 9528ab949d8SShri Abhyankar z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7; 9538ab949d8SShri 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; 9548ab949d8SShri Abhyankar 9558ab949d8SShri Abhyankar if (!usecprow) z += 15; 9568ab949d8SShri Abhyankar } 9578ab949d8SShri Abhyankar 9588ab949d8SShri Abhyankar ierr = VecRestoreArray(xx,(PetscScalar**)&x);CHKERRQ(ierr); 9598ab949d8SShri Abhyankar ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr); 9608ab949d8SShri Abhyankar ierr = PetscLogFlops(450.0*a->nz - 15.0*nonzerorow);CHKERRQ(ierr); 9618ab949d8SShri Abhyankar PetscFunctionReturn(0); 9628ab949d8SShri Abhyankar } 9638ab949d8SShri Abhyankar 9648ab949d8SShri Abhyankar 9653f1db9ecSBarry Smith /* 9663f1db9ecSBarry Smith This will not work with MatScalar == float because it calls the BLAS 9673f1db9ecSBarry Smith */ 9684a2ae208SSatish Balay #undef __FUNCT__ 9694a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_N" 970dfbe8321SBarry Smith PetscErrorCode MatMult_SeqBAIJ_N(Mat A,Vec xx,Vec zz) 9712d61bbb3SSatish Balay { 9722d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 973dcf5cc72SBarry Smith PetscScalar *x,*z = 0,*xb,*work,*workt,*zarray; 9743f1db9ecSBarry Smith MatScalar *v; 975dfbe8321SBarry Smith PetscErrorCode ierr; 976d0f46423SBarry Smith PetscInt mbs=a->mbs,i,*idx,*ii,bs=A->rmap->bs,j,n,bs2=a->bs2; 97798c9bda7SSatish Balay PetscInt ncols,k,*ridx=PETSC_NULL,nonzerorow=0; 97826e093fcSHong Zhang PetscTruth usecprow=a->compressedrow.use; 9792d61bbb3SSatish Balay 9802d61bbb3SSatish Balay PetscFunctionBegin; 9811ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 98226e093fcSHong Zhang ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr); 9832d61bbb3SSatish Balay 9842d61bbb3SSatish Balay idx = a->j; 9852d61bbb3SSatish Balay v = a->a; 98626e093fcSHong Zhang if (usecprow){ 98726e093fcSHong Zhang mbs = a->compressedrow.nrows; 98826e093fcSHong Zhang ii = a->compressedrow.i; 9897b2bb3b9SHong Zhang ridx = a->compressedrow.rindex; 99026e093fcSHong Zhang } else { 99126e093fcSHong Zhang mbs = a->mbs; 9922d61bbb3SSatish Balay ii = a->i; 99326e093fcSHong Zhang z = zarray; 99426e093fcSHong Zhang } 995218c64b6SSatish Balay 9962d61bbb3SSatish Balay if (!a->mult_work) { 997d0f46423SBarry Smith k = PetscMax(A->rmap->n,A->cmap->n); 99887828ca2SBarry Smith ierr = PetscMalloc((k+1)*sizeof(PetscScalar),&a->mult_work);CHKERRQ(ierr); 9992d61bbb3SSatish Balay } 10002d61bbb3SSatish Balay work = a->mult_work; 10012d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 10022d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 10032d61bbb3SSatish Balay ncols = n*bs; 10042d61bbb3SSatish Balay workt = work; 100598c9bda7SSatish Balay nonzerorow += (n>0); 10062d61bbb3SSatish Balay for (j=0; j<n; j++) { 10072d61bbb3SSatish Balay xb = x + bs*(*idx++); 10082d61bbb3SSatish Balay for (k=0; k<bs; k++) workt[k] = xb[k]; 10092d61bbb3SSatish Balay workt += bs; 10102d61bbb3SSatish Balay } 10117b2bb3b9SHong Zhang if (usecprow) z = zarray + bs*ridx[i]; 1012737d121aSSatish Balay Kernel_w_gets_Ar_times_v(bs,ncols,work,v,z); 101371044d3cSBarry Smith /* BLASgemv_("N",&bs,&ncols,&_DOne,v,&bs,work,&_One,&_DZero,z,&_One); */ 10142d61bbb3SSatish Balay v += n*bs2; 101526e093fcSHong Zhang if (!usecprow) z += bs; 10162d61bbb3SSatish Balay } 10171ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 101826e093fcSHong Zhang ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr); 1019dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*a->nz*bs2 - bs*nonzerorow);CHKERRQ(ierr); 10202d61bbb3SSatish Balay PetscFunctionReturn(0); 10212d61bbb3SSatish Balay } 10222d61bbb3SSatish Balay 10236a65c766SHong Zhang extern PetscErrorCode VecCopy_Seq(Vec,Vec); 10244a2ae208SSatish Balay #undef __FUNCT__ 10254a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_1" 1026dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_1(Mat A,Vec xx,Vec yy,Vec zz) 10272d61bbb3SSatish Balay { 10282d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1029122f12eaSBarry Smith const PetscScalar *x; 1030122f12eaSBarry Smith PetscScalar *y,*z,sum; 1031122f12eaSBarry Smith const MatScalar *v; 1032dfbe8321SBarry Smith PetscErrorCode ierr; 1033122f12eaSBarry Smith PetscInt mbs=a->mbs,i,n,*ridx=PETSC_NULL,nonzerorow=0; 1034122f12eaSBarry Smith const PetscInt *idx,*ii; 103526e093fcSHong Zhang PetscTruth usecprow=a->compressedrow.use; 10362d61bbb3SSatish Balay 10372d61bbb3SSatish Balay PetscFunctionBegin; 1038122f12eaSBarry Smith ierr = VecGetArray(xx,(PetscScalar**)&x);CHKERRQ(ierr); 1039122f12eaSBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 10402e8a6d31SBarry Smith if (zz != yy) { 1041122f12eaSBarry Smith ierr = VecGetArray(zz,&z);CHKERRQ(ierr); 10422e8a6d31SBarry Smith } else { 1043122f12eaSBarry Smith z = y; 10442e8a6d31SBarry Smith } 10452d61bbb3SSatish Balay 10462d61bbb3SSatish Balay idx = a->j; 10472d61bbb3SSatish Balay v = a->a; 104826e093fcSHong Zhang if (usecprow){ 10494eb6d288SHong Zhang if (zz != yy){ 1050122f12eaSBarry Smith ierr = PetscMemcpy(z,y,mbs*sizeof(PetscScalar));CHKERRQ(ierr); 10514eb6d288SHong Zhang } 105226e093fcSHong Zhang mbs = a->compressedrow.nrows; 105326e093fcSHong Zhang ii = a->compressedrow.i; 10547b2bb3b9SHong Zhang ridx = a->compressedrow.rindex; 105526e093fcSHong Zhang } else { 10562d61bbb3SSatish Balay ii = a->i; 105726e093fcSHong Zhang } 10582d61bbb3SSatish Balay 10592d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 1060122f12eaSBarry Smith n = ii[1] - ii[0]; 1061122f12eaSBarry Smith ii++; 106226e093fcSHong Zhang if (!usecprow){ 1063122f12eaSBarry Smith nonzerorow += (n>0); 1064122f12eaSBarry Smith sum = y[i]; 1065122f12eaSBarry Smith } else { 1066122f12eaSBarry Smith sum = y[ridx[i]]; 1067122f12eaSBarry Smith } 1068122f12eaSBarry Smith PetscSparseDensePlusDot(sum,x,v,idx,n); 1069122f12eaSBarry Smith v += n; 1070122f12eaSBarry Smith idx += n; 1071122f12eaSBarry Smith if (usecprow){ 1072122f12eaSBarry Smith z[ridx[i]] = sum; 1073122f12eaSBarry Smith } else { 1074122f12eaSBarry Smith z[i] = sum; 107526e093fcSHong Zhang } 10762d61bbb3SSatish Balay } 1077122f12eaSBarry Smith ierr = VecRestoreArray(xx,(PetscScalar**)&x);CHKERRQ(ierr); 1078122f12eaSBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 10792e8a6d31SBarry Smith if (zz != yy) { 1080122f12eaSBarry Smith ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr); 10812e8a6d31SBarry Smith } 1082122f12eaSBarry Smith ierr = PetscLogFlops(2.0*a->nz - nonzerorow);CHKERRQ(ierr); 10832d61bbb3SSatish Balay PetscFunctionReturn(0); 10842d61bbb3SSatish Balay } 10852d61bbb3SSatish Balay 10864a2ae208SSatish Balay #undef __FUNCT__ 10874a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_2" 1088dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_2(Mat A,Vec xx,Vec yy,Vec zz) 10892d61bbb3SSatish Balay { 10902d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1091dcf5cc72SBarry Smith PetscScalar *x,*y = 0,*z = 0,*xb,sum1,sum2; 109226e093fcSHong Zhang PetscScalar x1,x2,*yarray,*zarray; 10933f1db9ecSBarry Smith MatScalar *v; 1094dfbe8321SBarry Smith PetscErrorCode ierr; 10954eb6d288SHong Zhang PetscInt mbs=a->mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL; 109626e093fcSHong Zhang PetscTruth usecprow=a->compressedrow.use; 10972d61bbb3SSatish Balay 10982d61bbb3SSatish Balay PetscFunctionBegin; 10991ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 110026e093fcSHong Zhang ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr); 11012e8a6d31SBarry Smith if (zz != yy) { 110226e093fcSHong Zhang ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr); 11032e8a6d31SBarry Smith } else { 110426e093fcSHong Zhang zarray = yarray; 11052e8a6d31SBarry Smith } 11062d61bbb3SSatish Balay 11072d61bbb3SSatish Balay idx = a->j; 11082d61bbb3SSatish Balay v = a->a; 110926e093fcSHong Zhang if (usecprow){ 11104eb6d288SHong Zhang if (zz != yy){ 11114eb6d288SHong Zhang ierr = PetscMemcpy(zarray,yarray,2*mbs*sizeof(PetscScalar));CHKERRQ(ierr); 11124eb6d288SHong Zhang } 111326e093fcSHong Zhang mbs = a->compressedrow.nrows; 111426e093fcSHong Zhang ii = a->compressedrow.i; 11157b2bb3b9SHong Zhang ridx = a->compressedrow.rindex; 11164eb6d288SHong Zhang if (zz != yy){ 11174eb6d288SHong Zhang ierr = PetscMemcpy(zarray,yarray,a->mbs*sizeof(PetscScalar));CHKERRQ(ierr); 11184eb6d288SHong Zhang } 111926e093fcSHong Zhang } else { 11202d61bbb3SSatish Balay ii = a->i; 112126e093fcSHong Zhang y = yarray; 112226e093fcSHong Zhang z = zarray; 112326e093fcSHong Zhang } 11242d61bbb3SSatish Balay 11252d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 11262d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 112726e093fcSHong Zhang if (usecprow){ 11287b2bb3b9SHong Zhang z = zarray + 2*ridx[i]; 11297b2bb3b9SHong Zhang y = yarray + 2*ridx[i]; 113026e093fcSHong Zhang } 11312d61bbb3SSatish Balay sum1 = y[0]; sum2 = y[1]; 11322d61bbb3SSatish Balay for (j=0; j<n; j++) { 11332d61bbb3SSatish Balay xb = x + 2*(*idx++); x1 = xb[0]; x2 = xb[1]; 11342d61bbb3SSatish Balay sum1 += v[0]*x1 + v[2]*x2; 11352d61bbb3SSatish Balay sum2 += v[1]*x1 + v[3]*x2; 11362d61bbb3SSatish Balay v += 4; 11372d61bbb3SSatish Balay } 11382d61bbb3SSatish Balay z[0] = sum1; z[1] = sum2; 113926e093fcSHong Zhang if (!usecprow){ 11402d61bbb3SSatish Balay z += 2; y += 2; 11412d61bbb3SSatish Balay } 114226e093fcSHong Zhang } 11431ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 114426e093fcSHong Zhang ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr); 11452e8a6d31SBarry Smith if (zz != yy) { 114626e093fcSHong Zhang ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr); 11472e8a6d31SBarry Smith } 1148dc0b31edSSatish Balay ierr = PetscLogFlops(4.0*a->nz);CHKERRQ(ierr); 11492d61bbb3SSatish Balay PetscFunctionReturn(0); 11502d61bbb3SSatish Balay } 11512d61bbb3SSatish Balay 11524a2ae208SSatish Balay #undef __FUNCT__ 11534a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_3" 1154dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_3(Mat A,Vec xx,Vec yy,Vec zz) 11552d61bbb3SSatish Balay { 11562d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1157dcf5cc72SBarry Smith PetscScalar *x,*y = 0,*z = 0,*xb,sum1,sum2,sum3,x1,x2,x3,*yarray,*zarray; 11583f1db9ecSBarry Smith MatScalar *v; 1159dfbe8321SBarry Smith PetscErrorCode ierr; 11604eb6d288SHong Zhang PetscInt mbs=a->mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL; 116126e093fcSHong Zhang PetscTruth usecprow=a->compressedrow.use; 11622d61bbb3SSatish Balay 11632d61bbb3SSatish Balay PetscFunctionBegin; 11641ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 116526e093fcSHong Zhang ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr); 11662e8a6d31SBarry Smith if (zz != yy) { 116726e093fcSHong Zhang ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr); 11682e8a6d31SBarry Smith } else { 116926e093fcSHong Zhang zarray = yarray; 11702e8a6d31SBarry Smith } 11712d61bbb3SSatish Balay 11722d61bbb3SSatish Balay idx = a->j; 11732d61bbb3SSatish Balay v = a->a; 117426e093fcSHong Zhang if (usecprow){ 11754eb6d288SHong Zhang if (zz != yy){ 11764eb6d288SHong Zhang ierr = PetscMemcpy(zarray,yarray,3*mbs*sizeof(PetscScalar));CHKERRQ(ierr); 11774eb6d288SHong Zhang } 117826e093fcSHong Zhang mbs = a->compressedrow.nrows; 117926e093fcSHong Zhang ii = a->compressedrow.i; 11807b2bb3b9SHong Zhang ridx = a->compressedrow.rindex; 118126e093fcSHong Zhang } else { 11822d61bbb3SSatish Balay ii = a->i; 118326e093fcSHong Zhang y = yarray; 118426e093fcSHong Zhang z = zarray; 118526e093fcSHong Zhang } 11862d61bbb3SSatish Balay 11872d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 11882d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 118926e093fcSHong Zhang if (usecprow){ 11907b2bb3b9SHong Zhang z = zarray + 3*ridx[i]; 11917b2bb3b9SHong Zhang y = yarray + 3*ridx[i]; 119226e093fcSHong Zhang } 11932d61bbb3SSatish Balay sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; 11942d61bbb3SSatish Balay for (j=0; j<n; j++) { 11952d61bbb3SSatish Balay xb = x + 3*(*idx++); x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; 11962d61bbb3SSatish Balay sum1 += v[0]*x1 + v[3]*x2 + v[6]*x3; 11972d61bbb3SSatish Balay sum2 += v[1]*x1 + v[4]*x2 + v[7]*x3; 11982d61bbb3SSatish Balay sum3 += v[2]*x1 + v[5]*x2 + v[8]*x3; 11992d61bbb3SSatish Balay v += 9; 12002d61bbb3SSatish Balay } 12012d61bbb3SSatish Balay z[0] = sum1; z[1] = sum2; z[2] = sum3; 120226e093fcSHong Zhang if (!usecprow){ 12032d61bbb3SSatish Balay z += 3; y += 3; 12042d61bbb3SSatish Balay } 120526e093fcSHong Zhang } 12061ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 120726e093fcSHong Zhang ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr); 12082e8a6d31SBarry Smith if (zz != yy) { 120926e093fcSHong Zhang ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr); 12102e8a6d31SBarry Smith } 1211dc0b31edSSatish Balay ierr = PetscLogFlops(18.0*a->nz);CHKERRQ(ierr); 12122d61bbb3SSatish Balay PetscFunctionReturn(0); 12132d61bbb3SSatish Balay } 12142d61bbb3SSatish Balay 12154a2ae208SSatish Balay #undef __FUNCT__ 12164a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_4" 1217dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_4(Mat A,Vec xx,Vec yy,Vec zz) 12182d61bbb3SSatish Balay { 12192d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1220dcf5cc72SBarry Smith PetscScalar *x,*y = 0,*z = 0,*xb,sum1,sum2,sum3,sum4,x1,x2,x3,x4,*yarray,*zarray; 12213f1db9ecSBarry Smith MatScalar *v; 1222dfbe8321SBarry Smith PetscErrorCode ierr; 12234eb6d288SHong Zhang PetscInt mbs=a->mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL; 122426e093fcSHong Zhang PetscTruth usecprow=a->compressedrow.use; 12252d61bbb3SSatish Balay 12262d61bbb3SSatish Balay PetscFunctionBegin; 12271ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 122826e093fcSHong Zhang ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr); 12292e8a6d31SBarry Smith if (zz != yy) { 123026e093fcSHong Zhang ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr); 12312e8a6d31SBarry Smith } else { 123226e093fcSHong Zhang zarray = yarray; 12332e8a6d31SBarry Smith } 12342d61bbb3SSatish Balay 12352d61bbb3SSatish Balay idx = a->j; 12362d61bbb3SSatish Balay v = a->a; 123726e093fcSHong Zhang if (usecprow){ 12384eb6d288SHong Zhang if (zz != yy){ 12394eb6d288SHong Zhang ierr = PetscMemcpy(zarray,yarray,4*mbs*sizeof(PetscScalar));CHKERRQ(ierr); 12404eb6d288SHong Zhang } 124126e093fcSHong Zhang mbs = a->compressedrow.nrows; 124226e093fcSHong Zhang ii = a->compressedrow.i; 12437b2bb3b9SHong Zhang ridx = a->compressedrow.rindex; 124426e093fcSHong Zhang } else { 12452d61bbb3SSatish Balay ii = a->i; 124626e093fcSHong Zhang y = yarray; 124726e093fcSHong Zhang z = zarray; 124826e093fcSHong Zhang } 12492d61bbb3SSatish Balay 12502d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 12512d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 125226e093fcSHong Zhang if (usecprow){ 12537b2bb3b9SHong Zhang z = zarray + 4*ridx[i]; 12547b2bb3b9SHong Zhang y = yarray + 4*ridx[i]; 125526e093fcSHong Zhang } 12562d61bbb3SSatish Balay sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; 12572d61bbb3SSatish Balay for (j=0; j<n; j++) { 12582d61bbb3SSatish Balay xb = x + 4*(*idx++); 12592d61bbb3SSatish Balay x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; 12602d61bbb3SSatish Balay sum1 += v[0]*x1 + v[4]*x2 + v[8]*x3 + v[12]*x4; 12612d61bbb3SSatish Balay sum2 += v[1]*x1 + v[5]*x2 + v[9]*x3 + v[13]*x4; 12622d61bbb3SSatish Balay sum3 += v[2]*x1 + v[6]*x2 + v[10]*x3 + v[14]*x4; 12632d61bbb3SSatish Balay sum4 += v[3]*x1 + v[7]*x2 + v[11]*x3 + v[15]*x4; 12642d61bbb3SSatish Balay v += 16; 12652d61bbb3SSatish Balay } 12662d61bbb3SSatish Balay z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; 126726e093fcSHong Zhang if (!usecprow){ 12682d61bbb3SSatish Balay z += 4; y += 4; 12692d61bbb3SSatish Balay } 127026e093fcSHong Zhang } 12711ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 127226e093fcSHong Zhang ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr); 12732e8a6d31SBarry Smith if (zz != yy) { 127426e093fcSHong Zhang ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr); 12752e8a6d31SBarry Smith } 1276dc0b31edSSatish Balay ierr = PetscLogFlops(32.0*a->nz);CHKERRQ(ierr); 12772d61bbb3SSatish Balay PetscFunctionReturn(0); 12782d61bbb3SSatish Balay } 12792d61bbb3SSatish Balay 12804a2ae208SSatish Balay #undef __FUNCT__ 12814a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_5" 1282dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_5(Mat A,Vec xx,Vec yy,Vec zz) 12832d61bbb3SSatish Balay { 12842d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1285dcf5cc72SBarry Smith PetscScalar *x,*y = 0,*z = 0,*xb,sum1,sum2,sum3,sum4,sum5,x1,x2,x3,x4,x5; 128626e093fcSHong Zhang PetscScalar *yarray,*zarray; 12873f1db9ecSBarry Smith MatScalar *v; 1288dfbe8321SBarry Smith PetscErrorCode ierr; 12894eb6d288SHong Zhang PetscInt mbs=a->mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL; 129026e093fcSHong Zhang PetscTruth usecprow=a->compressedrow.use; 12912d61bbb3SSatish Balay 12922d61bbb3SSatish Balay PetscFunctionBegin; 12931ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 129426e093fcSHong Zhang ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr); 12952e8a6d31SBarry Smith if (zz != yy) { 129626e093fcSHong Zhang ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr); 12972e8a6d31SBarry Smith } else { 129826e093fcSHong Zhang zarray = yarray; 12992e8a6d31SBarry Smith } 13002d61bbb3SSatish Balay 13012d61bbb3SSatish Balay idx = a->j; 13022d61bbb3SSatish Balay v = a->a; 130326e093fcSHong Zhang if (usecprow){ 13044eb6d288SHong Zhang if (zz != yy){ 13054eb6d288SHong Zhang ierr = PetscMemcpy(zarray,yarray,5*mbs*sizeof(PetscScalar));CHKERRQ(ierr); 13064eb6d288SHong Zhang } 130726e093fcSHong Zhang mbs = a->compressedrow.nrows; 130826e093fcSHong Zhang ii = a->compressedrow.i; 13097b2bb3b9SHong Zhang ridx = a->compressedrow.rindex; 131026e093fcSHong Zhang } else { 13112d61bbb3SSatish Balay ii = a->i; 131226e093fcSHong Zhang y = yarray; 131326e093fcSHong Zhang z = zarray; 131426e093fcSHong Zhang } 13152d61bbb3SSatish Balay 13162d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 13172d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 131826e093fcSHong Zhang if (usecprow){ 13197b2bb3b9SHong Zhang z = zarray + 5*ridx[i]; 13207b2bb3b9SHong Zhang y = yarray + 5*ridx[i]; 132126e093fcSHong Zhang } 13222d61bbb3SSatish Balay sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4]; 13232d61bbb3SSatish Balay for (j=0; j<n; j++) { 13242d61bbb3SSatish Balay xb = x + 5*(*idx++); 13252d61bbb3SSatish Balay x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; 13262d61bbb3SSatish Balay sum1 += v[0]*x1 + v[5]*x2 + v[10]*x3 + v[15]*x4 + v[20]*x5; 13272d61bbb3SSatish Balay sum2 += v[1]*x1 + v[6]*x2 + v[11]*x3 + v[16]*x4 + v[21]*x5; 13282d61bbb3SSatish Balay sum3 += v[2]*x1 + v[7]*x2 + v[12]*x3 + v[17]*x4 + v[22]*x5; 13292d61bbb3SSatish Balay sum4 += v[3]*x1 + v[8]*x2 + v[13]*x3 + v[18]*x4 + v[23]*x5; 13302d61bbb3SSatish Balay sum5 += v[4]*x1 + v[9]*x2 + v[14]*x3 + v[19]*x4 + v[24]*x5; 13312d61bbb3SSatish Balay v += 25; 13322d61bbb3SSatish Balay } 13332d61bbb3SSatish Balay z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; 133426e093fcSHong Zhang if (!usecprow){ 13352d61bbb3SSatish Balay z += 5; y += 5; 13362d61bbb3SSatish Balay } 133726e093fcSHong Zhang } 13381ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 133926e093fcSHong Zhang ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr); 13402e8a6d31SBarry Smith if (zz != yy) { 134126e093fcSHong Zhang ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr); 13422e8a6d31SBarry Smith } 1343dc0b31edSSatish Balay ierr = PetscLogFlops(50.0*a->nz);CHKERRQ(ierr); 13442d61bbb3SSatish Balay PetscFunctionReturn(0); 13452d61bbb3SSatish Balay } 13464a2ae208SSatish Balay #undef __FUNCT__ 13474a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_6" 1348dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_6(Mat A,Vec xx,Vec yy,Vec zz) 134915091d37SBarry Smith { 135015091d37SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1351dcf5cc72SBarry Smith PetscScalar *x,*y = 0,*z = 0,*xb,sum1,sum2,sum3,sum4,sum5,sum6; 135226e093fcSHong Zhang PetscScalar x1,x2,x3,x4,x5,x6,*yarray,*zarray; 135315091d37SBarry Smith MatScalar *v; 1354dfbe8321SBarry Smith PetscErrorCode ierr; 13554eb6d288SHong Zhang PetscInt mbs=a->mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL; 135626e093fcSHong Zhang PetscTruth usecprow=a->compressedrow.use; 135715091d37SBarry Smith 135815091d37SBarry Smith PetscFunctionBegin; 13591ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 136026e093fcSHong Zhang ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr); 136115091d37SBarry Smith if (zz != yy) { 136226e093fcSHong Zhang ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr); 136315091d37SBarry Smith } else { 136426e093fcSHong Zhang zarray = yarray; 136515091d37SBarry Smith } 136615091d37SBarry Smith 136715091d37SBarry Smith idx = a->j; 136815091d37SBarry Smith v = a->a; 136926e093fcSHong Zhang if (usecprow){ 13704eb6d288SHong Zhang if (zz != yy){ 13714eb6d288SHong Zhang ierr = PetscMemcpy(zarray,yarray,6*mbs*sizeof(PetscScalar));CHKERRQ(ierr); 13724eb6d288SHong Zhang } 137326e093fcSHong Zhang mbs = a->compressedrow.nrows; 137426e093fcSHong Zhang ii = a->compressedrow.i; 13757b2bb3b9SHong Zhang ridx = a->compressedrow.rindex; 137626e093fcSHong Zhang } else { 137715091d37SBarry Smith ii = a->i; 137826e093fcSHong Zhang y = yarray; 137926e093fcSHong Zhang z = zarray; 138026e093fcSHong Zhang } 138115091d37SBarry Smith 138215091d37SBarry Smith for (i=0; i<mbs; i++) { 138315091d37SBarry Smith n = ii[1] - ii[0]; ii++; 138426e093fcSHong Zhang if (usecprow){ 13857b2bb3b9SHong Zhang z = zarray + 6*ridx[i]; 13867b2bb3b9SHong Zhang y = yarray + 6*ridx[i]; 138726e093fcSHong Zhang } 138815091d37SBarry Smith sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4]; sum6 = y[5]; 138915091d37SBarry Smith for (j=0; j<n; j++) { 13903b95cb0eSSatish Balay xb = x + 6*(*idx++); 139115091d37SBarry Smith x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; 139215091d37SBarry Smith sum1 += v[0]*x1 + v[6]*x2 + v[12]*x3 + v[18]*x4 + v[24]*x5 + v[30]*x6; 139315091d37SBarry Smith sum2 += v[1]*x1 + v[7]*x2 + v[13]*x3 + v[19]*x4 + v[25]*x5 + v[31]*x6; 139415091d37SBarry Smith sum3 += v[2]*x1 + v[8]*x2 + v[14]*x3 + v[20]*x4 + v[26]*x5 + v[32]*x6; 139515091d37SBarry Smith sum4 += v[3]*x1 + v[9]*x2 + v[15]*x3 + v[21]*x4 + v[27]*x5 + v[33]*x6; 139615091d37SBarry Smith sum5 += v[4]*x1 + v[10]*x2 + v[16]*x3 + v[22]*x4 + v[28]*x5 + v[34]*x6; 139715091d37SBarry Smith sum6 += v[5]*x1 + v[11]*x2 + v[17]*x3 + v[23]*x4 + v[29]*x5 + v[35]*x6; 139815091d37SBarry Smith v += 36; 139915091d37SBarry Smith } 140015091d37SBarry Smith z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; 140126e093fcSHong Zhang if (!usecprow){ 140215091d37SBarry Smith z += 6; y += 6; 140315091d37SBarry Smith } 140426e093fcSHong Zhang } 14051ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 140626e093fcSHong Zhang ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr); 140715091d37SBarry Smith if (zz != yy) { 140826e093fcSHong Zhang ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr); 140915091d37SBarry Smith } 1410dc0b31edSSatish Balay ierr = PetscLogFlops(72.0*a->nz);CHKERRQ(ierr); 141115091d37SBarry Smith PetscFunctionReturn(0); 141215091d37SBarry Smith } 14132d61bbb3SSatish Balay 14144a2ae208SSatish Balay #undef __FUNCT__ 14154a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_7" 1416dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_7(Mat A,Vec xx,Vec yy,Vec zz) 14172d61bbb3SSatish Balay { 14182d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1419dcf5cc72SBarry Smith PetscScalar *x,*y = 0,*z = 0,*xb,sum1,sum2,sum3,sum4,sum5,sum6,sum7; 142026e093fcSHong Zhang PetscScalar x1,x2,x3,x4,x5,x6,x7,*yarray,*zarray; 14213f1db9ecSBarry Smith MatScalar *v; 1422dfbe8321SBarry Smith PetscErrorCode ierr; 14237b2bb3b9SHong Zhang PetscInt mbs=a->mbs,i,*idx,*ii,j,n,*ridx=PETSC_NULL; 142426e093fcSHong Zhang PetscTruth usecprow=a->compressedrow.use; 14252d61bbb3SSatish Balay 14262d61bbb3SSatish Balay PetscFunctionBegin; 14271ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 142826e093fcSHong Zhang ierr = VecGetArray(yy,&yarray);CHKERRQ(ierr); 14292e8a6d31SBarry Smith if (zz != yy) { 143026e093fcSHong Zhang ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr); 14312e8a6d31SBarry Smith } else { 143226e093fcSHong Zhang zarray = yarray; 14332e8a6d31SBarry Smith } 14342d61bbb3SSatish Balay 14352d61bbb3SSatish Balay idx = a->j; 14362d61bbb3SSatish Balay v = a->a; 143726e093fcSHong Zhang if (usecprow){ 14384eb6d288SHong Zhang if (zz != yy){ 14394eb6d288SHong Zhang ierr = PetscMemcpy(zarray,yarray,7*mbs*sizeof(PetscScalar));CHKERRQ(ierr); 14404eb6d288SHong Zhang } 144126e093fcSHong Zhang mbs = a->compressedrow.nrows; 144226e093fcSHong Zhang ii = a->compressedrow.i; 14437b2bb3b9SHong Zhang ridx = a->compressedrow.rindex; 144426e093fcSHong Zhang } else { 14452d61bbb3SSatish Balay ii = a->i; 144626e093fcSHong Zhang y = yarray; 144726e093fcSHong Zhang z = zarray; 144826e093fcSHong Zhang } 14492d61bbb3SSatish Balay 14502d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 14512d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 145226e093fcSHong Zhang if (usecprow){ 14537b2bb3b9SHong Zhang z = zarray + 7*ridx[i]; 14547b2bb3b9SHong Zhang y = yarray + 7*ridx[i]; 145526e093fcSHong Zhang } 14562d61bbb3SSatish Balay sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4]; sum6 = y[5]; sum7 = y[6]; 14572d61bbb3SSatish Balay for (j=0; j<n; j++) { 14582d61bbb3SSatish Balay xb = x + 7*(*idx++); 14592d61bbb3SSatish Balay x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6]; 14602d61bbb3SSatish Balay sum1 += v[0]*x1 + v[7]*x2 + v[14]*x3 + v[21]*x4 + v[28]*x5 + v[35]*x6 + v[42]*x7; 14612d61bbb3SSatish Balay sum2 += v[1]*x1 + v[8]*x2 + v[15]*x3 + v[22]*x4 + v[29]*x5 + v[36]*x6 + v[43]*x7; 14622d61bbb3SSatish Balay sum3 += v[2]*x1 + v[9]*x2 + v[16]*x3 + v[23]*x4 + v[30]*x5 + v[37]*x6 + v[44]*x7; 14632d61bbb3SSatish Balay sum4 += v[3]*x1 + v[10]*x2 + v[17]*x3 + v[24]*x4 + v[31]*x5 + v[38]*x6 + v[45]*x7; 14642d61bbb3SSatish Balay sum5 += v[4]*x1 + v[11]*x2 + v[18]*x3 + v[25]*x4 + v[32]*x5 + v[39]*x6 + v[46]*x7; 14652d61bbb3SSatish Balay sum6 += v[5]*x1 + v[12]*x2 + v[19]*x3 + v[26]*x4 + v[33]*x5 + v[40]*x6 + v[47]*x7; 14662d61bbb3SSatish Balay sum7 += v[6]*x1 + v[13]*x2 + v[20]*x3 + v[27]*x4 + v[34]*x5 + v[41]*x6 + v[48]*x7; 14672d61bbb3SSatish Balay v += 49; 14682d61bbb3SSatish Balay } 14692d61bbb3SSatish Balay z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7; 147026e093fcSHong Zhang if (!usecprow){ 14712d61bbb3SSatish Balay z += 7; y += 7; 14722d61bbb3SSatish Balay } 147326e093fcSHong Zhang } 14741ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 147526e093fcSHong Zhang ierr = VecRestoreArray(yy,&yarray);CHKERRQ(ierr); 14762e8a6d31SBarry Smith if (zz != yy) { 147726e093fcSHong Zhang ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr); 14782e8a6d31SBarry Smith } 1479dc0b31edSSatish Balay ierr = PetscLogFlops(98.0*a->nz);CHKERRQ(ierr); 14802d61bbb3SSatish Balay PetscFunctionReturn(0); 14812d61bbb3SSatish Balay } 1482218c64b6SSatish Balay 14834a2ae208SSatish Balay #undef __FUNCT__ 14844a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_N" 1485dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqBAIJ_N(Mat A,Vec xx,Vec yy,Vec zz) 14862d61bbb3SSatish Balay { 14872d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1488bba805e6SBarry Smith PetscScalar *x,*z = 0,*xb,*work,*workt,*zarray; 14893f1db9ecSBarry Smith MatScalar *v; 14906849ba73SBarry Smith PetscErrorCode ierr; 1491d0f46423SBarry Smith PetscInt mbs,i,*idx,*ii,bs=A->rmap->bs,j,n,bs2=a->bs2; 14927b2bb3b9SHong Zhang PetscInt ncols,k,*ridx=PETSC_NULL; 149326e093fcSHong Zhang PetscTruth usecprow=a->compressedrow.use; 1494218c64b6SSatish Balay 14952d61bbb3SSatish Balay PetscFunctionBegin; 14966a65c766SHong Zhang ierr = VecCopy_Seq(yy,zz);CHKERRQ(ierr); 14971ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 149826e093fcSHong Zhang ierr = VecGetArray(zz,&zarray);CHKERRQ(ierr); 14992d61bbb3SSatish Balay 15002d61bbb3SSatish Balay idx = a->j; 15012d61bbb3SSatish Balay v = a->a; 150226e093fcSHong Zhang if (usecprow){ 150326e093fcSHong Zhang mbs = a->compressedrow.nrows; 150426e093fcSHong Zhang ii = a->compressedrow.i; 15057b2bb3b9SHong Zhang ridx = a->compressedrow.rindex; 150626e093fcSHong Zhang } else { 150726e093fcSHong Zhang mbs = a->mbs; 15082d61bbb3SSatish Balay ii = a->i; 150926e093fcSHong Zhang z = zarray; 151026e093fcSHong Zhang } 15112d61bbb3SSatish Balay 15122d61bbb3SSatish Balay if (!a->mult_work) { 1513d0f46423SBarry Smith k = PetscMax(A->rmap->n,A->cmap->n); 151487828ca2SBarry Smith ierr = PetscMalloc((k+1)*sizeof(PetscScalar),&a->mult_work);CHKERRQ(ierr); 15152d61bbb3SSatish Balay } 15162d61bbb3SSatish Balay work = a->mult_work; 15172d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 15182d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 15192d61bbb3SSatish Balay ncols = n*bs; 15202d61bbb3SSatish Balay workt = work; 15212d61bbb3SSatish Balay for (j=0; j<n; j++) { 15222d61bbb3SSatish Balay xb = x + bs*(*idx++); 15232d61bbb3SSatish Balay for (k=0; k<bs; k++) workt[k] = xb[k]; 15242d61bbb3SSatish Balay workt += bs; 15252d61bbb3SSatish Balay } 15267b2bb3b9SHong Zhang if (usecprow) z = zarray + bs*ridx[i]; 15273f1db9ecSBarry Smith Kernel_w_gets_w_plus_Ar_times_v(bs,ncols,work,v,z); 152871044d3cSBarry Smith /* BLASgemv_("N",&bs,&ncols,&_DOne,v,&bs,work,&_One,&_DOne,z,&_One); */ 15292d61bbb3SSatish Balay v += n*bs2; 153026e093fcSHong Zhang if (!usecprow){ 15312d61bbb3SSatish Balay z += bs; 15322d61bbb3SSatish Balay } 153326e093fcSHong Zhang } 15341ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 153526e093fcSHong Zhang ierr = VecRestoreArray(zz,&zarray);CHKERRQ(ierr); 1536dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*a->nz*bs2);CHKERRQ(ierr); 15372d61bbb3SSatish Balay PetscFunctionReturn(0); 15382d61bbb3SSatish Balay } 15392d61bbb3SSatish Balay 15404a2ae208SSatish Balay #undef __FUNCT__ 1541547795f9SHong Zhang #define __FUNCT__ "MatMultHermitianTranspose_SeqBAIJ" 1542547795f9SHong Zhang PetscErrorCode MatMultHermitianTranspose_SeqBAIJ(Mat A,Vec xx,Vec zz) 1543547795f9SHong Zhang { 1544547795f9SHong Zhang PetscScalar zero = 0.0; 1545547795f9SHong Zhang PetscErrorCode ierr; 1546547795f9SHong Zhang 1547547795f9SHong Zhang PetscFunctionBegin; 1548547795f9SHong Zhang ierr = VecSet(zz,zero);CHKERRQ(ierr); 1549547795f9SHong Zhang ierr = MatMultHermitianTransposeAdd_SeqBAIJ(A,xx,zz,zz);CHKERRQ(ierr); 1550547795f9SHong Zhang PetscFunctionReturn(0); 1551547795f9SHong Zhang } 1552547795f9SHong Zhang 1553547795f9SHong Zhang #undef __FUNCT__ 15544a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_SeqBAIJ" 1555dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_SeqBAIJ(Mat A,Vec xx,Vec zz) 15562d61bbb3SSatish Balay { 15573447b6efSHong Zhang PetscScalar zero = 0.0; 15586849ba73SBarry Smith PetscErrorCode ierr; 15592d61bbb3SSatish Balay 15602d61bbb3SSatish Balay PetscFunctionBegin; 15612dcb1b2aSMatthew Knepley ierr = VecSet(zz,zero);CHKERRQ(ierr); 15623447b6efSHong Zhang ierr = MatMultTransposeAdd_SeqBAIJ(A,xx,zz,zz);CHKERRQ(ierr); 15632d61bbb3SSatish Balay PetscFunctionReturn(0); 15642d61bbb3SSatish Balay } 15652d61bbb3SSatish Balay 15664a2ae208SSatish Balay #undef __FUNCT__ 1567547795f9SHong Zhang #define __FUNCT__ "MatMultHermitianTransposeAdd_SeqBAIJ" 1568547795f9SHong Zhang PetscErrorCode MatMultHermitianTransposeAdd_SeqBAIJ(Mat A,Vec xx,Vec yy,Vec zz) 1569547795f9SHong Zhang 1570547795f9SHong Zhang { 1571547795f9SHong Zhang Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1572547795f9SHong Zhang PetscScalar *zb,*x,*z,*xb = 0,x1,x2,x3,x4,x5; 1573547795f9SHong Zhang MatScalar *v; 1574547795f9SHong Zhang PetscErrorCode ierr; 1575547795f9SHong Zhang PetscInt mbs,i,*idx,*ii,rval,bs=A->rmap->bs,j,n,bs2=a->bs2,*ib,*ridx=PETSC_NULL; 1576547795f9SHong Zhang Mat_CompressedRow cprow = a->compressedrow; 1577547795f9SHong Zhang PetscTruth usecprow=cprow.use; 1578547795f9SHong Zhang 1579547795f9SHong Zhang PetscFunctionBegin; 1580547795f9SHong Zhang if (yy != zz) { ierr = VecCopy_Seq(yy,zz);CHKERRQ(ierr); } 1581547795f9SHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 1582547795f9SHong Zhang ierr = VecGetArray(zz,&z);CHKERRQ(ierr); 1583547795f9SHong Zhang 1584547795f9SHong Zhang idx = a->j; 1585547795f9SHong Zhang v = a->a; 1586547795f9SHong Zhang if (usecprow){ 1587547795f9SHong Zhang mbs = cprow.nrows; 1588547795f9SHong Zhang ii = cprow.i; 1589547795f9SHong Zhang ridx = cprow.rindex; 1590547795f9SHong Zhang } else { 1591547795f9SHong Zhang mbs=a->mbs; 1592547795f9SHong Zhang ii = a->i; 1593547795f9SHong Zhang xb = x; 1594547795f9SHong Zhang } 1595547795f9SHong Zhang 1596547795f9SHong Zhang switch (bs) { 1597547795f9SHong Zhang case 1: 1598547795f9SHong Zhang for (i=0; i<mbs; i++) { 1599547795f9SHong Zhang if (usecprow) xb = x + ridx[i]; 1600547795f9SHong Zhang x1 = xb[0]; 1601547795f9SHong Zhang ib = idx + ii[0]; 1602547795f9SHong Zhang n = ii[1] - ii[0]; ii++; 1603547795f9SHong Zhang for (j=0; j<n; j++) { 1604547795f9SHong Zhang rval = ib[j]; 1605547795f9SHong Zhang z[rval] += PetscConj(*v) * x1; 1606547795f9SHong Zhang v++; 1607547795f9SHong Zhang } 1608547795f9SHong Zhang if (!usecprow) xb++; 1609547795f9SHong Zhang } 1610547795f9SHong Zhang break; 1611547795f9SHong Zhang case 2: 1612547795f9SHong Zhang for (i=0; i<mbs; i++) { 1613547795f9SHong Zhang if (usecprow) xb = x + 2*ridx[i]; 1614547795f9SHong Zhang x1 = xb[0]; x2 = xb[1]; 1615547795f9SHong Zhang ib = idx + ii[0]; 1616547795f9SHong Zhang n = ii[1] - ii[0]; ii++; 1617547795f9SHong Zhang for (j=0; j<n; j++) { 1618547795f9SHong Zhang rval = ib[j]*2; 1619547795f9SHong Zhang z[rval++] += PetscConj(v[0])*x1 + PetscConj(v[1])*x2; 1620547795f9SHong Zhang z[rval++] += PetscConj(v[2])*x1 + PetscConj(v[3])*x2; 1621547795f9SHong Zhang v += 4; 1622547795f9SHong Zhang } 1623547795f9SHong Zhang if (!usecprow) xb += 2; 1624547795f9SHong Zhang } 1625547795f9SHong Zhang break; 1626547795f9SHong Zhang case 3: 1627547795f9SHong Zhang for (i=0; i<mbs; i++) { 1628547795f9SHong Zhang if (usecprow) xb = x + 3*ridx[i]; 1629547795f9SHong Zhang x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; 1630547795f9SHong Zhang ib = idx + ii[0]; 1631547795f9SHong Zhang n = ii[1] - ii[0]; ii++; 1632547795f9SHong Zhang for (j=0; j<n; j++) { 1633547795f9SHong Zhang rval = ib[j]*3; 1634547795f9SHong Zhang z[rval++] += PetscConj(v[0])*x1 + PetscConj(v[1])*x2 + PetscConj(v[2])*x3; 1635547795f9SHong Zhang z[rval++] += PetscConj(v[3])*x1 + PetscConj(v[4])*x2 + PetscConj(v[5])*x3; 1636547795f9SHong Zhang z[rval++] += PetscConj(v[6])*x1 + PetscConj(v[7])*x2 + PetscConj(v[8])*x3; 1637547795f9SHong Zhang v += 9; 1638547795f9SHong Zhang } 1639547795f9SHong Zhang if (!usecprow) xb += 3; 1640547795f9SHong Zhang } 1641547795f9SHong Zhang break; 1642547795f9SHong Zhang case 4: 1643547795f9SHong Zhang for (i=0; i<mbs; i++) { 1644547795f9SHong Zhang if (usecprow) xb = x + 4*ridx[i]; 1645547795f9SHong Zhang x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; 1646547795f9SHong Zhang ib = idx + ii[0]; 1647547795f9SHong Zhang n = ii[1] - ii[0]; ii++; 1648547795f9SHong Zhang for (j=0; j<n; j++) { 1649547795f9SHong Zhang rval = ib[j]*4; 1650547795f9SHong Zhang z[rval++] += PetscConj(v[0])*x1 + PetscConj(v[1])*x2 + PetscConj(v[2])*x3 + PetscConj(v[3])*x4; 1651547795f9SHong Zhang z[rval++] += PetscConj(v[4])*x1 + PetscConj(v[5])*x2 + PetscConj(v[6])*x3 + PetscConj(v[7])*x4; 1652547795f9SHong Zhang z[rval++] += PetscConj(v[8])*x1 + PetscConj(v[9])*x2 + PetscConj(v[10])*x3 + PetscConj(v[11])*x4; 1653547795f9SHong Zhang z[rval++] += PetscConj(v[12])*x1 + PetscConj(v[13])*x2 + PetscConj(v[14])*x3 + PetscConj(v[15])*x4; 1654547795f9SHong Zhang v += 16; 1655547795f9SHong Zhang } 1656547795f9SHong Zhang if (!usecprow) xb += 4; 1657547795f9SHong Zhang } 1658547795f9SHong Zhang break; 1659547795f9SHong Zhang case 5: 1660547795f9SHong Zhang for (i=0; i<mbs; i++) { 1661547795f9SHong Zhang if (usecprow) xb = x + 5*ridx[i]; 1662547795f9SHong Zhang x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; 1663547795f9SHong Zhang x4 = xb[3]; x5 = xb[4]; 1664547795f9SHong Zhang ib = idx + ii[0]; 1665547795f9SHong Zhang n = ii[1] - ii[0]; ii++; 1666547795f9SHong Zhang for (j=0; j<n; j++) { 1667547795f9SHong Zhang rval = ib[j]*5; 1668547795f9SHong Zhang z[rval++] += PetscConj(v[0])*x1 + PetscConj(v[1])*x2 + PetscConj(v[2])*x3 + PetscConj(v[3])*x4 + PetscConj(v[4])*x5; 1669547795f9SHong Zhang z[rval++] += PetscConj(v[5])*x1 + PetscConj(v[6])*x2 + PetscConj(v[7])*x3 + PetscConj(v[8])*x4 + PetscConj(v[9])*x5; 1670547795f9SHong Zhang z[rval++] += PetscConj(v[10])*x1 + PetscConj(v[11])*x2 + PetscConj(v[12])*x3 + PetscConj(v[13])*x4 + PetscConj(v[14])*x5; 1671547795f9SHong Zhang z[rval++] += PetscConj(v[15])*x1 + PetscConj(v[16])*x2 + PetscConj(v[17])*x3 + PetscConj(v[18])*x4 + PetscConj(v[19])*x5; 1672547795f9SHong Zhang z[rval++] += PetscConj(v[20])*x1 + PetscConj(v[21])*x2 + PetscConj(v[22])*x3 + PetscConj(v[23])*x4 + PetscConj(v[24])*x5; 1673547795f9SHong Zhang v += 25; 1674547795f9SHong Zhang } 1675547795f9SHong Zhang if (!usecprow) xb += 5; 1676547795f9SHong Zhang } 1677547795f9SHong Zhang break; 1678547795f9SHong Zhang default: { /* block sizes larger than 5 by 5 are handled by BLAS */ 1679547795f9SHong Zhang PetscInt ncols,k; 1680547795f9SHong Zhang PetscScalar *work,*workt,*xtmp; 1681547795f9SHong Zhang 1682547795f9SHong Zhang SETERRQ(PETSC_ERR_SUP,"block size larger than 5 is not supported yet"); 1683547795f9SHong Zhang if (!a->mult_work) { 1684547795f9SHong Zhang k = PetscMax(A->rmap->n,A->cmap->n); 1685547795f9SHong Zhang ierr = PetscMalloc((k+1)*sizeof(PetscScalar),&a->mult_work);CHKERRQ(ierr); 1686547795f9SHong Zhang } 1687547795f9SHong Zhang work = a->mult_work; 1688547795f9SHong Zhang xtmp = x; 1689547795f9SHong Zhang for (i=0; i<mbs; i++) { 1690547795f9SHong Zhang n = ii[1] - ii[0]; ii++; 1691547795f9SHong Zhang ncols = n*bs; 1692547795f9SHong Zhang ierr = PetscMemzero(work,ncols*sizeof(PetscScalar));CHKERRQ(ierr); 1693547795f9SHong Zhang if (usecprow) { 1694547795f9SHong Zhang xtmp = x + bs*ridx[i]; 1695547795f9SHong Zhang } 1696547795f9SHong Zhang Kernel_w_gets_w_plus_trans_Ar_times_v(bs,ncols,xtmp,v,work); 1697547795f9SHong Zhang /* BLASgemv_("T",&bs,&ncols,&_DOne,v,&bs,xtmp,&_One,&_DOne,work,&_One); */ 1698547795f9SHong Zhang v += n*bs2; 1699547795f9SHong Zhang if (!usecprow) xtmp += bs; 1700547795f9SHong Zhang workt = work; 1701547795f9SHong Zhang for (j=0; j<n; j++) { 1702547795f9SHong Zhang zb = z + bs*(*idx++); 1703547795f9SHong Zhang for (k=0; k<bs; k++) zb[k] += workt[k] ; 1704547795f9SHong Zhang workt += bs; 1705547795f9SHong Zhang } 1706547795f9SHong Zhang } 1707547795f9SHong Zhang } 1708547795f9SHong Zhang } 1709547795f9SHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 1710547795f9SHong Zhang ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr); 1711547795f9SHong Zhang ierr = PetscLogFlops(2.0*a->nz*a->bs2);CHKERRQ(ierr); 1712547795f9SHong Zhang PetscFunctionReturn(0); 1713547795f9SHong Zhang } 1714547795f9SHong Zhang 1715547795f9SHong Zhang #undef __FUNCT__ 17164a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_SeqBAIJ" 1717dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_SeqBAIJ(Mat A,Vec xx,Vec yy,Vec zz) 17182d61bbb3SSatish Balay { 17192d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1720dcf5cc72SBarry Smith PetscScalar *zb,*x,*z,*xb = 0,x1,x2,x3,x4,x5; 17213f1db9ecSBarry Smith MatScalar *v; 17226849ba73SBarry Smith PetscErrorCode ierr; 1723d0f46423SBarry Smith PetscInt mbs,i,*idx,*ii,rval,bs=A->rmap->bs,j,n,bs2=a->bs2,*ib,*ridx=PETSC_NULL; 17243447b6efSHong Zhang Mat_CompressedRow cprow = a->compressedrow; 17253447b6efSHong Zhang PetscTruth usecprow=cprow.use; 17262d61bbb3SSatish Balay 17272d61bbb3SSatish Balay PetscFunctionBegin; 17286a65c766SHong Zhang if (yy != zz) { ierr = VecCopy_Seq(yy,zz);CHKERRQ(ierr); } 17293447b6efSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 17303447b6efSHong Zhang ierr = VecGetArray(zz,&z);CHKERRQ(ierr); 17312d61bbb3SSatish Balay 17322d61bbb3SSatish Balay idx = a->j; 17332d61bbb3SSatish Balay v = a->a; 17343447b6efSHong Zhang if (usecprow){ 17353447b6efSHong Zhang mbs = cprow.nrows; 17363447b6efSHong Zhang ii = cprow.i; 17377b2bb3b9SHong Zhang ridx = cprow.rindex; 17383447b6efSHong Zhang } else { 17393447b6efSHong Zhang mbs=a->mbs; 17402d61bbb3SSatish Balay ii = a->i; 1741f1af5d2fSBarry Smith xb = x; 17423447b6efSHong Zhang } 17432d61bbb3SSatish Balay 17442d61bbb3SSatish Balay switch (bs) { 17452d61bbb3SSatish Balay case 1: 17462d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 17477b2bb3b9SHong Zhang if (usecprow) xb = x + ridx[i]; 1748f1af5d2fSBarry Smith x1 = xb[0]; 17493447b6efSHong Zhang ib = idx + ii[0]; 17503447b6efSHong Zhang n = ii[1] - ii[0]; ii++; 17512d61bbb3SSatish Balay for (j=0; j<n; j++) { 17522d61bbb3SSatish Balay rval = ib[j]; 1753f1af5d2fSBarry Smith z[rval] += *v * x1; 1754f1af5d2fSBarry Smith v++; 17552d61bbb3SSatish Balay } 17563447b6efSHong Zhang if (!usecprow) xb++; 17572d61bbb3SSatish Balay } 17582d61bbb3SSatish Balay break; 17592d61bbb3SSatish Balay case 2: 17602d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 17617b2bb3b9SHong Zhang if (usecprow) xb = x + 2*ridx[i]; 1762f1af5d2fSBarry Smith x1 = xb[0]; x2 = xb[1]; 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]*2; 17672d61bbb3SSatish Balay z[rval++] += v[0]*x1 + v[1]*x2; 17682d61bbb3SSatish Balay z[rval++] += v[2]*x1 + v[3]*x2; 17692d61bbb3SSatish Balay v += 4; 17702d61bbb3SSatish Balay } 17713447b6efSHong Zhang if (!usecprow) xb += 2; 17722d61bbb3SSatish Balay } 17732d61bbb3SSatish Balay break; 17742d61bbb3SSatish Balay case 3: 17752d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 17767b2bb3b9SHong Zhang if (usecprow) xb = x + 3*ridx[i]; 1777f1af5d2fSBarry Smith x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; 17783447b6efSHong Zhang ib = idx + ii[0]; 17793447b6efSHong Zhang n = ii[1] - ii[0]; ii++; 17802d61bbb3SSatish Balay for (j=0; j<n; j++) { 17812d61bbb3SSatish Balay rval = ib[j]*3; 17822d61bbb3SSatish Balay z[rval++] += v[0]*x1 + v[1]*x2 + v[2]*x3; 17832d61bbb3SSatish Balay z[rval++] += v[3]*x1 + v[4]*x2 + v[5]*x3; 17842d61bbb3SSatish Balay z[rval++] += v[6]*x1 + v[7]*x2 + v[8]*x3; 17852d61bbb3SSatish Balay v += 9; 17862d61bbb3SSatish Balay } 17873447b6efSHong Zhang if (!usecprow) xb += 3; 17882d61bbb3SSatish Balay } 17892d61bbb3SSatish Balay break; 17902d61bbb3SSatish Balay case 4: 17912d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 17927b2bb3b9SHong Zhang if (usecprow) xb = x + 4*ridx[i]; 1793f1af5d2fSBarry Smith x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; 17943447b6efSHong Zhang ib = idx + ii[0]; 17953447b6efSHong Zhang n = ii[1] - ii[0]; ii++; 17962d61bbb3SSatish Balay for (j=0; j<n; j++) { 17972d61bbb3SSatish Balay rval = ib[j]*4; 17982d61bbb3SSatish Balay z[rval++] += v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4; 17992d61bbb3SSatish Balay z[rval++] += v[4]*x1 + v[5]*x2 + v[6]*x3 + v[7]*x4; 18002d61bbb3SSatish Balay z[rval++] += v[8]*x1 + v[9]*x2 + v[10]*x3 + v[11]*x4; 18012d61bbb3SSatish Balay z[rval++] += v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4; 18022d61bbb3SSatish Balay v += 16; 18032d61bbb3SSatish Balay } 18043447b6efSHong Zhang if (!usecprow) xb += 4; 18052d61bbb3SSatish Balay } 18062d61bbb3SSatish Balay break; 18072d61bbb3SSatish Balay case 5: 18082d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 18097b2bb3b9SHong Zhang if (usecprow) xb = x + 5*ridx[i]; 1810f1af5d2fSBarry Smith x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; 18112d61bbb3SSatish Balay x4 = xb[3]; x5 = xb[4]; 18123447b6efSHong Zhang ib = idx + ii[0]; 18133447b6efSHong Zhang n = ii[1] - ii[0]; ii++; 18142d61bbb3SSatish Balay for (j=0; j<n; j++) { 18152d61bbb3SSatish Balay rval = ib[j]*5; 18162d61bbb3SSatish Balay z[rval++] += v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4 + v[4]*x5; 18172d61bbb3SSatish Balay z[rval++] += v[5]*x1 + v[6]*x2 + v[7]*x3 + v[8]*x4 + v[9]*x5; 18182d61bbb3SSatish Balay z[rval++] += v[10]*x1 + v[11]*x2 + v[12]*x3 + v[13]*x4 + v[14]*x5; 18192d61bbb3SSatish Balay z[rval++] += v[15]*x1 + v[16]*x2 + v[17]*x3 + v[18]*x4 + v[19]*x5; 18202d61bbb3SSatish Balay z[rval++] += v[20]*x1 + v[21]*x2 + v[22]*x3 + v[23]*x4 + v[24]*x5; 18212d61bbb3SSatish Balay v += 25; 18222d61bbb3SSatish Balay } 18233447b6efSHong Zhang if (!usecprow) xb += 5; 18242d61bbb3SSatish Balay } 18252d61bbb3SSatish Balay break; 1826f1af5d2fSBarry Smith default: { /* block sizes larger then 5 by 5 are handled by BLAS */ 1827690b6cddSBarry Smith PetscInt ncols,k; 18283447b6efSHong Zhang PetscScalar *work,*workt,*xtmp; 18293f1db9ecSBarry Smith 18302d61bbb3SSatish Balay if (!a->mult_work) { 1831d0f46423SBarry Smith k = PetscMax(A->rmap->n,A->cmap->n); 183287828ca2SBarry Smith ierr = PetscMalloc((k+1)*sizeof(PetscScalar),&a->mult_work);CHKERRQ(ierr); 18332d61bbb3SSatish Balay } 18342d61bbb3SSatish Balay work = a->mult_work; 18353447b6efSHong Zhang xtmp = x; 18362d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 18372d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 18382d61bbb3SSatish Balay ncols = n*bs; 183987828ca2SBarry Smith ierr = PetscMemzero(work,ncols*sizeof(PetscScalar));CHKERRQ(ierr); 18403447b6efSHong Zhang if (usecprow) { 18417b2bb3b9SHong Zhang xtmp = x + bs*ridx[i]; 18423447b6efSHong Zhang } 18433447b6efSHong Zhang Kernel_w_gets_w_plus_trans_Ar_times_v(bs,ncols,xtmp,v,work); 184471044d3cSBarry Smith /* BLASgemv_("T",&bs,&ncols,&_DOne,v,&bs,xtmp,&_One,&_DOne,work,&_One); */ 18452d61bbb3SSatish Balay v += n*bs2; 18463447b6efSHong Zhang if (!usecprow) xtmp += bs; 18472d61bbb3SSatish Balay workt = work; 18482d61bbb3SSatish Balay for (j=0; j<n; j++) { 18492d61bbb3SSatish Balay zb = z + bs*(*idx++); 18502d61bbb3SSatish Balay for (k=0; k<bs; k++) zb[k] += workt[k] ; 18512d61bbb3SSatish Balay workt += bs; 18522d61bbb3SSatish Balay } 18532d61bbb3SSatish Balay } 18542d61bbb3SSatish Balay } 18552d61bbb3SSatish Balay } 18563447b6efSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 18573447b6efSHong Zhang ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr); 1858dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*a->nz*a->bs2);CHKERRQ(ierr); 18592d61bbb3SSatish Balay PetscFunctionReturn(0); 18602d61bbb3SSatish Balay } 18612d61bbb3SSatish Balay 18624a2ae208SSatish Balay #undef __FUNCT__ 18634a2ae208SSatish Balay #define __FUNCT__ "MatScale_SeqBAIJ" 1864f4df32b1SMatthew Knepley PetscErrorCode MatScale_SeqBAIJ(Mat inA,PetscScalar alpha) 18652d61bbb3SSatish Balay { 18662d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)inA->data; 1867690b6cddSBarry Smith PetscInt totalnz = a->bs2*a->nz; 1868f4df32b1SMatthew Knepley PetscScalar oalpha = alpha; 1869efee365bSSatish Balay PetscErrorCode ierr; 18700805154bSBarry Smith PetscBLASInt one = 1,tnz = PetscBLASIntCast(totalnz); 18712d61bbb3SSatish Balay 18722d61bbb3SSatish Balay PetscFunctionBegin; 1873f4df32b1SMatthew Knepley BLASscal_(&tnz,&oalpha,a->a,&one); 1874efee365bSSatish Balay ierr = PetscLogFlops(totalnz);CHKERRQ(ierr); 18752d61bbb3SSatish Balay PetscFunctionReturn(0); 18762d61bbb3SSatish Balay } 18772d61bbb3SSatish Balay 18784a2ae208SSatish Balay #undef __FUNCT__ 18794a2ae208SSatish Balay #define __FUNCT__ "MatNorm_SeqBAIJ" 1880dfbe8321SBarry Smith PetscErrorCode MatNorm_SeqBAIJ(Mat A,NormType type,PetscReal *norm) 18812d61bbb3SSatish Balay { 18828a62d963SHong Zhang PetscErrorCode ierr; 18832d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 18843f1db9ecSBarry Smith MatScalar *v = a->a; 1885329f5518SBarry Smith PetscReal sum = 0.0; 1886d0f46423SBarry Smith PetscInt i,j,k,bs=A->rmap->bs,nz=a->nz,bs2=a->bs2,k1; 18872d61bbb3SSatish Balay 18882d61bbb3SSatish Balay PetscFunctionBegin; 18892d61bbb3SSatish Balay if (type == NORM_FROBENIUS) { 18902d61bbb3SSatish Balay for (i=0; i< bs2*nz; i++) { 1891aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 1892329f5518SBarry Smith sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 18932d61bbb3SSatish Balay #else 18942d61bbb3SSatish Balay sum += (*v)*(*v); v++; 18952d61bbb3SSatish Balay #endif 18962d61bbb3SSatish Balay } 18972d61bbb3SSatish Balay *norm = sqrt(sum); 18988a62d963SHong Zhang } else if (type == NORM_1) { /* maximum column sum */ 18998a62d963SHong Zhang PetscReal *tmp; 19008a62d963SHong Zhang PetscInt *bcol = a->j; 1901d0f46423SBarry Smith ierr = PetscMalloc((A->cmap->n+1)*sizeof(PetscReal),&tmp);CHKERRQ(ierr); 1902d0f46423SBarry Smith ierr = PetscMemzero(tmp,A->cmap->n*sizeof(PetscReal));CHKERRQ(ierr); 19038a62d963SHong Zhang for (i=0; i<nz; i++){ 19048a62d963SHong Zhang for (j=0; j<bs; j++){ 19058a62d963SHong Zhang k1 = bs*(*bcol) + j; /* column index */ 19068a62d963SHong Zhang for (k=0; k<bs; k++){ 19078a62d963SHong Zhang tmp[k1] += PetscAbsScalar(*v); v++; 19088a62d963SHong Zhang } 19098a62d963SHong Zhang } 19108a62d963SHong Zhang bcol++; 19118a62d963SHong Zhang } 19128a62d963SHong Zhang *norm = 0.0; 1913d0f46423SBarry Smith for (j=0; j<A->cmap->n; j++) { 19148a62d963SHong Zhang if (tmp[j] > *norm) *norm = tmp[j]; 19158a62d963SHong Zhang } 19168a62d963SHong Zhang ierr = PetscFree(tmp);CHKERRQ(ierr); 1917596552b5SBarry Smith } else if (type == NORM_INFINITY) { /* maximum row sum */ 1918596552b5SBarry Smith *norm = 0.0; 1919596552b5SBarry Smith for (k=0; k<bs; k++) { 192074f84c7bSSatish Balay for (j=0; j<a->mbs; j++) { 1921596552b5SBarry Smith v = a->a + bs2*a->i[j] + k; 1922596552b5SBarry Smith sum = 0.0; 1923596552b5SBarry Smith for (i=0; i<a->i[j+1]-a->i[j]; i++) { 19240e90e235SBarry Smith for (k1=0; k1<bs; k1++){ 1925596552b5SBarry Smith sum += PetscAbsScalar(*v); 1926596552b5SBarry Smith v += bs; 19272d61bbb3SSatish Balay } 19280e90e235SBarry Smith } 1929596552b5SBarry Smith if (sum > *norm) *norm = sum; 1930596552b5SBarry Smith } 1931596552b5SBarry Smith } 1932596552b5SBarry Smith } else { 193329bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"No support for this norm yet"); 19342d61bbb3SSatish Balay } 19352d61bbb3SSatish Balay PetscFunctionReturn(0); 19362d61bbb3SSatish Balay } 19372d61bbb3SSatish Balay 19382d61bbb3SSatish Balay 19394a2ae208SSatish Balay #undef __FUNCT__ 19404a2ae208SSatish Balay #define __FUNCT__ "MatEqual_SeqBAIJ" 1941dfbe8321SBarry Smith PetscErrorCode MatEqual_SeqBAIJ(Mat A,Mat B,PetscTruth* flg) 19422d61bbb3SSatish Balay { 19432d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data,*b = (Mat_SeqBAIJ *)B->data; 1944dfbe8321SBarry Smith PetscErrorCode ierr; 19452d61bbb3SSatish Balay 19462d61bbb3SSatish Balay PetscFunctionBegin; 19472d61bbb3SSatish Balay /* If the matrix/block dimensions are not equal, or no of nonzeros or shift */ 1948d0f46423SBarry Smith if ((A->rmap->N != B->rmap->N) || (A->cmap->n != B->cmap->n) || (A->rmap->bs != B->rmap->bs)|| (a->nz != b->nz)) { 1949273d9f13SBarry Smith *flg = PETSC_FALSE; 1950273d9f13SBarry Smith PetscFunctionReturn(0); 19512d61bbb3SSatish Balay } 19522d61bbb3SSatish Balay 19532d61bbb3SSatish Balay /* if the a->i are the same */ 1954690b6cddSBarry Smith ierr = PetscMemcmp(a->i,b->i,(a->mbs+1)*sizeof(PetscInt),flg);CHKERRQ(ierr); 1955abc0a331SBarry Smith if (!*flg) { 19560f5bd95cSBarry Smith PetscFunctionReturn(0); 19572d61bbb3SSatish Balay } 19582d61bbb3SSatish Balay 19592d61bbb3SSatish Balay /* if a->j are the same */ 1960690b6cddSBarry Smith ierr = PetscMemcmp(a->j,b->j,(a->nz)*sizeof(PetscInt),flg);CHKERRQ(ierr); 1961abc0a331SBarry Smith if (!*flg) { 19620f5bd95cSBarry Smith PetscFunctionReturn(0); 19632d61bbb3SSatish Balay } 19642d61bbb3SSatish Balay /* if a->a are the same */ 1965d0f46423SBarry Smith ierr = PetscMemcmp(a->a,b->a,(a->nz)*(A->rmap->bs)*(B->rmap->bs)*sizeof(PetscScalar),flg);CHKERRQ(ierr); 19662d61bbb3SSatish Balay PetscFunctionReturn(0); 19672d61bbb3SSatish Balay 19682d61bbb3SSatish Balay } 19692d61bbb3SSatish Balay 19704a2ae208SSatish Balay #undef __FUNCT__ 19714a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_SeqBAIJ" 1972dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_SeqBAIJ(Mat A,Vec v) 19732d61bbb3SSatish Balay { 19742d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1975dfbe8321SBarry Smith PetscErrorCode ierr; 1976690b6cddSBarry Smith PetscInt i,j,k,n,row,bs,*ai,*aj,ambs,bs2; 197787828ca2SBarry Smith PetscScalar *x,zero = 0.0; 19783f1db9ecSBarry Smith MatScalar *aa,*aa_j; 19792d61bbb3SSatish Balay 19802d61bbb3SSatish Balay PetscFunctionBegin; 198129bbc08cSBarry Smith if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 1982d0f46423SBarry Smith bs = A->rmap->bs; 19832d61bbb3SSatish Balay aa = a->a; 19842d61bbb3SSatish Balay ai = a->i; 19852d61bbb3SSatish Balay aj = a->j; 19862d61bbb3SSatish Balay ambs = a->mbs; 19872d61bbb3SSatish Balay bs2 = a->bs2; 19882d61bbb3SSatish Balay 19892dcb1b2aSMatthew Knepley ierr = VecSet(v,zero);CHKERRQ(ierr); 19901ebc52fbSHong Zhang ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1991e1311b90SBarry Smith ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr); 1992d0f46423SBarry Smith if (n != A->rmap->N) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector"); 19932d61bbb3SSatish Balay for (i=0; i<ambs; i++) { 19942d61bbb3SSatish Balay for (j=ai[i]; j<ai[i+1]; j++) { 19952d61bbb3SSatish Balay if (aj[j] == i) { 19962d61bbb3SSatish Balay row = i*bs; 19972d61bbb3SSatish Balay aa_j = aa+j*bs2; 19982d61bbb3SSatish Balay for (k=0; k<bs2; k+=(bs+1),row++) x[row] = aa_j[k]; 19992d61bbb3SSatish Balay break; 20002d61bbb3SSatish Balay } 20012d61bbb3SSatish Balay } 20022d61bbb3SSatish Balay } 20031ebc52fbSHong Zhang ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 20042d61bbb3SSatish Balay PetscFunctionReturn(0); 20052d61bbb3SSatish Balay } 20062d61bbb3SSatish Balay 20074a2ae208SSatish Balay #undef __FUNCT__ 20084a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_SeqBAIJ" 2009dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_SeqBAIJ(Mat A,Vec ll,Vec rr) 20102d61bbb3SSatish Balay { 20112d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 201287828ca2SBarry Smith PetscScalar *l,*r,x,*li,*ri; 20133f1db9ecSBarry Smith MatScalar *aa,*v; 2014dfbe8321SBarry Smith PetscErrorCode ierr; 2015690b6cddSBarry Smith PetscInt i,j,k,lm,rn,M,m,n,*ai,*aj,mbs,tmp,bs,bs2; 20162d61bbb3SSatish Balay 20172d61bbb3SSatish Balay PetscFunctionBegin; 20182d61bbb3SSatish Balay ai = a->i; 20192d61bbb3SSatish Balay aj = a->j; 20202d61bbb3SSatish Balay aa = a->a; 2021d0f46423SBarry Smith m = A->rmap->n; 2022d0f46423SBarry Smith n = A->cmap->n; 2023d0f46423SBarry Smith bs = A->rmap->bs; 20242d61bbb3SSatish Balay mbs = a->mbs; 20252d61bbb3SSatish Balay bs2 = a->bs2; 20262d61bbb3SSatish Balay if (ll) { 20271ebc52fbSHong Zhang ierr = VecGetArray(ll,&l);CHKERRQ(ierr); 2028e1311b90SBarry Smith ierr = VecGetLocalSize(ll,&lm);CHKERRQ(ierr); 202929bbc08cSBarry Smith if (lm != m) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vector wrong length"); 20302d61bbb3SSatish Balay for (i=0; i<mbs; i++) { /* for each block row */ 20312d61bbb3SSatish Balay M = ai[i+1] - ai[i]; 20322d61bbb3SSatish Balay li = l + i*bs; 20332d61bbb3SSatish Balay v = aa + bs2*ai[i]; 20342d61bbb3SSatish Balay for (j=0; j<M; j++) { /* for each block */ 20352d61bbb3SSatish Balay for (k=0; k<bs2; k++) { 20362d61bbb3SSatish Balay (*v++) *= li[k%bs]; 20372d61bbb3SSatish Balay } 20382d61bbb3SSatish Balay } 20392d61bbb3SSatish Balay } 20401ebc52fbSHong Zhang ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr); 2041efee365bSSatish Balay ierr = PetscLogFlops(a->nz);CHKERRQ(ierr); 20422d61bbb3SSatish Balay } 20432d61bbb3SSatish Balay 20442d61bbb3SSatish Balay if (rr) { 20451ebc52fbSHong Zhang ierr = VecGetArray(rr,&r);CHKERRQ(ierr); 2046e1311b90SBarry Smith ierr = VecGetLocalSize(rr,&rn);CHKERRQ(ierr); 204729bbc08cSBarry Smith if (rn != n) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vector wrong length"); 20482d61bbb3SSatish Balay for (i=0; i<mbs; i++) { /* for each block row */ 20492d61bbb3SSatish Balay M = ai[i+1] - ai[i]; 20502d61bbb3SSatish Balay v = aa + bs2*ai[i]; 20512d61bbb3SSatish Balay for (j=0; j<M; j++) { /* for each block */ 20522d61bbb3SSatish Balay ri = r + bs*aj[ai[i]+j]; 20532d61bbb3SSatish Balay for (k=0; k<bs; k++) { 20542d61bbb3SSatish Balay x = ri[k]; 20552d61bbb3SSatish Balay for (tmp=0; tmp<bs; tmp++) (*v++) *= x; 20562d61bbb3SSatish Balay } 20572d61bbb3SSatish Balay } 20582d61bbb3SSatish Balay } 20591ebc52fbSHong Zhang ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr); 2060efee365bSSatish Balay ierr = PetscLogFlops(a->nz);CHKERRQ(ierr); 20612d61bbb3SSatish Balay } 20622d61bbb3SSatish Balay PetscFunctionReturn(0); 20632d61bbb3SSatish Balay } 20642d61bbb3SSatish Balay 20652d61bbb3SSatish Balay 20664a2ae208SSatish Balay #undef __FUNCT__ 20674a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_SeqBAIJ" 2068dfbe8321SBarry Smith PetscErrorCode MatGetInfo_SeqBAIJ(Mat A,MatInfoType flag,MatInfo *info) 20692d61bbb3SSatish Balay { 20702d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 20712d61bbb3SSatish Balay 20722d61bbb3SSatish Balay PetscFunctionBegin; 20732d61bbb3SSatish Balay info->block_size = a->bs2; 20742d61bbb3SSatish Balay info->nz_allocated = a->maxnz; 20752d61bbb3SSatish Balay info->nz_used = a->bs2*a->nz; 20762d61bbb3SSatish Balay info->nz_unneeded = (double)(info->nz_allocated - info->nz_used); 20772d61bbb3SSatish Balay info->assemblies = A->num_ass; 20782d61bbb3SSatish Balay info->mallocs = a->reallocs; 20797adad957SLisandro Dalcin info->memory = ((PetscObject)A)->mem; 20802d61bbb3SSatish Balay if (A->factor) { 20812d61bbb3SSatish Balay info->fill_ratio_given = A->info.fill_ratio_given; 20822d61bbb3SSatish Balay info->fill_ratio_needed = A->info.fill_ratio_needed; 20832d61bbb3SSatish Balay info->factor_mallocs = A->info.factor_mallocs; 20842d61bbb3SSatish Balay } else { 20852d61bbb3SSatish Balay info->fill_ratio_given = 0; 20862d61bbb3SSatish Balay info->fill_ratio_needed = 0; 20872d61bbb3SSatish Balay info->factor_mallocs = 0; 20882d61bbb3SSatish Balay } 20892d61bbb3SSatish Balay PetscFunctionReturn(0); 20902d61bbb3SSatish Balay } 20912d61bbb3SSatish Balay 20922d61bbb3SSatish Balay 20934a2ae208SSatish Balay #undef __FUNCT__ 20944a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_SeqBAIJ" 2095dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_SeqBAIJ(Mat A) 20962d61bbb3SSatish Balay { 20972d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 2098dfbe8321SBarry Smith PetscErrorCode ierr; 20992d61bbb3SSatish Balay 21002d61bbb3SSatish Balay PetscFunctionBegin; 2101549d3d68SSatish Balay ierr = PetscMemzero(a->a,a->bs2*a->i[a->mbs]*sizeof(MatScalar));CHKERRQ(ierr); 21022d61bbb3SSatish Balay PetscFunctionReturn(0); 21032d61bbb3SSatish Balay } 2104