173f4d377SMatthew Knepley /*$Id: baij2.c,v 1.75 2001/09/07 20:09:49 bsmith Exp $*/ 2cac129eeSSatish Balay 370f55243SBarry Smith #include "src/mat/impls/baij/seq/baij.h" 42d61bbb3SSatish Balay #include "src/vec/vecimpl.h" 52d61bbb3SSatish Balay #include "src/inline/spops.h" 63f1db9ecSBarry Smith #include "src/inline/ilu.h" 70a835dfdSSatish Balay #include "petscbt.h" 8cac129eeSSatish Balay 94a2ae208SSatish Balay #undef __FUNCT__ 104a2ae208SSatish Balay #define __FUNCT__ "MatIncreaseOverlap_SeqBAIJ" 11736121d4SSatish Balay int MatIncreaseOverlap_SeqBAIJ(Mat A,int is_max,IS *is,int ov) 12a3192f15SSatish Balay { 13a3192f15SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 14218c64b6SSatish Balay int row,i,j,k,l,m,n,*idx,ierr,*nidx,isz,val,ival; 15218c64b6SSatish Balay int 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; 22218c64b6SSatish Balay bs = a->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); 27b0a32e0cSBarry Smith ierr = PetscMalloc((m+1)*sizeof(int),&nidx);CHKERRQ(ierr); 28b0a32e0cSBarry Smith ierr = PetscMalloc((A->m+1)*sizeof(int),&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 */ 4229bbc08cSBarry Smith 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" 767b2a1423SBarry Smith int MatGetSubMatrix_SeqBAIJ_Private(Mat A,IS isrow,IS iscol,int cs,MatReuse scall,Mat *B) 77736121d4SSatish Balay { 78736121d4SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*c; 7940011551SBarry Smith int *smap,i,k,kstart,kend,ierr,oldcols = a->nbs,*lens; 80218c64b6SSatish Balay int row,mat_i,*mat_j,tcol,*mat_ilen; 81736121d4SSatish Balay int *irow,*icol,nrows,ncols,*ssmap,bs=a->bs,bs2=a->bs2; 82218c64b6SSatish Balay int *aj = a->j,*ai = a->i; 833f1db9ecSBarry Smith MatScalar *mat_a; 84736121d4SSatish Balay Mat C; 850f5bd95cSBarry Smith PetscTruth flag; 86736121d4SSatish Balay 873a40ed3dSBarry Smith PetscFunctionBegin; 88888f2ed8SSatish Balay ierr = ISSorted(iscol,(PetscTruth*)&i);CHKERRQ(ierr); 8929bbc08cSBarry Smith if (!i) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"IS is not sorted"); 90736121d4SSatish Balay 91736121d4SSatish Balay ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr); 92218c64b6SSatish Balay ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr); 93b9b97703SBarry Smith ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr); 94b9b97703SBarry Smith ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr); 95736121d4SSatish Balay 96b0a32e0cSBarry Smith ierr = PetscMalloc((1+oldcols)*sizeof(int),&smap);CHKERRQ(ierr); 97736121d4SSatish Balay ssmap = smap; 98b0a32e0cSBarry Smith ierr = PetscMalloc((1+nrows)*sizeof(int),&lens);CHKERRQ(ierr); 99549d3d68SSatish Balay ierr = PetscMemzero(smap,oldcols*sizeof(int));CHKERRQ(ierr); 100736121d4SSatish Balay for (i=0; i<ncols; i++) smap[icol[i]] = i+1; 101736121d4SSatish Balay /* determine lens of each row */ 102736121d4SSatish Balay for (i=0; i<nrows; i++) { 103736121d4SSatish Balay kstart = ai[irow[i]]; 104736121d4SSatish Balay kend = kstart + a->ilen[irow[i]]; 105736121d4SSatish Balay lens[i] = 0; 106736121d4SSatish Balay for (k=kstart; k<kend; k++) { 107736121d4SSatish Balay if (ssmap[aj[k]]) { 108736121d4SSatish Balay lens[i]++; 109736121d4SSatish Balay } 110736121d4SSatish Balay } 111736121d4SSatish Balay } 112736121d4SSatish Balay /* Create and fill new matrix */ 113736121d4SSatish Balay if (scall == MAT_REUSE_MATRIX) { 114736121d4SSatish Balay c = (Mat_SeqBAIJ *)((*B)->data); 115736121d4SSatish Balay 11629bbc08cSBarry Smith if (c->mbs!=nrows || c->nbs!=ncols || c->bs!=bs) SETERRQ(PETSC_ERR_ARG_SIZ,"Submatrix wrong size"); 1170f5bd95cSBarry Smith ierr = PetscMemcmp(c->ilen,lens,c->mbs *sizeof(int),&flag);CHKERRQ(ierr); 1180f5bd95cSBarry Smith if (flag == PETSC_FALSE) { 11929bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong no of nonzeros"); 120736121d4SSatish Balay } 121549d3d68SSatish Balay ierr = PetscMemzero(c->ilen,c->mbs*sizeof(int));CHKERRQ(ierr); 122736121d4SSatish Balay C = *B; 1233a40ed3dSBarry Smith } else { 124736121d4SSatish Balay ierr = MatCreateSeqBAIJ(A->comm,bs,nrows*bs,ncols*bs,0,lens,&C);CHKERRQ(ierr); 125736121d4SSatish Balay } 126736121d4SSatish Balay c = (Mat_SeqBAIJ *)(C->data); 127736121d4SSatish Balay for (i=0; i<nrows; i++) { 128736121d4SSatish Balay row = irow[i]; 129736121d4SSatish Balay kstart = ai[row]; 130736121d4SSatish Balay kend = kstart + a->ilen[row]; 131736121d4SSatish Balay mat_i = c->i[i]; 132736121d4SSatish Balay mat_j = c->j + mat_i; 133218c64b6SSatish Balay mat_a = c->a + mat_i*bs2; 134736121d4SSatish Balay mat_ilen = c->ilen + i; 135736121d4SSatish Balay for (k=kstart; k<kend; k++) { 136736121d4SSatish Balay if ((tcol=ssmap[a->j[k]])) { 137736121d4SSatish Balay *mat_j++ = tcol - 1; 138549d3d68SSatish Balay ierr = PetscMemcpy(mat_a,a->a+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 139549d3d68SSatish Balay mat_a += bs2; 140736121d4SSatish Balay (*mat_ilen)++; 141736121d4SSatish Balay } 142736121d4SSatish Balay } 143736121d4SSatish Balay } 144218c64b6SSatish Balay 145736121d4SSatish Balay /* Free work space */ 146736121d4SSatish Balay ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr); 147606d414cSSatish Balay ierr = PetscFree(smap);CHKERRQ(ierr); 148606d414cSSatish Balay ierr = PetscFree(lens);CHKERRQ(ierr); 1496d4a8577SBarry Smith ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1506d4a8577SBarry Smith ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 151736121d4SSatish Balay 152736121d4SSatish Balay ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr); 153736121d4SSatish Balay *B = C; 1543a40ed3dSBarry Smith PetscFunctionReturn(0); 155736121d4SSatish Balay } 156736121d4SSatish Balay 1574a2ae208SSatish Balay #undef __FUNCT__ 1584a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqBAIJ" 1597b2a1423SBarry Smith int MatGetSubMatrix_SeqBAIJ(Mat A,IS isrow,IS iscol,int cs,MatReuse scall,Mat *B) 160218c64b6SSatish Balay { 161218c64b6SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 162218c64b6SSatish Balay IS is1,is2; 163218c64b6SSatish Balay int *vary,*iary,*irow,*icol,nrows,ncols,i,ierr,bs=a->bs,count; 164218c64b6SSatish Balay 1653a40ed3dSBarry Smith PetscFunctionBegin; 166218c64b6SSatish Balay ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr); 167218c64b6SSatish Balay ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr); 168b9b97703SBarry Smith ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr); 169b9b97703SBarry Smith ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr); 170218c64b6SSatish Balay 171218c64b6SSatish Balay /* Verify if the indices corespond to each element in a block 172218c64b6SSatish Balay and form the IS with compressed IS */ 173b0a32e0cSBarry Smith ierr = PetscMalloc(2*(a->mbs+1)*sizeof(int),&vary);CHKERRQ(ierr); 174218c64b6SSatish Balay iary = vary + a->mbs; 175549d3d68SSatish Balay ierr = PetscMemzero(vary,(a->mbs)*sizeof(int));CHKERRQ(ierr); 176218c64b6SSatish Balay for (i=0; i<nrows; i++) vary[irow[i]/bs]++; 177218c64b6SSatish Balay count = 0; 178218c64b6SSatish Balay for (i=0; i<a->mbs; i++) { 179ac355199SBarry Smith if (vary[i]!=0 && vary[i]!=bs) SETERRQ(1,"Index set does not match blocks"); 180218c64b6SSatish Balay if (vary[i]==bs) iary[count++] = i; 181218c64b6SSatish Balay } 182029af93fSBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,count,iary,&is1);CHKERRQ(ierr); 183218c64b6SSatish Balay 184549d3d68SSatish Balay ierr = PetscMemzero(vary,(a->mbs)*sizeof(int));CHKERRQ(ierr); 185218c64b6SSatish Balay for (i=0; i<ncols; i++) vary[icol[i]/bs]++; 186218c64b6SSatish Balay count = 0; 187218c64b6SSatish Balay for (i=0; i<a->mbs; i++) { 188ac355199SBarry Smith if (vary[i]!=0 && vary[i]!=bs) SETERRQ(1,"Internal error in PETSc"); 189218c64b6SSatish Balay if (vary[i]==bs) iary[count++] = i; 190218c64b6SSatish Balay } 191029af93fSBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,count,iary,&is2);CHKERRQ(ierr); 192218c64b6SSatish Balay ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr); 193218c64b6SSatish Balay ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr); 194606d414cSSatish Balay ierr = PetscFree(vary);CHKERRQ(ierr); 195218c64b6SSatish Balay 1966a6a5d1dSBarry Smith ierr = MatGetSubMatrix_SeqBAIJ_Private(A,is1,is2,cs,scall,B);CHKERRQ(ierr); 197218c64b6SSatish Balay ISDestroy(is1); 198218c64b6SSatish Balay ISDestroy(is2); 1993a40ed3dSBarry Smith PetscFunctionReturn(0); 200218c64b6SSatish Balay } 201218c64b6SSatish Balay 2024a2ae208SSatish Balay #undef __FUNCT__ 2034a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_SeqBAIJ" 2047b2a1423SBarry Smith int MatGetSubMatrices_SeqBAIJ(Mat A,int n,IS *irow,IS *icol,MatReuse scall,Mat **B) 205736121d4SSatish Balay { 206736121d4SSatish Balay int ierr,i; 207736121d4SSatish Balay 2083a40ed3dSBarry Smith PetscFunctionBegin; 209736121d4SSatish Balay if (scall == MAT_INITIAL_MATRIX) { 21082502324SSatish Balay ierr = PetscMalloc((n+1)*sizeof(Mat),B);CHKERRQ(ierr); 211736121d4SSatish Balay } 212736121d4SSatish Balay 213736121d4SSatish Balay for (i=0; i<n; i++) { 2146a6a5d1dSBarry Smith ierr = MatGetSubMatrix_SeqBAIJ(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr); 215736121d4SSatish Balay } 2163a40ed3dSBarry Smith PetscFunctionReturn(0); 217736121d4SSatish Balay } 218218c64b6SSatish Balay 219218c64b6SSatish Balay 2202d61bbb3SSatish Balay /* -------------------------------------------------------*/ 2212d61bbb3SSatish Balay /* Should check that shapes of vectors and matrices match */ 2222d61bbb3SSatish Balay /* -------------------------------------------------------*/ 223d9eff348SSatish Balay #include "petscblaslapack.h" 2242d61bbb3SSatish Balay 2254a2ae208SSatish Balay #undef __FUNCT__ 2264a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_1" 2272d61bbb3SSatish Balay int MatMult_SeqBAIJ_1(Mat A,Vec xx,Vec zz) 2282d61bbb3SSatish Balay { 2292d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 23087828ca2SBarry Smith PetscScalar *x,*z,sum; 2313f1db9ecSBarry Smith MatScalar *v; 232e1311b90SBarry Smith int mbs=a->mbs,i,*idx,*ii,n,ierr; 2332d61bbb3SSatish Balay 2342d61bbb3SSatish Balay PetscFunctionBegin; 235e1311b90SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 236e1311b90SBarry Smith ierr = VecGetArray(zz,&z);CHKERRQ(ierr); 2372d61bbb3SSatish Balay 2382d61bbb3SSatish Balay idx = a->j; 2392d61bbb3SSatish Balay v = a->a; 2402d61bbb3SSatish Balay ii = a->i; 2412d61bbb3SSatish Balay 2422d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 2432d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 2442d61bbb3SSatish Balay sum = 0.0; 2452d61bbb3SSatish Balay while (n--) sum += *v++ * x[*idx++]; 2462d61bbb3SSatish Balay z[i] = sum; 2472d61bbb3SSatish Balay } 248e1311b90SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 249e1311b90SBarry Smith ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr); 250b0a32e0cSBarry Smith PetscLogFlops(2*a->nz - A->m); 2512d61bbb3SSatish Balay PetscFunctionReturn(0); 2522d61bbb3SSatish Balay } 2532d61bbb3SSatish Balay 2544a2ae208SSatish Balay #undef __FUNCT__ 2554a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_2" 2562d61bbb3SSatish Balay int MatMult_SeqBAIJ_2(Mat A,Vec xx,Vec zz) 2572d61bbb3SSatish Balay { 2582d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 25987828ca2SBarry Smith PetscScalar *x,*z,*xb,sum1,sum2; 26087828ca2SBarry Smith PetscScalar x1,x2; 2613f1db9ecSBarry Smith MatScalar *v; 262e1311b90SBarry Smith int ierr,mbs=a->mbs,i,*idx,*ii,j,n; 2632d61bbb3SSatish Balay 2642d61bbb3SSatish Balay PetscFunctionBegin; 265e1311b90SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 266e1311b90SBarry Smith ierr = VecGetArray(zz,&z);CHKERRQ(ierr); 2672d61bbb3SSatish Balay 2682d61bbb3SSatish Balay idx = a->j; 2692d61bbb3SSatish Balay v = a->a; 2702d61bbb3SSatish Balay ii = a->i; 2712d61bbb3SSatish Balay 2722d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 2732d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 2742d61bbb3SSatish Balay sum1 = 0.0; sum2 = 0.0; 2752d61bbb3SSatish Balay for (j=0; j<n; j++) { 2762d61bbb3SSatish Balay xb = x + 2*(*idx++); x1 = xb[0]; x2 = xb[1]; 2772d61bbb3SSatish Balay sum1 += v[0]*x1 + v[2]*x2; 2782d61bbb3SSatish Balay sum2 += v[1]*x1 + v[3]*x2; 2792d61bbb3SSatish Balay v += 4; 2802d61bbb3SSatish Balay } 2812d61bbb3SSatish Balay z[0] = sum1; z[1] = sum2; 2822d61bbb3SSatish Balay z += 2; 2832d61bbb3SSatish Balay } 284e1311b90SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 285e1311b90SBarry Smith ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr); 286b0a32e0cSBarry Smith PetscLogFlops(8*a->nz - A->m); 2872d61bbb3SSatish Balay PetscFunctionReturn(0); 2882d61bbb3SSatish Balay } 2892d61bbb3SSatish Balay 2904a2ae208SSatish Balay #undef __FUNCT__ 2914a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_3" 2922d61bbb3SSatish Balay int MatMult_SeqBAIJ_3(Mat A,Vec xx,Vec zz) 2932d61bbb3SSatish Balay { 2942d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 29587828ca2SBarry Smith PetscScalar *x,*z,*xb,sum1,sum2,sum3,x1,x2,x3; 2963f1db9ecSBarry Smith MatScalar *v; 297e1311b90SBarry Smith int ierr,mbs=a->mbs,i,*idx,*ii,j,n; 2982d61bbb3SSatish Balay 299b6410449SSatish Balay #if defined(PETSC_HAVE_PRAGMA_DISJOINT) 300fee21e36SBarry Smith #pragma disjoint(*v,*z,*xb) 301fee21e36SBarry Smith #endif 302fee21e36SBarry Smith 3032d61bbb3SSatish Balay PetscFunctionBegin; 304e1311b90SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 305e1311b90SBarry Smith ierr = VecGetArray(zz,&z);CHKERRQ(ierr); 3062d61bbb3SSatish Balay 3072d61bbb3SSatish Balay idx = a->j; 3082d61bbb3SSatish Balay v = a->a; 3092d61bbb3SSatish Balay ii = a->i; 3102d61bbb3SSatish Balay 3112d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 3122d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 3132d61bbb3SSatish Balay sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; 3142d61bbb3SSatish Balay for (j=0; j<n; j++) { 3152d61bbb3SSatish Balay xb = x + 3*(*idx++); x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; 3162d61bbb3SSatish Balay sum1 += v[0]*x1 + v[3]*x2 + v[6]*x3; 3172d61bbb3SSatish Balay sum2 += v[1]*x1 + v[4]*x2 + v[7]*x3; 3182d61bbb3SSatish Balay sum3 += v[2]*x1 + v[5]*x2 + v[8]*x3; 3192d61bbb3SSatish Balay v += 9; 3202d61bbb3SSatish Balay } 3212d61bbb3SSatish Balay z[0] = sum1; z[1] = sum2; z[2] = sum3; 3222d61bbb3SSatish Balay z += 3; 3232d61bbb3SSatish Balay } 324e1311b90SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 325e1311b90SBarry Smith ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr); 326b0a32e0cSBarry Smith PetscLogFlops(18*a->nz - A->m); 3272d61bbb3SSatish Balay PetscFunctionReturn(0); 3282d61bbb3SSatish Balay } 3292d61bbb3SSatish Balay 3304a2ae208SSatish Balay #undef __FUNCT__ 3314a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_4" 3322d61bbb3SSatish Balay int MatMult_SeqBAIJ_4(Mat A,Vec xx,Vec zz) 3332d61bbb3SSatish Balay { 3342d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 33587828ca2SBarry Smith PetscScalar *x,*z,*xb,sum1,sum2,sum3,sum4,x1,x2,x3,x4; 3363f1db9ecSBarry Smith MatScalar *v; 337e1311b90SBarry Smith int ierr,mbs=a->mbs,i,*idx,*ii,j,n; 3382d61bbb3SSatish Balay 3392d61bbb3SSatish Balay PetscFunctionBegin; 340e1311b90SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 341e1311b90SBarry Smith ierr = VecGetArray(zz,&z);CHKERRQ(ierr); 3422d61bbb3SSatish Balay 3432d61bbb3SSatish Balay idx = a->j; 3442d61bbb3SSatish Balay v = a->a; 3452d61bbb3SSatish Balay ii = a->i; 3462d61bbb3SSatish Balay 3472d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 3482d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 3492d61bbb3SSatish Balay sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; 3502d61bbb3SSatish Balay for (j=0; j<n; j++) { 3512d61bbb3SSatish Balay xb = x + 4*(*idx++); 3522d61bbb3SSatish Balay x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; 3532d61bbb3SSatish Balay sum1 += v[0]*x1 + v[4]*x2 + v[8]*x3 + v[12]*x4; 3542d61bbb3SSatish Balay sum2 += v[1]*x1 + v[5]*x2 + v[9]*x3 + v[13]*x4; 3552d61bbb3SSatish Balay sum3 += v[2]*x1 + v[6]*x2 + v[10]*x3 + v[14]*x4; 3562d61bbb3SSatish Balay sum4 += v[3]*x1 + v[7]*x2 + v[11]*x3 + v[15]*x4; 3572d61bbb3SSatish Balay v += 16; 3582d61bbb3SSatish Balay } 3592d61bbb3SSatish Balay z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; 3602d61bbb3SSatish Balay z += 4; 3612d61bbb3SSatish Balay } 362e1311b90SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 363e1311b90SBarry Smith ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr); 364b0a32e0cSBarry Smith PetscLogFlops(32*a->nz - A->m); 3652d61bbb3SSatish Balay PetscFunctionReturn(0); 3662d61bbb3SSatish Balay } 3672d61bbb3SSatish Balay 3684a2ae208SSatish Balay #undef __FUNCT__ 3694a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_5" 3702d61bbb3SSatish Balay int MatMult_SeqBAIJ_5(Mat A,Vec xx,Vec zz) 3712d61bbb3SSatish Balay { 3722d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 37387828ca2SBarry Smith PetscScalar sum1,sum2,sum3,sum4,sum5,x1,x2,x3,x4,x5,*xb,*z,*x; 3743f1db9ecSBarry Smith MatScalar *v; 375e1311b90SBarry Smith int ierr,mbs=a->mbs,i,*idx,*ii,j,n; 3762d61bbb3SSatish Balay 377433994e6SBarry Smith PetscFunctionBegin; 378e1311b90SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 379e1311b90SBarry Smith ierr = VecGetArray(zz,&z);CHKERRQ(ierr); 3802d61bbb3SSatish Balay 3812d61bbb3SSatish Balay idx = a->j; 3822d61bbb3SSatish Balay v = a->a; 3832d61bbb3SSatish Balay ii = a->i; 3842d61bbb3SSatish Balay 3852d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 3862d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 3872d61bbb3SSatish Balay sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; 3882d61bbb3SSatish Balay for (j=0; j<n; j++) { 3892d61bbb3SSatish Balay xb = x + 5*(*idx++); 3902d61bbb3SSatish Balay x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; 3912d61bbb3SSatish Balay sum1 += v[0]*x1 + v[5]*x2 + v[10]*x3 + v[15]*x4 + v[20]*x5; 3922d61bbb3SSatish Balay sum2 += v[1]*x1 + v[6]*x2 + v[11]*x3 + v[16]*x4 + v[21]*x5; 3932d61bbb3SSatish Balay sum3 += v[2]*x1 + v[7]*x2 + v[12]*x3 + v[17]*x4 + v[22]*x5; 3942d61bbb3SSatish Balay sum4 += v[3]*x1 + v[8]*x2 + v[13]*x3 + v[18]*x4 + v[23]*x5; 3952d61bbb3SSatish Balay sum5 += v[4]*x1 + v[9]*x2 + v[14]*x3 + v[19]*x4 + v[24]*x5; 3962d61bbb3SSatish Balay v += 25; 3972d61bbb3SSatish Balay } 3982d61bbb3SSatish Balay z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; 3992d61bbb3SSatish Balay z += 5; 4002d61bbb3SSatish Balay } 401e1311b90SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 402e1311b90SBarry Smith ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr); 403b0a32e0cSBarry Smith PetscLogFlops(50*a->nz - A->m); 4042d61bbb3SSatish Balay PetscFunctionReturn(0); 4052d61bbb3SSatish Balay } 4062d61bbb3SSatish Balay 40715091d37SBarry Smith 4084a2ae208SSatish Balay #undef __FUNCT__ 4094a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_6" 41015091d37SBarry Smith int MatMult_SeqBAIJ_6(Mat A,Vec xx,Vec zz) 41115091d37SBarry Smith { 41215091d37SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 41387828ca2SBarry Smith PetscScalar *x,*z,*xb,sum1,sum2,sum3,sum4,sum5,sum6; 41487828ca2SBarry Smith PetscScalar x1,x2,x3,x4,x5,x6; 41515091d37SBarry Smith MatScalar *v; 41615091d37SBarry Smith int ierr,mbs=a->mbs,i,*idx,*ii,j,n; 41715091d37SBarry Smith 418433994e6SBarry Smith PetscFunctionBegin; 41915091d37SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 42015091d37SBarry Smith ierr = VecGetArray(zz,&z);CHKERRQ(ierr); 42115091d37SBarry Smith 42215091d37SBarry Smith idx = a->j; 42315091d37SBarry Smith v = a->a; 42415091d37SBarry Smith ii = a->i; 42515091d37SBarry Smith 42615091d37SBarry Smith for (i=0; i<mbs; i++) { 42715091d37SBarry Smith n = ii[1] - ii[0]; ii++; 42815091d37SBarry Smith sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0; 42915091d37SBarry Smith for (j=0; j<n; j++) { 43015091d37SBarry Smith xb = x + 6*(*idx++); 43115091d37SBarry Smith x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; 43215091d37SBarry Smith sum1 += v[0]*x1 + v[6]*x2 + v[12]*x3 + v[18]*x4 + v[24]*x5 + v[30]*x6; 43315091d37SBarry Smith sum2 += v[1]*x1 + v[7]*x2 + v[13]*x3 + v[19]*x4 + v[25]*x5 + v[31]*x6; 43415091d37SBarry Smith sum3 += v[2]*x1 + v[8]*x2 + v[14]*x3 + v[20]*x4 + v[26]*x5 + v[32]*x6; 43515091d37SBarry Smith sum4 += v[3]*x1 + v[9]*x2 + v[15]*x3 + v[21]*x4 + v[27]*x5 + v[33]*x6; 43615091d37SBarry Smith sum5 += v[4]*x1 + v[10]*x2 + v[16]*x3 + v[22]*x4 + v[28]*x5 + v[34]*x6; 43715091d37SBarry Smith sum6 += v[5]*x1 + v[11]*x2 + v[17]*x3 + v[23]*x4 + v[29]*x5 + v[35]*x6; 43815091d37SBarry Smith v += 36; 43915091d37SBarry Smith } 44015091d37SBarry Smith z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; 44115091d37SBarry Smith z += 6; 44215091d37SBarry Smith } 44315091d37SBarry Smith 44415091d37SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 44515091d37SBarry Smith ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr); 446b0a32e0cSBarry Smith PetscLogFlops(72*a->nz - A->m); 44715091d37SBarry Smith PetscFunctionReturn(0); 44815091d37SBarry Smith } 4494a2ae208SSatish Balay #undef __FUNCT__ 4504a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_7" 4512d61bbb3SSatish Balay int MatMult_SeqBAIJ_7(Mat A,Vec xx,Vec zz) 4522d61bbb3SSatish Balay { 4532d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 45487828ca2SBarry Smith PetscScalar *x,*z,*xb,sum1,sum2,sum3,sum4,sum5,sum6,sum7; 45587828ca2SBarry Smith PetscScalar x1,x2,x3,x4,x5,x6,x7; 4563f1db9ecSBarry Smith MatScalar *v; 457e1311b90SBarry Smith int ierr,mbs=a->mbs,i,*idx,*ii,j,n; 4582d61bbb3SSatish Balay 459433994e6SBarry Smith PetscFunctionBegin; 460e1311b90SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 461e1311b90SBarry Smith ierr = VecGetArray(zz,&z);CHKERRQ(ierr); 4622d61bbb3SSatish Balay 4632d61bbb3SSatish Balay idx = a->j; 4642d61bbb3SSatish Balay v = a->a; 4652d61bbb3SSatish Balay ii = a->i; 4662d61bbb3SSatish Balay 4672d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 4682d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 4692d61bbb3SSatish Balay sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0; sum7 = 0.0; 4702d61bbb3SSatish Balay for (j=0; j<n; j++) { 4712d61bbb3SSatish Balay xb = x + 7*(*idx++); 4722d61bbb3SSatish Balay x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6]; 4732d61bbb3SSatish Balay sum1 += v[0]*x1 + v[7]*x2 + v[14]*x3 + v[21]*x4 + v[28]*x5 + v[35]*x6 + v[42]*x7; 4742d61bbb3SSatish Balay sum2 += v[1]*x1 + v[8]*x2 + v[15]*x3 + v[22]*x4 + v[29]*x5 + v[36]*x6 + v[43]*x7; 4752d61bbb3SSatish Balay sum3 += v[2]*x1 + v[9]*x2 + v[16]*x3 + v[23]*x4 + v[30]*x5 + v[37]*x6 + v[44]*x7; 4762d61bbb3SSatish Balay sum4 += v[3]*x1 + v[10]*x2 + v[17]*x3 + v[24]*x4 + v[31]*x5 + v[38]*x6 + v[45]*x7; 4772d61bbb3SSatish Balay sum5 += v[4]*x1 + v[11]*x2 + v[18]*x3 + v[25]*x4 + v[32]*x5 + v[39]*x6 + v[46]*x7; 4782d61bbb3SSatish Balay sum6 += v[5]*x1 + v[12]*x2 + v[19]*x3 + v[26]*x4 + v[33]*x5 + v[40]*x6 + v[47]*x7; 4792d61bbb3SSatish Balay sum7 += v[6]*x1 + v[13]*x2 + v[20]*x3 + v[27]*x4 + v[34]*x5 + v[41]*x6 + v[48]*x7; 4802d61bbb3SSatish Balay v += 49; 4812d61bbb3SSatish Balay } 4822d61bbb3SSatish Balay z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7; 4832d61bbb3SSatish Balay z += 7; 4842d61bbb3SSatish Balay } 4852d61bbb3SSatish Balay 486e1311b90SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 487e1311b90SBarry Smith ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr); 488b0a32e0cSBarry Smith PetscLogFlops(98*a->nz - A->m); 4892d61bbb3SSatish Balay PetscFunctionReturn(0); 4902d61bbb3SSatish Balay } 4912d61bbb3SSatish Balay 4923f1db9ecSBarry Smith /* 4933f1db9ecSBarry Smith This will not work with MatScalar == float because it calls the BLAS 4943f1db9ecSBarry Smith */ 4954a2ae208SSatish Balay #undef __FUNCT__ 4964a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_N" 4972d61bbb3SSatish Balay int MatMult_SeqBAIJ_N(Mat A,Vec xx,Vec zz) 4982d61bbb3SSatish Balay { 4992d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 50087828ca2SBarry Smith PetscScalar *x,*z,*xb,*work,*workt; 5013f1db9ecSBarry Smith MatScalar *v; 502e1311b90SBarry Smith int ierr,mbs=a->mbs,i,*idx,*ii,bs=a->bs,j,n,bs2=a->bs2; 5033f1db9ecSBarry Smith int ncols,k; 5042d61bbb3SSatish Balay 5052d61bbb3SSatish Balay PetscFunctionBegin; 506e1311b90SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 507e1311b90SBarry Smith ierr = VecGetArray(zz,&z);CHKERRQ(ierr); 5082d61bbb3SSatish Balay 5092d61bbb3SSatish Balay idx = a->j; 5102d61bbb3SSatish Balay v = a->a; 5112d61bbb3SSatish Balay ii = a->i; 512218c64b6SSatish Balay 513218c64b6SSatish Balay 5142d61bbb3SSatish Balay if (!a->mult_work) { 515273d9f13SBarry Smith k = PetscMax(A->m,A->n); 51687828ca2SBarry Smith ierr = PetscMalloc((k+1)*sizeof(PetscScalar),&a->mult_work);CHKERRQ(ierr); 5172d61bbb3SSatish Balay } 5182d61bbb3SSatish Balay work = a->mult_work; 5192d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 5202d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 5212d61bbb3SSatish Balay ncols = n*bs; 5222d61bbb3SSatish Balay workt = work; 5232d61bbb3SSatish Balay for (j=0; j<n; j++) { 5242d61bbb3SSatish Balay xb = x + bs*(*idx++); 5252d61bbb3SSatish Balay for (k=0; k<bs; k++) workt[k] = xb[k]; 5262d61bbb3SSatish Balay workt += bs; 5272d61bbb3SSatish Balay } 5283f1db9ecSBarry Smith Kernel_w_gets_Ar_times_v(bs,ncols,work,v,z); 5293f1db9ecSBarry Smith /* LAgemv_("N",&bs,&ncols,&_DOne,v,&bs,work,&_One,&_DZero,z,&_One); */ 5302d61bbb3SSatish Balay v += n*bs2; 5312d61bbb3SSatish Balay z += bs; 5322d61bbb3SSatish Balay } 533e1311b90SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 534e1311b90SBarry Smith ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr); 535b0a32e0cSBarry Smith PetscLogFlops(2*a->nz*bs2 - A->m); 5362d61bbb3SSatish Balay PetscFunctionReturn(0); 5372d61bbb3SSatish Balay } 5382d61bbb3SSatish Balay 5394a2ae208SSatish Balay #undef __FUNCT__ 5404a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_1" 5412d61bbb3SSatish Balay int MatMultAdd_SeqBAIJ_1(Mat A,Vec xx,Vec yy,Vec zz) 5422d61bbb3SSatish Balay { 5432d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 54487828ca2SBarry Smith PetscScalar *x,*y,*z,sum; 5453f1db9ecSBarry Smith MatScalar *v; 546e1311b90SBarry Smith int ierr,mbs=a->mbs,i,*idx,*ii,n; 5472d61bbb3SSatish Balay 5482d61bbb3SSatish Balay PetscFunctionBegin; 549e1311b90SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 550e1311b90SBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 5512e8a6d31SBarry Smith if (zz != yy) { 552e1311b90SBarry Smith ierr = VecGetArray(zz,&z);CHKERRQ(ierr); 5532e8a6d31SBarry Smith } else { 5542e8a6d31SBarry Smith z = y; 5552e8a6d31SBarry Smith } 5562d61bbb3SSatish Balay 5572d61bbb3SSatish Balay idx = a->j; 5582d61bbb3SSatish Balay v = a->a; 5592d61bbb3SSatish Balay ii = a->i; 5602d61bbb3SSatish Balay 5612d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 5622d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 5632d61bbb3SSatish Balay sum = y[i]; 5642d61bbb3SSatish Balay while (n--) sum += *v++ * x[*idx++]; 5652d61bbb3SSatish Balay z[i] = sum; 5662d61bbb3SSatish Balay } 567e1311b90SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 568e1311b90SBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 5692e8a6d31SBarry Smith if (zz != yy) { 570e1311b90SBarry Smith ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr); 5712e8a6d31SBarry Smith } 572b0a32e0cSBarry Smith PetscLogFlops(2*a->nz); 5732d61bbb3SSatish Balay PetscFunctionReturn(0); 5742d61bbb3SSatish Balay } 5752d61bbb3SSatish Balay 5764a2ae208SSatish Balay #undef __FUNCT__ 5774a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_2" 5782d61bbb3SSatish Balay int MatMultAdd_SeqBAIJ_2(Mat A,Vec xx,Vec yy,Vec zz) 5792d61bbb3SSatish Balay { 5802d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 58187828ca2SBarry Smith PetscScalar *x,*y,*z,*xb,sum1,sum2; 58287828ca2SBarry Smith PetscScalar x1,x2; 5833f1db9ecSBarry Smith MatScalar *v; 584e1311b90SBarry Smith int ierr,mbs=a->mbs,i,*idx,*ii,j,n; 5852d61bbb3SSatish Balay 5862d61bbb3SSatish Balay PetscFunctionBegin; 587e1311b90SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 588e1311b90SBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 5892e8a6d31SBarry Smith if (zz != yy) { 590e1311b90SBarry Smith ierr = VecGetArray(zz,&z);CHKERRQ(ierr); 5912e8a6d31SBarry Smith } else { 5922e8a6d31SBarry Smith z = y; 5932e8a6d31SBarry Smith } 5942d61bbb3SSatish Balay 5952d61bbb3SSatish Balay idx = a->j; 5962d61bbb3SSatish Balay v = a->a; 5972d61bbb3SSatish Balay ii = a->i; 5982d61bbb3SSatish Balay 5992d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 6002d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 6012d61bbb3SSatish Balay sum1 = y[0]; sum2 = y[1]; 6022d61bbb3SSatish Balay for (j=0; j<n; j++) { 6032d61bbb3SSatish Balay xb = x + 2*(*idx++); x1 = xb[0]; x2 = xb[1]; 6042d61bbb3SSatish Balay sum1 += v[0]*x1 + v[2]*x2; 6052d61bbb3SSatish Balay sum2 += v[1]*x1 + v[3]*x2; 6062d61bbb3SSatish Balay v += 4; 6072d61bbb3SSatish Balay } 6082d61bbb3SSatish Balay z[0] = sum1; z[1] = sum2; 6092d61bbb3SSatish Balay z += 2; y += 2; 6102d61bbb3SSatish Balay } 611e1311b90SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 612e1311b90SBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 6132e8a6d31SBarry Smith if (zz != yy) { 614e1311b90SBarry Smith ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr); 6152e8a6d31SBarry Smith } 616b0a32e0cSBarry Smith PetscLogFlops(4*a->nz); 6172d61bbb3SSatish Balay PetscFunctionReturn(0); 6182d61bbb3SSatish Balay } 6192d61bbb3SSatish Balay 6204a2ae208SSatish Balay #undef __FUNCT__ 6214a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_3" 6222d61bbb3SSatish Balay int MatMultAdd_SeqBAIJ_3(Mat A,Vec xx,Vec yy,Vec zz) 6232d61bbb3SSatish Balay { 6242d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 62587828ca2SBarry Smith PetscScalar *x,*y,*z,*xb,sum1,sum2,sum3,x1,x2,x3; 6263f1db9ecSBarry Smith MatScalar *v; 627e1311b90SBarry Smith int ierr,mbs=a->mbs,i,*idx,*ii,j,n; 6282d61bbb3SSatish Balay 6292d61bbb3SSatish Balay PetscFunctionBegin; 630e1311b90SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 631e1311b90SBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 6322e8a6d31SBarry Smith if (zz != yy) { 633e1311b90SBarry Smith ierr = VecGetArray(zz,&z);CHKERRQ(ierr); 6342e8a6d31SBarry Smith } else { 6352e8a6d31SBarry Smith z = y; 6362e8a6d31SBarry Smith } 6372d61bbb3SSatish Balay 6382d61bbb3SSatish Balay idx = a->j; 6392d61bbb3SSatish Balay v = a->a; 6402d61bbb3SSatish Balay ii = a->i; 6412d61bbb3SSatish Balay 6422d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 6432d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 6442d61bbb3SSatish Balay sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; 6452d61bbb3SSatish Balay for (j=0; j<n; j++) { 6462d61bbb3SSatish Balay xb = x + 3*(*idx++); x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; 6472d61bbb3SSatish Balay sum1 += v[0]*x1 + v[3]*x2 + v[6]*x3; 6482d61bbb3SSatish Balay sum2 += v[1]*x1 + v[4]*x2 + v[7]*x3; 6492d61bbb3SSatish Balay sum3 += v[2]*x1 + v[5]*x2 + v[8]*x3; 6502d61bbb3SSatish Balay v += 9; 6512d61bbb3SSatish Balay } 6522d61bbb3SSatish Balay z[0] = sum1; z[1] = sum2; z[2] = sum3; 6532d61bbb3SSatish Balay z += 3; y += 3; 6542d61bbb3SSatish Balay } 655e1311b90SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 656e1311b90SBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 6572e8a6d31SBarry Smith if (zz != yy) { 658e1311b90SBarry Smith ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr); 6592e8a6d31SBarry Smith } 660b0a32e0cSBarry Smith PetscLogFlops(18*a->nz); 6612d61bbb3SSatish Balay PetscFunctionReturn(0); 6622d61bbb3SSatish Balay } 6632d61bbb3SSatish Balay 6644a2ae208SSatish Balay #undef __FUNCT__ 6654a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_4" 6662d61bbb3SSatish Balay int MatMultAdd_SeqBAIJ_4(Mat A,Vec xx,Vec yy,Vec zz) 6672d61bbb3SSatish Balay { 6682d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 66987828ca2SBarry Smith PetscScalar *x,*y,*z,*xb,sum1,sum2,sum3,sum4,x1,x2,x3,x4; 6703f1db9ecSBarry Smith MatScalar *v; 671e1311b90SBarry Smith int ierr,mbs=a->mbs,i,*idx,*ii; 6722d61bbb3SSatish Balay int j,n; 6732d61bbb3SSatish Balay 6742d61bbb3SSatish Balay PetscFunctionBegin; 675e1311b90SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 676e1311b90SBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 6772e8a6d31SBarry Smith if (zz != yy) { 678e1311b90SBarry Smith ierr = VecGetArray(zz,&z);CHKERRQ(ierr); 6792e8a6d31SBarry Smith } else { 6802e8a6d31SBarry Smith z = y; 6812e8a6d31SBarry Smith } 6822d61bbb3SSatish Balay 6832d61bbb3SSatish Balay idx = a->j; 6842d61bbb3SSatish Balay v = a->a; 6852d61bbb3SSatish Balay ii = a->i; 6862d61bbb3SSatish Balay 6872d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 6882d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 6892d61bbb3SSatish Balay sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; 6902d61bbb3SSatish Balay for (j=0; j<n; j++) { 6912d61bbb3SSatish Balay xb = x + 4*(*idx++); 6922d61bbb3SSatish Balay x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; 6932d61bbb3SSatish Balay sum1 += v[0]*x1 + v[4]*x2 + v[8]*x3 + v[12]*x4; 6942d61bbb3SSatish Balay sum2 += v[1]*x1 + v[5]*x2 + v[9]*x3 + v[13]*x4; 6952d61bbb3SSatish Balay sum3 += v[2]*x1 + v[6]*x2 + v[10]*x3 + v[14]*x4; 6962d61bbb3SSatish Balay sum4 += v[3]*x1 + v[7]*x2 + v[11]*x3 + v[15]*x4; 6972d61bbb3SSatish Balay v += 16; 6982d61bbb3SSatish Balay } 6992d61bbb3SSatish Balay z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; 7002d61bbb3SSatish Balay z += 4; y += 4; 7012d61bbb3SSatish Balay } 702e1311b90SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 703e1311b90SBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 7042e8a6d31SBarry Smith if (zz != yy) { 705e1311b90SBarry Smith ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr); 7062e8a6d31SBarry Smith } 707b0a32e0cSBarry Smith PetscLogFlops(32*a->nz); 7082d61bbb3SSatish Balay PetscFunctionReturn(0); 7092d61bbb3SSatish Balay } 7102d61bbb3SSatish Balay 7114a2ae208SSatish Balay #undef __FUNCT__ 7124a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_5" 7132d61bbb3SSatish Balay int MatMultAdd_SeqBAIJ_5(Mat A,Vec xx,Vec yy,Vec zz) 7142d61bbb3SSatish Balay { 7152d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 71687828ca2SBarry Smith PetscScalar *x,*y,*z,*xb,sum1,sum2,sum3,sum4,sum5,x1,x2,x3,x4,x5; 7173f1db9ecSBarry Smith MatScalar *v; 718e1311b90SBarry Smith int ierr,mbs=a->mbs,i,*idx,*ii,j,n; 7192d61bbb3SSatish Balay 7202d61bbb3SSatish Balay PetscFunctionBegin; 721e1311b90SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 722e1311b90SBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 7232e8a6d31SBarry Smith if (zz != yy) { 724e1311b90SBarry Smith ierr = VecGetArray(zz,&z);CHKERRQ(ierr); 7252e8a6d31SBarry Smith } else { 7262e8a6d31SBarry Smith z = y; 7272e8a6d31SBarry Smith } 7282d61bbb3SSatish Balay 7292d61bbb3SSatish Balay idx = a->j; 7302d61bbb3SSatish Balay v = a->a; 7312d61bbb3SSatish Balay ii = a->i; 7322d61bbb3SSatish Balay 7332d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 7342d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 7352d61bbb3SSatish Balay sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4]; 7362d61bbb3SSatish Balay for (j=0; j<n; j++) { 7372d61bbb3SSatish Balay xb = x + 5*(*idx++); 7382d61bbb3SSatish Balay x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; 7392d61bbb3SSatish Balay sum1 += v[0]*x1 + v[5]*x2 + v[10]*x3 + v[15]*x4 + v[20]*x5; 7402d61bbb3SSatish Balay sum2 += v[1]*x1 + v[6]*x2 + v[11]*x3 + v[16]*x4 + v[21]*x5; 7412d61bbb3SSatish Balay sum3 += v[2]*x1 + v[7]*x2 + v[12]*x3 + v[17]*x4 + v[22]*x5; 7422d61bbb3SSatish Balay sum4 += v[3]*x1 + v[8]*x2 + v[13]*x3 + v[18]*x4 + v[23]*x5; 7432d61bbb3SSatish Balay sum5 += v[4]*x1 + v[9]*x2 + v[14]*x3 + v[19]*x4 + v[24]*x5; 7442d61bbb3SSatish Balay v += 25; 7452d61bbb3SSatish Balay } 7462d61bbb3SSatish Balay z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; 7472d61bbb3SSatish Balay z += 5; y += 5; 7482d61bbb3SSatish Balay } 749e1311b90SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 750e1311b90SBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 7512e8a6d31SBarry Smith if (zz != yy) { 752e1311b90SBarry Smith ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr); 7532e8a6d31SBarry Smith } 754b0a32e0cSBarry Smith PetscLogFlops(50*a->nz); 7552d61bbb3SSatish Balay PetscFunctionReturn(0); 7562d61bbb3SSatish Balay } 7574a2ae208SSatish Balay #undef __FUNCT__ 7584a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_6" 75915091d37SBarry Smith int MatMultAdd_SeqBAIJ_6(Mat A,Vec xx,Vec yy,Vec zz) 76015091d37SBarry Smith { 76115091d37SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 76287828ca2SBarry Smith PetscScalar *x,*y,*z,*xb,sum1,sum2,sum3,sum4,sum5,sum6; 76387828ca2SBarry Smith PetscScalar x1,x2,x3,x4,x5,x6; 76415091d37SBarry Smith MatScalar *v; 76515091d37SBarry Smith int ierr,mbs=a->mbs,i,*idx,*ii,j,n; 76615091d37SBarry Smith 76715091d37SBarry Smith PetscFunctionBegin; 76815091d37SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 76915091d37SBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 77015091d37SBarry Smith if (zz != yy) { 77115091d37SBarry Smith ierr = VecGetArray(zz,&z);CHKERRQ(ierr); 77215091d37SBarry Smith } else { 77315091d37SBarry Smith z = y; 77415091d37SBarry Smith } 77515091d37SBarry Smith 77615091d37SBarry Smith idx = a->j; 77715091d37SBarry Smith v = a->a; 77815091d37SBarry Smith ii = a->i; 77915091d37SBarry Smith 78015091d37SBarry Smith for (i=0; i<mbs; i++) { 78115091d37SBarry Smith n = ii[1] - ii[0]; ii++; 78215091d37SBarry Smith sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4]; sum6 = y[5]; 78315091d37SBarry Smith for (j=0; j<n; j++) { 7843b95cb0eSSatish Balay xb = x + 6*(*idx++); 78515091d37SBarry Smith x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; 78615091d37SBarry Smith sum1 += v[0]*x1 + v[6]*x2 + v[12]*x3 + v[18]*x4 + v[24]*x5 + v[30]*x6; 78715091d37SBarry Smith sum2 += v[1]*x1 + v[7]*x2 + v[13]*x3 + v[19]*x4 + v[25]*x5 + v[31]*x6; 78815091d37SBarry Smith sum3 += v[2]*x1 + v[8]*x2 + v[14]*x3 + v[20]*x4 + v[26]*x5 + v[32]*x6; 78915091d37SBarry Smith sum4 += v[3]*x1 + v[9]*x2 + v[15]*x3 + v[21]*x4 + v[27]*x5 + v[33]*x6; 79015091d37SBarry Smith sum5 += v[4]*x1 + v[10]*x2 + v[16]*x3 + v[22]*x4 + v[28]*x5 + v[34]*x6; 79115091d37SBarry Smith sum6 += v[5]*x1 + v[11]*x2 + v[17]*x3 + v[23]*x4 + v[29]*x5 + v[35]*x6; 79215091d37SBarry Smith v += 36; 79315091d37SBarry Smith } 79415091d37SBarry Smith z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; 79515091d37SBarry Smith z += 6; y += 6; 79615091d37SBarry Smith } 79715091d37SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 79815091d37SBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 79915091d37SBarry Smith if (zz != yy) { 80015091d37SBarry Smith ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr); 80115091d37SBarry Smith } 802b0a32e0cSBarry Smith PetscLogFlops(72*a->nz); 80315091d37SBarry Smith PetscFunctionReturn(0); 80415091d37SBarry Smith } 8052d61bbb3SSatish Balay 8064a2ae208SSatish Balay #undef __FUNCT__ 8074a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_7" 8082d61bbb3SSatish Balay int MatMultAdd_SeqBAIJ_7(Mat A,Vec xx,Vec yy,Vec zz) 8092d61bbb3SSatish Balay { 8102d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 81187828ca2SBarry Smith PetscScalar *x,*y,*z,*xb,sum1,sum2,sum3,sum4,sum5,sum6,sum7; 81287828ca2SBarry Smith PetscScalar x1,x2,x3,x4,x5,x6,x7; 8133f1db9ecSBarry Smith MatScalar *v; 814e1311b90SBarry Smith int ierr,mbs=a->mbs,i,*idx,*ii,j,n; 8152d61bbb3SSatish Balay 8162d61bbb3SSatish Balay PetscFunctionBegin; 817e1311b90SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 818e1311b90SBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 8192e8a6d31SBarry Smith if (zz != yy) { 820e1311b90SBarry Smith ierr = VecGetArray(zz,&z);CHKERRQ(ierr); 8212e8a6d31SBarry Smith } else { 8222e8a6d31SBarry Smith z = y; 8232e8a6d31SBarry Smith } 8242d61bbb3SSatish Balay 8252d61bbb3SSatish Balay idx = a->j; 8262d61bbb3SSatish Balay v = a->a; 8272d61bbb3SSatish Balay ii = a->i; 8282d61bbb3SSatish Balay 8292d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 8302d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 8312d61bbb3SSatish Balay sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4]; sum6 = y[5]; sum7 = y[6]; 8322d61bbb3SSatish Balay for (j=0; j<n; j++) { 8332d61bbb3SSatish Balay xb = x + 7*(*idx++); 8342d61bbb3SSatish Balay x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6]; 8352d61bbb3SSatish Balay sum1 += v[0]*x1 + v[7]*x2 + v[14]*x3 + v[21]*x4 + v[28]*x5 + v[35]*x6 + v[42]*x7; 8362d61bbb3SSatish Balay sum2 += v[1]*x1 + v[8]*x2 + v[15]*x3 + v[22]*x4 + v[29]*x5 + v[36]*x6 + v[43]*x7; 8372d61bbb3SSatish Balay sum3 += v[2]*x1 + v[9]*x2 + v[16]*x3 + v[23]*x4 + v[30]*x5 + v[37]*x6 + v[44]*x7; 8382d61bbb3SSatish Balay sum4 += v[3]*x1 + v[10]*x2 + v[17]*x3 + v[24]*x4 + v[31]*x5 + v[38]*x6 + v[45]*x7; 8392d61bbb3SSatish Balay sum5 += v[4]*x1 + v[11]*x2 + v[18]*x3 + v[25]*x4 + v[32]*x5 + v[39]*x6 + v[46]*x7; 8402d61bbb3SSatish Balay sum6 += v[5]*x1 + v[12]*x2 + v[19]*x3 + v[26]*x4 + v[33]*x5 + v[40]*x6 + v[47]*x7; 8412d61bbb3SSatish Balay sum7 += v[6]*x1 + v[13]*x2 + v[20]*x3 + v[27]*x4 + v[34]*x5 + v[41]*x6 + v[48]*x7; 8422d61bbb3SSatish Balay v += 49; 8432d61bbb3SSatish Balay } 8442d61bbb3SSatish Balay z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7; 8452d61bbb3SSatish Balay z += 7; y += 7; 8462d61bbb3SSatish Balay } 847e1311b90SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 848e1311b90SBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 8492e8a6d31SBarry Smith if (zz != yy) { 850e1311b90SBarry Smith ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr); 8512e8a6d31SBarry Smith } 852b0a32e0cSBarry Smith PetscLogFlops(98*a->nz); 8532d61bbb3SSatish Balay PetscFunctionReturn(0); 8542d61bbb3SSatish Balay } 855218c64b6SSatish Balay 8564a2ae208SSatish Balay #undef __FUNCT__ 8574a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_N" 8582d61bbb3SSatish Balay int MatMultAdd_SeqBAIJ_N(Mat A,Vec xx,Vec yy,Vec zz) 8592d61bbb3SSatish Balay { 8602d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 861*b2e7d493SHong Zhang PetscScalar *x,*z,*xb,*work,*workt,*y; 8623f1db9ecSBarry Smith MatScalar *v; 8632d61bbb3SSatish Balay int mbs=a->mbs,i,*idx,*ii,bs=a->bs,j,n,bs2=a->bs2,ierr; 8643f1db9ecSBarry Smith int ncols,k; 865218c64b6SSatish Balay 8662d61bbb3SSatish Balay PetscFunctionBegin; 867e1311b90SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 868e1311b90SBarry Smith ierr = VecGetArray(zz,&z);CHKERRQ(ierr); 869*b2e7d493SHong Zhang if (zz != yy) { 870*b2e7d493SHong Zhang ierr = VecGetArrayFast(yy,&y);CHKERRQ(ierr); 871*b2e7d493SHong Zhang ierr = PetscMemcpy(z,y,yy->n*sizeof(PetscScalar));CHKERRQ(ierr); 872*b2e7d493SHong Zhang ierr = VecRestoreArrayFast(yy,&y);CHKERRQ(ierr); 873*b2e7d493SHong Zhang } 8742d61bbb3SSatish Balay 8752d61bbb3SSatish Balay idx = a->j; 8762d61bbb3SSatish Balay v = a->a; 8772d61bbb3SSatish Balay ii = a->i; 8782d61bbb3SSatish Balay 8792d61bbb3SSatish Balay 8802d61bbb3SSatish Balay if (!a->mult_work) { 881273d9f13SBarry Smith k = PetscMax(A->m,A->n); 88287828ca2SBarry Smith ierr = PetscMalloc((k+1)*sizeof(PetscScalar),&a->mult_work);CHKERRQ(ierr); 8832d61bbb3SSatish Balay } 8842d61bbb3SSatish Balay work = a->mult_work; 8852d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 8862d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 8872d61bbb3SSatish Balay ncols = n*bs; 8882d61bbb3SSatish Balay workt = work; 8892d61bbb3SSatish Balay for (j=0; j<n; j++) { 8902d61bbb3SSatish Balay xb = x + bs*(*idx++); 8912d61bbb3SSatish Balay for (k=0; k<bs; k++) workt[k] = xb[k]; 8922d61bbb3SSatish Balay workt += bs; 8932d61bbb3SSatish Balay } 8943f1db9ecSBarry Smith Kernel_w_gets_w_plus_Ar_times_v(bs,ncols,work,v,z); 8953f1db9ecSBarry Smith /* LAgemv_("N",&bs,&ncols,&_DOne,v,&bs,work,&_One,&_DOne,z,&_One); */ 8962d61bbb3SSatish Balay v += n*bs2; 8972d61bbb3SSatish Balay z += bs; 8982d61bbb3SSatish Balay } 899e1311b90SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 900e1311b90SBarry Smith ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr); 901b0a32e0cSBarry Smith PetscLogFlops(2*a->nz*bs2); 9022d61bbb3SSatish Balay PetscFunctionReturn(0); 9032d61bbb3SSatish Balay } 9042d61bbb3SSatish Balay 9054a2ae208SSatish Balay #undef __FUNCT__ 9064a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_SeqBAIJ" 9077c922b88SBarry Smith int MatMultTranspose_SeqBAIJ(Mat A,Vec xx,Vec zz) 9082d61bbb3SSatish Balay { 9092d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 91087828ca2SBarry Smith PetscScalar *xg,*zg,*zb,zero = 0.0; 91187828ca2SBarry Smith PetscScalar *x,*z,*xb,x1,x2,x3,x4,x5,x6,x7; 9123f1db9ecSBarry Smith MatScalar *v; 913f1af5d2fSBarry Smith int mbs=a->mbs,i,*idx,*ii,*ai=a->i,rval; 9142d61bbb3SSatish Balay int bs=a->bs,j,n,bs2=a->bs2,*ib,ierr; 9152d61bbb3SSatish Balay 9162d61bbb3SSatish Balay PetscFunctionBegin; 917f1af5d2fSBarry Smith ierr = VecSet(&zero,zz);CHKERRQ(ierr); 918e1311b90SBarry Smith ierr = VecGetArray(xx,&xg);CHKERRQ(ierr); x = xg; 919e1311b90SBarry Smith ierr = VecGetArray(zz,&zg);CHKERRQ(ierr); z = zg; 9202d61bbb3SSatish Balay 9212d61bbb3SSatish Balay idx = a->j; 9222d61bbb3SSatish Balay v = a->a; 9232d61bbb3SSatish Balay ii = a->i; 924f1af5d2fSBarry Smith xb = x; 9252d61bbb3SSatish Balay switch (bs) { 9262d61bbb3SSatish Balay case 1: 9272d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 9282d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 929f1af5d2fSBarry Smith x1 = xb[0]; 9302d61bbb3SSatish Balay ib = idx + ai[i]; 9312d61bbb3SSatish Balay for (j=0; j<n; j++) { 9322d61bbb3SSatish Balay rval = ib[j]; 933f1af5d2fSBarry Smith z[rval] += *v * x1; 934f1af5d2fSBarry Smith v++; 9352d61bbb3SSatish Balay } 936f1af5d2fSBarry Smith xb++; 9372d61bbb3SSatish Balay } 9382d61bbb3SSatish Balay break; 9392d61bbb3SSatish Balay case 2: 9402d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 9412d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 942f1af5d2fSBarry Smith x1 = xb[0]; x2 = xb[1]; 9432d61bbb3SSatish Balay ib = idx + ai[i]; 9442d61bbb3SSatish Balay for (j=0; j<n; j++) { 9452d61bbb3SSatish Balay rval = ib[j]*2; 9462d61bbb3SSatish Balay z[rval++] += v[0]*x1 + v[1]*x2; 947f1af5d2fSBarry Smith z[rval] += v[2]*x1 + v[3]*x2; 9482d61bbb3SSatish Balay v += 4; 9492d61bbb3SSatish Balay } 950f1af5d2fSBarry Smith xb += 2; 9512d61bbb3SSatish Balay } 9522d61bbb3SSatish Balay break; 9532d61bbb3SSatish Balay case 3: 9542d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 9552d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 956f1af5d2fSBarry Smith x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; 9572d61bbb3SSatish Balay ib = idx + ai[i]; 9582d61bbb3SSatish Balay for (j=0; j<n; j++) { 9592d61bbb3SSatish Balay rval = ib[j]*3; 9602d61bbb3SSatish Balay z[rval++] += v[0]*x1 + v[1]*x2 + v[2]*x3; 9612d61bbb3SSatish Balay z[rval++] += v[3]*x1 + v[4]*x2 + v[5]*x3; 962f1af5d2fSBarry Smith z[rval] += v[6]*x1 + v[7]*x2 + v[8]*x3; 9632d61bbb3SSatish Balay v += 9; 9642d61bbb3SSatish Balay } 965f1af5d2fSBarry Smith xb += 3; 9662d61bbb3SSatish Balay } 9672d61bbb3SSatish Balay break; 9682d61bbb3SSatish Balay case 4: 9692d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 9702d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 971f1af5d2fSBarry Smith x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; 9722d61bbb3SSatish Balay ib = idx + ai[i]; 9732d61bbb3SSatish Balay for (j=0; j<n; j++) { 9742d61bbb3SSatish Balay rval = ib[j]*4; 9752d61bbb3SSatish Balay z[rval++] += v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4; 9762d61bbb3SSatish Balay z[rval++] += v[4]*x1 + v[5]*x2 + v[6]*x3 + v[7]*x4; 9772d61bbb3SSatish Balay z[rval++] += v[8]*x1 + v[9]*x2 + v[10]*x3 + v[11]*x4; 978f1af5d2fSBarry Smith z[rval] += v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4; 9792d61bbb3SSatish Balay v += 16; 9802d61bbb3SSatish Balay } 981f1af5d2fSBarry Smith xb += 4; 9822d61bbb3SSatish Balay } 9832d61bbb3SSatish Balay break; 9842d61bbb3SSatish Balay case 5: 9852d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 9862d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 987f1af5d2fSBarry Smith x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; 9882d61bbb3SSatish Balay x4 = xb[3]; x5 = xb[4]; 9892d61bbb3SSatish Balay ib = idx + ai[i]; 9902d61bbb3SSatish Balay for (j=0; j<n; j++) { 9912d61bbb3SSatish Balay rval = ib[j]*5; 9922d61bbb3SSatish Balay z[rval++] += v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4 + v[4]*x5; 9932d61bbb3SSatish Balay z[rval++] += v[5]*x1 + v[6]*x2 + v[7]*x3 + v[8]*x4 + v[9]*x5; 9942d61bbb3SSatish Balay z[rval++] += v[10]*x1 + v[11]*x2 + v[12]*x3 + v[13]*x4 + v[14]*x5; 9952d61bbb3SSatish Balay z[rval++] += v[15]*x1 + v[16]*x2 + v[17]*x3 + v[18]*x4 + v[19]*x5; 996f1af5d2fSBarry Smith z[rval] += v[20]*x1 + v[21]*x2 + v[22]*x3 + v[23]*x4 + v[24]*x5; 9972d61bbb3SSatish Balay v += 25; 9982d61bbb3SSatish Balay } 999f1af5d2fSBarry Smith xb += 5; 10002d61bbb3SSatish Balay } 10012d61bbb3SSatish Balay break; 1002f1af5d2fSBarry Smith case 6: 1003f1af5d2fSBarry Smith for (i=0; i<mbs; i++) { 1004f1af5d2fSBarry Smith n = ii[1] - ii[0]; ii++; 1005f1af5d2fSBarry Smith x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; 1006f1af5d2fSBarry Smith x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; 1007f1af5d2fSBarry Smith ib = idx + ai[i]; 1008f1af5d2fSBarry Smith for (j=0; j<n; j++) { 1009f1af5d2fSBarry Smith rval = ib[j]*6; 1010f1af5d2fSBarry Smith z[rval++] += v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4 + v[4]*x5 + v[5]*x6; 1011f1af5d2fSBarry Smith z[rval++] += v[6]*x1 + v[7]*x2 + v[8]*x3 + v[9]*x4 + v[10]*x5 + v[11]*x6; 1012f1af5d2fSBarry Smith z[rval++] += v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4 + v[16]*x5 + v[17]*x6; 1013f1af5d2fSBarry Smith z[rval++] += v[18]*x1 + v[19]*x2 + v[20]*x3 + v[21]*x4 + v[22]*x5 + v[23]*x6; 1014f1af5d2fSBarry Smith z[rval++] += v[24]*x1 + v[25]*x2 + v[26]*x3 + v[27]*x4 + v[28]*x5 + v[29]*x6; 1015f1af5d2fSBarry Smith z[rval] += v[30]*x1 + v[31]*x2 + v[32]*x3 + v[33]*x4 + v[34]*x5 + v[35]*x6; 1016f1af5d2fSBarry Smith v += 36; 1017f1af5d2fSBarry Smith } 1018f1af5d2fSBarry Smith xb += 6; 1019f1af5d2fSBarry Smith } 1020f1af5d2fSBarry Smith break; 1021f1af5d2fSBarry Smith case 7: 1022f1af5d2fSBarry Smith for (i=0; i<mbs; i++) { 1023f1af5d2fSBarry Smith n = ii[1] - ii[0]; ii++; 1024f1af5d2fSBarry Smith x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; 1025f1af5d2fSBarry Smith x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6]; 1026f1af5d2fSBarry Smith ib = idx + ai[i]; 1027f1af5d2fSBarry Smith for (j=0; j<n; j++) { 1028f1af5d2fSBarry Smith rval = ib[j]*7; 1029f1af5d2fSBarry Smith z[rval++] += v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4 + v[4]*x5 + v[5]*x6 + v[6]*x7; 1030f1af5d2fSBarry Smith z[rval++] += v[7]*x1 + v[8]*x2 + v[9]*x3 + v[10]*x4 + v[11]*x5 + v[12]*x6 + v[13]*x7; 1031f1af5d2fSBarry Smith z[rval++] += v[14]*x1 + v[15]*x2 + v[16]*x3 + v[17]*x4 + v[18]*x5 + v[19]*x6 + v[20]*x7; 1032f1af5d2fSBarry Smith z[rval++] += v[21]*x1 + v[22]*x2 + v[23]*x3 + v[24]*x4 + v[25]*x5 + v[26]*x6 + v[27]*x7; 1033f1af5d2fSBarry Smith z[rval++] += v[28]*x1 + v[29]*x2 + v[30]*x3 + v[31]*x4 + v[32]*x5 + v[33]*x6 + v[34]*x7; 1034f1af5d2fSBarry Smith z[rval++] += v[35]*x1 + v[36]*x2 + v[37]*x3 + v[38]*x4 + v[39]*x5 + v[40]*x6 + v[41]*x7; 1035f1af5d2fSBarry Smith z[rval] += v[42]*x1 + v[43]*x2 + v[44]*x3 + v[45]*x4 + v[46]*x5 + v[47]*x6 + v[48]*x7; 1036f1af5d2fSBarry Smith v += 49; 1037f1af5d2fSBarry Smith } 1038f1af5d2fSBarry Smith xb += 7; 1039f1af5d2fSBarry Smith } 1040f1af5d2fSBarry Smith break; 1041f1af5d2fSBarry Smith default: { /* block sizes larger then 7 by 7 are handled by BLAS */ 10423f1db9ecSBarry Smith int ncols,k; 104387828ca2SBarry Smith PetscScalar *work,*workt; 10443f1db9ecSBarry Smith 10452d61bbb3SSatish Balay if (!a->mult_work) { 1046273d9f13SBarry Smith k = PetscMax(A->m,A->n); 104787828ca2SBarry Smith ierr = PetscMalloc((k+1)*sizeof(PetscScalar),&a->mult_work);CHKERRQ(ierr); 10482d61bbb3SSatish Balay } 10492d61bbb3SSatish Balay work = a->mult_work; 10502d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 10512d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 10522d61bbb3SSatish Balay ncols = n*bs; 105387828ca2SBarry Smith ierr = PetscMemzero(work,ncols*sizeof(PetscScalar));CHKERRQ(ierr); 1054d824769bSBarry Smith Kernel_w_gets_w_plus_trans_Ar_times_v(bs,ncols,x,v,work); 10553f1db9ecSBarry Smith /* LAgemv_("T",&bs,&ncols,&_DOne,v,&bs,x,&_One,&_DOne,work,&_One); */ 10562d61bbb3SSatish Balay v += n*bs2; 10572d61bbb3SSatish Balay x += bs; 10582d61bbb3SSatish Balay workt = work; 10592d61bbb3SSatish Balay for (j=0; j<n; j++) { 10602d61bbb3SSatish Balay zb = z + bs*(*idx++); 10612d61bbb3SSatish Balay for (k=0; k<bs; k++) zb[k] += workt[k] ; 10622d61bbb3SSatish Balay workt += bs; 10632d61bbb3SSatish Balay } 10642d61bbb3SSatish Balay } 10652d61bbb3SSatish Balay } 10662d61bbb3SSatish Balay } 10672d61bbb3SSatish Balay ierr = VecRestoreArray(xx,&xg);CHKERRQ(ierr); 10682d61bbb3SSatish Balay ierr = VecRestoreArray(zz,&zg);CHKERRQ(ierr); 1069b0a32e0cSBarry Smith PetscLogFlops(2*a->nz*a->bs2 - A->n); 10702d61bbb3SSatish Balay PetscFunctionReturn(0); 10712d61bbb3SSatish Balay } 10722d61bbb3SSatish Balay 10734a2ae208SSatish Balay #undef __FUNCT__ 10744a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_SeqBAIJ" 10757c922b88SBarry Smith int MatMultTransposeAdd_SeqBAIJ(Mat A,Vec xx,Vec yy,Vec zz) 10762d61bbb3SSatish Balay 10772d61bbb3SSatish Balay { 10782d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 107987828ca2SBarry Smith PetscScalar *xg,*zg,*zb,*x,*z,*xb,x1,x2,x3,x4,x5; 10803f1db9ecSBarry Smith MatScalar *v; 1081f1af5d2fSBarry Smith int mbs=a->mbs,i,*idx,*ii,*ai=a->i,rval,bs=a->bs,j,n,bs2=a->bs2,*ib,ierr; 10822d61bbb3SSatish Balay 10832d61bbb3SSatish Balay PetscFunctionBegin; 1084ef66eb69SBarry Smith if (yy != zz) { ierr = VecCopy(yy,zz);CHKERRQ(ierr); } 1085e1311b90SBarry Smith ierr = VecGetArray(xx,&xg);CHKERRQ(ierr); x = xg; 1086e1311b90SBarry Smith ierr = VecGetArray(zz,&zg);CHKERRQ(ierr); z = zg; 10872d61bbb3SSatish Balay 10882d61bbb3SSatish Balay 10892d61bbb3SSatish Balay idx = a->j; 10902d61bbb3SSatish Balay v = a->a; 10912d61bbb3SSatish Balay ii = a->i; 1092f1af5d2fSBarry Smith xb = x; 10932d61bbb3SSatish Balay 10942d61bbb3SSatish Balay switch (bs) { 10952d61bbb3SSatish Balay case 1: 10962d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 10972d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 1098f1af5d2fSBarry Smith x1 = xb[0]; 10992d61bbb3SSatish Balay ib = idx + ai[i]; 11002d61bbb3SSatish Balay for (j=0; j<n; j++) { 11012d61bbb3SSatish Balay rval = ib[j]; 1102f1af5d2fSBarry Smith z[rval] += *v * x1; 1103f1af5d2fSBarry Smith v++; 11042d61bbb3SSatish Balay } 1105f1af5d2fSBarry Smith xb++; 11062d61bbb3SSatish Balay } 11072d61bbb3SSatish Balay break; 11082d61bbb3SSatish Balay case 2: 11092d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 11102d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 1111f1af5d2fSBarry Smith x1 = xb[0]; x2 = xb[1]; 11122d61bbb3SSatish Balay ib = idx + ai[i]; 11132d61bbb3SSatish Balay for (j=0; j<n; j++) { 11142d61bbb3SSatish Balay rval = ib[j]*2; 11152d61bbb3SSatish Balay z[rval++] += v[0]*x1 + v[1]*x2; 11162d61bbb3SSatish Balay z[rval++] += v[2]*x1 + v[3]*x2; 11172d61bbb3SSatish Balay v += 4; 11182d61bbb3SSatish Balay } 1119f1af5d2fSBarry Smith xb += 2; 11202d61bbb3SSatish Balay } 11212d61bbb3SSatish Balay break; 11222d61bbb3SSatish Balay case 3: 11232d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 11242d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 1125f1af5d2fSBarry Smith x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; 11262d61bbb3SSatish Balay ib = idx + ai[i]; 11272d61bbb3SSatish Balay for (j=0; j<n; j++) { 11282d61bbb3SSatish Balay rval = ib[j]*3; 11292d61bbb3SSatish Balay z[rval++] += v[0]*x1 + v[1]*x2 + v[2]*x3; 11302d61bbb3SSatish Balay z[rval++] += v[3]*x1 + v[4]*x2 + v[5]*x3; 11312d61bbb3SSatish Balay z[rval++] += v[6]*x1 + v[7]*x2 + v[8]*x3; 11322d61bbb3SSatish Balay v += 9; 11332d61bbb3SSatish Balay } 1134f1af5d2fSBarry Smith xb += 3; 11352d61bbb3SSatish Balay } 11362d61bbb3SSatish Balay break; 11372d61bbb3SSatish Balay case 4: 11382d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 11392d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 1140f1af5d2fSBarry Smith x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; 11412d61bbb3SSatish Balay ib = idx + ai[i]; 11422d61bbb3SSatish Balay for (j=0; j<n; j++) { 11432d61bbb3SSatish Balay rval = ib[j]*4; 11442d61bbb3SSatish Balay z[rval++] += v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4; 11452d61bbb3SSatish Balay z[rval++] += v[4]*x1 + v[5]*x2 + v[6]*x3 + v[7]*x4; 11462d61bbb3SSatish Balay z[rval++] += v[8]*x1 + v[9]*x2 + v[10]*x3 + v[11]*x4; 11472d61bbb3SSatish Balay z[rval++] += v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4; 11482d61bbb3SSatish Balay v += 16; 11492d61bbb3SSatish Balay } 1150f1af5d2fSBarry Smith xb += 4; 11512d61bbb3SSatish Balay } 11522d61bbb3SSatish Balay break; 11532d61bbb3SSatish Balay case 5: 11542d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 11552d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 1156f1af5d2fSBarry Smith x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; 11572d61bbb3SSatish Balay x4 = xb[3]; x5 = xb[4]; 11582d61bbb3SSatish Balay ib = idx + ai[i]; 11592d61bbb3SSatish Balay for (j=0; j<n; j++) { 11602d61bbb3SSatish Balay rval = ib[j]*5; 11612d61bbb3SSatish Balay z[rval++] += v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4 + v[4]*x5; 11622d61bbb3SSatish Balay z[rval++] += v[5]*x1 + v[6]*x2 + v[7]*x3 + v[8]*x4 + v[9]*x5; 11632d61bbb3SSatish Balay z[rval++] += v[10]*x1 + v[11]*x2 + v[12]*x3 + v[13]*x4 + v[14]*x5; 11642d61bbb3SSatish Balay z[rval++] += v[15]*x1 + v[16]*x2 + v[17]*x3 + v[18]*x4 + v[19]*x5; 11652d61bbb3SSatish Balay z[rval++] += v[20]*x1 + v[21]*x2 + v[22]*x3 + v[23]*x4 + v[24]*x5; 11662d61bbb3SSatish Balay v += 25; 11672d61bbb3SSatish Balay } 1168f1af5d2fSBarry Smith xb += 5; 11692d61bbb3SSatish Balay } 11702d61bbb3SSatish Balay break; 1171f1af5d2fSBarry Smith default: { /* block sizes larger then 5 by 5 are handled by BLAS */ 11723f1db9ecSBarry Smith int ncols,k; 117387828ca2SBarry Smith PetscScalar *work,*workt; 11743f1db9ecSBarry Smith 11752d61bbb3SSatish Balay if (!a->mult_work) { 1176273d9f13SBarry Smith k = PetscMax(A->m,A->n); 117787828ca2SBarry Smith ierr = PetscMalloc((k+1)*sizeof(PetscScalar),&a->mult_work);CHKERRQ(ierr); 11782d61bbb3SSatish Balay } 11792d61bbb3SSatish Balay work = a->mult_work; 11802d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 11812d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 11822d61bbb3SSatish Balay ncols = n*bs; 118387828ca2SBarry Smith ierr = PetscMemzero(work,ncols*sizeof(PetscScalar));CHKERRQ(ierr); 11843f1db9ecSBarry Smith Kernel_w_gets_w_plus_trans_Ar_times_v(bs,ncols,x,v,work); 11853f1db9ecSBarry Smith /* LAgemv_("T",&bs,&ncols,&_DOne,v,&bs,x,&_One,&_DOne,work,&_One); */ 11862d61bbb3SSatish Balay v += n*bs2; 11872d61bbb3SSatish Balay x += bs; 11882d61bbb3SSatish Balay workt = work; 11892d61bbb3SSatish Balay for (j=0; j<n; j++) { 11902d61bbb3SSatish Balay zb = z + bs*(*idx++); 11912d61bbb3SSatish Balay for (k=0; k<bs; k++) zb[k] += workt[k] ; 11922d61bbb3SSatish Balay workt += bs; 11932d61bbb3SSatish Balay } 11942d61bbb3SSatish Balay } 11952d61bbb3SSatish Balay } 11962d61bbb3SSatish Balay } 11972d61bbb3SSatish Balay ierr = VecRestoreArray(xx,&xg);CHKERRQ(ierr); 11982d61bbb3SSatish Balay ierr = VecRestoreArray(zz,&zg);CHKERRQ(ierr); 1199b0a32e0cSBarry Smith PetscLogFlops(2*a->nz*a->bs2); 12002d61bbb3SSatish Balay PetscFunctionReturn(0); 12012d61bbb3SSatish Balay } 12022d61bbb3SSatish Balay 12034a2ae208SSatish Balay #undef __FUNCT__ 12044a2ae208SSatish Balay #define __FUNCT__ "MatScale_SeqBAIJ" 120587828ca2SBarry Smith int MatScale_SeqBAIJ(PetscScalar *alpha,Mat inA) 12062d61bbb3SSatish Balay { 12072d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)inA->data; 120887052302SDinesh Kaushik int totalnz = a->bs2*a->nz; 12093eda8832SBarry Smith #if defined(PETSC_USE_MAT_SINGLE) 12103eda8832SBarry Smith int i; 121187052302SDinesh Kaushik #else 121287052302SDinesh Kaushik int one = 1; 12133eda8832SBarry Smith #endif 12142d61bbb3SSatish Balay 12152d61bbb3SSatish Balay PetscFunctionBegin; 12163eda8832SBarry Smith #if defined(PETSC_USE_MAT_SINGLE) 12173eda8832SBarry Smith for (i=0; i<totalnz; i++) a->a[i] *= *alpha; 12183eda8832SBarry Smith #else 12192d61bbb3SSatish Balay BLscal_(&totalnz,alpha,a->a,&one); 12203eda8832SBarry Smith #endif 1221b0a32e0cSBarry Smith PetscLogFlops(totalnz); 12222d61bbb3SSatish Balay PetscFunctionReturn(0); 12232d61bbb3SSatish Balay } 12242d61bbb3SSatish Balay 12254a2ae208SSatish Balay #undef __FUNCT__ 12264a2ae208SSatish Balay #define __FUNCT__ "MatNorm_SeqBAIJ" 1227329f5518SBarry Smith int MatNorm_SeqBAIJ(Mat A,NormType type,PetscReal *norm) 12282d61bbb3SSatish Balay { 12292d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 12303f1db9ecSBarry Smith MatScalar *v = a->a; 1231329f5518SBarry Smith PetscReal sum = 0.0; 12320e90e235SBarry Smith int i,j,k,bs = a->bs,nz=a->nz,bs2=a->bs2,k1; 12332d61bbb3SSatish Balay 12342d61bbb3SSatish Balay PetscFunctionBegin; 12352d61bbb3SSatish Balay if (type == NORM_FROBENIUS) { 12362d61bbb3SSatish Balay for (i=0; i< bs2*nz; i++) { 1237aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 1238329f5518SBarry Smith sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 12392d61bbb3SSatish Balay #else 12402d61bbb3SSatish Balay sum += (*v)*(*v); v++; 12412d61bbb3SSatish Balay #endif 12422d61bbb3SSatish Balay } 12432d61bbb3SSatish Balay *norm = sqrt(sum); 1244596552b5SBarry Smith } else if (type == NORM_INFINITY) { /* maximum row sum */ 1245596552b5SBarry Smith *norm = 0.0; 1246596552b5SBarry Smith for (k=0; k<bs; k++) { 124774f84c7bSSatish Balay for (j=0; j<a->mbs; j++) { 1248596552b5SBarry Smith v = a->a + bs2*a->i[j] + k; 1249596552b5SBarry Smith sum = 0.0; 1250596552b5SBarry Smith for (i=0; i<a->i[j+1]-a->i[j]; i++) { 12510e90e235SBarry Smith for (k1=0; k1<bs; k1++){ 1252596552b5SBarry Smith sum += PetscAbsScalar(*v); 1253596552b5SBarry Smith v += bs; 12542d61bbb3SSatish Balay } 12550e90e235SBarry Smith } 1256596552b5SBarry Smith if (sum > *norm) *norm = sum; 1257596552b5SBarry Smith } 1258596552b5SBarry Smith } 1259596552b5SBarry Smith } else { 126029bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"No support for this norm yet"); 12612d61bbb3SSatish Balay } 12622d61bbb3SSatish Balay PetscFunctionReturn(0); 12632d61bbb3SSatish Balay } 12642d61bbb3SSatish Balay 12652d61bbb3SSatish Balay 12664a2ae208SSatish Balay #undef __FUNCT__ 12674a2ae208SSatish Balay #define __FUNCT__ "MatEqual_SeqBAIJ" 12682d61bbb3SSatish Balay int MatEqual_SeqBAIJ(Mat A,Mat B,PetscTruth* flg) 12692d61bbb3SSatish Balay { 12702d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data,*b = (Mat_SeqBAIJ *)B->data; 12710f5bd95cSBarry Smith int ierr; 1272273d9f13SBarry Smith PetscTruth flag; 12732d61bbb3SSatish Balay 12742d61bbb3SSatish Balay PetscFunctionBegin; 1275273d9f13SBarry Smith ierr = PetscTypeCompare((PetscObject)B,MATSEQBAIJ,&flag);CHKERRQ(ierr); 1276273d9f13SBarry Smith if (!flag) SETERRQ(PETSC_ERR_ARG_INCOMP,"Matrices must be same type"); 12772d61bbb3SSatish Balay 12782d61bbb3SSatish Balay /* If the matrix/block dimensions are not equal, or no of nonzeros or shift */ 1279273d9f13SBarry Smith if ((A->m != B->m) || (A->n != B->n) || (a->bs != b->bs)|| (a->nz != b->nz)) { 1280273d9f13SBarry Smith *flg = PETSC_FALSE; 1281273d9f13SBarry Smith PetscFunctionReturn(0); 12822d61bbb3SSatish Balay } 12832d61bbb3SSatish Balay 12842d61bbb3SSatish Balay /* if the a->i are the same */ 12850f5bd95cSBarry Smith ierr = PetscMemcmp(a->i,b->i,(a->mbs+1)*sizeof(int),flg);CHKERRQ(ierr); 12860f5bd95cSBarry Smith if (*flg == PETSC_FALSE) { 12870f5bd95cSBarry Smith PetscFunctionReturn(0); 12882d61bbb3SSatish Balay } 12892d61bbb3SSatish Balay 12902d61bbb3SSatish Balay /* if a->j are the same */ 12910f5bd95cSBarry Smith ierr = PetscMemcmp(a->j,b->j,(a->nz)*sizeof(int),flg);CHKERRQ(ierr); 12920f5bd95cSBarry Smith if (*flg == PETSC_FALSE) { 12930f5bd95cSBarry Smith PetscFunctionReturn(0); 12942d61bbb3SSatish Balay } 12952d61bbb3SSatish Balay /* if a->a are the same */ 129687828ca2SBarry Smith ierr = PetscMemcmp(a->a,b->a,(a->nz)*(a->bs)*(a->bs)*sizeof(PetscScalar),flg);CHKERRQ(ierr); 12972d61bbb3SSatish Balay PetscFunctionReturn(0); 12982d61bbb3SSatish Balay 12992d61bbb3SSatish Balay } 13002d61bbb3SSatish Balay 13014a2ae208SSatish Balay #undef __FUNCT__ 13024a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_SeqBAIJ" 13032d61bbb3SSatish Balay int MatGetDiagonal_SeqBAIJ(Mat A,Vec v) 13042d61bbb3SSatish Balay { 13052d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1306e1311b90SBarry Smith int ierr,i,j,k,n,row,bs,*ai,*aj,ambs,bs2; 130787828ca2SBarry Smith PetscScalar *x,zero = 0.0; 13083f1db9ecSBarry Smith MatScalar *aa,*aa_j; 13092d61bbb3SSatish Balay 13102d61bbb3SSatish Balay PetscFunctionBegin; 131129bbc08cSBarry Smith if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 13122d61bbb3SSatish Balay bs = a->bs; 13132d61bbb3SSatish Balay aa = a->a; 13142d61bbb3SSatish Balay ai = a->i; 13152d61bbb3SSatish Balay aj = a->j; 13162d61bbb3SSatish Balay ambs = a->mbs; 13172d61bbb3SSatish Balay bs2 = a->bs2; 13182d61bbb3SSatish Balay 1319e1311b90SBarry Smith ierr = VecSet(&zero,v);CHKERRQ(ierr); 1320e1311b90SBarry Smith ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1321e1311b90SBarry Smith ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr); 1322273d9f13SBarry Smith if (n != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector"); 13232d61bbb3SSatish Balay for (i=0; i<ambs; i++) { 13242d61bbb3SSatish Balay for (j=ai[i]; j<ai[i+1]; j++) { 13252d61bbb3SSatish Balay if (aj[j] == i) { 13262d61bbb3SSatish Balay row = i*bs; 13272d61bbb3SSatish Balay aa_j = aa+j*bs2; 13282d61bbb3SSatish Balay for (k=0; k<bs2; k+=(bs+1),row++) x[row] = aa_j[k]; 13292d61bbb3SSatish Balay break; 13302d61bbb3SSatish Balay } 13312d61bbb3SSatish Balay } 13322d61bbb3SSatish Balay } 133315091d37SBarry Smith ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 13342d61bbb3SSatish Balay PetscFunctionReturn(0); 13352d61bbb3SSatish Balay } 13362d61bbb3SSatish Balay 13374a2ae208SSatish Balay #undef __FUNCT__ 13384a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_SeqBAIJ" 13392d61bbb3SSatish Balay int MatDiagonalScale_SeqBAIJ(Mat A,Vec ll,Vec rr) 13402d61bbb3SSatish Balay { 13412d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 134287828ca2SBarry Smith PetscScalar *l,*r,x,*li,*ri; 13433f1db9ecSBarry Smith MatScalar *aa,*v; 1344e1311b90SBarry Smith int ierr,i,j,k,lm,rn,M,m,n,*ai,*aj,mbs,tmp,bs,bs2; 13452d61bbb3SSatish Balay 13462d61bbb3SSatish Balay PetscFunctionBegin; 13472d61bbb3SSatish Balay ai = a->i; 13482d61bbb3SSatish Balay aj = a->j; 13492d61bbb3SSatish Balay aa = a->a; 1350273d9f13SBarry Smith m = A->m; 1351273d9f13SBarry Smith n = A->n; 13522d61bbb3SSatish Balay bs = a->bs; 13532d61bbb3SSatish Balay mbs = a->mbs; 13542d61bbb3SSatish Balay bs2 = a->bs2; 13552d61bbb3SSatish Balay if (ll) { 1356e1311b90SBarry Smith ierr = VecGetArray(ll,&l);CHKERRQ(ierr); 1357e1311b90SBarry Smith ierr = VecGetLocalSize(ll,&lm);CHKERRQ(ierr); 135829bbc08cSBarry Smith if (lm != m) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vector wrong length"); 13592d61bbb3SSatish Balay for (i=0; i<mbs; i++) { /* for each block row */ 13602d61bbb3SSatish Balay M = ai[i+1] - ai[i]; 13612d61bbb3SSatish Balay li = l + i*bs; 13622d61bbb3SSatish Balay v = aa + bs2*ai[i]; 13632d61bbb3SSatish Balay for (j=0; j<M; j++) { /* for each block */ 13642d61bbb3SSatish Balay for (k=0; k<bs2; k++) { 13652d61bbb3SSatish Balay (*v++) *= li[k%bs]; 13662d61bbb3SSatish Balay } 13672d61bbb3SSatish Balay } 13682d61bbb3SSatish Balay } 136960d550acSSatish Balay ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr); 1370b0a32e0cSBarry Smith PetscLogFlops(a->nz); 13712d61bbb3SSatish Balay } 13722d61bbb3SSatish Balay 13732d61bbb3SSatish Balay if (rr) { 1374e1311b90SBarry Smith ierr = VecGetArray(rr,&r);CHKERRQ(ierr); 1375e1311b90SBarry Smith ierr = VecGetLocalSize(rr,&rn);CHKERRQ(ierr); 137629bbc08cSBarry Smith if (rn != n) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vector wrong length"); 13772d61bbb3SSatish Balay for (i=0; i<mbs; i++) { /* for each block row */ 13782d61bbb3SSatish Balay M = ai[i+1] - ai[i]; 13792d61bbb3SSatish Balay v = aa + bs2*ai[i]; 13802d61bbb3SSatish Balay for (j=0; j<M; j++) { /* for each block */ 13812d61bbb3SSatish Balay ri = r + bs*aj[ai[i]+j]; 13822d61bbb3SSatish Balay for (k=0; k<bs; k++) { 13832d61bbb3SSatish Balay x = ri[k]; 13842d61bbb3SSatish Balay for (tmp=0; tmp<bs; tmp++) (*v++) *= x; 13852d61bbb3SSatish Balay } 13862d61bbb3SSatish Balay } 13872d61bbb3SSatish Balay } 138860d550acSSatish Balay ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr); 1389b0a32e0cSBarry Smith PetscLogFlops(a->nz); 13902d61bbb3SSatish Balay } 13912d61bbb3SSatish Balay PetscFunctionReturn(0); 13922d61bbb3SSatish Balay } 13932d61bbb3SSatish Balay 13942d61bbb3SSatish Balay 13954a2ae208SSatish Balay #undef __FUNCT__ 13964a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_SeqBAIJ" 13972d61bbb3SSatish Balay int MatGetInfo_SeqBAIJ(Mat A,MatInfoType flag,MatInfo *info) 13982d61bbb3SSatish Balay { 13992d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 14002d61bbb3SSatish Balay 14012d61bbb3SSatish Balay PetscFunctionBegin; 1402273d9f13SBarry Smith info->rows_global = (double)A->m; 1403273d9f13SBarry Smith info->columns_global = (double)A->n; 1404273d9f13SBarry Smith info->rows_local = (double)A->m; 1405273d9f13SBarry Smith info->columns_local = (double)A->n; 14062d61bbb3SSatish Balay info->block_size = a->bs2; 14072d61bbb3SSatish Balay info->nz_allocated = a->maxnz; 14082d61bbb3SSatish Balay info->nz_used = a->bs2*a->nz; 14092d61bbb3SSatish Balay info->nz_unneeded = (double)(info->nz_allocated - info->nz_used); 14102d61bbb3SSatish Balay info->assemblies = A->num_ass; 14112d61bbb3SSatish Balay info->mallocs = a->reallocs; 14122d61bbb3SSatish Balay info->memory = A->mem; 14132d61bbb3SSatish Balay if (A->factor) { 14142d61bbb3SSatish Balay info->fill_ratio_given = A->info.fill_ratio_given; 14152d61bbb3SSatish Balay info->fill_ratio_needed = A->info.fill_ratio_needed; 14162d61bbb3SSatish Balay info->factor_mallocs = A->info.factor_mallocs; 14172d61bbb3SSatish Balay } else { 14182d61bbb3SSatish Balay info->fill_ratio_given = 0; 14192d61bbb3SSatish Balay info->fill_ratio_needed = 0; 14202d61bbb3SSatish Balay info->factor_mallocs = 0; 14212d61bbb3SSatish Balay } 14222d61bbb3SSatish Balay PetscFunctionReturn(0); 14232d61bbb3SSatish Balay } 14242d61bbb3SSatish Balay 14252d61bbb3SSatish Balay 14264a2ae208SSatish Balay #undef __FUNCT__ 14274a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_SeqBAIJ" 14282d61bbb3SSatish Balay int MatZeroEntries_SeqBAIJ(Mat A) 14292d61bbb3SSatish Balay { 14302d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1431549d3d68SSatish Balay int ierr; 14322d61bbb3SSatish Balay 14332d61bbb3SSatish Balay PetscFunctionBegin; 1434549d3d68SSatish Balay ierr = PetscMemzero(a->a,a->bs2*a->i[a->mbs]*sizeof(MatScalar));CHKERRQ(ierr); 14352d61bbb3SSatish Balay PetscFunctionReturn(0); 14362d61bbb3SSatish Balay } 1437