xref: /petsc/src/mat/impls/aij/seq/matmatmult.c (revision 1fa1209c6352ae9e10f5ed13f1a182607a450bef)
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;
196*1fa1209cSHong Zhang   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
197*1fa1209cSHong Zhang   Mat_SeqAIJ         *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
198*1fa1209cSHong Zhang   PetscInt           *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*ci,*cj,*acol,*bcol;
199*1fa1209cSHong Zhang   PetscInt           am=A->rmap->N,bm=B->rmap->N;
200*1fa1209cSHong Zhang   PetscInt           i,j,anzi,bnzj,cnzi,nlnk,*lnk,nspacedouble=0,ka,kb,index[1];
201*1fa1209cSHong Zhang   MatScalar          *ca;
202*1fa1209cSHong Zhang   PetscBT            lnkbt;
203*1fa1209cSHong Zhang   PetscReal          afill;
2045df89d91SHong Zhang 
205*1fa1209cSHong Zhang   PetscFunctionBegin;
206*1fa1209cSHong Zhang   /* Allocate row pointer array ci  */
207*1fa1209cSHong Zhang   ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
208*1fa1209cSHong Zhang   ci[0] = 0;
209*1fa1209cSHong Zhang 
210*1fa1209cSHong Zhang   /* Create and initialize a linked list for C columns */
211*1fa1209cSHong Zhang   nlnk = bm+1;
212*1fa1209cSHong Zhang   ierr = PetscLLCreate(bm,bm,nlnk,lnk,lnkbt);CHKERRQ(ierr);
213*1fa1209cSHong Zhang 
214*1fa1209cSHong Zhang   /* Initial FreeSpace with size fill*(nnz(A)+nnz(B)) */
215*1fa1209cSHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr);
216*1fa1209cSHong Zhang   current_space = free_space;
217*1fa1209cSHong Zhang 
218*1fa1209cSHong Zhang   /* Determine symbolic info for each row of the product A*B^T: */
219*1fa1209cSHong Zhang   for (i=0; i<am; i++) {
220*1fa1209cSHong Zhang     anzi = ai[i+1] - ai[i];
221*1fa1209cSHong Zhang     cnzi = 0;
222*1fa1209cSHong Zhang     acol = aj + ai[i];
223*1fa1209cSHong Zhang     for (j=0; j<bm; j++){
224*1fa1209cSHong Zhang       bnzj = bi[j+1] - bi[j];
225*1fa1209cSHong Zhang       bcol= bj + bi[j];
226*1fa1209cSHong Zhang       /* sparse inner-product A[i,:] * B[j,:]^T */
227*1fa1209cSHong Zhang       ka = 0; kb = 0;
228*1fa1209cSHong Zhang       while (ka < anzi && kb < bnzj){
229*1fa1209cSHong Zhang         while (acol[ka] < bcol[kb] && ka < anzi){
230*1fa1209cSHong Zhang           ka++;
231*1fa1209cSHong Zhang         }
232*1fa1209cSHong Zhang         while (acol[ka] > bcol[kb] && kb < bnzj){
233*1fa1209cSHong Zhang           kb++;
234*1fa1209cSHong Zhang         }
235*1fa1209cSHong Zhang         if (acol[ka] == bcol[kb]){ /* add nonzero c(i,j) to lnk */
236*1fa1209cSHong Zhang           index[0] = j;
237*1fa1209cSHong Zhang           ierr = PetscLLAdd(1,index,bm,nlnk,lnk,lnkbt);CHKERRQ(ierr);
238*1fa1209cSHong Zhang           cnzi++;
239*1fa1209cSHong Zhang           break;
240*1fa1209cSHong Zhang         }
241*1fa1209cSHong Zhang       }
242*1fa1209cSHong Zhang     }
243*1fa1209cSHong Zhang 
244*1fa1209cSHong Zhang     /* If free space is not available, make more free space */
245*1fa1209cSHong Zhang     /* Double the amount of total space in the list */
246*1fa1209cSHong Zhang     if (current_space->local_remaining<cnzi) {
247*1fa1209cSHong Zhang       ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,&current_space);CHKERRQ(ierr);
248*1fa1209cSHong Zhang       nspacedouble++;
249*1fa1209cSHong Zhang     }
250*1fa1209cSHong Zhang 
251*1fa1209cSHong Zhang     /* Copy data into free space, then initialize lnk */
252*1fa1209cSHong Zhang     ierr = PetscLLClean(bm,bm,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
253*1fa1209cSHong Zhang     current_space->array           += cnzi;
254*1fa1209cSHong Zhang     current_space->local_used      += cnzi;
255*1fa1209cSHong Zhang     current_space->local_remaining -= cnzi;
256*1fa1209cSHong Zhang 
257*1fa1209cSHong Zhang     ci[i+1] = ci[i] + cnzi;
258*1fa1209cSHong Zhang   }
259*1fa1209cSHong Zhang 
260*1fa1209cSHong Zhang 
261*1fa1209cSHong Zhang   /* Column indices are in the list of free space.
262*1fa1209cSHong Zhang      Allocate array cj, copy column indices to cj, and destroy list of free space */
263*1fa1209cSHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
264*1fa1209cSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
265*1fa1209cSHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
266*1fa1209cSHong Zhang 
267*1fa1209cSHong Zhang   /* Allocate space for ca */
268*1fa1209cSHong Zhang   ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
269*1fa1209cSHong Zhang   ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr);
270*1fa1209cSHong Zhang 
271*1fa1209cSHong Zhang   /* put together the new symbolic matrix */
272*1fa1209cSHong Zhang   ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bm,ci,cj,ca,C);CHKERRQ(ierr);
273*1fa1209cSHong Zhang 
274*1fa1209cSHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
275*1fa1209cSHong Zhang   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
276*1fa1209cSHong Zhang   c = (Mat_SeqAIJ *)((*C)->data);
277*1fa1209cSHong Zhang   c->free_a   = PETSC_TRUE;
278*1fa1209cSHong Zhang   c->free_ij  = PETSC_TRUE;
279*1fa1209cSHong Zhang   c->nonew    = 0;
280*1fa1209cSHong Zhang 
281*1fa1209cSHong Zhang   /* set MatInfo */
282*1fa1209cSHong Zhang   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
283*1fa1209cSHong Zhang   if (afill < 1.0) afill = 1.0;
284*1fa1209cSHong Zhang   c->maxnz                     = ci[am];
285*1fa1209cSHong Zhang   c->nz                        = ci[am];
286*1fa1209cSHong Zhang   (*C)->info.mallocs           = nspacedouble;
287*1fa1209cSHong Zhang   (*C)->info.fill_ratio_given  = fill;
288*1fa1209cSHong Zhang   (*C)->info.fill_ratio_needed = afill;
289*1fa1209cSHong Zhang 
290*1fa1209cSHong Zhang #if defined(PETSC_USE_INFO)
291*1fa1209cSHong Zhang   if (ci[am]) {
292*1fa1209cSHong Zhang     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
293*1fa1209cSHong Zhang     ierr = PetscInfo1((*C),"Use MatMatMultTranspose(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
294*1fa1209cSHong Zhang   } else {
295*1fa1209cSHong Zhang     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
296*1fa1209cSHong Zhang   }
297*1fa1209cSHong Zhang #endif
2985df89d91SHong Zhang   PetscFunctionReturn(0);
2995df89d91SHong Zhang }
3005df89d91SHong Zhang 
3015df89d91SHong Zhang #undef __FUNCT__
3025df89d91SHong Zhang #define __FUNCT__ "MatMatMultTransposeNumeric_SeqAIJ_SeqAIJ"
3035df89d91SHong Zhang PetscErrorCode MatMatMultTransposeNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
3045df89d91SHong Zhang {
305*1fa1209cSHong Zhang #if defined(TMP)
3065df89d91SHong Zhang   PetscErrorCode ierr;
3075df89d91SHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data;
3085df89d91SHong Zhang   PetscInt       am=A->rmap->n,anzi,*ai=a->i,*aj=a->j,*bi=b->i,*bj,bnzi,nextb;
309*1fa1209cSHong Zhang   PetscInt       cm=C->rmap->n,*ci=c->i,*cj=c->j,i,j,k,cnzi,*ccol;
3105df89d91SHong Zhang   PetscLogDouble flops=0.0;
3115df89d91SHong Zhang   MatScalar      *aa=a->a,*ba,*ca=c->a,*caj;
312*1fa1209cSHong Zhang #endif
3135df89d91SHong Zhang 
3145df89d91SHong Zhang   PetscFunctionBegin;
315*1fa1209cSHong Zhang #if defined(TMP)
316*1fa1209cSHong Zhang   ierr = MatView(A, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
317*1fa1209cSHong Zhang   ierr = MatView(B, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
3185df89d91SHong Zhang 
319*1fa1209cSHong Zhang   for (i=0; i<cm; i++) {
320*1fa1209cSHong Zhang     printf("C row %d\n",i);
321*1fa1209cSHong Zhang     cnzi = ci[i+1] - ci[i];
322*1fa1209cSHong Zhang     ccol = cj + ci[i];
323*1fa1209cSHong Zhang     for (j=0; j<cnzi; j++){
324*1fa1209cSHong Zhang       printf(" %d, ",ccol[j]);
325*1fa1209cSHong Zhang     }
326*1fa1209cSHong Zhang     printf("   cnzi %d\n",cnzi);
327*1fa1209cSHong Zhang   }
328*1fa1209cSHong Zhang #endif
329*1fa1209cSHong Zhang   SETERRQ(PETSC_COMM_SELF,0,"MatMatMultTransposeNumeric_SeqAIJ_SeqAIJ Not done yet");
3305df89d91SHong Zhang   PetscFunctionReturn(0);
3315df89d91SHong Zhang }
3325df89d91SHong Zhang 
3335df89d91SHong Zhang #undef __FUNCT__
3345df89d91SHong Zhang #define __FUNCT__ "MatMatTransposeMult_SeqAIJ_SeqAIJ"
3355df89d91SHong Zhang PetscErrorCode MatMatTransposeMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) {
3365df89d91SHong Zhang   PetscErrorCode ierr;
3375df89d91SHong Zhang 
3385df89d91SHong Zhang   PetscFunctionBegin;
3395df89d91SHong Zhang   if (scall == MAT_INITIAL_MATRIX){
3405df89d91SHong Zhang     ierr = MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
3415df89d91SHong Zhang   }
3425df89d91SHong Zhang   ierr = MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr);
3435df89d91SHong Zhang   PetscFunctionReturn(0);
3445df89d91SHong Zhang }
3455df89d91SHong Zhang 
3465df89d91SHong Zhang #undef __FUNCT__
3475df89d91SHong Zhang #define __FUNCT__ "MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ"
3485df89d91SHong Zhang PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
3495df89d91SHong Zhang {
350bc011b1eSHong Zhang   PetscErrorCode ierr;
351bc011b1eSHong Zhang   Mat            At;
35238baddfdSBarry Smith   PetscInt       *ati,*atj;
353bc011b1eSHong Zhang 
354bc011b1eSHong Zhang   PetscFunctionBegin;
355bc011b1eSHong Zhang   /* create symbolic At */
356bc011b1eSHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr);
357d0f46423SBarry Smith   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,A->cmap->n,A->rmap->n,ati,atj,PETSC_NULL,&At);CHKERRQ(ierr);
358bc011b1eSHong Zhang 
359bc011b1eSHong Zhang   /* get symbolic C=At*B */
360bc011b1eSHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(At,B,fill,C);CHKERRQ(ierr);
361bc011b1eSHong Zhang 
362bc011b1eSHong Zhang   /* clean up */
3636bf464f9SBarry Smith   ierr = MatDestroy(&At);CHKERRQ(ierr);
364bc011b1eSHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr);
365bc011b1eSHong Zhang   PetscFunctionReturn(0);
366bc011b1eSHong Zhang }
367bc011b1eSHong Zhang 
368bc011b1eSHong Zhang #undef __FUNCT__
3695df89d91SHong Zhang #define __FUNCT__ "MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ"
3705df89d91SHong Zhang PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
371bc011b1eSHong Zhang {
372bc011b1eSHong Zhang   PetscErrorCode ierr;
3730fbc74f4SHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data;
374d0f46423SBarry Smith   PetscInt       am=A->rmap->n,anzi,*ai=a->i,*aj=a->j,*bi=b->i,*bj,bnzi,nextb;
375d0f46423SBarry Smith   PetscInt       cm=C->rmap->n,*ci=c->i,*cj=c->j,crow,*cjj,i,j,k;
376d13dce4bSSatish Balay   PetscLogDouble flops=0.0;
3770fbc74f4SHong Zhang   MatScalar      *aa=a->a,*ba,*ca=c->a,*caj;
378bc011b1eSHong Zhang 
379bc011b1eSHong Zhang   PetscFunctionBegin;
380bc011b1eSHong Zhang   /* clear old values in C */
381bc011b1eSHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
382bc011b1eSHong Zhang 
383bc011b1eSHong Zhang   /* compute A^T*B using outer product (A^T)[:,i]*B[i,:] */
384bc011b1eSHong Zhang   for (i=0;i<am;i++) {
385bc011b1eSHong Zhang     bj   = b->j + bi[i];
386bc011b1eSHong Zhang     ba   = b->a + bi[i];
387bc011b1eSHong Zhang     bnzi = bi[i+1] - bi[i];
388bc011b1eSHong Zhang     anzi = ai[i+1] - ai[i];
389bc011b1eSHong Zhang     for (j=0; j<anzi; j++) {
390bc011b1eSHong Zhang       nextb = 0;
3910fbc74f4SHong Zhang       crow  = *aj++;
392bc011b1eSHong Zhang       cjj   = cj + ci[crow];
393bc011b1eSHong Zhang       caj   = ca + ci[crow];
394bc011b1eSHong Zhang       /* perform sparse axpy operation.  Note cjj includes bj. */
395bc011b1eSHong Zhang       for (k=0; nextb<bnzi; k++) {
3960fbc74f4SHong Zhang         if (cjj[k] == *(bj+nextb)) { /* ccol == bcol */
3970fbc74f4SHong Zhang           caj[k] += (*aa)*(*(ba+nextb));
398bc011b1eSHong Zhang           nextb++;
399bc011b1eSHong Zhang         }
400bc011b1eSHong Zhang       }
401bc011b1eSHong Zhang       flops += 2*bnzi;
4020fbc74f4SHong Zhang       aa++;
403bc011b1eSHong Zhang     }
404bc011b1eSHong Zhang   }
405bc011b1eSHong Zhang 
406bc011b1eSHong Zhang   /* Assemble the final matrix and clean up */
407bc011b1eSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
408bc011b1eSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
409bc011b1eSHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
410bc011b1eSHong Zhang   PetscFunctionReturn(0);
411bc011b1eSHong Zhang }
4129123193aSHong Zhang 
413ec01deb9SMatthew Knepley EXTERN_C_BEGIN
4149123193aSHong Zhang #undef __FUNCT__
4159123193aSHong Zhang #define __FUNCT__ "MatMatMult_SeqAIJ_SeqDense"
4169123193aSHong Zhang PetscErrorCode MatMatMult_SeqAIJ_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
4179123193aSHong Zhang {
4189123193aSHong Zhang   PetscErrorCode ierr;
4199123193aSHong Zhang 
4209123193aSHong Zhang   PetscFunctionBegin;
4219123193aSHong Zhang   if (scall == MAT_INITIAL_MATRIX){
4229123193aSHong Zhang     ierr = MatMatMultSymbolic_SeqAIJ_SeqDense(A,B,fill,C);CHKERRQ(ierr);
4239123193aSHong Zhang   }
4249123193aSHong Zhang   ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,B,*C);CHKERRQ(ierr);
4259123193aSHong Zhang   PetscFunctionReturn(0);
4269123193aSHong Zhang }
427ec01deb9SMatthew Knepley EXTERN_C_END
428edd81eecSMatthew Knepley 
4299123193aSHong Zhang #undef __FUNCT__
4309123193aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqDense"
4319123193aSHong Zhang PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
4329123193aSHong Zhang {
4339123193aSHong Zhang   PetscErrorCode ierr;
4349123193aSHong Zhang 
4359123193aSHong Zhang   PetscFunctionBegin;
4365a586d82SBarry Smith   ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,0.0,C);CHKERRQ(ierr);
4379123193aSHong Zhang   PetscFunctionReturn(0);
4389123193aSHong Zhang }
4399123193aSHong Zhang 
4409123193aSHong Zhang #undef __FUNCT__
4419123193aSHong Zhang #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqDense"
4429123193aSHong Zhang PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
4439123193aSHong Zhang {
444f32d5d43SBarry Smith   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
4459123193aSHong Zhang   PetscErrorCode ierr;
446dd6ea824SBarry Smith   PetscScalar    *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4;
447dd6ea824SBarry Smith   MatScalar      *aa;
448d0f46423SBarry Smith   PetscInt       cm=C->rmap->n, cn=B->cmap->n, bm=B->rmap->n, col, i,j,n,*aj, am = A->rmap->n;
449f32d5d43SBarry Smith   PetscInt       am2 = 2*am, am3 = 3*am,  bm4 = 4*bm,colam;
4509123193aSHong Zhang 
4519123193aSHong Zhang   PetscFunctionBegin;
452f32d5d43SBarry Smith   if (!cm || !cn) PetscFunctionReturn(0);
453e32f2f54SBarry 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);
454e32f2f54SBarry 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);
455e32f2f54SBarry 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);
456f32d5d43SBarry Smith   ierr = MatGetArray(B,&b);CHKERRQ(ierr);
457f32d5d43SBarry Smith   ierr = MatGetArray(C,&c);CHKERRQ(ierr);
458f32d5d43SBarry Smith   b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
459f32d5d43SBarry Smith   for (col=0; col<cn-4; col += 4){  /* over columns of C */
460f32d5d43SBarry Smith     colam = col*am;
461f32d5d43SBarry Smith     for (i=0; i<am; i++) {        /* over rows of C in those columns */
462f32d5d43SBarry Smith       r1 = r2 = r3 = r4 = 0.0;
463f32d5d43SBarry Smith       n   = a->i[i+1] - a->i[i];
464f32d5d43SBarry Smith       aj  = a->j + a->i[i];
465f32d5d43SBarry Smith       aa  = a->a + a->i[i];
466f32d5d43SBarry Smith       for (j=0; j<n; j++) {
467f32d5d43SBarry Smith         r1 += (*aa)*b1[*aj];
468f32d5d43SBarry Smith         r2 += (*aa)*b2[*aj];
469f32d5d43SBarry Smith         r3 += (*aa)*b3[*aj];
470f32d5d43SBarry Smith         r4 += (*aa++)*b4[*aj++];
4719123193aSHong Zhang       }
472f32d5d43SBarry Smith       c[colam + i]       = r1;
473f32d5d43SBarry Smith       c[colam + am + i]  = r2;
474f32d5d43SBarry Smith       c[colam + am2 + i] = r3;
475f32d5d43SBarry Smith       c[colam + am3 + i] = r4;
476f32d5d43SBarry Smith     }
477f32d5d43SBarry Smith     b1 += bm4;
478f32d5d43SBarry Smith     b2 += bm4;
479f32d5d43SBarry Smith     b3 += bm4;
480f32d5d43SBarry Smith     b4 += bm4;
481f32d5d43SBarry Smith   }
482f32d5d43SBarry Smith   for (;col<cn; col++){     /* over extra columns of C */
483f32d5d43SBarry Smith     for (i=0; i<am; i++) {  /* over rows of C in those columns */
484f32d5d43SBarry Smith       r1 = 0.0;
485f32d5d43SBarry Smith       n   = a->i[i+1] - a->i[i];
486f32d5d43SBarry Smith       aj  = a->j + a->i[i];
487f32d5d43SBarry Smith       aa  = a->a + a->i[i];
488f32d5d43SBarry Smith 
489f32d5d43SBarry Smith       for (j=0; j<n; j++) {
490f32d5d43SBarry Smith         r1 += (*aa++)*b1[*aj++];
491f32d5d43SBarry Smith       }
492f32d5d43SBarry Smith       c[col*am + i]     = r1;
493f32d5d43SBarry Smith     }
494f32d5d43SBarry Smith     b1 += bm;
495f32d5d43SBarry Smith   }
496dc0b31edSSatish Balay   ierr = PetscLogFlops(cn*(2.0*a->nz));CHKERRQ(ierr);
497f32d5d43SBarry Smith   ierr = MatRestoreArray(B,&b);CHKERRQ(ierr);
498f32d5d43SBarry Smith   ierr = MatRestoreArray(C,&c);CHKERRQ(ierr);
499f32d5d43SBarry Smith   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
500f32d5d43SBarry Smith   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
501f32d5d43SBarry Smith   PetscFunctionReturn(0);
502f32d5d43SBarry Smith }
503f32d5d43SBarry Smith 
504f32d5d43SBarry Smith /*
5054324174eSBarry Smith    Note very similar to MatMult_SeqAIJ(), should generate both codes from same base
506f32d5d43SBarry Smith */
507f32d5d43SBarry Smith #undef __FUNCT__
508f32d5d43SBarry Smith #define __FUNCT__ "MatMatMultNumericAdd_SeqAIJ_SeqDense"
509f32d5d43SBarry Smith PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
510f32d5d43SBarry Smith {
511f32d5d43SBarry Smith   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
512f32d5d43SBarry Smith   PetscErrorCode ierr;
513dd6ea824SBarry Smith   PetscScalar    *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4;
514dd6ea824SBarry Smith   MatScalar      *aa;
515d0f46423SBarry 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;
5164324174eSBarry Smith   PetscInt       am2 = 2*am, am3 = 3*am,  bm4 = 4*bm,colam,*ridx;
517f32d5d43SBarry Smith 
518f32d5d43SBarry Smith   PetscFunctionBegin;
519f32d5d43SBarry Smith   if (!cm || !cn) PetscFunctionReturn(0);
520f32d5d43SBarry Smith   ierr = MatGetArray(B,&b);CHKERRQ(ierr);
521f32d5d43SBarry Smith   ierr = MatGetArray(C,&c);CHKERRQ(ierr);
522f32d5d43SBarry Smith   b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
5234324174eSBarry Smith 
5244324174eSBarry Smith   if (a->compressedrow.use){ /* use compressed row format */
5254324174eSBarry Smith     for (col=0; col<cn-4; col += 4){  /* over columns of C */
5264324174eSBarry Smith       colam = col*am;
5274324174eSBarry Smith       arm   = a->compressedrow.nrows;
5284324174eSBarry Smith       ii    = a->compressedrow.i;
5294324174eSBarry Smith       ridx  = a->compressedrow.rindex;
5304324174eSBarry Smith       for (i=0; i<arm; i++) {        /* over rows of C in those columns */
5314324174eSBarry Smith 	r1 = r2 = r3 = r4 = 0.0;
5324324174eSBarry Smith 	n   = ii[i+1] - ii[i];
5334324174eSBarry Smith 	aj  = a->j + ii[i];
5344324174eSBarry Smith 	aa  = a->a + ii[i];
5354324174eSBarry Smith 	for (j=0; j<n; j++) {
5364324174eSBarry Smith 	  r1 += (*aa)*b1[*aj];
5374324174eSBarry Smith 	  r2 += (*aa)*b2[*aj];
5384324174eSBarry Smith 	  r3 += (*aa)*b3[*aj];
5394324174eSBarry Smith 	  r4 += (*aa++)*b4[*aj++];
5404324174eSBarry Smith 	}
5414324174eSBarry Smith 	c[colam       + ridx[i]] += r1;
5424324174eSBarry Smith 	c[colam + am  + ridx[i]] += r2;
5434324174eSBarry Smith 	c[colam + am2 + ridx[i]] += r3;
5444324174eSBarry Smith 	c[colam + am3 + ridx[i]] += r4;
5454324174eSBarry Smith       }
5464324174eSBarry Smith       b1 += bm4;
5474324174eSBarry Smith       b2 += bm4;
5484324174eSBarry Smith       b3 += bm4;
5494324174eSBarry Smith       b4 += bm4;
5504324174eSBarry Smith     }
5514324174eSBarry Smith     for (;col<cn; col++){     /* over extra columns of C */
5524324174eSBarry Smith       colam = col*am;
5534324174eSBarry Smith       arm   = a->compressedrow.nrows;
5544324174eSBarry Smith       ii    = a->compressedrow.i;
5554324174eSBarry Smith       ridx  = a->compressedrow.rindex;
5564324174eSBarry Smith       for (i=0; i<arm; i++) {  /* over rows of C in those columns */
5574324174eSBarry Smith 	r1 = 0.0;
5584324174eSBarry Smith 	n   = ii[i+1] - ii[i];
5594324174eSBarry Smith 	aj  = a->j + ii[i];
5604324174eSBarry Smith 	aa  = a->a + ii[i];
5614324174eSBarry Smith 
5624324174eSBarry Smith 	for (j=0; j<n; j++) {
5634324174eSBarry Smith 	  r1 += (*aa++)*b1[*aj++];
5644324174eSBarry Smith 	}
5654324174eSBarry Smith 	c[col*am + ridx[i]] += r1;
5664324174eSBarry Smith       }
5674324174eSBarry Smith       b1 += bm;
5684324174eSBarry Smith     }
5694324174eSBarry Smith   } else {
570f32d5d43SBarry Smith     for (col=0; col<cn-4; col += 4){  /* over columns of C */
571f32d5d43SBarry Smith       colam = col*am;
572f32d5d43SBarry Smith       for (i=0; i<am; i++) {        /* over rows of C in those columns */
573f32d5d43SBarry Smith 	r1 = r2 = r3 = r4 = 0.0;
574f32d5d43SBarry Smith 	n   = a->i[i+1] - a->i[i];
575f32d5d43SBarry Smith 	aj  = a->j + a->i[i];
576f32d5d43SBarry Smith 	aa  = a->a + a->i[i];
577f32d5d43SBarry Smith 	for (j=0; j<n; j++) {
578f32d5d43SBarry Smith 	  r1 += (*aa)*b1[*aj];
579f32d5d43SBarry Smith 	  r2 += (*aa)*b2[*aj];
580f32d5d43SBarry Smith 	  r3 += (*aa)*b3[*aj];
581f32d5d43SBarry Smith 	  r4 += (*aa++)*b4[*aj++];
582f32d5d43SBarry Smith 	}
583f32d5d43SBarry Smith 	c[colam + i]       += r1;
584f32d5d43SBarry Smith 	c[colam + am + i]  += r2;
585f32d5d43SBarry Smith 	c[colam + am2 + i] += r3;
586f32d5d43SBarry Smith 	c[colam + am3 + i] += r4;
587f32d5d43SBarry Smith       }
588f32d5d43SBarry Smith       b1 += bm4;
589f32d5d43SBarry Smith       b2 += bm4;
590f32d5d43SBarry Smith       b3 += bm4;
591f32d5d43SBarry Smith       b4 += bm4;
592f32d5d43SBarry Smith     }
593f32d5d43SBarry Smith     for (;col<cn; col++){     /* over extra columns of C */
594f32d5d43SBarry Smith       for (i=0; i<am; i++) {  /* over rows of C in those columns */
595f32d5d43SBarry Smith 	r1 = 0.0;
596f32d5d43SBarry Smith 	n   = a->i[i+1] - a->i[i];
597f32d5d43SBarry Smith 	aj  = a->j + a->i[i];
598f32d5d43SBarry Smith 	aa  = a->a + a->i[i];
599f32d5d43SBarry Smith 
600f32d5d43SBarry Smith 	for (j=0; j<n; j++) {
601f32d5d43SBarry Smith 	  r1 += (*aa++)*b1[*aj++];
602f32d5d43SBarry Smith 	}
603f32d5d43SBarry Smith 	c[col*am + i]     += r1;
604f32d5d43SBarry Smith       }
605f32d5d43SBarry Smith       b1 += bm;
606f32d5d43SBarry Smith     }
6074324174eSBarry Smith   }
608dc0b31edSSatish Balay   ierr = PetscLogFlops(cn*2.0*a->nz);CHKERRQ(ierr);
609f32d5d43SBarry Smith   ierr = MatRestoreArray(B,&b);CHKERRQ(ierr);
610f32d5d43SBarry Smith   ierr = MatRestoreArray(C,&c);CHKERRQ(ierr);
6119123193aSHong Zhang   PetscFunctionReturn(0);
6129123193aSHong Zhang }
613