xref: /petsc/src/mat/impls/aij/seq/matptap.c (revision 4c627768fec2a8f2c2ebd97828224a20479036a2)
1eb9c0419SKris Buschelman /*
29af31e4aSHong Zhang   Defines projective product routines where A is a AIJ matrix
3eb9c0419SKris Buschelman           C = P^T * A * P
4eb9c0419SKris Buschelman */
5eb9c0419SKris Buschelman 
6231952e2SKris Buschelman #include "src/mat/impls/aij/seq/aij.h"   /*I "petscmat.h" I*/
7eb9c0419SKris Buschelman #include "src/mat/utils/freespace.h"
89af31e4aSHong Zhang #include "src/mat/impls/aij/mpi/mpiaij.h"
9eb9c0419SKris Buschelman 
10eb9c0419SKris Buschelman #undef __FUNCT__
119af31e4aSHong Zhang #define __FUNCT__ "MatPtAP"
124d3841fdSKris Buschelman /*@
139af31e4aSHong Zhang    MatPtAP - Creates the matrix projection C = P^T * A * P
144d3841fdSKris Buschelman 
154d3841fdSKris Buschelman    Collective on Mat
164d3841fdSKris Buschelman 
174d3841fdSKris Buschelman    Input Parameters:
184d3841fdSKris Buschelman +  A - the matrix
19f747e1aeSHong Zhang .  P - the projection matrix
20f747e1aeSHong Zhang .  scall - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX
21f747e1aeSHong Zhang -  fill - expected fill as ratio of nnz(C)/(nnz(A) + nnz(P))
224d3841fdSKris Buschelman 
234d3841fdSKris Buschelman    Output Parameters:
244d3841fdSKris Buschelman .  C - the product matrix
254d3841fdSKris Buschelman 
264d3841fdSKris Buschelman    Notes:
274d3841fdSKris Buschelman    C will be created and must be destroyed by the user with MatDestroy().
284d3841fdSKris Buschelman 
294d3841fdSKris Buschelman    This routine is currently only implemented for pairs of SeqAIJ matrices and classes
304d3841fdSKris Buschelman    which inherit from SeqAIJ.  C will be of type MATSEQAIJ.
314d3841fdSKris Buschelman 
324d3841fdSKris Buschelman    Level: intermediate
334d3841fdSKris Buschelman 
349af31e4aSHong Zhang .seealso: MatPtAPSymbolic(),MatPtAPNumeric(),MatMatMult()
354d3841fdSKris Buschelman @*/
36dfbe8321SBarry Smith PetscErrorCode MatPtAP(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C) {
37dfbe8321SBarry Smith   PetscErrorCode ierr;
38534c1384SKris Buschelman   PetscErrorCode (*fA)(Mat,Mat,MatReuse,PetscReal,Mat *);
39534c1384SKris Buschelman   PetscErrorCode (*fP)(Mat,Mat,MatReuse,PetscReal,Mat *);
40eb9c0419SKris Buschelman 
41eb9c0419SKris Buschelman   PetscFunctionBegin;
429af31e4aSHong Zhang   PetscValidHeaderSpecific(A,MAT_COOKIE,1);
439af31e4aSHong Zhang   PetscValidType(A,1);
449af31e4aSHong Zhang   MatPreallocated(A);
459af31e4aSHong Zhang   if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
469af31e4aSHong Zhang   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
479af31e4aSHong Zhang   PetscValidHeaderSpecific(P,MAT_COOKIE,2);
489af31e4aSHong Zhang   PetscValidType(P,2);
499af31e4aSHong Zhang   MatPreallocated(P);
509af31e4aSHong Zhang   if (!P->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
519af31e4aSHong Zhang   if (P->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
529af31e4aSHong Zhang   PetscValidPointer(C,3);
53534c1384SKris Buschelman 
549af31e4aSHong Zhang   if (P->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->M,A->N);
55eb9c0419SKris Buschelman 
569af31e4aSHong Zhang   if (fill <=0.0) SETERRQ1(PETSC_ERR_ARG_SIZ,"fill=%g must be > 0.0",fill);
57eb9c0419SKris Buschelman 
58534c1384SKris Buschelman   /* For now, we do not dispatch based on the type of A and P */
59534c1384SKris Buschelman   /* When implementations like _SeqAIJ_MAIJ exist, attack the multiple dispatch problem. */
60534c1384SKris Buschelman   fA = A->ops->ptap;
61534c1384SKris Buschelman   if (!fA) SETERRQ1(PETSC_ERR_SUP,"MatPtAP not supported for A of type %s",A->type_name);
62534c1384SKris Buschelman   fP = P->ops->ptap;
63534c1384SKris Buschelman   if (!fP) SETERRQ1(PETSC_ERR_SUP,"MatPtAP not supported for P of type %s",P->type_name);
64534c1384SKris Buschelman   if (fP!=fA) SETERRQ2(PETSC_ERR_ARG_INCOMP,"MatPtAP requires A, %s, to be compatible with P, %s",A->type_name,P->type_name);
65534c1384SKris Buschelman 
669af31e4aSHong Zhang   ierr = PetscLogEventBegin(MAT_PtAP,A,P,0,0);CHKERRQ(ierr);
67534c1384SKris Buschelman   ierr = (*fA)(A,P,scall,fill,C);CHKERRQ(ierr);
689af31e4aSHong Zhang   ierr = PetscLogEventEnd(MAT_PtAP,A,P,0,0);CHKERRQ(ierr);
69eb9c0419SKris Buschelman   PetscFunctionReturn(0);
70eb9c0419SKris Buschelman }
71eb9c0419SKris Buschelman 
720e36024fSHong Zhang EXTERN PetscErrorCode MatPtAP_SeqAIJ_SeqAIJ_ReducedPt(Mat,Mat,MatReuse,PetscReal,int,int,Mat*);
730e36024fSHong Zhang EXTERN PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_ReducedPt(Mat,Mat,PetscReal,int,int,Mat*);
740e36024fSHong Zhang EXTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_ReducedPt(Mat,Mat,int,int,Mat);
75b90dcfe3SHong Zhang 
76eb9c0419SKris Buschelman #undef __FUNCT__
77ff134f7aSHong Zhang #define __FUNCT__ "MatPtAP_MPIAIJ_MPIAIJ"
78ff134f7aSHong Zhang PetscErrorCode MatPtAP_MPIAIJ_MPIAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
79ff134f7aSHong Zhang {
80ff134f7aSHong Zhang   PetscErrorCode    ierr;
81b90dcfe3SHong Zhang 
82b90dcfe3SHong Zhang   PetscFunctionBegin;
83b90dcfe3SHong Zhang   if (scall == MAT_INITIAL_MATRIX){
84*4c627768SHong Zhang     ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
85b90dcfe3SHong Zhang     ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ(A,P,fill,C);CHKERRQ(ierr);/* numeric product is computed as well */
86*4c627768SHong Zhang     ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
87b90dcfe3SHong Zhang   } else if (scall == MAT_REUSE_MATRIX){
88*4c627768SHong Zhang     ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
89b90dcfe3SHong Zhang     ierr = MatPtAPNumeric_MPIAIJ_MPIAIJ(A,P,*C);CHKERRQ(ierr);
90*4c627768SHong Zhang     ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
91b90dcfe3SHong Zhang   } else {
92b90dcfe3SHong Zhang     SETERRQ1(PETSC_ERR_ARG_WRONG,"Invalid MatReuse %d",scall);
93b90dcfe3SHong Zhang   }
94b90dcfe3SHong Zhang   PetscFunctionReturn(0);
95b90dcfe3SHong Zhang }
96b90dcfe3SHong Zhang 
97b90dcfe3SHong Zhang #undef __FUNCT__
98b90dcfe3SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ"
99b90dcfe3SHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
100b90dcfe3SHong Zhang {
101b90dcfe3SHong Zhang   PetscErrorCode    ierr;
10225616d81SHong Zhang   Mat               P_seq,A_loc,C_seq;
1030e36024fSHong Zhang   int               prstart,prend,m=P->m;
10425616d81SHong Zhang   IS                isrowp,iscolp;
105ff134f7aSHong Zhang 
106ff134f7aSHong Zhang   PetscFunctionBegin;
10725616d81SHong Zhang   /* get P_seq = submatrix of P by taking rows of P that equal to nonzero col of A */
108b90dcfe3SHong Zhang   ierr = MatGetBrowsOfAcols(A,P,MAT_INITIAL_MATRIX,&isrowp,&iscolp,&prstart,&P_seq);CHKERRQ(ierr);
10925616d81SHong Zhang   ierr = ISDestroy(iscolp);CHKERRQ(ierr);
110ff134f7aSHong Zhang 
11125616d81SHong Zhang   /* get A_loc = submatrix of A by taking all local rows of A */
112b90dcfe3SHong Zhang   ierr = MatGetLocalMat(A,MAT_INITIAL_MATRIX,PETSC_NULL,&isrowp,&A_loc);CHKERRQ(ierr);
11325616d81SHong Zhang   ierr = ISDestroy(isrowp);CHKERRQ(ierr);
1140e36024fSHong Zhang 
11525616d81SHong Zhang   /* compute C_seq = P_loc^T * A_loc * P_seq */
116ff134f7aSHong Zhang   prend   = prstart + m;
117b90dcfe3SHong Zhang   ierr = MatPtAP_SeqAIJ_SeqAIJ_ReducedPt(A_loc,P_seq,MAT_INITIAL_MATRIX,fill,prstart,prend,&C_seq);CHKERRQ(ierr);
11825616d81SHong Zhang   ierr = MatDestroy(P_seq);CHKERRQ(ierr);
11925616d81SHong Zhang   ierr = MatDestroy(A_loc);CHKERRQ(ierr);
120b90dcfe3SHong Zhang 
121b90dcfe3SHong Zhang   /* add C_seq into mpi C */
122b90dcfe3SHong Zhang   ierr = MatMerge_SeqsToMPI(A->comm,C_seq,P->n,P->n,MAT_INITIAL_MATRIX,C);CHKERRQ(ierr);
123b90dcfe3SHong Zhang 
124ff134f7aSHong Zhang   PetscFunctionReturn(0);
125ff134f7aSHong Zhang }
126ff134f7aSHong Zhang 
127ff134f7aSHong Zhang #undef __FUNCT__
128ff134f7aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ"
129b90dcfe3SHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
130ff134f7aSHong Zhang {
131b90dcfe3SHong Zhang   PetscErrorCode       ierr;
132b90dcfe3SHong Zhang   Mat                  P_seq,A_loc,C_seq;
133b90dcfe3SHong Zhang   int                  prstart,prend,m=P->m;
134b90dcfe3SHong Zhang   IS                   isrowp,iscolp;
135671beff6SHong Zhang   Mat_Merge_SeqsToMPI  *merge;
136671beff6SHong Zhang   PetscObjectContainer container;
137ff134f7aSHong Zhang 
138ff134f7aSHong Zhang   PetscFunctionBegin;
139671beff6SHong Zhang   ierr = PetscObjectQuery((PetscObject)C,"MatMergeSeqsToMPI",(PetscObject *)&container);CHKERRQ(ierr);
140671beff6SHong Zhang   if (container) {
1417f79fd58SMatthew Knepley     ierr  = PetscObjectContainerGetPointer(container,(void **)&merge);CHKERRQ(ierr);
1427f79fd58SMatthew Knepley   } else {
1437f79fd58SMatthew Knepley     SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Matrix C does not posses and object container");
144671beff6SHong Zhang   }
145671beff6SHong Zhang 
146b90dcfe3SHong Zhang   /* get P_seq = submatrix of P by taking rows of P that equal to nonzero col of A */
147b90dcfe3SHong Zhang   ierr = MatGetBrowsOfAcols(A,P,MAT_INITIAL_MATRIX,&isrowp,&iscolp,&prstart,&P_seq);CHKERRQ(ierr);
148b90dcfe3SHong Zhang   ierr = ISDestroy(iscolp);CHKERRQ(ierr);
149ff134f7aSHong Zhang 
150b90dcfe3SHong Zhang   /* get A_loc = submatrix of A by taking all local rows of A */
151b90dcfe3SHong Zhang   ierr = MatGetLocalMat(A,MAT_INITIAL_MATRIX,PETSC_NULL,&isrowp,&A_loc);CHKERRQ(ierr);
152b90dcfe3SHong Zhang   ierr = ISDestroy(isrowp);CHKERRQ(ierr);
153ff134f7aSHong Zhang 
154b90dcfe3SHong Zhang   /* compute C_seq = P_loc^T * A_loc * P_seq */
155b90dcfe3SHong Zhang   prend = prstart + m;
156b90dcfe3SHong Zhang   C_seq = merge->C_seq;
157b90dcfe3SHong Zhang   ierr = MatPtAP_SeqAIJ_SeqAIJ_ReducedPt(A_loc,P_seq,MAT_REUSE_MATRIX,1.0,prstart,prend,&C_seq);CHKERRQ(ierr);
158b90dcfe3SHong Zhang   ierr = MatDestroy(P_seq);CHKERRQ(ierr);
159b90dcfe3SHong Zhang   ierr = MatDestroy(A_loc);CHKERRQ(ierr);
160b90dcfe3SHong Zhang 
161b90dcfe3SHong Zhang   /* add C_seq into mpi C */
162b90dcfe3SHong Zhang   ierr = MatMerge_SeqsToMPI(A->comm,C_seq,P->n,P->n,MAT_REUSE_MATRIX,&C);CHKERRQ(ierr);
163b90dcfe3SHong Zhang 
164ff134f7aSHong Zhang   PetscFunctionReturn(0);
165ff134f7aSHong Zhang }
166ff134f7aSHong Zhang 
167ff134f7aSHong Zhang #undef __FUNCT__
1689af31e4aSHong Zhang #define __FUNCT__ "MatPtAP_SeqAIJ_SeqAIJ"
169dfbe8321SBarry Smith PetscErrorCode MatPtAP_SeqAIJ_SeqAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
1709af31e4aSHong Zhang {
171dfbe8321SBarry Smith   PetscErrorCode ierr;
1729af31e4aSHong Zhang   PetscFunctionBegin;
1739af31e4aSHong Zhang   if (scall == MAT_INITIAL_MATRIX){
174d20bfe6fSHong Zhang     ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
1759af31e4aSHong Zhang     ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ(A,P,fill,C);CHKERRQ(ierr);
176d20bfe6fSHong Zhang     ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
1779af31e4aSHong Zhang   }
178d20bfe6fSHong Zhang   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
1799af31e4aSHong Zhang   ierr = MatPtAPNumeric_SeqAIJ_SeqAIJ(A,P,*C);CHKERRQ(ierr);
180d20bfe6fSHong Zhang   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
1819af31e4aSHong Zhang   PetscFunctionReturn(0);
1829af31e4aSHong Zhang }
1839af31e4aSHong Zhang 
1849af31e4aSHong Zhang #undef __FUNCT__
1859af31e4aSHong Zhang #define __FUNCT__ "MatPtAPSymbolic"
1866849ba73SBarry Smith /*
1879af31e4aSHong Zhang    MatPtAPSymbolic - Creates the (i,j) structure of the matrix projection C = P^T * A * P
1884d3841fdSKris Buschelman 
1894d3841fdSKris Buschelman    Collective on Mat
1904d3841fdSKris Buschelman 
1914d3841fdSKris Buschelman    Input Parameters:
1924d3841fdSKris Buschelman +  A - the matrix
1934d3841fdSKris Buschelman -  P - the projection matrix
1944d3841fdSKris Buschelman 
1954d3841fdSKris Buschelman    Output Parameters:
1964d3841fdSKris Buschelman .  C - the (i,j) structure of the product matrix
1974d3841fdSKris Buschelman 
1984d3841fdSKris Buschelman    Notes:
1994d3841fdSKris Buschelman    C will be created and must be destroyed by the user with MatDestroy().
2004d3841fdSKris Buschelman 
2014d3841fdSKris Buschelman    This routine is currently only implemented for pairs of SeqAIJ matrices and classes
2024d3841fdSKris Buschelman    which inherit from SeqAIJ.  C will be of type MATSEQAIJ.  The product is computed using
2039af31e4aSHong Zhang    this (i,j) structure by calling MatPtAPNumeric().
2044d3841fdSKris Buschelman 
2054d3841fdSKris Buschelman    Level: intermediate
2064d3841fdSKris Buschelman 
2079af31e4aSHong Zhang .seealso: MatPtAP(),MatPtAPNumeric(),MatMatMultSymbolic()
2086849ba73SBarry Smith */
209dfbe8321SBarry Smith PetscErrorCode MatPtAPSymbolic(Mat A,Mat P,PetscReal fill,Mat *C) {
210dfbe8321SBarry Smith   PetscErrorCode ierr;
211534c1384SKris Buschelman   PetscErrorCode (*fA)(Mat,Mat,PetscReal,Mat*);
212534c1384SKris Buschelman   PetscErrorCode (*fP)(Mat,Mat,PetscReal,Mat*);
213eb9c0419SKris Buschelman 
214eb9c0419SKris Buschelman   PetscFunctionBegin;
215eb9c0419SKris Buschelman 
2164482741eSBarry Smith   PetscValidHeaderSpecific(A,MAT_COOKIE,1);
217c9780b6fSBarry Smith   PetscValidType(A,1);
218eb9c0419SKris Buschelman   MatPreallocated(A);
219eb9c0419SKris Buschelman   if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
220eb9c0419SKris Buschelman   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
221eb9c0419SKris Buschelman 
2224482741eSBarry Smith   PetscValidHeaderSpecific(P,MAT_COOKIE,2);
223c9780b6fSBarry Smith   PetscValidType(P,2);
224eb9c0419SKris Buschelman   MatPreallocated(P);
225eb9c0419SKris Buschelman   if (!P->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
226eb9c0419SKris Buschelman   if (P->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
227eb9c0419SKris Buschelman 
2284482741eSBarry Smith   PetscValidPointer(C,3);
2294482741eSBarry Smith 
230eb9c0419SKris Buschelman   if (P->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->M,A->N);
231eb9c0419SKris Buschelman   if (A->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix 'A' must be square, %d != %d",A->M,A->N);
232eb9c0419SKris Buschelman 
233534c1384SKris Buschelman   /* For now, we do not dispatch based on the type of A and P */
234534c1384SKris Buschelman   /* When implementations like _SeqAIJ_MAIJ exist, attack the multiple dispatch problem. */
235534c1384SKris Buschelman   fA = A->ops->ptapsymbolic;
236534c1384SKris Buschelman   if (!fA) SETERRQ1(PETSC_ERR_SUP,"MatPtAPSymbolic not supported for A of type %s",A->type_name);
237534c1384SKris Buschelman   fP = P->ops->ptapsymbolic;
238534c1384SKris Buschelman   if (!fP) SETERRQ1(PETSC_ERR_SUP,"MatPtAPSymbolic not supported for P of type %s",P->type_name);
239534c1384SKris Buschelman   if (fP!=fA) SETERRQ2(PETSC_ERR_ARG_INCOMP,"MatPtAPSymbolic requires A, %s, to be compatible with P, %s",A->type_name,P->type_name);
2404d3841fdSKris Buschelman 
241534c1384SKris Buschelman   ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
242534c1384SKris Buschelman   ierr = (*fA)(A,P,fill,C);CHKERRQ(ierr);
243534c1384SKris Buschelman   ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
244eb9c0419SKris Buschelman 
245eb9c0419SKris Buschelman   PetscFunctionReturn(0);
246eb9c0419SKris Buschelman }
247eb9c0419SKris Buschelman 
248f747e1aeSHong Zhang typedef struct {
249f747e1aeSHong Zhang   Mat    symAP;
250f747e1aeSHong Zhang } Mat_PtAPstruct;
251f747e1aeSHong Zhang 
25278a80504SBarry Smith EXTERN PetscErrorCode MatDestroy_SeqAIJ(Mat);
25378a80504SBarry Smith 
254f747e1aeSHong Zhang #undef __FUNCT__
255f747e1aeSHong Zhang #define __FUNCT__ "MatDestroy_SeqAIJ_PtAP"
256f4a850bbSBarry Smith PetscErrorCode MatDestroy_SeqAIJ_PtAP(Mat A)
257f747e1aeSHong Zhang {
258f4a850bbSBarry Smith   PetscErrorCode    ierr;
259f747e1aeSHong Zhang   Mat_PtAPstruct    *ptap=(Mat_PtAPstruct*)A->spptr;
260f747e1aeSHong Zhang 
261f747e1aeSHong Zhang   PetscFunctionBegin;
262f747e1aeSHong Zhang   ierr = MatDestroy(ptap->symAP);CHKERRQ(ierr);
263f747e1aeSHong Zhang   ierr = PetscFree(ptap);CHKERRQ(ierr);
26478a80504SBarry Smith   ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr);
265f747e1aeSHong Zhang   PetscFunctionReturn(0);
266f747e1aeSHong Zhang }
267f747e1aeSHong Zhang 
268eb9c0419SKris Buschelman #undef __FUNCT__
2699af31e4aSHong Zhang #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqAIJ"
270dfbe8321SBarry Smith PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat P,PetscReal fill,Mat *C) {
271dfbe8321SBarry Smith   PetscErrorCode ierr;
272d20bfe6fSHong Zhang   FreeSpaceList  free_space=PETSC_NULL,current_space=PETSC_NULL;
273d20bfe6fSHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*p=(Mat_SeqAIJ*)P->data,*c;
274d20bfe6fSHong Zhang   int            *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj;
275d20bfe6fSHong Zhang   int            *ci,*cj,*denserow,*sparserow,*ptadenserow,*ptasparserow,*ptaj;
276d20bfe6fSHong Zhang   int            an=A->N,am=A->M,pn=P->N,pm=P->M;
277d20bfe6fSHong Zhang   int            i,j,k,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi;
278d20bfe6fSHong Zhang   MatScalar      *ca;
279eb9c0419SKris Buschelman 
280eb9c0419SKris Buschelman   PetscFunctionBegin;
281d20bfe6fSHong Zhang   /* Get ij structure of P^T */
282eb9c0419SKris Buschelman   ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
283d20bfe6fSHong Zhang   ptJ=ptj;
284eb9c0419SKris Buschelman 
285d20bfe6fSHong Zhang   /* Allocate ci array, arrays for fill computation and */
286d20bfe6fSHong Zhang   /* free space for accumulating nonzero column info */
287d20bfe6fSHong Zhang   ierr = PetscMalloc((pn+1)*sizeof(int),&ci);CHKERRQ(ierr);
288d20bfe6fSHong Zhang   ci[0] = 0;
289eb9c0419SKris Buschelman 
290d20bfe6fSHong Zhang   ierr = PetscMalloc((2*pn+2*an+1)*sizeof(int),&ptadenserow);CHKERRQ(ierr);
291d20bfe6fSHong Zhang   ierr = PetscMemzero(ptadenserow,(2*pn+2*an+1)*sizeof(int));CHKERRQ(ierr);
292d20bfe6fSHong Zhang   ptasparserow = ptadenserow  + an;
293d20bfe6fSHong Zhang   denserow     = ptasparserow + an;
294d20bfe6fSHong Zhang   sparserow    = denserow     + pn;
295eb9c0419SKris Buschelman 
296d20bfe6fSHong Zhang   /* Set initial free space to be nnz(A) scaled by aspect ratio of P. */
297d20bfe6fSHong Zhang   /* This should be reasonable if sparsity of PtAP is similar to that of A. */
298d20bfe6fSHong Zhang   ierr          = GetMoreSpace((ai[am]/pm)*pn,&free_space);
299d20bfe6fSHong Zhang   current_space = free_space;
300d20bfe6fSHong Zhang 
301d20bfe6fSHong Zhang   /* Determine symbolic info for each row of C: */
302d20bfe6fSHong Zhang   for (i=0;i<pn;i++) {
303d20bfe6fSHong Zhang     ptnzi  = pti[i+1] - pti[i];
304d20bfe6fSHong Zhang     ptanzi = 0;
305d20bfe6fSHong Zhang     /* Determine symbolic row of PtA: */
306d20bfe6fSHong Zhang     for (j=0;j<ptnzi;j++) {
307d20bfe6fSHong Zhang       arow = *ptJ++;
308d20bfe6fSHong Zhang       anzj = ai[arow+1] - ai[arow];
309d20bfe6fSHong Zhang       ajj  = aj + ai[arow];
310d20bfe6fSHong Zhang       for (k=0;k<anzj;k++) {
311d20bfe6fSHong Zhang         if (!ptadenserow[ajj[k]]) {
312d20bfe6fSHong Zhang           ptadenserow[ajj[k]]    = -1;
313d20bfe6fSHong Zhang           ptasparserow[ptanzi++] = ajj[k];
314d20bfe6fSHong Zhang         }
315d20bfe6fSHong Zhang       }
316d20bfe6fSHong Zhang     }
317d20bfe6fSHong Zhang       /* Using symbolic info for row of PtA, determine symbolic info for row of C: */
318d20bfe6fSHong Zhang     ptaj = ptasparserow;
319d20bfe6fSHong Zhang     cnzi   = 0;
320d20bfe6fSHong Zhang     for (j=0;j<ptanzi;j++) {
321d20bfe6fSHong Zhang       prow = *ptaj++;
322d20bfe6fSHong Zhang       pnzj = pi[prow+1] - pi[prow];
323d20bfe6fSHong Zhang       pjj  = pj + pi[prow];
324d20bfe6fSHong Zhang       for (k=0;k<pnzj;k++) {
325d20bfe6fSHong Zhang         if (!denserow[pjj[k]]) {
326d20bfe6fSHong Zhang             denserow[pjj[k]]  = -1;
327d20bfe6fSHong Zhang             sparserow[cnzi++] = pjj[k];
328d20bfe6fSHong Zhang         }
329d20bfe6fSHong Zhang       }
330d20bfe6fSHong Zhang     }
331d20bfe6fSHong Zhang 
332d20bfe6fSHong Zhang     /* sort sparserow */
333d20bfe6fSHong Zhang     ierr = PetscSortInt(cnzi,sparserow);CHKERRQ(ierr);
334d20bfe6fSHong Zhang 
335d20bfe6fSHong Zhang     /* If free space is not available, make more free space */
336d20bfe6fSHong Zhang     /* Double the amount of total space in the list */
337d20bfe6fSHong Zhang     if (current_space->local_remaining<cnzi) {
338d20bfe6fSHong Zhang       ierr = GetMoreSpace(current_space->total_array_size,&current_space);CHKERRQ(ierr);
339d20bfe6fSHong Zhang     }
340d20bfe6fSHong Zhang 
341d20bfe6fSHong Zhang     /* Copy data into free space, and zero out denserows */
342d20bfe6fSHong Zhang     ierr = PetscMemcpy(current_space->array,sparserow,cnzi*sizeof(int));CHKERRQ(ierr);
343d20bfe6fSHong Zhang     current_space->array           += cnzi;
344d20bfe6fSHong Zhang     current_space->local_used      += cnzi;
345d20bfe6fSHong Zhang     current_space->local_remaining -= cnzi;
346d20bfe6fSHong Zhang 
347d20bfe6fSHong Zhang     for (j=0;j<ptanzi;j++) {
348d20bfe6fSHong Zhang       ptadenserow[ptasparserow[j]] = 0;
349d20bfe6fSHong Zhang     }
350d20bfe6fSHong Zhang     for (j=0;j<cnzi;j++) {
351d20bfe6fSHong Zhang       denserow[sparserow[j]] = 0;
352d20bfe6fSHong Zhang     }
353d20bfe6fSHong Zhang       /* Aside: Perhaps we should save the pta info for the numerical factorization. */
354d20bfe6fSHong Zhang       /*        For now, we will recompute what is needed. */
355d20bfe6fSHong Zhang     ci[i+1] = ci[i] + cnzi;
356d20bfe6fSHong Zhang   }
357d20bfe6fSHong Zhang   /* nnz is now stored in ci[ptm], column indices are in the list of free space */
358d20bfe6fSHong Zhang   /* Allocate space for cj, initialize cj, and */
359d20bfe6fSHong Zhang   /* destroy list of free space and other temporary array(s) */
360d20bfe6fSHong Zhang   ierr = PetscMalloc((ci[pn]+1)*sizeof(int),&cj);CHKERRQ(ierr);
361d20bfe6fSHong Zhang   ierr = MakeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
362d20bfe6fSHong Zhang   ierr = PetscFree(ptadenserow);CHKERRQ(ierr);
363d20bfe6fSHong Zhang 
364d20bfe6fSHong Zhang   /* Allocate space for ca */
365d20bfe6fSHong Zhang   ierr = PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
366d20bfe6fSHong Zhang   ierr = PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));CHKERRQ(ierr);
367d20bfe6fSHong Zhang 
368d20bfe6fSHong Zhang   /* put together the new matrix */
369d20bfe6fSHong Zhang   ierr = MatCreateSeqAIJWithArrays(A->comm,pn,pn,ci,cj,ca,C);CHKERRQ(ierr);
370d20bfe6fSHong Zhang 
371d20bfe6fSHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
372d20bfe6fSHong Zhang   /* Since these are PETSc arrays, change flags to free them as necessary. */
373d20bfe6fSHong Zhang   c = (Mat_SeqAIJ *)((*C)->data);
374d20bfe6fSHong Zhang   c->freedata = PETSC_TRUE;
375d20bfe6fSHong Zhang   c->nonew    = 0;
376d20bfe6fSHong Zhang 
377d20bfe6fSHong Zhang   /* Clean up. */
378d20bfe6fSHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
379eb9c0419SKris Buschelman 
380eb9c0419SKris Buschelman   PetscFunctionReturn(0);
381eb9c0419SKris Buschelman }
382eb9c0419SKris Buschelman 
3833985e5eaSKris Buschelman #include "src/mat/impls/maij/maij.h"
3843985e5eaSKris Buschelman EXTERN_C_BEGIN
3853985e5eaSKris Buschelman #undef __FUNCT__
3869af31e4aSHong Zhang #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqMAIJ"
387dfbe8321SBarry Smith PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqMAIJ(Mat A,Mat PP,Mat *C) {
3885c66b693SKris Buschelman   /* This routine requires testing -- I don't think it works. */
389dfbe8321SBarry Smith   PetscErrorCode ierr;
3903985e5eaSKris Buschelman   FreeSpaceList  free_space=PETSC_NULL,current_space=PETSC_NULL;
3913985e5eaSKris Buschelman   Mat_SeqMAIJ    *pp=(Mat_SeqMAIJ*)PP->data;
3923985e5eaSKris Buschelman   Mat            P=pp->AIJ;
3933985e5eaSKris Buschelman   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*p=(Mat_SeqAIJ*)P->data,*c;
3943985e5eaSKris Buschelman   int            *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj;
3953985e5eaSKris Buschelman   int            *ci,*cj,*denserow,*sparserow,*ptadenserow,*ptasparserow,*ptaj;
3963985e5eaSKris Buschelman   int            an=A->N,am=A->M,pn=P->N,pm=P->M,ppdof=pp->dof;
397fe05a634SKris Buschelman   int            i,j,k,dof,pdof,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi;
3983985e5eaSKris Buschelman   MatScalar      *ca;
3993985e5eaSKris Buschelman 
4003985e5eaSKris Buschelman   PetscFunctionBegin;
4013985e5eaSKris Buschelman   /* Start timer */
4029af31e4aSHong Zhang   ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,PP,0,0);CHKERRQ(ierr);
4033985e5eaSKris Buschelman 
4043985e5eaSKris Buschelman   /* Get ij structure of P^T */
4053985e5eaSKris Buschelman   ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
4063985e5eaSKris Buschelman 
4073985e5eaSKris Buschelman   /* Allocate ci array, arrays for fill computation and */
4083985e5eaSKris Buschelman   /* free space for accumulating nonzero column info */
4093985e5eaSKris Buschelman   ierr = PetscMalloc((pn+1)*sizeof(int),&ci);CHKERRQ(ierr);
4103985e5eaSKris Buschelman   ci[0] = 0;
4113985e5eaSKris Buschelman 
4123985e5eaSKris Buschelman   ierr = PetscMalloc((2*pn+2*an+1)*sizeof(int),&ptadenserow);CHKERRQ(ierr);
4133985e5eaSKris Buschelman   ierr = PetscMemzero(ptadenserow,(2*pn+2*an+1)*sizeof(int));CHKERRQ(ierr);
4143985e5eaSKris Buschelman   ptasparserow = ptadenserow  + an;
4153985e5eaSKris Buschelman   denserow     = ptasparserow + an;
4163985e5eaSKris Buschelman   sparserow    = denserow     + pn;
4173985e5eaSKris Buschelman 
4183985e5eaSKris Buschelman   /* Set initial free space to be nnz(A) scaled by aspect ratio of P. */
4193985e5eaSKris Buschelman   /* This should be reasonable if sparsity of PtAP is similar to that of A. */
4203985e5eaSKris Buschelman   ierr          = GetMoreSpace((ai[am]/pm)*pn,&free_space);
4213985e5eaSKris Buschelman   current_space = free_space;
4223985e5eaSKris Buschelman 
4233985e5eaSKris Buschelman   /* Determine symbolic info for each row of C: */
4243985e5eaSKris Buschelman   for (i=0;i<pn/ppdof;i++) {
4253985e5eaSKris Buschelman     ptnzi  = pti[i+1] - pti[i];
4263985e5eaSKris Buschelman     ptanzi = 0;
4273985e5eaSKris Buschelman     ptJ    = ptj + pti[i];
4283985e5eaSKris Buschelman     for (dof=0;dof<ppdof;dof++) {
4293985e5eaSKris Buschelman     /* Determine symbolic row of PtA: */
4303985e5eaSKris Buschelman       for (j=0;j<ptnzi;j++) {
4313985e5eaSKris Buschelman         arow = ptJ[j] + dof;
4323985e5eaSKris Buschelman         anzj = ai[arow+1] - ai[arow];
4333985e5eaSKris Buschelman         ajj  = aj + ai[arow];
4343985e5eaSKris Buschelman         for (k=0;k<anzj;k++) {
4353985e5eaSKris Buschelman           if (!ptadenserow[ajj[k]]) {
4363985e5eaSKris Buschelman             ptadenserow[ajj[k]]    = -1;
4373985e5eaSKris Buschelman             ptasparserow[ptanzi++] = ajj[k];
4383985e5eaSKris Buschelman           }
4393985e5eaSKris Buschelman         }
4403985e5eaSKris Buschelman       }
4413985e5eaSKris Buschelman       /* Using symbolic info for row of PtA, determine symbolic info for row of C: */
4423985e5eaSKris Buschelman       ptaj = ptasparserow;
4433985e5eaSKris Buschelman       cnzi   = 0;
4443985e5eaSKris Buschelman       for (j=0;j<ptanzi;j++) {
445fe05a634SKris Buschelman         pdof = *ptaj%dof;
4463985e5eaSKris Buschelman         prow = (*ptaj++)/dof;
4473985e5eaSKris Buschelman         pnzj = pi[prow+1] - pi[prow];
4483985e5eaSKris Buschelman         pjj  = pj + pi[prow];
4493985e5eaSKris Buschelman         for (k=0;k<pnzj;k++) {
450fe05a634SKris Buschelman           if (!denserow[pjj[k]+pdof]) {
451fe05a634SKris Buschelman             denserow[pjj[k]+pdof] = -1;
452fe05a634SKris Buschelman             sparserow[cnzi++]     = pjj[k]+pdof;
4533985e5eaSKris Buschelman           }
4543985e5eaSKris Buschelman         }
4553985e5eaSKris Buschelman       }
4563985e5eaSKris Buschelman 
4573985e5eaSKris Buschelman       /* sort sparserow */
4583985e5eaSKris Buschelman       ierr = PetscSortInt(cnzi,sparserow);CHKERRQ(ierr);
4593985e5eaSKris Buschelman 
4603985e5eaSKris Buschelman       /* If free space is not available, make more free space */
4613985e5eaSKris Buschelman       /* Double the amount of total space in the list */
4623985e5eaSKris Buschelman       if (current_space->local_remaining<cnzi) {
4633985e5eaSKris Buschelman         ierr = GetMoreSpace(current_space->total_array_size,&current_space);CHKERRQ(ierr);
4643985e5eaSKris Buschelman       }
4653985e5eaSKris Buschelman 
4663985e5eaSKris Buschelman       /* Copy data into free space, and zero out denserows */
4673985e5eaSKris Buschelman       ierr = PetscMemcpy(current_space->array,sparserow,cnzi*sizeof(int));CHKERRQ(ierr);
4683985e5eaSKris Buschelman       current_space->array           += cnzi;
4693985e5eaSKris Buschelman       current_space->local_used      += cnzi;
4703985e5eaSKris Buschelman       current_space->local_remaining -= cnzi;
4713985e5eaSKris Buschelman 
4723985e5eaSKris Buschelman       for (j=0;j<ptanzi;j++) {
4733985e5eaSKris Buschelman         ptadenserow[ptasparserow[j]] = 0;
4743985e5eaSKris Buschelman       }
4753985e5eaSKris Buschelman       for (j=0;j<cnzi;j++) {
4763985e5eaSKris Buschelman         denserow[sparserow[j]] = 0;
4773985e5eaSKris Buschelman       }
4783985e5eaSKris Buschelman       /* Aside: Perhaps we should save the pta info for the numerical factorization. */
4793985e5eaSKris Buschelman       /*        For now, we will recompute what is needed. */
4803985e5eaSKris Buschelman       ci[i+1+dof] = ci[i+dof] + cnzi;
4813985e5eaSKris Buschelman     }
4823985e5eaSKris Buschelman   }
4833985e5eaSKris Buschelman   /* nnz is now stored in ci[ptm], column indices are in the list of free space */
4843985e5eaSKris Buschelman   /* Allocate space for cj, initialize cj, and */
4853985e5eaSKris Buschelman   /* destroy list of free space and other temporary array(s) */
4863985e5eaSKris Buschelman   ierr = PetscMalloc((ci[pn]+1)*sizeof(int),&cj);CHKERRQ(ierr);
4873985e5eaSKris Buschelman   ierr = MakeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
4883985e5eaSKris Buschelman   ierr = PetscFree(ptadenserow);CHKERRQ(ierr);
4893985e5eaSKris Buschelman 
4903985e5eaSKris Buschelman   /* Allocate space for ca */
4913985e5eaSKris Buschelman   ierr = PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
4923985e5eaSKris Buschelman   ierr = PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));CHKERRQ(ierr);
4933985e5eaSKris Buschelman 
4943985e5eaSKris Buschelman   /* put together the new matrix */
4953985e5eaSKris Buschelman   ierr = MatCreateSeqAIJWithArrays(A->comm,pn,pn,ci,cj,ca,C);CHKERRQ(ierr);
4963985e5eaSKris Buschelman 
4973985e5eaSKris Buschelman   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
4983985e5eaSKris Buschelman   /* Since these are PETSc arrays, change flags to free them as necessary. */
4993985e5eaSKris Buschelman   c = (Mat_SeqAIJ *)((*C)->data);
5003985e5eaSKris Buschelman   c->freedata = PETSC_TRUE;
5013985e5eaSKris Buschelman   c->nonew    = 0;
5023985e5eaSKris Buschelman 
5033985e5eaSKris Buschelman   /* Clean up. */
5043985e5eaSKris Buschelman   ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
5053985e5eaSKris Buschelman 
5069af31e4aSHong Zhang   ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,PP,0,0);CHKERRQ(ierr);
5073985e5eaSKris Buschelman   PetscFunctionReturn(0);
5083985e5eaSKris Buschelman }
5093985e5eaSKris Buschelman EXTERN_C_END
5103985e5eaSKris Buschelman 
511eb9c0419SKris Buschelman #undef __FUNCT__
5129af31e4aSHong Zhang #define __FUNCT__ "MatPtAPNumeric"
5136849ba73SBarry Smith /*
5149af31e4aSHong Zhang    MatPtAPNumeric - Computes the matrix projection C = P^T * A * P
5154d3841fdSKris Buschelman 
5164d3841fdSKris Buschelman    Collective on Mat
5174d3841fdSKris Buschelman 
5184d3841fdSKris Buschelman    Input Parameters:
5194d3841fdSKris Buschelman +  A - the matrix
5204d3841fdSKris Buschelman -  P - the projection matrix
5214d3841fdSKris Buschelman 
5224d3841fdSKris Buschelman    Output Parameters:
5234d3841fdSKris Buschelman .  C - the product matrix
5244d3841fdSKris Buschelman 
5254d3841fdSKris Buschelman    Notes:
5269af31e4aSHong Zhang    C must have been created by calling MatPtAPSymbolic and must be destroyed by
5274d3841fdSKris Buschelman    the user using MatDeatroy().
5284d3841fdSKris Buschelman 
529170ef064SHong Zhang    This routine is currently only implemented for pairs of AIJ matrices and classes
530170ef064SHong Zhang    which inherit from AIJ.  C will be of type MATAIJ.
5314d3841fdSKris Buschelman 
5324d3841fdSKris Buschelman    Level: intermediate
5334d3841fdSKris Buschelman 
5349af31e4aSHong Zhang .seealso: MatPtAP(),MatPtAPSymbolic(),MatMatMultNumeric()
5356849ba73SBarry Smith */
536dfbe8321SBarry Smith PetscErrorCode MatPtAPNumeric(Mat A,Mat P,Mat C) {
537dfbe8321SBarry Smith   PetscErrorCode ierr;
538534c1384SKris Buschelman   PetscErrorCode (*fA)(Mat,Mat,Mat);
539534c1384SKris Buschelman   PetscErrorCode (*fP)(Mat,Mat,Mat);
540eb9c0419SKris Buschelman 
541eb9c0419SKris Buschelman   PetscFunctionBegin;
542eb9c0419SKris Buschelman 
5434482741eSBarry Smith   PetscValidHeaderSpecific(A,MAT_COOKIE,1);
544c9780b6fSBarry Smith   PetscValidType(A,1);
545eb9c0419SKris Buschelman   MatPreallocated(A);
546eb9c0419SKris Buschelman   if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
547eb9c0419SKris Buschelman   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
548eb9c0419SKris Buschelman 
5494482741eSBarry Smith   PetscValidHeaderSpecific(P,MAT_COOKIE,2);
550c9780b6fSBarry Smith   PetscValidType(P,2);
551eb9c0419SKris Buschelman   MatPreallocated(P);
552eb9c0419SKris Buschelman   if (!P->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
553eb9c0419SKris Buschelman   if (P->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
554eb9c0419SKris Buschelman 
5554482741eSBarry Smith   PetscValidHeaderSpecific(C,MAT_COOKIE,3);
556c9780b6fSBarry Smith   PetscValidType(C,3);
557eb9c0419SKris Buschelman   MatPreallocated(C);
558eb9c0419SKris Buschelman   if (!C->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
559eb9c0419SKris Buschelman   if (C->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
560eb9c0419SKris Buschelman 
561eb9c0419SKris Buschelman   if (P->N!=C->M) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->N,C->M);
562eb9c0419SKris Buschelman   if (P->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->M,A->N);
563eb9c0419SKris Buschelman   if (A->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix 'A' must be square, %d != %d",A->M,A->N);
564eb9c0419SKris Buschelman   if (P->N!=C->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->N,C->N);
565eb9c0419SKris Buschelman 
566534c1384SKris Buschelman   /* For now, we do not dispatch based on the type of A and P */
567534c1384SKris Buschelman   /* When implementations like _SeqAIJ_MAIJ exist, attack the multiple dispatch problem. */
568534c1384SKris Buschelman   fA = A->ops->ptapnumeric;
569534c1384SKris Buschelman   if (!fA) SETERRQ1(PETSC_ERR_SUP,"MatPtAPNumeric not supported for A of type %s",A->type_name);
570534c1384SKris Buschelman   fP = P->ops->ptapnumeric;
571534c1384SKris Buschelman   if (!fP) SETERRQ1(PETSC_ERR_SUP,"MatPtAPNumeric not supported for P of type %s",P->type_name);
572534c1384SKris Buschelman   if (fP!=fA) SETERRQ2(PETSC_ERR_ARG_INCOMP,"MatPtAPNumeric requires A, %s, to be compatible with P, %s",A->type_name,P->type_name);
5734d3841fdSKris Buschelman 
574534c1384SKris Buschelman   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
575534c1384SKris Buschelman   ierr = (*fA)(A,P,C);CHKERRQ(ierr);
576534c1384SKris Buschelman   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
577eb9c0419SKris Buschelman 
578eb9c0419SKris Buschelman   PetscFunctionReturn(0);
579eb9c0419SKris Buschelman }
580eb9c0419SKris Buschelman 
581eb9c0419SKris Buschelman #undef __FUNCT__
5829af31e4aSHong Zhang #define __FUNCT__ "MatPtAPNumeric_SeqAIJ_SeqAIJ"
583dfbe8321SBarry Smith PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat A,Mat P,Mat C)
584dfbe8321SBarry Smith {
585dfbe8321SBarry Smith   PetscErrorCode ierr;
586d20bfe6fSHong Zhang   int            flops=0;
587d20bfe6fSHong Zhang   Mat_SeqAIJ     *a  = (Mat_SeqAIJ *) A->data;
588d20bfe6fSHong Zhang   Mat_SeqAIJ     *p  = (Mat_SeqAIJ *) P->data;
589d20bfe6fSHong Zhang   Mat_SeqAIJ     *c  = (Mat_SeqAIJ *) C->data;
590d20bfe6fSHong Zhang   int            *ai=a->i,*aj=a->j,*apj,*apjdense,*pi=p->i,*pj=p->j,*pJ=p->j,*pjj;
591d20bfe6fSHong Zhang   int            *ci=c->i,*cj=c->j,*cjj;
592d20bfe6fSHong Zhang   int            am=A->M,cn=C->N,cm=C->M;
593d20bfe6fSHong Zhang   int            i,j,k,anzi,pnzi,apnzj,nextap,pnzj,prow,crow;
594d20bfe6fSHong Zhang   MatScalar      *aa=a->a,*apa,*pa=p->a,*pA=p->a,*paj,*ca=c->a,*caj;
595eb9c0419SKris Buschelman 
596eb9c0419SKris Buschelman   PetscFunctionBegin;
597d20bfe6fSHong Zhang   /* Allocate temporary array for storage of one row of A*P */
598d20bfe6fSHong Zhang   ierr = PetscMalloc(cn*(sizeof(MatScalar)+2*sizeof(int)),&apa);CHKERRQ(ierr);
599d20bfe6fSHong Zhang   ierr = PetscMemzero(apa,cn*(sizeof(MatScalar)+2*sizeof(int)));CHKERRQ(ierr);
600eb9c0419SKris Buschelman 
601d20bfe6fSHong Zhang   apj      = (int *)(apa + cn);
602d20bfe6fSHong Zhang   apjdense = apj + cn;
603d20bfe6fSHong Zhang 
604d20bfe6fSHong Zhang   /* Clear old values in C */
605d20bfe6fSHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
606d20bfe6fSHong Zhang 
607d20bfe6fSHong Zhang   for (i=0;i<am;i++) {
608d20bfe6fSHong Zhang     /* Form sparse row of A*P */
609d20bfe6fSHong Zhang     anzi  = ai[i+1] - ai[i];
610d20bfe6fSHong Zhang     apnzj = 0;
611d20bfe6fSHong Zhang     for (j=0;j<anzi;j++) {
612d20bfe6fSHong Zhang       prow = *aj++;
613d20bfe6fSHong Zhang       pnzj = pi[prow+1] - pi[prow];
614d20bfe6fSHong Zhang       pjj  = pj + pi[prow];
615d20bfe6fSHong Zhang       paj  = pa + pi[prow];
616d20bfe6fSHong Zhang       for (k=0;k<pnzj;k++) {
617d20bfe6fSHong Zhang         if (!apjdense[pjj[k]]) {
618d20bfe6fSHong Zhang           apjdense[pjj[k]] = -1;
619d20bfe6fSHong Zhang           apj[apnzj++]     = pjj[k];
620d20bfe6fSHong Zhang         }
621d20bfe6fSHong Zhang         apa[pjj[k]] += (*aa)*paj[k];
622d20bfe6fSHong Zhang       }
623d20bfe6fSHong Zhang       flops += 2*pnzj;
624d20bfe6fSHong Zhang       aa++;
625d20bfe6fSHong Zhang     }
626d20bfe6fSHong Zhang 
627d20bfe6fSHong Zhang     /* Sort the j index array for quick sparse axpy. */
628d20bfe6fSHong Zhang     ierr = PetscSortInt(apnzj,apj);CHKERRQ(ierr);
629d20bfe6fSHong Zhang 
630d20bfe6fSHong Zhang     /* Compute P^T*A*P using outer product (P^T)[:,j]*(A*P)[j,:]. */
631d20bfe6fSHong Zhang     pnzi = pi[i+1] - pi[i];
632d20bfe6fSHong Zhang     for (j=0;j<pnzi;j++) {
633d20bfe6fSHong Zhang       nextap = 0;
634d20bfe6fSHong Zhang       crow   = *pJ++;
635d20bfe6fSHong Zhang       cjj    = cj + ci[crow];
636d20bfe6fSHong Zhang       caj    = ca + ci[crow];
637d20bfe6fSHong Zhang       /* Perform sparse axpy operation.  Note cjj includes apj. */
638d20bfe6fSHong Zhang       for (k=0;nextap<apnzj;k++) {
639d20bfe6fSHong Zhang         if (cjj[k]==apj[nextap]) {
640d20bfe6fSHong Zhang           caj[k] += (*pA)*apa[apj[nextap++]];
641d20bfe6fSHong Zhang         }
642d20bfe6fSHong Zhang       }
643d20bfe6fSHong Zhang       flops += 2*apnzj;
644d20bfe6fSHong Zhang       pA++;
645d20bfe6fSHong Zhang     }
646d20bfe6fSHong Zhang 
647d20bfe6fSHong Zhang     /* Zero the current row info for A*P */
648d20bfe6fSHong Zhang     for (j=0;j<apnzj;j++) {
649d20bfe6fSHong Zhang       apa[apj[j]]      = 0.;
650d20bfe6fSHong Zhang       apjdense[apj[j]] = 0;
651d20bfe6fSHong Zhang     }
652d20bfe6fSHong Zhang   }
653d20bfe6fSHong Zhang 
654d20bfe6fSHong Zhang   /* Assemble the final matrix and clean up */
655d20bfe6fSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
656d20bfe6fSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
657d20bfe6fSHong Zhang   ierr = PetscFree(apa);CHKERRQ(ierr);
658d20bfe6fSHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
659d20bfe6fSHong Zhang 
660eb9c0419SKris Buschelman   PetscFunctionReturn(0);
661eb9c0419SKris Buschelman }
6620e36024fSHong Zhang 
6630e36024fSHong Zhang /* Compute C = P[rstart:rend,:]^T * A * P of seqaij matrices - used by MatPtAP_MPIAIJ_MPIAIJ() */
6640e36024fSHong Zhang 
6650e36024fSHong Zhang #undef __FUNCT__
6660e36024fSHong Zhang #define __FUNCT__ "MatPtAPNumeric_SeqAIJ_SeqAIJ_ReducedPt"
6677f79fd58SMatthew Knepley /*@C
668e9b43d0fSSatish Balay    MatPtAPNumeric_SeqAIJ_SeqAIJ_ReducedPt - Creates C = P[rstart:rend,:]^T * A * P of seqaij matrices,
669b90dcfe3SHong Zhang                 used by MatPtAP_MPIAIJ_MPIAIJ()
670b90dcfe3SHong Zhang 
671b90dcfe3SHong Zhang    Collective on Mat
672b90dcfe3SHong Zhang 
673b90dcfe3SHong Zhang    Input Parameters:
674b90dcfe3SHong Zhang +  A - the matrix in seqaij format
675b90dcfe3SHong Zhang .  P - the projection matrix in seqaij format
676b90dcfe3SHong Zhang .  scall - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX
677b90dcfe3SHong Zhang .  fill - expected fill (not being used when scall==MAT_REUSE_MATRIX)
678b90dcfe3SHong Zhang .  prstart, prend - the starting and ending-1-th row of the matrix P to be used for transpose
679b90dcfe3SHong Zhang 
680b90dcfe3SHong Zhang    Output Parameters:
681b90dcfe3SHong Zhang .  C - the product matrix in seqaij format
682b90dcfe3SHong Zhang 
683b90dcfe3SHong Zhang    Level: developer
684b90dcfe3SHong Zhang 
685b90dcfe3SHong Zhang .seealso: MatPtAPSymbolic(),MatPtAPNumeric(),MatMatMult()
686b90dcfe3SHong Zhang @*/
6870e36024fSHong Zhang PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_ReducedPt(Mat A,Mat P,int prstart,int prend,Mat C)
6880e36024fSHong Zhang {
6890e36024fSHong Zhang   PetscErrorCode ierr;
6900e36024fSHong Zhang   int            flops=0;
6910e36024fSHong Zhang   Mat_SeqAIJ     *a  = (Mat_SeqAIJ *) A->data;
6920e36024fSHong Zhang   Mat_SeqAIJ     *p  = (Mat_SeqAIJ *) P->data;
6930e36024fSHong Zhang   Mat_SeqAIJ     *c  = (Mat_SeqAIJ *) C->data;
69420d4747cSHong Zhang   int            *ai=a->i,*aj=a->j,*apj,*apjdense;
69520d4747cSHong Zhang   int            *pi=p->i,*pj=p->j,*pJ=p->j+pi[prstart],*pjj;
6960e36024fSHong Zhang   int            *ci=c->i,*cj=c->j,*cjj;
6970e36024fSHong Zhang   int            am=A->M,cn=C->N,cm=C->M;
6980e36024fSHong Zhang   int            i,j,k,anzi,pnzi,apnzj,nextap,pnzj,prow,crow;
6990e36024fSHong Zhang   MatScalar      *aa=a->a,*apa,*pa=p->a,*pA=p->a,*paj,*ca=c->a,*caj;
7000e36024fSHong Zhang 
7010e36024fSHong Zhang   PetscFunctionBegin;
7020e36024fSHong Zhang   pA=p->a+pi[prstart];
7030e36024fSHong Zhang   /* Allocate temporary array for storage of one row of A*P */
7040e36024fSHong Zhang   ierr = PetscMalloc(cn*(sizeof(MatScalar)+2*sizeof(int)),&apa);CHKERRQ(ierr);
7050e36024fSHong Zhang   ierr = PetscMemzero(apa,cn*(sizeof(MatScalar)+2*sizeof(int)));CHKERRQ(ierr);
7060e36024fSHong Zhang 
7070e36024fSHong Zhang   apj      = (int *)(apa + cn);
7080e36024fSHong Zhang   apjdense = apj + cn;
7090e36024fSHong Zhang 
7100e36024fSHong Zhang   /* Clear old values in C */
7110e36024fSHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
7120e36024fSHong Zhang 
7130e36024fSHong Zhang   for (i=0;i<am;i++) {
7140e36024fSHong Zhang     /* Form sparse row of A*P */
7150e36024fSHong Zhang     anzi  = ai[i+1] - ai[i];
7160e36024fSHong Zhang     apnzj = 0;
7170e36024fSHong Zhang     for (j=0;j<anzi;j++) {
7180e36024fSHong Zhang       prow = *aj++;
7190e36024fSHong Zhang       pnzj = pi[prow+1] - pi[prow];
7200e36024fSHong Zhang       pjj  = pj + pi[prow];
7210e36024fSHong Zhang       paj  = pa + pi[prow];
7220e36024fSHong Zhang       for (k=0;k<pnzj;k++) {
7230e36024fSHong Zhang         if (!apjdense[pjj[k]]) {
7240e36024fSHong Zhang           apjdense[pjj[k]] = -1;
7250e36024fSHong Zhang           apj[apnzj++]     = pjj[k];
7260e36024fSHong Zhang         }
7270e36024fSHong Zhang         apa[pjj[k]] += (*aa)*paj[k];
7280e36024fSHong Zhang       }
7290e36024fSHong Zhang       flops += 2*pnzj;
7300e36024fSHong Zhang       aa++;
7310e36024fSHong Zhang     }
7320e36024fSHong Zhang 
7330e36024fSHong Zhang     /* Sort the j index array for quick sparse axpy. */
7340e36024fSHong Zhang     ierr = PetscSortInt(apnzj,apj);CHKERRQ(ierr);
7350e36024fSHong Zhang 
7360e36024fSHong Zhang     /* Compute P[[prstart:prend,:]^T*A*P using outer product (P^T)[:,j+prstart]*(A*P)[j,:]. */
7370e36024fSHong Zhang     pnzi = pi[i+1+prstart] - pi[i+prstart];
7380e36024fSHong Zhang     for (j=0;j<pnzi;j++) {
7390e36024fSHong Zhang       nextap = 0;
7400e36024fSHong Zhang       crow   = *pJ++;
7410e36024fSHong Zhang       cjj    = cj + ci[crow];
7420e36024fSHong Zhang       caj    = ca + ci[crow];
7430e36024fSHong Zhang       /* Perform sparse axpy operation.  Note cjj includes apj. */
7440e36024fSHong Zhang       for (k=0;nextap<apnzj;k++) {
7450e36024fSHong Zhang         if (cjj[k]==apj[nextap]) {
7460e36024fSHong Zhang           caj[k] += (*pA)*apa[apj[nextap++]];
7470e36024fSHong Zhang         }
7480e36024fSHong Zhang       }
7490e36024fSHong Zhang       flops += 2*apnzj;
7500e36024fSHong Zhang       pA++;
7510e36024fSHong Zhang     }
7520e36024fSHong Zhang 
7530e36024fSHong Zhang     /* Zero the current row info for A*P */
7540e36024fSHong Zhang     for (j=0;j<apnzj;j++) {
7550e36024fSHong Zhang       apa[apj[j]]      = 0.;
7560e36024fSHong Zhang       apjdense[apj[j]] = 0;
7570e36024fSHong Zhang     }
7580e36024fSHong Zhang   }
7590e36024fSHong Zhang 
7600e36024fSHong Zhang   /* Assemble the final matrix and clean up */
7610e36024fSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7620e36024fSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7630e36024fSHong Zhang   ierr = PetscFree(apa);CHKERRQ(ierr);
7640e36024fSHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
7650e36024fSHong Zhang 
7660e36024fSHong Zhang   PetscFunctionReturn(0);
7670e36024fSHong Zhang }
7680e36024fSHong Zhang 
7690e36024fSHong Zhang #undef __FUNCT__
7700e36024fSHong Zhang #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqAIJ_ReducedPt"
7710e36024fSHong Zhang PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_ReducedPt(Mat A,Mat P,PetscReal fill,int prstart,int prend,Mat *C) {
7720e36024fSHong Zhang   PetscErrorCode ierr;
7730e36024fSHong Zhang   FreeSpaceList  free_space=PETSC_NULL,current_space=PETSC_NULL;
7740e36024fSHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*p=(Mat_SeqAIJ*)P->data,*c;
7750e36024fSHong Zhang   int            *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj;
7760e36024fSHong Zhang   int            *ci,*cj,*denserow,*sparserow,*ptadenserow,*ptasparserow,*ptaj;
7770e36024fSHong Zhang   int            an=A->N,am=A->M,pn=P->N,pm=P->M;
7780e36024fSHong Zhang   int            i,j,k,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi;
7790e36024fSHong Zhang   MatScalar      *ca;
7800e36024fSHong Zhang   Mat            *psub,P_sub;
7810e36024fSHong Zhang   IS             isrow,iscol;
7820e36024fSHong Zhang   int            m = prend - prstart;
7830b89d903Svictorle 
7840b89d903Svictorle   PetscFunctionBegin;
7850b89d903Svictorle   /* Get ij structure of P[rstart:rend,:]^T */
7860e36024fSHong Zhang   ierr = ISCreateStride(PETSC_COMM_SELF,m,prstart,1,&isrow);CHKERRQ(ierr);
7870e36024fSHong Zhang   ierr = ISCreateStride(PETSC_COMM_SELF,P->n,0,1,&iscol);CHKERRQ(ierr);
7880e36024fSHong Zhang   ierr = MatGetSubMatrices(P,1,&isrow,&iscol,MAT_INITIAL_MATRIX,&psub);CHKERRQ(ierr);
7890e36024fSHong Zhang   ierr = ISDestroy(isrow);CHKERRQ(ierr);
7900e36024fSHong Zhang   ierr = ISDestroy(iscol);CHKERRQ(ierr);
7910e36024fSHong Zhang   P_sub = psub[0];
7920e36024fSHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(P_sub,&pti,&ptj);CHKERRQ(ierr);
7930e36024fSHong Zhang   ierr = MatDestroyMatrices(1,&psub);CHKERRQ(ierr);
7940e36024fSHong Zhang   ptJ=ptj;
7950e36024fSHong Zhang 
7960e36024fSHong Zhang   /* Allocate ci array, arrays for fill computation and */
7970e36024fSHong Zhang   /* free space for accumulating nonzero column info */
7980e36024fSHong Zhang   ierr = PetscMalloc((pn+1)*sizeof(int),&ci);CHKERRQ(ierr);
7990e36024fSHong Zhang   ci[0] = 0;
8000e36024fSHong Zhang 
8010e36024fSHong Zhang   ierr = PetscMalloc((2*pn+2*an+1)*sizeof(int),&ptadenserow);CHKERRQ(ierr);
8020e36024fSHong Zhang   ierr = PetscMemzero(ptadenserow,(2*pn+2*an+1)*sizeof(int));CHKERRQ(ierr);
8030e36024fSHong Zhang   ptasparserow = ptadenserow  + an;
8040e36024fSHong Zhang   denserow     = ptasparserow + an;
8050e36024fSHong Zhang   sparserow    = denserow     + pn;
8060e36024fSHong Zhang 
8070e36024fSHong Zhang   /* Set initial free space to be nnz(A) scaled by aspect ratio of P. */
8080e36024fSHong Zhang   /* This should be reasonable if sparsity of PtAP is similar to that of A. */
809b90dcfe3SHong Zhang   ierr          = GetMoreSpace((int)(fill*ai[am]/pm)*pn,&free_space);
8100e36024fSHong Zhang   current_space = free_space;
8110e36024fSHong Zhang 
8120e36024fSHong Zhang   /* Determine symbolic info for each row of C: */
8130e36024fSHong Zhang   for (i=0;i<pn;i++) {
8140e36024fSHong Zhang     ptnzi  = pti[i+1] - pti[i];
8150e36024fSHong Zhang     ptanzi = 0;
8160e36024fSHong Zhang     /* Determine symbolic row of PtA_reduced: */
8170e36024fSHong Zhang     for (j=0;j<ptnzi;j++) {
8180e36024fSHong Zhang       arow = *ptJ++;
8190e36024fSHong Zhang       anzj = ai[arow+1] - ai[arow];
8200e36024fSHong Zhang       ajj  = aj + ai[arow];
8210e36024fSHong Zhang       for (k=0;k<anzj;k++) {
8220e36024fSHong Zhang         if (!ptadenserow[ajj[k]]) {
8230e36024fSHong Zhang           ptadenserow[ajj[k]]    = -1;
8240e36024fSHong Zhang           ptasparserow[ptanzi++] = ajj[k];
8250e36024fSHong Zhang         }
8260e36024fSHong Zhang       }
8270e36024fSHong Zhang     }
8280e36024fSHong Zhang       /* Using symbolic info for row of PtA, determine symbolic info for row of C: */
8290e36024fSHong Zhang     ptaj = ptasparserow;
8300e36024fSHong Zhang     cnzi   = 0;
8310e36024fSHong Zhang     for (j=0;j<ptanzi;j++) {
8320e36024fSHong Zhang       prow = *ptaj++;
8330e36024fSHong Zhang       pnzj = pi[prow+1] - pi[prow];
8340e36024fSHong Zhang       pjj  = pj + pi[prow];
8350e36024fSHong Zhang       for (k=0;k<pnzj;k++) {
8360e36024fSHong Zhang         if (!denserow[pjj[k]]) {
8370e36024fSHong Zhang             denserow[pjj[k]]  = -1;
8380e36024fSHong Zhang             sparserow[cnzi++] = pjj[k];
8390e36024fSHong Zhang         }
8400e36024fSHong Zhang       }
8410e36024fSHong Zhang     }
8420e36024fSHong Zhang 
8430e36024fSHong Zhang     /* sort sparserow */
8440e36024fSHong Zhang     ierr = PetscSortInt(cnzi,sparserow);CHKERRQ(ierr);
8450e36024fSHong Zhang 
8460e36024fSHong Zhang     /* If free space is not available, make more free space */
8470e36024fSHong Zhang     /* Double the amount of total space in the list */
8480e36024fSHong Zhang     if (current_space->local_remaining<cnzi) {
8490e36024fSHong Zhang       ierr = GetMoreSpace(current_space->total_array_size,&current_space);CHKERRQ(ierr);
8500e36024fSHong Zhang     }
8510e36024fSHong Zhang 
8520e36024fSHong Zhang     /* Copy data into free space, and zero out denserows */
8530e36024fSHong Zhang     ierr = PetscMemcpy(current_space->array,sparserow,cnzi*sizeof(int));CHKERRQ(ierr);
8540e36024fSHong Zhang     current_space->array           += cnzi;
8550e36024fSHong Zhang     current_space->local_used      += cnzi;
8560e36024fSHong Zhang     current_space->local_remaining -= cnzi;
8570e36024fSHong Zhang 
8580e36024fSHong Zhang     for (j=0;j<ptanzi;j++) {
8590e36024fSHong Zhang       ptadenserow[ptasparserow[j]] = 0;
8600e36024fSHong Zhang     }
8610e36024fSHong Zhang     for (j=0;j<cnzi;j++) {
8620e36024fSHong Zhang       denserow[sparserow[j]] = 0;
8630e36024fSHong Zhang     }
8640e36024fSHong Zhang       /* Aside: Perhaps we should save the pta info for the numerical factorization. */
8650e36024fSHong Zhang       /*        For now, we will recompute what is needed. */
8660e36024fSHong Zhang     ci[i+1] = ci[i] + cnzi;
8670e36024fSHong Zhang   }
8680e36024fSHong Zhang   /* nnz is now stored in ci[ptm], column indices are in the list of free space */
8690e36024fSHong Zhang   /* Allocate space for cj, initialize cj, and */
8700e36024fSHong Zhang   /* destroy list of free space and other temporary array(s) */
8710e36024fSHong Zhang   ierr = PetscMalloc((ci[pn]+1)*sizeof(int),&cj);CHKERRQ(ierr);
8720e36024fSHong Zhang   ierr = MakeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
8730e36024fSHong Zhang   ierr = PetscFree(ptadenserow);CHKERRQ(ierr);
8740e36024fSHong Zhang 
8750e36024fSHong Zhang   /* Allocate space for ca */
8760e36024fSHong Zhang   ierr = PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
8770e36024fSHong Zhang   ierr = PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));CHKERRQ(ierr);
8780e36024fSHong Zhang 
8790e36024fSHong Zhang   /* put together the new matrix */
8800e36024fSHong Zhang   ierr = MatCreateSeqAIJWithArrays(A->comm,pn,pn,ci,cj,ca,C);CHKERRQ(ierr);
8810e36024fSHong Zhang 
8820e36024fSHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
8830e36024fSHong Zhang   /* Since these are PETSc arrays, change flags to free them as necessary. */
8840e36024fSHong Zhang   c = (Mat_SeqAIJ *)((*C)->data);
8850e36024fSHong Zhang   c->freedata = PETSC_TRUE;
8860e36024fSHong Zhang   c->nonew    = 0;
8870e36024fSHong Zhang 
8880e36024fSHong Zhang   /* Clean up. */
8890e36024fSHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
8900e36024fSHong Zhang 
8910e36024fSHong Zhang   PetscFunctionReturn(0);
8920e36024fSHong Zhang }
8930e36024fSHong Zhang 
8940e36024fSHong Zhang #undef __FUNCT__
8950e36024fSHong Zhang #define __FUNCT__ "MatPtAP_SeqAIJ_SeqAIJ_ReducedPt"
8960e36024fSHong Zhang PetscErrorCode MatPtAP_SeqAIJ_SeqAIJ_ReducedPt(Mat A,Mat P,MatReuse scall,PetscReal fill,int prstart,int prend,Mat *C)
8970e36024fSHong Zhang {
8980e36024fSHong Zhang   PetscErrorCode ierr;
8990e36024fSHong Zhang   PetscFunctionBegin;
9000e36024fSHong Zhang   if (A->m != prend-prstart) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",A->m,prend-prstart);
9010e36024fSHong Zhang   if (prend-prstart > P->m) SETERRQ2(PETSC_ERR_ARG_SIZ," prend-prstart %d cannot be larger than P->m %d",prend-prstart,P->m);
9020e36024fSHong Zhang   if (scall == MAT_INITIAL_MATRIX){
9030e36024fSHong Zhang     ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ_ReducedPt(A,P,fill,prstart,prend,C);CHKERRQ(ierr);
9040e36024fSHong Zhang   }
9050e36024fSHong Zhang 
9060e36024fSHong Zhang   ierr = MatPtAPNumeric_SeqAIJ_SeqAIJ_ReducedPt(A,P,prstart,prend,*C);CHKERRQ(ierr);
9070e36024fSHong Zhang 
9080e36024fSHong Zhang   PetscFunctionReturn(0);
9090e36024fSHong Zhang }
9100e36024fSHong Zhang 
911