xref: /petsc/src/mat/impls/aij/seq/matptap.c (revision 12ce3465735c01efdf86887542af4c5dcade2f95)
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 #include "petscbt.h"
10 
11 #undef __FUNCT__
12 #define __FUNCT__ "MatPtAP"
13 /*@
14    MatPtAP - Creates the matrix projection C = P^T * A * P
15 
16    Collective on Mat
17 
18    Input Parameters:
19 +  A - the matrix
20 .  P - the projection matrix
21 .  scall - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX
22 -  fill - expected fill as ratio of nnz(C)/(nnz(A) + nnz(P))
23 
24    Output Parameters:
25 .  C - the product matrix
26 
27    Notes:
28    C will be created and must be destroyed by the user with MatDestroy().
29 
30    This routine is currently only implemented for pairs of SeqAIJ matrices and classes
31    which inherit from SeqAIJ.  C will be of type MATSEQAIJ.
32 
33    Level: intermediate
34 
35 .seealso: MatPtAPSymbolic(),MatPtAPNumeric(),MatMatMult()
36 @*/
37 PetscErrorCode MatPtAP(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
38 {
39   PetscErrorCode ierr;
40   PetscErrorCode (*fA)(Mat,Mat,MatReuse,PetscReal,Mat *);
41   PetscErrorCode (*fP)(Mat,Mat,MatReuse,PetscReal,Mat *);
42 
43   PetscFunctionBegin;
44   PetscValidHeaderSpecific(A,MAT_COOKIE,1);
45   PetscValidType(A,1);
46   MatPreallocated(A);
47   if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
48   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
49   PetscValidHeaderSpecific(P,MAT_COOKIE,2);
50   PetscValidType(P,2);
51   MatPreallocated(P);
52   if (!P->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
53   if (P->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
54   PetscValidPointer(C,3);
55 
56   if (P->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->M,A->N);
57 
58   if (fill <=0.0) SETERRQ1(PETSC_ERR_ARG_SIZ,"fill=%g must be > 0.0",fill);
59 
60   /* For now, we do not dispatch based on the type of A and P */
61   /* When implementations like _SeqAIJ_MAIJ exist, attack the multiple dispatch problem. */
62   fA = A->ops->ptap;
63   if (!fA) SETERRQ1(PETSC_ERR_SUP,"MatPtAP not supported for A of type %s",A->type_name);
64   fP = P->ops->ptap;
65   if (!fP) SETERRQ1(PETSC_ERR_SUP,"MatPtAP not supported for P of type %s",P->type_name);
66   if (fP!=fA) SETERRQ2(PETSC_ERR_ARG_INCOMP,"MatPtAP requires A, %s, to be compatible with P, %s",A->type_name,P->type_name);
67 
68   ierr = PetscLogEventBegin(MAT_PtAP,A,P,0,0);CHKERRQ(ierr);
69   ierr = (*fA)(A,P,scall,fill,C);CHKERRQ(ierr);
70   ierr = PetscLogEventEnd(MAT_PtAP,A,P,0,0);CHKERRQ(ierr);
71   PetscFunctionReturn(0);
72 }
73 
74 EXTERN PetscErrorCode MatPtAPReducedPt_MPIAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,PetscInt,PetscInt,Mat*);
75 EXTERN PetscErrorCode MatPtAPReducedPtSymbolic_MPIAIJ_SeqAIJ(Mat,Mat,PetscReal,PetscInt,PetscInt,Mat*);
76 EXTERN PetscErrorCode MatPtAPReducedPtNumeric_MPIAIJ_SeqAIJ(Mat,Mat,PetscInt,PetscInt,Mat);
77 
78 
79 EXTERN PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_Local(Mat,Mat,Mat,PetscInt,PetscReal,Mat*);
80 EXTERN PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_Local(Mat,Mat,Mat,PetscInt,Mat C);
81 
82 #undef __FUNCT__
83 #define __FUNCT__ "MatPtAP_MPIAIJ_MPIAIJ"
84 PetscErrorCode MatPtAP_MPIAIJ_MPIAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
85 {
86   PetscErrorCode    ierr;
87 
88   PetscFunctionBegin;
89   if (scall == MAT_INITIAL_MATRIX){
90     ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
91     ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ(A,P,fill,C);CHKERRQ(ierr);/* numeric product is computed as well */
92     ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
93   } else if (scall == MAT_REUSE_MATRIX){ /* not done yet! */
94     ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
95     ierr = MatDestroy(*C);CHKERRQ(ierr);
96     ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ(A,P,fill,C);CHKERRQ(ierr);
97 
98     /*
99     ierr = MatPtAPNumeric_MPIAIJ_MPIAIJ(A,P,*C);CHKERRQ(ierr);
100     */
101     ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
102   } else {
103     SETERRQ1(PETSC_ERR_ARG_WRONG,"Invalid MatReuse %d",scall);
104   }
105   PetscFunctionReturn(0);
106 }
107 #define NEW
108 #undef __FUNCT__
109 #define __FUNCT__ "MatPtAPSymbolic_MPIAIJ_MPIAIJ"
110 PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
111 {
112   PetscErrorCode    ierr;
113   Mat               P_seq,C_seq;
114   Mat               P_loc,P_oth,*ploc;
115   PetscInt          prstart,prend,m=P->m,low;
116   IS                isrow,iscol;
117 
118   PetscFunctionBegin;
119   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
120   ierr = MatGetBrowsOfAoCols_Private(A,P,MAT_INITIAL_MATRIX,PETSC_NULL,PETSC_NULL,&prstart,&P_oth);CHKERRQ(ierr);
121   prend = prstart + m;
122   /* get P_loc by taking all local rows of P */
123   ierr = MatGetOwnershipRange(P,&low,PETSC_NULL);CHKERRQ(ierr);
124   ierr = ISCreateStride(PETSC_COMM_SELF,P->m,low,1,&isrow);CHKERRQ(ierr);
125   ierr = ISCreateStride(PETSC_COMM_SELF,P->N,0,1,&iscol);CHKERRQ(ierr);
126   ierr = MatGetSubMatrices(P,1,&isrow,&iscol,MAT_INITIAL_MATRIX,&ploc);CHKERRQ(ierr);
127   P_loc = ploc[0];
128   ierr = ISDestroy(isrow);CHKERRQ(ierr);
129   ierr = ISDestroy(iscol);CHKERRQ(ierr);
130 
131   /* compute C_seq = P_loc^T * A_loc * P_seq */
132   ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ_Local(A,P_loc,P_oth,prstart,fill,&C_seq);CHKERRQ(ierr);
133   ierr = MatPtAPNumeric_MPIAIJ_MPIAIJ_Local(A,P_loc,P_oth,prstart,C_seq);CHKERRQ(ierr);
134 
135   ierr = MatDestroyMatrices(1,&ploc);CHKERRQ(ierr);
136   ierr = MatDestroy(P_oth);CHKERRQ(ierr);
137 #ifdef OLD
138   /* get P_seq = submatrix of P by taking rows of P that equal to nonzero col of A */
139   ierr = MatGetBrowsOfAcols(A,P,MAT_INITIAL_MATRIX,PETSC_NULL,PETSC_NULL,&prstart,&P_seq);CHKERRQ(ierr);
140   /* compute C_seq = P_loc^T * A_loc * P_seq */
141   prend  = prstart + m;
142   ierr = MatPtAPReducedPt_MPIAIJ_SeqAIJ(A,P_seq,MAT_INITIAL_MATRIX,fill,prstart,prend,&C_seq);CHKERRQ(ierr);
143   ierr = MatDestroy(P_seq);CHKERRQ(ierr);
144 #endif
145   /* add C_seq into mpi C */
146   ierr = MatMerge_SeqsToMPI(A->comm,C_seq,P->n,P->n,MAT_INITIAL_MATRIX,C);CHKERRQ(ierr);
147 
148   PetscFunctionReturn(0);
149 }
150 
151 #undef __FUNCT__
152 #define __FUNCT__ "MatPtAPNumeric_MPIAIJ_MPIAIJ"
153 PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
154 {
155   PetscErrorCode       ierr;
156   Mat                  P_seq,C_seq;
157   PetscInt             prstart,prend,m=P->m;
158   Mat_Merge_SeqsToMPI  *merge;
159   PetscObjectContainer container;
160 
161   PetscFunctionBegin;
162   ierr = PetscObjectQuery((PetscObject)C,"MatMergeSeqsToMPI",(PetscObject *)&container);CHKERRQ(ierr);
163   if (container) {
164     ierr  = PetscObjectContainerGetPointer(container,(void **)&merge);CHKERRQ(ierr);
165   } else {
166     SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Matrix C does not posses and object container");
167   }
168 
169   /* get P_seq = submatrix of P by taking rows of P that equal to nonzero col of A */
170   ierr = MatGetBrowsOfAcols(A,P,MAT_INITIAL_MATRIX,PETSC_NULL,PETSC_NULL,&prstart,&P_seq);CHKERRQ(ierr);
171 
172   /* compute C_seq = P_loc^T * A_loc * P_seq */
173   prend = prstart + m;
174   C_seq = merge->C_seq;
175   ierr = MatPtAPReducedPt_MPIAIJ_SeqAIJ(A,P_seq,MAT_REUSE_MATRIX,1.0,prstart,prend,&C_seq);CHKERRQ(ierr);
176   ierr = MatDestroy(P_seq);CHKERRQ(ierr);
177 
178   /* add C_seq into mpi C */
179   ierr = MatMerge_SeqsToMPI(A->comm,C_seq,P->n,P->n,MAT_REUSE_MATRIX,&C);CHKERRQ(ierr);
180 
181   PetscFunctionReturn(0);
182 }
183 
184 #undef __FUNCT__
185 #define __FUNCT__ "MatPtAP_SeqAIJ_SeqAIJ"
186 PetscErrorCode MatPtAP_SeqAIJ_SeqAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
187 {
188   PetscErrorCode ierr;
189 
190   PetscFunctionBegin;
191   if (scall == MAT_INITIAL_MATRIX){
192     ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
193     ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ(A,P,fill,C);CHKERRQ(ierr);
194     ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
195   }
196   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
197   ierr = MatPtAPNumeric_SeqAIJ_SeqAIJ(A,P,*C);CHKERRQ(ierr);
198   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
199   PetscFunctionReturn(0);
200 }
201 
202 /*
203    MatPtAPSymbolic - Creates the (i,j) structure of the matrix projection C = P^T * A * P
204 
205    Collective on Mat
206 
207    Input Parameters:
208 +  A - the matrix
209 -  P - the projection matrix
210 
211    Output Parameters:
212 .  C - the (i,j) structure of the product matrix
213 
214    Notes:
215    C will be created and must be destroyed by the user with MatDestroy().
216 
217    This routine is currently only implemented for pairs of SeqAIJ matrices and classes
218    which inherit from SeqAIJ.  C will be of type MATSEQAIJ.  The product is computed using
219    this (i,j) structure by calling MatPtAPNumeric().
220 
221    Level: intermediate
222 
223 .seealso: MatPtAP(),MatPtAPNumeric(),MatMatMultSymbolic()
224 */
225 #undef __FUNCT__
226 #define __FUNCT__ "MatPtAPSymbolic"
227 PetscErrorCode MatPtAPSymbolic(Mat A,Mat P,PetscReal fill,Mat *C)
228 {
229   PetscErrorCode ierr;
230   PetscErrorCode (*fA)(Mat,Mat,PetscReal,Mat*);
231   PetscErrorCode (*fP)(Mat,Mat,PetscReal,Mat*);
232 
233   PetscFunctionBegin;
234 
235   PetscValidHeaderSpecific(A,MAT_COOKIE,1);
236   PetscValidType(A,1);
237   MatPreallocated(A);
238   if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
239   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
240 
241   PetscValidHeaderSpecific(P,MAT_COOKIE,2);
242   PetscValidType(P,2);
243   MatPreallocated(P);
244   if (!P->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
245   if (P->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
246 
247   PetscValidPointer(C,3);
248 
249   if (P->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->M,A->N);
250   if (A->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix 'A' must be square, %d != %d",A->M,A->N);
251 
252   /* For now, we do not dispatch based on the type of A and P */
253   /* When implementations like _SeqAIJ_MAIJ exist, attack the multiple dispatch problem. */
254   fA = A->ops->ptapsymbolic;
255   if (!fA) SETERRQ1(PETSC_ERR_SUP,"MatPtAPSymbolic not supported for A of type %s",A->type_name);
256   fP = P->ops->ptapsymbolic;
257   if (!fP) SETERRQ1(PETSC_ERR_SUP,"MatPtAPSymbolic not supported for P of type %s",P->type_name);
258   if (fP!=fA) SETERRQ2(PETSC_ERR_ARG_INCOMP,"MatPtAPSymbolic requires A, %s, to be compatible with P, %s",A->type_name,P->type_name);
259 
260   ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
261   ierr = (*fA)(A,P,fill,C);CHKERRQ(ierr);
262   ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
263 
264   PetscFunctionReturn(0);
265 }
266 
267 typedef struct {
268   Mat    symAP;
269 } Mat_PtAPstruct;
270 
271 EXTERN PetscErrorCode MatDestroy_SeqAIJ(Mat);
272 
273 #undef __FUNCT__
274 #define __FUNCT__ "MatDestroy_SeqAIJ_PtAP"
275 PetscErrorCode MatDestroy_SeqAIJ_PtAP(Mat A)
276 {
277   PetscErrorCode    ierr;
278   Mat_PtAPstruct    *ptap=(Mat_PtAPstruct*)A->spptr;
279 
280   PetscFunctionBegin;
281   ierr = MatDestroy(ptap->symAP);CHKERRQ(ierr);
282   ierr = PetscFree(ptap);CHKERRQ(ierr);
283   ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr);
284   PetscFunctionReturn(0);
285 }
286 
287 #undef __FUNCT__
288 #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqAIJ"
289 PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
290 {
291   PetscErrorCode ierr;
292   FreeSpaceList  free_space=PETSC_NULL,current_space=PETSC_NULL;
293   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data,*p = (Mat_SeqAIJ*)P->data,*c;
294   PetscInt       *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj;
295   PetscInt       *ci,*cj,*ptadenserow,*ptasparserow,*ptaj;
296   PetscInt       an=A->N,am=A->M,pn=P->N,pm=P->M;
297   PetscInt       i,j,k,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi,nlnk,*lnk;
298   MatScalar      *ca;
299   PetscBT        lnkbt;
300 
301   PetscFunctionBegin;
302   /* Get ij structure of P^T */
303   ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
304   ptJ=ptj;
305 
306   /* Allocate ci array, arrays for fill computation and */
307   /* free space for accumulating nonzero column info */
308   ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
309   ci[0] = 0;
310 
311   ierr = PetscMalloc((2*an+1)*sizeof(PetscInt),&ptadenserow);CHKERRQ(ierr);
312   ierr = PetscMemzero(ptadenserow,(2*an+1)*sizeof(PetscInt));CHKERRQ(ierr);
313   ptasparserow = ptadenserow  + an;
314 
315   /* create and initialize a linked list */
316   nlnk = pn+1;
317   ierr = PetscLLCreate(pn,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr);
318 
319   /* Set initial free space to be nnz(A) scaled by aspect ratio of P. */
320   /* This should be reasonable if sparsity of PtAP is similar to that of A. */
321   ierr          = GetMoreSpace((ai[am]/pm)*pn,&free_space);
322   current_space = free_space;
323 
324   /* Determine symbolic info for each row of C: */
325   for (i=0;i<pn;i++) {
326     ptnzi  = pti[i+1] - pti[i];
327     ptanzi = 0;
328     /* Determine symbolic row of PtA: */
329     for (j=0;j<ptnzi;j++) {
330       arow = *ptJ++;
331       anzj = ai[arow+1] - ai[arow];
332       ajj  = aj + ai[arow];
333       for (k=0;k<anzj;k++) {
334         if (!ptadenserow[ajj[k]]) {
335           ptadenserow[ajj[k]]    = -1;
336           ptasparserow[ptanzi++] = ajj[k];
337         }
338       }
339     }
340       /* Using symbolic info for row of PtA, determine symbolic info for row of C: */
341     ptaj = ptasparserow;
342     cnzi   = 0;
343     for (j=0;j<ptanzi;j++) {
344       prow = *ptaj++;
345       pnzj = pi[prow+1] - pi[prow];
346       pjj  = pj + pi[prow];
347       /* add non-zero cols of P into the sorted linked list lnk */
348       ierr = PetscLLAdd(pnzj,pjj,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr);
349       cnzi += nlnk;
350     }
351 
352     /* If free space is not available, make more free space */
353     /* Double the amount of total space in the list */
354     if (current_space->local_remaining<cnzi) {
355       ierr = GetMoreSpace(current_space->total_array_size,&current_space);CHKERRQ(ierr);
356     }
357 
358     /* Copy data into free space, and zero out denserows */
359     ierr = PetscLLClean(pn,pn,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
360     current_space->array           += cnzi;
361     current_space->local_used      += cnzi;
362     current_space->local_remaining -= cnzi;
363 
364     for (j=0;j<ptanzi;j++) {
365       ptadenserow[ptasparserow[j]] = 0;
366     }
367     /* Aside: Perhaps we should save the pta info for the numerical factorization. */
368     /*        For now, we will recompute what is needed. */
369     ci[i+1] = ci[i] + cnzi;
370   }
371   /* nnz is now stored in ci[ptm], column indices are in the list of free space */
372   /* Allocate space for cj, initialize cj, and */
373   /* destroy list of free space and other temporary array(s) */
374   ierr = PetscMalloc((ci[pn]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
375   ierr = MakeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
376   ierr = PetscFree(ptadenserow);CHKERRQ(ierr);
377   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
378 
379   /* Allocate space for ca */
380   ierr = PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
381   ierr = PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));CHKERRQ(ierr);
382 
383   /* put together the new matrix */
384   ierr = MatCreateSeqAIJWithArrays(A->comm,pn,pn,ci,cj,ca,C);CHKERRQ(ierr);
385 
386   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
387   /* Since these are PETSc arrays, change flags to free them as necessary. */
388   c = (Mat_SeqAIJ *)((*C)->data);
389   c->freedata = PETSC_TRUE;
390   c->nonew    = 0;
391 
392   /* Clean up. */
393   ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
394 
395   PetscFunctionReturn(0);
396 }
397 
398 #include "src/mat/impls/maij/maij.h"
399 EXTERN_C_BEGIN
400 #undef __FUNCT__
401 #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqMAIJ"
402 PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqMAIJ(Mat A,Mat PP,Mat *C)
403 {
404   /* This routine requires testing -- I don't think it works. */
405   PetscErrorCode ierr;
406   FreeSpaceList  free_space=PETSC_NULL,current_space=PETSC_NULL;
407   Mat_SeqMAIJ    *pp=(Mat_SeqMAIJ*)PP->data;
408   Mat            P=pp->AIJ;
409   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*p=(Mat_SeqAIJ*)P->data,*c;
410   PetscInt       *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj;
411   PetscInt       *ci,*cj,*denserow,*sparserow,*ptadenserow,*ptasparserow,*ptaj;
412   PetscInt       an=A->N,am=A->M,pn=P->N,pm=P->M,ppdof=pp->dof;
413   PetscInt       i,j,k,dof,pdof,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi;
414   MatScalar      *ca;
415 
416   PetscFunctionBegin;
417   /* Start timer */
418   ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,PP,0,0);CHKERRQ(ierr);
419 
420   /* Get ij structure of P^T */
421   ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
422 
423   /* Allocate ci array, arrays for fill computation and */
424   /* free space for accumulating nonzero column info */
425   ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
426   ci[0] = 0;
427 
428   ierr = PetscMalloc((2*pn+2*an+1)*sizeof(PetscInt),&ptadenserow);CHKERRQ(ierr);
429   ierr = PetscMemzero(ptadenserow,(2*pn+2*an+1)*sizeof(PetscInt));CHKERRQ(ierr);
430   ptasparserow = ptadenserow  + an;
431   denserow     = ptasparserow + an;
432   sparserow    = denserow     + pn;
433 
434   /* Set initial free space to be nnz(A) scaled by aspect ratio of P. */
435   /* This should be reasonable if sparsity of PtAP is similar to that of A. */
436   ierr          = GetMoreSpace((ai[am]/pm)*pn,&free_space);
437   current_space = free_space;
438 
439   /* Determine symbolic info for each row of C: */
440   for (i=0;i<pn/ppdof;i++) {
441     ptnzi  = pti[i+1] - pti[i];
442     ptanzi = 0;
443     ptJ    = ptj + pti[i];
444     for (dof=0;dof<ppdof;dof++) {
445     /* Determine symbolic row of PtA: */
446       for (j=0;j<ptnzi;j++) {
447         arow = ptJ[j] + dof;
448         anzj = ai[arow+1] - ai[arow];
449         ajj  = aj + ai[arow];
450         for (k=0;k<anzj;k++) {
451           if (!ptadenserow[ajj[k]]) {
452             ptadenserow[ajj[k]]    = -1;
453             ptasparserow[ptanzi++] = ajj[k];
454           }
455         }
456       }
457       /* Using symbolic info for row of PtA, determine symbolic info for row of C: */
458       ptaj = ptasparserow;
459       cnzi   = 0;
460       for (j=0;j<ptanzi;j++) {
461         pdof = *ptaj%dof;
462         prow = (*ptaj++)/dof;
463         pnzj = pi[prow+1] - pi[prow];
464         pjj  = pj + pi[prow];
465         for (k=0;k<pnzj;k++) {
466           if (!denserow[pjj[k]+pdof]) {
467             denserow[pjj[k]+pdof] = -1;
468             sparserow[cnzi++]     = pjj[k]+pdof;
469           }
470         }
471       }
472 
473       /* sort sparserow */
474       ierr = PetscSortInt(cnzi,sparserow);CHKERRQ(ierr);
475 
476       /* If free space is not available, make more free space */
477       /* Double the amount of total space in the list */
478       if (current_space->local_remaining<cnzi) {
479         ierr = GetMoreSpace(current_space->total_array_size,&current_space);CHKERRQ(ierr);
480       }
481 
482       /* Copy data into free space, and zero out denserows */
483       ierr = PetscMemcpy(current_space->array,sparserow,cnzi*sizeof(PetscInt));CHKERRQ(ierr);
484       current_space->array           += cnzi;
485       current_space->local_used      += cnzi;
486       current_space->local_remaining -= cnzi;
487 
488       for (j=0;j<ptanzi;j++) {
489         ptadenserow[ptasparserow[j]] = 0;
490       }
491       for (j=0;j<cnzi;j++) {
492         denserow[sparserow[j]] = 0;
493       }
494       /* Aside: Perhaps we should save the pta info for the numerical factorization. */
495       /*        For now, we will recompute what is needed. */
496       ci[i+1+dof] = ci[i+dof] + cnzi;
497     }
498   }
499   /* nnz is now stored in ci[ptm], column indices are in the list of free space */
500   /* Allocate space for cj, initialize cj, and */
501   /* destroy list of free space and other temporary array(s) */
502   ierr = PetscMalloc((ci[pn]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
503   ierr = MakeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
504   ierr = PetscFree(ptadenserow);CHKERRQ(ierr);
505 
506   /* Allocate space for ca */
507   ierr = PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
508   ierr = PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));CHKERRQ(ierr);
509 
510   /* put together the new matrix */
511   ierr = MatCreateSeqAIJWithArrays(A->comm,pn,pn,ci,cj,ca,C);CHKERRQ(ierr);
512 
513   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
514   /* Since these are PETSc arrays, change flags to free them as necessary. */
515   c = (Mat_SeqAIJ *)((*C)->data);
516   c->freedata = PETSC_TRUE;
517   c->nonew    = 0;
518 
519   /* Clean up. */
520   ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
521 
522   ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,PP,0,0);CHKERRQ(ierr);
523   PetscFunctionReturn(0);
524 }
525 EXTERN_C_END
526 
527 /*
528    MatPtAPNumeric - Computes the matrix projection C = P^T * A * P
529 
530    Collective on Mat
531 
532    Input Parameters:
533 +  A - the matrix
534 -  P - the projection matrix
535 
536    Output Parameters:
537 .  C - the product matrix
538 
539    Notes:
540    C must have been created by calling MatPtAPSymbolic and must be destroyed by
541    the user using MatDeatroy().
542 
543    This routine is currently only implemented for pairs of AIJ matrices and classes
544    which inherit from AIJ.  C will be of type MATAIJ.
545 
546    Level: intermediate
547 
548 .seealso: MatPtAP(),MatPtAPSymbolic(),MatMatMultNumeric()
549 */
550 #undef __FUNCT__
551 #define __FUNCT__ "MatPtAPNumeric"
552 PetscErrorCode MatPtAPNumeric(Mat A,Mat P,Mat C)
553 {
554   PetscErrorCode ierr;
555   PetscErrorCode (*fA)(Mat,Mat,Mat);
556   PetscErrorCode (*fP)(Mat,Mat,Mat);
557 
558   PetscFunctionBegin;
559 
560   PetscValidHeaderSpecific(A,MAT_COOKIE,1);
561   PetscValidType(A,1);
562   MatPreallocated(A);
563   if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
564   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
565 
566   PetscValidHeaderSpecific(P,MAT_COOKIE,2);
567   PetscValidType(P,2);
568   MatPreallocated(P);
569   if (!P->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
570   if (P->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
571 
572   PetscValidHeaderSpecific(C,MAT_COOKIE,3);
573   PetscValidType(C,3);
574   MatPreallocated(C);
575   if (!C->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
576   if (C->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
577 
578   if (P->N!=C->M) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->N,C->M);
579   if (P->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->M,A->N);
580   if (A->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix 'A' must be square, %d != %d",A->M,A->N);
581   if (P->N!=C->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->N,C->N);
582 
583   /* For now, we do not dispatch based on the type of A and P */
584   /* When implementations like _SeqAIJ_MAIJ exist, attack the multiple dispatch problem. */
585   fA = A->ops->ptapnumeric;
586   if (!fA) SETERRQ1(PETSC_ERR_SUP,"MatPtAPNumeric not supported for A of type %s",A->type_name);
587   fP = P->ops->ptapnumeric;
588   if (!fP) SETERRQ1(PETSC_ERR_SUP,"MatPtAPNumeric not supported for P of type %s",P->type_name);
589   if (fP!=fA) SETERRQ2(PETSC_ERR_ARG_INCOMP,"MatPtAPNumeric requires A, %s, to be compatible with P, %s",A->type_name,P->type_name);
590 
591   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
592   ierr = (*fA)(A,P,C);CHKERRQ(ierr);
593   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
594 
595   PetscFunctionReturn(0);
596 }
597 
598 #undef __FUNCT__
599 #define __FUNCT__ "MatPtAPNumeric_SeqAIJ_SeqAIJ"
600 PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat A,Mat P,Mat C)
601 {
602   PetscErrorCode ierr;
603   PetscInt       flops=0;
604   Mat_SeqAIJ     *a  = (Mat_SeqAIJ *) A->data;
605   Mat_SeqAIJ     *p  = (Mat_SeqAIJ *) P->data;
606   Mat_SeqAIJ     *c  = (Mat_SeqAIJ *) C->data;
607   PetscInt       *ai=a->i,*aj=a->j,*apj,*apjdense,*pi=p->i,*pj=p->j,*pJ=p->j,*pjj;
608   PetscInt       *ci=c->i,*cj=c->j,*cjj;
609   PetscInt       am=A->M,cn=C->N,cm=C->M;
610   PetscInt       i,j,k,anzi,pnzi,apnzj,nextap,pnzj,prow,crow;
611   MatScalar      *aa=a->a,*apa,*pa=p->a,*pA=p->a,*paj,*ca=c->a,*caj;
612 
613   PetscFunctionBegin;
614   /* Allocate temporary array for storage of one row of A*P */
615   ierr = PetscMalloc(cn*(sizeof(MatScalar)+2*sizeof(PetscInt)),&apa);CHKERRQ(ierr);
616   ierr = PetscMemzero(apa,cn*(sizeof(MatScalar)+2*sizeof(PetscInt)));CHKERRQ(ierr);
617 
618   apj      = (PetscInt *)(apa + cn);
619   apjdense = apj + cn;
620 
621   /* Clear old values in C */
622   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
623 
624   for (i=0;i<am;i++) {
625     /* Form sparse row of A*P */
626     anzi  = ai[i+1] - ai[i];
627     apnzj = 0;
628     for (j=0;j<anzi;j++) {
629       prow = *aj++;
630       pnzj = pi[prow+1] - pi[prow];
631       pjj  = pj + pi[prow];
632       paj  = pa + pi[prow];
633       for (k=0;k<pnzj;k++) {
634         if (!apjdense[pjj[k]]) {
635           apjdense[pjj[k]] = -1;
636           apj[apnzj++]     = pjj[k];
637         }
638         apa[pjj[k]] += (*aa)*paj[k];
639       }
640       flops += 2*pnzj;
641       aa++;
642     }
643 
644     /* Sort the j index array for quick sparse axpy. */
645     ierr = PetscSortInt(apnzj,apj);CHKERRQ(ierr);
646 
647     /* Compute P^T*A*P using outer product (P^T)[:,j]*(A*P)[j,:]. */
648     pnzi = pi[i+1] - pi[i];
649     for (j=0;j<pnzi;j++) {
650       nextap = 0;
651       crow   = *pJ++;
652       cjj    = cj + ci[crow];
653       caj    = ca + ci[crow];
654       /* Perform sparse axpy operation.  Note cjj includes apj. */
655       for (k=0;nextap<apnzj;k++) {
656         if (cjj[k]==apj[nextap]) {
657           caj[k] += (*pA)*apa[apj[nextap++]];
658         }
659       }
660       flops += 2*apnzj;
661       pA++;
662     }
663 
664     /* Zero the current row info for A*P */
665     for (j=0;j<apnzj;j++) {
666       apa[apj[j]]      = 0.;
667       apjdense[apj[j]] = 0;
668     }
669   }
670 
671   /* Assemble the final matrix and clean up */
672   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
673   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
674   ierr = PetscFree(apa);CHKERRQ(ierr);
675   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
676 
677   PetscFunctionReturn(0);
678 }
679 
680 /* Compute local C = P[rstart:rend,:]^T * A * P - used by MatPtAP_MPIAIJ_MPIAIJ() */
681 
682 #undef __FUNCT__
683 #define __FUNCT__ "MatPtAPReducedPtNumeric_MPIAIJ_SeqAIJ"
684 PetscErrorCode MatPtAPReducedPtNumeric_MPIAIJ_SeqAIJ(Mat A,Mat P,PetscInt prstart,PetscInt prend,Mat C)
685 {
686   PetscErrorCode ierr;
687   PetscInt       flops=0;
688   Mat_MPIAIJ     *a=(Mat_MPIAIJ*)A->data;
689   Mat_SeqAIJ     *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
690   Mat_SeqAIJ     *p  = (Mat_SeqAIJ *) P->data;
691   Mat_SeqAIJ     *c  = (Mat_SeqAIJ *) C->data;
692   PetscInt       *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,*apj,*apjdense,cstart=a->cstart,cend=a->cend,col;
693   PetscInt       *pi=p->i,*pj=p->j,*pJ=p->j+pi[prstart],*pjj;
694   PetscInt       *ci=c->i,*cj=c->j,*cjj;
695   PetscInt       am=A->m,cn=C->N,cm=C->M;
696   PetscInt       i,j,k,anzi,pnzi,apnzj,nextap,pnzj,prow,crow;
697   MatScalar      *ada=ad->a,*aoa=ao->a,*apa,*pa=p->a,*pA=p->a,*paj,*ca=c->a,*caj;
698 
699   PetscFunctionBegin;
700   pA=p->a+pi[prstart];
701   /* Allocate temporary array for storage of one row of A*P */
702   ierr = PetscMalloc(cn*(sizeof(MatScalar)+2*sizeof(PetscInt)),&apa);CHKERRQ(ierr);
703   ierr = PetscMemzero(apa,cn*(sizeof(MatScalar)+2*sizeof(PetscInt)));CHKERRQ(ierr);
704 
705   apj      = (PetscInt *)(apa + cn);
706   apjdense = apj + cn;
707 
708   /* Clear old values in C */
709   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
710 
711   for (i=0;i<am;i++) {
712     /* Form sparse row of A*P */
713     /* diagonal portion of A */
714     anzi  = adi[i+1] - adi[i];
715     apnzj = 0;
716     for (j=0;j<anzi;j++) {
717       prow = *adj + prstart;
718       adj++;
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]] += (*ada)*paj[k];
728       }
729       flops += 2*pnzj;
730       ada++;
731     }
732     /* off-diagonal portion of A */
733     anzi  = aoi[i+1] - aoi[i];
734     for (j=0;j<anzi;j++) {
735       col = a->garray[*aoj];
736       if (col < cstart){
737         prow = *aoj;
738       } else if (col >= cend){
739         prow = *aoj + prend-prstart;
740       } else {
741         SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Off-diagonal portion of A has wrong column map");
742       }
743       aoj++;
744       pnzj = pi[prow+1] - pi[prow];
745       pjj  = pj + pi[prow];
746       paj  = pa + pi[prow];
747       for (k=0;k<pnzj;k++) {
748         if (!apjdense[pjj[k]]) {
749           apjdense[pjj[k]] = -1;
750           apj[apnzj++]     = pjj[k];
751         }
752         apa[pjj[k]] += (*aoa)*paj[k];
753       }
754       flops += 2*pnzj;
755       aoa++;
756     }
757 
758     /* Sort the j index array for quick sparse axpy. */
759     ierr = PetscSortInt(apnzj,apj);CHKERRQ(ierr);
760 
761     /* Compute P[[prstart:prend,:]^T*A*P using outer product (P^T)[:,j+prstart]*(A*P)[j,:]. */
762     pnzi = pi[i+1+prstart] - pi[i+prstart];
763     for (j=0;j<pnzi;j++) {
764       nextap = 0;
765       crow   = *pJ++;
766       cjj    = cj + ci[crow];
767       caj    = ca + ci[crow];
768       /* Perform sparse axpy operation.  Note cjj includes apj. */
769       for (k=0;nextap<apnzj;k++) {
770         if (cjj[k]==apj[nextap]) {
771           caj[k] += (*pA)*apa[apj[nextap++]];
772         }
773       }
774       flops += 2*apnzj;
775       pA++;
776     }
777 
778     /* Zero the current row info for A*P */
779     for (j=0;j<apnzj;j++) {
780       apa[apj[j]]      = 0.;
781       apjdense[apj[j]] = 0;
782     }
783   }
784 
785   /* Assemble the final matrix and clean up */
786   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
787   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
788   ierr = PetscFree(apa);CHKERRQ(ierr);
789   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
790 
791   PetscFunctionReturn(0);
792 }
793 
794 #undef __FUNCT__
795 #define __FUNCT__ "MatPtAPReducedPtSymbolic_MPIAIJ_SeqAIJ"
796 PetscErrorCode MatPtAPReducedPtSymbolic_MPIAIJ_SeqAIJ(Mat A,Mat P,PetscReal fill,PetscInt prstart,PetscInt prend,Mat *C)
797 {
798   PetscErrorCode ierr;
799   FreeSpaceList  free_space=PETSC_NULL,current_space=PETSC_NULL;
800   Mat_MPIAIJ     *a=(Mat_MPIAIJ*)A->data;
801   Mat_SeqAIJ     *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
802   Mat_SeqAIJ     *p=(Mat_SeqAIJ*)P->data,*c;
803   PetscInt       *pti,*ptj,*ptJ,*ajj,*pi=p->i,*pj=p->j,*pjj;
804   PetscInt       *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,nnz,cstart=a->cstart,cend=a->cend,col;
805   PetscInt       *ci,*cj,*ptadenserow,*ptasparserow,*ptaj;
806   PetscInt       an=A->N,am=A->m,pn=P->N,pm=P->M;
807   PetscInt       i,j,k,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi;
808   PetscInt       m=prend-prstart,nlnk,*lnk;
809   MatScalar      *ca;
810   PetscBT        lnkbt;
811 
812   PetscFunctionBegin;
813   int rank;
814   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
815 
816   /* Get ij structure of P[rstart:rend,:]^T */
817   ierr = MatGetSymbolicTransposeReduced_SeqAIJ(P,prstart,prend,&pti,&ptj);CHKERRQ(ierr);
818   ptJ=ptj;
819 
820   /* Allocate ci array, arrays for fill computation and */
821   /* free space for accumulating nonzero column info */
822   ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
823   ci[0] = 0;
824 
825   ierr = PetscMalloc((2*an+1)*sizeof(PetscInt),&ptadenserow);CHKERRQ(ierr);
826   ierr = PetscMemzero(ptadenserow,(2*an+1)*sizeof(PetscInt));CHKERRQ(ierr);
827   ptasparserow = ptadenserow  + an;
828 
829   /* create and initialize a linked list */
830   nlnk = pn+1;
831   ierr = PetscLLCreate(pn,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr);
832 
833   /* Set initial free space to be nnz(A) scaled by aspect ratio of P. */
834   /* This should be reasonable if sparsity of PtAP is similar to that of A. */
835   nnz           = adi[am] + aoi[am];
836   ierr          = GetMoreSpace((PetscInt)(fill*nnz/pm)*pn,&free_space);
837   current_space = free_space;
838 
839   /* Determine symbolic info for each row of C: */
840   for (i=0;i<pn;i++) {
841     ptnzi  = pti[i+1] - pti[i];
842     ptanzi = 0;
843     /* Determine symbolic row of PtA_reduced: */
844     for (j=0;j<ptnzi;j++) {
845       arow = *ptJ++;
846       if (rank==1) printf("i: %d, ptnzi: %d, arow: %d\n",i,ptnzi,arow);
847       /* diagonal portion of A */
848       anzj = adi[arow+1] - adi[arow];
849       ajj  = adj + adi[arow];
850       for (k=0;k<anzj;k++) {
851         col = ajj[k]+prstart;
852         if (!ptadenserow[col]) {
853           ptadenserow[col]    = -1;
854           ptasparserow[ptanzi++] = col;
855         }
856       }
857       /* off-diagonal portion of A */
858       anzj = aoi[arow+1] - aoi[arow];
859       ajj  = aoj + aoi[arow];
860       for (k=0;k<anzj;k++) {
861         col = a->garray[ajj[k]];  /* global col */
862         if (col < cstart){
863           col = ajj[k];
864         } else if (col >= cend){
865           col = ajj[k] + m;
866         } else {
867           SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Off-diagonal portion of A has wrong column map");
868         }
869         if (!ptadenserow[col]) {
870           ptadenserow[col]    = -1;
871           ptasparserow[ptanzi++] = col;
872         }
873       }
874     }
875 
876     printf(" [%d] i: %d, ptanzi: %d \n",rank,i,ptanzi);
877     /* Using symbolic info for row of PtA, determine symbolic info for row of C: */
878     ptaj = ptasparserow;
879     cnzi   = 0;
880     for (j=0;j<ptanzi;j++) {
881       prow = *ptaj++;
882       pnzj = pi[prow+1] - pi[prow];
883       pjj  = pj + pi[prow];
884       /* if (rank==1) printf(" prow: %d, pnzj: %d\n",prow,pnzj); */
885       /* add non-zero cols of P into the sorted linked list lnk */
886       ierr = PetscLLAdd(pnzj,pjj,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr);
887       cnzi += nlnk;
888     }
889 
890     /* If free space is not available, make more free space */
891     /* Double the amount of total space in the list */
892     if (current_space->local_remaining<cnzi) {
893       ierr = GetMoreSpace(current_space->total_array_size,&current_space);CHKERRQ(ierr);
894     }
895 
896     /* Copy data into free space, and zero out denserows */
897     ierr = PetscLLClean(pn,pn,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
898     current_space->array           += cnzi;
899     current_space->local_used      += cnzi;
900     current_space->local_remaining -= cnzi;
901 
902     for (j=0;j<ptanzi;j++) {
903       ptadenserow[ptasparserow[j]] = 0;
904     }
905 
906     /* Aside: Perhaps we should save the pta info for the numerical factorization. */
907     /*        For now, we will recompute what is needed. */
908     ci[i+1] = ci[i] + cnzi;
909   }
910   /* nnz is now stored in ci[ptm], column indices are in the list of free space */
911   /* Allocate space for cj, initialize cj, and */
912   /* destroy list of free space and other temporary array(s) */
913   ierr = PetscMalloc((ci[pn]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
914   ierr = MakeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
915   ierr = PetscFree(ptadenserow);CHKERRQ(ierr);
916   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
917 
918   /* Allocate space for ca */
919   ierr = PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
920   ierr = PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));CHKERRQ(ierr);
921 
922   /* put together the new matrix */
923   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,pn,pn,ci,cj,ca,C);CHKERRQ(ierr);
924 
925   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
926   /* Since these are PETSc arrays, change flags to free them as necessary. */
927   c = (Mat_SeqAIJ *)((*C)->data);
928   c->freedata = PETSC_TRUE;
929   c->nonew    = 0;
930 
931   /* Clean up. */
932   ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
933 
934   PetscFunctionReturn(0);
935 }
936 
937 /*@C
938    MatPtAPReducedPt - Creates C = P[rstart:rend,:]^T * A * P of seqaij matrices,
939                 used by MatPtAP_MPIAIJ_MPIAIJ()
940 
941    Collective on Mat
942 
943    Input Parameters:
944 +  A - the matrix in mpiaij format
945 .  P - the projection matrix in seqaij format
946 .  scall - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX
947 .  fill - expected fill (not being used when scall==MAT_REUSE_MATRIX)
948 .  prstart, prend - the starting and ending-1-th row of the matrix P to be used for transpose
949 
950    Output Parameters:
951 .  C - the product matrix in seqaij format
952 
953    Level: developer
954 
955 .seealso: MatPtAPSymbolic(),MatPtAPNumeric(),MatMatMult()
956 @*/
957 static PetscEvent logkey_matptapreducedpt = 0;
958 #undef __FUNCT__
959 #define __FUNCT__ "MatPtAPReducedPt_MPIAIJ_SeqAIJ"
960 PetscErrorCode MatPtAPReducedPt_MPIAIJ_SeqAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,PetscInt prstart,PetscInt prend,Mat *C)
961 {
962   PetscErrorCode ierr;
963 
964   PetscFunctionBegin;
965   if (A->m != prend-prstart) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",A->m,prend-prstart);
966   if (prend-prstart > P->m) SETERRQ2(PETSC_ERR_ARG_SIZ," prend-prstart %d cannot be larger than P->m %d",prend-prstart,P->m);
967   if (!logkey_matptapreducedpt) {
968     ierr = PetscLogEventRegister(&logkey_matptapreducedpt,"MatPtAP_ReducedPt",MAT_COOKIE);
969   }
970   PetscLogEventBegin(logkey_matptapreducedpt,A,P,0,0);CHKERRQ(ierr);
971 
972   if (scall == MAT_INITIAL_MATRIX){
973     ierr = MatPtAPReducedPtSymbolic_MPIAIJ_SeqAIJ(A,P,fill,prstart,prend,C);CHKERRQ(ierr);
974   }
975   ierr = MatPtAPReducedPtNumeric_MPIAIJ_SeqAIJ(A,P,prstart,prend,*C);CHKERRQ(ierr);
976   ierr = PetscLogEventEnd(logkey_matptapreducedpt,A,P,0,0);CHKERRQ(ierr);
977 
978   PetscFunctionReturn(0);
979 }
980 
981 /*@C
982     MatGetBrowsOfAoCols - Creates a SeqAIJ matrix by taking rows of B that equal to nonzero columns
983      of the off-diagonal portion of A
984 
985     Collective on Mat
986 
987    Input Parameters:
988 +    A,B - the matrices in mpiaij format
989 .    scall - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX
990 -    rowb, colb - index sets of rows and columns of B to extract (or PETSC_NULL)
991 
992    Output Parameter:
993 +    rowb, colb - index sets of rows and columns of B to extract
994 .    brstart - row index of B_seq from which next B->m rows are taken from B's local rows
995 -    B_seq - the sequential matrix generated
996 
997     Level: developer
998 
999 @*/
1000 static PetscEvent logkey_GetBrowsOfAocols = 0;
1001 #undef __FUNCT__
1002 #define __FUNCT__ "MatGetBrowsOfAoCols_Private"
1003 PetscErrorCode MatGetBrowsOfAoCols_Private(Mat A,Mat B,MatReuse scall,IS *rowb,IS *colb,PetscInt *brstart,Mat *B_seq)
1004 {
1005   Mat_MPIAIJ        *a=(Mat_MPIAIJ*)A->data,*b=(Mat_MPIAIJ*)B->data;
1006   PetscErrorCode    ierr;
1007   PetscInt          *idx,i,start,ncols,nzB,*cmap,imark;
1008   IS                isrowb,iscolb;
1009   Mat               *bseq;
1010 
1011   PetscFunctionBegin;
1012   if (a->cstart != b->rstart || a->cend != b->rend){
1013     SETERRQ4(PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, (%d, %d) != (%d,%d)",a->cstart,a->cend,b->rstart,b->rend);
1014   }
1015   if (!logkey_GetBrowsOfAocols) {
1016     ierr = PetscLogEventRegister(&logkey_GetBrowsOfAocols,"MatGetBrowsOfAoCols",MAT_COOKIE);
1017   }
1018   ierr = PetscLogEventBegin(logkey_GetBrowsOfAocols,A,B,0,0);CHKERRQ(ierr);
1019 
1020   if (scall == MAT_INITIAL_MATRIX){
1021     start = a->cstart;
1022     cmap  = a->garray;
1023     nzB   = a->B->n;
1024     if (!nzB){ /* output a null matrix */
1025       B_seq = PETSC_NULL;
1026       ierr = PetscLogEventEnd(logkey_GetBrowsOfAocols,A,B,0,0);CHKERRQ(ierr);
1027       PetscFunctionReturn(0);
1028     }
1029     ierr  = PetscMalloc(nzB*sizeof(PetscInt), &idx);CHKERRQ(ierr);
1030     ncols = 0;
1031     for (i=0; i<nzB; i++) {  /* row < local row index */
1032       if (cmap[i] < start) idx[ncols++] = cmap[i];
1033       else break;
1034     }
1035     imark = i;
1036     /* for (i=0; i<nzA; i++) idx[ncols++] = start + i; */ /* local rows */
1037     for (i=imark; i<nzB; i++) idx[ncols++] = cmap[i]; /* row > local row index */
1038     ierr = ISCreateGeneral(PETSC_COMM_SELF,ncols,idx,&isrowb);CHKERRQ(ierr);
1039     ierr = PetscFree(idx);CHKERRQ(ierr);
1040     *brstart = imark;
1041     ierr = ISCreateStride(PETSC_COMM_SELF,B->N,0,1,&iscolb);CHKERRQ(ierr);
1042   } else {
1043     if (!rowb || !colb) SETERRQ(PETSC_ERR_SUP,"IS rowb and colb must be provided for MAT_REUSE_MATRIX");
1044     isrowb = *rowb; iscolb = *colb;
1045     ierr = PetscMalloc(sizeof(Mat),&bseq);CHKERRQ(ierr);
1046     bseq[0] = *B_seq;
1047   }
1048   ierr = MatGetSubMatrices(B,1,&isrowb,&iscolb,scall,&bseq);CHKERRQ(ierr);
1049   *B_seq = bseq[0];
1050   ierr = PetscFree(bseq);CHKERRQ(ierr);
1051   if (!rowb){
1052     ierr = ISDestroy(isrowb);CHKERRQ(ierr);
1053   } else {
1054     *rowb = isrowb;
1055   }
1056   if (!colb){
1057     ierr = ISDestroy(iscolb);CHKERRQ(ierr);
1058   } else {
1059     *colb = iscolb;
1060   }
1061   ierr = PetscLogEventEnd(logkey_GetBrowsOfAocols,A,B,0,0);CHKERRQ(ierr);
1062   PetscFunctionReturn(0);
1063 }
1064 
1065 /* ------------- New --------------------*/
1066 static PetscEvent logkey_matptapnumeric_local = 0;
1067 #undef __FUNCT__
1068 #define __FUNCT__ "MatPtAPNumeric_MPIAIJ_MPIAIJ_Local"
1069 PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_Local(Mat A,Mat P_loc,Mat P_oth,PetscInt prstart,Mat C)
1070 {
1071   PetscErrorCode ierr;
1072   PetscInt       flops=0;
1073   Mat_MPIAIJ     *a=(Mat_MPIAIJ*)A->data;
1074   Mat_SeqAIJ     *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
1075   Mat_SeqAIJ     *c  = (Mat_SeqAIJ *) C->data;
1076   Mat_SeqAIJ     *p_loc=(Mat_SeqAIJ*)P_loc->data,*p_oth=(Mat_SeqAIJ*)P_oth->data;
1077   PetscInt       *pi_loc=p_loc->i,*pj_loc=p_loc->j,*pi_oth=p_oth->i,*pj_oth=p_oth->j;
1078   PetscInt       *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,*apj,*apjdense,cstart=a->cstart,cend=a->cend,col;
1079   PetscInt       *pJ=pj_loc,*pjj;
1080   PetscInt       *ci=c->i,*cj=c->j,*cjj;
1081   PetscInt       am=A->m,cn=C->N,cm=C->M;
1082   PetscInt       i,j,k,anzi,pnzi,apnzj,nextap,pnzj,prow,crow;
1083   MatScalar      *ada=ad->a,*aoa=ao->a,*apa,*paj,*ca=c->a,*caj;
1084   MatScalar      *pa_loc=p_loc->a,*pA=pa_loc,*pa_oth=p_oth->a;
1085 
1086   PetscFunctionBegin;
1087     if (!logkey_matptapnumeric_local) {
1088     ierr = PetscLogEventRegister(&logkey_matptapnumeric_local,"MatPtAPNumeric_MPIAIJ_MPIAIJ_Local",MAT_COOKIE);
1089   }
1090   PetscLogEventBegin(logkey_matptapnumeric_local,A,P_loc,P_oth,0);CHKERRQ(ierr);
1091 
1092   /* Allocate temporary array for storage of one row of A*P */
1093   ierr = PetscMalloc(cn*(sizeof(MatScalar)+2*sizeof(PetscInt)),&apa);CHKERRQ(ierr);
1094   ierr = PetscMemzero(apa,cn*(sizeof(MatScalar)+2*sizeof(PetscInt)));CHKERRQ(ierr);
1095   apj      = (PetscInt *)(apa + cn);
1096   apjdense = apj + cn;
1097 
1098   /* Clear old values in C */
1099   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
1100 
1101   for (i=0;i<am;i++) {
1102     /* Form i-th sparse row of A*P */
1103      apnzj = 0;
1104     /* diagonal portion of A */
1105     anzi  = adi[i+1] - adi[i];
1106     for (j=0;j<anzi;j++) {
1107       prow = *adj;
1108       adj++;
1109       pnzj = pi_loc[prow+1] - pi_loc[prow];
1110       pjj  = pj_loc + pi_loc[prow];
1111       paj  = pa_loc + pi_loc[prow];
1112       for (k=0;k<pnzj;k++) {
1113         if (!apjdense[pjj[k]]) {
1114           apjdense[pjj[k]] = -1;
1115           apj[apnzj++]     = pjj[k];
1116         }
1117         apa[pjj[k]] += (*ada)*paj[k];
1118       }
1119       flops += 2*pnzj;
1120       ada++;
1121     }
1122     /* off-diagonal portion of A */
1123     anzi  = aoi[i+1] - aoi[i];
1124     for (j=0;j<anzi;j++) {
1125       col = a->garray[*aoj];
1126       if (col < cstart){
1127         prow = *aoj;
1128       } else if (col >= cend){
1129         prow = *aoj;
1130       } else {
1131         SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Off-diagonal portion of A has wrong column map");
1132       }
1133       aoj++;
1134       pnzj = pi_oth[prow+1] - pi_oth[prow];
1135       pjj  = pj_oth + pi_oth[prow];
1136       paj  = pa_oth + pi_oth[prow];
1137       for (k=0;k<pnzj;k++) {
1138         if (!apjdense[pjj[k]]) {
1139           apjdense[pjj[k]] = -1;
1140           apj[apnzj++]     = pjj[k];
1141         }
1142         apa[pjj[k]] += (*aoa)*paj[k];
1143       }
1144       flops += 2*pnzj;
1145       aoa++;
1146     }
1147     /* Sort the j index array for quick sparse axpy. */
1148     ierr = PetscSortInt(apnzj,apj);CHKERRQ(ierr);
1149 
1150     /* Compute P_loc[i,:]^T*AP[i,:] using outer product */
1151     pnzi = pi_loc[i+1] - pi_loc[i];
1152     for (j=0;j<pnzi;j++) {
1153       nextap = 0;
1154       crow   = *pJ++;
1155       cjj    = cj + ci[crow];
1156       caj    = ca + ci[crow];
1157       /* Perform sparse axpy operation.  Note cjj includes apj. */
1158       for (k=0;nextap<apnzj;k++) {
1159         if (cjj[k]==apj[nextap]) {
1160           caj[k] += (*pA)*apa[apj[nextap++]];
1161         }
1162       }
1163       flops += 2*apnzj;
1164       pA++;
1165     }
1166 
1167     /* Zero the current row info for A*P */
1168     for (j=0;j<apnzj;j++) {
1169       apa[apj[j]]      = 0.;
1170       apjdense[apj[j]] = 0;
1171     }
1172   }
1173 
1174   /* Assemble the final matrix and clean up */
1175   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1176   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1177   ierr = PetscFree(apa);CHKERRQ(ierr);
1178   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
1179   PetscLogEventEnd(logkey_matptapnumeric_local,A,P_loc,P_oth,0);CHKERRQ(ierr);
1180   PetscFunctionReturn(0);
1181 }
1182 static PetscEvent logkey_matptapsymbolic_local = 0;
1183 #undef __FUNCT__
1184 #define __FUNCT__ "MatPtAPSymbolic_MPIAIJ_MPIAIJ_Local"
1185 PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_Local(Mat A,Mat P_loc,Mat P_oth,PetscInt prstart,PetscReal fill,Mat *C)
1186 {
1187   PetscErrorCode ierr;
1188   FreeSpaceList  free_space=PETSC_NULL,current_space=PETSC_NULL;
1189   Mat_MPIAIJ     *a=(Mat_MPIAIJ*)A->data;
1190   Mat_SeqAIJ     *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*c;
1191   PetscInt       *pti,*ptj,*ptJ,*ajj,*pjj;
1192   PetscInt       *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,nnz,cstart=a->cstart,cend=a->cend,col;
1193   PetscInt       *ci,*cj,*ptaj;
1194   PetscInt       an=A->N,am=A->m,pN=P_loc->N;
1195   PetscInt       i,j,k,ptnzi,arow,anzj,prow,pnzj,cnzi;
1196   PetscInt       pm=P_loc->m,nlnk,*lnk;
1197   MatScalar      *ca;
1198   PetscBT        lnkbt;
1199   PetscInt       prend,nprow_loc,nprow_oth;
1200   PetscInt       *ptadenserow_loc,*ptadenserow_oth,*ptasparserow_loc,*ptasparserow_oth;
1201   Mat_SeqAIJ     *p_loc=(Mat_SeqAIJ*)P_loc->data,*p_oth=(Mat_SeqAIJ*)P_oth->data;
1202   PetscInt       *pi_loc=p_loc->i,*pj_loc=p_loc->j,*pi_oth=p_oth->i,*pj_oth=p_oth->j;
1203   /* PetscInt       rank; */
1204 
1205   PetscFunctionBegin;
1206   if (!logkey_matptapsymbolic_local) {
1207     ierr = PetscLogEventRegister(&logkey_matptapsymbolic_local,"MatPtAPSymbolic_MPIAIJ_MPIAIJ_Local",MAT_COOKIE);
1208   }
1209   ierr = PetscLogEventBegin(logkey_matptapsymbolic_local,A,P_loc,P_oth,0);CHKERRQ(ierr);
1210 
1211   /* ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr); */
1212   prend = prstart + pm;
1213 
1214   /* get ij structure of P_loc^T */
1215   ierr = MatGetSymbolicTranspose_SeqAIJ(P_loc,&pti,&ptj);CHKERRQ(ierr);
1216   ptJ=ptj;
1217 
1218   /* Allocate ci array, arrays for fill computation and */
1219   /* free space for accumulating nonzero column info */
1220   ierr = PetscMalloc((pN+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
1221   ci[0] = 0;
1222 
1223   ierr = PetscMalloc((4*an+1)*sizeof(PetscInt),&ptadenserow_loc);CHKERRQ(ierr);
1224   ierr = PetscMemzero(ptadenserow_loc,(4*an+1)*sizeof(PetscInt));CHKERRQ(ierr);
1225   ptasparserow_loc = ptadenserow_loc + an;
1226   ptadenserow_oth  = ptasparserow_loc + an;
1227   ptasparserow_oth = ptadenserow_oth + an;
1228 
1229   /* create and initialize a linked list */
1230   nlnk = pN+1;
1231   ierr = PetscLLCreate(pN,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
1232 
1233   /* Set initial free space to be nnz(A) scaled by aspect ratio of P (P->M=A->N). */
1234   /* This should be reasonable if sparsity of PtAP is similar to that of A. */
1235   nnz           = adi[am] + aoi[am];
1236   ierr          = GetMoreSpace((PetscInt)(fill*nnz/A->N)*pN,&free_space);
1237   current_space = free_space;
1238 
1239   /* determine symbolic info for each row of C: */
1240   for (i=0;i<pN;i++) {
1241     ptnzi  = pti[i+1] - pti[i];
1242     nprow_loc = 0; nprow_oth = 0;
1243     /* i-th row of symbolic P_loc^T*A_loc: */
1244     for (j=0;j<ptnzi;j++) {
1245       arow = *ptJ++;
1246       /* if (rank==1) printf("i: %d, ptnzi: %d, arow: %d\n",i,ptnzi,arow); */
1247       /* diagonal portion of A */
1248       anzj = adi[arow+1] - adi[arow];
1249       ajj  = adj + adi[arow];
1250       for (k=0;k<anzj;k++) {
1251         col = ajj[k]+prstart;
1252         if (!ptadenserow_loc[col]) {
1253           ptadenserow_loc[col]    = -1;
1254           ptasparserow_loc[nprow_loc++] = col;
1255         }
1256       }
1257       /* off-diagonal portion of A */
1258       anzj = aoi[arow+1] - aoi[arow];
1259       ajj  = aoj + aoi[arow];
1260       for (k=0;k<anzj;k++) {
1261         col = a->garray[ajj[k]];  /* global col */
1262         if (col < cstart){
1263           col = ajj[k];
1264         } else if (col >= cend){
1265           col = ajj[k] + pm;
1266         } else {
1267           SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Off-diagonal portion of A has wrong column map");
1268         }
1269         if (!ptadenserow_oth[col]) {
1270           ptadenserow_oth[col]    = -1;
1271           ptasparserow_oth[nprow_oth++] = col;
1272         }
1273       }
1274     }
1275 
1276     /* printf(" [%d] i: %d, nprow_loc: %d, nprow_oth: %d\n",rank,i,nprow_loc,nprow_oth); */
1277     /* using symbolic info of local PtA, determine symbolic info for row of C: */
1278     cnzi   = 0;
1279     /* rows of P_loc */
1280     ptaj = ptasparserow_loc;
1281     for (j=0; j<nprow_loc; j++) {
1282       prow = *ptaj++;
1283       prow -= prstart; /* rm */
1284       pnzj = pi_loc[prow+1] - pi_loc[prow];
1285       pjj  = pj_loc + pi_loc[prow];
1286       /* add non-zero cols of P into the sorted linked list lnk */
1287       ierr = PetscLLAdd(pnzj,pjj,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
1288       cnzi += nlnk;
1289     }
1290     /* rows of P_oth */
1291     ptaj = ptasparserow_oth;
1292     for (j=0; j<nprow_oth; j++) {
1293       prow = *ptaj++;
1294       if (prow >= prend) prow -= pm; /* rm */
1295       pnzj = pi_oth[prow+1] - pi_oth[prow];
1296       pjj  = pj_oth + pi_oth[prow];
1297       /* add non-zero cols of P into the sorted linked list lnk */
1298       ierr = PetscLLAdd(pnzj,pjj,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
1299       cnzi += nlnk;
1300     }
1301 
1302     /* If free space is not available, make more free space */
1303     /* Double the amount of total space in the list */
1304     if (current_space->local_remaining<cnzi) {
1305       ierr = GetMoreSpace(current_space->total_array_size,&current_space);CHKERRQ(ierr);
1306     }
1307 
1308     /* Copy data into free space, and zero out denserows */
1309     ierr = PetscLLClean(pN,pN,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
1310     current_space->array           += cnzi;
1311     current_space->local_used      += cnzi;
1312     current_space->local_remaining -= cnzi;
1313 
1314     for (j=0;j<nprow_loc; j++) {
1315       ptadenserow_loc[ptasparserow_loc[j]] = 0;
1316     }
1317     for (j=0;j<nprow_oth; j++) {
1318       ptadenserow_oth[ptasparserow_oth[j]] = 0;
1319     }
1320 
1321     /* Aside: Perhaps we should save the pta info for the numerical factorization. */
1322     /*        For now, we will recompute what is needed. */
1323     ci[i+1] = ci[i] + cnzi;
1324   }
1325   /* nnz is now stored in ci[ptm], column indices are in the list of free space */
1326   /* Allocate space for cj, initialize cj, and */
1327   /* destroy list of free space and other temporary array(s) */
1328   ierr = PetscMalloc((ci[pN]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
1329   ierr = MakeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
1330   ierr = PetscFree(ptadenserow_loc);CHKERRQ(ierr);
1331   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
1332 
1333   /* Allocate space for ca */
1334   ierr = PetscMalloc((ci[pN]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
1335   ierr = PetscMemzero(ca,(ci[pN]+1)*sizeof(MatScalar));CHKERRQ(ierr);
1336 
1337   /* put together the new matrix */
1338   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,pN,pN,ci,cj,ca,C);CHKERRQ(ierr);
1339 
1340   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
1341   /* Since these are PETSc arrays, change flags to free them as necessary. */
1342   c = (Mat_SeqAIJ *)((*C)->data);
1343   c->freedata = PETSC_TRUE;
1344   c->nonew    = 0;
1345 
1346   /* Clean up. */
1347   ierr = MatRestoreSymbolicTranspose_SeqAIJ(P_loc,&pti,&ptj);CHKERRQ(ierr);
1348   /*
1349   ierr = MPI_Barrier(PETSC_COMM_WORLD);
1350   if (rank == 1){
1351     printf("[%d] C: %d, %d\n",rank,(*C)->M,(*C)->N);
1352     ierr = MatView(*C, PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
1353     } */
1354   ierr = PetscLogEventEnd(logkey_matptapsymbolic_local,A,P_loc,P_oth,0);CHKERRQ(ierr);
1355   PetscFunctionReturn(0);
1356 }
1357