xref: /petsc/src/mat/impls/aij/seq/matptap.c (revision b1d57f1590ff1d2f11fead0759f3cc66a5edcec8)
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 @*/
36*b1d57f15SBarry Smith PetscErrorCode MatPtAP(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
37*b1d57f15SBarry Smith {
38dfbe8321SBarry Smith   PetscErrorCode ierr;
39534c1384SKris Buschelman   PetscErrorCode (*fA)(Mat,Mat,MatReuse,PetscReal,Mat *);
40534c1384SKris Buschelman   PetscErrorCode (*fP)(Mat,Mat,MatReuse,PetscReal,Mat *);
41eb9c0419SKris Buschelman 
42eb9c0419SKris Buschelman   PetscFunctionBegin;
439af31e4aSHong Zhang   PetscValidHeaderSpecific(A,MAT_COOKIE,1);
449af31e4aSHong Zhang   PetscValidType(A,1);
459af31e4aSHong Zhang   MatPreallocated(A);
469af31e4aSHong Zhang   if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
479af31e4aSHong Zhang   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
489af31e4aSHong Zhang   PetscValidHeaderSpecific(P,MAT_COOKIE,2);
499af31e4aSHong Zhang   PetscValidType(P,2);
509af31e4aSHong Zhang   MatPreallocated(P);
519af31e4aSHong Zhang   if (!P->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
529af31e4aSHong Zhang   if (P->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
539af31e4aSHong Zhang   PetscValidPointer(C,3);
54534c1384SKris Buschelman 
559af31e4aSHong Zhang   if (P->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->M,A->N);
56eb9c0419SKris Buschelman 
579af31e4aSHong Zhang   if (fill <=0.0) SETERRQ1(PETSC_ERR_ARG_SIZ,"fill=%g must be > 0.0",fill);
58eb9c0419SKris Buschelman 
59534c1384SKris Buschelman   /* For now, we do not dispatch based on the type of A and P */
60534c1384SKris Buschelman   /* When implementations like _SeqAIJ_MAIJ exist, attack the multiple dispatch problem. */
61534c1384SKris Buschelman   fA = A->ops->ptap;
62534c1384SKris Buschelman   if (!fA) SETERRQ1(PETSC_ERR_SUP,"MatPtAP not supported for A of type %s",A->type_name);
63534c1384SKris Buschelman   fP = P->ops->ptap;
64534c1384SKris Buschelman   if (!fP) SETERRQ1(PETSC_ERR_SUP,"MatPtAP not supported for P of type %s",P->type_name);
65534c1384SKris 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);
66534c1384SKris Buschelman 
679af31e4aSHong Zhang   ierr = PetscLogEventBegin(MAT_PtAP,A,P,0,0);CHKERRQ(ierr);
68534c1384SKris Buschelman   ierr = (*fA)(A,P,scall,fill,C);CHKERRQ(ierr);
699af31e4aSHong Zhang   ierr = PetscLogEventEnd(MAT_PtAP,A,P,0,0);CHKERRQ(ierr);
70eb9c0419SKris Buschelman   PetscFunctionReturn(0);
71eb9c0419SKris Buschelman }
72eb9c0419SKris Buschelman 
73*b1d57f15SBarry Smith EXTERN PetscErrorCode MatPtAP_SeqAIJ_SeqAIJ_ReducedPt(Mat,Mat,MatReuse,PetscReal,PetscInt,PetscInt,Mat*);
74*b1d57f15SBarry Smith EXTERN PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_ReducedPt(Mat,Mat,PetscReal,PetscInt,PetscInt,Mat*);
75*b1d57f15SBarry Smith EXTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_ReducedPt(Mat,Mat,PetscInt,PetscInt,Mat);
76b90dcfe3SHong Zhang 
77eb9c0419SKris Buschelman #undef __FUNCT__
78ff134f7aSHong Zhang #define __FUNCT__ "MatPtAP_MPIAIJ_MPIAIJ"
79ff134f7aSHong Zhang PetscErrorCode MatPtAP_MPIAIJ_MPIAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
80ff134f7aSHong Zhang {
81ff134f7aSHong Zhang   PetscErrorCode    ierr;
82b90dcfe3SHong Zhang 
83b90dcfe3SHong Zhang   PetscFunctionBegin;
84b90dcfe3SHong Zhang   if (scall == MAT_INITIAL_MATRIX){
85b90dcfe3SHong Zhang     ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ(A,P,fill,C);CHKERRQ(ierr);/* numeric product is computed as well */
86b90dcfe3SHong Zhang   } else if (scall == MAT_REUSE_MATRIX){
87b90dcfe3SHong Zhang     ierr = MatPtAPNumeric_MPIAIJ_MPIAIJ(A,P,*C);CHKERRQ(ierr);
88b90dcfe3SHong Zhang   } else {
89b90dcfe3SHong Zhang     SETERRQ1(PETSC_ERR_ARG_WRONG,"Invalid MatReuse %d",scall);
90b90dcfe3SHong Zhang   }
91b90dcfe3SHong Zhang   PetscFunctionReturn(0);
92b90dcfe3SHong Zhang }
93b90dcfe3SHong Zhang 
94b90dcfe3SHong Zhang #undef __FUNCT__
95b90dcfe3SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ"
96b90dcfe3SHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
97b90dcfe3SHong Zhang {
98b90dcfe3SHong Zhang   PetscErrorCode    ierr;
9925616d81SHong Zhang   Mat               P_seq,A_loc,C_seq;
100*b1d57f15SBarry Smith   PetscInt          prstart,prend,m=P->m;
10125616d81SHong Zhang   IS                isrowp,iscolp;
102ff134f7aSHong Zhang 
103ff134f7aSHong Zhang   PetscFunctionBegin;
10425616d81SHong Zhang   /* get P_seq = submatrix of P by taking rows of P that equal to nonzero col of A */
105b90dcfe3SHong Zhang   ierr = MatGetBrowsOfAcols(A,P,MAT_INITIAL_MATRIX,&isrowp,&iscolp,&prstart,&P_seq);CHKERRQ(ierr);
10625616d81SHong Zhang   ierr = ISDestroy(iscolp);CHKERRQ(ierr);
107ff134f7aSHong Zhang 
10825616d81SHong Zhang   /* get A_loc = submatrix of A by taking all local rows of A */
109b90dcfe3SHong Zhang   ierr = MatGetLocalMat(A,MAT_INITIAL_MATRIX,PETSC_NULL,&isrowp,&A_loc);CHKERRQ(ierr);
11025616d81SHong Zhang   ierr = ISDestroy(isrowp);CHKERRQ(ierr);
1110e36024fSHong Zhang 
11225616d81SHong Zhang   /* compute C_seq = P_loc^T * A_loc * P_seq */
113ff134f7aSHong Zhang   prend   = prstart + m;
114b90dcfe3SHong Zhang   ierr = MatPtAP_SeqAIJ_SeqAIJ_ReducedPt(A_loc,P_seq,MAT_INITIAL_MATRIX,fill,prstart,prend,&C_seq);CHKERRQ(ierr);
11525616d81SHong Zhang   ierr = MatDestroy(P_seq);CHKERRQ(ierr);
11625616d81SHong Zhang   ierr = MatDestroy(A_loc);CHKERRQ(ierr);
117b90dcfe3SHong Zhang 
118b90dcfe3SHong Zhang   /* add C_seq into mpi C */
119b90dcfe3SHong Zhang   ierr = MatMerge_SeqsToMPI(A->comm,C_seq,P->n,P->n,MAT_INITIAL_MATRIX,C);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;
130*b1d57f15SBarry Smith   PetscInt             prstart,prend,m=P->m;
131b90dcfe3SHong Zhang   IS                   isrowp,iscolp;
132671beff6SHong Zhang   Mat_Merge_SeqsToMPI  *merge;
133671beff6SHong Zhang   PetscObjectContainer container;
134ff134f7aSHong Zhang 
135ff134f7aSHong Zhang   PetscFunctionBegin;
136671beff6SHong Zhang   ierr = PetscObjectQuery((PetscObject)C,"MatMergeSeqsToMPI",(PetscObject *)&container);CHKERRQ(ierr);
137671beff6SHong Zhang   if (container) {
1387f79fd58SMatthew Knepley     ierr  = PetscObjectContainerGetPointer(container,(void **)&merge);CHKERRQ(ierr);
1397f79fd58SMatthew Knepley   } else {
1407f79fd58SMatthew Knepley     SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Matrix C does not posses and object container");
141671beff6SHong Zhang   }
142671beff6SHong Zhang 
143b90dcfe3SHong Zhang   /* get P_seq = submatrix of P by taking rows of P that equal to nonzero col of A */
144b90dcfe3SHong Zhang   ierr = MatGetBrowsOfAcols(A,P,MAT_INITIAL_MATRIX,&isrowp,&iscolp,&prstart,&P_seq);CHKERRQ(ierr);
145b90dcfe3SHong Zhang   ierr = ISDestroy(iscolp);CHKERRQ(ierr);
146ff134f7aSHong Zhang 
147b90dcfe3SHong Zhang   /* get A_loc = submatrix of A by taking all local rows of A */
148b90dcfe3SHong Zhang   ierr = MatGetLocalMat(A,MAT_INITIAL_MATRIX,PETSC_NULL,&isrowp,&A_loc);CHKERRQ(ierr);
149b90dcfe3SHong Zhang   ierr = ISDestroy(isrowp);CHKERRQ(ierr);
150ff134f7aSHong Zhang 
151b90dcfe3SHong Zhang   /* compute C_seq = P_loc^T * A_loc * P_seq */
152b90dcfe3SHong Zhang   prend = prstart + m;
153b90dcfe3SHong Zhang   C_seq = merge->C_seq;
154b90dcfe3SHong Zhang   ierr = MatPtAP_SeqAIJ_SeqAIJ_ReducedPt(A_loc,P_seq,MAT_REUSE_MATRIX,1.0,prstart,prend,&C_seq);CHKERRQ(ierr);
155b90dcfe3SHong Zhang   ierr = MatDestroy(P_seq);CHKERRQ(ierr);
156b90dcfe3SHong Zhang   ierr = MatDestroy(A_loc);CHKERRQ(ierr);
157b90dcfe3SHong Zhang 
158b90dcfe3SHong Zhang   /* add C_seq into mpi C */
159b90dcfe3SHong Zhang   ierr = MatMerge_SeqsToMPI(A->comm,C_seq,P->n,P->n,MAT_REUSE_MATRIX,&C);CHKERRQ(ierr);
160b90dcfe3SHong Zhang 
161ff134f7aSHong Zhang   PetscFunctionReturn(0);
162ff134f7aSHong Zhang }
163ff134f7aSHong Zhang 
164ff134f7aSHong Zhang #undef __FUNCT__
1659af31e4aSHong Zhang #define __FUNCT__ "MatPtAP_SeqAIJ_SeqAIJ"
166dfbe8321SBarry Smith PetscErrorCode MatPtAP_SeqAIJ_SeqAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
1679af31e4aSHong Zhang {
168dfbe8321SBarry Smith   PetscErrorCode ierr;
169*b1d57f15SBarry Smith 
1709af31e4aSHong Zhang   PetscFunctionBegin;
1719af31e4aSHong Zhang   if (scall == MAT_INITIAL_MATRIX){
172d20bfe6fSHong Zhang     ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
1739af31e4aSHong Zhang     ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ(A,P,fill,C);CHKERRQ(ierr);
174d20bfe6fSHong Zhang     ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
1759af31e4aSHong Zhang   }
176d20bfe6fSHong Zhang   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
1779af31e4aSHong Zhang   ierr = MatPtAPNumeric_SeqAIJ_SeqAIJ(A,P,*C);CHKERRQ(ierr);
178d20bfe6fSHong Zhang   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
1799af31e4aSHong Zhang   PetscFunctionReturn(0);
1809af31e4aSHong Zhang }
1819af31e4aSHong Zhang 
1829af31e4aSHong Zhang #undef __FUNCT__
1839af31e4aSHong Zhang #define __FUNCT__ "MatPtAPSymbolic"
1846849ba73SBarry Smith /*
1859af31e4aSHong Zhang    MatPtAPSymbolic - Creates the (i,j) structure of the matrix projection C = P^T * A * P
1864d3841fdSKris Buschelman 
1874d3841fdSKris Buschelman    Collective on Mat
1884d3841fdSKris Buschelman 
1894d3841fdSKris Buschelman    Input Parameters:
1904d3841fdSKris Buschelman +  A - the matrix
1914d3841fdSKris Buschelman -  P - the projection matrix
1924d3841fdSKris Buschelman 
1934d3841fdSKris Buschelman    Output Parameters:
1944d3841fdSKris Buschelman .  C - the (i,j) structure of the product matrix
1954d3841fdSKris Buschelman 
1964d3841fdSKris Buschelman    Notes:
1974d3841fdSKris Buschelman    C will be created and must be destroyed by the user with MatDestroy().
1984d3841fdSKris Buschelman 
1994d3841fdSKris Buschelman    This routine is currently only implemented for pairs of SeqAIJ matrices and classes
2004d3841fdSKris Buschelman    which inherit from SeqAIJ.  C will be of type MATSEQAIJ.  The product is computed using
2019af31e4aSHong Zhang    this (i,j) structure by calling MatPtAPNumeric().
2024d3841fdSKris Buschelman 
2034d3841fdSKris Buschelman    Level: intermediate
2044d3841fdSKris Buschelman 
2059af31e4aSHong Zhang .seealso: MatPtAP(),MatPtAPNumeric(),MatMatMultSymbolic()
2066849ba73SBarry Smith */
207*b1d57f15SBarry Smith PetscErrorCode MatPtAPSymbolic(Mat A,Mat P,PetscReal fill,Mat *C)
208*b1d57f15SBarry Smith {
209dfbe8321SBarry Smith   PetscErrorCode ierr;
210534c1384SKris Buschelman   PetscErrorCode (*fA)(Mat,Mat,PetscReal,Mat*);
211534c1384SKris Buschelman   PetscErrorCode (*fP)(Mat,Mat,PetscReal,Mat*);
212eb9c0419SKris Buschelman 
213eb9c0419SKris Buschelman   PetscFunctionBegin;
214eb9c0419SKris Buschelman 
2154482741eSBarry Smith   PetscValidHeaderSpecific(A,MAT_COOKIE,1);
216c9780b6fSBarry Smith   PetscValidType(A,1);
217eb9c0419SKris Buschelman   MatPreallocated(A);
218eb9c0419SKris Buschelman   if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
219eb9c0419SKris Buschelman   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
220eb9c0419SKris Buschelman 
2214482741eSBarry Smith   PetscValidHeaderSpecific(P,MAT_COOKIE,2);
222c9780b6fSBarry Smith   PetscValidType(P,2);
223eb9c0419SKris Buschelman   MatPreallocated(P);
224eb9c0419SKris Buschelman   if (!P->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
225eb9c0419SKris Buschelman   if (P->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
226eb9c0419SKris Buschelman 
2274482741eSBarry Smith   PetscValidPointer(C,3);
2284482741eSBarry Smith 
229eb9c0419SKris Buschelman   if (P->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->M,A->N);
230eb9c0419SKris Buschelman   if (A->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix 'A' must be square, %d != %d",A->M,A->N);
231eb9c0419SKris Buschelman 
232534c1384SKris Buschelman   /* For now, we do not dispatch based on the type of A and P */
233534c1384SKris Buschelman   /* When implementations like _SeqAIJ_MAIJ exist, attack the multiple dispatch problem. */
234534c1384SKris Buschelman   fA = A->ops->ptapsymbolic;
235534c1384SKris Buschelman   if (!fA) SETERRQ1(PETSC_ERR_SUP,"MatPtAPSymbolic not supported for A of type %s",A->type_name);
236534c1384SKris Buschelman   fP = P->ops->ptapsymbolic;
237534c1384SKris Buschelman   if (!fP) SETERRQ1(PETSC_ERR_SUP,"MatPtAPSymbolic not supported for P of type %s",P->type_name);
238534c1384SKris 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);
2394d3841fdSKris Buschelman 
240534c1384SKris Buschelman   ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
241534c1384SKris Buschelman   ierr = (*fA)(A,P,fill,C);CHKERRQ(ierr);
242534c1384SKris Buschelman   ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
243eb9c0419SKris Buschelman 
244eb9c0419SKris Buschelman   PetscFunctionReturn(0);
245eb9c0419SKris Buschelman }
246eb9c0419SKris Buschelman 
247f747e1aeSHong Zhang typedef struct {
248f747e1aeSHong Zhang   Mat    symAP;
249f747e1aeSHong Zhang } Mat_PtAPstruct;
250f747e1aeSHong Zhang 
25178a80504SBarry Smith EXTERN PetscErrorCode MatDestroy_SeqAIJ(Mat);
25278a80504SBarry Smith 
253f747e1aeSHong Zhang #undef __FUNCT__
254f747e1aeSHong Zhang #define __FUNCT__ "MatDestroy_SeqAIJ_PtAP"
255f4a850bbSBarry Smith PetscErrorCode MatDestroy_SeqAIJ_PtAP(Mat A)
256f747e1aeSHong Zhang {
257f4a850bbSBarry Smith   PetscErrorCode    ierr;
258f747e1aeSHong Zhang   Mat_PtAPstruct    *ptap=(Mat_PtAPstruct*)A->spptr;
259f747e1aeSHong Zhang 
260f747e1aeSHong Zhang   PetscFunctionBegin;
261f747e1aeSHong Zhang   ierr = MatDestroy(ptap->symAP);CHKERRQ(ierr);
262f747e1aeSHong Zhang   ierr = PetscFree(ptap);CHKERRQ(ierr);
26378a80504SBarry Smith   ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr);
264f747e1aeSHong Zhang   PetscFunctionReturn(0);
265f747e1aeSHong Zhang }
266f747e1aeSHong Zhang 
267eb9c0419SKris Buschelman #undef __FUNCT__
2689af31e4aSHong Zhang #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqAIJ"
269*b1d57f15SBarry Smith PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
270*b1d57f15SBarry Smith {
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;
274*b1d57f15SBarry Smith   PetscInt       *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj;
275*b1d57f15SBarry Smith   PetscInt       *ci,*cj,*denserow,*sparserow,*ptadenserow,*ptasparserow,*ptaj;
276*b1d57f15SBarry Smith   PetscInt       an=A->N,am=A->M,pn=P->N,pm=P->M;
277*b1d57f15SBarry Smith   PetscInt       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 */
287*b1d57f15SBarry Smith   ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
288d20bfe6fSHong Zhang   ci[0] = 0;
289eb9c0419SKris Buschelman 
290*b1d57f15SBarry Smith   ierr = PetscMalloc((2*pn+2*an+1)*sizeof(PetscInt),&ptadenserow);CHKERRQ(ierr);
291*b1d57f15SBarry Smith   ierr = PetscMemzero(ptadenserow,(2*pn+2*an+1)*sizeof(PetscInt));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 */
342*b1d57f15SBarry Smith     ierr = PetscMemcpy(current_space->array,sparserow,cnzi*sizeof(PetscInt));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) */
360*b1d57f15SBarry Smith   ierr = PetscMalloc((ci[pn]+1)*sizeof(PetscInt),&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"
387*b1d57f15SBarry Smith PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqMAIJ(Mat A,Mat PP,Mat *C)
388*b1d57f15SBarry Smith {
3895c66b693SKris Buschelman   /* This routine requires testing -- I don't think it works. */
390dfbe8321SBarry Smith   PetscErrorCode ierr;
3913985e5eaSKris Buschelman   FreeSpaceList  free_space=PETSC_NULL,current_space=PETSC_NULL;
3923985e5eaSKris Buschelman   Mat_SeqMAIJ    *pp=(Mat_SeqMAIJ*)PP->data;
3933985e5eaSKris Buschelman   Mat            P=pp->AIJ;
3943985e5eaSKris Buschelman   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*p=(Mat_SeqAIJ*)P->data,*c;
395*b1d57f15SBarry Smith   PetscInt       *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj;
396*b1d57f15SBarry Smith   PetscInt       *ci,*cj,*denserow,*sparserow,*ptadenserow,*ptasparserow,*ptaj;
397*b1d57f15SBarry Smith   PetscInt       an=A->N,am=A->M,pn=P->N,pm=P->M,ppdof=pp->dof;
398*b1d57f15SBarry Smith   PetscInt       i,j,k,dof,pdof,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi;
3993985e5eaSKris Buschelman   MatScalar      *ca;
4003985e5eaSKris Buschelman 
4013985e5eaSKris Buschelman   PetscFunctionBegin;
4023985e5eaSKris Buschelman   /* Start timer */
4039af31e4aSHong Zhang   ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,PP,0,0);CHKERRQ(ierr);
4043985e5eaSKris Buschelman 
4053985e5eaSKris Buschelman   /* Get ij structure of P^T */
4063985e5eaSKris Buschelman   ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
4073985e5eaSKris Buschelman 
4083985e5eaSKris Buschelman   /* Allocate ci array, arrays for fill computation and */
4093985e5eaSKris Buschelman   /* free space for accumulating nonzero column info */
410*b1d57f15SBarry Smith   ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
4113985e5eaSKris Buschelman   ci[0] = 0;
4123985e5eaSKris Buschelman 
413*b1d57f15SBarry Smith   ierr = PetscMalloc((2*pn+2*an+1)*sizeof(PetscInt),&ptadenserow);CHKERRQ(ierr);
414*b1d57f15SBarry Smith   ierr = PetscMemzero(ptadenserow,(2*pn+2*an+1)*sizeof(PetscInt));CHKERRQ(ierr);
4153985e5eaSKris Buschelman   ptasparserow = ptadenserow  + an;
4163985e5eaSKris Buschelman   denserow     = ptasparserow + an;
4173985e5eaSKris Buschelman   sparserow    = denserow     + pn;
4183985e5eaSKris Buschelman 
4193985e5eaSKris Buschelman   /* Set initial free space to be nnz(A) scaled by aspect ratio of P. */
4203985e5eaSKris Buschelman   /* This should be reasonable if sparsity of PtAP is similar to that of A. */
4213985e5eaSKris Buschelman   ierr          = GetMoreSpace((ai[am]/pm)*pn,&free_space);
4223985e5eaSKris Buschelman   current_space = free_space;
4233985e5eaSKris Buschelman 
4243985e5eaSKris Buschelman   /* Determine symbolic info for each row of C: */
4253985e5eaSKris Buschelman   for (i=0;i<pn/ppdof;i++) {
4263985e5eaSKris Buschelman     ptnzi  = pti[i+1] - pti[i];
4273985e5eaSKris Buschelman     ptanzi = 0;
4283985e5eaSKris Buschelman     ptJ    = ptj + pti[i];
4293985e5eaSKris Buschelman     for (dof=0;dof<ppdof;dof++) {
4303985e5eaSKris Buschelman     /* Determine symbolic row of PtA: */
4313985e5eaSKris Buschelman       for (j=0;j<ptnzi;j++) {
4323985e5eaSKris Buschelman         arow = ptJ[j] + dof;
4333985e5eaSKris Buschelman         anzj = ai[arow+1] - ai[arow];
4343985e5eaSKris Buschelman         ajj  = aj + ai[arow];
4353985e5eaSKris Buschelman         for (k=0;k<anzj;k++) {
4363985e5eaSKris Buschelman           if (!ptadenserow[ajj[k]]) {
4373985e5eaSKris Buschelman             ptadenserow[ajj[k]]    = -1;
4383985e5eaSKris Buschelman             ptasparserow[ptanzi++] = ajj[k];
4393985e5eaSKris Buschelman           }
4403985e5eaSKris Buschelman         }
4413985e5eaSKris Buschelman       }
4423985e5eaSKris Buschelman       /* Using symbolic info for row of PtA, determine symbolic info for row of C: */
4433985e5eaSKris Buschelman       ptaj = ptasparserow;
4443985e5eaSKris Buschelman       cnzi   = 0;
4453985e5eaSKris Buschelman       for (j=0;j<ptanzi;j++) {
446fe05a634SKris Buschelman         pdof = *ptaj%dof;
4473985e5eaSKris Buschelman         prow = (*ptaj++)/dof;
4483985e5eaSKris Buschelman         pnzj = pi[prow+1] - pi[prow];
4493985e5eaSKris Buschelman         pjj  = pj + pi[prow];
4503985e5eaSKris Buschelman         for (k=0;k<pnzj;k++) {
451fe05a634SKris Buschelman           if (!denserow[pjj[k]+pdof]) {
452fe05a634SKris Buschelman             denserow[pjj[k]+pdof] = -1;
453fe05a634SKris Buschelman             sparserow[cnzi++]     = pjj[k]+pdof;
4543985e5eaSKris Buschelman           }
4553985e5eaSKris Buschelman         }
4563985e5eaSKris Buschelman       }
4573985e5eaSKris Buschelman 
4583985e5eaSKris Buschelman       /* sort sparserow */
4593985e5eaSKris Buschelman       ierr = PetscSortInt(cnzi,sparserow);CHKERRQ(ierr);
4603985e5eaSKris Buschelman 
4613985e5eaSKris Buschelman       /* If free space is not available, make more free space */
4623985e5eaSKris Buschelman       /* Double the amount of total space in the list */
4633985e5eaSKris Buschelman       if (current_space->local_remaining<cnzi) {
4643985e5eaSKris Buschelman         ierr = GetMoreSpace(current_space->total_array_size,&current_space);CHKERRQ(ierr);
4653985e5eaSKris Buschelman       }
4663985e5eaSKris Buschelman 
4673985e5eaSKris Buschelman       /* Copy data into free space, and zero out denserows */
468*b1d57f15SBarry Smith       ierr = PetscMemcpy(current_space->array,sparserow,cnzi*sizeof(PetscInt));CHKERRQ(ierr);
4693985e5eaSKris Buschelman       current_space->array           += cnzi;
4703985e5eaSKris Buschelman       current_space->local_used      += cnzi;
4713985e5eaSKris Buschelman       current_space->local_remaining -= cnzi;
4723985e5eaSKris Buschelman 
4733985e5eaSKris Buschelman       for (j=0;j<ptanzi;j++) {
4743985e5eaSKris Buschelman         ptadenserow[ptasparserow[j]] = 0;
4753985e5eaSKris Buschelman       }
4763985e5eaSKris Buschelman       for (j=0;j<cnzi;j++) {
4773985e5eaSKris Buschelman         denserow[sparserow[j]] = 0;
4783985e5eaSKris Buschelman       }
4793985e5eaSKris Buschelman       /* Aside: Perhaps we should save the pta info for the numerical factorization. */
4803985e5eaSKris Buschelman       /*        For now, we will recompute what is needed. */
4813985e5eaSKris Buschelman       ci[i+1+dof] = ci[i+dof] + cnzi;
4823985e5eaSKris Buschelman     }
4833985e5eaSKris Buschelman   }
4843985e5eaSKris Buschelman   /* nnz is now stored in ci[ptm], column indices are in the list of free space */
4853985e5eaSKris Buschelman   /* Allocate space for cj, initialize cj, and */
4863985e5eaSKris Buschelman   /* destroy list of free space and other temporary array(s) */
487*b1d57f15SBarry Smith   ierr = PetscMalloc((ci[pn]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
4883985e5eaSKris Buschelman   ierr = MakeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
4893985e5eaSKris Buschelman   ierr = PetscFree(ptadenserow);CHKERRQ(ierr);
4903985e5eaSKris Buschelman 
4913985e5eaSKris Buschelman   /* Allocate space for ca */
4923985e5eaSKris Buschelman   ierr = PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
4933985e5eaSKris Buschelman   ierr = PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));CHKERRQ(ierr);
4943985e5eaSKris Buschelman 
4953985e5eaSKris Buschelman   /* put together the new matrix */
4963985e5eaSKris Buschelman   ierr = MatCreateSeqAIJWithArrays(A->comm,pn,pn,ci,cj,ca,C);CHKERRQ(ierr);
4973985e5eaSKris Buschelman 
4983985e5eaSKris Buschelman   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
4993985e5eaSKris Buschelman   /* Since these are PETSc arrays, change flags to free them as necessary. */
5003985e5eaSKris Buschelman   c = (Mat_SeqAIJ *)((*C)->data);
5013985e5eaSKris Buschelman   c->freedata = PETSC_TRUE;
5023985e5eaSKris Buschelman   c->nonew    = 0;
5033985e5eaSKris Buschelman 
5043985e5eaSKris Buschelman   /* Clean up. */
5053985e5eaSKris Buschelman   ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
5063985e5eaSKris Buschelman 
5079af31e4aSHong Zhang   ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,PP,0,0);CHKERRQ(ierr);
5083985e5eaSKris Buschelman   PetscFunctionReturn(0);
5093985e5eaSKris Buschelman }
5103985e5eaSKris Buschelman EXTERN_C_END
5113985e5eaSKris Buschelman 
512eb9c0419SKris Buschelman #undef __FUNCT__
5139af31e4aSHong Zhang #define __FUNCT__ "MatPtAPNumeric"
5146849ba73SBarry Smith /*
5159af31e4aSHong Zhang    MatPtAPNumeric - Computes the matrix projection C = P^T * A * P
5164d3841fdSKris Buschelman 
5174d3841fdSKris Buschelman    Collective on Mat
5184d3841fdSKris Buschelman 
5194d3841fdSKris Buschelman    Input Parameters:
5204d3841fdSKris Buschelman +  A - the matrix
5214d3841fdSKris Buschelman -  P - the projection matrix
5224d3841fdSKris Buschelman 
5234d3841fdSKris Buschelman    Output Parameters:
5244d3841fdSKris Buschelman .  C - the product matrix
5254d3841fdSKris Buschelman 
5264d3841fdSKris Buschelman    Notes:
5279af31e4aSHong Zhang    C must have been created by calling MatPtAPSymbolic and must be destroyed by
5284d3841fdSKris Buschelman    the user using MatDeatroy().
5294d3841fdSKris Buschelman 
530170ef064SHong Zhang    This routine is currently only implemented for pairs of AIJ matrices and classes
531170ef064SHong Zhang    which inherit from AIJ.  C will be of type MATAIJ.
5324d3841fdSKris Buschelman 
5334d3841fdSKris Buschelman    Level: intermediate
5344d3841fdSKris Buschelman 
5359af31e4aSHong Zhang .seealso: MatPtAP(),MatPtAPSymbolic(),MatMatMultNumeric()
5366849ba73SBarry Smith */
537*b1d57f15SBarry Smith PetscErrorCode MatPtAPNumeric(Mat A,Mat P,Mat C)
538*b1d57f15SBarry Smith {
539dfbe8321SBarry Smith   PetscErrorCode ierr;
540534c1384SKris Buschelman   PetscErrorCode (*fA)(Mat,Mat,Mat);
541534c1384SKris Buschelman   PetscErrorCode (*fP)(Mat,Mat,Mat);
542eb9c0419SKris Buschelman 
543eb9c0419SKris Buschelman   PetscFunctionBegin;
544eb9c0419SKris Buschelman 
5454482741eSBarry Smith   PetscValidHeaderSpecific(A,MAT_COOKIE,1);
546c9780b6fSBarry Smith   PetscValidType(A,1);
547eb9c0419SKris Buschelman   MatPreallocated(A);
548eb9c0419SKris Buschelman   if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
549eb9c0419SKris Buschelman   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
550eb9c0419SKris Buschelman 
5514482741eSBarry Smith   PetscValidHeaderSpecific(P,MAT_COOKIE,2);
552c9780b6fSBarry Smith   PetscValidType(P,2);
553eb9c0419SKris Buschelman   MatPreallocated(P);
554eb9c0419SKris Buschelman   if (!P->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
555eb9c0419SKris Buschelman   if (P->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
556eb9c0419SKris Buschelman 
5574482741eSBarry Smith   PetscValidHeaderSpecific(C,MAT_COOKIE,3);
558c9780b6fSBarry Smith   PetscValidType(C,3);
559eb9c0419SKris Buschelman   MatPreallocated(C);
560eb9c0419SKris Buschelman   if (!C->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
561eb9c0419SKris Buschelman   if (C->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
562eb9c0419SKris Buschelman 
563eb9c0419SKris Buschelman   if (P->N!=C->M) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->N,C->M);
564eb9c0419SKris Buschelman   if (P->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->M,A->N);
565eb9c0419SKris Buschelman   if (A->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix 'A' must be square, %d != %d",A->M,A->N);
566eb9c0419SKris Buschelman   if (P->N!=C->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->N,C->N);
567eb9c0419SKris Buschelman 
568534c1384SKris Buschelman   /* For now, we do not dispatch based on the type of A and P */
569534c1384SKris Buschelman   /* When implementations like _SeqAIJ_MAIJ exist, attack the multiple dispatch problem. */
570534c1384SKris Buschelman   fA = A->ops->ptapnumeric;
571534c1384SKris Buschelman   if (!fA) SETERRQ1(PETSC_ERR_SUP,"MatPtAPNumeric not supported for A of type %s",A->type_name);
572534c1384SKris Buschelman   fP = P->ops->ptapnumeric;
573534c1384SKris Buschelman   if (!fP) SETERRQ1(PETSC_ERR_SUP,"MatPtAPNumeric not supported for P of type %s",P->type_name);
574534c1384SKris 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);
5754d3841fdSKris Buschelman 
576534c1384SKris Buschelman   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
577534c1384SKris Buschelman   ierr = (*fA)(A,P,C);CHKERRQ(ierr);
578534c1384SKris Buschelman   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
579eb9c0419SKris Buschelman 
580eb9c0419SKris Buschelman   PetscFunctionReturn(0);
581eb9c0419SKris Buschelman }
582eb9c0419SKris Buschelman 
583eb9c0419SKris Buschelman #undef __FUNCT__
5849af31e4aSHong Zhang #define __FUNCT__ "MatPtAPNumeric_SeqAIJ_SeqAIJ"
585dfbe8321SBarry Smith PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat A,Mat P,Mat C)
586dfbe8321SBarry Smith {
587dfbe8321SBarry Smith   PetscErrorCode ierr;
588*b1d57f15SBarry Smith   PetscInt       flops=0;
589d20bfe6fSHong Zhang   Mat_SeqAIJ     *a  = (Mat_SeqAIJ *) A->data;
590d20bfe6fSHong Zhang   Mat_SeqAIJ     *p  = (Mat_SeqAIJ *) P->data;
591d20bfe6fSHong Zhang   Mat_SeqAIJ     *c  = (Mat_SeqAIJ *) C->data;
592*b1d57f15SBarry Smith   PetscInt       *ai=a->i,*aj=a->j,*apj,*apjdense,*pi=p->i,*pj=p->j,*pJ=p->j,*pjj;
593*b1d57f15SBarry Smith   PetscInt       *ci=c->i,*cj=c->j,*cjj;
594*b1d57f15SBarry Smith   PetscInt       am=A->M,cn=C->N,cm=C->M;
595*b1d57f15SBarry Smith   PetscInt       i,j,k,anzi,pnzi,apnzj,nextap,pnzj,prow,crow;
596d20bfe6fSHong Zhang   MatScalar      *aa=a->a,*apa,*pa=p->a,*pA=p->a,*paj,*ca=c->a,*caj;
597eb9c0419SKris Buschelman 
598eb9c0419SKris Buschelman   PetscFunctionBegin;
599d20bfe6fSHong Zhang   /* Allocate temporary array for storage of one row of A*P */
600*b1d57f15SBarry Smith   ierr = PetscMalloc(cn*(sizeof(MatScalar)+2*sizeof(PetscInt)),&apa);CHKERRQ(ierr);
601*b1d57f15SBarry Smith   ierr = PetscMemzero(apa,cn*(sizeof(MatScalar)+2*sizeof(PetscInt)));CHKERRQ(ierr);
602eb9c0419SKris Buschelman 
603*b1d57f15SBarry Smith   apj      = (PetscInt *)(apa + cn);
604d20bfe6fSHong Zhang   apjdense = apj + cn;
605d20bfe6fSHong Zhang 
606d20bfe6fSHong Zhang   /* Clear old values in C */
607d20bfe6fSHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
608d20bfe6fSHong Zhang 
609d20bfe6fSHong Zhang   for (i=0;i<am;i++) {
610d20bfe6fSHong Zhang     /* Form sparse row of A*P */
611d20bfe6fSHong Zhang     anzi  = ai[i+1] - ai[i];
612d20bfe6fSHong Zhang     apnzj = 0;
613d20bfe6fSHong Zhang     for (j=0;j<anzi;j++) {
614d20bfe6fSHong Zhang       prow = *aj++;
615d20bfe6fSHong Zhang       pnzj = pi[prow+1] - pi[prow];
616d20bfe6fSHong Zhang       pjj  = pj + pi[prow];
617d20bfe6fSHong Zhang       paj  = pa + pi[prow];
618d20bfe6fSHong Zhang       for (k=0;k<pnzj;k++) {
619d20bfe6fSHong Zhang         if (!apjdense[pjj[k]]) {
620d20bfe6fSHong Zhang           apjdense[pjj[k]] = -1;
621d20bfe6fSHong Zhang           apj[apnzj++]     = pjj[k];
622d20bfe6fSHong Zhang         }
623d20bfe6fSHong Zhang         apa[pjj[k]] += (*aa)*paj[k];
624d20bfe6fSHong Zhang       }
625d20bfe6fSHong Zhang       flops += 2*pnzj;
626d20bfe6fSHong Zhang       aa++;
627d20bfe6fSHong Zhang     }
628d20bfe6fSHong Zhang 
629d20bfe6fSHong Zhang     /* Sort the j index array for quick sparse axpy. */
630d20bfe6fSHong Zhang     ierr = PetscSortInt(apnzj,apj);CHKERRQ(ierr);
631d20bfe6fSHong Zhang 
632d20bfe6fSHong Zhang     /* Compute P^T*A*P using outer product (P^T)[:,j]*(A*P)[j,:]. */
633d20bfe6fSHong Zhang     pnzi = pi[i+1] - pi[i];
634d20bfe6fSHong Zhang     for (j=0;j<pnzi;j++) {
635d20bfe6fSHong Zhang       nextap = 0;
636d20bfe6fSHong Zhang       crow   = *pJ++;
637d20bfe6fSHong Zhang       cjj    = cj + ci[crow];
638d20bfe6fSHong Zhang       caj    = ca + ci[crow];
639d20bfe6fSHong Zhang       /* Perform sparse axpy operation.  Note cjj includes apj. */
640d20bfe6fSHong Zhang       for (k=0;nextap<apnzj;k++) {
641d20bfe6fSHong Zhang         if (cjj[k]==apj[nextap]) {
642d20bfe6fSHong Zhang           caj[k] += (*pA)*apa[apj[nextap++]];
643d20bfe6fSHong Zhang         }
644d20bfe6fSHong Zhang       }
645d20bfe6fSHong Zhang       flops += 2*apnzj;
646d20bfe6fSHong Zhang       pA++;
647d20bfe6fSHong Zhang     }
648d20bfe6fSHong Zhang 
649d20bfe6fSHong Zhang     /* Zero the current row info for A*P */
650d20bfe6fSHong Zhang     for (j=0;j<apnzj;j++) {
651d20bfe6fSHong Zhang       apa[apj[j]]      = 0.;
652d20bfe6fSHong Zhang       apjdense[apj[j]] = 0;
653d20bfe6fSHong Zhang     }
654d20bfe6fSHong Zhang   }
655d20bfe6fSHong Zhang 
656d20bfe6fSHong Zhang   /* Assemble the final matrix and clean up */
657d20bfe6fSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
658d20bfe6fSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
659d20bfe6fSHong Zhang   ierr = PetscFree(apa);CHKERRQ(ierr);
660d20bfe6fSHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
661d20bfe6fSHong Zhang 
662eb9c0419SKris Buschelman   PetscFunctionReturn(0);
663eb9c0419SKris Buschelman }
6640e36024fSHong Zhang 
6650e36024fSHong Zhang /* Compute C = P[rstart:rend,:]^T * A * P of seqaij matrices - used by MatPtAP_MPIAIJ_MPIAIJ() */
6660e36024fSHong Zhang 
6670e36024fSHong Zhang #undef __FUNCT__
6680e36024fSHong Zhang #define __FUNCT__ "MatPtAPNumeric_SeqAIJ_SeqAIJ_ReducedPt"
6697f79fd58SMatthew Knepley /*@C
670e9b43d0fSSatish Balay    MatPtAPNumeric_SeqAIJ_SeqAIJ_ReducedPt - Creates C = P[rstart:rend,:]^T * A * P of seqaij matrices,
671b90dcfe3SHong Zhang                 used by MatPtAP_MPIAIJ_MPIAIJ()
672b90dcfe3SHong Zhang 
673b90dcfe3SHong Zhang    Collective on Mat
674b90dcfe3SHong Zhang 
675b90dcfe3SHong Zhang    Input Parameters:
676b90dcfe3SHong Zhang +  A - the matrix in seqaij format
677b90dcfe3SHong Zhang .  P - the projection matrix in seqaij format
678b90dcfe3SHong Zhang .  scall - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX
679b90dcfe3SHong Zhang .  fill - expected fill (not being used when scall==MAT_REUSE_MATRIX)
680b90dcfe3SHong Zhang .  prstart, prend - the starting and ending-1-th row of the matrix P to be used for transpose
681b90dcfe3SHong Zhang 
682b90dcfe3SHong Zhang    Output Parameters:
683b90dcfe3SHong Zhang .  C - the product matrix in seqaij format
684b90dcfe3SHong Zhang 
685b90dcfe3SHong Zhang    Level: developer
686b90dcfe3SHong Zhang 
687b90dcfe3SHong Zhang .seealso: MatPtAPSymbolic(),MatPtAPNumeric(),MatMatMult()
688b90dcfe3SHong Zhang @*/
689*b1d57f15SBarry Smith PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_ReducedPt(Mat A,Mat P,PetscInt prstart,PetscInt prend,Mat C)
6900e36024fSHong Zhang {
6910e36024fSHong Zhang   PetscErrorCode ierr;
692*b1d57f15SBarry Smith   PetscInt       flops=0;
6930e36024fSHong Zhang   Mat_SeqAIJ     *a  = (Mat_SeqAIJ *) A->data;
6940e36024fSHong Zhang   Mat_SeqAIJ     *p  = (Mat_SeqAIJ *) P->data;
6950e36024fSHong Zhang   Mat_SeqAIJ     *c  = (Mat_SeqAIJ *) C->data;
696*b1d57f15SBarry Smith   PetscInt       *ai=a->i,*aj=a->j,*apj,*apjdense;
697*b1d57f15SBarry Smith   PetscInt       *pi=p->i,*pj=p->j,*pJ=p->j+pi[prstart],*pjj;
698*b1d57f15SBarry Smith   PetscInt       *ci=c->i,*cj=c->j,*cjj;
699*b1d57f15SBarry Smith   PetscInt       am=A->M,cn=C->N,cm=C->M;
700*b1d57f15SBarry Smith   PetscInt       i,j,k,anzi,pnzi,apnzj,nextap,pnzj,prow,crow;
7010e36024fSHong Zhang   MatScalar      *aa=a->a,*apa,*pa=p->a,*pA=p->a,*paj,*ca=c->a,*caj;
7020e36024fSHong Zhang 
7030e36024fSHong Zhang   PetscFunctionBegin;
7040e36024fSHong Zhang   pA=p->a+pi[prstart];
7050e36024fSHong Zhang   /* Allocate temporary array for storage of one row of A*P */
706*b1d57f15SBarry Smith   ierr = PetscMalloc(cn*(sizeof(MatScalar)+2*sizeof(PetscInt)),&apa);CHKERRQ(ierr);
707*b1d57f15SBarry Smith   ierr = PetscMemzero(apa,cn*(sizeof(MatScalar)+2*sizeof(PetscInt)));CHKERRQ(ierr);
7080e36024fSHong Zhang 
709*b1d57f15SBarry Smith   apj      = (PetscInt *)(apa + cn);
7100e36024fSHong Zhang   apjdense = apj + cn;
7110e36024fSHong Zhang 
7120e36024fSHong Zhang   /* Clear old values in C */
7130e36024fSHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
7140e36024fSHong Zhang 
7150e36024fSHong Zhang   for (i=0;i<am;i++) {
7160e36024fSHong Zhang     /* Form sparse row of A*P */
7170e36024fSHong Zhang     anzi  = ai[i+1] - ai[i];
7180e36024fSHong Zhang     apnzj = 0;
7190e36024fSHong Zhang     for (j=0;j<anzi;j++) {
7200e36024fSHong Zhang       prow = *aj++;
7210e36024fSHong Zhang       pnzj = pi[prow+1] - pi[prow];
7220e36024fSHong Zhang       pjj  = pj + pi[prow];
7230e36024fSHong Zhang       paj  = pa + pi[prow];
7240e36024fSHong Zhang       for (k=0;k<pnzj;k++) {
7250e36024fSHong Zhang         if (!apjdense[pjj[k]]) {
7260e36024fSHong Zhang           apjdense[pjj[k]] = -1;
7270e36024fSHong Zhang           apj[apnzj++]     = pjj[k];
7280e36024fSHong Zhang         }
7290e36024fSHong Zhang         apa[pjj[k]] += (*aa)*paj[k];
7300e36024fSHong Zhang       }
7310e36024fSHong Zhang       flops += 2*pnzj;
7320e36024fSHong Zhang       aa++;
7330e36024fSHong Zhang     }
7340e36024fSHong Zhang 
7350e36024fSHong Zhang     /* Sort the j index array for quick sparse axpy. */
7360e36024fSHong Zhang     ierr = PetscSortInt(apnzj,apj);CHKERRQ(ierr);
7370e36024fSHong Zhang 
7380e36024fSHong Zhang     /* Compute P[[prstart:prend,:]^T*A*P using outer product (P^T)[:,j+prstart]*(A*P)[j,:]. */
7390e36024fSHong Zhang     pnzi = pi[i+1+prstart] - pi[i+prstart];
7400e36024fSHong Zhang     for (j=0;j<pnzi;j++) {
7410e36024fSHong Zhang       nextap = 0;
7420e36024fSHong Zhang       crow   = *pJ++;
7430e36024fSHong Zhang       cjj    = cj + ci[crow];
7440e36024fSHong Zhang       caj    = ca + ci[crow];
7450e36024fSHong Zhang       /* Perform sparse axpy operation.  Note cjj includes apj. */
7460e36024fSHong Zhang       for (k=0;nextap<apnzj;k++) {
7470e36024fSHong Zhang         if (cjj[k]==apj[nextap]) {
7480e36024fSHong Zhang           caj[k] += (*pA)*apa[apj[nextap++]];
7490e36024fSHong Zhang         }
7500e36024fSHong Zhang       }
7510e36024fSHong Zhang       flops += 2*apnzj;
7520e36024fSHong Zhang       pA++;
7530e36024fSHong Zhang     }
7540e36024fSHong Zhang 
7550e36024fSHong Zhang     /* Zero the current row info for A*P */
7560e36024fSHong Zhang     for (j=0;j<apnzj;j++) {
7570e36024fSHong Zhang       apa[apj[j]]      = 0.;
7580e36024fSHong Zhang       apjdense[apj[j]] = 0;
7590e36024fSHong Zhang     }
7600e36024fSHong Zhang   }
7610e36024fSHong Zhang 
7620e36024fSHong Zhang   /* Assemble the final matrix and clean up */
7630e36024fSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7640e36024fSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7650e36024fSHong Zhang   ierr = PetscFree(apa);CHKERRQ(ierr);
7660e36024fSHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
7670e36024fSHong Zhang 
7680e36024fSHong Zhang   PetscFunctionReturn(0);
7690e36024fSHong Zhang }
7700e36024fSHong Zhang 
7710e36024fSHong Zhang #undef __FUNCT__
7720e36024fSHong Zhang #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqAIJ_ReducedPt"
773*b1d57f15SBarry Smith PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_ReducedPt(Mat A,Mat P,PetscReal fill,PetscInt prstart,PetscInt prend,Mat *C)
774*b1d57f15SBarry Smith {
7750e36024fSHong Zhang   PetscErrorCode ierr;
7760e36024fSHong Zhang   FreeSpaceList  free_space=PETSC_NULL,current_space=PETSC_NULL;
7770e36024fSHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*p=(Mat_SeqAIJ*)P->data,*c;
778*b1d57f15SBarry Smith   PetscInt       *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj;
779*b1d57f15SBarry Smith   PetscInt       *ci,*cj,*denserow,*sparserow,*ptadenserow,*ptasparserow,*ptaj;
780*b1d57f15SBarry Smith   PetscInt       an=A->N,am=A->M,pn=P->N,pm=P->M;
781*b1d57f15SBarry Smith   PetscInt       i,j,k,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi;
782*b1d57f15SBarry Smith   PetscInt       m = prend - prstart;
7830e36024fSHong Zhang   MatScalar      *ca;
7840e36024fSHong Zhang   Mat            *psub,P_sub;
7850e36024fSHong Zhang   IS             isrow,iscol;
7860b89d903Svictorle 
7870b89d903Svictorle   PetscFunctionBegin;
7880b89d903Svictorle   /* Get ij structure of P[rstart:rend,:]^T */
7890e36024fSHong Zhang   ierr = ISCreateStride(PETSC_COMM_SELF,m,prstart,1,&isrow);CHKERRQ(ierr);
7900e36024fSHong Zhang   ierr = ISCreateStride(PETSC_COMM_SELF,P->n,0,1,&iscol);CHKERRQ(ierr);
7910e36024fSHong Zhang   ierr = MatGetSubMatrices(P,1,&isrow,&iscol,MAT_INITIAL_MATRIX,&psub);CHKERRQ(ierr);
7920e36024fSHong Zhang   ierr = ISDestroy(isrow);CHKERRQ(ierr);
7930e36024fSHong Zhang   ierr = ISDestroy(iscol);CHKERRQ(ierr);
7940e36024fSHong Zhang   P_sub = psub[0];
7950e36024fSHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(P_sub,&pti,&ptj);CHKERRQ(ierr);
7960e36024fSHong Zhang   ierr = MatDestroyMatrices(1,&psub);CHKERRQ(ierr);
7970e36024fSHong Zhang   ptJ=ptj;
7980e36024fSHong Zhang 
7990e36024fSHong Zhang   /* Allocate ci array, arrays for fill computation and */
8000e36024fSHong Zhang   /* free space for accumulating nonzero column info */
801*b1d57f15SBarry Smith   ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
8020e36024fSHong Zhang   ci[0] = 0;
8030e36024fSHong Zhang 
804*b1d57f15SBarry Smith   ierr = PetscMalloc((2*pn+2*an+1)*sizeof(PetscInt),&ptadenserow);CHKERRQ(ierr);
805*b1d57f15SBarry Smith   ierr = PetscMemzero(ptadenserow,(2*pn+2*an+1)*sizeof(PetscInt));CHKERRQ(ierr);
8060e36024fSHong Zhang   ptasparserow = ptadenserow  + an;
8070e36024fSHong Zhang   denserow     = ptasparserow + an;
8080e36024fSHong Zhang   sparserow    = denserow     + pn;
8090e36024fSHong Zhang 
8100e36024fSHong Zhang   /* Set initial free space to be nnz(A) scaled by aspect ratio of P. */
8110e36024fSHong Zhang   /* This should be reasonable if sparsity of PtAP is similar to that of A. */
812*b1d57f15SBarry Smith   ierr          = GetMoreSpace((PetscInt)(fill*ai[am]/pm)*pn,&free_space);
8130e36024fSHong Zhang   current_space = free_space;
8140e36024fSHong Zhang 
8150e36024fSHong Zhang   /* Determine symbolic info for each row of C: */
8160e36024fSHong Zhang   for (i=0;i<pn;i++) {
8170e36024fSHong Zhang     ptnzi  = pti[i+1] - pti[i];
8180e36024fSHong Zhang     ptanzi = 0;
8190e36024fSHong Zhang     /* Determine symbolic row of PtA_reduced: */
8200e36024fSHong Zhang     for (j=0;j<ptnzi;j++) {
8210e36024fSHong Zhang       arow = *ptJ++;
8220e36024fSHong Zhang       anzj = ai[arow+1] - ai[arow];
8230e36024fSHong Zhang       ajj  = aj + ai[arow];
8240e36024fSHong Zhang       for (k=0;k<anzj;k++) {
8250e36024fSHong Zhang         if (!ptadenserow[ajj[k]]) {
8260e36024fSHong Zhang           ptadenserow[ajj[k]]    = -1;
8270e36024fSHong Zhang           ptasparserow[ptanzi++] = ajj[k];
8280e36024fSHong Zhang         }
8290e36024fSHong Zhang       }
8300e36024fSHong Zhang     }
8310e36024fSHong Zhang       /* Using symbolic info for row of PtA, determine symbolic info for row of C: */
8320e36024fSHong Zhang     ptaj = ptasparserow;
8330e36024fSHong Zhang     cnzi   = 0;
8340e36024fSHong Zhang     for (j=0;j<ptanzi;j++) {
8350e36024fSHong Zhang       prow = *ptaj++;
8360e36024fSHong Zhang       pnzj = pi[prow+1] - pi[prow];
8370e36024fSHong Zhang       pjj  = pj + pi[prow];
8380e36024fSHong Zhang       for (k=0;k<pnzj;k++) {
8390e36024fSHong Zhang         if (!denserow[pjj[k]]) {
8400e36024fSHong Zhang             denserow[pjj[k]]  = -1;
8410e36024fSHong Zhang             sparserow[cnzi++] = pjj[k];
8420e36024fSHong Zhang         }
8430e36024fSHong Zhang       }
8440e36024fSHong Zhang     }
8450e36024fSHong Zhang 
8460e36024fSHong Zhang     /* sort sparserow */
8470e36024fSHong Zhang     ierr = PetscSortInt(cnzi,sparserow);CHKERRQ(ierr);
8480e36024fSHong Zhang 
8490e36024fSHong Zhang     /* If free space is not available, make more free space */
8500e36024fSHong Zhang     /* Double the amount of total space in the list */
8510e36024fSHong Zhang     if (current_space->local_remaining<cnzi) {
8520e36024fSHong Zhang       ierr = GetMoreSpace(current_space->total_array_size,&current_space);CHKERRQ(ierr);
8530e36024fSHong Zhang     }
8540e36024fSHong Zhang 
8550e36024fSHong Zhang     /* Copy data into free space, and zero out denserows */
856*b1d57f15SBarry Smith     ierr = PetscMemcpy(current_space->array,sparserow,cnzi*sizeof(PetscInt));CHKERRQ(ierr);
8570e36024fSHong Zhang     current_space->array           += cnzi;
8580e36024fSHong Zhang     current_space->local_used      += cnzi;
8590e36024fSHong Zhang     current_space->local_remaining -= cnzi;
8600e36024fSHong Zhang 
8610e36024fSHong Zhang     for (j=0;j<ptanzi;j++) {
8620e36024fSHong Zhang       ptadenserow[ptasparserow[j]] = 0;
8630e36024fSHong Zhang     }
8640e36024fSHong Zhang     for (j=0;j<cnzi;j++) {
8650e36024fSHong Zhang       denserow[sparserow[j]] = 0;
8660e36024fSHong Zhang     }
8670e36024fSHong Zhang       /* Aside: Perhaps we should save the pta info for the numerical factorization. */
8680e36024fSHong Zhang       /*        For now, we will recompute what is needed. */
8690e36024fSHong Zhang     ci[i+1] = ci[i] + cnzi;
8700e36024fSHong Zhang   }
8710e36024fSHong Zhang   /* nnz is now stored in ci[ptm], column indices are in the list of free space */
8720e36024fSHong Zhang   /* Allocate space for cj, initialize cj, and */
8730e36024fSHong Zhang   /* destroy list of free space and other temporary array(s) */
874*b1d57f15SBarry Smith   ierr = PetscMalloc((ci[pn]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
8750e36024fSHong Zhang   ierr = MakeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
8760e36024fSHong Zhang   ierr = PetscFree(ptadenserow);CHKERRQ(ierr);
8770e36024fSHong Zhang 
8780e36024fSHong Zhang   /* Allocate space for ca */
8790e36024fSHong Zhang   ierr = PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
8800e36024fSHong Zhang   ierr = PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));CHKERRQ(ierr);
8810e36024fSHong Zhang 
8820e36024fSHong Zhang   /* put together the new matrix */
8830e36024fSHong Zhang   ierr = MatCreateSeqAIJWithArrays(A->comm,pn,pn,ci,cj,ca,C);CHKERRQ(ierr);
8840e36024fSHong Zhang 
8850e36024fSHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
8860e36024fSHong Zhang   /* Since these are PETSc arrays, change flags to free them as necessary. */
8870e36024fSHong Zhang   c = (Mat_SeqAIJ *)((*C)->data);
8880e36024fSHong Zhang   c->freedata = PETSC_TRUE;
8890e36024fSHong Zhang   c->nonew    = 0;
8900e36024fSHong Zhang 
8910e36024fSHong Zhang   /* Clean up. */
8920e36024fSHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
8930e36024fSHong Zhang 
8940e36024fSHong Zhang   PetscFunctionReturn(0);
8950e36024fSHong Zhang }
8960e36024fSHong Zhang 
8970e36024fSHong Zhang #undef __FUNCT__
8980e36024fSHong Zhang #define __FUNCT__ "MatPtAP_SeqAIJ_SeqAIJ_ReducedPt"
899*b1d57f15SBarry Smith PetscErrorCode MatPtAP_SeqAIJ_SeqAIJ_ReducedPt(Mat A,Mat P,MatReuse scall,PetscReal fill,PetscInt prstart,PetscInt prend,Mat *C)
9000e36024fSHong Zhang {
9010e36024fSHong Zhang   PetscErrorCode ierr;
902*b1d57f15SBarry Smith 
9030e36024fSHong Zhang   PetscFunctionBegin;
9040e36024fSHong Zhang   if (A->m != prend-prstart) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",A->m,prend-prstart);
9050e36024fSHong 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);
9060e36024fSHong Zhang   if (scall == MAT_INITIAL_MATRIX){
9070e36024fSHong Zhang     ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ_ReducedPt(A,P,fill,prstart,prend,C);CHKERRQ(ierr);
9080e36024fSHong Zhang   }
9090e36024fSHong Zhang 
9100e36024fSHong Zhang   ierr = MatPtAPNumeric_SeqAIJ_SeqAIJ_ReducedPt(A,P,prstart,prend,*C);CHKERRQ(ierr);
9110e36024fSHong Zhang 
9120e36024fSHong Zhang   PetscFunctionReturn(0);
9130e36024fSHong Zhang }
9140e36024fSHong Zhang 
915