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