1be1d678aSKris Buschelman #define PETSCMAT_DLL 2be1d678aSKris Buschelman 3d50806bdSBarry Smith /* 42499ec78SHong Zhang Defines matrix-matrix product routines for pairs of SeqAIJ matrices 5d50806bdSBarry Smith C = A * B 6d50806bdSBarry Smith */ 7d50806bdSBarry Smith 87c4f633dSBarry Smith #include "../src/mat/impls/aij/seq/aij.h" /*I "petscmat.h" I*/ 97c4f633dSBarry Smith #include "../src/mat/utils/freespace.h" 10be0fcf8dSHong Zhang #include "petscbt.h" 117c4f633dSBarry Smith #include "../src/mat/impls/dense/seq/dense.h" /*I "petscmat.h" I*/ 126284ec50SHong Zhang 13ec01deb9SMatthew Knepley EXTERN_C_BEGIN 146284ec50SHong Zhang #undef __FUNCT__ 155c66b693SKris Buschelman #define __FUNCT__ "MatMatMult_SeqAIJ_SeqAIJ" 1638baddfdSBarry Smith PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 1738baddfdSBarry Smith { 18dfbe8321SBarry Smith PetscErrorCode ierr; 195c66b693SKris Buschelman 205c66b693SKris Buschelman PetscFunctionBegin; 2126be0446SHong Zhang if (scall == MAT_INITIAL_MATRIX){ 2226be0446SHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr); 2326be0446SHong Zhang } 2426be0446SHong Zhang ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr); 255c66b693SKris Buschelman PetscFunctionReturn(0); 265c66b693SKris Buschelman } 27ec01deb9SMatthew Knepley EXTERN_C_END 281c24bd37SHong Zhang 2926be0446SHong Zhang #undef __FUNCT__ 3026be0446SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ" 31dfbe8321SBarry Smith PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C) 3258c24d83SHong Zhang { 33dfbe8321SBarry Smith PetscErrorCode ierr; 34a1a86e44SBarry Smith PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 3558c24d83SHong Zhang Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c; 3638baddfdSBarry Smith PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci,*cj; 37d0f46423SBarry Smith PetscInt am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N; 3838baddfdSBarry Smith PetscInt i,j,anzi,brow,bnzj,cnzi,nlnk,*lnk,nspacedouble=0; 3958c24d83SHong Zhang MatScalar *ca; 40be0fcf8dSHong Zhang PetscBT lnkbt; 4158c24d83SHong Zhang 4258c24d83SHong Zhang PetscFunctionBegin; 4358c24d83SHong Zhang /* Set up */ 4458c24d83SHong Zhang /* Allocate ci array, arrays for fill computation and */ 4558c24d83SHong Zhang /* free space for accumulating nonzero column info */ 4638baddfdSBarry Smith ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 4758c24d83SHong Zhang ci[0] = 0; 4858c24d83SHong Zhang 49be0fcf8dSHong Zhang /* create and initialize a linked list */ 50be0fcf8dSHong Zhang nlnk = bn+1; 51be0fcf8dSHong Zhang ierr = PetscLLCreate(bn,bn,nlnk,lnk,lnkbt);CHKERRQ(ierr); 5258c24d83SHong Zhang 53c5db241fSHong Zhang /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */ 54a1a86e44SBarry Smith ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr); 5558c24d83SHong Zhang current_space = free_space; 5658c24d83SHong Zhang 5758c24d83SHong Zhang /* Determine symbolic info for each row of the product: */ 5858c24d83SHong Zhang for (i=0;i<am;i++) { 5958c24d83SHong Zhang anzi = ai[i+1] - ai[i]; 6058c24d83SHong Zhang cnzi = 0; 612d09714cSHong Zhang j = anzi; 622d09714cSHong Zhang aj = a->j + ai[i]; 632d09714cSHong Zhang while (j){/* assume cols are almost in increasing order, starting from its end saves computation */ 642d09714cSHong Zhang j--; 652d09714cSHong Zhang brow = *(aj + j); 6658c24d83SHong Zhang bnzj = bi[brow+1] - bi[brow]; 6758c24d83SHong Zhang bjj = bj + bi[brow]; 681c239cc6SHong Zhang /* add non-zero cols of B into the sorted linked list lnk */ 69be0fcf8dSHong Zhang ierr = PetscLLAdd(bnzj,bjj,bn,nlnk,lnk,lnkbt);CHKERRQ(ierr); 701c239cc6SHong Zhang cnzi += nlnk; 7158c24d83SHong Zhang } 7258c24d83SHong Zhang 7358c24d83SHong Zhang /* If free space is not available, make more free space */ 7458c24d83SHong Zhang /* Double the amount of total space in the list */ 7558c24d83SHong Zhang if (current_space->local_remaining<cnzi) { 764238b7adSHong Zhang ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 77c5db241fSHong Zhang nspacedouble++; 7858c24d83SHong Zhang } 7958c24d83SHong Zhang 80c5db241fSHong Zhang /* Copy data into free space, then initialize lnk */ 81be0fcf8dSHong Zhang ierr = PetscLLClean(bn,bn,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 82c5db241fSHong Zhang current_space->array += cnzi; 8358c24d83SHong Zhang current_space->local_used += cnzi; 8458c24d83SHong Zhang current_space->local_remaining -= cnzi; 8558c24d83SHong Zhang 8658c24d83SHong Zhang ci[i+1] = ci[i] + cnzi; 8758c24d83SHong Zhang } 8858c24d83SHong Zhang 8958c24d83SHong Zhang /* Column indices are in the list of free space */ 9058c24d83SHong Zhang /* Allocate space for cj, initialize cj, and */ 9158c24d83SHong Zhang /* destroy list of free space and other temporary array(s) */ 9238baddfdSBarry Smith ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr); 93a1a86e44SBarry Smith ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 94be0fcf8dSHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 9558c24d83SHong Zhang 9658c24d83SHong Zhang /* Allocate space for ca */ 9758c24d83SHong Zhang ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 9858c24d83SHong Zhang ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr); 9958c24d83SHong Zhang 10026be0446SHong Zhang /* put together the new symbolic matrix */ 1017adad957SLisandro Dalcin ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bn,ci,cj,ca,C);CHKERRQ(ierr); 10258c24d83SHong Zhang 10358c24d83SHong Zhang /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 10458c24d83SHong Zhang /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */ 10558c24d83SHong Zhang c = (Mat_SeqAIJ *)((*C)->data); 106e6b907acSBarry Smith c->free_a = PETSC_TRUE; 107e6b907acSBarry Smith c->free_ij = PETSC_TRUE; 10858c24d83SHong Zhang c->nonew = 0; 10958c24d83SHong Zhang 110f2b054eeSHong Zhang #if defined(PETSC_USE_INFO) 111f2b054eeSHong Zhang if (ci[am] != 0) { 112f2b054eeSHong Zhang PetscReal afill = ((PetscReal)ci[am])/(ai[am]+bi[bm]); 113f2b054eeSHong Zhang if (afill < 1.0) afill = 1.0; 114f2b054eeSHong Zhang ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr); 115f2b054eeSHong Zhang ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr); 116f2b054eeSHong Zhang } else { 117f2b054eeSHong Zhang ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr); 118be0fcf8dSHong Zhang } 119f2b054eeSHong Zhang #endif 12058c24d83SHong Zhang PetscFunctionReturn(0); 12158c24d83SHong Zhang } 122d50806bdSBarry Smith 1235c66b693SKris Buschelman 12426be0446SHong Zhang #undef __FUNCT__ 12526be0446SHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqAIJ" 126dfbe8321SBarry Smith PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C) 127d50806bdSBarry Smith { 128dfbe8321SBarry Smith PetscErrorCode ierr; 129d13dce4bSSatish Balay PetscLogDouble flops=0.0; 130d50806bdSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 131d50806bdSBarry Smith Mat_SeqAIJ *b = (Mat_SeqAIJ *)B->data; 132d50806bdSBarry Smith Mat_SeqAIJ *c = (Mat_SeqAIJ *)C->data; 13338baddfdSBarry Smith PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j; 134d0f46423SBarry Smith PetscInt am=A->rmap->N,cm=C->rmap->N; 135c124e916SHong Zhang PetscInt i,j,k,anzi,bnzi,cnzi,brow,nextb; 136c124e916SHong Zhang MatScalar *aa=a->a,*ba=b->a,*baj,*ca=c->a; 137d50806bdSBarry Smith 138d50806bdSBarry Smith PetscFunctionBegin; 139c124e916SHong Zhang /* clean old values in C */ 140c124e916SHong Zhang ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 141d50806bdSBarry Smith /* Traverse A row-wise. */ 142d50806bdSBarry Smith /* Build the ith row in C by summing over nonzero columns in A, */ 143d50806bdSBarry Smith /* the rows of B corresponding to nonzeros of A. */ 144d50806bdSBarry Smith for (i=0;i<am;i++) { 145d50806bdSBarry Smith anzi = ai[i+1] - ai[i]; 146d50806bdSBarry Smith for (j=0;j<anzi;j++) { 147d50806bdSBarry Smith brow = *aj++; 148d50806bdSBarry Smith bnzi = bi[brow+1] - bi[brow]; 149d50806bdSBarry Smith bjj = bj + bi[brow]; 150d50806bdSBarry Smith baj = ba + bi[brow]; 151c124e916SHong Zhang nextb = 0; 152c124e916SHong Zhang for (k=0; nextb<bnzi; k++) { 153c124e916SHong Zhang if (cj[k] == bjj[nextb]){ /* ccol == bcol */ 154c124e916SHong Zhang ca[k] += (*aa)*baj[nextb++]; 155c124e916SHong Zhang } 156d50806bdSBarry Smith } 157d50806bdSBarry Smith flops += 2*bnzi; 158d50806bdSBarry Smith aa++; 159d50806bdSBarry Smith } 160d50806bdSBarry Smith cnzi = ci[i+1] - ci[i]; 161d50806bdSBarry Smith ca += cnzi; 162d50806bdSBarry Smith cj += cnzi; 163d50806bdSBarry Smith } 164716bacf3SKris Buschelman ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 165716bacf3SKris Buschelman ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 166716bacf3SKris Buschelman 167d50806bdSBarry Smith ierr = PetscLogFlops(flops);CHKERRQ(ierr); 168d50806bdSBarry Smith PetscFunctionReturn(0); 169d50806bdSBarry Smith } 170bc011b1eSHong Zhang 171bc011b1eSHong Zhang 172bc011b1eSHong Zhang #undef __FUNCT__ 173bc011b1eSHong Zhang #define __FUNCT__ "MatMatMultTranspose_SeqAIJ_SeqAIJ" 174bc011b1eSHong Zhang PetscErrorCode MatMatMultTranspose_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) { 175bc011b1eSHong Zhang PetscErrorCode ierr; 176bc011b1eSHong Zhang 177bc011b1eSHong Zhang PetscFunctionBegin; 178bc011b1eSHong Zhang if (scall == MAT_INITIAL_MATRIX){ 179bc011b1eSHong Zhang ierr = MatMatMultTransposeSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr); 180bc011b1eSHong Zhang } 181bc011b1eSHong Zhang ierr = MatMatMultTransposeNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr); 182bc011b1eSHong Zhang PetscFunctionReturn(0); 183bc011b1eSHong Zhang } 184bc011b1eSHong Zhang 185bc011b1eSHong Zhang #undef __FUNCT__ 186bc011b1eSHong Zhang #define __FUNCT__ "MatMatMultTransposeSymbolic_SeqAIJ_SeqAIJ" 187bc011b1eSHong Zhang PetscErrorCode MatMatMultTransposeSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C) 188bc011b1eSHong Zhang { 189bc011b1eSHong Zhang PetscErrorCode ierr; 190bc011b1eSHong Zhang Mat At; 19138baddfdSBarry Smith PetscInt *ati,*atj; 192bc011b1eSHong Zhang 193bc011b1eSHong Zhang PetscFunctionBegin; 194bc011b1eSHong Zhang /* create symbolic At */ 195bc011b1eSHong Zhang ierr = MatGetSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr); 196d0f46423SBarry Smith ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,A->cmap->n,A->rmap->n,ati,atj,PETSC_NULL,&At);CHKERRQ(ierr); 197bc011b1eSHong Zhang 198bc011b1eSHong Zhang /* get symbolic C=At*B */ 199bc011b1eSHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(At,B,fill,C);CHKERRQ(ierr); 200bc011b1eSHong Zhang 201bc011b1eSHong Zhang /* clean up */ 202bc011b1eSHong Zhang ierr = MatDestroy(At);CHKERRQ(ierr); 203bc011b1eSHong Zhang ierr = MatRestoreSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr); 204bc011b1eSHong Zhang 205bc011b1eSHong Zhang PetscFunctionReturn(0); 206bc011b1eSHong Zhang } 207bc011b1eSHong Zhang 208bc011b1eSHong Zhang #undef __FUNCT__ 209bc011b1eSHong Zhang #define __FUNCT__ "MatMatMultTransposeNumeric_SeqAIJ_SeqAIJ" 210bc011b1eSHong Zhang PetscErrorCode MatMatMultTransposeNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C) 211bc011b1eSHong Zhang { 212bc011b1eSHong Zhang PetscErrorCode ierr; 2130fbc74f4SHong Zhang Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data; 214d0f46423SBarry Smith PetscInt am=A->rmap->n,anzi,*ai=a->i,*aj=a->j,*bi=b->i,*bj,bnzi,nextb; 215d0f46423SBarry Smith PetscInt cm=C->rmap->n,*ci=c->i,*cj=c->j,crow,*cjj,i,j,k; 216d13dce4bSSatish Balay PetscLogDouble flops=0.0; 2170fbc74f4SHong Zhang MatScalar *aa=a->a,*ba,*ca=c->a,*caj; 218bc011b1eSHong Zhang 219bc011b1eSHong Zhang PetscFunctionBegin; 220bc011b1eSHong Zhang /* clear old values in C */ 221bc011b1eSHong Zhang ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 222bc011b1eSHong Zhang 223bc011b1eSHong Zhang /* compute A^T*B using outer product (A^T)[:,i]*B[i,:] */ 224bc011b1eSHong Zhang for (i=0;i<am;i++) { 225bc011b1eSHong Zhang bj = b->j + bi[i]; 226bc011b1eSHong Zhang ba = b->a + bi[i]; 227bc011b1eSHong Zhang bnzi = bi[i+1] - bi[i]; 228bc011b1eSHong Zhang anzi = ai[i+1] - ai[i]; 229bc011b1eSHong Zhang for (j=0; j<anzi; j++) { 230bc011b1eSHong Zhang nextb = 0; 2310fbc74f4SHong Zhang crow = *aj++; 232bc011b1eSHong Zhang cjj = cj + ci[crow]; 233bc011b1eSHong Zhang caj = ca + ci[crow]; 234bc011b1eSHong Zhang /* perform sparse axpy operation. Note cjj includes bj. */ 235bc011b1eSHong Zhang for (k=0; nextb<bnzi; k++) { 2360fbc74f4SHong Zhang if (cjj[k] == *(bj+nextb)) { /* ccol == bcol */ 2370fbc74f4SHong Zhang caj[k] += (*aa)*(*(ba+nextb)); 238bc011b1eSHong Zhang nextb++; 239bc011b1eSHong Zhang } 240bc011b1eSHong Zhang } 241bc011b1eSHong Zhang flops += 2*bnzi; 2420fbc74f4SHong Zhang aa++; 243bc011b1eSHong Zhang } 244bc011b1eSHong Zhang } 245bc011b1eSHong Zhang 246bc011b1eSHong Zhang /* Assemble the final matrix and clean up */ 247bc011b1eSHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 248bc011b1eSHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 249bc011b1eSHong Zhang ierr = PetscLogFlops(flops);CHKERRQ(ierr); 250bc011b1eSHong Zhang PetscFunctionReturn(0); 251bc011b1eSHong Zhang } 2529123193aSHong Zhang 253ec01deb9SMatthew Knepley EXTERN_C_BEGIN 2549123193aSHong Zhang #undef __FUNCT__ 2559123193aSHong Zhang #define __FUNCT__ "MatMatMult_SeqAIJ_SeqDense" 2569123193aSHong Zhang PetscErrorCode MatMatMult_SeqAIJ_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 2579123193aSHong Zhang { 2589123193aSHong Zhang PetscErrorCode ierr; 2599123193aSHong Zhang 2609123193aSHong Zhang PetscFunctionBegin; 2619123193aSHong Zhang if (scall == MAT_INITIAL_MATRIX){ 2629123193aSHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqDense(A,B,fill,C);CHKERRQ(ierr); 2639123193aSHong Zhang } 2649123193aSHong Zhang ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,B,*C);CHKERRQ(ierr); 2659123193aSHong Zhang PetscFunctionReturn(0); 2669123193aSHong Zhang } 267ec01deb9SMatthew Knepley EXTERN_C_END 268edd81eecSMatthew Knepley 2699123193aSHong Zhang #undef __FUNCT__ 2709123193aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqDense" 2719123193aSHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C) 2729123193aSHong Zhang { 2739123193aSHong Zhang PetscErrorCode ierr; 2749123193aSHong Zhang 2759123193aSHong Zhang PetscFunctionBegin; 2769123193aSHong Zhang ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,0.0,C); 2779123193aSHong Zhang PetscFunctionReturn(0); 2789123193aSHong Zhang } 2799123193aSHong Zhang 2809123193aSHong Zhang #undef __FUNCT__ 2819123193aSHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqDense" 2829123193aSHong Zhang PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C) 2839123193aSHong Zhang { 284f32d5d43SBarry Smith Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data; 2859123193aSHong Zhang PetscErrorCode ierr; 286dd6ea824SBarry Smith PetscScalar *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4; 287dd6ea824SBarry Smith MatScalar *aa; 288d0f46423SBarry Smith PetscInt cm=C->rmap->n, cn=B->cmap->n, bm=B->rmap->n, col, i,j,n,*aj, am = A->rmap->n; 289f32d5d43SBarry Smith PetscInt am2 = 2*am, am3 = 3*am, bm4 = 4*bm,colam; 2909123193aSHong Zhang 2919123193aSHong Zhang PetscFunctionBegin; 292f32d5d43SBarry Smith if (!cm || !cn) PetscFunctionReturn(0); 293*e32f2f54SBarry Smith if (bm != A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Number columns in A %D not equal rows in B %D\n",A->cmap->n,bm); 294*e32f2f54SBarry Smith if (A->rmap->n != C->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Number rows in C %D not equal rows in A %D\n",C->rmap->n,A->rmap->n); 295*e32f2f54SBarry Smith if (B->cmap->n != C->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Number columns in B %D not equal columns in C %D\n",B->cmap->n,C->cmap->n); 296f32d5d43SBarry Smith ierr = MatGetArray(B,&b);CHKERRQ(ierr); 297f32d5d43SBarry Smith ierr = MatGetArray(C,&c);CHKERRQ(ierr); 298f32d5d43SBarry Smith b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm; 299f32d5d43SBarry Smith for (col=0; col<cn-4; col += 4){ /* over columns of C */ 300f32d5d43SBarry Smith colam = col*am; 301f32d5d43SBarry Smith for (i=0; i<am; i++) { /* over rows of C in those columns */ 302f32d5d43SBarry Smith r1 = r2 = r3 = r4 = 0.0; 303f32d5d43SBarry Smith n = a->i[i+1] - a->i[i]; 304f32d5d43SBarry Smith aj = a->j + a->i[i]; 305f32d5d43SBarry Smith aa = a->a + a->i[i]; 306f32d5d43SBarry Smith for (j=0; j<n; j++) { 307f32d5d43SBarry Smith r1 += (*aa)*b1[*aj]; 308f32d5d43SBarry Smith r2 += (*aa)*b2[*aj]; 309f32d5d43SBarry Smith r3 += (*aa)*b3[*aj]; 310f32d5d43SBarry Smith r4 += (*aa++)*b4[*aj++]; 3119123193aSHong Zhang } 312f32d5d43SBarry Smith c[colam + i] = r1; 313f32d5d43SBarry Smith c[colam + am + i] = r2; 314f32d5d43SBarry Smith c[colam + am2 + i] = r3; 315f32d5d43SBarry Smith c[colam + am3 + i] = r4; 316f32d5d43SBarry Smith } 317f32d5d43SBarry Smith b1 += bm4; 318f32d5d43SBarry Smith b2 += bm4; 319f32d5d43SBarry Smith b3 += bm4; 320f32d5d43SBarry Smith b4 += bm4; 321f32d5d43SBarry Smith } 322f32d5d43SBarry Smith for (;col<cn; col++){ /* over extra columns of C */ 323f32d5d43SBarry Smith for (i=0; i<am; i++) { /* over rows of C in those columns */ 324f32d5d43SBarry Smith r1 = 0.0; 325f32d5d43SBarry Smith n = a->i[i+1] - a->i[i]; 326f32d5d43SBarry Smith aj = a->j + a->i[i]; 327f32d5d43SBarry Smith aa = a->a + a->i[i]; 328f32d5d43SBarry Smith 329f32d5d43SBarry Smith for (j=0; j<n; j++) { 330f32d5d43SBarry Smith r1 += (*aa++)*b1[*aj++]; 331f32d5d43SBarry Smith } 332f32d5d43SBarry Smith c[col*am + i] = r1; 333f32d5d43SBarry Smith } 334f32d5d43SBarry Smith b1 += bm; 335f32d5d43SBarry Smith } 336dc0b31edSSatish Balay ierr = PetscLogFlops(cn*(2.0*a->nz));CHKERRQ(ierr); 337f32d5d43SBarry Smith ierr = MatRestoreArray(B,&b);CHKERRQ(ierr); 338f32d5d43SBarry Smith ierr = MatRestoreArray(C,&c);CHKERRQ(ierr); 339f32d5d43SBarry Smith ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 340f32d5d43SBarry Smith ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 341f32d5d43SBarry Smith PetscFunctionReturn(0); 342f32d5d43SBarry Smith } 343f32d5d43SBarry Smith 344f32d5d43SBarry Smith /* 3454324174eSBarry Smith Note very similar to MatMult_SeqAIJ(), should generate both codes from same base 346f32d5d43SBarry Smith */ 347f32d5d43SBarry Smith #undef __FUNCT__ 348f32d5d43SBarry Smith #define __FUNCT__ "MatMatMultNumericAdd_SeqAIJ_SeqDense" 349f32d5d43SBarry Smith PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat A,Mat B,Mat C) 350f32d5d43SBarry Smith { 351f32d5d43SBarry Smith Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data; 352f32d5d43SBarry Smith PetscErrorCode ierr; 353dd6ea824SBarry Smith PetscScalar *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4; 354dd6ea824SBarry Smith MatScalar *aa; 355d0f46423SBarry Smith PetscInt cm=C->rmap->n, cn=B->cmap->n, bm=B->rmap->n, col, i,j,n,*aj, am = A->rmap->n,*ii,arm; 3564324174eSBarry Smith PetscInt am2 = 2*am, am3 = 3*am, bm4 = 4*bm,colam,*ridx; 357f32d5d43SBarry Smith 358f32d5d43SBarry Smith PetscFunctionBegin; 359f32d5d43SBarry Smith if (!cm || !cn) PetscFunctionReturn(0); 360f32d5d43SBarry Smith ierr = MatGetArray(B,&b);CHKERRQ(ierr); 361f32d5d43SBarry Smith ierr = MatGetArray(C,&c);CHKERRQ(ierr); 362f32d5d43SBarry Smith b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm; 3634324174eSBarry Smith 3644324174eSBarry Smith if (a->compressedrow.use){ /* use compressed row format */ 3654324174eSBarry Smith for (col=0; col<cn-4; col += 4){ /* over columns of C */ 3664324174eSBarry Smith colam = col*am; 3674324174eSBarry Smith arm = a->compressedrow.nrows; 3684324174eSBarry Smith ii = a->compressedrow.i; 3694324174eSBarry Smith ridx = a->compressedrow.rindex; 3704324174eSBarry Smith for (i=0; i<arm; i++) { /* over rows of C in those columns */ 3714324174eSBarry Smith r1 = r2 = r3 = r4 = 0.0; 3724324174eSBarry Smith n = ii[i+1] - ii[i]; 3734324174eSBarry Smith aj = a->j + ii[i]; 3744324174eSBarry Smith aa = a->a + ii[i]; 3754324174eSBarry Smith for (j=0; j<n; j++) { 3764324174eSBarry Smith r1 += (*aa)*b1[*aj]; 3774324174eSBarry Smith r2 += (*aa)*b2[*aj]; 3784324174eSBarry Smith r3 += (*aa)*b3[*aj]; 3794324174eSBarry Smith r4 += (*aa++)*b4[*aj++]; 3804324174eSBarry Smith } 3814324174eSBarry Smith c[colam + ridx[i]] += r1; 3824324174eSBarry Smith c[colam + am + ridx[i]] += r2; 3834324174eSBarry Smith c[colam + am2 + ridx[i]] += r3; 3844324174eSBarry Smith c[colam + am3 + ridx[i]] += r4; 3854324174eSBarry Smith } 3864324174eSBarry Smith b1 += bm4; 3874324174eSBarry Smith b2 += bm4; 3884324174eSBarry Smith b3 += bm4; 3894324174eSBarry Smith b4 += bm4; 3904324174eSBarry Smith } 3914324174eSBarry Smith for (;col<cn; col++){ /* over extra columns of C */ 3924324174eSBarry Smith colam = col*am; 3934324174eSBarry Smith arm = a->compressedrow.nrows; 3944324174eSBarry Smith ii = a->compressedrow.i; 3954324174eSBarry Smith ridx = a->compressedrow.rindex; 3964324174eSBarry Smith for (i=0; i<arm; i++) { /* over rows of C in those columns */ 3974324174eSBarry Smith r1 = 0.0; 3984324174eSBarry Smith n = ii[i+1] - ii[i]; 3994324174eSBarry Smith aj = a->j + ii[i]; 4004324174eSBarry Smith aa = a->a + ii[i]; 4014324174eSBarry Smith 4024324174eSBarry Smith for (j=0; j<n; j++) { 4034324174eSBarry Smith r1 += (*aa++)*b1[*aj++]; 4044324174eSBarry Smith } 4054324174eSBarry Smith c[col*am + ridx[i]] += r1; 4064324174eSBarry Smith } 4074324174eSBarry Smith b1 += bm; 4084324174eSBarry Smith } 4094324174eSBarry Smith } else { 410f32d5d43SBarry Smith for (col=0; col<cn-4; col += 4){ /* over columns of C */ 411f32d5d43SBarry Smith colam = col*am; 412f32d5d43SBarry Smith for (i=0; i<am; i++) { /* over rows of C in those columns */ 413f32d5d43SBarry Smith r1 = r2 = r3 = r4 = 0.0; 414f32d5d43SBarry Smith n = a->i[i+1] - a->i[i]; 415f32d5d43SBarry Smith aj = a->j + a->i[i]; 416f32d5d43SBarry Smith aa = a->a + a->i[i]; 417f32d5d43SBarry Smith for (j=0; j<n; j++) { 418f32d5d43SBarry Smith r1 += (*aa)*b1[*aj]; 419f32d5d43SBarry Smith r2 += (*aa)*b2[*aj]; 420f32d5d43SBarry Smith r3 += (*aa)*b3[*aj]; 421f32d5d43SBarry Smith r4 += (*aa++)*b4[*aj++]; 422f32d5d43SBarry Smith } 423f32d5d43SBarry Smith c[colam + i] += r1; 424f32d5d43SBarry Smith c[colam + am + i] += r2; 425f32d5d43SBarry Smith c[colam + am2 + i] += r3; 426f32d5d43SBarry Smith c[colam + am3 + i] += r4; 427f32d5d43SBarry Smith } 428f32d5d43SBarry Smith b1 += bm4; 429f32d5d43SBarry Smith b2 += bm4; 430f32d5d43SBarry Smith b3 += bm4; 431f32d5d43SBarry Smith b4 += bm4; 432f32d5d43SBarry Smith } 433f32d5d43SBarry Smith for (;col<cn; col++){ /* over extra columns of C */ 434f32d5d43SBarry Smith for (i=0; i<am; i++) { /* over rows of C in those columns */ 435f32d5d43SBarry Smith r1 = 0.0; 436f32d5d43SBarry Smith n = a->i[i+1] - a->i[i]; 437f32d5d43SBarry Smith aj = a->j + a->i[i]; 438f32d5d43SBarry Smith aa = a->a + a->i[i]; 439f32d5d43SBarry Smith 440f32d5d43SBarry Smith for (j=0; j<n; j++) { 441f32d5d43SBarry Smith r1 += (*aa++)*b1[*aj++]; 442f32d5d43SBarry Smith } 443f32d5d43SBarry Smith c[col*am + i] += r1; 444f32d5d43SBarry Smith } 445f32d5d43SBarry Smith b1 += bm; 446f32d5d43SBarry Smith } 4474324174eSBarry Smith } 448dc0b31edSSatish Balay ierr = PetscLogFlops(cn*2.0*a->nz);CHKERRQ(ierr); 449f32d5d43SBarry Smith ierr = MatRestoreArray(B,&b);CHKERRQ(ierr); 450f32d5d43SBarry Smith ierr = MatRestoreArray(C,&c);CHKERRQ(ierr); 4519123193aSHong Zhang PetscFunctionReturn(0); 4529123193aSHong Zhang } 453