xref: /petsc/src/mat/impls/aij/seq/matmatmult.c (revision e2cac8caf559aad8d0a261cc5db22e4573b50cbf)
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,&current_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;
19681d82fe4SHong Zhang   Mat            Bt;
19781d82fe4SHong Zhang   PetscInt       *bti,*btj;
19881d82fe4SHong Zhang 
19981d82fe4SHong Zhang   PetscFunctionBegin;
20081d82fe4SHong Zhang    /* create symbolic Bt */
20181d82fe4SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr);
20281d82fe4SHong Zhang   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,B->cmap->n,B->rmap->n,bti,btj,PETSC_NULL,&Bt);CHKERRQ(ierr);
20381d82fe4SHong Zhang 
20481d82fe4SHong Zhang   /* get symbolic C=A*Bt */
20581d82fe4SHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,Bt,fill,C);CHKERRQ(ierr);
20681d82fe4SHong Zhang 
20781d82fe4SHong Zhang   /* clean up */
20881d82fe4SHong Zhang   ierr = MatDestroy(&Bt);CHKERRQ(ierr);
20981d82fe4SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr);
21081d82fe4SHong Zhang 
21181d82fe4SHong Zhang #if defined(INEFFICIENT_ALGORITHM)
21281d82fe4SHong 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];
24281d82fe4SHong 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){
24581d82fe4SHong Zhang         while (acol[ka] < bcol[kb] && ka < anzi) ka++;
24681d82fe4SHong Zhang         if (ka == anzi) break;
24781d82fe4SHong Zhang         while (acol[ka] > bcol[kb] && kb < bnzj) kb++;
24881d82fe4SHong 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,&current_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
31281d82fe4SHong Zhang #endif
3135df89d91SHong Zhang   PetscFunctionReturn(0);
3145df89d91SHong Zhang }
3155df89d91SHong Zhang 
3166973769bSHong Zhang /* #define USE_ARRAY - for sparse dot product. Slower than !USE_ARRAY */
3175df89d91SHong Zhang #undef __FUNCT__
3185df89d91SHong Zhang #define __FUNCT__ "MatMatMultTransposeNumeric_SeqAIJ_SeqAIJ"
3195df89d91SHong Zhang PetscErrorCode MatMatMultTransposeNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
3205df89d91SHong Zhang {
3215df89d91SHong Zhang   PetscErrorCode ierr;
3225df89d91SHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data;
323*e2cac8caSJed Brown   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,anzi,bnzj,nexta,nextb,*acol,*bcol,brow;
32481d82fe4SHong Zhang   PetscInt       cm=C->rmap->n,*ci=c->i,*cj=c->j,i,j,cnzi,*ccol;
3255df89d91SHong Zhang   PetscLogDouble flops=0.0;
3266973769bSHong Zhang   MatScalar      *aa=a->a,*aval,*ba=b->a,*bval,*ca=c->a,*cval;
3276973769bSHong Zhang #if defined(USE_ARRAY)
3286973769bSHong Zhang   MatScalar *spdot;
3296973769bSHong Zhang #endif
3305df89d91SHong Zhang 
3315df89d91SHong Zhang   PetscFunctionBegin;
3326973769bSHong Zhang #if defined(USE_ARRAY)
3336973769bSHong Zhang   /* allocate an array for implementing sparse inner-product */
334*e2cac8caSJed Brown   ierr = PetscMalloc((A->cmap->n+1)*sizeof(MatScalar),&spdot);CHKERRQ(ierr);
335*e2cac8caSJed Brown   ierr = PetscMemzero(spdot,(A->cmap->n+1)*sizeof(MatScalar));CHKERRQ(ierr);
3366973769bSHong Zhang #endif
3376973769bSHong Zhang 
33881d82fe4SHong Zhang   /* clear old values in C */
33981d82fe4SHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
3405df89d91SHong Zhang 
3411fa1209cSHong Zhang   for (i=0; i<cm; i++) {
34281d82fe4SHong Zhang     anzi = ai[i+1] - ai[i];
34381d82fe4SHong Zhang     acol = aj + ai[i];
3446973769bSHong Zhang     aval = aa + ai[i];
3451fa1209cSHong Zhang     cnzi = ci[i+1] - ci[i];
3461fa1209cSHong Zhang     ccol = cj + ci[i];
3476973769bSHong Zhang     cval = ca + ci[i];
3481fa1209cSHong Zhang     for (j=0; j<cnzi; j++){
34981d82fe4SHong Zhang       brow = ccol[j];
35081d82fe4SHong Zhang       bnzj = bi[brow+1] - bi[brow];
35181d82fe4SHong Zhang       bcol = bj + bi[brow];
3526973769bSHong Zhang       bval = ba + bi[brow];
3536973769bSHong Zhang 
35481d82fe4SHong Zhang       /* perform sparse inner-product c(i,j)=A[i,:]*B[j,:]^T */
3556973769bSHong Zhang #if defined(USE_ARRAY)
3566973769bSHong Zhang       /* put ba to spdot array */
3576973769bSHong Zhang       for (nextb=0; nextb<bnzj; nextb++) spdot[bcol[nextb]] = bval[nextb];
3586973769bSHong Zhang       /* c(i,j)=A[i,:]*B[j,:]^T */
3596973769bSHong Zhang       for (nexta=0; nexta<anzi; nexta++){
3606973769bSHong Zhang         cval[j] += spdot[acol[nexta]]*aval[nexta];
3616973769bSHong Zhang       }
3626973769bSHong Zhang       /* zero spdot array */
3636973769bSHong Zhang       for (nextb=0; nextb<bnzj; nextb++) spdot[bcol[nextb]] = 0.0;
3646973769bSHong Zhang #else
36581d82fe4SHong Zhang       nexta = 0; nextb = 0;
36681d82fe4SHong Zhang       while (nexta<anzi && nextb<bnzj){
36781d82fe4SHong Zhang         while (acol[nexta] < bcol[nextb] && nexta < anzi) nexta++;
36881d82fe4SHong Zhang         if (nexta == anzi) break;
36981d82fe4SHong Zhang         while (acol[nexta] > bcol[nextb] && nextb < bnzj) nextb++;
37081d82fe4SHong Zhang         if (nextb == bnzj) break;
37181d82fe4SHong Zhang         if (acol[nexta] == bcol[nextb]){
3726973769bSHong Zhang           cval[j] += aval[nexta]*bval[nextb];
37381d82fe4SHong Zhang           nexta++; nextb++;
37481d82fe4SHong Zhang           flops += 2;
3751fa1209cSHong Zhang         }
3761fa1209cSHong Zhang       }
3776973769bSHong Zhang #endif
37881d82fe4SHong Zhang     }
37981d82fe4SHong Zhang   }
38081d82fe4SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
38181d82fe4SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
38281d82fe4SHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
3836973769bSHong Zhang #if defined(USE_ARRAY)
3846973769bSHong Zhang   ierr = PetscFree(spdot);
3856973769bSHong Zhang #endif
3865df89d91SHong Zhang   PetscFunctionReturn(0);
3875df89d91SHong Zhang }
3885df89d91SHong Zhang 
3895df89d91SHong Zhang #undef __FUNCT__
3905df89d91SHong Zhang #define __FUNCT__ "MatMatTransposeMult_SeqAIJ_SeqAIJ"
3915df89d91SHong Zhang PetscErrorCode MatMatTransposeMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) {
3925df89d91SHong Zhang   PetscErrorCode ierr;
3935df89d91SHong Zhang 
3945df89d91SHong Zhang   PetscFunctionBegin;
3955df89d91SHong Zhang   if (scall == MAT_INITIAL_MATRIX){
3965df89d91SHong Zhang     ierr = MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
3975df89d91SHong Zhang   }
3985df89d91SHong Zhang   ierr = MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr);
3995df89d91SHong Zhang   PetscFunctionReturn(0);
4005df89d91SHong Zhang }
4015df89d91SHong Zhang 
4025df89d91SHong Zhang #undef __FUNCT__
4035df89d91SHong Zhang #define __FUNCT__ "MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ"
4045df89d91SHong Zhang PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
4055df89d91SHong Zhang {
406bc011b1eSHong Zhang   PetscErrorCode ierr;
407bc011b1eSHong Zhang   Mat            At;
40838baddfdSBarry Smith   PetscInt       *ati,*atj;
409bc011b1eSHong Zhang 
410bc011b1eSHong Zhang   PetscFunctionBegin;
411bc011b1eSHong Zhang   /* create symbolic At */
412bc011b1eSHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr);
413d0f46423SBarry Smith   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,A->cmap->n,A->rmap->n,ati,atj,PETSC_NULL,&At);CHKERRQ(ierr);
414bc011b1eSHong Zhang 
415bc011b1eSHong Zhang   /* get symbolic C=At*B */
416bc011b1eSHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(At,B,fill,C);CHKERRQ(ierr);
417bc011b1eSHong Zhang 
418bc011b1eSHong Zhang   /* clean up */
4196bf464f9SBarry Smith   ierr = MatDestroy(&At);CHKERRQ(ierr);
420bc011b1eSHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr);
421bc011b1eSHong Zhang   PetscFunctionReturn(0);
422bc011b1eSHong Zhang }
423bc011b1eSHong Zhang 
424bc011b1eSHong Zhang #undef __FUNCT__
4255df89d91SHong Zhang #define __FUNCT__ "MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ"
4265df89d91SHong Zhang PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
427bc011b1eSHong Zhang {
428bc011b1eSHong Zhang   PetscErrorCode ierr;
4290fbc74f4SHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data;
430d0f46423SBarry Smith   PetscInt       am=A->rmap->n,anzi,*ai=a->i,*aj=a->j,*bi=b->i,*bj,bnzi,nextb;
431d0f46423SBarry Smith   PetscInt       cm=C->rmap->n,*ci=c->i,*cj=c->j,crow,*cjj,i,j,k;
432d13dce4bSSatish Balay   PetscLogDouble flops=0.0;
4330fbc74f4SHong Zhang   MatScalar      *aa=a->a,*ba,*ca=c->a,*caj;
434bc011b1eSHong Zhang 
435bc011b1eSHong Zhang   PetscFunctionBegin;
436bc011b1eSHong Zhang   /* clear old values in C */
437bc011b1eSHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
438bc011b1eSHong Zhang 
439bc011b1eSHong Zhang   /* compute A^T*B using outer product (A^T)[:,i]*B[i,:] */
440bc011b1eSHong Zhang   for (i=0;i<am;i++) {
441bc011b1eSHong Zhang     bj   = b->j + bi[i];
442bc011b1eSHong Zhang     ba   = b->a + bi[i];
443bc011b1eSHong Zhang     bnzi = bi[i+1] - bi[i];
444bc011b1eSHong Zhang     anzi = ai[i+1] - ai[i];
445bc011b1eSHong Zhang     for (j=0; j<anzi; j++) {
446bc011b1eSHong Zhang       nextb = 0;
4470fbc74f4SHong Zhang       crow  = *aj++;
448bc011b1eSHong Zhang       cjj   = cj + ci[crow];
449bc011b1eSHong Zhang       caj   = ca + ci[crow];
450bc011b1eSHong Zhang       /* perform sparse axpy operation.  Note cjj includes bj. */
451bc011b1eSHong Zhang       for (k=0; nextb<bnzi; k++) {
4520fbc74f4SHong Zhang         if (cjj[k] == *(bj+nextb)) { /* ccol == bcol */
4530fbc74f4SHong Zhang           caj[k] += (*aa)*(*(ba+nextb));
454bc011b1eSHong Zhang           nextb++;
455bc011b1eSHong Zhang         }
456bc011b1eSHong Zhang       }
457bc011b1eSHong Zhang       flops += 2*bnzi;
4580fbc74f4SHong Zhang       aa++;
459bc011b1eSHong Zhang     }
460bc011b1eSHong Zhang   }
461bc011b1eSHong Zhang 
462bc011b1eSHong Zhang   /* Assemble the final matrix and clean up */
463bc011b1eSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
464bc011b1eSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
465bc011b1eSHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
466bc011b1eSHong Zhang   PetscFunctionReturn(0);
467bc011b1eSHong Zhang }
4689123193aSHong Zhang 
469ec01deb9SMatthew Knepley EXTERN_C_BEGIN
4709123193aSHong Zhang #undef __FUNCT__
4719123193aSHong Zhang #define __FUNCT__ "MatMatMult_SeqAIJ_SeqDense"
4729123193aSHong Zhang PetscErrorCode MatMatMult_SeqAIJ_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
4739123193aSHong Zhang {
4749123193aSHong Zhang   PetscErrorCode ierr;
4759123193aSHong Zhang 
4769123193aSHong Zhang   PetscFunctionBegin;
4779123193aSHong Zhang   if (scall == MAT_INITIAL_MATRIX){
4789123193aSHong Zhang     ierr = MatMatMultSymbolic_SeqAIJ_SeqDense(A,B,fill,C);CHKERRQ(ierr);
4799123193aSHong Zhang   }
4809123193aSHong Zhang   ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,B,*C);CHKERRQ(ierr);
4819123193aSHong Zhang   PetscFunctionReturn(0);
4829123193aSHong Zhang }
483ec01deb9SMatthew Knepley EXTERN_C_END
484edd81eecSMatthew Knepley 
4859123193aSHong Zhang #undef __FUNCT__
4869123193aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqDense"
4879123193aSHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
4889123193aSHong Zhang {
4899123193aSHong Zhang   PetscErrorCode ierr;
4909123193aSHong Zhang 
4919123193aSHong Zhang   PetscFunctionBegin;
4925a586d82SBarry Smith   ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,0.0,C);CHKERRQ(ierr);
4939123193aSHong Zhang   PetscFunctionReturn(0);
4949123193aSHong Zhang }
4959123193aSHong Zhang 
4969123193aSHong Zhang #undef __FUNCT__
4979123193aSHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqDense"
4989123193aSHong Zhang PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
4999123193aSHong Zhang {
500f32d5d43SBarry Smith   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
5019123193aSHong Zhang   PetscErrorCode ierr;
502dd6ea824SBarry Smith   PetscScalar    *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4;
503dd6ea824SBarry Smith   MatScalar      *aa;
504d0f46423SBarry Smith   PetscInt       cm=C->rmap->n, cn=B->cmap->n, bm=B->rmap->n, col, i,j,n,*aj, am = A->rmap->n;
505f32d5d43SBarry Smith   PetscInt       am2 = 2*am, am3 = 3*am,  bm4 = 4*bm,colam;
5069123193aSHong Zhang 
5079123193aSHong Zhang   PetscFunctionBegin;
508f32d5d43SBarry Smith   if (!cm || !cn) PetscFunctionReturn(0);
509e32f2f54SBarry 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);
510e32f2f54SBarry 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);
511e32f2f54SBarry 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);
512f32d5d43SBarry Smith   ierr = MatGetArray(B,&b);CHKERRQ(ierr);
513f32d5d43SBarry Smith   ierr = MatGetArray(C,&c);CHKERRQ(ierr);
514f32d5d43SBarry Smith   b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
515f32d5d43SBarry Smith   for (col=0; col<cn-4; col += 4){  /* over columns of C */
516f32d5d43SBarry Smith     colam = col*am;
517f32d5d43SBarry Smith     for (i=0; i<am; i++) {        /* over rows of C in those columns */
518f32d5d43SBarry Smith       r1 = r2 = r3 = r4 = 0.0;
519f32d5d43SBarry Smith       n   = a->i[i+1] - a->i[i];
520f32d5d43SBarry Smith       aj  = a->j + a->i[i];
521f32d5d43SBarry Smith       aa  = a->a + a->i[i];
522f32d5d43SBarry Smith       for (j=0; j<n; j++) {
523f32d5d43SBarry Smith         r1 += (*aa)*b1[*aj];
524f32d5d43SBarry Smith         r2 += (*aa)*b2[*aj];
525f32d5d43SBarry Smith         r3 += (*aa)*b3[*aj];
526f32d5d43SBarry Smith         r4 += (*aa++)*b4[*aj++];
5279123193aSHong Zhang       }
528f32d5d43SBarry Smith       c[colam + i]       = r1;
529f32d5d43SBarry Smith       c[colam + am + i]  = r2;
530f32d5d43SBarry Smith       c[colam + am2 + i] = r3;
531f32d5d43SBarry Smith       c[colam + am3 + i] = r4;
532f32d5d43SBarry Smith     }
533f32d5d43SBarry Smith     b1 += bm4;
534f32d5d43SBarry Smith     b2 += bm4;
535f32d5d43SBarry Smith     b3 += bm4;
536f32d5d43SBarry Smith     b4 += bm4;
537f32d5d43SBarry Smith   }
538f32d5d43SBarry Smith   for (;col<cn; col++){     /* over extra columns of C */
539f32d5d43SBarry Smith     for (i=0; i<am; i++) {  /* over rows of C in those columns */
540f32d5d43SBarry Smith       r1 = 0.0;
541f32d5d43SBarry Smith       n   = a->i[i+1] - a->i[i];
542f32d5d43SBarry Smith       aj  = a->j + a->i[i];
543f32d5d43SBarry Smith       aa  = a->a + a->i[i];
544f32d5d43SBarry Smith 
545f32d5d43SBarry Smith       for (j=0; j<n; j++) {
546f32d5d43SBarry Smith         r1 += (*aa++)*b1[*aj++];
547f32d5d43SBarry Smith       }
548f32d5d43SBarry Smith       c[col*am + i]     = r1;
549f32d5d43SBarry Smith     }
550f32d5d43SBarry Smith     b1 += bm;
551f32d5d43SBarry Smith   }
552dc0b31edSSatish Balay   ierr = PetscLogFlops(cn*(2.0*a->nz));CHKERRQ(ierr);
553f32d5d43SBarry Smith   ierr = MatRestoreArray(B,&b);CHKERRQ(ierr);
554f32d5d43SBarry Smith   ierr = MatRestoreArray(C,&c);CHKERRQ(ierr);
555f32d5d43SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
556f32d5d43SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
557f32d5d43SBarry Smith   PetscFunctionReturn(0);
558f32d5d43SBarry Smith }
559f32d5d43SBarry Smith 
560f32d5d43SBarry Smith /*
5614324174eSBarry Smith    Note very similar to MatMult_SeqAIJ(), should generate both codes from same base
562f32d5d43SBarry Smith */
563f32d5d43SBarry Smith #undef __FUNCT__
564f32d5d43SBarry Smith #define __FUNCT__ "MatMatMultNumericAdd_SeqAIJ_SeqDense"
565f32d5d43SBarry Smith PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
566f32d5d43SBarry Smith {
567f32d5d43SBarry Smith   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
568f32d5d43SBarry Smith   PetscErrorCode ierr;
569dd6ea824SBarry Smith   PetscScalar    *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4;
570dd6ea824SBarry Smith   MatScalar      *aa;
571d0f46423SBarry 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;
5724324174eSBarry Smith   PetscInt       am2 = 2*am, am3 = 3*am,  bm4 = 4*bm,colam,*ridx;
573f32d5d43SBarry Smith 
574f32d5d43SBarry Smith   PetscFunctionBegin;
575f32d5d43SBarry Smith   if (!cm || !cn) PetscFunctionReturn(0);
576f32d5d43SBarry Smith   ierr = MatGetArray(B,&b);CHKERRQ(ierr);
577f32d5d43SBarry Smith   ierr = MatGetArray(C,&c);CHKERRQ(ierr);
578f32d5d43SBarry Smith   b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
5794324174eSBarry Smith 
5804324174eSBarry Smith   if (a->compressedrow.use){ /* use compressed row format */
5814324174eSBarry Smith     for (col=0; col<cn-4; col += 4){  /* over columns of C */
5824324174eSBarry Smith       colam = col*am;
5834324174eSBarry Smith       arm   = a->compressedrow.nrows;
5844324174eSBarry Smith       ii    = a->compressedrow.i;
5854324174eSBarry Smith       ridx  = a->compressedrow.rindex;
5864324174eSBarry Smith       for (i=0; i<arm; i++) {        /* over rows of C in those columns */
5874324174eSBarry Smith 	r1 = r2 = r3 = r4 = 0.0;
5884324174eSBarry Smith 	n   = ii[i+1] - ii[i];
5894324174eSBarry Smith 	aj  = a->j + ii[i];
5904324174eSBarry Smith 	aa  = a->a + ii[i];
5914324174eSBarry Smith 	for (j=0; j<n; j++) {
5924324174eSBarry Smith 	  r1 += (*aa)*b1[*aj];
5934324174eSBarry Smith 	  r2 += (*aa)*b2[*aj];
5944324174eSBarry Smith 	  r3 += (*aa)*b3[*aj];
5954324174eSBarry Smith 	  r4 += (*aa++)*b4[*aj++];
5964324174eSBarry Smith 	}
5974324174eSBarry Smith 	c[colam       + ridx[i]] += r1;
5984324174eSBarry Smith 	c[colam + am  + ridx[i]] += r2;
5994324174eSBarry Smith 	c[colam + am2 + ridx[i]] += r3;
6004324174eSBarry Smith 	c[colam + am3 + ridx[i]] += r4;
6014324174eSBarry Smith       }
6024324174eSBarry Smith       b1 += bm4;
6034324174eSBarry Smith       b2 += bm4;
6044324174eSBarry Smith       b3 += bm4;
6054324174eSBarry Smith       b4 += bm4;
6064324174eSBarry Smith     }
6074324174eSBarry Smith     for (;col<cn; col++){     /* over extra columns of C */
6084324174eSBarry Smith       colam = col*am;
6094324174eSBarry Smith       arm   = a->compressedrow.nrows;
6104324174eSBarry Smith       ii    = a->compressedrow.i;
6114324174eSBarry Smith       ridx  = a->compressedrow.rindex;
6124324174eSBarry Smith       for (i=0; i<arm; i++) {  /* over rows of C in those columns */
6134324174eSBarry Smith 	r1 = 0.0;
6144324174eSBarry Smith 	n   = ii[i+1] - ii[i];
6154324174eSBarry Smith 	aj  = a->j + ii[i];
6164324174eSBarry Smith 	aa  = a->a + ii[i];
6174324174eSBarry Smith 
6184324174eSBarry Smith 	for (j=0; j<n; j++) {
6194324174eSBarry Smith 	  r1 += (*aa++)*b1[*aj++];
6204324174eSBarry Smith 	}
6214324174eSBarry Smith 	c[col*am + ridx[i]] += r1;
6224324174eSBarry Smith       }
6234324174eSBarry Smith       b1 += bm;
6244324174eSBarry Smith     }
6254324174eSBarry Smith   } else {
626f32d5d43SBarry Smith     for (col=0; col<cn-4; col += 4){  /* over columns of C */
627f32d5d43SBarry Smith       colam = col*am;
628f32d5d43SBarry Smith       for (i=0; i<am; i++) {        /* over rows of C in those columns */
629f32d5d43SBarry Smith 	r1 = r2 = r3 = r4 = 0.0;
630f32d5d43SBarry Smith 	n   = a->i[i+1] - a->i[i];
631f32d5d43SBarry Smith 	aj  = a->j + a->i[i];
632f32d5d43SBarry Smith 	aa  = a->a + a->i[i];
633f32d5d43SBarry Smith 	for (j=0; j<n; j++) {
634f32d5d43SBarry Smith 	  r1 += (*aa)*b1[*aj];
635f32d5d43SBarry Smith 	  r2 += (*aa)*b2[*aj];
636f32d5d43SBarry Smith 	  r3 += (*aa)*b3[*aj];
637f32d5d43SBarry Smith 	  r4 += (*aa++)*b4[*aj++];
638f32d5d43SBarry Smith 	}
639f32d5d43SBarry Smith 	c[colam + i]       += r1;
640f32d5d43SBarry Smith 	c[colam + am + i]  += r2;
641f32d5d43SBarry Smith 	c[colam + am2 + i] += r3;
642f32d5d43SBarry Smith 	c[colam + am3 + i] += r4;
643f32d5d43SBarry Smith       }
644f32d5d43SBarry Smith       b1 += bm4;
645f32d5d43SBarry Smith       b2 += bm4;
646f32d5d43SBarry Smith       b3 += bm4;
647f32d5d43SBarry Smith       b4 += bm4;
648f32d5d43SBarry Smith     }
649f32d5d43SBarry Smith     for (;col<cn; col++){     /* over extra columns of C */
650f32d5d43SBarry Smith       for (i=0; i<am; i++) {  /* over rows of C in those columns */
651f32d5d43SBarry Smith 	r1 = 0.0;
652f32d5d43SBarry Smith 	n   = a->i[i+1] - a->i[i];
653f32d5d43SBarry Smith 	aj  = a->j + a->i[i];
654f32d5d43SBarry Smith 	aa  = a->a + a->i[i];
655f32d5d43SBarry Smith 
656f32d5d43SBarry Smith 	for (j=0; j<n; j++) {
657f32d5d43SBarry Smith 	  r1 += (*aa++)*b1[*aj++];
658f32d5d43SBarry Smith 	}
659f32d5d43SBarry Smith 	c[col*am + i]     += r1;
660f32d5d43SBarry Smith       }
661f32d5d43SBarry Smith       b1 += bm;
662f32d5d43SBarry Smith     }
6634324174eSBarry Smith   }
664dc0b31edSSatish Balay   ierr = PetscLogFlops(cn*2.0*a->nz);CHKERRQ(ierr);
665f32d5d43SBarry Smith   ierr = MatRestoreArray(B,&b);CHKERRQ(ierr);
666f32d5d43SBarry Smith   ierr = MatRestoreArray(C,&c);CHKERRQ(ierr);
6679123193aSHong Zhang   PetscFunctionReturn(0);
6689123193aSHong Zhang }
669