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/inline/spops.h" 53f1db9ecSBarry Smith #include "src/inline/ilu.h" 60a835dfdSSatish Balay #include "petscbt.h" 7cac129eeSSatish Balay 84a2ae208SSatish Balay #undef __FUNCT__ 94a2ae208SSatish Balay #define __FUNCT__ "MatIncreaseOverlap_SeqBAIJ" 10268466fbSBarry Smith int MatIncreaseOverlap_SeqBAIJ(Mat A,int is_max,IS is[],int ov) 11a3192f15SSatish Balay { 12a3192f15SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 13218c64b6SSatish Balay int row,i,j,k,l,m,n,*idx,ierr,*nidx,isz,val,ival; 14218c64b6SSatish Balay int start,end,*ai,*aj,bs,*nidx2; 15f1af5d2fSBarry Smith PetscBT table; 16a3192f15SSatish Balay 173a40ed3dSBarry Smith PetscFunctionBegin; 18a3192f15SSatish Balay m = a->mbs; 19a3192f15SSatish Balay ai = a->i; 20a3192f15SSatish Balay aj = a->j; 21218c64b6SSatish Balay bs = a->bs; 22a3192f15SSatish Balay 2329bbc08cSBarry Smith if (ov < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap specified"); 24a3192f15SSatish Balay 256831982aSBarry Smith ierr = PetscBTCreate(m,table);CHKERRQ(ierr); 26b0a32e0cSBarry Smith ierr = PetscMalloc((m+1)*sizeof(int),&nidx);CHKERRQ(ierr); 27b0a32e0cSBarry Smith ierr = PetscMalloc((A->m+1)*sizeof(int),&nidx2);CHKERRQ(ierr); 28a3192f15SSatish Balay 29a3192f15SSatish Balay for (i=0; i<is_max; i++) { 30a3192f15SSatish Balay /* Initialise the two local arrays */ 31a3192f15SSatish Balay isz = 0; 326831982aSBarry Smith ierr = PetscBTMemzero(m,table);CHKERRQ(ierr); 33a3192f15SSatish Balay 34a3192f15SSatish Balay /* Extract the indices, assume there can be duplicate entries */ 35a3192f15SSatish Balay ierr = ISGetIndices(is[i],&idx);CHKERRQ(ierr); 36b9b97703SBarry Smith ierr = ISGetLocalSize(is[i],&n);CHKERRQ(ierr); 37a3192f15SSatish Balay 38a3192f15SSatish Balay /* Enter these into the temp arrays i.e mark table[row], enter row into new index */ 39a3192f15SSatish Balay for (j=0; j<n ; ++j){ 40218c64b6SSatish Balay ival = idx[j]/bs; /* convert the indices into block indices */ 4129bbc08cSBarry Smith if (ival>m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"index greater than mat-dim"); 42f1af5d2fSBarry Smith if(!PetscBTLookupSet(table,ival)) { nidx[isz++] = ival;} 43a3192f15SSatish Balay } 44a3192f15SSatish Balay ierr = ISRestoreIndices(is[i],&idx);CHKERRQ(ierr); 45a3192f15SSatish Balay ierr = ISDestroy(is[i]);CHKERRQ(ierr); 46a3192f15SSatish Balay 47a3192f15SSatish Balay k = 0; 48a3192f15SSatish Balay for (j=0; j<ov; j++){ /* for each overlap*/ 49a3192f15SSatish Balay n = isz; 50a3192f15SSatish Balay for (; k<n ; k++){ /* do only those rows in nidx[k], which are not done yet */ 51a3192f15SSatish Balay row = nidx[k]; 52a3192f15SSatish Balay start = ai[row]; 53a3192f15SSatish Balay end = ai[row+1]; 54a3192f15SSatish Balay for (l = start; l<end ; l++){ 55a3192f15SSatish Balay val = aj[l]; 56f1af5d2fSBarry Smith if (!PetscBTLookupSet(table,val)) {nidx[isz++] = val;} 57a3192f15SSatish Balay } 58a3192f15SSatish Balay } 59a3192f15SSatish Balay } 60218c64b6SSatish Balay /* expand the Index Set */ 61218c64b6SSatish Balay for (j=0; j<isz; j++) { 62218c64b6SSatish Balay for (k=0; k<bs; k++) 63218c64b6SSatish Balay nidx2[j*bs+k] = nidx[j]*bs+k; 64218c64b6SSatish Balay } 65329f5518SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,isz*bs,nidx2,is+i);CHKERRQ(ierr); 66a3192f15SSatish Balay } 676831982aSBarry Smith ierr = PetscBTDestroy(table);CHKERRQ(ierr); 68606d414cSSatish Balay ierr = PetscFree(nidx);CHKERRQ(ierr); 69606d414cSSatish Balay ierr = PetscFree(nidx2);CHKERRQ(ierr); 703a40ed3dSBarry Smith PetscFunctionReturn(0); 71a3192f15SSatish Balay } 721c351548SSatish Balay 734a2ae208SSatish Balay #undef __FUNCT__ 744a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqBAIJ_Private" 757b2a1423SBarry Smith int MatGetSubMatrix_SeqBAIJ_Private(Mat A,IS isrow,IS iscol,int cs,MatReuse scall,Mat *B) 76736121d4SSatish Balay { 77736121d4SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*c; 7840011551SBarry Smith int *smap,i,k,kstart,kend,ierr,oldcols = a->nbs,*lens; 79218c64b6SSatish Balay int row,mat_i,*mat_j,tcol,*mat_ilen; 80736121d4SSatish Balay int *irow,*icol,nrows,ncols,*ssmap,bs=a->bs,bs2=a->bs2; 81218c64b6SSatish Balay int *aj = a->j,*ai = a->i; 823f1db9ecSBarry Smith MatScalar *mat_a; 83736121d4SSatish Balay Mat C; 840f5bd95cSBarry Smith PetscTruth flag; 85736121d4SSatish Balay 863a40ed3dSBarry Smith PetscFunctionBegin; 87888f2ed8SSatish Balay ierr = ISSorted(iscol,(PetscTruth*)&i);CHKERRQ(ierr); 8829bbc08cSBarry Smith if (!i) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"IS is not sorted"); 89736121d4SSatish Balay 90736121d4SSatish Balay ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr); 91218c64b6SSatish Balay ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr); 92b9b97703SBarry Smith ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr); 93b9b97703SBarry Smith ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr); 94736121d4SSatish Balay 95b0a32e0cSBarry Smith ierr = PetscMalloc((1+oldcols)*sizeof(int),&smap);CHKERRQ(ierr); 96736121d4SSatish Balay ssmap = smap; 97b0a32e0cSBarry Smith ierr = PetscMalloc((1+nrows)*sizeof(int),&lens);CHKERRQ(ierr); 98549d3d68SSatish Balay ierr = PetscMemzero(smap,oldcols*sizeof(int));CHKERRQ(ierr); 99736121d4SSatish Balay for (i=0; i<ncols; i++) smap[icol[i]] = i+1; 100736121d4SSatish Balay /* determine lens of each row */ 101736121d4SSatish Balay for (i=0; i<nrows; i++) { 102736121d4SSatish Balay kstart = ai[irow[i]]; 103736121d4SSatish Balay kend = kstart + a->ilen[irow[i]]; 104736121d4SSatish Balay lens[i] = 0; 105736121d4SSatish Balay for (k=kstart; k<kend; k++) { 106736121d4SSatish Balay if (ssmap[aj[k]]) { 107736121d4SSatish Balay lens[i]++; 108736121d4SSatish Balay } 109736121d4SSatish Balay } 110736121d4SSatish Balay } 111736121d4SSatish Balay /* Create and fill new matrix */ 112736121d4SSatish Balay if (scall == MAT_REUSE_MATRIX) { 113736121d4SSatish Balay c = (Mat_SeqBAIJ *)((*B)->data); 114736121d4SSatish Balay 11529bbc08cSBarry Smith if (c->mbs!=nrows || c->nbs!=ncols || c->bs!=bs) SETERRQ(PETSC_ERR_ARG_SIZ,"Submatrix wrong size"); 1160f5bd95cSBarry Smith ierr = PetscMemcmp(c->ilen,lens,c->mbs *sizeof(int),&flag);CHKERRQ(ierr); 1170f5bd95cSBarry Smith if (flag == PETSC_FALSE) { 11829bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong no of nonzeros"); 119736121d4SSatish Balay } 120549d3d68SSatish Balay ierr = PetscMemzero(c->ilen,c->mbs*sizeof(int));CHKERRQ(ierr); 121736121d4SSatish Balay C = *B; 1223a40ed3dSBarry Smith } else { 123*e2d9671bSKris Buschelman ierr = MatCreate(A->comm,nrows*bs,ncols*bs,PETSC_DETERMINE,PETSC_DETERMINE,&C);CHKERRQ(ierr); 124*e2d9671bSKris Buschelman ierr = MatSetType(C,A->type_name);CHKERRQ(ierr); 125*e2d9671bSKris Buschelman ierr = MatSeqBAIJSetPreallocation(C,bs,0,lens);CHKERRQ(ierr); 126736121d4SSatish Balay } 127736121d4SSatish Balay c = (Mat_SeqBAIJ *)(C->data); 128736121d4SSatish Balay for (i=0; i<nrows; i++) { 129736121d4SSatish Balay row = irow[i]; 130736121d4SSatish Balay kstart = ai[row]; 131736121d4SSatish Balay kend = kstart + a->ilen[row]; 132736121d4SSatish Balay mat_i = c->i[i]; 133736121d4SSatish Balay mat_j = c->j + mat_i; 134218c64b6SSatish Balay mat_a = c->a + mat_i*bs2; 135736121d4SSatish Balay mat_ilen = c->ilen + i; 136736121d4SSatish Balay for (k=kstart; k<kend; k++) { 137736121d4SSatish Balay if ((tcol=ssmap[a->j[k]])) { 138736121d4SSatish Balay *mat_j++ = tcol - 1; 139549d3d68SSatish Balay ierr = PetscMemcpy(mat_a,a->a+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 140549d3d68SSatish Balay mat_a += bs2; 141736121d4SSatish Balay (*mat_ilen)++; 142736121d4SSatish Balay } 143736121d4SSatish Balay } 144736121d4SSatish Balay } 145218c64b6SSatish Balay 146736121d4SSatish Balay /* Free work space */ 147736121d4SSatish Balay ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr); 148606d414cSSatish Balay ierr = PetscFree(smap);CHKERRQ(ierr); 149606d414cSSatish Balay ierr = PetscFree(lens);CHKERRQ(ierr); 1506d4a8577SBarry Smith ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1516d4a8577SBarry Smith ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 152736121d4SSatish Balay 153736121d4SSatish Balay ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr); 154736121d4SSatish Balay *B = C; 1553a40ed3dSBarry Smith PetscFunctionReturn(0); 156736121d4SSatish Balay } 157736121d4SSatish Balay 1584a2ae208SSatish Balay #undef __FUNCT__ 1594a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqBAIJ" 1607b2a1423SBarry Smith int MatGetSubMatrix_SeqBAIJ(Mat A,IS isrow,IS iscol,int cs,MatReuse scall,Mat *B) 161218c64b6SSatish Balay { 162218c64b6SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 163218c64b6SSatish Balay IS is1,is2; 164218c64b6SSatish Balay int *vary,*iary,*irow,*icol,nrows,ncols,i,ierr,bs=a->bs,count; 165218c64b6SSatish Balay 1663a40ed3dSBarry Smith PetscFunctionBegin; 167218c64b6SSatish Balay ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr); 168218c64b6SSatish Balay ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr); 169b9b97703SBarry Smith ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr); 170b9b97703SBarry Smith ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr); 171218c64b6SSatish Balay 172218c64b6SSatish Balay /* Verify if the indices corespond to each element in a block 173218c64b6SSatish Balay and form the IS with compressed IS */ 174b0a32e0cSBarry Smith ierr = PetscMalloc(2*(a->mbs+1)*sizeof(int),&vary);CHKERRQ(ierr); 175218c64b6SSatish Balay iary = vary + a->mbs; 176549d3d68SSatish Balay ierr = PetscMemzero(vary,(a->mbs)*sizeof(int));CHKERRQ(ierr); 177218c64b6SSatish Balay for (i=0; i<nrows; i++) vary[irow[i]/bs]++; 178218c64b6SSatish Balay count = 0; 179218c64b6SSatish Balay for (i=0; i<a->mbs; i++) { 180ac355199SBarry Smith if (vary[i]!=0 && vary[i]!=bs) SETERRQ(1,"Index set does not match blocks"); 181218c64b6SSatish Balay if (vary[i]==bs) iary[count++] = i; 182218c64b6SSatish Balay } 183029af93fSBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,count,iary,&is1);CHKERRQ(ierr); 184218c64b6SSatish Balay 185549d3d68SSatish Balay ierr = PetscMemzero(vary,(a->mbs)*sizeof(int));CHKERRQ(ierr); 186218c64b6SSatish Balay for (i=0; i<ncols; i++) vary[icol[i]/bs]++; 187218c64b6SSatish Balay count = 0; 188218c64b6SSatish Balay for (i=0; i<a->mbs; i++) { 189ac355199SBarry Smith if (vary[i]!=0 && vary[i]!=bs) SETERRQ(1,"Internal error in PETSc"); 190218c64b6SSatish Balay if (vary[i]==bs) iary[count++] = i; 191218c64b6SSatish Balay } 192029af93fSBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,count,iary,&is2);CHKERRQ(ierr); 193218c64b6SSatish Balay ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr); 194218c64b6SSatish Balay ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr); 195606d414cSSatish Balay ierr = PetscFree(vary);CHKERRQ(ierr); 196218c64b6SSatish Balay 1976a6a5d1dSBarry Smith ierr = MatGetSubMatrix_SeqBAIJ_Private(A,is1,is2,cs,scall,B);CHKERRQ(ierr); 198218c64b6SSatish Balay ISDestroy(is1); 199218c64b6SSatish Balay ISDestroy(is2); 2003a40ed3dSBarry Smith PetscFunctionReturn(0); 201218c64b6SSatish Balay } 202218c64b6SSatish Balay 2034a2ae208SSatish Balay #undef __FUNCT__ 2044a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_SeqBAIJ" 205268466fbSBarry Smith int MatGetSubMatrices_SeqBAIJ(Mat A,int n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[]) 206736121d4SSatish Balay { 207736121d4SSatish Balay int ierr,i; 208736121d4SSatish Balay 2093a40ed3dSBarry Smith PetscFunctionBegin; 210736121d4SSatish Balay if (scall == MAT_INITIAL_MATRIX) { 21182502324SSatish Balay ierr = PetscMalloc((n+1)*sizeof(Mat),B);CHKERRQ(ierr); 212736121d4SSatish Balay } 213736121d4SSatish Balay 214736121d4SSatish Balay for (i=0; i<n; i++) { 2156a6a5d1dSBarry Smith ierr = MatGetSubMatrix_SeqBAIJ(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr); 216736121d4SSatish Balay } 2173a40ed3dSBarry Smith PetscFunctionReturn(0); 218736121d4SSatish Balay } 219218c64b6SSatish Balay 220218c64b6SSatish Balay 2212d61bbb3SSatish Balay /* -------------------------------------------------------*/ 2222d61bbb3SSatish Balay /* Should check that shapes of vectors and matrices match */ 2232d61bbb3SSatish Balay /* -------------------------------------------------------*/ 224d9eff348SSatish Balay #include "petscblaslapack.h" 2252d61bbb3SSatish Balay 2264a2ae208SSatish Balay #undef __FUNCT__ 2274a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_1" 2282d61bbb3SSatish Balay int MatMult_SeqBAIJ_1(Mat A,Vec xx,Vec zz) 2292d61bbb3SSatish Balay { 2302d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 23187828ca2SBarry Smith PetscScalar *x,*z,sum; 2323f1db9ecSBarry Smith MatScalar *v; 233e1311b90SBarry Smith int mbs=a->mbs,i,*idx,*ii,n,ierr; 2342d61bbb3SSatish Balay 2352d61bbb3SSatish Balay PetscFunctionBegin; 236b1d4fb26SBarry Smith ierr = VecGetArrayFast(xx,&x);CHKERRQ(ierr); 237b1d4fb26SBarry Smith ierr = VecGetArrayFast(zz,&z);CHKERRQ(ierr); 2382d61bbb3SSatish Balay 2392d61bbb3SSatish Balay idx = a->j; 2402d61bbb3SSatish Balay v = a->a; 2412d61bbb3SSatish Balay ii = a->i; 2422d61bbb3SSatish Balay 2432d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 2442d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 2452d61bbb3SSatish Balay sum = 0.0; 2462d61bbb3SSatish Balay while (n--) sum += *v++ * x[*idx++]; 2472d61bbb3SSatish Balay z[i] = sum; 2482d61bbb3SSatish Balay } 249b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(xx,&x);CHKERRQ(ierr); 250b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(zz,&z);CHKERRQ(ierr); 251b0a32e0cSBarry Smith PetscLogFlops(2*a->nz - A->m); 2522d61bbb3SSatish Balay PetscFunctionReturn(0); 2532d61bbb3SSatish Balay } 2542d61bbb3SSatish Balay 2554a2ae208SSatish Balay #undef __FUNCT__ 2564a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_2" 2572d61bbb3SSatish Balay int MatMult_SeqBAIJ_2(Mat A,Vec xx,Vec zz) 2582d61bbb3SSatish Balay { 2592d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 26087828ca2SBarry Smith PetscScalar *x,*z,*xb,sum1,sum2; 26187828ca2SBarry Smith PetscScalar x1,x2; 2623f1db9ecSBarry Smith MatScalar *v; 263e1311b90SBarry Smith int ierr,mbs=a->mbs,i,*idx,*ii,j,n; 2642d61bbb3SSatish Balay 2652d61bbb3SSatish Balay PetscFunctionBegin; 266b1d4fb26SBarry Smith ierr = VecGetArrayFast(xx,&x);CHKERRQ(ierr); 267b1d4fb26SBarry Smith ierr = VecGetArrayFast(zz,&z);CHKERRQ(ierr); 2682d61bbb3SSatish Balay 2692d61bbb3SSatish Balay idx = a->j; 2702d61bbb3SSatish Balay v = a->a; 2712d61bbb3SSatish Balay ii = a->i; 2722d61bbb3SSatish Balay 2732d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 2742d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 2752d61bbb3SSatish Balay sum1 = 0.0; sum2 = 0.0; 2762d61bbb3SSatish Balay for (j=0; j<n; j++) { 2772d61bbb3SSatish Balay xb = x + 2*(*idx++); x1 = xb[0]; x2 = xb[1]; 2782d61bbb3SSatish Balay sum1 += v[0]*x1 + v[2]*x2; 2792d61bbb3SSatish Balay sum2 += v[1]*x1 + v[3]*x2; 2802d61bbb3SSatish Balay v += 4; 2812d61bbb3SSatish Balay } 2822d61bbb3SSatish Balay z[0] = sum1; z[1] = sum2; 2832d61bbb3SSatish Balay z += 2; 2842d61bbb3SSatish Balay } 285b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(xx,&x);CHKERRQ(ierr); 286b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(zz,&z);CHKERRQ(ierr); 287b0a32e0cSBarry Smith PetscLogFlops(8*a->nz - A->m); 2882d61bbb3SSatish Balay PetscFunctionReturn(0); 2892d61bbb3SSatish Balay } 2902d61bbb3SSatish Balay 2914a2ae208SSatish Balay #undef __FUNCT__ 2924a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_3" 2932d61bbb3SSatish Balay int MatMult_SeqBAIJ_3(Mat A,Vec xx,Vec zz) 2942d61bbb3SSatish Balay { 2952d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 29687828ca2SBarry Smith PetscScalar *x,*z,*xb,sum1,sum2,sum3,x1,x2,x3; 2973f1db9ecSBarry Smith MatScalar *v; 298e1311b90SBarry Smith int ierr,mbs=a->mbs,i,*idx,*ii,j,n; 2992d61bbb3SSatish Balay 300b6410449SSatish Balay #if defined(PETSC_HAVE_PRAGMA_DISJOINT) 301fee21e36SBarry Smith #pragma disjoint(*v,*z,*xb) 302fee21e36SBarry Smith #endif 303fee21e36SBarry Smith 3042d61bbb3SSatish Balay PetscFunctionBegin; 305b1d4fb26SBarry Smith ierr = VecGetArrayFast(xx,&x);CHKERRQ(ierr); 306b1d4fb26SBarry Smith ierr = VecGetArrayFast(zz,&z);CHKERRQ(ierr); 3072d61bbb3SSatish Balay 3082d61bbb3SSatish Balay idx = a->j; 3092d61bbb3SSatish Balay v = a->a; 3102d61bbb3SSatish Balay ii = a->i; 3112d61bbb3SSatish Balay 3122d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 3132d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 3142d61bbb3SSatish Balay sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; 3152d61bbb3SSatish Balay for (j=0; j<n; j++) { 3162d61bbb3SSatish Balay xb = x + 3*(*idx++); x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; 3172d61bbb3SSatish Balay sum1 += v[0]*x1 + v[3]*x2 + v[6]*x3; 3182d61bbb3SSatish Balay sum2 += v[1]*x1 + v[4]*x2 + v[7]*x3; 3192d61bbb3SSatish Balay sum3 += v[2]*x1 + v[5]*x2 + v[8]*x3; 3202d61bbb3SSatish Balay v += 9; 3212d61bbb3SSatish Balay } 3222d61bbb3SSatish Balay z[0] = sum1; z[1] = sum2; z[2] = sum3; 3232d61bbb3SSatish Balay z += 3; 3242d61bbb3SSatish Balay } 325b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(xx,&x);CHKERRQ(ierr); 326b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(zz,&z);CHKERRQ(ierr); 327b0a32e0cSBarry Smith PetscLogFlops(18*a->nz - A->m); 3282d61bbb3SSatish Balay PetscFunctionReturn(0); 3292d61bbb3SSatish Balay } 3302d61bbb3SSatish Balay 3314a2ae208SSatish Balay #undef __FUNCT__ 3324a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_4" 3332d61bbb3SSatish Balay int MatMult_SeqBAIJ_4(Mat A,Vec xx,Vec zz) 3342d61bbb3SSatish Balay { 3352d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 33687828ca2SBarry Smith PetscScalar *x,*z,*xb,sum1,sum2,sum3,sum4,x1,x2,x3,x4; 3373f1db9ecSBarry Smith MatScalar *v; 338e1311b90SBarry Smith int ierr,mbs=a->mbs,i,*idx,*ii,j,n; 3392d61bbb3SSatish Balay 3402d61bbb3SSatish Balay PetscFunctionBegin; 341b1d4fb26SBarry Smith ierr = VecGetArrayFast(xx,&x);CHKERRQ(ierr); 342b1d4fb26SBarry Smith ierr = VecGetArrayFast(zz,&z);CHKERRQ(ierr); 3432d61bbb3SSatish Balay 3442d61bbb3SSatish Balay idx = a->j; 3452d61bbb3SSatish Balay v = a->a; 3462d61bbb3SSatish Balay ii = a->i; 3472d61bbb3SSatish Balay 3482d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 3492d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 3502d61bbb3SSatish Balay sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; 3512d61bbb3SSatish Balay for (j=0; j<n; j++) { 3522d61bbb3SSatish Balay xb = x + 4*(*idx++); 3532d61bbb3SSatish Balay x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; 3542d61bbb3SSatish Balay sum1 += v[0]*x1 + v[4]*x2 + v[8]*x3 + v[12]*x4; 3552d61bbb3SSatish Balay sum2 += v[1]*x1 + v[5]*x2 + v[9]*x3 + v[13]*x4; 3562d61bbb3SSatish Balay sum3 += v[2]*x1 + v[6]*x2 + v[10]*x3 + v[14]*x4; 3572d61bbb3SSatish Balay sum4 += v[3]*x1 + v[7]*x2 + v[11]*x3 + v[15]*x4; 3582d61bbb3SSatish Balay v += 16; 3592d61bbb3SSatish Balay } 3602d61bbb3SSatish Balay z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; 3612d61bbb3SSatish Balay z += 4; 3622d61bbb3SSatish Balay } 363b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(xx,&x);CHKERRQ(ierr); 364b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(zz,&z);CHKERRQ(ierr); 365b0a32e0cSBarry Smith PetscLogFlops(32*a->nz - A->m); 3662d61bbb3SSatish Balay PetscFunctionReturn(0); 3672d61bbb3SSatish Balay } 3682d61bbb3SSatish Balay 3694a2ae208SSatish Balay #undef __FUNCT__ 3704a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_5" 3712d61bbb3SSatish Balay int MatMult_SeqBAIJ_5(Mat A,Vec xx,Vec zz) 3722d61bbb3SSatish Balay { 3732d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 37487828ca2SBarry Smith PetscScalar sum1,sum2,sum3,sum4,sum5,x1,x2,x3,x4,x5,*xb,*z,*x; 3753f1db9ecSBarry Smith MatScalar *v; 376e1311b90SBarry Smith int ierr,mbs=a->mbs,i,*idx,*ii,j,n; 3772d61bbb3SSatish Balay 378433994e6SBarry Smith PetscFunctionBegin; 379b1d4fb26SBarry Smith ierr = VecGetArrayFast(xx,&x);CHKERRQ(ierr); 380b1d4fb26SBarry Smith ierr = VecGetArrayFast(zz,&z);CHKERRQ(ierr); 3812d61bbb3SSatish Balay 3822d61bbb3SSatish Balay idx = a->j; 3832d61bbb3SSatish Balay v = a->a; 3842d61bbb3SSatish Balay ii = a->i; 3852d61bbb3SSatish Balay 3862d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 3872d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 3882d61bbb3SSatish Balay sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; 3892d61bbb3SSatish Balay for (j=0; j<n; j++) { 3902d61bbb3SSatish Balay xb = x + 5*(*idx++); 3912d61bbb3SSatish Balay x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; 3922d61bbb3SSatish Balay sum1 += v[0]*x1 + v[5]*x2 + v[10]*x3 + v[15]*x4 + v[20]*x5; 3932d61bbb3SSatish Balay sum2 += v[1]*x1 + v[6]*x2 + v[11]*x3 + v[16]*x4 + v[21]*x5; 3942d61bbb3SSatish Balay sum3 += v[2]*x1 + v[7]*x2 + v[12]*x3 + v[17]*x4 + v[22]*x5; 3952d61bbb3SSatish Balay sum4 += v[3]*x1 + v[8]*x2 + v[13]*x3 + v[18]*x4 + v[23]*x5; 3962d61bbb3SSatish Balay sum5 += v[4]*x1 + v[9]*x2 + v[14]*x3 + v[19]*x4 + v[24]*x5; 3972d61bbb3SSatish Balay v += 25; 3982d61bbb3SSatish Balay } 3992d61bbb3SSatish Balay z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; 4002d61bbb3SSatish Balay z += 5; 4012d61bbb3SSatish Balay } 402b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(xx,&x);CHKERRQ(ierr); 403b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(zz,&z);CHKERRQ(ierr); 404b0a32e0cSBarry Smith PetscLogFlops(50*a->nz - A->m); 4052d61bbb3SSatish Balay PetscFunctionReturn(0); 4062d61bbb3SSatish Balay } 4072d61bbb3SSatish Balay 40815091d37SBarry Smith 4094a2ae208SSatish Balay #undef __FUNCT__ 4104a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_6" 41115091d37SBarry Smith int MatMult_SeqBAIJ_6(Mat A,Vec xx,Vec zz) 41215091d37SBarry Smith { 41315091d37SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 41487828ca2SBarry Smith PetscScalar *x,*z,*xb,sum1,sum2,sum3,sum4,sum5,sum6; 41587828ca2SBarry Smith PetscScalar x1,x2,x3,x4,x5,x6; 41615091d37SBarry Smith MatScalar *v; 41715091d37SBarry Smith int ierr,mbs=a->mbs,i,*idx,*ii,j,n; 41815091d37SBarry Smith 419433994e6SBarry Smith PetscFunctionBegin; 420b1d4fb26SBarry Smith ierr = VecGetArrayFast(xx,&x);CHKERRQ(ierr); 421b1d4fb26SBarry Smith ierr = VecGetArrayFast(zz,&z);CHKERRQ(ierr); 42215091d37SBarry Smith 42315091d37SBarry Smith idx = a->j; 42415091d37SBarry Smith v = a->a; 42515091d37SBarry Smith ii = a->i; 42615091d37SBarry Smith 42715091d37SBarry Smith for (i=0; i<mbs; i++) { 42815091d37SBarry Smith n = ii[1] - ii[0]; ii++; 42915091d37SBarry Smith sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0; 43015091d37SBarry Smith for (j=0; j<n; j++) { 43115091d37SBarry Smith xb = x + 6*(*idx++); 43215091d37SBarry Smith x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; 43315091d37SBarry Smith sum1 += v[0]*x1 + v[6]*x2 + v[12]*x3 + v[18]*x4 + v[24]*x5 + v[30]*x6; 43415091d37SBarry Smith sum2 += v[1]*x1 + v[7]*x2 + v[13]*x3 + v[19]*x4 + v[25]*x5 + v[31]*x6; 43515091d37SBarry Smith sum3 += v[2]*x1 + v[8]*x2 + v[14]*x3 + v[20]*x4 + v[26]*x5 + v[32]*x6; 43615091d37SBarry Smith sum4 += v[3]*x1 + v[9]*x2 + v[15]*x3 + v[21]*x4 + v[27]*x5 + v[33]*x6; 43715091d37SBarry Smith sum5 += v[4]*x1 + v[10]*x2 + v[16]*x3 + v[22]*x4 + v[28]*x5 + v[34]*x6; 43815091d37SBarry Smith sum6 += v[5]*x1 + v[11]*x2 + v[17]*x3 + v[23]*x4 + v[29]*x5 + v[35]*x6; 43915091d37SBarry Smith v += 36; 44015091d37SBarry Smith } 44115091d37SBarry Smith z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; 44215091d37SBarry Smith z += 6; 44315091d37SBarry Smith } 44415091d37SBarry Smith 445b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(xx,&x);CHKERRQ(ierr); 446b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(zz,&z);CHKERRQ(ierr); 447b0a32e0cSBarry Smith PetscLogFlops(72*a->nz - A->m); 44815091d37SBarry Smith PetscFunctionReturn(0); 44915091d37SBarry Smith } 4504a2ae208SSatish Balay #undef __FUNCT__ 4514a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_7" 4522d61bbb3SSatish Balay int MatMult_SeqBAIJ_7(Mat A,Vec xx,Vec zz) 4532d61bbb3SSatish Balay { 4542d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 45587828ca2SBarry Smith PetscScalar *x,*z,*xb,sum1,sum2,sum3,sum4,sum5,sum6,sum7; 45687828ca2SBarry Smith PetscScalar x1,x2,x3,x4,x5,x6,x7; 4573f1db9ecSBarry Smith MatScalar *v; 458e1311b90SBarry Smith int ierr,mbs=a->mbs,i,*idx,*ii,j,n; 4592d61bbb3SSatish Balay 460433994e6SBarry Smith PetscFunctionBegin; 461b1d4fb26SBarry Smith ierr = VecGetArrayFast(xx,&x);CHKERRQ(ierr); 462b1d4fb26SBarry Smith ierr = VecGetArrayFast(zz,&z);CHKERRQ(ierr); 4632d61bbb3SSatish Balay 4642d61bbb3SSatish Balay idx = a->j; 4652d61bbb3SSatish Balay v = a->a; 4662d61bbb3SSatish Balay ii = a->i; 4672d61bbb3SSatish Balay 4682d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 4692d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 4702d61bbb3SSatish Balay sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0; sum7 = 0.0; 4712d61bbb3SSatish Balay for (j=0; j<n; j++) { 4722d61bbb3SSatish Balay xb = x + 7*(*idx++); 4732d61bbb3SSatish Balay x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6]; 4742d61bbb3SSatish Balay sum1 += v[0]*x1 + v[7]*x2 + v[14]*x3 + v[21]*x4 + v[28]*x5 + v[35]*x6 + v[42]*x7; 4752d61bbb3SSatish Balay sum2 += v[1]*x1 + v[8]*x2 + v[15]*x3 + v[22]*x4 + v[29]*x5 + v[36]*x6 + v[43]*x7; 4762d61bbb3SSatish Balay sum3 += v[2]*x1 + v[9]*x2 + v[16]*x3 + v[23]*x4 + v[30]*x5 + v[37]*x6 + v[44]*x7; 4772d61bbb3SSatish Balay sum4 += v[3]*x1 + v[10]*x2 + v[17]*x3 + v[24]*x4 + v[31]*x5 + v[38]*x6 + v[45]*x7; 4782d61bbb3SSatish Balay sum5 += v[4]*x1 + v[11]*x2 + v[18]*x3 + v[25]*x4 + v[32]*x5 + v[39]*x6 + v[46]*x7; 4792d61bbb3SSatish Balay sum6 += v[5]*x1 + v[12]*x2 + v[19]*x3 + v[26]*x4 + v[33]*x5 + v[40]*x6 + v[47]*x7; 4802d61bbb3SSatish Balay sum7 += v[6]*x1 + v[13]*x2 + v[20]*x3 + v[27]*x4 + v[34]*x5 + v[41]*x6 + v[48]*x7; 4812d61bbb3SSatish Balay v += 49; 4822d61bbb3SSatish Balay } 4832d61bbb3SSatish Balay z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7; 4842d61bbb3SSatish Balay z += 7; 4852d61bbb3SSatish Balay } 4862d61bbb3SSatish Balay 487b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(xx,&x);CHKERRQ(ierr); 488b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(zz,&z);CHKERRQ(ierr); 489b0a32e0cSBarry Smith PetscLogFlops(98*a->nz - A->m); 4902d61bbb3SSatish Balay PetscFunctionReturn(0); 4912d61bbb3SSatish Balay } 4922d61bbb3SSatish Balay 4933f1db9ecSBarry Smith /* 4943f1db9ecSBarry Smith This will not work with MatScalar == float because it calls the BLAS 4953f1db9ecSBarry Smith */ 4964a2ae208SSatish Balay #undef __FUNCT__ 4974a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqBAIJ_N" 4982d61bbb3SSatish Balay int MatMult_SeqBAIJ_N(Mat A,Vec xx,Vec zz) 4992d61bbb3SSatish Balay { 5002d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 50187828ca2SBarry Smith PetscScalar *x,*z,*xb,*work,*workt; 5023f1db9ecSBarry Smith MatScalar *v; 503e1311b90SBarry Smith int ierr,mbs=a->mbs,i,*idx,*ii,bs=a->bs,j,n,bs2=a->bs2; 5043f1db9ecSBarry Smith int ncols,k; 5052d61bbb3SSatish Balay 5062d61bbb3SSatish Balay PetscFunctionBegin; 507b1d4fb26SBarry Smith ierr = VecGetArrayFast(xx,&x);CHKERRQ(ierr); 508b1d4fb26SBarry Smith ierr = VecGetArrayFast(zz,&z);CHKERRQ(ierr); 5092d61bbb3SSatish Balay 5102d61bbb3SSatish Balay idx = a->j; 5112d61bbb3SSatish Balay v = a->a; 5122d61bbb3SSatish Balay ii = a->i; 513218c64b6SSatish Balay 514218c64b6SSatish Balay 5152d61bbb3SSatish Balay if (!a->mult_work) { 516273d9f13SBarry Smith k = PetscMax(A->m,A->n); 51787828ca2SBarry Smith ierr = PetscMalloc((k+1)*sizeof(PetscScalar),&a->mult_work);CHKERRQ(ierr); 5182d61bbb3SSatish Balay } 5192d61bbb3SSatish Balay work = a->mult_work; 5202d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 5212d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 5222d61bbb3SSatish Balay ncols = n*bs; 5232d61bbb3SSatish Balay workt = work; 5242d61bbb3SSatish Balay for (j=0; j<n; j++) { 5252d61bbb3SSatish Balay xb = x + bs*(*idx++); 5262d61bbb3SSatish Balay for (k=0; k<bs; k++) workt[k] = xb[k]; 5272d61bbb3SSatish Balay workt += bs; 5282d61bbb3SSatish Balay } 5293f1db9ecSBarry Smith Kernel_w_gets_Ar_times_v(bs,ncols,work,v,z); 5303f1db9ecSBarry Smith /* LAgemv_("N",&bs,&ncols,&_DOne,v,&bs,work,&_One,&_DZero,z,&_One); */ 5312d61bbb3SSatish Balay v += n*bs2; 5322d61bbb3SSatish Balay z += bs; 5332d61bbb3SSatish Balay } 534b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(xx,&x);CHKERRQ(ierr); 535b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(zz,&z);CHKERRQ(ierr); 536b0a32e0cSBarry Smith PetscLogFlops(2*a->nz*bs2 - A->m); 5372d61bbb3SSatish Balay PetscFunctionReturn(0); 5382d61bbb3SSatish Balay } 5392d61bbb3SSatish Balay 5404a2ae208SSatish Balay #undef __FUNCT__ 5414a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_1" 5422d61bbb3SSatish Balay int MatMultAdd_SeqBAIJ_1(Mat A,Vec xx,Vec yy,Vec zz) 5432d61bbb3SSatish Balay { 5442d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 54587828ca2SBarry Smith PetscScalar *x,*y,*z,sum; 5463f1db9ecSBarry Smith MatScalar *v; 547e1311b90SBarry Smith int ierr,mbs=a->mbs,i,*idx,*ii,n; 5482d61bbb3SSatish Balay 5492d61bbb3SSatish Balay PetscFunctionBegin; 550b1d4fb26SBarry Smith ierr = VecGetArrayFast(xx,&x);CHKERRQ(ierr); 551b1d4fb26SBarry Smith ierr = VecGetArrayFast(yy,&y);CHKERRQ(ierr); 5522e8a6d31SBarry Smith if (zz != yy) { 553b1d4fb26SBarry Smith ierr = VecGetArrayFast(zz,&z);CHKERRQ(ierr); 5542e8a6d31SBarry Smith } else { 5552e8a6d31SBarry Smith z = y; 5562e8a6d31SBarry Smith } 5572d61bbb3SSatish Balay 5582d61bbb3SSatish Balay idx = a->j; 5592d61bbb3SSatish Balay v = a->a; 5602d61bbb3SSatish Balay ii = a->i; 5612d61bbb3SSatish Balay 5622d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 5632d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 5642d61bbb3SSatish Balay sum = y[i]; 5652d61bbb3SSatish Balay while (n--) sum += *v++ * x[*idx++]; 5662d61bbb3SSatish Balay z[i] = sum; 5672d61bbb3SSatish Balay } 568b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(xx,&x);CHKERRQ(ierr); 569b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(yy,&y);CHKERRQ(ierr); 5702e8a6d31SBarry Smith if (zz != yy) { 571b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(zz,&z);CHKERRQ(ierr); 5722e8a6d31SBarry Smith } 573b0a32e0cSBarry Smith PetscLogFlops(2*a->nz); 5742d61bbb3SSatish Balay PetscFunctionReturn(0); 5752d61bbb3SSatish Balay } 5762d61bbb3SSatish Balay 5774a2ae208SSatish Balay #undef __FUNCT__ 5784a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_2" 5792d61bbb3SSatish Balay int MatMultAdd_SeqBAIJ_2(Mat A,Vec xx,Vec yy,Vec zz) 5802d61bbb3SSatish Balay { 5812d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 58287828ca2SBarry Smith PetscScalar *x,*y,*z,*xb,sum1,sum2; 58387828ca2SBarry Smith PetscScalar x1,x2; 5843f1db9ecSBarry Smith MatScalar *v; 585e1311b90SBarry Smith int ierr,mbs=a->mbs,i,*idx,*ii,j,n; 5862d61bbb3SSatish Balay 5872d61bbb3SSatish Balay PetscFunctionBegin; 588b1d4fb26SBarry Smith ierr = VecGetArrayFast(xx,&x);CHKERRQ(ierr); 589b1d4fb26SBarry Smith ierr = VecGetArrayFast(yy,&y);CHKERRQ(ierr); 5902e8a6d31SBarry Smith if (zz != yy) { 591b1d4fb26SBarry Smith ierr = VecGetArrayFast(zz,&z);CHKERRQ(ierr); 5922e8a6d31SBarry Smith } else { 5932e8a6d31SBarry Smith z = y; 5942e8a6d31SBarry Smith } 5952d61bbb3SSatish Balay 5962d61bbb3SSatish Balay idx = a->j; 5972d61bbb3SSatish Balay v = a->a; 5982d61bbb3SSatish Balay ii = a->i; 5992d61bbb3SSatish Balay 6002d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 6012d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 6022d61bbb3SSatish Balay sum1 = y[0]; sum2 = y[1]; 6032d61bbb3SSatish Balay for (j=0; j<n; j++) { 6042d61bbb3SSatish Balay xb = x + 2*(*idx++); x1 = xb[0]; x2 = xb[1]; 6052d61bbb3SSatish Balay sum1 += v[0]*x1 + v[2]*x2; 6062d61bbb3SSatish Balay sum2 += v[1]*x1 + v[3]*x2; 6072d61bbb3SSatish Balay v += 4; 6082d61bbb3SSatish Balay } 6092d61bbb3SSatish Balay z[0] = sum1; z[1] = sum2; 6102d61bbb3SSatish Balay z += 2; y += 2; 6112d61bbb3SSatish Balay } 612b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(xx,&x);CHKERRQ(ierr); 613b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(yy,&y);CHKERRQ(ierr); 6142e8a6d31SBarry Smith if (zz != yy) { 615b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(zz,&z);CHKERRQ(ierr); 6162e8a6d31SBarry Smith } 617b0a32e0cSBarry Smith PetscLogFlops(4*a->nz); 6182d61bbb3SSatish Balay PetscFunctionReturn(0); 6192d61bbb3SSatish Balay } 6202d61bbb3SSatish Balay 6214a2ae208SSatish Balay #undef __FUNCT__ 6224a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_3" 6232d61bbb3SSatish Balay int MatMultAdd_SeqBAIJ_3(Mat A,Vec xx,Vec yy,Vec zz) 6242d61bbb3SSatish Balay { 6252d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 62687828ca2SBarry Smith PetscScalar *x,*y,*z,*xb,sum1,sum2,sum3,x1,x2,x3; 6273f1db9ecSBarry Smith MatScalar *v; 628e1311b90SBarry Smith int ierr,mbs=a->mbs,i,*idx,*ii,j,n; 6292d61bbb3SSatish Balay 6302d61bbb3SSatish Balay PetscFunctionBegin; 631b1d4fb26SBarry Smith ierr = VecGetArrayFast(xx,&x);CHKERRQ(ierr); 632b1d4fb26SBarry Smith ierr = VecGetArrayFast(yy,&y);CHKERRQ(ierr); 6332e8a6d31SBarry Smith if (zz != yy) { 634b1d4fb26SBarry Smith ierr = VecGetArrayFast(zz,&z);CHKERRQ(ierr); 6352e8a6d31SBarry Smith } else { 6362e8a6d31SBarry Smith z = y; 6372e8a6d31SBarry Smith } 6382d61bbb3SSatish Balay 6392d61bbb3SSatish Balay idx = a->j; 6402d61bbb3SSatish Balay v = a->a; 6412d61bbb3SSatish Balay ii = a->i; 6422d61bbb3SSatish Balay 6432d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 6442d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 6452d61bbb3SSatish Balay sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; 6462d61bbb3SSatish Balay for (j=0; j<n; j++) { 6472d61bbb3SSatish Balay xb = x + 3*(*idx++); x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; 6482d61bbb3SSatish Balay sum1 += v[0]*x1 + v[3]*x2 + v[6]*x3; 6492d61bbb3SSatish Balay sum2 += v[1]*x1 + v[4]*x2 + v[7]*x3; 6502d61bbb3SSatish Balay sum3 += v[2]*x1 + v[5]*x2 + v[8]*x3; 6512d61bbb3SSatish Balay v += 9; 6522d61bbb3SSatish Balay } 6532d61bbb3SSatish Balay z[0] = sum1; z[1] = sum2; z[2] = sum3; 6542d61bbb3SSatish Balay z += 3; y += 3; 6552d61bbb3SSatish Balay } 656b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(xx,&x);CHKERRQ(ierr); 657b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(yy,&y);CHKERRQ(ierr); 6582e8a6d31SBarry Smith if (zz != yy) { 659b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(zz,&z);CHKERRQ(ierr); 6602e8a6d31SBarry Smith } 661b0a32e0cSBarry Smith PetscLogFlops(18*a->nz); 6622d61bbb3SSatish Balay PetscFunctionReturn(0); 6632d61bbb3SSatish Balay } 6642d61bbb3SSatish Balay 6654a2ae208SSatish Balay #undef __FUNCT__ 6664a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_4" 6672d61bbb3SSatish Balay int MatMultAdd_SeqBAIJ_4(Mat A,Vec xx,Vec yy,Vec zz) 6682d61bbb3SSatish Balay { 6692d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 67087828ca2SBarry Smith PetscScalar *x,*y,*z,*xb,sum1,sum2,sum3,sum4,x1,x2,x3,x4; 6713f1db9ecSBarry Smith MatScalar *v; 672e1311b90SBarry Smith int ierr,mbs=a->mbs,i,*idx,*ii; 6732d61bbb3SSatish Balay int j,n; 6742d61bbb3SSatish Balay 6752d61bbb3SSatish Balay PetscFunctionBegin; 676b1d4fb26SBarry Smith ierr = VecGetArrayFast(xx,&x);CHKERRQ(ierr); 677b1d4fb26SBarry Smith ierr = VecGetArrayFast(yy,&y);CHKERRQ(ierr); 6782e8a6d31SBarry Smith if (zz != yy) { 679b1d4fb26SBarry Smith ierr = VecGetArrayFast(zz,&z);CHKERRQ(ierr); 6802e8a6d31SBarry Smith } else { 6812e8a6d31SBarry Smith z = y; 6822e8a6d31SBarry Smith } 6832d61bbb3SSatish Balay 6842d61bbb3SSatish Balay idx = a->j; 6852d61bbb3SSatish Balay v = a->a; 6862d61bbb3SSatish Balay ii = a->i; 6872d61bbb3SSatish Balay 6882d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 6892d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 6902d61bbb3SSatish Balay sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; 6912d61bbb3SSatish Balay for (j=0; j<n; j++) { 6922d61bbb3SSatish Balay xb = x + 4*(*idx++); 6932d61bbb3SSatish Balay x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; 6942d61bbb3SSatish Balay sum1 += v[0]*x1 + v[4]*x2 + v[8]*x3 + v[12]*x4; 6952d61bbb3SSatish Balay sum2 += v[1]*x1 + v[5]*x2 + v[9]*x3 + v[13]*x4; 6962d61bbb3SSatish Balay sum3 += v[2]*x1 + v[6]*x2 + v[10]*x3 + v[14]*x4; 6972d61bbb3SSatish Balay sum4 += v[3]*x1 + v[7]*x2 + v[11]*x3 + v[15]*x4; 6982d61bbb3SSatish Balay v += 16; 6992d61bbb3SSatish Balay } 7002d61bbb3SSatish Balay z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; 7012d61bbb3SSatish Balay z += 4; y += 4; 7022d61bbb3SSatish Balay } 703b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(xx,&x);CHKERRQ(ierr); 704b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(yy,&y);CHKERRQ(ierr); 7052e8a6d31SBarry Smith if (zz != yy) { 706b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(zz,&z);CHKERRQ(ierr); 7072e8a6d31SBarry Smith } 708b0a32e0cSBarry Smith PetscLogFlops(32*a->nz); 7092d61bbb3SSatish Balay PetscFunctionReturn(0); 7102d61bbb3SSatish Balay } 7112d61bbb3SSatish Balay 7124a2ae208SSatish Balay #undef __FUNCT__ 7134a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_5" 7142d61bbb3SSatish Balay int MatMultAdd_SeqBAIJ_5(Mat A,Vec xx,Vec yy,Vec zz) 7152d61bbb3SSatish Balay { 7162d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 71787828ca2SBarry Smith PetscScalar *x,*y,*z,*xb,sum1,sum2,sum3,sum4,sum5,x1,x2,x3,x4,x5; 7183f1db9ecSBarry Smith MatScalar *v; 719e1311b90SBarry Smith int ierr,mbs=a->mbs,i,*idx,*ii,j,n; 7202d61bbb3SSatish Balay 7212d61bbb3SSatish Balay PetscFunctionBegin; 722b1d4fb26SBarry Smith ierr = VecGetArrayFast(xx,&x);CHKERRQ(ierr); 723b1d4fb26SBarry Smith ierr = VecGetArrayFast(yy,&y);CHKERRQ(ierr); 7242e8a6d31SBarry Smith if (zz != yy) { 725b1d4fb26SBarry Smith ierr = VecGetArrayFast(zz,&z);CHKERRQ(ierr); 7262e8a6d31SBarry Smith } else { 7272e8a6d31SBarry Smith z = y; 7282e8a6d31SBarry Smith } 7292d61bbb3SSatish Balay 7302d61bbb3SSatish Balay idx = a->j; 7312d61bbb3SSatish Balay v = a->a; 7322d61bbb3SSatish Balay ii = a->i; 7332d61bbb3SSatish Balay 7342d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 7352d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 7362d61bbb3SSatish Balay sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4]; 7372d61bbb3SSatish Balay for (j=0; j<n; j++) { 7382d61bbb3SSatish Balay xb = x + 5*(*idx++); 7392d61bbb3SSatish Balay x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; 7402d61bbb3SSatish Balay sum1 += v[0]*x1 + v[5]*x2 + v[10]*x3 + v[15]*x4 + v[20]*x5; 7412d61bbb3SSatish Balay sum2 += v[1]*x1 + v[6]*x2 + v[11]*x3 + v[16]*x4 + v[21]*x5; 7422d61bbb3SSatish Balay sum3 += v[2]*x1 + v[7]*x2 + v[12]*x3 + v[17]*x4 + v[22]*x5; 7432d61bbb3SSatish Balay sum4 += v[3]*x1 + v[8]*x2 + v[13]*x3 + v[18]*x4 + v[23]*x5; 7442d61bbb3SSatish Balay sum5 += v[4]*x1 + v[9]*x2 + v[14]*x3 + v[19]*x4 + v[24]*x5; 7452d61bbb3SSatish Balay v += 25; 7462d61bbb3SSatish Balay } 7472d61bbb3SSatish Balay z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; 7482d61bbb3SSatish Balay z += 5; y += 5; 7492d61bbb3SSatish Balay } 750b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(xx,&x);CHKERRQ(ierr); 751b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(yy,&y);CHKERRQ(ierr); 7522e8a6d31SBarry Smith if (zz != yy) { 753b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(zz,&z);CHKERRQ(ierr); 7542e8a6d31SBarry Smith } 755b0a32e0cSBarry Smith PetscLogFlops(50*a->nz); 7562d61bbb3SSatish Balay PetscFunctionReturn(0); 7572d61bbb3SSatish Balay } 7584a2ae208SSatish Balay #undef __FUNCT__ 7594a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_6" 76015091d37SBarry Smith int MatMultAdd_SeqBAIJ_6(Mat A,Vec xx,Vec yy,Vec zz) 76115091d37SBarry Smith { 76215091d37SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 76387828ca2SBarry Smith PetscScalar *x,*y,*z,*xb,sum1,sum2,sum3,sum4,sum5,sum6; 76487828ca2SBarry Smith PetscScalar x1,x2,x3,x4,x5,x6; 76515091d37SBarry Smith MatScalar *v; 76615091d37SBarry Smith int ierr,mbs=a->mbs,i,*idx,*ii,j,n; 76715091d37SBarry Smith 76815091d37SBarry Smith PetscFunctionBegin; 769b1d4fb26SBarry Smith ierr = VecGetArrayFast(xx,&x);CHKERRQ(ierr); 770b1d4fb26SBarry Smith ierr = VecGetArrayFast(yy,&y);CHKERRQ(ierr); 77115091d37SBarry Smith if (zz != yy) { 772b1d4fb26SBarry Smith ierr = VecGetArrayFast(zz,&z);CHKERRQ(ierr); 77315091d37SBarry Smith } else { 77415091d37SBarry Smith z = y; 77515091d37SBarry Smith } 77615091d37SBarry Smith 77715091d37SBarry Smith idx = a->j; 77815091d37SBarry Smith v = a->a; 77915091d37SBarry Smith ii = a->i; 78015091d37SBarry Smith 78115091d37SBarry Smith for (i=0; i<mbs; i++) { 78215091d37SBarry Smith n = ii[1] - ii[0]; ii++; 78315091d37SBarry Smith sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4]; sum6 = y[5]; 78415091d37SBarry Smith for (j=0; j<n; j++) { 7853b95cb0eSSatish Balay xb = x + 6*(*idx++); 78615091d37SBarry Smith x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; 78715091d37SBarry Smith sum1 += v[0]*x1 + v[6]*x2 + v[12]*x3 + v[18]*x4 + v[24]*x5 + v[30]*x6; 78815091d37SBarry Smith sum2 += v[1]*x1 + v[7]*x2 + v[13]*x3 + v[19]*x4 + v[25]*x5 + v[31]*x6; 78915091d37SBarry Smith sum3 += v[2]*x1 + v[8]*x2 + v[14]*x3 + v[20]*x4 + v[26]*x5 + v[32]*x6; 79015091d37SBarry Smith sum4 += v[3]*x1 + v[9]*x2 + v[15]*x3 + v[21]*x4 + v[27]*x5 + v[33]*x6; 79115091d37SBarry Smith sum5 += v[4]*x1 + v[10]*x2 + v[16]*x3 + v[22]*x4 + v[28]*x5 + v[34]*x6; 79215091d37SBarry Smith sum6 += v[5]*x1 + v[11]*x2 + v[17]*x3 + v[23]*x4 + v[29]*x5 + v[35]*x6; 79315091d37SBarry Smith v += 36; 79415091d37SBarry Smith } 79515091d37SBarry Smith z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; 79615091d37SBarry Smith z += 6; y += 6; 79715091d37SBarry Smith } 798b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(xx,&x);CHKERRQ(ierr); 799b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(yy,&y);CHKERRQ(ierr); 80015091d37SBarry Smith if (zz != yy) { 801b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(zz,&z);CHKERRQ(ierr); 80215091d37SBarry Smith } 803b0a32e0cSBarry Smith PetscLogFlops(72*a->nz); 80415091d37SBarry Smith PetscFunctionReturn(0); 80515091d37SBarry Smith } 8062d61bbb3SSatish Balay 8074a2ae208SSatish Balay #undef __FUNCT__ 8084a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_7" 8092d61bbb3SSatish Balay int MatMultAdd_SeqBAIJ_7(Mat A,Vec xx,Vec yy,Vec zz) 8102d61bbb3SSatish Balay { 8112d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 81287828ca2SBarry Smith PetscScalar *x,*y,*z,*xb,sum1,sum2,sum3,sum4,sum5,sum6,sum7; 81387828ca2SBarry Smith PetscScalar x1,x2,x3,x4,x5,x6,x7; 8143f1db9ecSBarry Smith MatScalar *v; 815e1311b90SBarry Smith int ierr,mbs=a->mbs,i,*idx,*ii,j,n; 8162d61bbb3SSatish Balay 8172d61bbb3SSatish Balay PetscFunctionBegin; 818b1d4fb26SBarry Smith ierr = VecGetArrayFast(xx,&x);CHKERRQ(ierr); 819b1d4fb26SBarry Smith ierr = VecGetArrayFast(yy,&y);CHKERRQ(ierr); 8202e8a6d31SBarry Smith if (zz != yy) { 821b1d4fb26SBarry Smith ierr = VecGetArrayFast(zz,&z);CHKERRQ(ierr); 8222e8a6d31SBarry Smith } else { 8232e8a6d31SBarry Smith z = y; 8242e8a6d31SBarry Smith } 8252d61bbb3SSatish Balay 8262d61bbb3SSatish Balay idx = a->j; 8272d61bbb3SSatish Balay v = a->a; 8282d61bbb3SSatish Balay ii = a->i; 8292d61bbb3SSatish Balay 8302d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 8312d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 8322d61bbb3SSatish Balay sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4]; sum6 = y[5]; sum7 = y[6]; 8332d61bbb3SSatish Balay for (j=0; j<n; j++) { 8342d61bbb3SSatish Balay xb = x + 7*(*idx++); 8352d61bbb3SSatish Balay x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6]; 8362d61bbb3SSatish Balay sum1 += v[0]*x1 + v[7]*x2 + v[14]*x3 + v[21]*x4 + v[28]*x5 + v[35]*x6 + v[42]*x7; 8372d61bbb3SSatish Balay sum2 += v[1]*x1 + v[8]*x2 + v[15]*x3 + v[22]*x4 + v[29]*x5 + v[36]*x6 + v[43]*x7; 8382d61bbb3SSatish Balay sum3 += v[2]*x1 + v[9]*x2 + v[16]*x3 + v[23]*x4 + v[30]*x5 + v[37]*x6 + v[44]*x7; 8392d61bbb3SSatish Balay sum4 += v[3]*x1 + v[10]*x2 + v[17]*x3 + v[24]*x4 + v[31]*x5 + v[38]*x6 + v[45]*x7; 8402d61bbb3SSatish Balay sum5 += v[4]*x1 + v[11]*x2 + v[18]*x3 + v[25]*x4 + v[32]*x5 + v[39]*x6 + v[46]*x7; 8412d61bbb3SSatish Balay sum6 += v[5]*x1 + v[12]*x2 + v[19]*x3 + v[26]*x4 + v[33]*x5 + v[40]*x6 + v[47]*x7; 8422d61bbb3SSatish Balay sum7 += v[6]*x1 + v[13]*x2 + v[20]*x3 + v[27]*x4 + v[34]*x5 + v[41]*x6 + v[48]*x7; 8432d61bbb3SSatish Balay v += 49; 8442d61bbb3SSatish Balay } 8452d61bbb3SSatish Balay z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7; 8462d61bbb3SSatish Balay z += 7; y += 7; 8472d61bbb3SSatish Balay } 848b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(xx,&x);CHKERRQ(ierr); 849b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(yy,&y);CHKERRQ(ierr); 8502e8a6d31SBarry Smith if (zz != yy) { 851b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(zz,&z);CHKERRQ(ierr); 8522e8a6d31SBarry Smith } 853b0a32e0cSBarry Smith PetscLogFlops(98*a->nz); 8542d61bbb3SSatish Balay PetscFunctionReturn(0); 8552d61bbb3SSatish Balay } 856218c64b6SSatish Balay 8574a2ae208SSatish Balay #undef __FUNCT__ 8584a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqBAIJ_N" 8592d61bbb3SSatish Balay int MatMultAdd_SeqBAIJ_N(Mat A,Vec xx,Vec yy,Vec zz) 8602d61bbb3SSatish Balay { 8612d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 862b2e7d493SHong Zhang PetscScalar *x,*z,*xb,*work,*workt,*y; 8633f1db9ecSBarry Smith MatScalar *v; 8642d61bbb3SSatish Balay int mbs=a->mbs,i,*idx,*ii,bs=a->bs,j,n,bs2=a->bs2,ierr; 8653f1db9ecSBarry Smith int ncols,k; 866218c64b6SSatish Balay 8672d61bbb3SSatish Balay PetscFunctionBegin; 868b1d4fb26SBarry Smith ierr = VecGetArrayFast(xx,&x);CHKERRQ(ierr); 869b1d4fb26SBarry Smith ierr = VecGetArrayFast(zz,&z);CHKERRQ(ierr); 870b2e7d493SHong Zhang if (zz != yy) { 871b2e7d493SHong Zhang ierr = VecGetArrayFast(yy,&y);CHKERRQ(ierr); 872b2e7d493SHong Zhang ierr = PetscMemcpy(z,y,yy->n*sizeof(PetscScalar));CHKERRQ(ierr); 873b2e7d493SHong Zhang ierr = VecRestoreArrayFast(yy,&y);CHKERRQ(ierr); 874b2e7d493SHong Zhang } 8752d61bbb3SSatish Balay 8762d61bbb3SSatish Balay idx = a->j; 8772d61bbb3SSatish Balay v = a->a; 8782d61bbb3SSatish Balay ii = a->i; 8792d61bbb3SSatish Balay 8802d61bbb3SSatish Balay 8812d61bbb3SSatish Balay if (!a->mult_work) { 882273d9f13SBarry Smith k = PetscMax(A->m,A->n); 88387828ca2SBarry Smith ierr = PetscMalloc((k+1)*sizeof(PetscScalar),&a->mult_work);CHKERRQ(ierr); 8842d61bbb3SSatish Balay } 8852d61bbb3SSatish Balay work = a->mult_work; 8862d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 8872d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 8882d61bbb3SSatish Balay ncols = n*bs; 8892d61bbb3SSatish Balay workt = work; 8902d61bbb3SSatish Balay for (j=0; j<n; j++) { 8912d61bbb3SSatish Balay xb = x + bs*(*idx++); 8922d61bbb3SSatish Balay for (k=0; k<bs; k++) workt[k] = xb[k]; 8932d61bbb3SSatish Balay workt += bs; 8942d61bbb3SSatish Balay } 8953f1db9ecSBarry Smith Kernel_w_gets_w_plus_Ar_times_v(bs,ncols,work,v,z); 8963f1db9ecSBarry Smith /* LAgemv_("N",&bs,&ncols,&_DOne,v,&bs,work,&_One,&_DOne,z,&_One); */ 8972d61bbb3SSatish Balay v += n*bs2; 8982d61bbb3SSatish Balay z += bs; 8992d61bbb3SSatish Balay } 900b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(xx,&x);CHKERRQ(ierr); 901b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(zz,&z);CHKERRQ(ierr); 902b0a32e0cSBarry Smith PetscLogFlops(2*a->nz*bs2); 9032d61bbb3SSatish Balay PetscFunctionReturn(0); 9042d61bbb3SSatish Balay } 9052d61bbb3SSatish Balay 9064a2ae208SSatish Balay #undef __FUNCT__ 9074a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_SeqBAIJ" 9087c922b88SBarry Smith int MatMultTranspose_SeqBAIJ(Mat A,Vec xx,Vec zz) 9092d61bbb3SSatish Balay { 9102d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 91187828ca2SBarry Smith PetscScalar *xg,*zg,*zb,zero = 0.0; 91287828ca2SBarry Smith PetscScalar *x,*z,*xb,x1,x2,x3,x4,x5,x6,x7; 9133f1db9ecSBarry Smith MatScalar *v; 914f1af5d2fSBarry Smith int mbs=a->mbs,i,*idx,*ii,*ai=a->i,rval; 9152d61bbb3SSatish Balay int bs=a->bs,j,n,bs2=a->bs2,*ib,ierr; 9162d61bbb3SSatish Balay 9172d61bbb3SSatish Balay PetscFunctionBegin; 918f1af5d2fSBarry Smith ierr = VecSet(&zero,zz);CHKERRQ(ierr); 919b1d4fb26SBarry Smith ierr = VecGetArrayFast(xx,&xg);CHKERRQ(ierr); x = xg; 920b1d4fb26SBarry Smith ierr = VecGetArrayFast(zz,&zg);CHKERRQ(ierr); z = zg; 9212d61bbb3SSatish Balay 9222d61bbb3SSatish Balay idx = a->j; 9232d61bbb3SSatish Balay v = a->a; 9242d61bbb3SSatish Balay ii = a->i; 925f1af5d2fSBarry Smith xb = x; 9262d61bbb3SSatish Balay switch (bs) { 9272d61bbb3SSatish Balay case 1: 9282d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 9292d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 930f1af5d2fSBarry Smith x1 = xb[0]; 9312d61bbb3SSatish Balay ib = idx + ai[i]; 9322d61bbb3SSatish Balay for (j=0; j<n; j++) { 9332d61bbb3SSatish Balay rval = ib[j]; 934f1af5d2fSBarry Smith z[rval] += *v * x1; 935f1af5d2fSBarry Smith v++; 9362d61bbb3SSatish Balay } 937f1af5d2fSBarry Smith xb++; 9382d61bbb3SSatish Balay } 9392d61bbb3SSatish Balay break; 9402d61bbb3SSatish Balay case 2: 9412d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 9422d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 943f1af5d2fSBarry Smith x1 = xb[0]; x2 = xb[1]; 9442d61bbb3SSatish Balay ib = idx + ai[i]; 9452d61bbb3SSatish Balay for (j=0; j<n; j++) { 9462d61bbb3SSatish Balay rval = ib[j]*2; 9472d61bbb3SSatish Balay z[rval++] += v[0]*x1 + v[1]*x2; 948f1af5d2fSBarry Smith z[rval] += v[2]*x1 + v[3]*x2; 9492d61bbb3SSatish Balay v += 4; 9502d61bbb3SSatish Balay } 951f1af5d2fSBarry Smith xb += 2; 9522d61bbb3SSatish Balay } 9532d61bbb3SSatish Balay break; 9542d61bbb3SSatish Balay case 3: 9552d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 9562d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 957f1af5d2fSBarry Smith x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; 9582d61bbb3SSatish Balay ib = idx + ai[i]; 9592d61bbb3SSatish Balay for (j=0; j<n; j++) { 9602d61bbb3SSatish Balay rval = ib[j]*3; 9612d61bbb3SSatish Balay z[rval++] += v[0]*x1 + v[1]*x2 + v[2]*x3; 9622d61bbb3SSatish Balay z[rval++] += v[3]*x1 + v[4]*x2 + v[5]*x3; 963f1af5d2fSBarry Smith z[rval] += v[6]*x1 + v[7]*x2 + v[8]*x3; 9642d61bbb3SSatish Balay v += 9; 9652d61bbb3SSatish Balay } 966f1af5d2fSBarry Smith xb += 3; 9672d61bbb3SSatish Balay } 9682d61bbb3SSatish Balay break; 9692d61bbb3SSatish Balay case 4: 9702d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 9712d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 972f1af5d2fSBarry Smith x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; 9732d61bbb3SSatish Balay ib = idx + ai[i]; 9742d61bbb3SSatish Balay for (j=0; j<n; j++) { 9752d61bbb3SSatish Balay rval = ib[j]*4; 9762d61bbb3SSatish Balay z[rval++] += v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4; 9772d61bbb3SSatish Balay z[rval++] += v[4]*x1 + v[5]*x2 + v[6]*x3 + v[7]*x4; 9782d61bbb3SSatish Balay z[rval++] += v[8]*x1 + v[9]*x2 + v[10]*x3 + v[11]*x4; 979f1af5d2fSBarry Smith z[rval] += v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4; 9802d61bbb3SSatish Balay v += 16; 9812d61bbb3SSatish Balay } 982f1af5d2fSBarry Smith xb += 4; 9832d61bbb3SSatish Balay } 9842d61bbb3SSatish Balay break; 9852d61bbb3SSatish Balay case 5: 9862d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 9872d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 988f1af5d2fSBarry Smith x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; 9892d61bbb3SSatish Balay x4 = xb[3]; x5 = xb[4]; 9902d61bbb3SSatish Balay ib = idx + ai[i]; 9912d61bbb3SSatish Balay for (j=0; j<n; j++) { 9922d61bbb3SSatish Balay rval = ib[j]*5; 9932d61bbb3SSatish Balay z[rval++] += v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4 + v[4]*x5; 9942d61bbb3SSatish Balay z[rval++] += v[5]*x1 + v[6]*x2 + v[7]*x3 + v[8]*x4 + v[9]*x5; 9952d61bbb3SSatish Balay z[rval++] += v[10]*x1 + v[11]*x2 + v[12]*x3 + v[13]*x4 + v[14]*x5; 9962d61bbb3SSatish Balay z[rval++] += v[15]*x1 + v[16]*x2 + v[17]*x3 + v[18]*x4 + v[19]*x5; 997f1af5d2fSBarry Smith z[rval] += v[20]*x1 + v[21]*x2 + v[22]*x3 + v[23]*x4 + v[24]*x5; 9982d61bbb3SSatish Balay v += 25; 9992d61bbb3SSatish Balay } 1000f1af5d2fSBarry Smith xb += 5; 10012d61bbb3SSatish Balay } 10022d61bbb3SSatish Balay break; 1003f1af5d2fSBarry Smith case 6: 1004f1af5d2fSBarry Smith for (i=0; i<mbs; i++) { 1005f1af5d2fSBarry Smith n = ii[1] - ii[0]; ii++; 1006f1af5d2fSBarry Smith x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; 1007f1af5d2fSBarry Smith x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; 1008f1af5d2fSBarry Smith ib = idx + ai[i]; 1009f1af5d2fSBarry Smith for (j=0; j<n; j++) { 1010f1af5d2fSBarry Smith rval = ib[j]*6; 1011f1af5d2fSBarry Smith z[rval++] += v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4 + v[4]*x5 + v[5]*x6; 1012f1af5d2fSBarry Smith z[rval++] += v[6]*x1 + v[7]*x2 + v[8]*x3 + v[9]*x4 + v[10]*x5 + v[11]*x6; 1013f1af5d2fSBarry Smith z[rval++] += v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4 + v[16]*x5 + v[17]*x6; 1014f1af5d2fSBarry Smith z[rval++] += v[18]*x1 + v[19]*x2 + v[20]*x3 + v[21]*x4 + v[22]*x5 + v[23]*x6; 1015f1af5d2fSBarry Smith z[rval++] += v[24]*x1 + v[25]*x2 + v[26]*x3 + v[27]*x4 + v[28]*x5 + v[29]*x6; 1016f1af5d2fSBarry Smith z[rval] += v[30]*x1 + v[31]*x2 + v[32]*x3 + v[33]*x4 + v[34]*x5 + v[35]*x6; 1017f1af5d2fSBarry Smith v += 36; 1018f1af5d2fSBarry Smith } 1019f1af5d2fSBarry Smith xb += 6; 1020f1af5d2fSBarry Smith } 1021f1af5d2fSBarry Smith break; 1022f1af5d2fSBarry Smith case 7: 1023f1af5d2fSBarry Smith for (i=0; i<mbs; i++) { 1024f1af5d2fSBarry Smith n = ii[1] - ii[0]; ii++; 1025f1af5d2fSBarry Smith x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; 1026f1af5d2fSBarry Smith x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6]; 1027f1af5d2fSBarry Smith ib = idx + ai[i]; 1028f1af5d2fSBarry Smith for (j=0; j<n; j++) { 1029f1af5d2fSBarry Smith rval = ib[j]*7; 1030f1af5d2fSBarry Smith z[rval++] += v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4 + v[4]*x5 + v[5]*x6 + v[6]*x7; 1031f1af5d2fSBarry Smith z[rval++] += v[7]*x1 + v[8]*x2 + v[9]*x3 + v[10]*x4 + v[11]*x5 + v[12]*x6 + v[13]*x7; 1032f1af5d2fSBarry Smith z[rval++] += v[14]*x1 + v[15]*x2 + v[16]*x3 + v[17]*x4 + v[18]*x5 + v[19]*x6 + v[20]*x7; 1033f1af5d2fSBarry Smith z[rval++] += v[21]*x1 + v[22]*x2 + v[23]*x3 + v[24]*x4 + v[25]*x5 + v[26]*x6 + v[27]*x7; 1034f1af5d2fSBarry Smith z[rval++] += v[28]*x1 + v[29]*x2 + v[30]*x3 + v[31]*x4 + v[32]*x5 + v[33]*x6 + v[34]*x7; 1035f1af5d2fSBarry Smith z[rval++] += v[35]*x1 + v[36]*x2 + v[37]*x3 + v[38]*x4 + v[39]*x5 + v[40]*x6 + v[41]*x7; 1036f1af5d2fSBarry Smith z[rval] += v[42]*x1 + v[43]*x2 + v[44]*x3 + v[45]*x4 + v[46]*x5 + v[47]*x6 + v[48]*x7; 1037f1af5d2fSBarry Smith v += 49; 1038f1af5d2fSBarry Smith } 1039f1af5d2fSBarry Smith xb += 7; 1040f1af5d2fSBarry Smith } 1041f1af5d2fSBarry Smith break; 1042f1af5d2fSBarry Smith default: { /* block sizes larger then 7 by 7 are handled by BLAS */ 10433f1db9ecSBarry Smith int ncols,k; 104487828ca2SBarry Smith PetscScalar *work,*workt; 10453f1db9ecSBarry Smith 10462d61bbb3SSatish Balay if (!a->mult_work) { 1047273d9f13SBarry Smith k = PetscMax(A->m,A->n); 104887828ca2SBarry Smith ierr = PetscMalloc((k+1)*sizeof(PetscScalar),&a->mult_work);CHKERRQ(ierr); 10492d61bbb3SSatish Balay } 10502d61bbb3SSatish Balay work = a->mult_work; 10512d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 10522d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 10532d61bbb3SSatish Balay ncols = n*bs; 105487828ca2SBarry Smith ierr = PetscMemzero(work,ncols*sizeof(PetscScalar));CHKERRQ(ierr); 1055d824769bSBarry Smith Kernel_w_gets_w_plus_trans_Ar_times_v(bs,ncols,x,v,work); 10563f1db9ecSBarry Smith /* LAgemv_("T",&bs,&ncols,&_DOne,v,&bs,x,&_One,&_DOne,work,&_One); */ 10572d61bbb3SSatish Balay v += n*bs2; 10582d61bbb3SSatish Balay x += bs; 10592d61bbb3SSatish Balay workt = work; 10602d61bbb3SSatish Balay for (j=0; j<n; j++) { 10612d61bbb3SSatish Balay zb = z + bs*(*idx++); 10622d61bbb3SSatish Balay for (k=0; k<bs; k++) zb[k] += workt[k] ; 10632d61bbb3SSatish Balay workt += bs; 10642d61bbb3SSatish Balay } 10652d61bbb3SSatish Balay } 10662d61bbb3SSatish Balay } 10672d61bbb3SSatish Balay } 1068b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(xx,&xg);CHKERRQ(ierr); 1069b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(zz,&zg);CHKERRQ(ierr); 1070b0a32e0cSBarry Smith PetscLogFlops(2*a->nz*a->bs2 - A->n); 10712d61bbb3SSatish Balay PetscFunctionReturn(0); 10722d61bbb3SSatish Balay } 10732d61bbb3SSatish Balay 10744a2ae208SSatish Balay #undef __FUNCT__ 10754a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_SeqBAIJ" 10767c922b88SBarry Smith int MatMultTransposeAdd_SeqBAIJ(Mat A,Vec xx,Vec yy,Vec zz) 10772d61bbb3SSatish Balay 10782d61bbb3SSatish Balay { 10792d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 108087828ca2SBarry Smith PetscScalar *xg,*zg,*zb,*x,*z,*xb,x1,x2,x3,x4,x5; 10813f1db9ecSBarry Smith MatScalar *v; 1082f1af5d2fSBarry Smith int mbs=a->mbs,i,*idx,*ii,*ai=a->i,rval,bs=a->bs,j,n,bs2=a->bs2,*ib,ierr; 10832d61bbb3SSatish Balay 10842d61bbb3SSatish Balay PetscFunctionBegin; 1085ef66eb69SBarry Smith if (yy != zz) { ierr = VecCopy(yy,zz);CHKERRQ(ierr); } 1086b1d4fb26SBarry Smith ierr = VecGetArrayFast(xx,&xg);CHKERRQ(ierr); x = xg; 1087b1d4fb26SBarry Smith ierr = VecGetArrayFast(zz,&zg);CHKERRQ(ierr); z = zg; 10882d61bbb3SSatish Balay 10892d61bbb3SSatish Balay 10902d61bbb3SSatish Balay idx = a->j; 10912d61bbb3SSatish Balay v = a->a; 10922d61bbb3SSatish Balay ii = a->i; 1093f1af5d2fSBarry Smith xb = x; 10942d61bbb3SSatish Balay 10952d61bbb3SSatish Balay switch (bs) { 10962d61bbb3SSatish Balay case 1: 10972d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 10982d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 1099f1af5d2fSBarry Smith x1 = xb[0]; 11002d61bbb3SSatish Balay ib = idx + ai[i]; 11012d61bbb3SSatish Balay for (j=0; j<n; j++) { 11022d61bbb3SSatish Balay rval = ib[j]; 1103f1af5d2fSBarry Smith z[rval] += *v * x1; 1104f1af5d2fSBarry Smith v++; 11052d61bbb3SSatish Balay } 1106f1af5d2fSBarry Smith xb++; 11072d61bbb3SSatish Balay } 11082d61bbb3SSatish Balay break; 11092d61bbb3SSatish Balay case 2: 11102d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 11112d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 1112f1af5d2fSBarry Smith x1 = xb[0]; x2 = xb[1]; 11132d61bbb3SSatish Balay ib = idx + ai[i]; 11142d61bbb3SSatish Balay for (j=0; j<n; j++) { 11152d61bbb3SSatish Balay rval = ib[j]*2; 11162d61bbb3SSatish Balay z[rval++] += v[0]*x1 + v[1]*x2; 11172d61bbb3SSatish Balay z[rval++] += v[2]*x1 + v[3]*x2; 11182d61bbb3SSatish Balay v += 4; 11192d61bbb3SSatish Balay } 1120f1af5d2fSBarry Smith xb += 2; 11212d61bbb3SSatish Balay } 11222d61bbb3SSatish Balay break; 11232d61bbb3SSatish Balay case 3: 11242d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 11252d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 1126f1af5d2fSBarry Smith x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; 11272d61bbb3SSatish Balay ib = idx + ai[i]; 11282d61bbb3SSatish Balay for (j=0; j<n; j++) { 11292d61bbb3SSatish Balay rval = ib[j]*3; 11302d61bbb3SSatish Balay z[rval++] += v[0]*x1 + v[1]*x2 + v[2]*x3; 11312d61bbb3SSatish Balay z[rval++] += v[3]*x1 + v[4]*x2 + v[5]*x3; 11322d61bbb3SSatish Balay z[rval++] += v[6]*x1 + v[7]*x2 + v[8]*x3; 11332d61bbb3SSatish Balay v += 9; 11342d61bbb3SSatish Balay } 1135f1af5d2fSBarry Smith xb += 3; 11362d61bbb3SSatish Balay } 11372d61bbb3SSatish Balay break; 11382d61bbb3SSatish Balay case 4: 11392d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 11402d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 1141f1af5d2fSBarry Smith x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; 11422d61bbb3SSatish Balay ib = idx + ai[i]; 11432d61bbb3SSatish Balay for (j=0; j<n; j++) { 11442d61bbb3SSatish Balay rval = ib[j]*4; 11452d61bbb3SSatish Balay z[rval++] += v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4; 11462d61bbb3SSatish Balay z[rval++] += v[4]*x1 + v[5]*x2 + v[6]*x3 + v[7]*x4; 11472d61bbb3SSatish Balay z[rval++] += v[8]*x1 + v[9]*x2 + v[10]*x3 + v[11]*x4; 11482d61bbb3SSatish Balay z[rval++] += v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4; 11492d61bbb3SSatish Balay v += 16; 11502d61bbb3SSatish Balay } 1151f1af5d2fSBarry Smith xb += 4; 11522d61bbb3SSatish Balay } 11532d61bbb3SSatish Balay break; 11542d61bbb3SSatish Balay case 5: 11552d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 11562d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 1157f1af5d2fSBarry Smith x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; 11582d61bbb3SSatish Balay x4 = xb[3]; x5 = xb[4]; 11592d61bbb3SSatish Balay ib = idx + ai[i]; 11602d61bbb3SSatish Balay for (j=0; j<n; j++) { 11612d61bbb3SSatish Balay rval = ib[j]*5; 11622d61bbb3SSatish Balay z[rval++] += v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4 + v[4]*x5; 11632d61bbb3SSatish Balay z[rval++] += v[5]*x1 + v[6]*x2 + v[7]*x3 + v[8]*x4 + v[9]*x5; 11642d61bbb3SSatish Balay z[rval++] += v[10]*x1 + v[11]*x2 + v[12]*x3 + v[13]*x4 + v[14]*x5; 11652d61bbb3SSatish Balay z[rval++] += v[15]*x1 + v[16]*x2 + v[17]*x3 + v[18]*x4 + v[19]*x5; 11662d61bbb3SSatish Balay z[rval++] += v[20]*x1 + v[21]*x2 + v[22]*x3 + v[23]*x4 + v[24]*x5; 11672d61bbb3SSatish Balay v += 25; 11682d61bbb3SSatish Balay } 1169f1af5d2fSBarry Smith xb += 5; 11702d61bbb3SSatish Balay } 11712d61bbb3SSatish Balay break; 1172f1af5d2fSBarry Smith default: { /* block sizes larger then 5 by 5 are handled by BLAS */ 11733f1db9ecSBarry Smith int ncols,k; 117487828ca2SBarry Smith PetscScalar *work,*workt; 11753f1db9ecSBarry Smith 11762d61bbb3SSatish Balay if (!a->mult_work) { 1177273d9f13SBarry Smith k = PetscMax(A->m,A->n); 117887828ca2SBarry Smith ierr = PetscMalloc((k+1)*sizeof(PetscScalar),&a->mult_work);CHKERRQ(ierr); 11792d61bbb3SSatish Balay } 11802d61bbb3SSatish Balay work = a->mult_work; 11812d61bbb3SSatish Balay for (i=0; i<mbs; i++) { 11822d61bbb3SSatish Balay n = ii[1] - ii[0]; ii++; 11832d61bbb3SSatish Balay ncols = n*bs; 118487828ca2SBarry Smith ierr = PetscMemzero(work,ncols*sizeof(PetscScalar));CHKERRQ(ierr); 11853f1db9ecSBarry Smith Kernel_w_gets_w_plus_trans_Ar_times_v(bs,ncols,x,v,work); 11863f1db9ecSBarry Smith /* LAgemv_("T",&bs,&ncols,&_DOne,v,&bs,x,&_One,&_DOne,work,&_One); */ 11872d61bbb3SSatish Balay v += n*bs2; 11882d61bbb3SSatish Balay x += bs; 11892d61bbb3SSatish Balay workt = work; 11902d61bbb3SSatish Balay for (j=0; j<n; j++) { 11912d61bbb3SSatish Balay zb = z + bs*(*idx++); 11922d61bbb3SSatish Balay for (k=0; k<bs; k++) zb[k] += workt[k] ; 11932d61bbb3SSatish Balay workt += bs; 11942d61bbb3SSatish Balay } 11952d61bbb3SSatish Balay } 11962d61bbb3SSatish Balay } 11972d61bbb3SSatish Balay } 1198b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(xx,&xg);CHKERRQ(ierr); 1199b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(zz,&zg);CHKERRQ(ierr); 1200b0a32e0cSBarry Smith PetscLogFlops(2*a->nz*a->bs2); 12012d61bbb3SSatish Balay PetscFunctionReturn(0); 12022d61bbb3SSatish Balay } 12032d61bbb3SSatish Balay 12044a2ae208SSatish Balay #undef __FUNCT__ 12054a2ae208SSatish Balay #define __FUNCT__ "MatScale_SeqBAIJ" 1206268466fbSBarry Smith int MatScale_SeqBAIJ(const PetscScalar *alpha,Mat inA) 12072d61bbb3SSatish Balay { 12082d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)inA->data; 120987052302SDinesh Kaushik int totalnz = a->bs2*a->nz; 12103eda8832SBarry Smith #if defined(PETSC_USE_MAT_SINGLE) 12113eda8832SBarry Smith int i; 121287052302SDinesh Kaushik #else 121387052302SDinesh Kaushik int one = 1; 12143eda8832SBarry Smith #endif 12152d61bbb3SSatish Balay 12162d61bbb3SSatish Balay PetscFunctionBegin; 12173eda8832SBarry Smith #if defined(PETSC_USE_MAT_SINGLE) 12183eda8832SBarry Smith for (i=0; i<totalnz; i++) a->a[i] *= *alpha; 12193eda8832SBarry Smith #else 1220268466fbSBarry Smith BLscal_(&totalnz,(PetscScalar*)alpha,a->a,&one); 12213eda8832SBarry Smith #endif 1222b0a32e0cSBarry Smith PetscLogFlops(totalnz); 12232d61bbb3SSatish Balay PetscFunctionReturn(0); 12242d61bbb3SSatish Balay } 12252d61bbb3SSatish Balay 12264a2ae208SSatish Balay #undef __FUNCT__ 12274a2ae208SSatish Balay #define __FUNCT__ "MatNorm_SeqBAIJ" 1228329f5518SBarry Smith int MatNorm_SeqBAIJ(Mat A,NormType type,PetscReal *norm) 12292d61bbb3SSatish Balay { 12302d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 12313f1db9ecSBarry Smith MatScalar *v = a->a; 1232329f5518SBarry Smith PetscReal sum = 0.0; 12330e90e235SBarry Smith int i,j,k,bs = a->bs,nz=a->nz,bs2=a->bs2,k1; 12342d61bbb3SSatish Balay 12352d61bbb3SSatish Balay PetscFunctionBegin; 12362d61bbb3SSatish Balay if (type == NORM_FROBENIUS) { 12372d61bbb3SSatish Balay for (i=0; i< bs2*nz; i++) { 1238aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 1239329f5518SBarry Smith sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 12402d61bbb3SSatish Balay #else 12412d61bbb3SSatish Balay sum += (*v)*(*v); v++; 12422d61bbb3SSatish Balay #endif 12432d61bbb3SSatish Balay } 12442d61bbb3SSatish Balay *norm = sqrt(sum); 1245596552b5SBarry Smith } else if (type == NORM_INFINITY) { /* maximum row sum */ 1246596552b5SBarry Smith *norm = 0.0; 1247596552b5SBarry Smith for (k=0; k<bs; k++) { 124874f84c7bSSatish Balay for (j=0; j<a->mbs; j++) { 1249596552b5SBarry Smith v = a->a + bs2*a->i[j] + k; 1250596552b5SBarry Smith sum = 0.0; 1251596552b5SBarry Smith for (i=0; i<a->i[j+1]-a->i[j]; i++) { 12520e90e235SBarry Smith for (k1=0; k1<bs; k1++){ 1253596552b5SBarry Smith sum += PetscAbsScalar(*v); 1254596552b5SBarry Smith v += bs; 12552d61bbb3SSatish Balay } 12560e90e235SBarry Smith } 1257596552b5SBarry Smith if (sum > *norm) *norm = sum; 1258596552b5SBarry Smith } 1259596552b5SBarry Smith } 1260596552b5SBarry Smith } else { 126129bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"No support for this norm yet"); 12622d61bbb3SSatish Balay } 12632d61bbb3SSatish Balay PetscFunctionReturn(0); 12642d61bbb3SSatish Balay } 12652d61bbb3SSatish Balay 12662d61bbb3SSatish Balay 12674a2ae208SSatish Balay #undef __FUNCT__ 12684a2ae208SSatish Balay #define __FUNCT__ "MatEqual_SeqBAIJ" 12692d61bbb3SSatish Balay int MatEqual_SeqBAIJ(Mat A,Mat B,PetscTruth* flg) 12702d61bbb3SSatish Balay { 12712d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data,*b = (Mat_SeqBAIJ *)B->data; 12720f5bd95cSBarry Smith int ierr; 12732d61bbb3SSatish Balay 12742d61bbb3SSatish Balay PetscFunctionBegin; 12752d61bbb3SSatish Balay /* If the matrix/block dimensions are not equal, or no of nonzeros or shift */ 1276273d9f13SBarry Smith if ((A->m != B->m) || (A->n != B->n) || (a->bs != b->bs)|| (a->nz != b->nz)) { 1277273d9f13SBarry Smith *flg = PETSC_FALSE; 1278273d9f13SBarry Smith PetscFunctionReturn(0); 12792d61bbb3SSatish Balay } 12802d61bbb3SSatish Balay 12812d61bbb3SSatish Balay /* if the a->i are the same */ 12820f5bd95cSBarry Smith ierr = PetscMemcmp(a->i,b->i,(a->mbs+1)*sizeof(int),flg);CHKERRQ(ierr); 12830f5bd95cSBarry Smith if (*flg == PETSC_FALSE) { 12840f5bd95cSBarry Smith PetscFunctionReturn(0); 12852d61bbb3SSatish Balay } 12862d61bbb3SSatish Balay 12872d61bbb3SSatish Balay /* if a->j are the same */ 12880f5bd95cSBarry Smith ierr = PetscMemcmp(a->j,b->j,(a->nz)*sizeof(int),flg);CHKERRQ(ierr); 12890f5bd95cSBarry Smith if (*flg == PETSC_FALSE) { 12900f5bd95cSBarry Smith PetscFunctionReturn(0); 12912d61bbb3SSatish Balay } 12922d61bbb3SSatish Balay /* if a->a are the same */ 129387828ca2SBarry Smith ierr = PetscMemcmp(a->a,b->a,(a->nz)*(a->bs)*(a->bs)*sizeof(PetscScalar),flg);CHKERRQ(ierr); 12942d61bbb3SSatish Balay PetscFunctionReturn(0); 12952d61bbb3SSatish Balay 12962d61bbb3SSatish Balay } 12972d61bbb3SSatish Balay 12984a2ae208SSatish Balay #undef __FUNCT__ 12994a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_SeqBAIJ" 13002d61bbb3SSatish Balay int MatGetDiagonal_SeqBAIJ(Mat A,Vec v) 13012d61bbb3SSatish Balay { 13022d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1303e1311b90SBarry Smith int ierr,i,j,k,n,row,bs,*ai,*aj,ambs,bs2; 130487828ca2SBarry Smith PetscScalar *x,zero = 0.0; 13053f1db9ecSBarry Smith MatScalar *aa,*aa_j; 13062d61bbb3SSatish Balay 13072d61bbb3SSatish Balay PetscFunctionBegin; 130829bbc08cSBarry Smith if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 13092d61bbb3SSatish Balay bs = a->bs; 13102d61bbb3SSatish Balay aa = a->a; 13112d61bbb3SSatish Balay ai = a->i; 13122d61bbb3SSatish Balay aj = a->j; 13132d61bbb3SSatish Balay ambs = a->mbs; 13142d61bbb3SSatish Balay bs2 = a->bs2; 13152d61bbb3SSatish Balay 1316e1311b90SBarry Smith ierr = VecSet(&zero,v);CHKERRQ(ierr); 1317b1d4fb26SBarry Smith ierr = VecGetArrayFast(v,&x);CHKERRQ(ierr); 1318e1311b90SBarry Smith ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr); 1319273d9f13SBarry Smith if (n != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector"); 13202d61bbb3SSatish Balay for (i=0; i<ambs; i++) { 13212d61bbb3SSatish Balay for (j=ai[i]; j<ai[i+1]; j++) { 13222d61bbb3SSatish Balay if (aj[j] == i) { 13232d61bbb3SSatish Balay row = i*bs; 13242d61bbb3SSatish Balay aa_j = aa+j*bs2; 13252d61bbb3SSatish Balay for (k=0; k<bs2; k+=(bs+1),row++) x[row] = aa_j[k]; 13262d61bbb3SSatish Balay break; 13272d61bbb3SSatish Balay } 13282d61bbb3SSatish Balay } 13292d61bbb3SSatish Balay } 1330b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(v,&x);CHKERRQ(ierr); 13312d61bbb3SSatish Balay PetscFunctionReturn(0); 13322d61bbb3SSatish Balay } 13332d61bbb3SSatish Balay 13344a2ae208SSatish Balay #undef __FUNCT__ 13354a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_SeqBAIJ" 13362d61bbb3SSatish Balay int MatDiagonalScale_SeqBAIJ(Mat A,Vec ll,Vec rr) 13372d61bbb3SSatish Balay { 13382d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 133987828ca2SBarry Smith PetscScalar *l,*r,x,*li,*ri; 13403f1db9ecSBarry Smith MatScalar *aa,*v; 1341e1311b90SBarry Smith int ierr,i,j,k,lm,rn,M,m,n,*ai,*aj,mbs,tmp,bs,bs2; 13422d61bbb3SSatish Balay 13432d61bbb3SSatish Balay PetscFunctionBegin; 13442d61bbb3SSatish Balay ai = a->i; 13452d61bbb3SSatish Balay aj = a->j; 13462d61bbb3SSatish Balay aa = a->a; 1347273d9f13SBarry Smith m = A->m; 1348273d9f13SBarry Smith n = A->n; 13492d61bbb3SSatish Balay bs = a->bs; 13502d61bbb3SSatish Balay mbs = a->mbs; 13512d61bbb3SSatish Balay bs2 = a->bs2; 13522d61bbb3SSatish Balay if (ll) { 1353b1d4fb26SBarry Smith ierr = VecGetArrayFast(ll,&l);CHKERRQ(ierr); 1354e1311b90SBarry Smith ierr = VecGetLocalSize(ll,&lm);CHKERRQ(ierr); 135529bbc08cSBarry Smith if (lm != m) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vector wrong length"); 13562d61bbb3SSatish Balay for (i=0; i<mbs; i++) { /* for each block row */ 13572d61bbb3SSatish Balay M = ai[i+1] - ai[i]; 13582d61bbb3SSatish Balay li = l + i*bs; 13592d61bbb3SSatish Balay v = aa + bs2*ai[i]; 13602d61bbb3SSatish Balay for (j=0; j<M; j++) { /* for each block */ 13612d61bbb3SSatish Balay for (k=0; k<bs2; k++) { 13622d61bbb3SSatish Balay (*v++) *= li[k%bs]; 13632d61bbb3SSatish Balay } 13642d61bbb3SSatish Balay } 13652d61bbb3SSatish Balay } 1366b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(ll,&l);CHKERRQ(ierr); 1367b0a32e0cSBarry Smith PetscLogFlops(a->nz); 13682d61bbb3SSatish Balay } 13692d61bbb3SSatish Balay 13702d61bbb3SSatish Balay if (rr) { 1371b1d4fb26SBarry Smith ierr = VecGetArrayFast(rr,&r);CHKERRQ(ierr); 1372e1311b90SBarry Smith ierr = VecGetLocalSize(rr,&rn);CHKERRQ(ierr); 137329bbc08cSBarry Smith if (rn != n) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vector wrong length"); 13742d61bbb3SSatish Balay for (i=0; i<mbs; i++) { /* for each block row */ 13752d61bbb3SSatish Balay M = ai[i+1] - ai[i]; 13762d61bbb3SSatish Balay v = aa + bs2*ai[i]; 13772d61bbb3SSatish Balay for (j=0; j<M; j++) { /* for each block */ 13782d61bbb3SSatish Balay ri = r + bs*aj[ai[i]+j]; 13792d61bbb3SSatish Balay for (k=0; k<bs; k++) { 13802d61bbb3SSatish Balay x = ri[k]; 13812d61bbb3SSatish Balay for (tmp=0; tmp<bs; tmp++) (*v++) *= x; 13822d61bbb3SSatish Balay } 13832d61bbb3SSatish Balay } 13842d61bbb3SSatish Balay } 1385b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(rr,&r);CHKERRQ(ierr); 1386b0a32e0cSBarry Smith PetscLogFlops(a->nz); 13872d61bbb3SSatish Balay } 13882d61bbb3SSatish Balay PetscFunctionReturn(0); 13892d61bbb3SSatish Balay } 13902d61bbb3SSatish Balay 13912d61bbb3SSatish Balay 13924a2ae208SSatish Balay #undef __FUNCT__ 13934a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_SeqBAIJ" 13942d61bbb3SSatish Balay int MatGetInfo_SeqBAIJ(Mat A,MatInfoType flag,MatInfo *info) 13952d61bbb3SSatish Balay { 13962d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 13972d61bbb3SSatish Balay 13982d61bbb3SSatish Balay PetscFunctionBegin; 1399273d9f13SBarry Smith info->rows_global = (double)A->m; 1400273d9f13SBarry Smith info->columns_global = (double)A->n; 1401273d9f13SBarry Smith info->rows_local = (double)A->m; 1402273d9f13SBarry Smith info->columns_local = (double)A->n; 14032d61bbb3SSatish Balay info->block_size = a->bs2; 14042d61bbb3SSatish Balay info->nz_allocated = a->maxnz; 14052d61bbb3SSatish Balay info->nz_used = a->bs2*a->nz; 14062d61bbb3SSatish Balay info->nz_unneeded = (double)(info->nz_allocated - info->nz_used); 14072d61bbb3SSatish Balay info->assemblies = A->num_ass; 14082d61bbb3SSatish Balay info->mallocs = a->reallocs; 14092d61bbb3SSatish Balay info->memory = A->mem; 14102d61bbb3SSatish Balay if (A->factor) { 14112d61bbb3SSatish Balay info->fill_ratio_given = A->info.fill_ratio_given; 14122d61bbb3SSatish Balay info->fill_ratio_needed = A->info.fill_ratio_needed; 14132d61bbb3SSatish Balay info->factor_mallocs = A->info.factor_mallocs; 14142d61bbb3SSatish Balay } else { 14152d61bbb3SSatish Balay info->fill_ratio_given = 0; 14162d61bbb3SSatish Balay info->fill_ratio_needed = 0; 14172d61bbb3SSatish Balay info->factor_mallocs = 0; 14182d61bbb3SSatish Balay } 14192d61bbb3SSatish Balay PetscFunctionReturn(0); 14202d61bbb3SSatish Balay } 14212d61bbb3SSatish Balay 14222d61bbb3SSatish Balay 14234a2ae208SSatish Balay #undef __FUNCT__ 14244a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_SeqBAIJ" 14252d61bbb3SSatish Balay int MatZeroEntries_SeqBAIJ(Mat A) 14262d61bbb3SSatish Balay { 14272d61bbb3SSatish Balay Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1428549d3d68SSatish Balay int ierr; 14292d61bbb3SSatish Balay 14302d61bbb3SSatish Balay PetscFunctionBegin; 1431549d3d68SSatish Balay ierr = PetscMemzero(a->a,a->bs2*a->i[a->mbs]*sizeof(MatScalar));CHKERRQ(ierr); 14322d61bbb3SSatish Balay PetscFunctionReturn(0); 14332d61bbb3SSatish Balay } 1434