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