xref: /petsc/src/mat/impls/aij/seq/matptap.c (revision 9c44eddc3abd5ec499e19db53d7a9ea05dfa3e9e)
1 
2 /*
3   Defines projective product routines where A is a SeqAIJ matrix
4           C = P^T * A * P
5 */
6 
7 #include <../src/mat/impls/aij/seq/aij.h>   /*I "petscmat.h" I*/
8 #include <../src/mat/utils/freespace.h>
9 #include <petscbt.h>
10 
11 #undef __FUNCT__
12 #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ"
13 PetscErrorCode MatPtAPSymbolic_SeqAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
14 {
15   PetscErrorCode ierr;
16 
17   PetscFunctionBegin;
18   if (!P->ops->ptapsymbolic_seqaij) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_SUP,"Not implemented for A=%s and P=%s",((PetscObject)A)->type_name,((PetscObject)P)->type_name);
19   ierr = (*P->ops->ptapsymbolic_seqaij)(A,P,fill,C);CHKERRQ(ierr);
20   PetscFunctionReturn(0);
21 }
22 
23 #undef __FUNCT__
24 #define __FUNCT__ "MatPtAPNumeric_SeqAIJ"
25 PetscErrorCode MatPtAPNumeric_SeqAIJ(Mat A,Mat P,Mat C)
26 {
27   PetscErrorCode ierr;
28 
29   PetscFunctionBegin;
30   if (!P->ops->ptapnumeric_seqaij) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_SUP,"Not implemented for A=%s and P=%s",((PetscObject)A)->type_name,((PetscObject)P)->type_name);
31   ierr = (*P->ops->ptapnumeric_seqaij)(A,P,C);CHKERRQ(ierr);
32   PetscFunctionReturn(0);
33 }
34 
35 #undef __FUNCT__
36 #define __FUNCT__ "MatDestroy_SeqAIJ_PtAP"
37 PetscErrorCode MatDestroy_SeqAIJ_PtAP(Mat A)
38 {
39   PetscErrorCode ierr;
40   Mat_SeqAIJ     *a = (Mat_SeqAIJ *)A->data;
41   Mat_PtAP       *ptap = a->ptap;
42 
43   PetscFunctionBegin;
44   /* free ptap, then A */
45   ierr = PetscFree(ptap->apa);CHKERRQ(ierr);
46   ierr = PetscFree(ptap->api);CHKERRQ(ierr);
47   ierr = PetscFree(ptap->apj);CHKERRQ(ierr);
48   ierr = (ptap->destroy)(A);CHKERRQ(ierr);
49   ierr = PetscFree(ptap);CHKERRQ(ierr);
50   PetscFunctionReturn(0);
51 }
52 
53 #undef __FUNCT__
54 #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy2"
55 PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy2(Mat A,Mat P,PetscReal fill,Mat *C)
56 {
57   PetscErrorCode     ierr;
58   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
59   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data,*p = (Mat_SeqAIJ*)P->data,*c;
60   PetscInt           *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj;
61   PetscInt           *ci,*cj,*ptadenserow,*ptasparserow,*ptaj,nspacedouble=0;
62   PetscInt           an=A->cmap->N,am=A->rmap->N,pn=P->cmap->N;
63   PetscInt           i,j,k,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi,nlnk,*lnk;
64   MatScalar          *ca;
65   PetscBT            lnkbt;
66 
67   PetscFunctionBegin;
68   /* Get ij structure of P^T */
69   ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
70   ptJ=ptj;
71 
72   /* Allocate ci array, arrays for fill computation and */
73   /* free space for accumulating nonzero column info */
74   ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
75   ci[0] = 0;
76 
77   ierr = PetscMalloc((2*an+1)*sizeof(PetscInt),&ptadenserow);CHKERRQ(ierr);
78   ierr = PetscMemzero(ptadenserow,(2*an+1)*sizeof(PetscInt));CHKERRQ(ierr);
79   ptasparserow = ptadenserow  + an;
80 
81   /* create and initialize a linked list */
82   nlnk = pn+1;
83   ierr = PetscLLCreate(pn,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr);
84 
85   /* Set initial free space to be fill*nnz(A). */
86   /* This should be reasonable if sparsity of PtAP is similar to that of A. */
87   ierr          = PetscFreeSpaceGet((PetscInt)(fill*ai[am]),&free_space);
88   current_space = free_space;
89 
90   /* Determine symbolic info for each row of C: */
91   for (i=0;i<pn;i++) {
92     ptnzi  = pti[i+1] - pti[i];
93     ptanzi = 0;
94     /* Determine symbolic row of PtA: */
95     for (j=0;j<ptnzi;j++) {
96       arow = *ptJ++;
97       anzj = ai[arow+1] - ai[arow];
98       ajj  = aj + ai[arow];
99       for (k=0;k<anzj;k++) {
100         if (!ptadenserow[ajj[k]]) {
101           ptadenserow[ajj[k]]    = -1;
102           ptasparserow[ptanzi++] = ajj[k];
103         }
104       }
105     }
106     /* Using symbolic info for row of PtA, determine symbolic info for row of C: */
107     ptaj = ptasparserow;
108     cnzi   = 0;
109     for (j=0;j<ptanzi;j++) {
110       prow = *ptaj++;
111       pnzj = pi[prow+1] - pi[prow];
112       pjj  = pj + pi[prow];
113       /* add non-zero cols of P into the sorted linked list lnk */
114       ierr = PetscLLAddSorted(pnzj,pjj,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr);
115       cnzi += nlnk;
116     }
117 
118     /* If free space is not available, make more free space */
119     /* Double the amount of total space in the list */
120     if (current_space->local_remaining<cnzi) {
121       ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,&current_space);CHKERRQ(ierr);
122       nspacedouble++;
123     }
124 
125     /* Copy data into free space, and zero out denserows */
126     ierr = PetscLLClean(pn,pn,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
127     current_space->array           += cnzi;
128     current_space->local_used      += cnzi;
129     current_space->local_remaining -= cnzi;
130 
131     for (j=0;j<ptanzi;j++) {
132       ptadenserow[ptasparserow[j]] = 0;
133     }
134     /* Aside: Perhaps we should save the pta info for the numerical factorization. */
135     /*        For now, we will recompute what is needed. */
136     ci[i+1] = ci[i] + cnzi;
137   }
138   /* nnz is now stored in ci[ptm], column indices are in the list of free space */
139   /* Allocate space for cj, initialize cj, and */
140   /* destroy list of free space and other temporary array(s) */
141   ierr = PetscMalloc((ci[pn]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
142   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
143   ierr = PetscFree(ptadenserow);CHKERRQ(ierr);
144   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
145 
146   /* Allocate space for ca */
147   ierr = PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
148   ierr = PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));CHKERRQ(ierr);
149 
150   /* put together the new matrix */
151   ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,pn,pn,ci,cj,ca,C);CHKERRQ(ierr);
152 
153   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
154   /* Since these are PETSc arrays, change flags to free them as necessary. */
155   c = (Mat_SeqAIJ *)((*C)->data);
156   c->free_a  = PETSC_TRUE;
157   c->free_ij = PETSC_TRUE;
158   c->nonew   = 0;
159   A->ops->ptapnumeric = MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy2; /* should use *C->ops until PtAP insterface is updated to double dispatch as MatMatMult() */
160 
161   /* Clean up. */
162   ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
163 #if defined(PETSC_USE_INFO)
164   if (ci[pn] != 0) {
165     PetscReal afill = ((PetscReal)ci[pn])/ai[am];
166     if (afill < 1.0) afill = 1.0;
167     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
168     ierr = PetscInfo1((*C),"Use MatPtAP(A,P,MatReuse,%G,&C) for best performance.\n",afill);CHKERRQ(ierr);
169   } else {
170     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
171   }
172 #endif
173   PetscFunctionReturn(0);
174 }
175 
176 #undef __FUNCT__
177 #define __FUNCT__ "MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy2"
178 PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy2(Mat A,Mat P,Mat C)
179 {
180   PetscErrorCode ierr;
181   Mat_SeqAIJ     *a  = (Mat_SeqAIJ *) A->data;
182   Mat_SeqAIJ     *p  = (Mat_SeqAIJ *) P->data;
183   Mat_SeqAIJ     *c  = (Mat_SeqAIJ *) C->data;
184   PetscInt       *ai=a->i,*aj=a->j,*apj,*apjdense,*pi=p->i,*pj=p->j,*pJ=p->j,*pjj;
185   PetscInt       *ci=c->i,*cj=c->j,*cjj;
186   PetscInt       am=A->rmap->N,cn=C->cmap->N,cm=C->rmap->N;
187   PetscInt       i,j,k,anzi,pnzi,apnzj,nextap,pnzj,prow,crow;
188   MatScalar      *aa=a->a,*apa,*pa=p->a,*pA=p->a,*paj,*ca=c->a,*caj;
189 
190   PetscFunctionBegin;
191   /* Allocate temporary array for storage of one row of A*P */
192   ierr = PetscMalloc(cn*(sizeof(MatScalar)+sizeof(PetscInt))+c->rmax*sizeof(PetscInt),&apa);CHKERRQ(ierr);
193 
194   apjdense = (PetscInt *)(apa + cn);
195   apj      = apjdense + cn;
196   ierr = PetscMemzero(apa,cn*(sizeof(MatScalar)+sizeof(PetscInt)));CHKERRQ(ierr);
197 
198   /* Clear old values in C */
199   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
200 
201   for (i=0; i<am; i++) {
202     /* Form sparse row of A*P */
203     anzi  = ai[i+1] - ai[i];
204     apnzj = 0;
205     for (j=0; j<anzi; j++) {
206       prow = *aj++;
207       pnzj = pi[prow+1] - pi[prow];
208       pjj  = pj + pi[prow];
209       paj  = pa + pi[prow];
210       for (k=0;k<pnzj;k++) {
211         if (!apjdense[pjj[k]]) {
212           apjdense[pjj[k]] = -1;
213           apj[apnzj++]     = pjj[k];
214         }
215         apa[pjj[k]] += (*aa)*paj[k];
216       }
217       ierr = PetscLogFlops(2.0*pnzj);CHKERRQ(ierr);
218       aa++;
219     }
220 
221     /* Sort the j index array for quick sparse axpy. */
222     /* Note: a array does not need sorting as it is in dense storage locations. */
223     ierr = PetscSortInt(apnzj,apj);CHKERRQ(ierr);
224 
225     /* Compute P^T*A*P using outer product (P^T)[:,j]*(A*P)[j,:]. */
226     pnzi = pi[i+1] - pi[i];
227     for (j=0; j<pnzi; j++) {
228       nextap = 0;
229       crow   = *pJ++;
230       cjj    = cj + ci[crow];
231       caj    = ca + ci[crow];
232       /* Perform sparse axpy operation.  Note cjj includes apj. */
233       for (k=0;nextap<apnzj;k++) {
234 #if defined(PETSC_USE_DEBUG)
235         if (k >= ci[crow+1] - ci[crow]) {
236           SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"k too large k %d, crow %d",k,crow);
237         }
238 #endif
239         if (cjj[k]==apj[nextap]) {
240           caj[k] += (*pA)*apa[apj[nextap++]];
241         }
242       }
243       ierr = PetscLogFlops(2.0*apnzj);CHKERRQ(ierr);
244       pA++;
245     }
246 
247     /* Zero the current row info for A*P */
248     for (j=0; j<apnzj; j++) {
249       apa[apj[j]]      = 0.;
250       apjdense[apj[j]] = 0;
251     }
252   }
253 
254   /* Assemble the final matrix and clean up */
255   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
256   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
257   ierr = PetscFree(apa);CHKERRQ(ierr);
258   PetscFunctionReturn(0);
259 }
260 
261 #undef __FUNCT__
262 #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqAIJ"
263 PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
264 {
265   PetscErrorCode     ierr;
266   Mat_SeqAIJ         *ap,*c;
267   PetscInt           *api,*apj,*ci,*cj,pn=P->cmap->N,sparse_axpy=0;
268   MatScalar          *ca;
269   Mat_PtAP           *ptap;
270   Mat                Pt,AP;
271 
272   PetscFunctionBegin;
273   /* flag 'sparse_axpy' determines which implementations to be used:
274        0: do dense axpy in MatPtAPNumeric() - fastest, but requires storage of struct A*P; (default)
275        1: do one sparse axpy - uses same memory as sparse_axpy=0 and might execute less flops
276           (apnz vs. cnz in the outerproduct), slower than case '0' when cnz is not too large than apnz;
277        2: do two sparse axpy in MatPtAPNumeric() - slowest, does not store structure of A*P. */
278   ierr = PetscOptionsGetInt(PETSC_NULL,"-matptap_sparseaxpy",&sparse_axpy,PETSC_NULL);CHKERRQ(ierr);
279   if (sparse_axpy == 2){
280     ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy2(A,P,fill,C);CHKERRQ(ierr);
281     PetscFunctionReturn(0);
282   }
283 
284   /* Get symbolic Pt = P^T */
285   ierr = MatTransposeSymbolic_SeqAIJ(P,&Pt);CHKERRQ(ierr);
286 
287   /* Get symbolic AP = A*P */
288   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,P,fill,&AP);CHKERRQ(ierr);
289   ap          = (Mat_SeqAIJ*)AP->data;
290   api         = ap->i;
291   apj         = ap->j;
292   ap->free_ij = PETSC_FALSE; /* api and apj are kept in struct ptap, cannot be destroyed with AP */
293 
294   /* Get C = Pt*AP */
295   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(Pt,AP,fill,C);CHKERRQ(ierr);
296   c  = (Mat_SeqAIJ*)(*C)->data;
297   ci = c->i;
298   cj = c->j;
299   ierr = PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
300   ierr = PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));CHKERRQ(ierr);
301   c->a       = ca;
302   c->free_a  = PETSC_TRUE;
303 
304   /* Create a supporting struct for reuse by MatPtAPNumeric() */
305   ierr = PetscNew(Mat_PtAP,&ptap);CHKERRQ(ierr);
306   c->ptap            = ptap;
307   ptap->destroy      = (*C)->ops->destroy;
308   (*C)->ops->destroy = MatDestroy_SeqAIJ_PtAP;
309 
310   /* Allocate temporary array for storage of one row of A*P */
311   ierr = PetscMalloc((pn+1)*sizeof(PetscScalar),&ptap->apa);CHKERRQ(ierr);
312   ierr = PetscMemzero(ptap->apa,(pn+1)*sizeof(PetscScalar));CHKERRQ(ierr);
313   if (sparse_axpy == 1){
314     A->ops->ptapnumeric = MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy;
315   } else {
316     A->ops->ptapnumeric = MatPtAPNumeric_SeqAIJ_SeqAIJ;
317   }
318   ptap->api = api;
319   ptap->apj = apj;
320 
321   /* Clean up. */
322   ierr = MatDestroy(&Pt);CHKERRQ(ierr);
323   ierr = MatDestroy(&AP);CHKERRQ(ierr);
324   PetscFunctionReturn(0);
325 }
326 
327 /* #define PROFILE_MatPtAPNumeric */
328 #undef __FUNCT__
329 #define __FUNCT__ "MatPtAPNumeric_SeqAIJ_SeqAIJ"
330 PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat A,Mat P,Mat C)
331 {
332   PetscErrorCode    ierr;
333   Mat_SeqAIJ        *a  = (Mat_SeqAIJ *) A->data;
334   Mat_SeqAIJ        *p  = (Mat_SeqAIJ *) P->data;
335   Mat_SeqAIJ        *c  = (Mat_SeqAIJ *) C->data;
336   const PetscInt    *ai=a->i,*aj=a->j,*pi=p->i,*pj=p->j,*ci=c->i,*cj=c->j;
337   const PetscScalar *aa=a->a,*pa=p->a,*pval;
338   const PetscInt    *apj,*pcol,*cjj;
339   const PetscInt    am=A->rmap->N,cm=C->rmap->N;
340   PetscInt          i,j,k,anz,apnz,pnz,prow,crow,cnz;
341   PetscScalar       *apa,*ca=c->a,*caj,pvalj;
342   Mat_PtAP          *ptap = c->ptap;
343 #if defined(PROFILE_MatPtAPNumeric)
344   PetscLogDouble    t0,tf,time_Cseq0=0.0,time_Cseq1=0.0;
345   PetscInt          flops0=0,flops1=0;
346 #endif
347 
348   PetscFunctionBegin;
349   /* Get temporary array for storage of one row of A*P */
350   apa = ptap->apa;
351 
352   /* Clear old values in C */
353   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
354 
355   for (i=0;i<am;i++) {
356     /* Form sparse row of AP[i,:] = A[i,:]*P */
357 #if defined(PROFILE_MatPtAPNumeric)
358     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
359 #endif
360     anz  = ai[i+1] - ai[i];
361     apnz = 0;
362     for (j=0; j<anz; j++) {
363       prow = aj[j];
364       pnz  = pi[prow+1] - pi[prow];
365       pcol = pj + pi[prow];
366       pval = pa + pi[prow];
367       for (k=0; k<pnz; k++) {
368         apa[pcol[k]] += aa[j]*pval[k];
369       }
370       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
371 #if defined(PROFILE_MatPtAPNumeric)
372       flops0 += 2.0*pnz;
373 #endif
374     }
375     aj += anz; aa += anz;
376 #if defined(PROFILE_MatPtAPNumeric)
377     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
378     time_Cseq0 += tf - t0;
379 #endif
380 
381     /* Compute P^T*A*P using outer product P[i,:]^T*AP[i,:]. */
382 #if defined(PROFILE_MatPtAPNumeric)
383     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
384 #endif
385     apj  = ptap->apj + ptap->api[i];
386     apnz = ptap->api[i+1] - ptap->api[i];
387     pnz  = pi[i+1] - pi[i];
388     pcol = pj + pi[i];
389     pval = pa + pi[i];
390 
391     /* Perform dense axpy */
392     for (j=0; j<pnz; j++) {
393       crow  = pcol[j];
394       cjj   = cj + ci[crow];
395       caj   = ca + ci[crow];
396       pvalj = pval[j];
397       cnz   = ci[crow+1] - ci[crow];
398       for (k=0; k<cnz; k++){
399         caj[k] += pvalj*apa[cjj[k]];
400       }
401       ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
402 #if defined(PROFILE_MatPtAPNumeric)
403       flops1 += 2.0*cnz;
404 #endif
405     }
406 #if defined(PROFILE_MatPtAPNumeric)
407     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
408     time_Cseq1 += tf - t0;
409 #endif
410 
411     /* Zero the current row info for A*P */
412     for (j=0; j<apnz; j++) apa[apj[j]] = 0.0;
413   }
414 
415   /* Assemble the final matrix and clean up */
416   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
417   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
418 #if defined(PROFILE_MatPtAPNumeric)
419   printf("PtAPNumeric_SeqAIJ time %g + %g, flops %d %d\n",time_Cseq0,time_Cseq1,flops0,flops1);
420 #endif
421   PetscFunctionReturn(0);
422 }
423 
424 #undef __FUNCT__
425 #define __FUNCT__ "MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy"
426 PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat A,Mat P,Mat C)
427 {
428   PetscErrorCode    ierr;
429   Mat_SeqAIJ        *a  = (Mat_SeqAIJ *) A->data;
430   Mat_SeqAIJ        *p  = (Mat_SeqAIJ *) P->data;
431   Mat_SeqAIJ        *c  = (Mat_SeqAIJ *) C->data;
432   const PetscInt    *ai=a->i,*aj=a->j,*pi=p->i,*pj=p->j,*ci=c->i,*cj=c->j;
433   const PetscScalar *aa=a->a,*pa=p->a,*pval;
434   const PetscInt    *apj,*pcol,*cjj;
435   const PetscInt    am=A->rmap->N,cm=C->rmap->N;
436   PetscInt          i,j,k,anz,apnz,pnz,prow,crow,apcol,nextap;
437   PetscScalar       *apa,*ca=c->a,*caj,pvalj;
438   Mat_PtAP          *ptap = c->ptap;
439 #if defined(PROFILE_MatPtAPNumeric)
440   PetscLogDouble    t0,tf,time_Cseq0=0.0,time_Cseq1=0.0;
441   PetscInt          flops0=0,flops1=0;
442 #endif
443 
444   PetscFunctionBegin;
445   /* Get temporary array for storage of one row of A*P */
446   apa = ptap->apa;
447 
448   /* Clear old values in C */
449   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
450 
451   for (i=0;i<am;i++) {
452     /* Form sparse row of AP[i,:] = A[i,:]*P */
453 #if defined(PROFILE_MatPtAPNumeric)
454     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
455 #endif
456     anz  = ai[i+1] - ai[i];
457     apnz = 0;
458     for (j=0; j<anz; j++) {
459       prow = aj[j];
460       pnz  = pi[prow+1] - pi[prow];
461       pcol = pj + pi[prow];
462       pval = pa + pi[prow];
463       for (k=0; k<pnz; k++) {
464         apa[pcol[k]] += aa[j]*pval[k];
465       }
466       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
467 #if defined(PROFILE_MatPtAPNumeric)
468       flops0 += 2.0*pnz;
469 #endif
470     }
471     aj += anz; aa += anz;
472 #if defined(PROFILE_MatPtAPNumeric)
473     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
474     time_Cseq0 += tf - t0;
475 #endif
476 
477     /* Compute P^T*A*P using outer product P[i,:]^T*AP[i,:]. */
478 #if defined(PROFILE_MatPtAPNumeric)
479     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
480 #endif
481     apj  = ptap->apj + ptap->api[i];
482     apnz = ptap->api[i+1] - ptap->api[i];
483     pnz  = pi[i+1] - pi[i];
484     pcol = pj + pi[i];
485     pval = pa + pi[i];
486 
487     /* Perform sparse axpy */
488     for (j=0; j<pnz; j++) {
489       crow   = pcol[j];
490       cjj    = cj + ci[crow];
491       caj    = ca + ci[crow];
492       pvalj  = pval[j];
493       nextap = 1;
494       apcol  = apj[nextap];
495       for (k=0; nextap<apnz; k++) {
496         if (cjj[k] == apcol) {
497           caj[k] += pvalj*apa[apcol];
498           apcol   = apj[nextap++];
499         }
500       }
501       ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
502 #if defined(PROFILE_MatPtAPNumeric)
503       flops1 += 2.0*apnz;
504 #endif
505     }
506 #if defined(PROFILE_MatPtAPNumeric)
507     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
508     time_Cseq1 += tf - t0;
509 #endif
510 
511     /* Zero the current row info for A*P */
512     for (j=0; j<apnz; j++) {
513       apcol      = apj[j];
514       apa[apcol] = 0.;
515     }
516   }
517 
518   /* Assemble the final matrix and clean up */
519   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
520   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
521 #if defined(PROFILE_MatPtAPNumeric)
522   printf("MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy time %g + %g, flops %d %d\n",time_Cseq0,time_Cseq1,flops0,flops1);
523 #endif
524   PetscFunctionReturn(0);
525 }
526