1be1d678aSKris Buschelman 2d50806bdSBarry Smith /* 32499ec78SHong Zhang Defines matrix-matrix product routines for pairs of SeqAIJ matrices 4d50806bdSBarry Smith C = A * B 5d50806bdSBarry Smith */ 6d50806bdSBarry Smith 7c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h> /*I "petscmat.h" I*/ 8c6db04a5SJed Brown #include <../src/mat/utils/freespace.h> 9c6db04a5SJed Brown #include <petscbt.h> 10c6db04a5SJed Brown #include <../src/mat/impls/dense/seq/dense.h> /*I "petscmat.h" I*/ 116284ec50SHong Zhang 12ec01deb9SMatthew Knepley EXTERN_C_BEGIN 136284ec50SHong Zhang #undef __FUNCT__ 145c66b693SKris Buschelman #define __FUNCT__ "MatMatMult_SeqAIJ_SeqAIJ" 1538baddfdSBarry Smith PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 1638baddfdSBarry Smith { 17dfbe8321SBarry Smith PetscErrorCode ierr; 185c66b693SKris Buschelman 195c66b693SKris Buschelman PetscFunctionBegin; 2026be0446SHong Zhang if (scall == MAT_INITIAL_MATRIX){ 2126be0446SHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr); 2226be0446SHong Zhang } 2326be0446SHong Zhang ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr); 245c66b693SKris Buschelman PetscFunctionReturn(0); 255c66b693SKris Buschelman } 26ec01deb9SMatthew Knepley EXTERN_C_END 271c24bd37SHong Zhang 2826be0446SHong Zhang #undef __FUNCT__ 2926be0446SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ" 30dfbe8321SBarry Smith PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C) 3158c24d83SHong Zhang { 32dfbe8321SBarry Smith PetscErrorCode ierr; 33a1a86e44SBarry Smith PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 3458c24d83SHong Zhang Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c; 3538baddfdSBarry Smith PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci,*cj; 36d0f46423SBarry Smith PetscInt am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N; 3738baddfdSBarry Smith PetscInt i,j,anzi,brow,bnzj,cnzi,nlnk,*lnk,nspacedouble=0; 3858c24d83SHong Zhang MatScalar *ca; 39be0fcf8dSHong Zhang PetscBT lnkbt; 407212da7cSHong Zhang PetscReal afill; 4158c24d83SHong Zhang 4258c24d83SHong Zhang PetscFunctionBegin; 4358c24d83SHong Zhang /* Allocate ci array, arrays for fill computation and */ 4458c24d83SHong Zhang /* free space for accumulating nonzero column info */ 4538baddfdSBarry Smith ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 4658c24d83SHong Zhang ci[0] = 0; 4758c24d83SHong Zhang 48be0fcf8dSHong Zhang /* create and initialize a linked list */ 49be0fcf8dSHong Zhang nlnk = bn+1; 50be0fcf8dSHong Zhang ierr = PetscLLCreate(bn,bn,nlnk,lnk,lnkbt);CHKERRQ(ierr); 5158c24d83SHong Zhang 52c5db241fSHong Zhang /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */ 53a1a86e44SBarry Smith ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr); 5458c24d83SHong Zhang current_space = free_space; 5558c24d83SHong Zhang 5658c24d83SHong Zhang /* Determine symbolic info for each row of the product: */ 5758c24d83SHong Zhang for (i=0;i<am;i++) { 5858c24d83SHong Zhang anzi = ai[i+1] - ai[i]; 5958c24d83SHong Zhang cnzi = 0; 602d09714cSHong Zhang j = anzi; 612d09714cSHong Zhang aj = a->j + ai[i]; 622d09714cSHong Zhang while (j){/* assume cols are almost in increasing order, starting from its end saves computation */ 632d09714cSHong Zhang j--; 642d09714cSHong Zhang brow = *(aj + j); 6558c24d83SHong Zhang bnzj = bi[brow+1] - bi[brow]; 6658c24d83SHong Zhang bjj = bj + bi[brow]; 671c239cc6SHong Zhang /* add non-zero cols of B into the sorted linked list lnk */ 68be0fcf8dSHong Zhang ierr = PetscLLAdd(bnzj,bjj,bn,nlnk,lnk,lnkbt);CHKERRQ(ierr); 691c239cc6SHong Zhang cnzi += nlnk; 7058c24d83SHong Zhang } 7158c24d83SHong Zhang 7258c24d83SHong Zhang /* If free space is not available, make more free space */ 7358c24d83SHong Zhang /* Double the amount of total space in the list */ 7458c24d83SHong Zhang if (current_space->local_remaining<cnzi) { 754238b7adSHong Zhang ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 76c5db241fSHong Zhang nspacedouble++; 7758c24d83SHong Zhang } 7858c24d83SHong Zhang 79c5db241fSHong Zhang /* Copy data into free space, then initialize lnk */ 80be0fcf8dSHong Zhang ierr = PetscLLClean(bn,bn,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 81c5db241fSHong Zhang current_space->array += cnzi; 8258c24d83SHong Zhang current_space->local_used += cnzi; 8358c24d83SHong Zhang current_space->local_remaining -= cnzi; 8458c24d83SHong Zhang 8558c24d83SHong Zhang ci[i+1] = ci[i] + cnzi; 8658c24d83SHong Zhang } 8758c24d83SHong Zhang 8858c24d83SHong Zhang /* Column indices are in the list of free space */ 8958c24d83SHong Zhang /* Allocate space for cj, initialize cj, and */ 9058c24d83SHong Zhang /* destroy list of free space and other temporary array(s) */ 9138baddfdSBarry Smith ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr); 92a1a86e44SBarry Smith ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 93be0fcf8dSHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 9458c24d83SHong Zhang 9558c24d83SHong Zhang /* Allocate space for ca */ 9658c24d83SHong Zhang ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 9758c24d83SHong Zhang ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr); 9858c24d83SHong Zhang 9926be0446SHong Zhang /* put together the new symbolic matrix */ 1007adad957SLisandro Dalcin ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bn,ci,cj,ca,C);CHKERRQ(ierr); 10158c24d83SHong Zhang 10258c24d83SHong Zhang /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 10358c24d83SHong Zhang /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */ 10458c24d83SHong Zhang c = (Mat_SeqAIJ *)((*C)->data); 105e6b907acSBarry Smith c->free_a = PETSC_TRUE; 106e6b907acSBarry Smith c->free_ij = PETSC_TRUE; 10758c24d83SHong Zhang c->nonew = 0; 10858c24d83SHong Zhang 1097212da7cSHong Zhang /* set MatInfo */ 1107212da7cSHong Zhang afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5; 111f2b054eeSHong Zhang if (afill < 1.0) afill = 1.0; 1127212da7cSHong Zhang c->maxnz = ci[am]; 1137212da7cSHong Zhang c->nz = ci[am]; 1147212da7cSHong Zhang (*C)->info.mallocs = nspacedouble; 1157212da7cSHong Zhang (*C)->info.fill_ratio_given = fill; 1167212da7cSHong Zhang (*C)->info.fill_ratio_needed = afill; 1177212da7cSHong Zhang 1187212da7cSHong Zhang #if defined(PETSC_USE_INFO) 1197212da7cSHong Zhang if (ci[am]) { 120f2b054eeSHong Zhang ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr); 121f2b054eeSHong Zhang ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr); 122f2b054eeSHong Zhang } else { 123f2b054eeSHong Zhang ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr); 124be0fcf8dSHong Zhang } 125f2b054eeSHong Zhang #endif 12658c24d83SHong Zhang PetscFunctionReturn(0); 12758c24d83SHong Zhang } 128d50806bdSBarry Smith 1295c66b693SKris Buschelman 13026be0446SHong Zhang #undef __FUNCT__ 13126be0446SHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqAIJ" 132dfbe8321SBarry Smith PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C) 133d50806bdSBarry Smith { 134dfbe8321SBarry Smith PetscErrorCode ierr; 135d13dce4bSSatish Balay PetscLogDouble flops=0.0; 136d50806bdSBarry Smith Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 137d50806bdSBarry Smith Mat_SeqAIJ *b = (Mat_SeqAIJ *)B->data; 138d50806bdSBarry Smith Mat_SeqAIJ *c = (Mat_SeqAIJ *)C->data; 13938baddfdSBarry Smith PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j; 140d0f46423SBarry Smith PetscInt am=A->rmap->N,cm=C->rmap->N; 141c124e916SHong Zhang PetscInt i,j,k,anzi,bnzi,cnzi,brow,nextb; 142c124e916SHong Zhang MatScalar *aa=a->a,*ba=b->a,*baj,*ca=c->a; 143d50806bdSBarry Smith 144d50806bdSBarry Smith PetscFunctionBegin; 145c124e916SHong Zhang /* clean old values in C */ 146c124e916SHong Zhang ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 147d50806bdSBarry Smith /* Traverse A row-wise. */ 148d50806bdSBarry Smith /* Build the ith row in C by summing over nonzero columns in A, */ 149d50806bdSBarry Smith /* the rows of B corresponding to nonzeros of A. */ 150d50806bdSBarry Smith for (i=0;i<am;i++) { 151d50806bdSBarry Smith anzi = ai[i+1] - ai[i]; 152d50806bdSBarry Smith for (j=0;j<anzi;j++) { 153d50806bdSBarry Smith brow = *aj++; 154d50806bdSBarry Smith bnzi = bi[brow+1] - bi[brow]; 155d50806bdSBarry Smith bjj = bj + bi[brow]; 156d50806bdSBarry Smith baj = ba + bi[brow]; 157c124e916SHong Zhang nextb = 0; 158c124e916SHong Zhang for (k=0; nextb<bnzi; k++) { 159c124e916SHong Zhang if (cj[k] == bjj[nextb]){ /* ccol == bcol */ 160c124e916SHong Zhang ca[k] += (*aa)*baj[nextb++]; 161c124e916SHong Zhang } 162d50806bdSBarry Smith } 163d50806bdSBarry Smith flops += 2*bnzi; 164d50806bdSBarry Smith aa++; 165d50806bdSBarry Smith } 166d50806bdSBarry Smith cnzi = ci[i+1] - ci[i]; 167d50806bdSBarry Smith ca += cnzi; 168d50806bdSBarry Smith cj += cnzi; 169d50806bdSBarry Smith } 170716bacf3SKris Buschelman ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 171716bacf3SKris Buschelman ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 172716bacf3SKris Buschelman 173d50806bdSBarry Smith ierr = PetscLogFlops(flops);CHKERRQ(ierr); 174d50806bdSBarry Smith PetscFunctionReturn(0); 175d50806bdSBarry Smith } 176bc011b1eSHong Zhang 177bc011b1eSHong Zhang #undef __FUNCT__ 178bc011b1eSHong Zhang #define __FUNCT__ "MatMatMultTranspose_SeqAIJ_SeqAIJ" 1795df89d91SHong Zhang PetscErrorCode MatMatMultTranspose_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 1805df89d91SHong Zhang { 181bc011b1eSHong Zhang PetscErrorCode ierr; 182bc011b1eSHong Zhang 183bc011b1eSHong Zhang PetscFunctionBegin; 184bc011b1eSHong Zhang if (scall == MAT_INITIAL_MATRIX){ 185bc011b1eSHong Zhang ierr = MatMatMultTransposeSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr); 186bc011b1eSHong Zhang } 187bc011b1eSHong Zhang ierr = MatMatMultTransposeNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr); 188bc011b1eSHong Zhang PetscFunctionReturn(0); 189bc011b1eSHong Zhang } 190bc011b1eSHong Zhang 191bc011b1eSHong Zhang #undef __FUNCT__ 192bc011b1eSHong Zhang #define __FUNCT__ "MatMatMultTransposeSymbolic_SeqAIJ_SeqAIJ" 193bc011b1eSHong Zhang PetscErrorCode MatMatMultTransposeSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C) 194bc011b1eSHong Zhang { 1955df89d91SHong Zhang PetscErrorCode ierr; 196*81d82fe4SHong Zhang Mat Bt; 197*81d82fe4SHong Zhang PetscInt *bti,*btj; 198*81d82fe4SHong Zhang 199*81d82fe4SHong Zhang PetscFunctionBegin; 200*81d82fe4SHong Zhang /* create symbolic Bt */ 201*81d82fe4SHong Zhang ierr = MatGetSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr); 202*81d82fe4SHong Zhang ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,B->cmap->n,B->rmap->n,bti,btj,PETSC_NULL,&Bt);CHKERRQ(ierr); 203*81d82fe4SHong Zhang 204*81d82fe4SHong Zhang /* get symbolic C=A*Bt */ 205*81d82fe4SHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,Bt,fill,C);CHKERRQ(ierr); 206*81d82fe4SHong Zhang 207*81d82fe4SHong Zhang /* clean up */ 208*81d82fe4SHong Zhang ierr = MatDestroy(&Bt);CHKERRQ(ierr); 209*81d82fe4SHong Zhang ierr = MatRestoreSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr); 210*81d82fe4SHong Zhang 211*81d82fe4SHong Zhang #if defined(INEFFICIENT_ALGORITHM) 212*81d82fe4SHong Zhang /* The algorithm below computes am*bm sparse inner-product - inefficient! It will be deleted later. */ 2131fa1209cSHong Zhang PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 2141fa1209cSHong Zhang Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c; 2151fa1209cSHong Zhang PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*ci,*cj,*acol,*bcol; 2161fa1209cSHong Zhang PetscInt am=A->rmap->N,bm=B->rmap->N; 2171fa1209cSHong Zhang PetscInt i,j,anzi,bnzj,cnzi,nlnk,*lnk,nspacedouble=0,ka,kb,index[1]; 2181fa1209cSHong Zhang MatScalar *ca; 2191fa1209cSHong Zhang PetscBT lnkbt; 2201fa1209cSHong Zhang PetscReal afill; 2215df89d91SHong Zhang 2221fa1209cSHong Zhang /* Allocate row pointer array ci */ 2231fa1209cSHong Zhang ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 2241fa1209cSHong Zhang ci[0] = 0; 2251fa1209cSHong Zhang 2261fa1209cSHong Zhang /* Create and initialize a linked list for C columns */ 2271fa1209cSHong Zhang nlnk = bm+1; 2281fa1209cSHong Zhang ierr = PetscLLCreate(bm,bm,nlnk,lnk,lnkbt);CHKERRQ(ierr); 2291fa1209cSHong Zhang 2301fa1209cSHong Zhang /* Initial FreeSpace with size fill*(nnz(A)+nnz(B)) */ 2311fa1209cSHong Zhang ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr); 2321fa1209cSHong Zhang current_space = free_space; 2331fa1209cSHong Zhang 2341fa1209cSHong Zhang /* Determine symbolic info for each row of the product A*B^T: */ 2351fa1209cSHong Zhang for (i=0; i<am; i++) { 2361fa1209cSHong Zhang anzi = ai[i+1] - ai[i]; 2371fa1209cSHong Zhang cnzi = 0; 2381fa1209cSHong Zhang acol = aj + ai[i]; 2391fa1209cSHong Zhang for (j=0; j<bm; j++){ 2401fa1209cSHong Zhang bnzj = bi[j+1] - bi[j]; 2411fa1209cSHong Zhang bcol= bj + bi[j]; 242*81d82fe4SHong Zhang /* sparse inner-product c(i,j)=A[i,:]*B[j,:]^T */ 2431fa1209cSHong Zhang ka = 0; kb = 0; 2441fa1209cSHong Zhang while (ka < anzi && kb < bnzj){ 245*81d82fe4SHong Zhang while (acol[ka] < bcol[kb] && ka < anzi) ka++; 246*81d82fe4SHong Zhang if (ka == anzi) break; 247*81d82fe4SHong Zhang while (acol[ka] > bcol[kb] && kb < bnzj) kb++; 248*81d82fe4SHong Zhang if (kb == bnzj) break; 2491fa1209cSHong Zhang if (acol[ka] == bcol[kb]){ /* add nonzero c(i,j) to lnk */ 2501fa1209cSHong Zhang index[0] = j; 2511fa1209cSHong Zhang ierr = PetscLLAdd(1,index,bm,nlnk,lnk,lnkbt);CHKERRQ(ierr); 2521fa1209cSHong Zhang cnzi++; 2531fa1209cSHong Zhang break; 2541fa1209cSHong Zhang } 2551fa1209cSHong Zhang } 2561fa1209cSHong Zhang } 2571fa1209cSHong Zhang 2581fa1209cSHong Zhang /* If free space is not available, make more free space */ 2591fa1209cSHong Zhang /* Double the amount of total space in the list */ 2601fa1209cSHong Zhang if (current_space->local_remaining<cnzi) { 2611fa1209cSHong Zhang ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 2621fa1209cSHong Zhang nspacedouble++; 2631fa1209cSHong Zhang } 2641fa1209cSHong Zhang 2651fa1209cSHong Zhang /* Copy data into free space, then initialize lnk */ 2661fa1209cSHong Zhang ierr = PetscLLClean(bm,bm,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 2671fa1209cSHong Zhang current_space->array += cnzi; 2681fa1209cSHong Zhang current_space->local_used += cnzi; 2691fa1209cSHong Zhang current_space->local_remaining -= cnzi; 2701fa1209cSHong Zhang 2711fa1209cSHong Zhang ci[i+1] = ci[i] + cnzi; 2721fa1209cSHong Zhang } 2731fa1209cSHong Zhang 2741fa1209cSHong Zhang 2751fa1209cSHong Zhang /* Column indices are in the list of free space. 2761fa1209cSHong Zhang Allocate array cj, copy column indices to cj, and destroy list of free space */ 2771fa1209cSHong Zhang ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr); 2781fa1209cSHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 2791fa1209cSHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 2801fa1209cSHong Zhang 2811fa1209cSHong Zhang /* Allocate space for ca */ 2821fa1209cSHong Zhang ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 2831fa1209cSHong Zhang ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr); 2841fa1209cSHong Zhang 2851fa1209cSHong Zhang /* put together the new symbolic matrix */ 2861fa1209cSHong Zhang ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bm,ci,cj,ca,C);CHKERRQ(ierr); 2871fa1209cSHong Zhang 2881fa1209cSHong Zhang /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 2891fa1209cSHong Zhang /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */ 2901fa1209cSHong Zhang c = (Mat_SeqAIJ *)((*C)->data); 2911fa1209cSHong Zhang c->free_a = PETSC_TRUE; 2921fa1209cSHong Zhang c->free_ij = PETSC_TRUE; 2931fa1209cSHong Zhang c->nonew = 0; 2941fa1209cSHong Zhang 2951fa1209cSHong Zhang /* set MatInfo */ 2961fa1209cSHong Zhang afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5; 2971fa1209cSHong Zhang if (afill < 1.0) afill = 1.0; 2981fa1209cSHong Zhang c->maxnz = ci[am]; 2991fa1209cSHong Zhang c->nz = ci[am]; 3001fa1209cSHong Zhang (*C)->info.mallocs = nspacedouble; 3011fa1209cSHong Zhang (*C)->info.fill_ratio_given = fill; 3021fa1209cSHong Zhang (*C)->info.fill_ratio_needed = afill; 3031fa1209cSHong Zhang 3041fa1209cSHong Zhang #if defined(PETSC_USE_INFO) 3051fa1209cSHong Zhang if (ci[am]) { 3061fa1209cSHong Zhang ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr); 3071fa1209cSHong Zhang ierr = PetscInfo1((*C),"Use MatMatMultTranspose(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr); 3081fa1209cSHong Zhang } else { 3091fa1209cSHong Zhang ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr); 3101fa1209cSHong Zhang } 3111fa1209cSHong Zhang #endif 312*81d82fe4SHong Zhang #endif 3135df89d91SHong Zhang PetscFunctionReturn(0); 3145df89d91SHong Zhang } 3155df89d91SHong Zhang 3165df89d91SHong Zhang #undef __FUNCT__ 3175df89d91SHong Zhang #define __FUNCT__ "MatMatMultTransposeNumeric_SeqAIJ_SeqAIJ" 3185df89d91SHong Zhang PetscErrorCode MatMatMultTransposeNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C) 3195df89d91SHong Zhang { 3205df89d91SHong Zhang PetscErrorCode ierr; 3215df89d91SHong Zhang Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data; 322*81d82fe4SHong Zhang PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,anzi,bnzj,nexta,nextb,*acol,*bcol,brow; 323*81d82fe4SHong Zhang PetscInt cm=C->rmap->n,*ci=c->i,*cj=c->j,i,j,cnzi,*ccol; 3245df89d91SHong Zhang PetscLogDouble flops=0.0; 325*81d82fe4SHong Zhang MatScalar *aa=a->a,*ba=b->a,*ca=c->a; 3265df89d91SHong Zhang 3275df89d91SHong Zhang PetscFunctionBegin; 328*81d82fe4SHong Zhang /* clear old values in C */ 329*81d82fe4SHong Zhang ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 3305df89d91SHong Zhang 3311fa1209cSHong Zhang for (i=0; i<cm; i++) { 332*81d82fe4SHong Zhang anzi = ai[i+1] - ai[i]; 333*81d82fe4SHong Zhang acol = aj + ai[i]; 3341fa1209cSHong Zhang cnzi = ci[i+1] - ci[i]; 3351fa1209cSHong Zhang ccol = cj + ci[i]; 3361fa1209cSHong Zhang for (j=0; j<cnzi; j++){ 337*81d82fe4SHong Zhang brow = ccol[j]; 338*81d82fe4SHong Zhang bnzj = bi[brow+1] - bi[brow]; 339*81d82fe4SHong Zhang bcol = bj + bi[brow]; 340*81d82fe4SHong Zhang /* perform sparse inner-product c(i,j)=A[i,:]*B[j,:]^T */ 341*81d82fe4SHong Zhang nexta = 0; nextb = 0; 342*81d82fe4SHong Zhang while (nexta<anzi && nextb<bnzj){ 343*81d82fe4SHong Zhang while (acol[nexta] < bcol[nextb] && nexta < anzi) nexta++; 344*81d82fe4SHong Zhang if (nexta == anzi) break; 345*81d82fe4SHong Zhang while (acol[nexta] > bcol[nextb] && nextb < bnzj) nextb++; 346*81d82fe4SHong Zhang if (nextb == bnzj) break; 347*81d82fe4SHong Zhang if (acol[nexta] == bcol[nextb]){ 348*81d82fe4SHong Zhang *(ca+ci[i]+j) += (*(aa+ai[i]+nexta))*(*(ba+bi[brow]+nextb)); 349*81d82fe4SHong Zhang nexta++; nextb++; 350*81d82fe4SHong Zhang flops += 2; 3511fa1209cSHong Zhang } 3521fa1209cSHong Zhang } 353*81d82fe4SHong Zhang } 354*81d82fe4SHong Zhang } 355*81d82fe4SHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 356*81d82fe4SHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 357*81d82fe4SHong Zhang ierr = PetscLogFlops(flops);CHKERRQ(ierr); 3585df89d91SHong Zhang PetscFunctionReturn(0); 3595df89d91SHong Zhang } 3605df89d91SHong Zhang 3615df89d91SHong Zhang #undef __FUNCT__ 3625df89d91SHong Zhang #define __FUNCT__ "MatMatTransposeMult_SeqAIJ_SeqAIJ" 3635df89d91SHong Zhang PetscErrorCode MatMatTransposeMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) { 3645df89d91SHong Zhang PetscErrorCode ierr; 3655df89d91SHong Zhang 3665df89d91SHong Zhang PetscFunctionBegin; 3675df89d91SHong Zhang if (scall == MAT_INITIAL_MATRIX){ 3685df89d91SHong Zhang ierr = MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr); 3695df89d91SHong Zhang } 3705df89d91SHong Zhang ierr = MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr); 3715df89d91SHong Zhang PetscFunctionReturn(0); 3725df89d91SHong Zhang } 3735df89d91SHong Zhang 3745df89d91SHong Zhang #undef __FUNCT__ 3755df89d91SHong Zhang #define __FUNCT__ "MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ" 3765df89d91SHong Zhang PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C) 3775df89d91SHong Zhang { 378bc011b1eSHong Zhang PetscErrorCode ierr; 379bc011b1eSHong Zhang Mat At; 38038baddfdSBarry Smith PetscInt *ati,*atj; 381bc011b1eSHong Zhang 382bc011b1eSHong Zhang PetscFunctionBegin; 383bc011b1eSHong Zhang /* create symbolic At */ 384bc011b1eSHong Zhang ierr = MatGetSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr); 385d0f46423SBarry Smith ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,A->cmap->n,A->rmap->n,ati,atj,PETSC_NULL,&At);CHKERRQ(ierr); 386bc011b1eSHong Zhang 387bc011b1eSHong Zhang /* get symbolic C=At*B */ 388bc011b1eSHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(At,B,fill,C);CHKERRQ(ierr); 389bc011b1eSHong Zhang 390bc011b1eSHong Zhang /* clean up */ 3916bf464f9SBarry Smith ierr = MatDestroy(&At);CHKERRQ(ierr); 392bc011b1eSHong Zhang ierr = MatRestoreSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr); 393bc011b1eSHong Zhang PetscFunctionReturn(0); 394bc011b1eSHong Zhang } 395bc011b1eSHong Zhang 396bc011b1eSHong Zhang #undef __FUNCT__ 3975df89d91SHong Zhang #define __FUNCT__ "MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ" 3985df89d91SHong Zhang PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C) 399bc011b1eSHong Zhang { 400bc011b1eSHong Zhang PetscErrorCode ierr; 4010fbc74f4SHong Zhang Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data; 402d0f46423SBarry Smith PetscInt am=A->rmap->n,anzi,*ai=a->i,*aj=a->j,*bi=b->i,*bj,bnzi,nextb; 403d0f46423SBarry Smith PetscInt cm=C->rmap->n,*ci=c->i,*cj=c->j,crow,*cjj,i,j,k; 404d13dce4bSSatish Balay PetscLogDouble flops=0.0; 4050fbc74f4SHong Zhang MatScalar *aa=a->a,*ba,*ca=c->a,*caj; 406bc011b1eSHong Zhang 407bc011b1eSHong Zhang PetscFunctionBegin; 408bc011b1eSHong Zhang /* clear old values in C */ 409bc011b1eSHong Zhang ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 410bc011b1eSHong Zhang 411bc011b1eSHong Zhang /* compute A^T*B using outer product (A^T)[:,i]*B[i,:] */ 412bc011b1eSHong Zhang for (i=0;i<am;i++) { 413bc011b1eSHong Zhang bj = b->j + bi[i]; 414bc011b1eSHong Zhang ba = b->a + bi[i]; 415bc011b1eSHong Zhang bnzi = bi[i+1] - bi[i]; 416bc011b1eSHong Zhang anzi = ai[i+1] - ai[i]; 417bc011b1eSHong Zhang for (j=0; j<anzi; j++) { 418bc011b1eSHong Zhang nextb = 0; 4190fbc74f4SHong Zhang crow = *aj++; 420bc011b1eSHong Zhang cjj = cj + ci[crow]; 421bc011b1eSHong Zhang caj = ca + ci[crow]; 422bc011b1eSHong Zhang /* perform sparse axpy operation. Note cjj includes bj. */ 423bc011b1eSHong Zhang for (k=0; nextb<bnzi; k++) { 4240fbc74f4SHong Zhang if (cjj[k] == *(bj+nextb)) { /* ccol == bcol */ 4250fbc74f4SHong Zhang caj[k] += (*aa)*(*(ba+nextb)); 426bc011b1eSHong Zhang nextb++; 427bc011b1eSHong Zhang } 428bc011b1eSHong Zhang } 429bc011b1eSHong Zhang flops += 2*bnzi; 4300fbc74f4SHong Zhang aa++; 431bc011b1eSHong Zhang } 432bc011b1eSHong Zhang } 433bc011b1eSHong Zhang 434bc011b1eSHong Zhang /* Assemble the final matrix and clean up */ 435bc011b1eSHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 436bc011b1eSHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 437bc011b1eSHong Zhang ierr = PetscLogFlops(flops);CHKERRQ(ierr); 438bc011b1eSHong Zhang PetscFunctionReturn(0); 439bc011b1eSHong Zhang } 4409123193aSHong Zhang 441ec01deb9SMatthew Knepley EXTERN_C_BEGIN 4429123193aSHong Zhang #undef __FUNCT__ 4439123193aSHong Zhang #define __FUNCT__ "MatMatMult_SeqAIJ_SeqDense" 4449123193aSHong Zhang PetscErrorCode MatMatMult_SeqAIJ_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 4459123193aSHong Zhang { 4469123193aSHong Zhang PetscErrorCode ierr; 4479123193aSHong Zhang 4489123193aSHong Zhang PetscFunctionBegin; 4499123193aSHong Zhang if (scall == MAT_INITIAL_MATRIX){ 4509123193aSHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqDense(A,B,fill,C);CHKERRQ(ierr); 4519123193aSHong Zhang } 4529123193aSHong Zhang ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,B,*C);CHKERRQ(ierr); 4539123193aSHong Zhang PetscFunctionReturn(0); 4549123193aSHong Zhang } 455ec01deb9SMatthew Knepley EXTERN_C_END 456edd81eecSMatthew Knepley 4579123193aSHong Zhang #undef __FUNCT__ 4589123193aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqDense" 4599123193aSHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C) 4609123193aSHong Zhang { 4619123193aSHong Zhang PetscErrorCode ierr; 4629123193aSHong Zhang 4639123193aSHong Zhang PetscFunctionBegin; 4645a586d82SBarry Smith ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,0.0,C);CHKERRQ(ierr); 4659123193aSHong Zhang PetscFunctionReturn(0); 4669123193aSHong Zhang } 4679123193aSHong Zhang 4689123193aSHong Zhang #undef __FUNCT__ 4699123193aSHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqDense" 4709123193aSHong Zhang PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C) 4719123193aSHong Zhang { 472f32d5d43SBarry Smith Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data; 4739123193aSHong Zhang PetscErrorCode ierr; 474dd6ea824SBarry Smith PetscScalar *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4; 475dd6ea824SBarry Smith MatScalar *aa; 476d0f46423SBarry Smith PetscInt cm=C->rmap->n, cn=B->cmap->n, bm=B->rmap->n, col, i,j,n,*aj, am = A->rmap->n; 477f32d5d43SBarry Smith PetscInt am2 = 2*am, am3 = 3*am, bm4 = 4*bm,colam; 4789123193aSHong Zhang 4799123193aSHong Zhang PetscFunctionBegin; 480f32d5d43SBarry Smith if (!cm || !cn) PetscFunctionReturn(0); 481e32f2f54SBarry 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); 482e32f2f54SBarry 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); 483e32f2f54SBarry 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); 484f32d5d43SBarry Smith ierr = MatGetArray(B,&b);CHKERRQ(ierr); 485f32d5d43SBarry Smith ierr = MatGetArray(C,&c);CHKERRQ(ierr); 486f32d5d43SBarry Smith b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm; 487f32d5d43SBarry Smith for (col=0; col<cn-4; col += 4){ /* over columns of C */ 488f32d5d43SBarry Smith colam = col*am; 489f32d5d43SBarry Smith for (i=0; i<am; i++) { /* over rows of C in those columns */ 490f32d5d43SBarry Smith r1 = r2 = r3 = r4 = 0.0; 491f32d5d43SBarry Smith n = a->i[i+1] - a->i[i]; 492f32d5d43SBarry Smith aj = a->j + a->i[i]; 493f32d5d43SBarry Smith aa = a->a + a->i[i]; 494f32d5d43SBarry Smith for (j=0; j<n; j++) { 495f32d5d43SBarry Smith r1 += (*aa)*b1[*aj]; 496f32d5d43SBarry Smith r2 += (*aa)*b2[*aj]; 497f32d5d43SBarry Smith r3 += (*aa)*b3[*aj]; 498f32d5d43SBarry Smith r4 += (*aa++)*b4[*aj++]; 4999123193aSHong Zhang } 500f32d5d43SBarry Smith c[colam + i] = r1; 501f32d5d43SBarry Smith c[colam + am + i] = r2; 502f32d5d43SBarry Smith c[colam + am2 + i] = r3; 503f32d5d43SBarry Smith c[colam + am3 + i] = r4; 504f32d5d43SBarry Smith } 505f32d5d43SBarry Smith b1 += bm4; 506f32d5d43SBarry Smith b2 += bm4; 507f32d5d43SBarry Smith b3 += bm4; 508f32d5d43SBarry Smith b4 += bm4; 509f32d5d43SBarry Smith } 510f32d5d43SBarry Smith for (;col<cn; col++){ /* over extra columns of C */ 511f32d5d43SBarry Smith for (i=0; i<am; i++) { /* over rows of C in those columns */ 512f32d5d43SBarry Smith r1 = 0.0; 513f32d5d43SBarry Smith n = a->i[i+1] - a->i[i]; 514f32d5d43SBarry Smith aj = a->j + a->i[i]; 515f32d5d43SBarry Smith aa = a->a + a->i[i]; 516f32d5d43SBarry Smith 517f32d5d43SBarry Smith for (j=0; j<n; j++) { 518f32d5d43SBarry Smith r1 += (*aa++)*b1[*aj++]; 519f32d5d43SBarry Smith } 520f32d5d43SBarry Smith c[col*am + i] = r1; 521f32d5d43SBarry Smith } 522f32d5d43SBarry Smith b1 += bm; 523f32d5d43SBarry Smith } 524dc0b31edSSatish Balay ierr = PetscLogFlops(cn*(2.0*a->nz));CHKERRQ(ierr); 525f32d5d43SBarry Smith ierr = MatRestoreArray(B,&b);CHKERRQ(ierr); 526f32d5d43SBarry Smith ierr = MatRestoreArray(C,&c);CHKERRQ(ierr); 527f32d5d43SBarry Smith ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 528f32d5d43SBarry Smith ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 529f32d5d43SBarry Smith PetscFunctionReturn(0); 530f32d5d43SBarry Smith } 531f32d5d43SBarry Smith 532f32d5d43SBarry Smith /* 5334324174eSBarry Smith Note very similar to MatMult_SeqAIJ(), should generate both codes from same base 534f32d5d43SBarry Smith */ 535f32d5d43SBarry Smith #undef __FUNCT__ 536f32d5d43SBarry Smith #define __FUNCT__ "MatMatMultNumericAdd_SeqAIJ_SeqDense" 537f32d5d43SBarry Smith PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat A,Mat B,Mat C) 538f32d5d43SBarry Smith { 539f32d5d43SBarry Smith Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data; 540f32d5d43SBarry Smith PetscErrorCode ierr; 541dd6ea824SBarry Smith PetscScalar *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4; 542dd6ea824SBarry Smith MatScalar *aa; 543d0f46423SBarry 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; 5444324174eSBarry Smith PetscInt am2 = 2*am, am3 = 3*am, bm4 = 4*bm,colam,*ridx; 545f32d5d43SBarry Smith 546f32d5d43SBarry Smith PetscFunctionBegin; 547f32d5d43SBarry Smith if (!cm || !cn) PetscFunctionReturn(0); 548f32d5d43SBarry Smith ierr = MatGetArray(B,&b);CHKERRQ(ierr); 549f32d5d43SBarry Smith ierr = MatGetArray(C,&c);CHKERRQ(ierr); 550f32d5d43SBarry Smith b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm; 5514324174eSBarry Smith 5524324174eSBarry Smith if (a->compressedrow.use){ /* use compressed row format */ 5534324174eSBarry Smith for (col=0; col<cn-4; col += 4){ /* over columns of C */ 5544324174eSBarry Smith colam = col*am; 5554324174eSBarry Smith arm = a->compressedrow.nrows; 5564324174eSBarry Smith ii = a->compressedrow.i; 5574324174eSBarry Smith ridx = a->compressedrow.rindex; 5584324174eSBarry Smith for (i=0; i<arm; i++) { /* over rows of C in those columns */ 5594324174eSBarry Smith r1 = r2 = r3 = r4 = 0.0; 5604324174eSBarry Smith n = ii[i+1] - ii[i]; 5614324174eSBarry Smith aj = a->j + ii[i]; 5624324174eSBarry Smith aa = a->a + ii[i]; 5634324174eSBarry Smith for (j=0; j<n; j++) { 5644324174eSBarry Smith r1 += (*aa)*b1[*aj]; 5654324174eSBarry Smith r2 += (*aa)*b2[*aj]; 5664324174eSBarry Smith r3 += (*aa)*b3[*aj]; 5674324174eSBarry Smith r4 += (*aa++)*b4[*aj++]; 5684324174eSBarry Smith } 5694324174eSBarry Smith c[colam + ridx[i]] += r1; 5704324174eSBarry Smith c[colam + am + ridx[i]] += r2; 5714324174eSBarry Smith c[colam + am2 + ridx[i]] += r3; 5724324174eSBarry Smith c[colam + am3 + ridx[i]] += r4; 5734324174eSBarry Smith } 5744324174eSBarry Smith b1 += bm4; 5754324174eSBarry Smith b2 += bm4; 5764324174eSBarry Smith b3 += bm4; 5774324174eSBarry Smith b4 += bm4; 5784324174eSBarry Smith } 5794324174eSBarry Smith for (;col<cn; col++){ /* over extra columns of C */ 5804324174eSBarry Smith colam = col*am; 5814324174eSBarry Smith arm = a->compressedrow.nrows; 5824324174eSBarry Smith ii = a->compressedrow.i; 5834324174eSBarry Smith ridx = a->compressedrow.rindex; 5844324174eSBarry Smith for (i=0; i<arm; i++) { /* over rows of C in those columns */ 5854324174eSBarry Smith r1 = 0.0; 5864324174eSBarry Smith n = ii[i+1] - ii[i]; 5874324174eSBarry Smith aj = a->j + ii[i]; 5884324174eSBarry Smith aa = a->a + ii[i]; 5894324174eSBarry Smith 5904324174eSBarry Smith for (j=0; j<n; j++) { 5914324174eSBarry Smith r1 += (*aa++)*b1[*aj++]; 5924324174eSBarry Smith } 5934324174eSBarry Smith c[col*am + ridx[i]] += r1; 5944324174eSBarry Smith } 5954324174eSBarry Smith b1 += bm; 5964324174eSBarry Smith } 5974324174eSBarry Smith } else { 598f32d5d43SBarry Smith for (col=0; col<cn-4; col += 4){ /* over columns of C */ 599f32d5d43SBarry Smith colam = col*am; 600f32d5d43SBarry Smith for (i=0; i<am; i++) { /* over rows of C in those columns */ 601f32d5d43SBarry Smith r1 = r2 = r3 = r4 = 0.0; 602f32d5d43SBarry Smith n = a->i[i+1] - a->i[i]; 603f32d5d43SBarry Smith aj = a->j + a->i[i]; 604f32d5d43SBarry Smith aa = a->a + a->i[i]; 605f32d5d43SBarry Smith for (j=0; j<n; j++) { 606f32d5d43SBarry Smith r1 += (*aa)*b1[*aj]; 607f32d5d43SBarry Smith r2 += (*aa)*b2[*aj]; 608f32d5d43SBarry Smith r3 += (*aa)*b3[*aj]; 609f32d5d43SBarry Smith r4 += (*aa++)*b4[*aj++]; 610f32d5d43SBarry Smith } 611f32d5d43SBarry Smith c[colam + i] += r1; 612f32d5d43SBarry Smith c[colam + am + i] += r2; 613f32d5d43SBarry Smith c[colam + am2 + i] += r3; 614f32d5d43SBarry Smith c[colam + am3 + i] += r4; 615f32d5d43SBarry Smith } 616f32d5d43SBarry Smith b1 += bm4; 617f32d5d43SBarry Smith b2 += bm4; 618f32d5d43SBarry Smith b3 += bm4; 619f32d5d43SBarry Smith b4 += bm4; 620f32d5d43SBarry Smith } 621f32d5d43SBarry Smith for (;col<cn; col++){ /* over extra columns of C */ 622f32d5d43SBarry Smith for (i=0; i<am; i++) { /* over rows of C in those columns */ 623f32d5d43SBarry Smith r1 = 0.0; 624f32d5d43SBarry Smith n = a->i[i+1] - a->i[i]; 625f32d5d43SBarry Smith aj = a->j + a->i[i]; 626f32d5d43SBarry Smith aa = a->a + a->i[i]; 627f32d5d43SBarry Smith 628f32d5d43SBarry Smith for (j=0; j<n; j++) { 629f32d5d43SBarry Smith r1 += (*aa++)*b1[*aj++]; 630f32d5d43SBarry Smith } 631f32d5d43SBarry Smith c[col*am + i] += r1; 632f32d5d43SBarry Smith } 633f32d5d43SBarry Smith b1 += bm; 634f32d5d43SBarry Smith } 6354324174eSBarry Smith } 636dc0b31edSSatish Balay ierr = PetscLogFlops(cn*2.0*a->nz);CHKERRQ(ierr); 637f32d5d43SBarry Smith ierr = MatRestoreArray(B,&b);CHKERRQ(ierr); 638f32d5d43SBarry Smith ierr = MatRestoreArray(C,&c);CHKERRQ(ierr); 6399123193aSHong Zhang PetscFunctionReturn(0); 6409123193aSHong Zhang } 641