xref: /petsc/src/mat/impls/aij/seq/matptap.c (revision 6849ba73f22fecb8f92ef896a42e4e8bd4cd6965)
1 /*
2   Defines projective product routines where A is a AIJ matrix
3           C = P^T * A * P
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 EXTERN PetscErrorCode RegisterMatMatMultRoutines_Private(Mat);
11 
12 static PetscEvent MAT_PtAPSymbolic = 0;
13 static PetscEvent MAT_PtAPNumeric  = 0;
14 
15 #undef __FUNCT__
16 #define __FUNCT__ "MatPtAP"
17 /*@
18    MatPtAP - Creates the matrix projection C = P^T * A * P
19 
20    Collective on Mat
21 
22    Input Parameters:
23 +  A - the matrix
24 .  P - the projection matrix
25 .  scall - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX
26 -  fill - expected fill as ratio of nnz(C)/(nnz(A) + nnz(P))
27 
28    Output Parameters:
29 .  C - the product matrix
30 
31    Notes:
32    C will be created and must be destroyed by the user with MatDestroy().
33 
34    This routine is currently only implemented for pairs of SeqAIJ matrices and classes
35    which inherit from SeqAIJ.  C will be of type MATSEQAIJ.
36 
37    Level: intermediate
38 
39 .seealso: MatPtAPSymbolic(),MatPtAPNumeric(),MatMatMult()
40 @*/
41 PetscErrorCode MatPtAP(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C) {
42   PetscErrorCode ierr;
43 
44   PetscFunctionBegin;
45   PetscValidHeaderSpecific(A,MAT_COOKIE,1);
46   PetscValidType(A,1);
47   MatPreallocated(A);
48   if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
49   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
50   PetscValidHeaderSpecific(P,MAT_COOKIE,2);
51   PetscValidType(P,2);
52   MatPreallocated(P);
53   if (!P->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
54   if (P->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
55   PetscValidPointer(C,3);
56   if (P->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->M,A->N);
57 
58   if (fill <=0.0) SETERRQ1(PETSC_ERR_ARG_SIZ,"fill=%g must be > 0.0",fill);
59 
60   ierr = PetscLogEventBegin(MAT_PtAP,A,P,0,0);CHKERRQ(ierr);
61   ierr = (*A->ops->matptap)(A,P,scall,fill,C);CHKERRQ(ierr);
62   ierr = PetscLogEventEnd(MAT_PtAP,A,P,0,0);CHKERRQ(ierr);
63   PetscFunctionReturn(0);
64 }
65 
66 #undef __FUNCT__
67 #define __FUNCT__ "MatPtAP_SeqAIJ_SeqAIJ"
68 PetscErrorCode MatPtAP_SeqAIJ_SeqAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
69 {
70   PetscErrorCode ierr;
71   PetscFunctionBegin;
72   if (scall == MAT_INITIAL_MATRIX){
73     ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ(A,P,fill,C);CHKERRQ(ierr);
74   }
75   ierr = MatPtAPNumeric_SeqAIJ_SeqAIJ(A,P,*C);CHKERRQ(ierr);
76   PetscFunctionReturn(0);
77 }
78 
79 #undef __FUNCT__
80 #define __FUNCT__ "MatPtAPSymbolic"
81 /*
82    MatPtAPSymbolic - Creates the (i,j) structure of the matrix projection C = P^T * A * P
83 
84    Collective on Mat
85 
86    Input Parameters:
87 +  A - the matrix
88 -  P - the projection matrix
89 
90    Output Parameters:
91 .  C - the (i,j) structure of the product matrix
92 
93    Notes:
94    C will be created and must be destroyed by the user with MatDestroy().
95 
96    This routine is currently only implemented for pairs of SeqAIJ matrices and classes
97    which inherit from SeqAIJ.  C will be of type MATSEQAIJ.  The product is computed using
98    this (i,j) structure by calling MatPtAPNumeric().
99 
100    Level: intermediate
101 
102 .seealso: MatPtAP(),MatPtAPNumeric(),MatMatMultSymbolic()
103 */
104 PetscErrorCode MatPtAPSymbolic(Mat A,Mat P,PetscReal fill,Mat *C) {
105   PetscErrorCode ierr;
106   char funct[80];
107   PetscErrorCode (*f)(Mat,Mat,Mat*);
108 
109   PetscFunctionBegin;
110 
111   PetscValidHeaderSpecific(A,MAT_COOKIE,1);
112   PetscValidType(A,1);
113   MatPreallocated(A);
114   if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
115   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
116 
117   PetscValidHeaderSpecific(P,MAT_COOKIE,2);
118   PetscValidType(P,2);
119   MatPreallocated(P);
120   if (!P->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
121   if (P->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
122 
123   PetscValidPointer(C,3);
124 
125   if (P->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->M,A->N);
126   if (A->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix 'A' must be square, %d != %d",A->M,A->N);
127 
128   /* Currently only _seqaij_seqaij is implemented, so just query for it. */
129   /* When other implementations exist, attack the multiple dispatch problem. */
130   ierr = PetscStrcpy(funct,"MatPtAPSymbolic_seqaij_seqaij");CHKERRQ(ierr);
131   ierr = PetscObjectQueryFunction((PetscObject)P,funct,(PetscVoidFunction)&f);CHKERRQ(ierr);
132   if (!f) SETERRQ1(PETSC_ERR_SUP,"MatPtAPSymbolic is not supported for P of type %s",P->type_name);
133   ierr = PetscObjectQueryFunction((PetscObject)A,funct,(PetscVoidFunction)&f);CHKERRQ(ierr);
134   if (!f) SETERRQ1(PETSC_ERR_SUP,"MatPtAPSymbolic is not supported for A of type %s",A->type_name);
135 
136   ierr = (*f)(A,P,C);CHKERRQ(ierr);
137 
138   PetscFunctionReturn(0);
139 }
140 
141 typedef struct {
142   Mat    symAP;
143 } Mat_PtAPstruct;
144 
145 #undef __FUNCT__
146 #define __FUNCT__ "MatDestroy_SeqAIJ_PtAP"
147 PetscErrorCode MatDestroy_SeqAIJ_PtAP(Mat A)
148 {
149   PetscErrorCode    ierr;
150   Mat_PtAPstruct    *ptap=(Mat_PtAPstruct*)A->spptr;
151 
152   PetscFunctionBegin;
153   ierr = MatDestroy(ptap->symAP);CHKERRQ(ierr);
154   ierr = PetscFree(ptap);CHKERRQ(ierr);
155 
156   ierr = MatDestroy(A);CHKERRQ(ierr);
157   PetscFunctionReturn(0);
158 }
159 
160 #undef __FUNCT__
161 #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqAIJ"
162 PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat P,PetscReal fill,Mat *C) {
163   PetscErrorCode ierr;
164   int            *pti,*ptj;
165   Mat            Pt,AP;
166   Mat_PtAPstruct *ptap;
167 
168   PetscFunctionBegin;
169   /* create symbolic Pt */
170   ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
171   ierr = MatCreateSeqAIJWithArrays(P->comm,P->N,P->M,pti,ptj,PETSC_NULL,&Pt);CHKERRQ(ierr);
172 
173   /* get symbolic AP=A*P and C=Pt*AP */
174   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,P,fill,&AP);CHKERRQ(ierr);
175   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(Pt,AP,fill,C);CHKERRQ(ierr);
176 
177   /* clean up */
178   ierr = MatRestoreSymbolicTranspose_SeqAIJ(Pt,&pti,&ptj);CHKERRQ(ierr);
179   ierr = MatDestroy(Pt);CHKERRQ(ierr);
180 
181   /* save symbolic AP - to be used by MatPtAPNumeric_SeqAIJ_SeqAIJ() */
182   ierr = PetscNew(Mat_PtAPstruct,&ptap);CHKERRQ(ierr);
183   ptap->symAP = AP;
184   (*C)->spptr = (void*)ptap;
185   (*C)->ops->destroy  = MatDestroy_SeqAIJ_PtAP;
186 
187   PetscFunctionReturn(0);
188 }
189 
190 #include "src/mat/impls/maij/maij.h"
191 EXTERN_C_BEGIN
192 #undef __FUNCT__
193 #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqMAIJ"
194 PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqMAIJ(Mat A,Mat PP,Mat *C) {
195   /* This routine requires testing -- I don't think it works. */
196   PetscErrorCode ierr;
197   FreeSpaceList  free_space=PETSC_NULL,current_space=PETSC_NULL;
198   Mat_SeqMAIJ    *pp=(Mat_SeqMAIJ*)PP->data;
199   Mat            P=pp->AIJ;
200   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*p=(Mat_SeqAIJ*)P->data,*c;
201   int            *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj;
202   int            *ci,*cj,*denserow,*sparserow,*ptadenserow,*ptasparserow,*ptaj;
203   int            an=A->N,am=A->M,pn=P->N,pm=P->M,ppdof=pp->dof;
204   int            i,j,k,dof,pdof,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi;
205   MatScalar      *ca;
206 
207   PetscFunctionBegin;
208   /* Start timer */
209   ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,PP,0,0);CHKERRQ(ierr);
210 
211   /* Get ij structure of P^T */
212   ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
213 
214   /* Allocate ci array, arrays for fill computation and */
215   /* free space for accumulating nonzero column info */
216   ierr = PetscMalloc((pn+1)*sizeof(int),&ci);CHKERRQ(ierr);
217   ci[0] = 0;
218 
219   ierr = PetscMalloc((2*pn+2*an+1)*sizeof(int),&ptadenserow);CHKERRQ(ierr);
220   ierr = PetscMemzero(ptadenserow,(2*pn+2*an+1)*sizeof(int));CHKERRQ(ierr);
221   ptasparserow = ptadenserow  + an;
222   denserow     = ptasparserow + an;
223   sparserow    = denserow     + pn;
224 
225   /* Set initial free space to be nnz(A) scaled by aspect ratio of P. */
226   /* This should be reasonable if sparsity of PtAP is similar to that of A. */
227   ierr          = GetMoreSpace((ai[am]/pm)*pn,&free_space);
228   current_space = free_space;
229 
230   /* Determine symbolic info for each row of C: */
231   for (i=0;i<pn/ppdof;i++) {
232     ptnzi  = pti[i+1] - pti[i];
233     ptanzi = 0;
234     ptJ    = ptj + pti[i];
235     for (dof=0;dof<ppdof;dof++) {
236     /* Determine symbolic row of PtA: */
237       for (j=0;j<ptnzi;j++) {
238         arow = ptJ[j] + dof;
239         anzj = ai[arow+1] - ai[arow];
240         ajj  = aj + ai[arow];
241         for (k=0;k<anzj;k++) {
242           if (!ptadenserow[ajj[k]]) {
243             ptadenserow[ajj[k]]    = -1;
244             ptasparserow[ptanzi++] = ajj[k];
245           }
246         }
247       }
248       /* Using symbolic info for row of PtA, determine symbolic info for row of C: */
249       ptaj = ptasparserow;
250       cnzi   = 0;
251       for (j=0;j<ptanzi;j++) {
252         pdof = *ptaj%dof;
253         prow = (*ptaj++)/dof;
254         pnzj = pi[prow+1] - pi[prow];
255         pjj  = pj + pi[prow];
256         for (k=0;k<pnzj;k++) {
257           if (!denserow[pjj[k]+pdof]) {
258             denserow[pjj[k]+pdof] = -1;
259             sparserow[cnzi++]     = pjj[k]+pdof;
260           }
261         }
262       }
263 
264       /* sort sparserow */
265       ierr = PetscSortInt(cnzi,sparserow);CHKERRQ(ierr);
266 
267       /* If free space is not available, make more free space */
268       /* Double the amount of total space in the list */
269       if (current_space->local_remaining<cnzi) {
270         ierr = GetMoreSpace(current_space->total_array_size,&current_space);CHKERRQ(ierr);
271       }
272 
273       /* Copy data into free space, and zero out denserows */
274       ierr = PetscMemcpy(current_space->array,sparserow,cnzi*sizeof(int));CHKERRQ(ierr);
275       current_space->array           += cnzi;
276       current_space->local_used      += cnzi;
277       current_space->local_remaining -= cnzi;
278 
279       for (j=0;j<ptanzi;j++) {
280         ptadenserow[ptasparserow[j]] = 0;
281       }
282       for (j=0;j<cnzi;j++) {
283         denserow[sparserow[j]] = 0;
284       }
285       /* Aside: Perhaps we should save the pta info for the numerical factorization. */
286       /*        For now, we will recompute what is needed. */
287       ci[i+1+dof] = ci[i+dof] + cnzi;
288     }
289   }
290   /* nnz is now stored in ci[ptm], column indices are in the list of free space */
291   /* Allocate space for cj, initialize cj, and */
292   /* destroy list of free space and other temporary array(s) */
293   ierr = PetscMalloc((ci[pn]+1)*sizeof(int),&cj);CHKERRQ(ierr);
294   ierr = MakeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
295   ierr = PetscFree(ptadenserow);CHKERRQ(ierr);
296 
297   /* Allocate space for ca */
298   ierr = PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
299   ierr = PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));CHKERRQ(ierr);
300 
301   /* put together the new matrix */
302   ierr = MatCreateSeqAIJWithArrays(A->comm,pn,pn,ci,cj,ca,C);CHKERRQ(ierr);
303 
304   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
305   /* Since these are PETSc arrays, change flags to free them as necessary. */
306   c = (Mat_SeqAIJ *)((*C)->data);
307   c->freedata = PETSC_TRUE;
308   c->nonew    = 0;
309 
310   /* Clean up. */
311   ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
312 
313   ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,PP,0,0);CHKERRQ(ierr);
314   PetscFunctionReturn(0);
315 }
316 EXTERN_C_END
317 
318 #undef __FUNCT__
319 #define __FUNCT__ "MatPtAPNumeric"
320 /*
321    MatPtAPNumeric - Computes the matrix projection C = P^T * A * P
322 
323    Collective on Mat
324 
325    Input Parameters:
326 +  A - the matrix
327 -  P - the projection matrix
328 
329    Output Parameters:
330 .  C - the product matrix
331 
332    Notes:
333    C must have been created by calling MatPtAPSymbolic and must be destroyed by
334    the user using MatDeatroy().
335 
336    This routine is currently only implemented for pairs of SeqAIJ matrices and classes
337    which inherit from SeqAIJ.  C will be of type MATSEQAIJ.
338 
339    Level: intermediate
340 
341 .seealso: MatPtAP(),MatPtAPSymbolic(),MatMatMultNumeric()
342 */
343 PetscErrorCode MatPtAPNumeric(Mat A,Mat P,Mat C) {
344   PetscErrorCode ierr;
345   char funct[80];
346   PetscErrorCode (*f)(Mat,Mat,Mat);
347 
348   PetscFunctionBegin;
349 
350   PetscValidHeaderSpecific(A,MAT_COOKIE,1);
351   PetscValidType(A,1);
352   MatPreallocated(A);
353   if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
354   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
355 
356   PetscValidHeaderSpecific(P,MAT_COOKIE,2);
357   PetscValidType(P,2);
358   MatPreallocated(P);
359   if (!P->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
360   if (P->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
361 
362   PetscValidHeaderSpecific(C,MAT_COOKIE,3);
363   PetscValidType(C,3);
364   MatPreallocated(C);
365   if (!C->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
366   if (C->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
367 
368   if (P->N!=C->M) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->N,C->M);
369   if (P->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->M,A->N);
370   if (A->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix 'A' must be square, %d != %d",A->M,A->N);
371   if (P->N!=C->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->N,C->N);
372 
373   /* Currently only _seqaij_seqaij is implemented, so just query for it. */
374   /* When other implementations exist, attack the multiple dispatch problem. */
375   ierr = PetscStrcpy(funct,"MatPtAPNumeric_seqaij_seqaij");CHKERRQ(ierr);
376   ierr = PetscObjectQueryFunction((PetscObject)P,funct,(PetscVoidFunction)&f);CHKERRQ(ierr);
377   if (!f) SETERRQ1(PETSC_ERR_SUP,"MatPtAPNumeric is not supported for P of type %s",P->type_name);
378   ierr = PetscObjectQueryFunction((PetscObject)A,funct,(PetscVoidFunction)&f);CHKERRQ(ierr);
379   if (!f) SETERRQ1(PETSC_ERR_SUP,"MatPtAPNumeric is not supported for A of type %s",A->type_name);
380 
381   ierr = (*f)(A,P,C);CHKERRQ(ierr);
382 
383   PetscFunctionReturn(0);
384 }
385 
386 #undef __FUNCT__
387 #define __FUNCT__ "MatPtAPNumeric_SeqAIJ_SeqAIJ"
388 PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat A,Mat P,Mat C)
389 {
390   PetscErrorCode ierr;
391   int        flops=0;
392   Mat_SeqAIJ *a  = (Mat_SeqAIJ *) A->data;
393   Mat_SeqAIJ *p  = (Mat_SeqAIJ *) P->data;
394   Mat_SeqAIJ *c  = (Mat_SeqAIJ *) C->data;
395   int        *ai=a->i,*aj=a->j,*pi=p->i,*pj=p->j,*pJ=p->j,*pjj;
396   int        *ci=c->i,*cj=c->j,*cjj;
397   int        am=A->M,cn=C->N,cm=C->M;
398   int        i,j,k,anzi,pnzi,apnzj,nextap,pnzj,prow,crow;
399   MatScalar  *aa=a->a,*apa,*pa=p->a,*pA=p->a,*paj,*ca=c->a,*caj;
400   Mat_PtAPstruct *ptap=(Mat_PtAPstruct*)C->spptr;
401   Mat_SeqAIJ     *ap = (Mat_SeqAIJ *)(ptap->symAP)->data;
402   int            *api=ap->i,*apj=ap->j,apj_nextap;
403 
404   PetscFunctionBegin;
405   /* Allocate temporary array for storage of one row of A*P */
406   ierr = PetscMalloc(cn*sizeof(MatScalar),&apa);CHKERRQ(ierr);
407   ierr = PetscMemzero(apa,cn*sizeof(MatScalar));CHKERRQ(ierr);
408 
409   /* Clear old values in C */
410   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
411 
412   for (i=0;i<am;i++) {
413     /* Get sparse values of A*P[i,:] */
414     anzi  = ai[i+1] - ai[i];
415     apnzj = 0;
416     for (j=0;j<anzi;j++) {
417       prow = *aj++;
418       pnzj = pi[prow+1] - pi[prow];
419       pjj  = pj + pi[prow];
420       paj  = pa + pi[prow];
421       for (k=0;k<pnzj;k++) {
422         apa[pjj[k]] += (*aa)*paj[k];
423       }
424       flops += 2*pnzj;
425       aa++;
426     }
427 
428     /* Compute P^T*A*P using outer product (P^T)[:,j]*(A*P)[j,:]. */
429     apj   = ap->j + api[i];
430     apnzj = api[i+1] - api[i];
431     pnzi  = pi[i+1] - pi[i];
432     for (j=0;j<pnzi;j++) {
433       nextap = 0;
434       crow   = *pJ++;
435       cjj    = cj + ci[crow];
436       caj    = ca + ci[crow];
437       /* Perform sparse axpy operation.  Note cjj includes apj. */
438       for (k=0; nextap<apnzj; k++) {
439         apj_nextap = *(apj+nextap);
440         if (cjj[k]==apj_nextap) {
441           caj[k] += (*pA)*apa[apj_nextap];
442           nextap++;
443         }
444       }
445       flops += 2*apnzj;
446       pA++;
447     }
448 
449     /* Zero the current row values for A*P */
450     for (j=0;j<apnzj;j++) apa[apj[j]] = 0.0;
451   }
452 
453   /* Assemble the final matrix and clean up */
454   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
455   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
456   ierr = PetscFree(apa);CHKERRQ(ierr);
457   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
458 
459   PetscFunctionReturn(0);
460 }
461 
462 #undef __FUNCT__
463 #define __FUNCT__ "RegisterPtAPRoutines_Private"
464 PetscErrorCode RegisterPtAPRoutines_Private(Mat A)
465 {
466   PetscErrorCode ierr;
467 
468   PetscFunctionBegin;
469   if (!MAT_PtAPSymbolic) {
470     ierr = PetscLogEventRegister(&MAT_PtAPSymbolic,"MatPtAPSymbolic",MAT_COOKIE);CHKERRQ(ierr);
471   }
472   if (!MAT_PtAPNumeric) {
473     ierr = PetscLogEventRegister(&MAT_PtAPNumeric,"MatPtAPNumeric",MAT_COOKIE);CHKERRQ(ierr);
474   }
475   ierr = RegisterMatMatMultRoutines_Private(A);CHKERRQ(ierr);
476 
477   PetscFunctionReturn(0);
478 }
479