xref: /petsc/src/mat/impls/aij/seq/matptap.c (revision e9b43d0f224bbe4aed6dc1cb9ef9a101ac74b063)
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){
84b90dcfe3SHong Zhang     ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ(A,P,fill,C);CHKERRQ(ierr);/* numeric product is computed as well */
85b90dcfe3SHong Zhang   } else if (scall == MAT_REUSE_MATRIX){
86b90dcfe3SHong Zhang     ierr = MatPtAPNumeric_MPIAIJ_MPIAIJ(A,P,*C);CHKERRQ(ierr);
87b90dcfe3SHong Zhang   } else {
88b90dcfe3SHong Zhang     SETERRQ1(PETSC_ERR_ARG_WRONG,"Invalid MatReuse %d",scall);
89b90dcfe3SHong Zhang   }
90b90dcfe3SHong Zhang   PetscFunctionReturn(0);
91b90dcfe3SHong Zhang }
92b90dcfe3SHong Zhang 
93b90dcfe3SHong Zhang #undef __FUNCT__
94b90dcfe3SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ"
95b90dcfe3SHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
96b90dcfe3SHong Zhang {
97b90dcfe3SHong Zhang   PetscErrorCode    ierr;
9825616d81SHong Zhang   Mat               P_seq,A_loc,C_seq;
990e36024fSHong Zhang   int               prstart,prend,m=P->m;
10025616d81SHong Zhang   IS                isrowp,iscolp;
101ff134f7aSHong Zhang 
102ff134f7aSHong Zhang   PetscFunctionBegin;
10325616d81SHong Zhang   /* get P_seq = submatrix of P by taking rows of P that equal to nonzero col of A */
104b90dcfe3SHong Zhang   ierr = MatGetBrowsOfAcols(A,P,MAT_INITIAL_MATRIX,&isrowp,&iscolp,&prstart,&P_seq);CHKERRQ(ierr);
10525616d81SHong Zhang   ierr = ISDestroy(iscolp);CHKERRQ(ierr);
106ff134f7aSHong Zhang 
10725616d81SHong Zhang   /* get A_loc = submatrix of A by taking all local rows of A */
108b90dcfe3SHong Zhang   ierr = MatGetLocalMat(A,MAT_INITIAL_MATRIX,PETSC_NULL,&isrowp,&A_loc);CHKERRQ(ierr);
10925616d81SHong Zhang   ierr = ISDestroy(isrowp);CHKERRQ(ierr);
1100e36024fSHong Zhang 
11125616d81SHong Zhang   /* compute C_seq = P_loc^T * A_loc * P_seq */
112ff134f7aSHong Zhang   prend   = prstart + m;
113b90dcfe3SHong Zhang   ierr = MatPtAP_SeqAIJ_SeqAIJ_ReducedPt(A_loc,P_seq,MAT_INITIAL_MATRIX,fill,prstart,prend,&C_seq);CHKERRQ(ierr);
11425616d81SHong Zhang   ierr = MatDestroy(P_seq);CHKERRQ(ierr);
11525616d81SHong Zhang   ierr = MatDestroy(A_loc);CHKERRQ(ierr);
116b90dcfe3SHong Zhang 
117b90dcfe3SHong Zhang   /* add C_seq into mpi C */
118b90dcfe3SHong Zhang   ierr = MatMerge_SeqsToMPI(A->comm,C_seq,P->n,P->n,MAT_INITIAL_MATRIX,C);CHKERRQ(ierr);
119b90dcfe3SHong Zhang   /* ierr = MatDestroy(C_seq);CHKERRQ(ierr); */
120b90dcfe3SHong Zhang 
121ff134f7aSHong Zhang   PetscFunctionReturn(0);
122ff134f7aSHong Zhang }
123ff134f7aSHong Zhang 
124ff134f7aSHong Zhang #undef __FUNCT__
125ff134f7aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ"
126b90dcfe3SHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
127ff134f7aSHong Zhang {
128b90dcfe3SHong Zhang   PetscErrorCode    ierr;
129b90dcfe3SHong Zhang   Mat               P_seq,A_loc,C_seq;
130b90dcfe3SHong Zhang   int               prstart,prend,m=P->m;
131b90dcfe3SHong Zhang   IS                isrowp,iscolp;
132b90dcfe3SHong Zhang   Mat_Merge_SeqsToMPI *merge=(Mat_Merge_SeqsToMPI*)C->spptr;
133ff134f7aSHong Zhang 
134ff134f7aSHong Zhang   PetscFunctionBegin;
135b90dcfe3SHong Zhang   /* get P_seq = submatrix of P by taking rows of P that equal to nonzero col of A */
136b90dcfe3SHong Zhang   ierr = MatGetBrowsOfAcols(A,P,MAT_INITIAL_MATRIX,&isrowp,&iscolp,&prstart,&P_seq);CHKERRQ(ierr);
137b90dcfe3SHong Zhang   ierr = ISDestroy(iscolp);CHKERRQ(ierr);
138ff134f7aSHong Zhang 
139b90dcfe3SHong Zhang   /* get A_loc = submatrix of A by taking all local rows of A */
140b90dcfe3SHong Zhang   ierr = MatGetLocalMat(A,MAT_INITIAL_MATRIX,PETSC_NULL,&isrowp,&A_loc);CHKERRQ(ierr);
141b90dcfe3SHong Zhang   ierr = ISDestroy(isrowp);CHKERRQ(ierr);
142ff134f7aSHong Zhang 
143b90dcfe3SHong Zhang   /* compute C_seq = P_loc^T * A_loc * P_seq */
144b90dcfe3SHong Zhang   prend = prstart + m;
145b90dcfe3SHong Zhang   C_seq = merge->C_seq;
146b90dcfe3SHong Zhang   ierr = MatPtAP_SeqAIJ_SeqAIJ_ReducedPt(A_loc,P_seq,MAT_REUSE_MATRIX,1.0,prstart,prend,&C_seq);CHKERRQ(ierr);
147b90dcfe3SHong Zhang   ierr = MatDestroy(P_seq);CHKERRQ(ierr);
148b90dcfe3SHong Zhang   ierr = MatDestroy(A_loc);CHKERRQ(ierr);
149b90dcfe3SHong Zhang 
150b90dcfe3SHong Zhang   /* add C_seq into mpi C */
151b90dcfe3SHong Zhang   ierr = MatMerge_SeqsToMPI(A->comm,C_seq,P->n,P->n,MAT_REUSE_MATRIX,&C);CHKERRQ(ierr);
152b90dcfe3SHong Zhang 
153ff134f7aSHong Zhang   PetscFunctionReturn(0);
154ff134f7aSHong Zhang }
155ff134f7aSHong Zhang 
156ff134f7aSHong Zhang #undef __FUNCT__
1579af31e4aSHong Zhang #define __FUNCT__ "MatPtAP_SeqAIJ_SeqAIJ"
158dfbe8321SBarry Smith PetscErrorCode MatPtAP_SeqAIJ_SeqAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
1599af31e4aSHong Zhang {
160dfbe8321SBarry Smith   PetscErrorCode ierr;
1619af31e4aSHong Zhang   PetscFunctionBegin;
1629af31e4aSHong Zhang   if (scall == MAT_INITIAL_MATRIX){
163d20bfe6fSHong Zhang     ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
1649af31e4aSHong Zhang     ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ(A,P,fill,C);CHKERRQ(ierr);
165d20bfe6fSHong Zhang     ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
1669af31e4aSHong Zhang   }
167d20bfe6fSHong Zhang   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
1689af31e4aSHong Zhang   ierr = MatPtAPNumeric_SeqAIJ_SeqAIJ(A,P,*C);CHKERRQ(ierr);
169d20bfe6fSHong Zhang   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
1709af31e4aSHong Zhang   PetscFunctionReturn(0);
1719af31e4aSHong Zhang }
1729af31e4aSHong Zhang 
1739af31e4aSHong Zhang #undef __FUNCT__
1749af31e4aSHong Zhang #define __FUNCT__ "MatPtAPSymbolic"
1756849ba73SBarry Smith /*
1769af31e4aSHong Zhang    MatPtAPSymbolic - Creates the (i,j) structure of the matrix projection C = P^T * A * P
1774d3841fdSKris Buschelman 
1784d3841fdSKris Buschelman    Collective on Mat
1794d3841fdSKris Buschelman 
1804d3841fdSKris Buschelman    Input Parameters:
1814d3841fdSKris Buschelman +  A - the matrix
1824d3841fdSKris Buschelman -  P - the projection matrix
1834d3841fdSKris Buschelman 
1844d3841fdSKris Buschelman    Output Parameters:
1854d3841fdSKris Buschelman .  C - the (i,j) structure of the product matrix
1864d3841fdSKris Buschelman 
1874d3841fdSKris Buschelman    Notes:
1884d3841fdSKris Buschelman    C will be created and must be destroyed by the user with MatDestroy().
1894d3841fdSKris Buschelman 
1904d3841fdSKris Buschelman    This routine is currently only implemented for pairs of SeqAIJ matrices and classes
1914d3841fdSKris Buschelman    which inherit from SeqAIJ.  C will be of type MATSEQAIJ.  The product is computed using
1929af31e4aSHong Zhang    this (i,j) structure by calling MatPtAPNumeric().
1934d3841fdSKris Buschelman 
1944d3841fdSKris Buschelman    Level: intermediate
1954d3841fdSKris Buschelman 
1969af31e4aSHong Zhang .seealso: MatPtAP(),MatPtAPNumeric(),MatMatMultSymbolic()
1976849ba73SBarry Smith */
198dfbe8321SBarry Smith PetscErrorCode MatPtAPSymbolic(Mat A,Mat P,PetscReal fill,Mat *C) {
199dfbe8321SBarry Smith   PetscErrorCode ierr;
200534c1384SKris Buschelman   PetscErrorCode (*fA)(Mat,Mat,PetscReal,Mat*);
201534c1384SKris Buschelman   PetscErrorCode (*fP)(Mat,Mat,PetscReal,Mat*);
202eb9c0419SKris Buschelman 
203eb9c0419SKris Buschelman   PetscFunctionBegin;
204eb9c0419SKris Buschelman 
2054482741eSBarry Smith   PetscValidHeaderSpecific(A,MAT_COOKIE,1);
206c9780b6fSBarry Smith   PetscValidType(A,1);
207eb9c0419SKris Buschelman   MatPreallocated(A);
208eb9c0419SKris Buschelman   if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
209eb9c0419SKris Buschelman   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
210eb9c0419SKris Buschelman 
2114482741eSBarry Smith   PetscValidHeaderSpecific(P,MAT_COOKIE,2);
212c9780b6fSBarry Smith   PetscValidType(P,2);
213eb9c0419SKris Buschelman   MatPreallocated(P);
214eb9c0419SKris Buschelman   if (!P->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
215eb9c0419SKris Buschelman   if (P->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
216eb9c0419SKris Buschelman 
2174482741eSBarry Smith   PetscValidPointer(C,3);
2184482741eSBarry Smith 
219eb9c0419SKris Buschelman   if (P->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->M,A->N);
220eb9c0419SKris Buschelman   if (A->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix 'A' must be square, %d != %d",A->M,A->N);
221eb9c0419SKris Buschelman 
222534c1384SKris Buschelman   /* For now, we do not dispatch based on the type of A and P */
223534c1384SKris Buschelman   /* When implementations like _SeqAIJ_MAIJ exist, attack the multiple dispatch problem. */
224534c1384SKris Buschelman   fA = A->ops->ptapsymbolic;
225534c1384SKris Buschelman   if (!fA) SETERRQ1(PETSC_ERR_SUP,"MatPtAPSymbolic not supported for A of type %s",A->type_name);
226534c1384SKris Buschelman   fP = P->ops->ptapsymbolic;
227534c1384SKris Buschelman   if (!fP) SETERRQ1(PETSC_ERR_SUP,"MatPtAPSymbolic not supported for P of type %s",P->type_name);
228534c1384SKris 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);
2294d3841fdSKris Buschelman 
230534c1384SKris Buschelman   ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
231534c1384SKris Buschelman   ierr = (*fA)(A,P,fill,C);CHKERRQ(ierr);
232534c1384SKris Buschelman   ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
233eb9c0419SKris Buschelman 
234eb9c0419SKris Buschelman   PetscFunctionReturn(0);
235eb9c0419SKris Buschelman }
236eb9c0419SKris Buschelman 
237f747e1aeSHong Zhang typedef struct {
238f747e1aeSHong Zhang   Mat    symAP;
239f747e1aeSHong Zhang } Mat_PtAPstruct;
240f747e1aeSHong Zhang 
24178a80504SBarry Smith EXTERN PetscErrorCode MatDestroy_SeqAIJ(Mat);
24278a80504SBarry Smith 
243f747e1aeSHong Zhang #undef __FUNCT__
244f747e1aeSHong Zhang #define __FUNCT__ "MatDestroy_SeqAIJ_PtAP"
245f4a850bbSBarry Smith PetscErrorCode MatDestroy_SeqAIJ_PtAP(Mat A)
246f747e1aeSHong Zhang {
247f4a850bbSBarry Smith   PetscErrorCode    ierr;
248f747e1aeSHong Zhang   Mat_PtAPstruct    *ptap=(Mat_PtAPstruct*)A->spptr;
249f747e1aeSHong Zhang 
250f747e1aeSHong Zhang   PetscFunctionBegin;
251f747e1aeSHong Zhang   ierr = MatDestroy(ptap->symAP);CHKERRQ(ierr);
252f747e1aeSHong Zhang   ierr = PetscFree(ptap);CHKERRQ(ierr);
25378a80504SBarry Smith   ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr);
254f747e1aeSHong Zhang   PetscFunctionReturn(0);
255f747e1aeSHong Zhang }
256f747e1aeSHong Zhang 
257eb9c0419SKris Buschelman #undef __FUNCT__
2589af31e4aSHong Zhang #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqAIJ"
259dfbe8321SBarry Smith PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat P,PetscReal fill,Mat *C) {
260dfbe8321SBarry Smith   PetscErrorCode ierr;
261d20bfe6fSHong Zhang   FreeSpaceList  free_space=PETSC_NULL,current_space=PETSC_NULL;
262d20bfe6fSHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*p=(Mat_SeqAIJ*)P->data,*c;
263d20bfe6fSHong Zhang   int            *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj;
264d20bfe6fSHong Zhang   int            *ci,*cj,*denserow,*sparserow,*ptadenserow,*ptasparserow,*ptaj;
265d20bfe6fSHong Zhang   int            an=A->N,am=A->M,pn=P->N,pm=P->M;
266d20bfe6fSHong Zhang   int            i,j,k,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi;
267d20bfe6fSHong Zhang   MatScalar      *ca;
268eb9c0419SKris Buschelman 
269eb9c0419SKris Buschelman   PetscFunctionBegin;
270d20bfe6fSHong Zhang   /* Get ij structure of P^T */
271eb9c0419SKris Buschelman   ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
272d20bfe6fSHong Zhang   ptJ=ptj;
273eb9c0419SKris Buschelman 
274d20bfe6fSHong Zhang   /* Allocate ci array, arrays for fill computation and */
275d20bfe6fSHong Zhang   /* free space for accumulating nonzero column info */
276d20bfe6fSHong Zhang   ierr = PetscMalloc((pn+1)*sizeof(int),&ci);CHKERRQ(ierr);
277d20bfe6fSHong Zhang   ci[0] = 0;
278eb9c0419SKris Buschelman 
279d20bfe6fSHong Zhang   ierr = PetscMalloc((2*pn+2*an+1)*sizeof(int),&ptadenserow);CHKERRQ(ierr);
280d20bfe6fSHong Zhang   ierr = PetscMemzero(ptadenserow,(2*pn+2*an+1)*sizeof(int));CHKERRQ(ierr);
281d20bfe6fSHong Zhang   ptasparserow = ptadenserow  + an;
282d20bfe6fSHong Zhang   denserow     = ptasparserow + an;
283d20bfe6fSHong Zhang   sparserow    = denserow     + pn;
284eb9c0419SKris Buschelman 
285d20bfe6fSHong Zhang   /* Set initial free space to be nnz(A) scaled by aspect ratio of P. */
286d20bfe6fSHong Zhang   /* This should be reasonable if sparsity of PtAP is similar to that of A. */
287d20bfe6fSHong Zhang   ierr          = GetMoreSpace((ai[am]/pm)*pn,&free_space);
288d20bfe6fSHong Zhang   current_space = free_space;
289d20bfe6fSHong Zhang 
290d20bfe6fSHong Zhang   /* Determine symbolic info for each row of C: */
291d20bfe6fSHong Zhang   for (i=0;i<pn;i++) {
292d20bfe6fSHong Zhang     ptnzi  = pti[i+1] - pti[i];
293d20bfe6fSHong Zhang     ptanzi = 0;
294d20bfe6fSHong Zhang     /* Determine symbolic row of PtA: */
295d20bfe6fSHong Zhang     for (j=0;j<ptnzi;j++) {
296d20bfe6fSHong Zhang       arow = *ptJ++;
297d20bfe6fSHong Zhang       anzj = ai[arow+1] - ai[arow];
298d20bfe6fSHong Zhang       ajj  = aj + ai[arow];
299d20bfe6fSHong Zhang       for (k=0;k<anzj;k++) {
300d20bfe6fSHong Zhang         if (!ptadenserow[ajj[k]]) {
301d20bfe6fSHong Zhang           ptadenserow[ajj[k]]    = -1;
302d20bfe6fSHong Zhang           ptasparserow[ptanzi++] = ajj[k];
303d20bfe6fSHong Zhang         }
304d20bfe6fSHong Zhang       }
305d20bfe6fSHong Zhang     }
306d20bfe6fSHong Zhang       /* Using symbolic info for row of PtA, determine symbolic info for row of C: */
307d20bfe6fSHong Zhang     ptaj = ptasparserow;
308d20bfe6fSHong Zhang     cnzi   = 0;
309d20bfe6fSHong Zhang     for (j=0;j<ptanzi;j++) {
310d20bfe6fSHong Zhang       prow = *ptaj++;
311d20bfe6fSHong Zhang       pnzj = pi[prow+1] - pi[prow];
312d20bfe6fSHong Zhang       pjj  = pj + pi[prow];
313d20bfe6fSHong Zhang       for (k=0;k<pnzj;k++) {
314d20bfe6fSHong Zhang         if (!denserow[pjj[k]]) {
315d20bfe6fSHong Zhang             denserow[pjj[k]]  = -1;
316d20bfe6fSHong Zhang             sparserow[cnzi++] = pjj[k];
317d20bfe6fSHong Zhang         }
318d20bfe6fSHong Zhang       }
319d20bfe6fSHong Zhang     }
320d20bfe6fSHong Zhang 
321d20bfe6fSHong Zhang     /* sort sparserow */
322d20bfe6fSHong Zhang     ierr = PetscSortInt(cnzi,sparserow);CHKERRQ(ierr);
323d20bfe6fSHong Zhang 
324d20bfe6fSHong Zhang     /* If free space is not available, make more free space */
325d20bfe6fSHong Zhang     /* Double the amount of total space in the list */
326d20bfe6fSHong Zhang     if (current_space->local_remaining<cnzi) {
327d20bfe6fSHong Zhang       ierr = GetMoreSpace(current_space->total_array_size,&current_space);CHKERRQ(ierr);
328d20bfe6fSHong Zhang     }
329d20bfe6fSHong Zhang 
330d20bfe6fSHong Zhang     /* Copy data into free space, and zero out denserows */
331d20bfe6fSHong Zhang     ierr = PetscMemcpy(current_space->array,sparserow,cnzi*sizeof(int));CHKERRQ(ierr);
332d20bfe6fSHong Zhang     current_space->array           += cnzi;
333d20bfe6fSHong Zhang     current_space->local_used      += cnzi;
334d20bfe6fSHong Zhang     current_space->local_remaining -= cnzi;
335d20bfe6fSHong Zhang 
336d20bfe6fSHong Zhang     for (j=0;j<ptanzi;j++) {
337d20bfe6fSHong Zhang       ptadenserow[ptasparserow[j]] = 0;
338d20bfe6fSHong Zhang     }
339d20bfe6fSHong Zhang     for (j=0;j<cnzi;j++) {
340d20bfe6fSHong Zhang       denserow[sparserow[j]] = 0;
341d20bfe6fSHong Zhang     }
342d20bfe6fSHong Zhang       /* Aside: Perhaps we should save the pta info for the numerical factorization. */
343d20bfe6fSHong Zhang       /*        For now, we will recompute what is needed. */
344d20bfe6fSHong Zhang     ci[i+1] = ci[i] + cnzi;
345d20bfe6fSHong Zhang   }
346d20bfe6fSHong Zhang   /* nnz is now stored in ci[ptm], column indices are in the list of free space */
347d20bfe6fSHong Zhang   /* Allocate space for cj, initialize cj, and */
348d20bfe6fSHong Zhang   /* destroy list of free space and other temporary array(s) */
349d20bfe6fSHong Zhang   ierr = PetscMalloc((ci[pn]+1)*sizeof(int),&cj);CHKERRQ(ierr);
350d20bfe6fSHong Zhang   ierr = MakeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
351d20bfe6fSHong Zhang   ierr = PetscFree(ptadenserow);CHKERRQ(ierr);
352d20bfe6fSHong Zhang 
353d20bfe6fSHong Zhang   /* Allocate space for ca */
354d20bfe6fSHong Zhang   ierr = PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
355d20bfe6fSHong Zhang   ierr = PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));CHKERRQ(ierr);
356d20bfe6fSHong Zhang 
357d20bfe6fSHong Zhang   /* put together the new matrix */
358d20bfe6fSHong Zhang   ierr = MatCreateSeqAIJWithArrays(A->comm,pn,pn,ci,cj,ca,C);CHKERRQ(ierr);
359d20bfe6fSHong Zhang 
360d20bfe6fSHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
361d20bfe6fSHong Zhang   /* Since these are PETSc arrays, change flags to free them as necessary. */
362d20bfe6fSHong Zhang   c = (Mat_SeqAIJ *)((*C)->data);
363d20bfe6fSHong Zhang   c->freedata = PETSC_TRUE;
364d20bfe6fSHong Zhang   c->nonew    = 0;
365d20bfe6fSHong Zhang 
366d20bfe6fSHong Zhang   /* Clean up. */
367d20bfe6fSHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
368eb9c0419SKris Buschelman 
369eb9c0419SKris Buschelman   PetscFunctionReturn(0);
370eb9c0419SKris Buschelman }
371eb9c0419SKris Buschelman 
3723985e5eaSKris Buschelman #include "src/mat/impls/maij/maij.h"
3733985e5eaSKris Buschelman EXTERN_C_BEGIN
3743985e5eaSKris Buschelman #undef __FUNCT__
3759af31e4aSHong Zhang #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqMAIJ"
376dfbe8321SBarry Smith PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqMAIJ(Mat A,Mat PP,Mat *C) {
3775c66b693SKris Buschelman   /* This routine requires testing -- I don't think it works. */
378dfbe8321SBarry Smith   PetscErrorCode ierr;
3793985e5eaSKris Buschelman   FreeSpaceList  free_space=PETSC_NULL,current_space=PETSC_NULL;
3803985e5eaSKris Buschelman   Mat_SeqMAIJ    *pp=(Mat_SeqMAIJ*)PP->data;
3813985e5eaSKris Buschelman   Mat            P=pp->AIJ;
3823985e5eaSKris Buschelman   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*p=(Mat_SeqAIJ*)P->data,*c;
3833985e5eaSKris Buschelman   int            *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj;
3843985e5eaSKris Buschelman   int            *ci,*cj,*denserow,*sparserow,*ptadenserow,*ptasparserow,*ptaj;
3853985e5eaSKris Buschelman   int            an=A->N,am=A->M,pn=P->N,pm=P->M,ppdof=pp->dof;
386fe05a634SKris Buschelman   int            i,j,k,dof,pdof,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi;
3873985e5eaSKris Buschelman   MatScalar      *ca;
3883985e5eaSKris Buschelman 
3893985e5eaSKris Buschelman   PetscFunctionBegin;
3903985e5eaSKris Buschelman   /* Start timer */
3919af31e4aSHong Zhang   ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,PP,0,0);CHKERRQ(ierr);
3923985e5eaSKris Buschelman 
3933985e5eaSKris Buschelman   /* Get ij structure of P^T */
3943985e5eaSKris Buschelman   ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
3953985e5eaSKris Buschelman 
3963985e5eaSKris Buschelman   /* Allocate ci array, arrays for fill computation and */
3973985e5eaSKris Buschelman   /* free space for accumulating nonzero column info */
3983985e5eaSKris Buschelman   ierr = PetscMalloc((pn+1)*sizeof(int),&ci);CHKERRQ(ierr);
3993985e5eaSKris Buschelman   ci[0] = 0;
4003985e5eaSKris Buschelman 
4013985e5eaSKris Buschelman   ierr = PetscMalloc((2*pn+2*an+1)*sizeof(int),&ptadenserow);CHKERRQ(ierr);
4023985e5eaSKris Buschelman   ierr = PetscMemzero(ptadenserow,(2*pn+2*an+1)*sizeof(int));CHKERRQ(ierr);
4033985e5eaSKris Buschelman   ptasparserow = ptadenserow  + an;
4043985e5eaSKris Buschelman   denserow     = ptasparserow + an;
4053985e5eaSKris Buschelman   sparserow    = denserow     + pn;
4063985e5eaSKris Buschelman 
4073985e5eaSKris Buschelman   /* Set initial free space to be nnz(A) scaled by aspect ratio of P. */
4083985e5eaSKris Buschelman   /* This should be reasonable if sparsity of PtAP is similar to that of A. */
4093985e5eaSKris Buschelman   ierr          = GetMoreSpace((ai[am]/pm)*pn,&free_space);
4103985e5eaSKris Buschelman   current_space = free_space;
4113985e5eaSKris Buschelman 
4123985e5eaSKris Buschelman   /* Determine symbolic info for each row of C: */
4133985e5eaSKris Buschelman   for (i=0;i<pn/ppdof;i++) {
4143985e5eaSKris Buschelman     ptnzi  = pti[i+1] - pti[i];
4153985e5eaSKris Buschelman     ptanzi = 0;
4163985e5eaSKris Buschelman     ptJ    = ptj + pti[i];
4173985e5eaSKris Buschelman     for (dof=0;dof<ppdof;dof++) {
4183985e5eaSKris Buschelman     /* Determine symbolic row of PtA: */
4193985e5eaSKris Buschelman       for (j=0;j<ptnzi;j++) {
4203985e5eaSKris Buschelman         arow = ptJ[j] + dof;
4213985e5eaSKris Buschelman         anzj = ai[arow+1] - ai[arow];
4223985e5eaSKris Buschelman         ajj  = aj + ai[arow];
4233985e5eaSKris Buschelman         for (k=0;k<anzj;k++) {
4243985e5eaSKris Buschelman           if (!ptadenserow[ajj[k]]) {
4253985e5eaSKris Buschelman             ptadenserow[ajj[k]]    = -1;
4263985e5eaSKris Buschelman             ptasparserow[ptanzi++] = ajj[k];
4273985e5eaSKris Buschelman           }
4283985e5eaSKris Buschelman         }
4293985e5eaSKris Buschelman       }
4303985e5eaSKris Buschelman       /* Using symbolic info for row of PtA, determine symbolic info for row of C: */
4313985e5eaSKris Buschelman       ptaj = ptasparserow;
4323985e5eaSKris Buschelman       cnzi   = 0;
4333985e5eaSKris Buschelman       for (j=0;j<ptanzi;j++) {
434fe05a634SKris Buschelman         pdof = *ptaj%dof;
4353985e5eaSKris Buschelman         prow = (*ptaj++)/dof;
4363985e5eaSKris Buschelman         pnzj = pi[prow+1] - pi[prow];
4373985e5eaSKris Buschelman         pjj  = pj + pi[prow];
4383985e5eaSKris Buschelman         for (k=0;k<pnzj;k++) {
439fe05a634SKris Buschelman           if (!denserow[pjj[k]+pdof]) {
440fe05a634SKris Buschelman             denserow[pjj[k]+pdof] = -1;
441fe05a634SKris Buschelman             sparserow[cnzi++]     = pjj[k]+pdof;
4423985e5eaSKris Buschelman           }
4433985e5eaSKris Buschelman         }
4443985e5eaSKris Buschelman       }
4453985e5eaSKris Buschelman 
4463985e5eaSKris Buschelman       /* sort sparserow */
4473985e5eaSKris Buschelman       ierr = PetscSortInt(cnzi,sparserow);CHKERRQ(ierr);
4483985e5eaSKris Buschelman 
4493985e5eaSKris Buschelman       /* If free space is not available, make more free space */
4503985e5eaSKris Buschelman       /* Double the amount of total space in the list */
4513985e5eaSKris Buschelman       if (current_space->local_remaining<cnzi) {
4523985e5eaSKris Buschelman         ierr = GetMoreSpace(current_space->total_array_size,&current_space);CHKERRQ(ierr);
4533985e5eaSKris Buschelman       }
4543985e5eaSKris Buschelman 
4553985e5eaSKris Buschelman       /* Copy data into free space, and zero out denserows */
4563985e5eaSKris Buschelman       ierr = PetscMemcpy(current_space->array,sparserow,cnzi*sizeof(int));CHKERRQ(ierr);
4573985e5eaSKris Buschelman       current_space->array           += cnzi;
4583985e5eaSKris Buschelman       current_space->local_used      += cnzi;
4593985e5eaSKris Buschelman       current_space->local_remaining -= cnzi;
4603985e5eaSKris Buschelman 
4613985e5eaSKris Buschelman       for (j=0;j<ptanzi;j++) {
4623985e5eaSKris Buschelman         ptadenserow[ptasparserow[j]] = 0;
4633985e5eaSKris Buschelman       }
4643985e5eaSKris Buschelman       for (j=0;j<cnzi;j++) {
4653985e5eaSKris Buschelman         denserow[sparserow[j]] = 0;
4663985e5eaSKris Buschelman       }
4673985e5eaSKris Buschelman       /* Aside: Perhaps we should save the pta info for the numerical factorization. */
4683985e5eaSKris Buschelman       /*        For now, we will recompute what is needed. */
4693985e5eaSKris Buschelman       ci[i+1+dof] = ci[i+dof] + cnzi;
4703985e5eaSKris Buschelman     }
4713985e5eaSKris Buschelman   }
4723985e5eaSKris Buschelman   /* nnz is now stored in ci[ptm], column indices are in the list of free space */
4733985e5eaSKris Buschelman   /* Allocate space for cj, initialize cj, and */
4743985e5eaSKris Buschelman   /* destroy list of free space and other temporary array(s) */
4753985e5eaSKris Buschelman   ierr = PetscMalloc((ci[pn]+1)*sizeof(int),&cj);CHKERRQ(ierr);
4763985e5eaSKris Buschelman   ierr = MakeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
4773985e5eaSKris Buschelman   ierr = PetscFree(ptadenserow);CHKERRQ(ierr);
4783985e5eaSKris Buschelman 
4793985e5eaSKris Buschelman   /* Allocate space for ca */
4803985e5eaSKris Buschelman   ierr = PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
4813985e5eaSKris Buschelman   ierr = PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));CHKERRQ(ierr);
4823985e5eaSKris Buschelman 
4833985e5eaSKris Buschelman   /* put together the new matrix */
4843985e5eaSKris Buschelman   ierr = MatCreateSeqAIJWithArrays(A->comm,pn,pn,ci,cj,ca,C);CHKERRQ(ierr);
4853985e5eaSKris Buschelman 
4863985e5eaSKris Buschelman   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
4873985e5eaSKris Buschelman   /* Since these are PETSc arrays, change flags to free them as necessary. */
4883985e5eaSKris Buschelman   c = (Mat_SeqAIJ *)((*C)->data);
4893985e5eaSKris Buschelman   c->freedata = PETSC_TRUE;
4903985e5eaSKris Buschelman   c->nonew    = 0;
4913985e5eaSKris Buschelman 
4923985e5eaSKris Buschelman   /* Clean up. */
4933985e5eaSKris Buschelman   ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
4943985e5eaSKris Buschelman 
4959af31e4aSHong Zhang   ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,PP,0,0);CHKERRQ(ierr);
4963985e5eaSKris Buschelman   PetscFunctionReturn(0);
4973985e5eaSKris Buschelman }
4983985e5eaSKris Buschelman EXTERN_C_END
4993985e5eaSKris Buschelman 
500eb9c0419SKris Buschelman #undef __FUNCT__
5019af31e4aSHong Zhang #define __FUNCT__ "MatPtAPNumeric"
5026849ba73SBarry Smith /*
5039af31e4aSHong Zhang    MatPtAPNumeric - Computes the matrix projection C = P^T * A * P
5044d3841fdSKris Buschelman 
5054d3841fdSKris Buschelman    Collective on Mat
5064d3841fdSKris Buschelman 
5074d3841fdSKris Buschelman    Input Parameters:
5084d3841fdSKris Buschelman +  A - the matrix
5094d3841fdSKris Buschelman -  P - the projection matrix
5104d3841fdSKris Buschelman 
5114d3841fdSKris Buschelman    Output Parameters:
5124d3841fdSKris Buschelman .  C - the product matrix
5134d3841fdSKris Buschelman 
5144d3841fdSKris Buschelman    Notes:
5159af31e4aSHong Zhang    C must have been created by calling MatPtAPSymbolic and must be destroyed by
5164d3841fdSKris Buschelman    the user using MatDeatroy().
5174d3841fdSKris Buschelman 
518170ef064SHong Zhang    This routine is currently only implemented for pairs of AIJ matrices and classes
519170ef064SHong Zhang    which inherit from AIJ.  C will be of type MATAIJ.
5204d3841fdSKris Buschelman 
5214d3841fdSKris Buschelman    Level: intermediate
5224d3841fdSKris Buschelman 
5239af31e4aSHong Zhang .seealso: MatPtAP(),MatPtAPSymbolic(),MatMatMultNumeric()
5246849ba73SBarry Smith */
525dfbe8321SBarry Smith PetscErrorCode MatPtAPNumeric(Mat A,Mat P,Mat C) {
526dfbe8321SBarry Smith   PetscErrorCode ierr;
527534c1384SKris Buschelman   PetscErrorCode (*fA)(Mat,Mat,Mat);
528534c1384SKris Buschelman   PetscErrorCode (*fP)(Mat,Mat,Mat);
529eb9c0419SKris Buschelman 
530eb9c0419SKris Buschelman   PetscFunctionBegin;
531eb9c0419SKris Buschelman 
5324482741eSBarry Smith   PetscValidHeaderSpecific(A,MAT_COOKIE,1);
533c9780b6fSBarry Smith   PetscValidType(A,1);
534eb9c0419SKris Buschelman   MatPreallocated(A);
535eb9c0419SKris Buschelman   if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
536eb9c0419SKris Buschelman   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
537eb9c0419SKris Buschelman 
5384482741eSBarry Smith   PetscValidHeaderSpecific(P,MAT_COOKIE,2);
539c9780b6fSBarry Smith   PetscValidType(P,2);
540eb9c0419SKris Buschelman   MatPreallocated(P);
541eb9c0419SKris Buschelman   if (!P->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
542eb9c0419SKris Buschelman   if (P->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
543eb9c0419SKris Buschelman 
5444482741eSBarry Smith   PetscValidHeaderSpecific(C,MAT_COOKIE,3);
545c9780b6fSBarry Smith   PetscValidType(C,3);
546eb9c0419SKris Buschelman   MatPreallocated(C);
547eb9c0419SKris Buschelman   if (!C->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
548eb9c0419SKris Buschelman   if (C->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
549eb9c0419SKris Buschelman 
550eb9c0419SKris Buschelman   if (P->N!=C->M) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->N,C->M);
551eb9c0419SKris Buschelman   if (P->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->M,A->N);
552eb9c0419SKris Buschelman   if (A->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix 'A' must be square, %d != %d",A->M,A->N);
553eb9c0419SKris Buschelman   if (P->N!=C->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->N,C->N);
554eb9c0419SKris Buschelman 
555534c1384SKris Buschelman   /* For now, we do not dispatch based on the type of A and P */
556534c1384SKris Buschelman   /* When implementations like _SeqAIJ_MAIJ exist, attack the multiple dispatch problem. */
557534c1384SKris Buschelman   fA = A->ops->ptapnumeric;
558534c1384SKris Buschelman   if (!fA) SETERRQ1(PETSC_ERR_SUP,"MatPtAPNumeric not supported for A of type %s",A->type_name);
559534c1384SKris Buschelman   fP = P->ops->ptapnumeric;
560534c1384SKris Buschelman   if (!fP) SETERRQ1(PETSC_ERR_SUP,"MatPtAPNumeric not supported for P of type %s",P->type_name);
561534c1384SKris 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);
5624d3841fdSKris Buschelman 
563534c1384SKris Buschelman   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
564534c1384SKris Buschelman   ierr = (*fA)(A,P,C);CHKERRQ(ierr);
565534c1384SKris Buschelman   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
566eb9c0419SKris Buschelman 
567eb9c0419SKris Buschelman   PetscFunctionReturn(0);
568eb9c0419SKris Buschelman }
569eb9c0419SKris Buschelman 
570eb9c0419SKris Buschelman #undef __FUNCT__
5719af31e4aSHong Zhang #define __FUNCT__ "MatPtAPNumeric_SeqAIJ_SeqAIJ"
572dfbe8321SBarry Smith PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat A,Mat P,Mat C)
573dfbe8321SBarry Smith {
574dfbe8321SBarry Smith   PetscErrorCode ierr;
575d20bfe6fSHong Zhang   int            flops=0;
576d20bfe6fSHong Zhang   Mat_SeqAIJ     *a  = (Mat_SeqAIJ *) A->data;
577d20bfe6fSHong Zhang   Mat_SeqAIJ     *p  = (Mat_SeqAIJ *) P->data;
578d20bfe6fSHong Zhang   Mat_SeqAIJ     *c  = (Mat_SeqAIJ *) C->data;
579d20bfe6fSHong Zhang   int            *ai=a->i,*aj=a->j,*apj,*apjdense,*pi=p->i,*pj=p->j,*pJ=p->j,*pjj;
580d20bfe6fSHong Zhang   int            *ci=c->i,*cj=c->j,*cjj;
581d20bfe6fSHong Zhang   int            am=A->M,cn=C->N,cm=C->M;
582d20bfe6fSHong Zhang   int            i,j,k,anzi,pnzi,apnzj,nextap,pnzj,prow,crow;
583d20bfe6fSHong Zhang   MatScalar      *aa=a->a,*apa,*pa=p->a,*pA=p->a,*paj,*ca=c->a,*caj;
584eb9c0419SKris Buschelman 
585eb9c0419SKris Buschelman   PetscFunctionBegin;
586d20bfe6fSHong Zhang   /* Allocate temporary array for storage of one row of A*P */
587d20bfe6fSHong Zhang   ierr = PetscMalloc(cn*(sizeof(MatScalar)+2*sizeof(int)),&apa);CHKERRQ(ierr);
588d20bfe6fSHong Zhang   ierr = PetscMemzero(apa,cn*(sizeof(MatScalar)+2*sizeof(int)));CHKERRQ(ierr);
589eb9c0419SKris Buschelman 
590d20bfe6fSHong Zhang   apj      = (int *)(apa + cn);
591d20bfe6fSHong Zhang   apjdense = apj + cn;
592d20bfe6fSHong Zhang 
593d20bfe6fSHong Zhang   /* Clear old values in C */
594d20bfe6fSHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
595d20bfe6fSHong Zhang 
596d20bfe6fSHong Zhang   for (i=0;i<am;i++) {
597d20bfe6fSHong Zhang     /* Form sparse row of A*P */
598d20bfe6fSHong Zhang     anzi  = ai[i+1] - ai[i];
599d20bfe6fSHong Zhang     apnzj = 0;
600d20bfe6fSHong Zhang     for (j=0;j<anzi;j++) {
601d20bfe6fSHong Zhang       prow = *aj++;
602d20bfe6fSHong Zhang       pnzj = pi[prow+1] - pi[prow];
603d20bfe6fSHong Zhang       pjj  = pj + pi[prow];
604d20bfe6fSHong Zhang       paj  = pa + pi[prow];
605d20bfe6fSHong Zhang       for (k=0;k<pnzj;k++) {
606d20bfe6fSHong Zhang         if (!apjdense[pjj[k]]) {
607d20bfe6fSHong Zhang           apjdense[pjj[k]] = -1;
608d20bfe6fSHong Zhang           apj[apnzj++]     = pjj[k];
609d20bfe6fSHong Zhang         }
610d20bfe6fSHong Zhang         apa[pjj[k]] += (*aa)*paj[k];
611d20bfe6fSHong Zhang       }
612d20bfe6fSHong Zhang       flops += 2*pnzj;
613d20bfe6fSHong Zhang       aa++;
614d20bfe6fSHong Zhang     }
615d20bfe6fSHong Zhang 
616d20bfe6fSHong Zhang     /* Sort the j index array for quick sparse axpy. */
617d20bfe6fSHong Zhang     ierr = PetscSortInt(apnzj,apj);CHKERRQ(ierr);
618d20bfe6fSHong Zhang 
619d20bfe6fSHong Zhang     /* Compute P^T*A*P using outer product (P^T)[:,j]*(A*P)[j,:]. */
620d20bfe6fSHong Zhang     pnzi = pi[i+1] - pi[i];
621d20bfe6fSHong Zhang     for (j=0;j<pnzi;j++) {
622d20bfe6fSHong Zhang       nextap = 0;
623d20bfe6fSHong Zhang       crow   = *pJ++;
624d20bfe6fSHong Zhang       cjj    = cj + ci[crow];
625d20bfe6fSHong Zhang       caj    = ca + ci[crow];
626d20bfe6fSHong Zhang       /* Perform sparse axpy operation.  Note cjj includes apj. */
627d20bfe6fSHong Zhang       for (k=0;nextap<apnzj;k++) {
628d20bfe6fSHong Zhang         if (cjj[k]==apj[nextap]) {
629d20bfe6fSHong Zhang           caj[k] += (*pA)*apa[apj[nextap++]];
630d20bfe6fSHong Zhang         }
631d20bfe6fSHong Zhang       }
632d20bfe6fSHong Zhang       flops += 2*apnzj;
633d20bfe6fSHong Zhang       pA++;
634d20bfe6fSHong Zhang     }
635d20bfe6fSHong Zhang 
636d20bfe6fSHong Zhang     /* Zero the current row info for A*P */
637d20bfe6fSHong Zhang     for (j=0;j<apnzj;j++) {
638d20bfe6fSHong Zhang       apa[apj[j]]      = 0.;
639d20bfe6fSHong Zhang       apjdense[apj[j]] = 0;
640d20bfe6fSHong Zhang     }
641d20bfe6fSHong Zhang   }
642d20bfe6fSHong Zhang 
643d20bfe6fSHong Zhang   /* Assemble the final matrix and clean up */
644d20bfe6fSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
645d20bfe6fSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
646d20bfe6fSHong Zhang   ierr = PetscFree(apa);CHKERRQ(ierr);
647d20bfe6fSHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
648d20bfe6fSHong Zhang 
649eb9c0419SKris Buschelman   PetscFunctionReturn(0);
650eb9c0419SKris Buschelman }
6510e36024fSHong Zhang 
6520e36024fSHong Zhang /* Compute C = P[rstart:rend,:]^T * A * P of seqaij matrices - used by MatPtAP_MPIAIJ_MPIAIJ() */
6530e36024fSHong Zhang 
6540e36024fSHong Zhang #undef __FUNCT__
6550e36024fSHong Zhang #define __FUNCT__ "MatPtAPNumeric_SeqAIJ_SeqAIJ_ReducedPt"
656b90dcfe3SHong Zhang /*@
657*e9b43d0fSSatish Balay    MatPtAPNumeric_SeqAIJ_SeqAIJ_ReducedPt - Creates C = P[rstart:rend,:]^T * A * P of seqaij matrices,
658b90dcfe3SHong Zhang                 used by MatPtAP_MPIAIJ_MPIAIJ()
659b90dcfe3SHong Zhang 
660b90dcfe3SHong Zhang    Collective on Mat
661b90dcfe3SHong Zhang 
662b90dcfe3SHong Zhang    Input Parameters:
663b90dcfe3SHong Zhang +  A - the matrix in seqaij format
664b90dcfe3SHong Zhang .  P - the projection matrix in seqaij format
665b90dcfe3SHong Zhang .  scall - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX
666b90dcfe3SHong Zhang .  fill - expected fill (not being used when scall==MAT_REUSE_MATRIX)
667b90dcfe3SHong Zhang .  prstart, prend - the starting and ending-1-th row of the matrix P to be used for transpose
668b90dcfe3SHong Zhang 
669b90dcfe3SHong Zhang    Output Parameters:
670b90dcfe3SHong Zhang .  C - the product matrix in seqaij format
671b90dcfe3SHong Zhang 
672b90dcfe3SHong Zhang    Level: developer
673b90dcfe3SHong Zhang 
674b90dcfe3SHong Zhang .seealso: MatPtAPSymbolic(),MatPtAPNumeric(),MatMatMult()
675b90dcfe3SHong Zhang @*/
6760e36024fSHong Zhang PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_ReducedPt(Mat A,Mat P,int prstart,int prend,Mat C)
6770e36024fSHong Zhang {
6780e36024fSHong Zhang   PetscErrorCode ierr;
6790e36024fSHong Zhang   int            flops=0;
6800e36024fSHong Zhang   Mat_SeqAIJ     *a  = (Mat_SeqAIJ *) A->data;
6810e36024fSHong Zhang   Mat_SeqAIJ     *p  = (Mat_SeqAIJ *) P->data;
6820e36024fSHong Zhang   Mat_SeqAIJ     *c  = (Mat_SeqAIJ *) C->data;
68320d4747cSHong Zhang   int            *ai=a->i,*aj=a->j,*apj,*apjdense;
68420d4747cSHong Zhang   int            *pi=p->i,*pj=p->j,*pJ=p->j+pi[prstart],*pjj;
6850e36024fSHong Zhang   int            *ci=c->i,*cj=c->j,*cjj;
6860e36024fSHong Zhang   int            am=A->M,cn=C->N,cm=C->M;
6870e36024fSHong Zhang   int            i,j,k,anzi,pnzi,apnzj,nextap,pnzj,prow,crow;
6880e36024fSHong Zhang   MatScalar      *aa=a->a,*apa,*pa=p->a,*pA=p->a,*paj,*ca=c->a,*caj;
6890e36024fSHong Zhang 
6900e36024fSHong Zhang   PetscFunctionBegin;
6910e36024fSHong Zhang   pA=p->a+pi[prstart];
6920e36024fSHong Zhang   /* Allocate temporary array for storage of one row of A*P */
6930e36024fSHong Zhang   ierr = PetscMalloc(cn*(sizeof(MatScalar)+2*sizeof(int)),&apa);CHKERRQ(ierr);
6940e36024fSHong Zhang   ierr = PetscMemzero(apa,cn*(sizeof(MatScalar)+2*sizeof(int)));CHKERRQ(ierr);
6950e36024fSHong Zhang 
6960e36024fSHong Zhang   apj      = (int *)(apa + cn);
6970e36024fSHong Zhang   apjdense = apj + cn;
6980e36024fSHong Zhang 
6990e36024fSHong Zhang   /* Clear old values in C */
7000e36024fSHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
7010e36024fSHong Zhang 
7020e36024fSHong Zhang   for (i=0;i<am;i++) {
7030e36024fSHong Zhang     /* Form sparse row of A*P */
7040e36024fSHong Zhang     anzi  = ai[i+1] - ai[i];
7050e36024fSHong Zhang     apnzj = 0;
7060e36024fSHong Zhang     for (j=0;j<anzi;j++) {
7070e36024fSHong Zhang       prow = *aj++;
7080e36024fSHong Zhang       pnzj = pi[prow+1] - pi[prow];
7090e36024fSHong Zhang       pjj  = pj + pi[prow];
7100e36024fSHong Zhang       paj  = pa + pi[prow];
7110e36024fSHong Zhang       for (k=0;k<pnzj;k++) {
7120e36024fSHong Zhang         if (!apjdense[pjj[k]]) {
7130e36024fSHong Zhang           apjdense[pjj[k]] = -1;
7140e36024fSHong Zhang           apj[apnzj++]     = pjj[k];
7150e36024fSHong Zhang         }
7160e36024fSHong Zhang         apa[pjj[k]] += (*aa)*paj[k];
7170e36024fSHong Zhang       }
7180e36024fSHong Zhang       flops += 2*pnzj;
7190e36024fSHong Zhang       aa++;
7200e36024fSHong Zhang     }
7210e36024fSHong Zhang 
7220e36024fSHong Zhang     /* Sort the j index array for quick sparse axpy. */
7230e36024fSHong Zhang     ierr = PetscSortInt(apnzj,apj);CHKERRQ(ierr);
7240e36024fSHong Zhang 
7250e36024fSHong Zhang     /* Compute P[[prstart:prend,:]^T*A*P using outer product (P^T)[:,j+prstart]*(A*P)[j,:]. */
7260e36024fSHong Zhang     pnzi = pi[i+1+prstart] - pi[i+prstart];
7270e36024fSHong Zhang     for (j=0;j<pnzi;j++) {
7280e36024fSHong Zhang       nextap = 0;
7290e36024fSHong Zhang       crow   = *pJ++;
7300e36024fSHong Zhang       cjj    = cj + ci[crow];
7310e36024fSHong Zhang       caj    = ca + ci[crow];
7320e36024fSHong Zhang       /* Perform sparse axpy operation.  Note cjj includes apj. */
7330e36024fSHong Zhang       for (k=0;nextap<apnzj;k++) {
7340e36024fSHong Zhang         if (cjj[k]==apj[nextap]) {
7350e36024fSHong Zhang           caj[k] += (*pA)*apa[apj[nextap++]];
7360e36024fSHong Zhang         }
7370e36024fSHong Zhang       }
7380e36024fSHong Zhang       flops += 2*apnzj;
7390e36024fSHong Zhang       pA++;
7400e36024fSHong Zhang     }
7410e36024fSHong Zhang 
7420e36024fSHong Zhang     /* Zero the current row info for A*P */
7430e36024fSHong Zhang     for (j=0;j<apnzj;j++) {
7440e36024fSHong Zhang       apa[apj[j]]      = 0.;
7450e36024fSHong Zhang       apjdense[apj[j]] = 0;
7460e36024fSHong Zhang     }
7470e36024fSHong Zhang   }
7480e36024fSHong Zhang 
7490e36024fSHong Zhang   /* Assemble the final matrix and clean up */
7500e36024fSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7510e36024fSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7520e36024fSHong Zhang   ierr = PetscFree(apa);CHKERRQ(ierr);
7530e36024fSHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
7540e36024fSHong Zhang 
7550e36024fSHong Zhang   PetscFunctionReturn(0);
7560e36024fSHong Zhang }
7570e36024fSHong Zhang 
7580e36024fSHong Zhang #undef __FUNCT__
7590e36024fSHong Zhang #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqAIJ_ReducedPt"
7600e36024fSHong Zhang PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_ReducedPt(Mat A,Mat P,PetscReal fill,int prstart,int prend,Mat *C) {
7610e36024fSHong Zhang   PetscErrorCode ierr;
7620e36024fSHong Zhang   FreeSpaceList  free_space=PETSC_NULL,current_space=PETSC_NULL;
7630e36024fSHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*p=(Mat_SeqAIJ*)P->data,*c;
7640e36024fSHong Zhang   int            *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj;
7650e36024fSHong Zhang   int            *ci,*cj,*denserow,*sparserow,*ptadenserow,*ptasparserow,*ptaj;
7660e36024fSHong Zhang   int            an=A->N,am=A->M,pn=P->N,pm=P->M;
7670e36024fSHong Zhang   int            i,j,k,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi;
7680e36024fSHong Zhang   MatScalar      *ca;
7690e36024fSHong Zhang   Mat            *psub,P_sub;
7700e36024fSHong Zhang   IS             isrow,iscol;
7710e36024fSHong Zhang   int            m = prend - prstart;
7720b89d903Svictorle 
7730b89d903Svictorle   PetscFunctionBegin;
7740b89d903Svictorle   /* Get ij structure of P[rstart:rend,:]^T */
7750e36024fSHong Zhang   ierr = ISCreateStride(PETSC_COMM_SELF,m,prstart,1,&isrow);CHKERRQ(ierr);
7760e36024fSHong Zhang   ierr = ISCreateStride(PETSC_COMM_SELF,P->n,0,1,&iscol);CHKERRQ(ierr);
7770e36024fSHong Zhang   ierr = MatGetSubMatrices(P,1,&isrow,&iscol,MAT_INITIAL_MATRIX,&psub);CHKERRQ(ierr);
7780e36024fSHong Zhang   ierr = ISDestroy(isrow);CHKERRQ(ierr);
7790e36024fSHong Zhang   ierr = ISDestroy(iscol);CHKERRQ(ierr);
7800e36024fSHong Zhang   P_sub = psub[0];
7810e36024fSHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(P_sub,&pti,&ptj);CHKERRQ(ierr);
7820e36024fSHong Zhang   ierr = MatDestroyMatrices(1,&psub);CHKERRQ(ierr);
7830e36024fSHong Zhang   ptJ=ptj;
7840e36024fSHong Zhang 
7850e36024fSHong Zhang   /* Allocate ci array, arrays for fill computation and */
7860e36024fSHong Zhang   /* free space for accumulating nonzero column info */
7870e36024fSHong Zhang   ierr = PetscMalloc((pn+1)*sizeof(int),&ci);CHKERRQ(ierr);
7880e36024fSHong Zhang   ci[0] = 0;
7890e36024fSHong Zhang 
7900e36024fSHong Zhang   ierr = PetscMalloc((2*pn+2*an+1)*sizeof(int),&ptadenserow);CHKERRQ(ierr);
7910e36024fSHong Zhang   ierr = PetscMemzero(ptadenserow,(2*pn+2*an+1)*sizeof(int));CHKERRQ(ierr);
7920e36024fSHong Zhang   ptasparserow = ptadenserow  + an;
7930e36024fSHong Zhang   denserow     = ptasparserow + an;
7940e36024fSHong Zhang   sparserow    = denserow     + pn;
7950e36024fSHong Zhang 
7960e36024fSHong Zhang   /* Set initial free space to be nnz(A) scaled by aspect ratio of P. */
7970e36024fSHong Zhang   /* This should be reasonable if sparsity of PtAP is similar to that of A. */
798b90dcfe3SHong Zhang   ierr          = GetMoreSpace((int)(fill*ai[am]/pm)*pn,&free_space);
7990e36024fSHong Zhang   current_space = free_space;
8000e36024fSHong Zhang 
8010e36024fSHong Zhang   /* Determine symbolic info for each row of C: */
8020e36024fSHong Zhang   for (i=0;i<pn;i++) {
8030e36024fSHong Zhang     ptnzi  = pti[i+1] - pti[i];
8040e36024fSHong Zhang     ptanzi = 0;
8050e36024fSHong Zhang     /* Determine symbolic row of PtA_reduced: */
8060e36024fSHong Zhang     for (j=0;j<ptnzi;j++) {
8070e36024fSHong Zhang       arow = *ptJ++;
8080e36024fSHong Zhang       anzj = ai[arow+1] - ai[arow];
8090e36024fSHong Zhang       ajj  = aj + ai[arow];
8100e36024fSHong Zhang       for (k=0;k<anzj;k++) {
8110e36024fSHong Zhang         if (!ptadenserow[ajj[k]]) {
8120e36024fSHong Zhang           ptadenserow[ajj[k]]    = -1;
8130e36024fSHong Zhang           ptasparserow[ptanzi++] = ajj[k];
8140e36024fSHong Zhang         }
8150e36024fSHong Zhang       }
8160e36024fSHong Zhang     }
8170e36024fSHong Zhang       /* Using symbolic info for row of PtA, determine symbolic info for row of C: */
8180e36024fSHong Zhang     ptaj = ptasparserow;
8190e36024fSHong Zhang     cnzi   = 0;
8200e36024fSHong Zhang     for (j=0;j<ptanzi;j++) {
8210e36024fSHong Zhang       prow = *ptaj++;
8220e36024fSHong Zhang       pnzj = pi[prow+1] - pi[prow];
8230e36024fSHong Zhang       pjj  = pj + pi[prow];
8240e36024fSHong Zhang       for (k=0;k<pnzj;k++) {
8250e36024fSHong Zhang         if (!denserow[pjj[k]]) {
8260e36024fSHong Zhang             denserow[pjj[k]]  = -1;
8270e36024fSHong Zhang             sparserow[cnzi++] = pjj[k];
8280e36024fSHong Zhang         }
8290e36024fSHong Zhang       }
8300e36024fSHong Zhang     }
8310e36024fSHong Zhang 
8320e36024fSHong Zhang     /* sort sparserow */
8330e36024fSHong Zhang     ierr = PetscSortInt(cnzi,sparserow);CHKERRQ(ierr);
8340e36024fSHong Zhang 
8350e36024fSHong Zhang     /* If free space is not available, make more free space */
8360e36024fSHong Zhang     /* Double the amount of total space in the list */
8370e36024fSHong Zhang     if (current_space->local_remaining<cnzi) {
8380e36024fSHong Zhang       ierr = GetMoreSpace(current_space->total_array_size,&current_space);CHKERRQ(ierr);
8390e36024fSHong Zhang     }
8400e36024fSHong Zhang 
8410e36024fSHong Zhang     /* Copy data into free space, and zero out denserows */
8420e36024fSHong Zhang     ierr = PetscMemcpy(current_space->array,sparserow,cnzi*sizeof(int));CHKERRQ(ierr);
8430e36024fSHong Zhang     current_space->array           += cnzi;
8440e36024fSHong Zhang     current_space->local_used      += cnzi;
8450e36024fSHong Zhang     current_space->local_remaining -= cnzi;
8460e36024fSHong Zhang 
8470e36024fSHong Zhang     for (j=0;j<ptanzi;j++) {
8480e36024fSHong Zhang       ptadenserow[ptasparserow[j]] = 0;
8490e36024fSHong Zhang     }
8500e36024fSHong Zhang     for (j=0;j<cnzi;j++) {
8510e36024fSHong Zhang       denserow[sparserow[j]] = 0;
8520e36024fSHong Zhang     }
8530e36024fSHong Zhang       /* Aside: Perhaps we should save the pta info for the numerical factorization. */
8540e36024fSHong Zhang       /*        For now, we will recompute what is needed. */
8550e36024fSHong Zhang     ci[i+1] = ci[i] + cnzi;
8560e36024fSHong Zhang   }
8570e36024fSHong Zhang   /* nnz is now stored in ci[ptm], column indices are in the list of free space */
8580e36024fSHong Zhang   /* Allocate space for cj, initialize cj, and */
8590e36024fSHong Zhang   /* destroy list of free space and other temporary array(s) */
8600e36024fSHong Zhang   ierr = PetscMalloc((ci[pn]+1)*sizeof(int),&cj);CHKERRQ(ierr);
8610e36024fSHong Zhang   ierr = MakeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
8620e36024fSHong Zhang   ierr = PetscFree(ptadenserow);CHKERRQ(ierr);
8630e36024fSHong Zhang 
8640e36024fSHong Zhang   /* Allocate space for ca */
8650e36024fSHong Zhang   ierr = PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
8660e36024fSHong Zhang   ierr = PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));CHKERRQ(ierr);
8670e36024fSHong Zhang 
8680e36024fSHong Zhang   /* put together the new matrix */
8690e36024fSHong Zhang   ierr = MatCreateSeqAIJWithArrays(A->comm,pn,pn,ci,cj,ca,C);CHKERRQ(ierr);
8700e36024fSHong Zhang 
8710e36024fSHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
8720e36024fSHong Zhang   /* Since these are PETSc arrays, change flags to free them as necessary. */
8730e36024fSHong Zhang   c = (Mat_SeqAIJ *)((*C)->data);
8740e36024fSHong Zhang   c->freedata = PETSC_TRUE;
8750e36024fSHong Zhang   c->nonew    = 0;
8760e36024fSHong Zhang 
8770e36024fSHong Zhang   /* Clean up. */
8780e36024fSHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
8790e36024fSHong Zhang 
8800e36024fSHong Zhang   PetscFunctionReturn(0);
8810e36024fSHong Zhang }
8820e36024fSHong Zhang 
8830e36024fSHong Zhang #undef __FUNCT__
8840e36024fSHong Zhang #define __FUNCT__ "MatPtAP_SeqAIJ_SeqAIJ_ReducedPt"
8850e36024fSHong Zhang PetscErrorCode MatPtAP_SeqAIJ_SeqAIJ_ReducedPt(Mat A,Mat P,MatReuse scall,PetscReal fill,int prstart,int prend,Mat *C)
8860e36024fSHong Zhang {
8870e36024fSHong Zhang   PetscErrorCode ierr;
8880e36024fSHong Zhang   PetscFunctionBegin;
8890e36024fSHong Zhang   if (A->m != prend-prstart) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",A->m,prend-prstart);
8900e36024fSHong 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);
8910e36024fSHong Zhang   if (scall == MAT_INITIAL_MATRIX){
8920e36024fSHong Zhang     ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ_ReducedPt(A,P,fill,prstart,prend,C);CHKERRQ(ierr);
8930e36024fSHong Zhang   }
8940e36024fSHong Zhang 
8950e36024fSHong Zhang   ierr = MatPtAPNumeric_SeqAIJ_SeqAIJ_ReducedPt(A,P,prstart,prend,*C);CHKERRQ(ierr);
8960e36024fSHong Zhang 
8970e36024fSHong Zhang   PetscFunctionReturn(0);
8980e36024fSHong Zhang }
8990e36024fSHong Zhang 
900