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