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