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