xref: /petsc/src/mat/impls/aij/seq/matptap.c (revision 4c627768fec2a8f2c2ebd97828224a20479036a2)
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 EXTERN PetscErrorCode MatPtAP_SeqAIJ_SeqAIJ_ReducedPt(Mat,Mat,MatReuse,PetscReal,int,int,Mat*);
73 EXTERN PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_ReducedPt(Mat,Mat,PetscReal,int,int,Mat*);
74 EXTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_ReducedPt(Mat,Mat,int,int,Mat);
75 
76 #undef __FUNCT__
77 #define __FUNCT__ "MatPtAP_MPIAIJ_MPIAIJ"
78 PetscErrorCode MatPtAP_MPIAIJ_MPIAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
79 {
80   PetscErrorCode    ierr;
81 
82   PetscFunctionBegin;
83   if (scall == MAT_INITIAL_MATRIX){
84     ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
85     ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ(A,P,fill,C);CHKERRQ(ierr);/* numeric product is computed as well */
86     ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
87   } else if (scall == MAT_REUSE_MATRIX){
88     ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
89     ierr = MatPtAPNumeric_MPIAIJ_MPIAIJ(A,P,*C);CHKERRQ(ierr);
90     ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
91   } else {
92     SETERRQ1(PETSC_ERR_ARG_WRONG,"Invalid MatReuse %d",scall);
93   }
94   PetscFunctionReturn(0);
95 }
96 
97 #undef __FUNCT__
98 #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ"
99 PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
100 {
101   PetscErrorCode    ierr;
102   Mat               P_seq,A_loc,C_seq;
103   int               prstart,prend,m=P->m;
104   IS                isrowp,iscolp;
105 
106   PetscFunctionBegin;
107   /* get P_seq = submatrix of P by taking rows of P that equal to nonzero col of A */
108   ierr = MatGetBrowsOfAcols(A,P,MAT_INITIAL_MATRIX,&isrowp,&iscolp,&prstart,&P_seq);CHKERRQ(ierr);
109   ierr = ISDestroy(iscolp);CHKERRQ(ierr);
110 
111   /* get A_loc = submatrix of A by taking all local rows of A */
112   ierr = MatGetLocalMat(A,MAT_INITIAL_MATRIX,PETSC_NULL,&isrowp,&A_loc);CHKERRQ(ierr);
113   ierr = ISDestroy(isrowp);CHKERRQ(ierr);
114 
115   /* compute C_seq = P_loc^T * A_loc * P_seq */
116   prend   = prstart + m;
117   ierr = MatPtAP_SeqAIJ_SeqAIJ_ReducedPt(A_loc,P_seq,MAT_INITIAL_MATRIX,fill,prstart,prend,&C_seq);CHKERRQ(ierr);
118   ierr = MatDestroy(P_seq);CHKERRQ(ierr);
119   ierr = MatDestroy(A_loc);CHKERRQ(ierr);
120 
121   /* add C_seq into mpi C */
122   ierr = MatMerge_SeqsToMPI(A->comm,C_seq,P->n,P->n,MAT_INITIAL_MATRIX,C);CHKERRQ(ierr);
123 
124   PetscFunctionReturn(0);
125 }
126 
127 #undef __FUNCT__
128 #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ"
129 PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
130 {
131   PetscErrorCode       ierr;
132   Mat                  P_seq,A_loc,C_seq;
133   int                  prstart,prend,m=P->m;
134   IS                   isrowp,iscolp;
135   Mat_Merge_SeqsToMPI  *merge;
136   PetscObjectContainer container;
137 
138   PetscFunctionBegin;
139   ierr = PetscObjectQuery((PetscObject)C,"MatMergeSeqsToMPI",(PetscObject *)&container);CHKERRQ(ierr);
140   if (container) {
141     ierr  = PetscObjectContainerGetPointer(container,(void **)&merge);CHKERRQ(ierr);
142   } else {
143     SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Matrix C does not posses and object container");
144   }
145 
146   /* get P_seq = submatrix of P by taking rows of P that equal to nonzero col of A */
147   ierr = MatGetBrowsOfAcols(A,P,MAT_INITIAL_MATRIX,&isrowp,&iscolp,&prstart,&P_seq);CHKERRQ(ierr);
148   ierr = ISDestroy(iscolp);CHKERRQ(ierr);
149 
150   /* get A_loc = submatrix of A by taking all local rows of A */
151   ierr = MatGetLocalMat(A,MAT_INITIAL_MATRIX,PETSC_NULL,&isrowp,&A_loc);CHKERRQ(ierr);
152   ierr = ISDestroy(isrowp);CHKERRQ(ierr);
153 
154   /* compute C_seq = P_loc^T * A_loc * P_seq */
155   prend = prstart + m;
156   C_seq = merge->C_seq;
157   ierr = MatPtAP_SeqAIJ_SeqAIJ_ReducedPt(A_loc,P_seq,MAT_REUSE_MATRIX,1.0,prstart,prend,&C_seq);CHKERRQ(ierr);
158   ierr = MatDestroy(P_seq);CHKERRQ(ierr);
159   ierr = MatDestroy(A_loc);CHKERRQ(ierr);
160 
161   /* add C_seq into mpi C */
162   ierr = MatMerge_SeqsToMPI(A->comm,C_seq,P->n,P->n,MAT_REUSE_MATRIX,&C);CHKERRQ(ierr);
163 
164   PetscFunctionReturn(0);
165 }
166 
167 #undef __FUNCT__
168 #define __FUNCT__ "MatPtAP_SeqAIJ_SeqAIJ"
169 PetscErrorCode MatPtAP_SeqAIJ_SeqAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
170 {
171   PetscErrorCode ierr;
172   PetscFunctionBegin;
173   if (scall == MAT_INITIAL_MATRIX){
174     ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
175     ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ(A,P,fill,C);CHKERRQ(ierr);
176     ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
177   }
178   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
179   ierr = MatPtAPNumeric_SeqAIJ_SeqAIJ(A,P,*C);CHKERRQ(ierr);
180   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
181   PetscFunctionReturn(0);
182 }
183 
184 #undef __FUNCT__
185 #define __FUNCT__ "MatPtAPSymbolic"
186 /*
187    MatPtAPSymbolic - Creates the (i,j) structure of the matrix projection C = P^T * A * P
188 
189    Collective on Mat
190 
191    Input Parameters:
192 +  A - the matrix
193 -  P - the projection matrix
194 
195    Output Parameters:
196 .  C - the (i,j) structure of the product matrix
197 
198    Notes:
199    C will be created and must be destroyed by the user with MatDestroy().
200 
201    This routine is currently only implemented for pairs of SeqAIJ matrices and classes
202    which inherit from SeqAIJ.  C will be of type MATSEQAIJ.  The product is computed using
203    this (i,j) structure by calling MatPtAPNumeric().
204 
205    Level: intermediate
206 
207 .seealso: MatPtAP(),MatPtAPNumeric(),MatMatMultSymbolic()
208 */
209 PetscErrorCode MatPtAPSymbolic(Mat A,Mat P,PetscReal fill,Mat *C) {
210   PetscErrorCode ierr;
211   PetscErrorCode (*fA)(Mat,Mat,PetscReal,Mat*);
212   PetscErrorCode (*fP)(Mat,Mat,PetscReal,Mat*);
213 
214   PetscFunctionBegin;
215 
216   PetscValidHeaderSpecific(A,MAT_COOKIE,1);
217   PetscValidType(A,1);
218   MatPreallocated(A);
219   if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
220   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
221 
222   PetscValidHeaderSpecific(P,MAT_COOKIE,2);
223   PetscValidType(P,2);
224   MatPreallocated(P);
225   if (!P->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
226   if (P->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
227 
228   PetscValidPointer(C,3);
229 
230   if (P->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->M,A->N);
231   if (A->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix 'A' must be square, %d != %d",A->M,A->N);
232 
233   /* For now, we do not dispatch based on the type of A and P */
234   /* When implementations like _SeqAIJ_MAIJ exist, attack the multiple dispatch problem. */
235   fA = A->ops->ptapsymbolic;
236   if (!fA) SETERRQ1(PETSC_ERR_SUP,"MatPtAPSymbolic not supported for A of type %s",A->type_name);
237   fP = P->ops->ptapsymbolic;
238   if (!fP) SETERRQ1(PETSC_ERR_SUP,"MatPtAPSymbolic not supported for P of type %s",P->type_name);
239   if (fP!=fA) SETERRQ2(PETSC_ERR_ARG_INCOMP,"MatPtAPSymbolic requires A, %s, to be compatible with P, %s",A->type_name,P->type_name);
240 
241   ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
242   ierr = (*fA)(A,P,fill,C);CHKERRQ(ierr);
243   ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
244 
245   PetscFunctionReturn(0);
246 }
247 
248 typedef struct {
249   Mat    symAP;
250 } Mat_PtAPstruct;
251 
252 EXTERN PetscErrorCode MatDestroy_SeqAIJ(Mat);
253 
254 #undef __FUNCT__
255 #define __FUNCT__ "MatDestroy_SeqAIJ_PtAP"
256 PetscErrorCode MatDestroy_SeqAIJ_PtAP(Mat A)
257 {
258   PetscErrorCode    ierr;
259   Mat_PtAPstruct    *ptap=(Mat_PtAPstruct*)A->spptr;
260 
261   PetscFunctionBegin;
262   ierr = MatDestroy(ptap->symAP);CHKERRQ(ierr);
263   ierr = PetscFree(ptap);CHKERRQ(ierr);
264   ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr);
265   PetscFunctionReturn(0);
266 }
267 
268 #undef __FUNCT__
269 #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqAIJ"
270 PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat P,PetscReal fill,Mat *C) {
271   PetscErrorCode ierr;
272   FreeSpaceList  free_space=PETSC_NULL,current_space=PETSC_NULL;
273   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*p=(Mat_SeqAIJ*)P->data,*c;
274   int            *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj;
275   int            *ci,*cj,*denserow,*sparserow,*ptadenserow,*ptasparserow,*ptaj;
276   int            an=A->N,am=A->M,pn=P->N,pm=P->M;
277   int            i,j,k,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi;
278   MatScalar      *ca;
279 
280   PetscFunctionBegin;
281   /* Get ij structure of P^T */
282   ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
283   ptJ=ptj;
284 
285   /* Allocate ci array, arrays for fill computation and */
286   /* free space for accumulating nonzero column info */
287   ierr = PetscMalloc((pn+1)*sizeof(int),&ci);CHKERRQ(ierr);
288   ci[0] = 0;
289 
290   ierr = PetscMalloc((2*pn+2*an+1)*sizeof(int),&ptadenserow);CHKERRQ(ierr);
291   ierr = PetscMemzero(ptadenserow,(2*pn+2*an+1)*sizeof(int));CHKERRQ(ierr);
292   ptasparserow = ptadenserow  + an;
293   denserow     = ptasparserow + an;
294   sparserow    = denserow     + pn;
295 
296   /* Set initial free space to be nnz(A) scaled by aspect ratio of P. */
297   /* This should be reasonable if sparsity of PtAP is similar to that of A. */
298   ierr          = GetMoreSpace((ai[am]/pm)*pn,&free_space);
299   current_space = free_space;
300 
301   /* Determine symbolic info for each row of C: */
302   for (i=0;i<pn;i++) {
303     ptnzi  = pti[i+1] - pti[i];
304     ptanzi = 0;
305     /* Determine symbolic row of PtA: */
306     for (j=0;j<ptnzi;j++) {
307       arow = *ptJ++;
308       anzj = ai[arow+1] - ai[arow];
309       ajj  = aj + ai[arow];
310       for (k=0;k<anzj;k++) {
311         if (!ptadenserow[ajj[k]]) {
312           ptadenserow[ajj[k]]    = -1;
313           ptasparserow[ptanzi++] = ajj[k];
314         }
315       }
316     }
317       /* Using symbolic info for row of PtA, determine symbolic info for row of C: */
318     ptaj = ptasparserow;
319     cnzi   = 0;
320     for (j=0;j<ptanzi;j++) {
321       prow = *ptaj++;
322       pnzj = pi[prow+1] - pi[prow];
323       pjj  = pj + pi[prow];
324       for (k=0;k<pnzj;k++) {
325         if (!denserow[pjj[k]]) {
326             denserow[pjj[k]]  = -1;
327             sparserow[cnzi++] = pjj[k];
328         }
329       }
330     }
331 
332     /* sort sparserow */
333     ierr = PetscSortInt(cnzi,sparserow);CHKERRQ(ierr);
334 
335     /* If free space is not available, make more free space */
336     /* Double the amount of total space in the list */
337     if (current_space->local_remaining<cnzi) {
338       ierr = GetMoreSpace(current_space->total_array_size,&current_space);CHKERRQ(ierr);
339     }
340 
341     /* Copy data into free space, and zero out denserows */
342     ierr = PetscMemcpy(current_space->array,sparserow,cnzi*sizeof(int));CHKERRQ(ierr);
343     current_space->array           += cnzi;
344     current_space->local_used      += cnzi;
345     current_space->local_remaining -= cnzi;
346 
347     for (j=0;j<ptanzi;j++) {
348       ptadenserow[ptasparserow[j]] = 0;
349     }
350     for (j=0;j<cnzi;j++) {
351       denserow[sparserow[j]] = 0;
352     }
353       /* Aside: Perhaps we should save the pta info for the numerical factorization. */
354       /*        For now, we will recompute what is needed. */
355     ci[i+1] = ci[i] + cnzi;
356   }
357   /* nnz is now stored in ci[ptm], column indices are in the list of free space */
358   /* Allocate space for cj, initialize cj, and */
359   /* destroy list of free space and other temporary array(s) */
360   ierr = PetscMalloc((ci[pn]+1)*sizeof(int),&cj);CHKERRQ(ierr);
361   ierr = MakeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
362   ierr = PetscFree(ptadenserow);CHKERRQ(ierr);
363 
364   /* Allocate space for ca */
365   ierr = PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
366   ierr = PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));CHKERRQ(ierr);
367 
368   /* put together the new matrix */
369   ierr = MatCreateSeqAIJWithArrays(A->comm,pn,pn,ci,cj,ca,C);CHKERRQ(ierr);
370 
371   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
372   /* Since these are PETSc arrays, change flags to free them as necessary. */
373   c = (Mat_SeqAIJ *)((*C)->data);
374   c->freedata = PETSC_TRUE;
375   c->nonew    = 0;
376 
377   /* Clean up. */
378   ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
379 
380   PetscFunctionReturn(0);
381 }
382 
383 #include "src/mat/impls/maij/maij.h"
384 EXTERN_C_BEGIN
385 #undef __FUNCT__
386 #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqMAIJ"
387 PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqMAIJ(Mat A,Mat PP,Mat *C) {
388   /* This routine requires testing -- I don't think it works. */
389   PetscErrorCode ierr;
390   FreeSpaceList  free_space=PETSC_NULL,current_space=PETSC_NULL;
391   Mat_SeqMAIJ    *pp=(Mat_SeqMAIJ*)PP->data;
392   Mat            P=pp->AIJ;
393   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*p=(Mat_SeqAIJ*)P->data,*c;
394   int            *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj;
395   int            *ci,*cj,*denserow,*sparserow,*ptadenserow,*ptasparserow,*ptaj;
396   int            an=A->N,am=A->M,pn=P->N,pm=P->M,ppdof=pp->dof;
397   int            i,j,k,dof,pdof,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi;
398   MatScalar      *ca;
399 
400   PetscFunctionBegin;
401   /* Start timer */
402   ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,PP,0,0);CHKERRQ(ierr);
403 
404   /* Get ij structure of P^T */
405   ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
406 
407   /* Allocate ci array, arrays for fill computation and */
408   /* free space for accumulating nonzero column info */
409   ierr = PetscMalloc((pn+1)*sizeof(int),&ci);CHKERRQ(ierr);
410   ci[0] = 0;
411 
412   ierr = PetscMalloc((2*pn+2*an+1)*sizeof(int),&ptadenserow);CHKERRQ(ierr);
413   ierr = PetscMemzero(ptadenserow,(2*pn+2*an+1)*sizeof(int));CHKERRQ(ierr);
414   ptasparserow = ptadenserow  + an;
415   denserow     = ptasparserow + an;
416   sparserow    = denserow     + pn;
417 
418   /* Set initial free space to be nnz(A) scaled by aspect ratio of P. */
419   /* This should be reasonable if sparsity of PtAP is similar to that of A. */
420   ierr          = GetMoreSpace((ai[am]/pm)*pn,&free_space);
421   current_space = free_space;
422 
423   /* Determine symbolic info for each row of C: */
424   for (i=0;i<pn/ppdof;i++) {
425     ptnzi  = pti[i+1] - pti[i];
426     ptanzi = 0;
427     ptJ    = ptj + pti[i];
428     for (dof=0;dof<ppdof;dof++) {
429     /* Determine symbolic row of PtA: */
430       for (j=0;j<ptnzi;j++) {
431         arow = ptJ[j] + dof;
432         anzj = ai[arow+1] - ai[arow];
433         ajj  = aj + ai[arow];
434         for (k=0;k<anzj;k++) {
435           if (!ptadenserow[ajj[k]]) {
436             ptadenserow[ajj[k]]    = -1;
437             ptasparserow[ptanzi++] = ajj[k];
438           }
439         }
440       }
441       /* Using symbolic info for row of PtA, determine symbolic info for row of C: */
442       ptaj = ptasparserow;
443       cnzi   = 0;
444       for (j=0;j<ptanzi;j++) {
445         pdof = *ptaj%dof;
446         prow = (*ptaj++)/dof;
447         pnzj = pi[prow+1] - pi[prow];
448         pjj  = pj + pi[prow];
449         for (k=0;k<pnzj;k++) {
450           if (!denserow[pjj[k]+pdof]) {
451             denserow[pjj[k]+pdof] = -1;
452             sparserow[cnzi++]     = pjj[k]+pdof;
453           }
454         }
455       }
456 
457       /* sort sparserow */
458       ierr = PetscSortInt(cnzi,sparserow);CHKERRQ(ierr);
459 
460       /* If free space is not available, make more free space */
461       /* Double the amount of total space in the list */
462       if (current_space->local_remaining<cnzi) {
463         ierr = GetMoreSpace(current_space->total_array_size,&current_space);CHKERRQ(ierr);
464       }
465 
466       /* Copy data into free space, and zero out denserows */
467       ierr = PetscMemcpy(current_space->array,sparserow,cnzi*sizeof(int));CHKERRQ(ierr);
468       current_space->array           += cnzi;
469       current_space->local_used      += cnzi;
470       current_space->local_remaining -= cnzi;
471 
472       for (j=0;j<ptanzi;j++) {
473         ptadenserow[ptasparserow[j]] = 0;
474       }
475       for (j=0;j<cnzi;j++) {
476         denserow[sparserow[j]] = 0;
477       }
478       /* Aside: Perhaps we should save the pta info for the numerical factorization. */
479       /*        For now, we will recompute what is needed. */
480       ci[i+1+dof] = ci[i+dof] + cnzi;
481     }
482   }
483   /* nnz is now stored in ci[ptm], column indices are in the list of free space */
484   /* Allocate space for cj, initialize cj, and */
485   /* destroy list of free space and other temporary array(s) */
486   ierr = PetscMalloc((ci[pn]+1)*sizeof(int),&cj);CHKERRQ(ierr);
487   ierr = MakeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
488   ierr = PetscFree(ptadenserow);CHKERRQ(ierr);
489 
490   /* Allocate space for ca */
491   ierr = PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
492   ierr = PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));CHKERRQ(ierr);
493 
494   /* put together the new matrix */
495   ierr = MatCreateSeqAIJWithArrays(A->comm,pn,pn,ci,cj,ca,C);CHKERRQ(ierr);
496 
497   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
498   /* Since these are PETSc arrays, change flags to free them as necessary. */
499   c = (Mat_SeqAIJ *)((*C)->data);
500   c->freedata = PETSC_TRUE;
501   c->nonew    = 0;
502 
503   /* Clean up. */
504   ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
505 
506   ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,PP,0,0);CHKERRQ(ierr);
507   PetscFunctionReturn(0);
508 }
509 EXTERN_C_END
510 
511 #undef __FUNCT__
512 #define __FUNCT__ "MatPtAPNumeric"
513 /*
514    MatPtAPNumeric - Computes the matrix projection C = P^T * A * P
515 
516    Collective on Mat
517 
518    Input Parameters:
519 +  A - the matrix
520 -  P - the projection matrix
521 
522    Output Parameters:
523 .  C - the product matrix
524 
525    Notes:
526    C must have been created by calling MatPtAPSymbolic and must be destroyed by
527    the user using MatDeatroy().
528 
529    This routine is currently only implemented for pairs of AIJ matrices and classes
530    which inherit from AIJ.  C will be of type MATAIJ.
531 
532    Level: intermediate
533 
534 .seealso: MatPtAP(),MatPtAPSymbolic(),MatMatMultNumeric()
535 */
536 PetscErrorCode MatPtAPNumeric(Mat A,Mat P,Mat C) {
537   PetscErrorCode ierr;
538   PetscErrorCode (*fA)(Mat,Mat,Mat);
539   PetscErrorCode (*fP)(Mat,Mat,Mat);
540 
541   PetscFunctionBegin;
542 
543   PetscValidHeaderSpecific(A,MAT_COOKIE,1);
544   PetscValidType(A,1);
545   MatPreallocated(A);
546   if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
547   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
548 
549   PetscValidHeaderSpecific(P,MAT_COOKIE,2);
550   PetscValidType(P,2);
551   MatPreallocated(P);
552   if (!P->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
553   if (P->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
554 
555   PetscValidHeaderSpecific(C,MAT_COOKIE,3);
556   PetscValidType(C,3);
557   MatPreallocated(C);
558   if (!C->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
559   if (C->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
560 
561   if (P->N!=C->M) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->N,C->M);
562   if (P->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->M,A->N);
563   if (A->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix 'A' must be square, %d != %d",A->M,A->N);
564   if (P->N!=C->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->N,C->N);
565 
566   /* For now, we do not dispatch based on the type of A and P */
567   /* When implementations like _SeqAIJ_MAIJ exist, attack the multiple dispatch problem. */
568   fA = A->ops->ptapnumeric;
569   if (!fA) SETERRQ1(PETSC_ERR_SUP,"MatPtAPNumeric not supported for A of type %s",A->type_name);
570   fP = P->ops->ptapnumeric;
571   if (!fP) SETERRQ1(PETSC_ERR_SUP,"MatPtAPNumeric not supported for P of type %s",P->type_name);
572   if (fP!=fA) SETERRQ2(PETSC_ERR_ARG_INCOMP,"MatPtAPNumeric requires A, %s, to be compatible with P, %s",A->type_name,P->type_name);
573 
574   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
575   ierr = (*fA)(A,P,C);CHKERRQ(ierr);
576   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
577 
578   PetscFunctionReturn(0);
579 }
580 
581 #undef __FUNCT__
582 #define __FUNCT__ "MatPtAPNumeric_SeqAIJ_SeqAIJ"
583 PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat A,Mat P,Mat C)
584 {
585   PetscErrorCode ierr;
586   int            flops=0;
587   Mat_SeqAIJ     *a  = (Mat_SeqAIJ *) A->data;
588   Mat_SeqAIJ     *p  = (Mat_SeqAIJ *) P->data;
589   Mat_SeqAIJ     *c  = (Mat_SeqAIJ *) C->data;
590   int            *ai=a->i,*aj=a->j,*apj,*apjdense,*pi=p->i,*pj=p->j,*pJ=p->j,*pjj;
591   int            *ci=c->i,*cj=c->j,*cjj;
592   int            am=A->M,cn=C->N,cm=C->M;
593   int            i,j,k,anzi,pnzi,apnzj,nextap,pnzj,prow,crow;
594   MatScalar      *aa=a->a,*apa,*pa=p->a,*pA=p->a,*paj,*ca=c->a,*caj;
595 
596   PetscFunctionBegin;
597   /* Allocate temporary array for storage of one row of A*P */
598   ierr = PetscMalloc(cn*(sizeof(MatScalar)+2*sizeof(int)),&apa);CHKERRQ(ierr);
599   ierr = PetscMemzero(apa,cn*(sizeof(MatScalar)+2*sizeof(int)));CHKERRQ(ierr);
600 
601   apj      = (int *)(apa + cn);
602   apjdense = apj + cn;
603 
604   /* Clear old values in C */
605   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
606 
607   for (i=0;i<am;i++) {
608     /* Form sparse row of A*P */
609     anzi  = ai[i+1] - ai[i];
610     apnzj = 0;
611     for (j=0;j<anzi;j++) {
612       prow = *aj++;
613       pnzj = pi[prow+1] - pi[prow];
614       pjj  = pj + pi[prow];
615       paj  = pa + pi[prow];
616       for (k=0;k<pnzj;k++) {
617         if (!apjdense[pjj[k]]) {
618           apjdense[pjj[k]] = -1;
619           apj[apnzj++]     = pjj[k];
620         }
621         apa[pjj[k]] += (*aa)*paj[k];
622       }
623       flops += 2*pnzj;
624       aa++;
625     }
626 
627     /* Sort the j index array for quick sparse axpy. */
628     ierr = PetscSortInt(apnzj,apj);CHKERRQ(ierr);
629 
630     /* Compute P^T*A*P using outer product (P^T)[:,j]*(A*P)[j,:]. */
631     pnzi = pi[i+1] - pi[i];
632     for (j=0;j<pnzi;j++) {
633       nextap = 0;
634       crow   = *pJ++;
635       cjj    = cj + ci[crow];
636       caj    = ca + ci[crow];
637       /* Perform sparse axpy operation.  Note cjj includes apj. */
638       for (k=0;nextap<apnzj;k++) {
639         if (cjj[k]==apj[nextap]) {
640           caj[k] += (*pA)*apa[apj[nextap++]];
641         }
642       }
643       flops += 2*apnzj;
644       pA++;
645     }
646 
647     /* Zero the current row info for A*P */
648     for (j=0;j<apnzj;j++) {
649       apa[apj[j]]      = 0.;
650       apjdense[apj[j]] = 0;
651     }
652   }
653 
654   /* Assemble the final matrix and clean up */
655   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
656   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
657   ierr = PetscFree(apa);CHKERRQ(ierr);
658   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
659 
660   PetscFunctionReturn(0);
661 }
662 
663 /* Compute C = P[rstart:rend,:]^T * A * P of seqaij matrices - used by MatPtAP_MPIAIJ_MPIAIJ() */
664 
665 #undef __FUNCT__
666 #define __FUNCT__ "MatPtAPNumeric_SeqAIJ_SeqAIJ_ReducedPt"
667 /*@C
668    MatPtAPNumeric_SeqAIJ_SeqAIJ_ReducedPt - Creates C = P[rstart:rend,:]^T * A * P of seqaij matrices,
669                 used by MatPtAP_MPIAIJ_MPIAIJ()
670 
671    Collective on Mat
672 
673    Input Parameters:
674 +  A - the matrix in seqaij format
675 .  P - the projection matrix in seqaij format
676 .  scall - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX
677 .  fill - expected fill (not being used when scall==MAT_REUSE_MATRIX)
678 .  prstart, prend - the starting and ending-1-th row of the matrix P to be used for transpose
679 
680    Output Parameters:
681 .  C - the product matrix in seqaij format
682 
683    Level: developer
684 
685 .seealso: MatPtAPSymbolic(),MatPtAPNumeric(),MatMatMult()
686 @*/
687 PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_ReducedPt(Mat A,Mat P,int prstart,int prend,Mat C)
688 {
689   PetscErrorCode ierr;
690   int            flops=0;
691   Mat_SeqAIJ     *a  = (Mat_SeqAIJ *) A->data;
692   Mat_SeqAIJ     *p  = (Mat_SeqAIJ *) P->data;
693   Mat_SeqAIJ     *c  = (Mat_SeqAIJ *) C->data;
694   int            *ai=a->i,*aj=a->j,*apj,*apjdense;
695   int            *pi=p->i,*pj=p->j,*pJ=p->j+pi[prstart],*pjj;
696   int            *ci=c->i,*cj=c->j,*cjj;
697   int            am=A->M,cn=C->N,cm=C->M;
698   int            i,j,k,anzi,pnzi,apnzj,nextap,pnzj,prow,crow;
699   MatScalar      *aa=a->a,*apa,*pa=p->a,*pA=p->a,*paj,*ca=c->a,*caj;
700 
701   PetscFunctionBegin;
702   pA=p->a+pi[prstart];
703   /* Allocate temporary array for storage of one row of A*P */
704   ierr = PetscMalloc(cn*(sizeof(MatScalar)+2*sizeof(int)),&apa);CHKERRQ(ierr);
705   ierr = PetscMemzero(apa,cn*(sizeof(MatScalar)+2*sizeof(int)));CHKERRQ(ierr);
706 
707   apj      = (int *)(apa + cn);
708   apjdense = apj + cn;
709 
710   /* Clear old values in C */
711   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
712 
713   for (i=0;i<am;i++) {
714     /* Form sparse row of A*P */
715     anzi  = ai[i+1] - ai[i];
716     apnzj = 0;
717     for (j=0;j<anzi;j++) {
718       prow = *aj++;
719       pnzj = pi[prow+1] - pi[prow];
720       pjj  = pj + pi[prow];
721       paj  = pa + pi[prow];
722       for (k=0;k<pnzj;k++) {
723         if (!apjdense[pjj[k]]) {
724           apjdense[pjj[k]] = -1;
725           apj[apnzj++]     = pjj[k];
726         }
727         apa[pjj[k]] += (*aa)*paj[k];
728       }
729       flops += 2*pnzj;
730       aa++;
731     }
732 
733     /* Sort the j index array for quick sparse axpy. */
734     ierr = PetscSortInt(apnzj,apj);CHKERRQ(ierr);
735 
736     /* Compute P[[prstart:prend,:]^T*A*P using outer product (P^T)[:,j+prstart]*(A*P)[j,:]. */
737     pnzi = pi[i+1+prstart] - pi[i+prstart];
738     for (j=0;j<pnzi;j++) {
739       nextap = 0;
740       crow   = *pJ++;
741       cjj    = cj + ci[crow];
742       caj    = ca + ci[crow];
743       /* Perform sparse axpy operation.  Note cjj includes apj. */
744       for (k=0;nextap<apnzj;k++) {
745         if (cjj[k]==apj[nextap]) {
746           caj[k] += (*pA)*apa[apj[nextap++]];
747         }
748       }
749       flops += 2*apnzj;
750       pA++;
751     }
752 
753     /* Zero the current row info for A*P */
754     for (j=0;j<apnzj;j++) {
755       apa[apj[j]]      = 0.;
756       apjdense[apj[j]] = 0;
757     }
758   }
759 
760   /* Assemble the final matrix and clean up */
761   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
762   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
763   ierr = PetscFree(apa);CHKERRQ(ierr);
764   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
765 
766   PetscFunctionReturn(0);
767 }
768 
769 #undef __FUNCT__
770 #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqAIJ_ReducedPt"
771 PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_ReducedPt(Mat A,Mat P,PetscReal fill,int prstart,int prend,Mat *C) {
772   PetscErrorCode ierr;
773   FreeSpaceList  free_space=PETSC_NULL,current_space=PETSC_NULL;
774   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*p=(Mat_SeqAIJ*)P->data,*c;
775   int            *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj;
776   int            *ci,*cj,*denserow,*sparserow,*ptadenserow,*ptasparserow,*ptaj;
777   int            an=A->N,am=A->M,pn=P->N,pm=P->M;
778   int            i,j,k,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi;
779   MatScalar      *ca;
780   Mat            *psub,P_sub;
781   IS             isrow,iscol;
782   int            m = prend - prstart;
783 
784   PetscFunctionBegin;
785   /* Get ij structure of P[rstart:rend,:]^T */
786   ierr = ISCreateStride(PETSC_COMM_SELF,m,prstart,1,&isrow);CHKERRQ(ierr);
787   ierr = ISCreateStride(PETSC_COMM_SELF,P->n,0,1,&iscol);CHKERRQ(ierr);
788   ierr = MatGetSubMatrices(P,1,&isrow,&iscol,MAT_INITIAL_MATRIX,&psub);CHKERRQ(ierr);
789   ierr = ISDestroy(isrow);CHKERRQ(ierr);
790   ierr = ISDestroy(iscol);CHKERRQ(ierr);
791   P_sub = psub[0];
792   ierr = MatGetSymbolicTranspose_SeqAIJ(P_sub,&pti,&ptj);CHKERRQ(ierr);
793   ierr = MatDestroyMatrices(1,&psub);CHKERRQ(ierr);
794   ptJ=ptj;
795 
796   /* Allocate ci array, arrays for fill computation and */
797   /* free space for accumulating nonzero column info */
798   ierr = PetscMalloc((pn+1)*sizeof(int),&ci);CHKERRQ(ierr);
799   ci[0] = 0;
800 
801   ierr = PetscMalloc((2*pn+2*an+1)*sizeof(int),&ptadenserow);CHKERRQ(ierr);
802   ierr = PetscMemzero(ptadenserow,(2*pn+2*an+1)*sizeof(int));CHKERRQ(ierr);
803   ptasparserow = ptadenserow  + an;
804   denserow     = ptasparserow + an;
805   sparserow    = denserow     + pn;
806 
807   /* Set initial free space to be nnz(A) scaled by aspect ratio of P. */
808   /* This should be reasonable if sparsity of PtAP is similar to that of A. */
809   ierr          = GetMoreSpace((int)(fill*ai[am]/pm)*pn,&free_space);
810   current_space = free_space;
811 
812   /* Determine symbolic info for each row of C: */
813   for (i=0;i<pn;i++) {
814     ptnzi  = pti[i+1] - pti[i];
815     ptanzi = 0;
816     /* Determine symbolic row of PtA_reduced: */
817     for (j=0;j<ptnzi;j++) {
818       arow = *ptJ++;
819       anzj = ai[arow+1] - ai[arow];
820       ajj  = aj + ai[arow];
821       for (k=0;k<anzj;k++) {
822         if (!ptadenserow[ajj[k]]) {
823           ptadenserow[ajj[k]]    = -1;
824           ptasparserow[ptanzi++] = ajj[k];
825         }
826       }
827     }
828       /* Using symbolic info for row of PtA, determine symbolic info for row of C: */
829     ptaj = ptasparserow;
830     cnzi   = 0;
831     for (j=0;j<ptanzi;j++) {
832       prow = *ptaj++;
833       pnzj = pi[prow+1] - pi[prow];
834       pjj  = pj + pi[prow];
835       for (k=0;k<pnzj;k++) {
836         if (!denserow[pjj[k]]) {
837             denserow[pjj[k]]  = -1;
838             sparserow[cnzi++] = pjj[k];
839         }
840       }
841     }
842 
843     /* sort sparserow */
844     ierr = PetscSortInt(cnzi,sparserow);CHKERRQ(ierr);
845 
846     /* If free space is not available, make more free space */
847     /* Double the amount of total space in the list */
848     if (current_space->local_remaining<cnzi) {
849       ierr = GetMoreSpace(current_space->total_array_size,&current_space);CHKERRQ(ierr);
850     }
851 
852     /* Copy data into free space, and zero out denserows */
853     ierr = PetscMemcpy(current_space->array,sparserow,cnzi*sizeof(int));CHKERRQ(ierr);
854     current_space->array           += cnzi;
855     current_space->local_used      += cnzi;
856     current_space->local_remaining -= cnzi;
857 
858     for (j=0;j<ptanzi;j++) {
859       ptadenserow[ptasparserow[j]] = 0;
860     }
861     for (j=0;j<cnzi;j++) {
862       denserow[sparserow[j]] = 0;
863     }
864       /* Aside: Perhaps we should save the pta info for the numerical factorization. */
865       /*        For now, we will recompute what is needed. */
866     ci[i+1] = ci[i] + cnzi;
867   }
868   /* nnz is now stored in ci[ptm], column indices are in the list of free space */
869   /* Allocate space for cj, initialize cj, and */
870   /* destroy list of free space and other temporary array(s) */
871   ierr = PetscMalloc((ci[pn]+1)*sizeof(int),&cj);CHKERRQ(ierr);
872   ierr = MakeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
873   ierr = PetscFree(ptadenserow);CHKERRQ(ierr);
874 
875   /* Allocate space for ca */
876   ierr = PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
877   ierr = PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));CHKERRQ(ierr);
878 
879   /* put together the new matrix */
880   ierr = MatCreateSeqAIJWithArrays(A->comm,pn,pn,ci,cj,ca,C);CHKERRQ(ierr);
881 
882   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
883   /* Since these are PETSc arrays, change flags to free them as necessary. */
884   c = (Mat_SeqAIJ *)((*C)->data);
885   c->freedata = PETSC_TRUE;
886   c->nonew    = 0;
887 
888   /* Clean up. */
889   ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
890 
891   PetscFunctionReturn(0);
892 }
893 
894 #undef __FUNCT__
895 #define __FUNCT__ "MatPtAP_SeqAIJ_SeqAIJ_ReducedPt"
896 PetscErrorCode MatPtAP_SeqAIJ_SeqAIJ_ReducedPt(Mat A,Mat P,MatReuse scall,PetscReal fill,int prstart,int prend,Mat *C)
897 {
898   PetscErrorCode ierr;
899   PetscFunctionBegin;
900   if (A->m != prend-prstart) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",A->m,prend-prstart);
901   if (prend-prstart > P->m) SETERRQ2(PETSC_ERR_ARG_SIZ," prend-prstart %d cannot be larger than P->m %d",prend-prstart,P->m);
902   if (scall == MAT_INITIAL_MATRIX){
903     ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ_ReducedPt(A,P,fill,prstart,prend,C);CHKERRQ(ierr);
904   }
905 
906   ierr = MatPtAPNumeric_SeqAIJ_SeqAIJ_ReducedPt(A,P,prstart,prend,*C);CHKERRQ(ierr);
907 
908   PetscFunctionReturn(0);
909 }
910 
911