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