xref: /petsc/src/mat/impls/aij/seq/matmatmult.c (revision 45a8bf62150c8074a4868c3b3823cd3ad0dac4d2)
1 /*
2   Defines matrix-matrix product routines for pairs of AIJ matrices
3           C = A * B
4 */
5 
6 #include "src/mat/impls/aij/seq/aij.h" /*I "petscmat.h" I*/
7 #include "src/mat/utils/freespace.h"
8 #include "src/mat/impls/aij/mpi/mpiaij.h"
9 
10 typedef struct { /* used by MatMatMult_MPIAIJ_MPIAIJ for reusing symbolic mat product */
11   IS     isrowa,isrowb,iscolb;
12   Mat    *aseq,*bseq,C_seq;
13 } Mat_MatMatMultMPI;
14 
15 static int logkey_matmatmult_symbolic = 0;
16 static int logkey_matmatmult_numeric  = 0;
17 
18 #undef __FUNCT__
19 #define __FUNCT__ "MatMatMult"
20 /*@
21    MatMatMult - Performs Matrix-Matrix Multiplication C=A*B.
22 
23    Collective on Mat
24 
25    Input Parameters:
26 +  A - the left matrix
27 .  B - the right matrix
28 .  scall - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX
29 -  fill - expected fill as ratio of nnz(C)/(nnz(A) + nnz(B))
30 
31    Output Parameters:
32 .  C - the product matrix
33 
34    Notes:
35    C will be created and must be destroyed by the user with MatDestroy().
36 
37    This routine is currently only implemented for pairs of SeqAIJ matrices and classes
38    which inherit from SeqAIJ.  C will be of type MATSEQAIJ.
39 
40    Level: intermediate
41 
42 .seealso: MatMatMultSymbolic(),MatMatMultNumeric()
43 @*/
44 int MatMatMult(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
45 {
46   int  ierr;
47 
48   PetscFunctionBegin;
49   PetscValidHeaderSpecific(A,MAT_COOKIE,1);
50   PetscValidType(A,1);
51   MatPreallocated(A);
52   if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
53   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
54   PetscValidHeaderSpecific(B,MAT_COOKIE,2);
55   PetscValidType(B,2);
56   MatPreallocated(B);
57   if (!B->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
58   if (B->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
59   PetscValidPointer(C,3);
60   if (B->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",B->M,A->N);
61 
62   if (fill <=0.0) SETERRQ1(PETSC_ERR_ARG_SIZ,"fill=%g must be > 0.0",fill);
63 
64   ierr = PetscLogEventBegin(MAT_MatMult,A,B,0,0);CHKERRQ(ierr);
65   ierr = (*A->ops->matmult)(A,B,scall,fill,C);CHKERRQ(ierr);
66   ierr = PetscLogEventEnd(MAT_MatMult,A,B,0,0);CHKERRQ(ierr);
67 
68   PetscFunctionReturn(0);
69 }
70 
71 #undef __FUNCT__
72 #define __FUNCT__ "MatMatMult_MPIAIJ_MPIAIJ"
73 int MatMatMult_MPIAIJ_MPIAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill, Mat *C)
74 {
75   int ierr;
76 
77   PetscFunctionBegin;
78   if (scall == MAT_INITIAL_MATRIX){
79     ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ(A,B,fill,C);CHKERRQ(ierr);/* numeric product is computed as well */
80   } else if (scall == MAT_REUSE_MATRIX){
81     ierr = MatMatMultNumeric_MPIAIJ_MPIAIJ(A,B,*C);CHKERRQ(ierr);
82   } else {
83     SETERRQ1(PETSC_ERR_ARG_WRONG,"Invalid MatReuse %d",scall);
84   }
85   PetscFunctionReturn(0);
86 }
87 
88 #undef __FUNCT__
89 #define __FUNCT__ "MatMatMult_SeqAIJ_SeqAIJ"
90 int MatMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) {
91   int ierr;
92 
93   PetscFunctionBegin;
94   if (scall == MAT_INITIAL_MATRIX){
95     ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
96   }
97   ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr);
98   PetscFunctionReturn(0);
99 }
100 
101 #undef __FUNCT__
102 #define __FUNCT__ "MatMatMultSymbolic"
103 /*@
104    MatMatMultSymbolic - Performs construction, preallocation, and computes the ij structure
105    of the matrix-matrix product C=A*B.  Call this routine before calling MatMatMultNumeric().
106 
107    Collective on Mat
108 
109    Input Parameters:
110 +  A - the left matrix
111 .  B - the right matrix
112 -  fill - expected fill as ratio of nnz(C)/(nnz(A) + nnz(B))
113 
114    Output Parameters:
115 .  C - the matrix containing the ij structure of product matrix
116 
117    Notes:
118    C will be created as a MATSEQAIJ matrix and must be destroyed by the user with MatDestroy().
119 
120    This routine is currently only implemented for SeqAIJ matrices and classes which inherit from SeqAIJ.
121 
122    Level: intermediate
123 
124 .seealso: MatMatMult(),MatMatMultNumeric()
125 @*/
126 int MatMatMultSymbolic(Mat A,Mat B,PetscReal fill,Mat *C) {
127   /* Perhaps this "interface" routine should be moved into the interface directory.*/
128   /* To facilitate implementations with varying types, QueryFunction is used.*/
129   /* It is assumed that implementations will be composed as "MatMatMultSymbolic_<type of A><type of B>". */
130   int  ierr;
131   char symfunct[80];
132   int  (*symbolic)(Mat,Mat,PetscReal,Mat *);
133 
134   PetscFunctionBegin;
135   PetscValidHeaderSpecific(A,MAT_COOKIE,1);
136   PetscValidType(A,1);
137   MatPreallocated(A);
138   if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
139   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
140 
141   PetscValidHeaderSpecific(B,MAT_COOKIE,2);
142   PetscValidType(B,2);
143   MatPreallocated(B);
144   if (!B->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
145   if (B->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
146   PetscValidPointer(C,3);
147 
148   if (B->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",B->M,A->N);
149   if (fill <=0.0) SETERRQ1(PETSC_ERR_ARG_SIZ,"fill=%g must be > 0.0",fill);
150 
151   /* Currently only _seqaijseqaij is implemented, so just query for it in A and B. */
152   /* When other implementations exist, attack the multiple dispatch problem. */
153   ierr = PetscStrcpy(symfunct,"MatMatMultSymbolic_seqaijseqaij");CHKERRQ(ierr);
154   ierr = PetscObjectQueryFunction((PetscObject)B,symfunct,(PetscVoidFunction)&symbolic);CHKERRQ(ierr);
155   if (!symbolic) SETERRQ1(PETSC_ERR_SUP,"C=A*B not implemented for B of type %s",B->type_name);
156   ierr = PetscObjectQueryFunction((PetscObject)A,symfunct,(PetscVoidFunction)&symbolic);CHKERRQ(ierr);
157   if (!symbolic) SETERRQ1(PETSC_ERR_SUP,"C=A*B not implemented for A of type %s",A->type_name);
158   ierr = (*symbolic)(A,B,fill,C);CHKERRQ(ierr);
159 
160   PetscFunctionReturn(0);
161 }
162 
163 EXTERN int MatDestroy_MPIAIJ(Mat);
164 #undef __FUNCT__
165 #define __FUNCT__ "MatDestroy_MPIAIJ_MatMatMult"
166 int MatDestroy_MPIAIJ_MatMatMult(Mat A)
167 {
168   int               ierr;
169   Mat_MatMatMultMPI *mult=(Mat_MatMatMultMPI*)A->spptr;
170 
171   PetscFunctionBegin;
172   ierr = ISDestroy(mult->isrowb);CHKERRQ(ierr);
173   ierr = ISDestroy(mult->iscolb);CHKERRQ(ierr);
174   ierr = ISDestroy(mult->isrowa);CHKERRQ(ierr);
175   ierr = MatDestroyMatrices(1,&mult->aseq);CHKERRQ(ierr);
176   ierr = MatDestroyMatrices(1,&mult->bseq);CHKERRQ(ierr);
177   ierr = MatDestroy(mult->C_seq);CHKERRQ(ierr);
178   ierr = PetscFree(mult);CHKERRQ(ierr);
179 
180   ierr = MatDestroy_MPIAIJ(A);CHKERRQ(ierr);
181 
182   PetscFunctionReturn(0);
183 }
184 
185 #undef __FUNCT__
186 #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ"
187 int MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
188 {
189   Mat_MPIAIJ        *a = (Mat_MPIAIJ*)A->data;
190   int               ierr,*idx,i,start,end,ncols,imark,nzA,nzB,*cmap;
191   Mat_MatMatMultMPI *mult;
192 
193   PetscFunctionBegin;
194   ierr = PetscNew(Mat_MatMatMultMPI,&mult);CHKERRQ(ierr);
195 
196   /* create a seq matrix B_seq = submatrix of B by taking rows of B that equal to nonzero col of A */
197   start = a->cstart;
198   cmap  = a->garray;
199   nzA   = a->A->n;
200   nzB   = a->B->n;
201   ierr  = PetscMalloc((nzA+nzB)*sizeof(int), &idx);CHKERRQ(ierr);
202   ncols = 0;
203   for (i=0; i<nzB; i++) {
204     if (cmap[i] < start) idx[ncols++] = cmap[i];
205     else break;
206   }
207   imark = i;
208   for (i=0; i<nzA; i++) idx[ncols++] = start + i;
209   for (i=imark; i<nzB; i++) idx[ncols++] = cmap[i];
210   ierr = ISCreateGeneral(PETSC_COMM_SELF,ncols,idx,&mult->isrowb);CHKERRQ(ierr);
211   ierr = PetscFree(idx);CHKERRQ(ierr);
212   ierr = ISCreateStride(PETSC_COMM_SELF,B->N,0,1,&mult->iscolb);CHKERRQ(ierr);
213   ierr = MatGetSubMatrices(B,1,&mult->isrowb,&mult->iscolb,MAT_INITIAL_MATRIX,&mult->bseq);CHKERRQ(ierr)
214 
215   /*  create a seq matrix A_seq = submatrix of A by taking all local rows of A */
216   start = a->rstart; end = a->rend;
217   ierr = ISCreateStride(PETSC_COMM_SELF,end-start,start,1,&mult->isrowa);CHKERRQ(ierr);
218   ierr = MatGetSubMatrices(A,1,&mult->isrowa,&mult->isrowb,MAT_INITIAL_MATRIX,&mult->aseq);CHKERRQ(ierr);
219 
220   /* compute C_seq = A_seq * B_seq */
221   ierr = MatMatMult_SeqAIJ_SeqAIJ(mult->aseq[0],mult->bseq[0],MAT_INITIAL_MATRIX,fill,&mult->C_seq);CHKERRQ(ierr);
222 
223   /* create mpi matrix C by concatinating C_seq */
224   ierr = PetscObjectReference((PetscObject)mult->C_seq);CHKERRQ(ierr); /* prevent C_seq being destroyed by MatMerge() */
225   ierr = MatMerge(A->comm,mult->C_seq,MAT_INITIAL_MATRIX,C);CHKERRQ(ierr);
226 
227   /* attach the supporting struct to C for reuse of symbolic C */
228   (*C)->spptr         = (void*)mult;
229   (*C)->ops->destroy  = MatDestroy_MPIAIJ_MatMatMult;
230 
231   PetscFunctionReturn(0);
232 }
233 
234 #undef __FUNCT__
235 #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ"
236 int MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
237 {
238   int            ierr;
239   FreeSpaceList  free_space=PETSC_NULL,current_space=PETSC_NULL;
240   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
241   int            *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj;
242   int            *ci,*cj,*lnk;
243   int            am=A->M,bn=B->N,bm=B->M;
244   int            i,j,anzi,brow,bnzj,cnzi,nlnk,lnk_init=-1,nspacedouble=0;
245   MatScalar      *ca;
246 
247   PetscFunctionBegin;
248   /* Start timers */
249   ierr = PetscLogEventBegin(logkey_matmatmult_symbolic,A,B,0,0);CHKERRQ(ierr);
250   /* Set up */
251   /* Allocate ci array, arrays for fill computation and */
252   /* free space for accumulating nonzero column info */
253   ierr = PetscMalloc(((am+1)+1)*sizeof(int),&ci);CHKERRQ(ierr);
254   ci[0] = 0;
255 
256   ierr = PetscMalloc((bn+1)*sizeof(int),&lnk);CHKERRQ(ierr);
257   ierr = PetscLLInitialize(lnk_init,bn,lnk);CHKERRQ(ierr);
258 
259   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
260   ierr = GetMoreSpace((int)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr);
261   current_space = free_space;
262 
263   /* Determine symbolic info for each row of the product: */
264   for (i=0;i<am;i++) {
265     anzi = ai[i+1] - ai[i];
266     cnzi = 0;
267     lnk[bn] = bn;
268     j       = anzi;
269     aj      = a->j + ai[i];
270     while (j){/* assume cols are almost in increasing order, starting from its end saves computation */
271       j--;
272       brow = *(aj + j);
273       bnzj = bi[brow+1] - bi[brow];
274       bjj  = bj + bi[brow];
275       /* add non-zero cols of B into the sorted linked list lnk */
276       ierr = PetscLLAdd(bnzj,bjj,bn,lnk_init,nlnk,lnk);CHKERRQ(ierr);
277       cnzi += nlnk;
278     }
279 
280     /* If free space is not available, make more free space */
281     /* Double the amount of total space in the list */
282     if (current_space->local_remaining<cnzi) {
283       /* printf("MatMatMultSymbolic_SeqAIJ_SeqAIJ()...%d -th row, double space ...\n",i);*/
284       ierr = GetMoreSpace(current_space->total_array_size,&current_space);CHKERRQ(ierr);
285       nspacedouble++;
286     }
287 
288     /* Copy data into free space, then initialize lnk */
289     ierr = PetscLLClear(bn,lnk_init,cnzi,lnk,current_space->array);CHKERRQ(ierr);
290     current_space->array += cnzi;
291 
292     current_space->local_used      += cnzi;
293     current_space->local_remaining -= cnzi;
294 
295     ci[i+1] = ci[i] + cnzi;
296   }
297 
298   /* Column indices are in the list of free space */
299   /* Allocate space for cj, initialize cj, and */
300   /* destroy list of free space and other temporary array(s) */
301   ierr = PetscMalloc((ci[am]+1)*sizeof(int),&cj);CHKERRQ(ierr);
302   ierr = MakeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
303   ierr = PetscFree(lnk);CHKERRQ(ierr);
304 
305   /* Allocate space for ca */
306   ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
307   ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr);
308 
309   /* put together the new symbolic matrix */
310   ierr = MatCreateSeqAIJWithArrays(A->comm,am,bn,ci,cj,ca,C);CHKERRQ(ierr);
311 
312   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
313   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
314   c = (Mat_SeqAIJ *)((*C)->data);
315   c->freedata = PETSC_TRUE;
316   c->nonew    = 0;
317 
318   ierr = PetscLogEventEnd(logkey_matmatmult_symbolic,A,B,0,0);CHKERRQ(ierr);
319   PetscLogInfo((PetscObject)(*C),"Number of calls to GetMoreSpace(): %d\n",nspacedouble);
320   PetscFunctionReturn(0);
321 }
322 
323 #undef __FUNCT__
324 #define __FUNCT__ "MatMatMultNumeric"
325 /*@
326    MatMatMultNumeric - Performs the numeric matrix-matrix product.
327    Call this routine after first calling MatMatMultSymbolic().
328 
329    Collective on Mat
330 
331    Input Parameters:
332 +  A - the left matrix
333 -  B - the right matrix
334 
335    Output Parameters:
336 .  C - the product matrix, whose ij structure was defined from MatMatMultSymbolic().
337 
338    Notes:
339    C must have been created with MatMatMultSymbolic.
340 
341    This routine is currently only implemented for SeqAIJ type matrices.
342 
343    Level: intermediate
344 
345 .seealso: MatMatMult(),MatMatMultSymbolic()
346 @*/
347 int MatMatMultNumeric(Mat A,Mat B,Mat C){
348   /* Perhaps this "interface" routine should be moved into the interface directory.*/
349   /* To facilitate implementations with varying types, QueryFunction is used.*/
350   /* It is assumed that implementations will be composed as "MatMatMultNumeric_<type of A><type of B>". */
351   int ierr;
352   char numfunct[80];
353   int (*numeric)(Mat,Mat,Mat);
354 
355   PetscFunctionBegin;
356 
357   PetscValidHeaderSpecific(A,MAT_COOKIE,1);
358   PetscValidType(A,1);
359   MatPreallocated(A);
360   if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
361   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
362 
363   PetscValidHeaderSpecific(B,MAT_COOKIE,2);
364   PetscValidType(B,2);
365   MatPreallocated(B);
366   if (!B->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
367   if (B->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
368 
369   PetscValidHeaderSpecific(C,MAT_COOKIE,3);
370   PetscValidType(C,3);
371   MatPreallocated(C);
372   if (!C->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
373   if (C->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
374 
375   if (B->N!=C->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",B->N,C->N);
376   if (B->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",B->M,A->N);
377   if (A->M!=C->M) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",A->M,C->M);
378 
379   /* Currently only _seqaijseqaij is implemented, so just query for it in A and B. */
380   /* When other implementations exist, attack the multiple dispatch problem. */
381   ierr = PetscStrcpy(numfunct,"MatMatMultNumeric_seqaijseqaij");CHKERRQ(ierr);
382   ierr = PetscObjectQueryFunction((PetscObject)A,numfunct,(PetscVoidFunction)&numeric);CHKERRQ(ierr);
383   if (!numeric) SETERRQ1(PETSC_ERR_SUP,"C=A*B not implemented for A of type %s",A->type_name);
384   ierr = PetscObjectQueryFunction((PetscObject)B,numfunct,(PetscVoidFunction)&numeric);CHKERRQ(ierr);
385   if (!numeric) SETERRQ1(PETSC_ERR_SUP,"C=A*B not implemented for B of type %s",B->type_name);
386 
387   ierr = (*numeric)(A,B,C);CHKERRQ(ierr);
388 
389   PetscFunctionReturn(0);
390 }
391 
392 /* This routine is called ONLY in the case of reusing previously computed symbolic C */
393 #undef __FUNCT__
394 #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ"
395 int MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat A,Mat B,Mat C)
396 {
397   int               ierr;
398   Mat_MatMatMultMPI *mult=(Mat_MatMatMultMPI*)C->spptr;
399 
400   PetscFunctionBegin;
401   ierr = MatGetSubMatrices(B,1,&mult->isrowb,&mult->iscolb,MAT_REUSE_MATRIX,&mult->bseq);CHKERRQ(ierr)
402   ierr = MatGetSubMatrices(A,1,&mult->isrowa,&mult->isrowb,MAT_REUSE_MATRIX,&mult->aseq);CHKERRQ(ierr);
403   ierr = MatMatMult_SeqAIJ_SeqAIJ(mult->aseq[0],mult->bseq[0],MAT_REUSE_MATRIX,0.0,&mult->C_seq);CHKERRQ(ierr);
404 
405   ierr = PetscObjectReference((PetscObject)mult->C_seq);CHKERRQ(ierr);
406   ierr = MatMerge(A->comm,mult->C_seq,MAT_REUSE_MATRIX,&C);CHKERRQ(ierr);
407 
408   PetscFunctionReturn(0);
409 }
410 
411 #undef __FUNCT__
412 #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqAIJ"
413 int MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
414 {
415   int        ierr,flops=0;
416   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;
417   Mat_SeqAIJ *b = (Mat_SeqAIJ *)B->data;
418   Mat_SeqAIJ *c = (Mat_SeqAIJ *)C->data;
419   int        *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j;
420   int        am=A->M,cn=C->N;
421   int        i,j,k,anzi,bnzi,cnzi,brow;
422   MatScalar  *aa=a->a,*ba=b->a,*baj,*ca=c->a,*temp;
423 
424   PetscFunctionBegin;
425 
426   /* Start timers */
427   ierr = PetscLogEventBegin(logkey_matmatmult_numeric,A,B,C,0);CHKERRQ(ierr);
428 
429   /* Allocate temp accumulation space to avoid searching for nonzero columns in C */
430   ierr = PetscMalloc((cn+1)*sizeof(MatScalar),&temp);CHKERRQ(ierr);
431   ierr = PetscMemzero(temp,cn*sizeof(MatScalar));CHKERRQ(ierr);
432   /* Traverse A row-wise. */
433   /* Build the ith row in C by summing over nonzero columns in A, */
434   /* the rows of B corresponding to nonzeros of A. */
435   for (i=0;i<am;i++) {
436     anzi = ai[i+1] - ai[i];
437     for (j=0;j<anzi;j++) {
438       brow = *aj++;
439       bnzi = bi[brow+1] - bi[brow];
440       bjj  = bj + bi[brow];
441       baj  = ba + bi[brow];
442       for (k=0;k<bnzi;k++) {
443         temp[bjj[k]] += (*aa)*baj[k];
444       }
445       flops += 2*bnzi;
446       aa++;
447     }
448     /* Store row back into C, and re-zero temp */
449     cnzi = ci[i+1] - ci[i];
450     for (j=0;j<cnzi;j++) {
451       ca[j] = temp[cj[j]];
452       temp[cj[j]] = 0.0;
453     }
454     ca += cnzi;
455     cj += cnzi;
456   }
457   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
458   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
459 
460   /* Free temp */
461   ierr = PetscFree(temp);CHKERRQ(ierr);
462   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
463   ierr = PetscLogEventEnd(logkey_matmatmult_numeric,A,B,C,0);CHKERRQ(ierr);
464   PetscFunctionReturn(0);
465 }
466 
467 #undef __FUNCT__
468 #define __FUNCT__ "RegisterMatMatMultRoutines_Private"
469 int RegisterMatMatMultRoutines_Private(Mat A) {
470   int ierr;
471 
472   PetscFunctionBegin;
473   if (!logkey_matmatmult_symbolic) {
474     ierr = PetscLogEventRegister(&logkey_matmatmult_symbolic,"MatMatMultSymbolic",MAT_COOKIE);CHKERRQ(ierr);
475   }
476   if (!logkey_matmatmult_numeric) {
477     ierr = PetscLogEventRegister(&logkey_matmatmult_numeric,"MatMatMultNumeric",MAT_COOKIE);CHKERRQ(ierr);
478   }
479 
480   PetscFunctionReturn(0);
481 }
482