xref: /petsc/src/mat/impls/aij/seq/matptap.c (revision 0b89d90399720d2b703ad352fcc8545943435fa9)
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);
75eb9c0419SKris Buschelman #undef __FUNCT__
76ff134f7aSHong Zhang #define __FUNCT__ "MatPtAP_MPIAIJ_MPIAIJ"
77ff134f7aSHong Zhang PetscErrorCode MatPtAP_MPIAIJ_MPIAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
78ff134f7aSHong Zhang {
79ff134f7aSHong Zhang   PetscErrorCode    ierr;
800e36024fSHong Zhang   Mat               AP,AP_seq,P_seq,A_loc,C_seq;
81ff134f7aSHong Zhang   Mat_MPIAIJ        *p = (Mat_MPIAIJ*)P->data;
82ff134f7aSHong Zhang   Mat_MatMatMultMPI *mult;
830e36024fSHong Zhang 
840e36024fSHong Zhang   int               prstart,prend,m=P->m;
850e36024fSHong Zhang   int               rank,prid=10;
86*0b89d903Svictorle   IS  isrow,iscol;
87*0b89d903Svictorle   Mat P_subseq,*psubseq;
88ff134f7aSHong Zhang 
89ff134f7aSHong Zhang   PetscFunctionBegin;
90ff134f7aSHong Zhang   ierr = MPI_Comm_rank(A->comm,&rank);CHKERRQ(ierr);
91ff134f7aSHong Zhang 
920e36024fSHong Zhang   /* compute symbolic and numeric AP = A*P */
930e36024fSHong Zhang   ierr = MatMatMult_MPIAIJ_MPIAIJ(A,P,scall,fill,&AP);CHKERRQ(ierr);
940e36024fSHong Zhang   mult = (Mat_MatMatMultMPI*)AP->spptr;
950e36024fSHong Zhang   P_seq = mult->bseq[0]; /* = submatrix of P by taking rows of P that equal to nonzero col of A */
960e36024fSHong Zhang   A_loc = mult->aseq[0]; /* = submatrix of A by taking all local rows of A */
97ff134f7aSHong Zhang   AP_seq  = mult->C_seq;
980e36024fSHong Zhang 
990e36024fSHong Zhang   if (rank == prid){
1000e36024fSHong Zhang     ierr = PetscPrintf(PETSC_COMM_SELF," [%d] A_loc: %d, %d\n",rank,A_loc->m,A_loc->n);CHKERRQ(ierr);
1010e36024fSHong Zhang     ierr = MatView(A_loc,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
1020e36024fSHong Zhang     ierr = PetscPrintf(PETSC_COMM_SELF," [%d] P_seq: %d, %d\n",rank,P_seq->m,P_seq->n);CHKERRQ(ierr);
1030e36024fSHong Zhang     ierr = MatView(P_seq,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
1040e36024fSHong Zhang     ierr = PetscPrintf(PETSC_COMM_SELF," [%d]  AP_seq: %d, %d\n",rank,AP_seq->m,AP_seq->n);
1050e36024fSHong Zhang     ierr = MatView(AP_seq,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
1060e36024fSHong Zhang   }
1070e36024fSHong Zhang 
108ff134f7aSHong Zhang   prstart = mult->brstart;
109ff134f7aSHong Zhang   prend   = prstart + m;
1100e36024fSHong Zhang   /*
111ff134f7aSHong Zhang   ierr = PetscPrintf(PETSC_COMM_SELF," [%d] prstart: %d, prend: %d, dim of P_seq: %d, %d\n",rank,prstart,prend,P_seq->m,P_seq->n);
1120e36024fSHong Zhang   */
113ff134f7aSHong Zhang 
114ff134f7aSHong Zhang     /* get seq matrix P_subseq by taking local rows of P */
115ff134f7aSHong Zhang   ierr = ISCreateStride(PETSC_COMM_SELF,m,prstart,1,&isrow);CHKERRQ(ierr);
116ff134f7aSHong Zhang   ierr = ISCreateStride(PETSC_COMM_SELF,P_seq->n,0,1,&iscol);CHKERRQ(ierr);
117ff134f7aSHong Zhang   ierr = MatGetSubMatrices(P_seq,1,&isrow,&iscol,MAT_INITIAL_MATRIX,&psubseq);CHKERRQ(ierr);
118ff134f7aSHong Zhang   P_subseq = psubseq[0];
1190e36024fSHong Zhang   ierr = ISDestroy(isrow);CHKERRQ(ierr);
1200e36024fSHong Zhang   ierr = ISDestroy(iscol);CHKERRQ(ierr);
121ff134f7aSHong Zhang 
1220e36024fSHong Zhang   /* compute P_subseq^T*AP_seq */
1230e36024fSHong Zhang   ierr = MatMatMultTranspose_SeqAIJ_SeqAIJ(P_subseq,AP_seq,scall,fill,&C_seq);CHKERRQ(ierr);
1240e36024fSHong Zhang   if (rank == prid){
1250e36024fSHong Zhang     ierr = PetscPrintf(PETSC_COMM_SELF," [%d] dim of P_subseq: %d, %d; AP_seq: %d, %d\n",rank,P_subseq->m,P_subseq->n,AP_seq->m,AP_seq->n);
1260e36024fSHong Zhang     ierr = PetscPrintf(PETSC_COMM_SELF," [%d] C_seq: %d, %d\n",rank,C_seq->m,C_seq->n);CHKERRQ(ierr);
1270e36024fSHong Zhang     ierr = MatView(C_seq,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
128ff134f7aSHong Zhang   }
1290e36024fSHong Zhang   ierr = MatDestroyMatrices(1,&psubseq);CHKERRQ(ierr);
130ff134f7aSHong Zhang 
1310e36024fSHong Zhang   /* ierr = MatPtAP_SeqAIJ_SeqAIJ_ReducedPt(A_loc,P_seq,scall,fill,prstart,prend,&C_seq);CHKERRQ(ierr); */
1320e36024fSHong Zhang   /* ierr = PetscPrintf(PETSC_COMM_SELF," [%d] C_seq dim: %d, %d\n",rank,C_seq->m,C_seq->n); */
1330e36024fSHong Zhang 
1340e36024fSHong Zhang   /* add C_seq into C */
1350e36024fSHong Zhang   ierr = MatMerge_SeqsToMPI(A->comm,C_seq,P->n,P->n,scall,C);CHKERRQ(ierr);
1360e36024fSHong Zhang 
1370e36024fSHong Zhang   /* clean up */
1380e36024fSHong Zhang   ierr = MatDestroy(AP);CHKERRQ(ierr);
1390e36024fSHong Zhang   ierr = MatDestroy(C_seq);CHKERRQ(ierr);
1400e36024fSHong Zhang 
141ff134f7aSHong Zhang   PetscFunctionReturn(0);
142ff134f7aSHong Zhang }
143ff134f7aSHong Zhang 
144ff134f7aSHong Zhang #undef __FUNCT__
145ff134f7aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ"
146ff134f7aSHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
147ff134f7aSHong Zhang {
148ff134f7aSHong Zhang 
149ff134f7aSHong Zhang   PetscFunctionBegin;
150ff134f7aSHong Zhang   PetscFunctionReturn(0);
151ff134f7aSHong Zhang }
152ff134f7aSHong Zhang 
153ff134f7aSHong Zhang #undef __FUNCT__
154ff134f7aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ"
155ff134f7aSHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat B,Mat C)
156ff134f7aSHong Zhang {
157ff134f7aSHong Zhang 
158ff134f7aSHong Zhang   PetscFunctionBegin;
159ff134f7aSHong Zhang   PetscFunctionReturn(0);
160ff134f7aSHong Zhang }
161ff134f7aSHong Zhang 
162ff134f7aSHong Zhang #undef __FUNCT__
1639af31e4aSHong Zhang #define __FUNCT__ "MatPtAP_SeqAIJ_SeqAIJ"
164dfbe8321SBarry Smith PetscErrorCode MatPtAP_SeqAIJ_SeqAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
1659af31e4aSHong Zhang {
166dfbe8321SBarry Smith   PetscErrorCode ierr;
1679af31e4aSHong Zhang   PetscFunctionBegin;
1689af31e4aSHong Zhang   if (scall == MAT_INITIAL_MATRIX){
169d20bfe6fSHong Zhang     ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
1709af31e4aSHong Zhang     ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ(A,P,fill,C);CHKERRQ(ierr);
171d20bfe6fSHong Zhang     ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
1729af31e4aSHong Zhang   }
173d20bfe6fSHong Zhang   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
1749af31e4aSHong Zhang   ierr = MatPtAPNumeric_SeqAIJ_SeqAIJ(A,P,*C);CHKERRQ(ierr);
175d20bfe6fSHong Zhang   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
1769af31e4aSHong Zhang   PetscFunctionReturn(0);
1779af31e4aSHong Zhang }
1789af31e4aSHong Zhang 
1799af31e4aSHong Zhang #undef __FUNCT__
1809af31e4aSHong Zhang #define __FUNCT__ "MatPtAPSymbolic"
1816849ba73SBarry Smith /*
1829af31e4aSHong Zhang    MatPtAPSymbolic - Creates the (i,j) structure of the matrix projection C = P^T * A * P
1834d3841fdSKris Buschelman 
1844d3841fdSKris Buschelman    Collective on Mat
1854d3841fdSKris Buschelman 
1864d3841fdSKris Buschelman    Input Parameters:
1874d3841fdSKris Buschelman +  A - the matrix
1884d3841fdSKris Buschelman -  P - the projection matrix
1894d3841fdSKris Buschelman 
1904d3841fdSKris Buschelman    Output Parameters:
1914d3841fdSKris Buschelman .  C - the (i,j) structure of the product matrix
1924d3841fdSKris Buschelman 
1934d3841fdSKris Buschelman    Notes:
1944d3841fdSKris Buschelman    C will be created and must be destroyed by the user with MatDestroy().
1954d3841fdSKris Buschelman 
1964d3841fdSKris Buschelman    This routine is currently only implemented for pairs of SeqAIJ matrices and classes
1974d3841fdSKris Buschelman    which inherit from SeqAIJ.  C will be of type MATSEQAIJ.  The product is computed using
1989af31e4aSHong Zhang    this (i,j) structure by calling MatPtAPNumeric().
1994d3841fdSKris Buschelman 
2004d3841fdSKris Buschelman    Level: intermediate
2014d3841fdSKris Buschelman 
2029af31e4aSHong Zhang .seealso: MatPtAP(),MatPtAPNumeric(),MatMatMultSymbolic()
2036849ba73SBarry Smith */
204dfbe8321SBarry Smith PetscErrorCode MatPtAPSymbolic(Mat A,Mat P,PetscReal fill,Mat *C) {
205dfbe8321SBarry Smith   PetscErrorCode ierr;
206534c1384SKris Buschelman   PetscErrorCode (*fA)(Mat,Mat,PetscReal,Mat*);
207534c1384SKris Buschelman   PetscErrorCode (*fP)(Mat,Mat,PetscReal,Mat*);
208eb9c0419SKris Buschelman 
209eb9c0419SKris Buschelman   PetscFunctionBegin;
210eb9c0419SKris Buschelman 
2114482741eSBarry Smith   PetscValidHeaderSpecific(A,MAT_COOKIE,1);
212c9780b6fSBarry Smith   PetscValidType(A,1);
213eb9c0419SKris Buschelman   MatPreallocated(A);
214eb9c0419SKris Buschelman   if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
215eb9c0419SKris Buschelman   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
216eb9c0419SKris Buschelman 
2174482741eSBarry Smith   PetscValidHeaderSpecific(P,MAT_COOKIE,2);
218c9780b6fSBarry Smith   PetscValidType(P,2);
219eb9c0419SKris Buschelman   MatPreallocated(P);
220eb9c0419SKris Buschelman   if (!P->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
221eb9c0419SKris Buschelman   if (P->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
222eb9c0419SKris Buschelman 
2234482741eSBarry Smith   PetscValidPointer(C,3);
2244482741eSBarry Smith 
225eb9c0419SKris Buschelman   if (P->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->M,A->N);
226eb9c0419SKris Buschelman   if (A->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix 'A' must be square, %d != %d",A->M,A->N);
227eb9c0419SKris Buschelman 
228534c1384SKris Buschelman   /* For now, we do not dispatch based on the type of A and P */
229534c1384SKris Buschelman   /* When implementations like _SeqAIJ_MAIJ exist, attack the multiple dispatch problem. */
230534c1384SKris Buschelman   fA = A->ops->ptapsymbolic;
231534c1384SKris Buschelman   if (!fA) SETERRQ1(PETSC_ERR_SUP,"MatPtAPSymbolic not supported for A of type %s",A->type_name);
232534c1384SKris Buschelman   fP = P->ops->ptapsymbolic;
233534c1384SKris Buschelman   if (!fP) SETERRQ1(PETSC_ERR_SUP,"MatPtAPSymbolic not supported for P of type %s",P->type_name);
234534c1384SKris 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);
2354d3841fdSKris Buschelman 
236534c1384SKris Buschelman   ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
237534c1384SKris Buschelman   ierr = (*fA)(A,P,fill,C);CHKERRQ(ierr);
238534c1384SKris Buschelman   ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
239eb9c0419SKris Buschelman 
240eb9c0419SKris Buschelman   PetscFunctionReturn(0);
241eb9c0419SKris Buschelman }
242eb9c0419SKris Buschelman 
243f747e1aeSHong Zhang typedef struct {
244f747e1aeSHong Zhang   Mat    symAP;
245f747e1aeSHong Zhang } Mat_PtAPstruct;
246f747e1aeSHong Zhang 
24778a80504SBarry Smith EXTERN PetscErrorCode MatDestroy_SeqAIJ(Mat);
24878a80504SBarry Smith 
249f747e1aeSHong Zhang #undef __FUNCT__
250f747e1aeSHong Zhang #define __FUNCT__ "MatDestroy_SeqAIJ_PtAP"
251f4a850bbSBarry Smith PetscErrorCode MatDestroy_SeqAIJ_PtAP(Mat A)
252f747e1aeSHong Zhang {
253f4a850bbSBarry Smith   PetscErrorCode    ierr;
254f747e1aeSHong Zhang   Mat_PtAPstruct    *ptap=(Mat_PtAPstruct*)A->spptr;
255f747e1aeSHong Zhang 
256f747e1aeSHong Zhang   PetscFunctionBegin;
257f747e1aeSHong Zhang   ierr = MatDestroy(ptap->symAP);CHKERRQ(ierr);
258f747e1aeSHong Zhang   ierr = PetscFree(ptap);CHKERRQ(ierr);
25978a80504SBarry Smith   ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr);
260f747e1aeSHong Zhang   PetscFunctionReturn(0);
261f747e1aeSHong Zhang }
262f747e1aeSHong Zhang 
263eb9c0419SKris Buschelman #undef __FUNCT__
2649af31e4aSHong Zhang #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqAIJ"
265dfbe8321SBarry Smith PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat P,PetscReal fill,Mat *C) {
266dfbe8321SBarry Smith   PetscErrorCode ierr;
267d20bfe6fSHong Zhang   FreeSpaceList  free_space=PETSC_NULL,current_space=PETSC_NULL;
268d20bfe6fSHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*p=(Mat_SeqAIJ*)P->data,*c;
269d20bfe6fSHong Zhang   int            *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj;
270d20bfe6fSHong Zhang   int            *ci,*cj,*denserow,*sparserow,*ptadenserow,*ptasparserow,*ptaj;
271d20bfe6fSHong Zhang   int            an=A->N,am=A->M,pn=P->N,pm=P->M;
272d20bfe6fSHong Zhang   int            i,j,k,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi;
273d20bfe6fSHong Zhang   MatScalar      *ca;
274eb9c0419SKris Buschelman 
275eb9c0419SKris Buschelman   PetscFunctionBegin;
276d20bfe6fSHong Zhang   /* Get ij structure of P^T */
277eb9c0419SKris Buschelman   ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
278d20bfe6fSHong Zhang   ptJ=ptj;
279eb9c0419SKris Buschelman 
280d20bfe6fSHong Zhang   /* Allocate ci array, arrays for fill computation and */
281d20bfe6fSHong Zhang   /* free space for accumulating nonzero column info */
282d20bfe6fSHong Zhang   ierr = PetscMalloc((pn+1)*sizeof(int),&ci);CHKERRQ(ierr);
283d20bfe6fSHong Zhang   ci[0] = 0;
284eb9c0419SKris Buschelman 
285d20bfe6fSHong Zhang   ierr = PetscMalloc((2*pn+2*an+1)*sizeof(int),&ptadenserow);CHKERRQ(ierr);
286d20bfe6fSHong Zhang   ierr = PetscMemzero(ptadenserow,(2*pn+2*an+1)*sizeof(int));CHKERRQ(ierr);
287d20bfe6fSHong Zhang   ptasparserow = ptadenserow  + an;
288d20bfe6fSHong Zhang   denserow     = ptasparserow + an;
289d20bfe6fSHong Zhang   sparserow    = denserow     + pn;
290eb9c0419SKris Buschelman 
291d20bfe6fSHong Zhang   /* Set initial free space to be nnz(A) scaled by aspect ratio of P. */
292d20bfe6fSHong Zhang   /* This should be reasonable if sparsity of PtAP is similar to that of A. */
293d20bfe6fSHong Zhang   ierr          = GetMoreSpace((ai[am]/pm)*pn,&free_space);
294d20bfe6fSHong Zhang   current_space = free_space;
295d20bfe6fSHong Zhang 
296d20bfe6fSHong Zhang   /* Determine symbolic info for each row of C: */
297d20bfe6fSHong Zhang   for (i=0;i<pn;i++) {
298d20bfe6fSHong Zhang     ptnzi  = pti[i+1] - pti[i];
299d20bfe6fSHong Zhang     ptanzi = 0;
300d20bfe6fSHong Zhang     /* Determine symbolic row of PtA: */
301d20bfe6fSHong Zhang     for (j=0;j<ptnzi;j++) {
302d20bfe6fSHong Zhang       arow = *ptJ++;
303d20bfe6fSHong Zhang       anzj = ai[arow+1] - ai[arow];
304d20bfe6fSHong Zhang       ajj  = aj + ai[arow];
305d20bfe6fSHong Zhang       for (k=0;k<anzj;k++) {
306d20bfe6fSHong Zhang         if (!ptadenserow[ajj[k]]) {
307d20bfe6fSHong Zhang           ptadenserow[ajj[k]]    = -1;
308d20bfe6fSHong Zhang           ptasparserow[ptanzi++] = ajj[k];
309d20bfe6fSHong Zhang         }
310d20bfe6fSHong Zhang       }
311d20bfe6fSHong Zhang     }
312d20bfe6fSHong Zhang       /* Using symbolic info for row of PtA, determine symbolic info for row of C: */
313d20bfe6fSHong Zhang     ptaj = ptasparserow;
314d20bfe6fSHong Zhang     cnzi   = 0;
315d20bfe6fSHong Zhang     for (j=0;j<ptanzi;j++) {
316d20bfe6fSHong Zhang       prow = *ptaj++;
317d20bfe6fSHong Zhang       pnzj = pi[prow+1] - pi[prow];
318d20bfe6fSHong Zhang       pjj  = pj + pi[prow];
319d20bfe6fSHong Zhang       for (k=0;k<pnzj;k++) {
320d20bfe6fSHong Zhang         if (!denserow[pjj[k]]) {
321d20bfe6fSHong Zhang             denserow[pjj[k]]  = -1;
322d20bfe6fSHong Zhang             sparserow[cnzi++] = pjj[k];
323d20bfe6fSHong Zhang         }
324d20bfe6fSHong Zhang       }
325d20bfe6fSHong Zhang     }
326d20bfe6fSHong Zhang 
327d20bfe6fSHong Zhang     /* sort sparserow */
328d20bfe6fSHong Zhang     ierr = PetscSortInt(cnzi,sparserow);CHKERRQ(ierr);
329d20bfe6fSHong Zhang 
330d20bfe6fSHong Zhang     /* If free space is not available, make more free space */
331d20bfe6fSHong Zhang     /* Double the amount of total space in the list */
332d20bfe6fSHong Zhang     if (current_space->local_remaining<cnzi) {
333d20bfe6fSHong Zhang       ierr = GetMoreSpace(current_space->total_array_size,&current_space);CHKERRQ(ierr);
334d20bfe6fSHong Zhang     }
335d20bfe6fSHong Zhang 
336d20bfe6fSHong Zhang     /* Copy data into free space, and zero out denserows */
337d20bfe6fSHong Zhang     ierr = PetscMemcpy(current_space->array,sparserow,cnzi*sizeof(int));CHKERRQ(ierr);
338d20bfe6fSHong Zhang     current_space->array           += cnzi;
339d20bfe6fSHong Zhang     current_space->local_used      += cnzi;
340d20bfe6fSHong Zhang     current_space->local_remaining -= cnzi;
341d20bfe6fSHong Zhang 
342d20bfe6fSHong Zhang     for (j=0;j<ptanzi;j++) {
343d20bfe6fSHong Zhang       ptadenserow[ptasparserow[j]] = 0;
344d20bfe6fSHong Zhang     }
345d20bfe6fSHong Zhang     for (j=0;j<cnzi;j++) {
346d20bfe6fSHong Zhang       denserow[sparserow[j]] = 0;
347d20bfe6fSHong Zhang     }
348d20bfe6fSHong Zhang       /* Aside: Perhaps we should save the pta info for the numerical factorization. */
349d20bfe6fSHong Zhang       /*        For now, we will recompute what is needed. */
350d20bfe6fSHong Zhang     ci[i+1] = ci[i] + cnzi;
351d20bfe6fSHong Zhang   }
352d20bfe6fSHong Zhang   /* nnz is now stored in ci[ptm], column indices are in the list of free space */
353d20bfe6fSHong Zhang   /* Allocate space for cj, initialize cj, and */
354d20bfe6fSHong Zhang   /* destroy list of free space and other temporary array(s) */
355d20bfe6fSHong Zhang   ierr = PetscMalloc((ci[pn]+1)*sizeof(int),&cj);CHKERRQ(ierr);
356d20bfe6fSHong Zhang   ierr = MakeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
357d20bfe6fSHong Zhang   ierr = PetscFree(ptadenserow);CHKERRQ(ierr);
358d20bfe6fSHong Zhang 
359d20bfe6fSHong Zhang   /* Allocate space for ca */
360d20bfe6fSHong Zhang   ierr = PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
361d20bfe6fSHong Zhang   ierr = PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));CHKERRQ(ierr);
362d20bfe6fSHong Zhang 
363d20bfe6fSHong Zhang   /* put together the new matrix */
364d20bfe6fSHong Zhang   ierr = MatCreateSeqAIJWithArrays(A->comm,pn,pn,ci,cj,ca,C);CHKERRQ(ierr);
365d20bfe6fSHong Zhang 
366d20bfe6fSHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
367d20bfe6fSHong Zhang   /* Since these are PETSc arrays, change flags to free them as necessary. */
368d20bfe6fSHong Zhang   c = (Mat_SeqAIJ *)((*C)->data);
369d20bfe6fSHong Zhang   c->freedata = PETSC_TRUE;
370d20bfe6fSHong Zhang   c->nonew    = 0;
371d20bfe6fSHong Zhang 
372d20bfe6fSHong Zhang   /* Clean up. */
373d20bfe6fSHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
374eb9c0419SKris Buschelman 
375eb9c0419SKris Buschelman   PetscFunctionReturn(0);
376eb9c0419SKris Buschelman }
377eb9c0419SKris Buschelman 
3783985e5eaSKris Buschelman #include "src/mat/impls/maij/maij.h"
3793985e5eaSKris Buschelman EXTERN_C_BEGIN
3803985e5eaSKris Buschelman #undef __FUNCT__
3819af31e4aSHong Zhang #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqMAIJ"
382dfbe8321SBarry Smith PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqMAIJ(Mat A,Mat PP,Mat *C) {
3835c66b693SKris Buschelman   /* This routine requires testing -- I don't think it works. */
384dfbe8321SBarry Smith   PetscErrorCode ierr;
3853985e5eaSKris Buschelman   FreeSpaceList  free_space=PETSC_NULL,current_space=PETSC_NULL;
3863985e5eaSKris Buschelman   Mat_SeqMAIJ    *pp=(Mat_SeqMAIJ*)PP->data;
3873985e5eaSKris Buschelman   Mat            P=pp->AIJ;
3883985e5eaSKris Buschelman   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*p=(Mat_SeqAIJ*)P->data,*c;
3893985e5eaSKris Buschelman   int            *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj;
3903985e5eaSKris Buschelman   int            *ci,*cj,*denserow,*sparserow,*ptadenserow,*ptasparserow,*ptaj;
3913985e5eaSKris Buschelman   int            an=A->N,am=A->M,pn=P->N,pm=P->M,ppdof=pp->dof;
392fe05a634SKris Buschelman   int            i,j,k,dof,pdof,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi;
3933985e5eaSKris Buschelman   MatScalar      *ca;
3943985e5eaSKris Buschelman 
3953985e5eaSKris Buschelman   PetscFunctionBegin;
3963985e5eaSKris Buschelman   /* Start timer */
3979af31e4aSHong Zhang   ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,PP,0,0);CHKERRQ(ierr);
3983985e5eaSKris Buschelman 
3993985e5eaSKris Buschelman   /* Get ij structure of P^T */
4003985e5eaSKris Buschelman   ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
4013985e5eaSKris Buschelman 
4023985e5eaSKris Buschelman   /* Allocate ci array, arrays for fill computation and */
4033985e5eaSKris Buschelman   /* free space for accumulating nonzero column info */
4043985e5eaSKris Buschelman   ierr = PetscMalloc((pn+1)*sizeof(int),&ci);CHKERRQ(ierr);
4053985e5eaSKris Buschelman   ci[0] = 0;
4063985e5eaSKris Buschelman 
4073985e5eaSKris Buschelman   ierr = PetscMalloc((2*pn+2*an+1)*sizeof(int),&ptadenserow);CHKERRQ(ierr);
4083985e5eaSKris Buschelman   ierr = PetscMemzero(ptadenserow,(2*pn+2*an+1)*sizeof(int));CHKERRQ(ierr);
4093985e5eaSKris Buschelman   ptasparserow = ptadenserow  + an;
4103985e5eaSKris Buschelman   denserow     = ptasparserow + an;
4113985e5eaSKris Buschelman   sparserow    = denserow     + pn;
4123985e5eaSKris Buschelman 
4133985e5eaSKris Buschelman   /* Set initial free space to be nnz(A) scaled by aspect ratio of P. */
4143985e5eaSKris Buschelman   /* This should be reasonable if sparsity of PtAP is similar to that of A. */
4153985e5eaSKris Buschelman   ierr          = GetMoreSpace((ai[am]/pm)*pn,&free_space);
4163985e5eaSKris Buschelman   current_space = free_space;
4173985e5eaSKris Buschelman 
4183985e5eaSKris Buschelman   /* Determine symbolic info for each row of C: */
4193985e5eaSKris Buschelman   for (i=0;i<pn/ppdof;i++) {
4203985e5eaSKris Buschelman     ptnzi  = pti[i+1] - pti[i];
4213985e5eaSKris Buschelman     ptanzi = 0;
4223985e5eaSKris Buschelman     ptJ    = ptj + pti[i];
4233985e5eaSKris Buschelman     for (dof=0;dof<ppdof;dof++) {
4243985e5eaSKris Buschelman     /* Determine symbolic row of PtA: */
4253985e5eaSKris Buschelman       for (j=0;j<ptnzi;j++) {
4263985e5eaSKris Buschelman         arow = ptJ[j] + dof;
4273985e5eaSKris Buschelman         anzj = ai[arow+1] - ai[arow];
4283985e5eaSKris Buschelman         ajj  = aj + ai[arow];
4293985e5eaSKris Buschelman         for (k=0;k<anzj;k++) {
4303985e5eaSKris Buschelman           if (!ptadenserow[ajj[k]]) {
4313985e5eaSKris Buschelman             ptadenserow[ajj[k]]    = -1;
4323985e5eaSKris Buschelman             ptasparserow[ptanzi++] = ajj[k];
4333985e5eaSKris Buschelman           }
4343985e5eaSKris Buschelman         }
4353985e5eaSKris Buschelman       }
4363985e5eaSKris Buschelman       /* Using symbolic info for row of PtA, determine symbolic info for row of C: */
4373985e5eaSKris Buschelman       ptaj = ptasparserow;
4383985e5eaSKris Buschelman       cnzi   = 0;
4393985e5eaSKris Buschelman       for (j=0;j<ptanzi;j++) {
440fe05a634SKris Buschelman         pdof = *ptaj%dof;
4413985e5eaSKris Buschelman         prow = (*ptaj++)/dof;
4423985e5eaSKris Buschelman         pnzj = pi[prow+1] - pi[prow];
4433985e5eaSKris Buschelman         pjj  = pj + pi[prow];
4443985e5eaSKris Buschelman         for (k=0;k<pnzj;k++) {
445fe05a634SKris Buschelman           if (!denserow[pjj[k]+pdof]) {
446fe05a634SKris Buschelman             denserow[pjj[k]+pdof] = -1;
447fe05a634SKris Buschelman             sparserow[cnzi++]     = pjj[k]+pdof;
4483985e5eaSKris Buschelman           }
4493985e5eaSKris Buschelman         }
4503985e5eaSKris Buschelman       }
4513985e5eaSKris Buschelman 
4523985e5eaSKris Buschelman       /* sort sparserow */
4533985e5eaSKris Buschelman       ierr = PetscSortInt(cnzi,sparserow);CHKERRQ(ierr);
4543985e5eaSKris Buschelman 
4553985e5eaSKris Buschelman       /* If free space is not available, make more free space */
4563985e5eaSKris Buschelman       /* Double the amount of total space in the list */
4573985e5eaSKris Buschelman       if (current_space->local_remaining<cnzi) {
4583985e5eaSKris Buschelman         ierr = GetMoreSpace(current_space->total_array_size,&current_space);CHKERRQ(ierr);
4593985e5eaSKris Buschelman       }
4603985e5eaSKris Buschelman 
4613985e5eaSKris Buschelman       /* Copy data into free space, and zero out denserows */
4623985e5eaSKris Buschelman       ierr = PetscMemcpy(current_space->array,sparserow,cnzi*sizeof(int));CHKERRQ(ierr);
4633985e5eaSKris Buschelman       current_space->array           += cnzi;
4643985e5eaSKris Buschelman       current_space->local_used      += cnzi;
4653985e5eaSKris Buschelman       current_space->local_remaining -= cnzi;
4663985e5eaSKris Buschelman 
4673985e5eaSKris Buschelman       for (j=0;j<ptanzi;j++) {
4683985e5eaSKris Buschelman         ptadenserow[ptasparserow[j]] = 0;
4693985e5eaSKris Buschelman       }
4703985e5eaSKris Buschelman       for (j=0;j<cnzi;j++) {
4713985e5eaSKris Buschelman         denserow[sparserow[j]] = 0;
4723985e5eaSKris Buschelman       }
4733985e5eaSKris Buschelman       /* Aside: Perhaps we should save the pta info for the numerical factorization. */
4743985e5eaSKris Buschelman       /*        For now, we will recompute what is needed. */
4753985e5eaSKris Buschelman       ci[i+1+dof] = ci[i+dof] + cnzi;
4763985e5eaSKris Buschelman     }
4773985e5eaSKris Buschelman   }
4783985e5eaSKris Buschelman   /* nnz is now stored in ci[ptm], column indices are in the list of free space */
4793985e5eaSKris Buschelman   /* Allocate space for cj, initialize cj, and */
4803985e5eaSKris Buschelman   /* destroy list of free space and other temporary array(s) */
4813985e5eaSKris Buschelman   ierr = PetscMalloc((ci[pn]+1)*sizeof(int),&cj);CHKERRQ(ierr);
4823985e5eaSKris Buschelman   ierr = MakeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
4833985e5eaSKris Buschelman   ierr = PetscFree(ptadenserow);CHKERRQ(ierr);
4843985e5eaSKris Buschelman 
4853985e5eaSKris Buschelman   /* Allocate space for ca */
4863985e5eaSKris Buschelman   ierr = PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
4873985e5eaSKris Buschelman   ierr = PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));CHKERRQ(ierr);
4883985e5eaSKris Buschelman 
4893985e5eaSKris Buschelman   /* put together the new matrix */
4903985e5eaSKris Buschelman   ierr = MatCreateSeqAIJWithArrays(A->comm,pn,pn,ci,cj,ca,C);CHKERRQ(ierr);
4913985e5eaSKris Buschelman 
4923985e5eaSKris Buschelman   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
4933985e5eaSKris Buschelman   /* Since these are PETSc arrays, change flags to free them as necessary. */
4943985e5eaSKris Buschelman   c = (Mat_SeqAIJ *)((*C)->data);
4953985e5eaSKris Buschelman   c->freedata = PETSC_TRUE;
4963985e5eaSKris Buschelman   c->nonew    = 0;
4973985e5eaSKris Buschelman 
4983985e5eaSKris Buschelman   /* Clean up. */
4993985e5eaSKris Buschelman   ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
5003985e5eaSKris Buschelman 
5019af31e4aSHong Zhang   ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,PP,0,0);CHKERRQ(ierr);
5023985e5eaSKris Buschelman   PetscFunctionReturn(0);
5033985e5eaSKris Buschelman }
5043985e5eaSKris Buschelman EXTERN_C_END
5053985e5eaSKris Buschelman 
506eb9c0419SKris Buschelman #undef __FUNCT__
5079af31e4aSHong Zhang #define __FUNCT__ "MatPtAPNumeric"
5086849ba73SBarry Smith /*
5099af31e4aSHong Zhang    MatPtAPNumeric - Computes the matrix projection C = P^T * A * P
5104d3841fdSKris Buschelman 
5114d3841fdSKris Buschelman    Collective on Mat
5124d3841fdSKris Buschelman 
5134d3841fdSKris Buschelman    Input Parameters:
5144d3841fdSKris Buschelman +  A - the matrix
5154d3841fdSKris Buschelman -  P - the projection matrix
5164d3841fdSKris Buschelman 
5174d3841fdSKris Buschelman    Output Parameters:
5184d3841fdSKris Buschelman .  C - the product matrix
5194d3841fdSKris Buschelman 
5204d3841fdSKris Buschelman    Notes:
5219af31e4aSHong Zhang    C must have been created by calling MatPtAPSymbolic and must be destroyed by
5224d3841fdSKris Buschelman    the user using MatDeatroy().
5234d3841fdSKris Buschelman 
524170ef064SHong Zhang    This routine is currently only implemented for pairs of AIJ matrices and classes
525170ef064SHong Zhang    which inherit from AIJ.  C will be of type MATAIJ.
5264d3841fdSKris Buschelman 
5274d3841fdSKris Buschelman    Level: intermediate
5284d3841fdSKris Buschelman 
5299af31e4aSHong Zhang .seealso: MatPtAP(),MatPtAPSymbolic(),MatMatMultNumeric()
5306849ba73SBarry Smith */
531dfbe8321SBarry Smith PetscErrorCode MatPtAPNumeric(Mat A,Mat P,Mat C) {
532dfbe8321SBarry Smith   PetscErrorCode ierr;
533534c1384SKris Buschelman   PetscErrorCode (*fA)(Mat,Mat,Mat);
534534c1384SKris Buschelman   PetscErrorCode (*fP)(Mat,Mat,Mat);
535eb9c0419SKris Buschelman 
536eb9c0419SKris Buschelman   PetscFunctionBegin;
537eb9c0419SKris Buschelman 
5384482741eSBarry Smith   PetscValidHeaderSpecific(A,MAT_COOKIE,1);
539c9780b6fSBarry Smith   PetscValidType(A,1);
540eb9c0419SKris Buschelman   MatPreallocated(A);
541eb9c0419SKris Buschelman   if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
542eb9c0419SKris Buschelman   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
543eb9c0419SKris Buschelman 
5444482741eSBarry Smith   PetscValidHeaderSpecific(P,MAT_COOKIE,2);
545c9780b6fSBarry Smith   PetscValidType(P,2);
546eb9c0419SKris Buschelman   MatPreallocated(P);
547eb9c0419SKris Buschelman   if (!P->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
548eb9c0419SKris Buschelman   if (P->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
549eb9c0419SKris Buschelman 
5504482741eSBarry Smith   PetscValidHeaderSpecific(C,MAT_COOKIE,3);
551c9780b6fSBarry Smith   PetscValidType(C,3);
552eb9c0419SKris Buschelman   MatPreallocated(C);
553eb9c0419SKris Buschelman   if (!C->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
554eb9c0419SKris Buschelman   if (C->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
555eb9c0419SKris Buschelman 
556eb9c0419SKris Buschelman   if (P->N!=C->M) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->N,C->M);
557eb9c0419SKris Buschelman   if (P->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->M,A->N);
558eb9c0419SKris Buschelman   if (A->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix 'A' must be square, %d != %d",A->M,A->N);
559eb9c0419SKris Buschelman   if (P->N!=C->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->N,C->N);
560eb9c0419SKris Buschelman 
561534c1384SKris Buschelman   /* For now, we do not dispatch based on the type of A and P */
562534c1384SKris Buschelman   /* When implementations like _SeqAIJ_MAIJ exist, attack the multiple dispatch problem. */
563534c1384SKris Buschelman   fA = A->ops->ptapnumeric;
564534c1384SKris Buschelman   if (!fA) SETERRQ1(PETSC_ERR_SUP,"MatPtAPNumeric not supported for A of type %s",A->type_name);
565534c1384SKris Buschelman   fP = P->ops->ptapnumeric;
566534c1384SKris Buschelman   if (!fP) SETERRQ1(PETSC_ERR_SUP,"MatPtAPNumeric not supported for P of type %s",P->type_name);
567534c1384SKris 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);
5684d3841fdSKris Buschelman 
569534c1384SKris Buschelman   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
570534c1384SKris Buschelman   ierr = (*fA)(A,P,C);CHKERRQ(ierr);
571534c1384SKris Buschelman   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
572eb9c0419SKris Buschelman 
573eb9c0419SKris Buschelman   PetscFunctionReturn(0);
574eb9c0419SKris Buschelman }
575eb9c0419SKris Buschelman 
576eb9c0419SKris Buschelman #undef __FUNCT__
5779af31e4aSHong Zhang #define __FUNCT__ "MatPtAPNumeric_SeqAIJ_SeqAIJ"
578dfbe8321SBarry Smith PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat A,Mat P,Mat C)
579dfbe8321SBarry Smith {
580dfbe8321SBarry Smith   PetscErrorCode ierr;
581d20bfe6fSHong Zhang   int            flops=0;
582d20bfe6fSHong Zhang   Mat_SeqAIJ     *a  = (Mat_SeqAIJ *) A->data;
583d20bfe6fSHong Zhang   Mat_SeqAIJ     *p  = (Mat_SeqAIJ *) P->data;
584d20bfe6fSHong Zhang   Mat_SeqAIJ     *c  = (Mat_SeqAIJ *) C->data;
585d20bfe6fSHong Zhang   int            *ai=a->i,*aj=a->j,*apj,*apjdense,*pi=p->i,*pj=p->j,*pJ=p->j,*pjj;
586d20bfe6fSHong Zhang   int            *ci=c->i,*cj=c->j,*cjj;
587d20bfe6fSHong Zhang   int            am=A->M,cn=C->N,cm=C->M;
588d20bfe6fSHong Zhang   int            i,j,k,anzi,pnzi,apnzj,nextap,pnzj,prow,crow;
589d20bfe6fSHong Zhang   MatScalar      *aa=a->a,*apa,*pa=p->a,*pA=p->a,*paj,*ca=c->a,*caj;
590eb9c0419SKris Buschelman 
591eb9c0419SKris Buschelman   PetscFunctionBegin;
592d20bfe6fSHong Zhang   /* Allocate temporary array for storage of one row of A*P */
593d20bfe6fSHong Zhang   ierr = PetscMalloc(cn*(sizeof(MatScalar)+2*sizeof(int)),&apa);CHKERRQ(ierr);
594d20bfe6fSHong Zhang   ierr = PetscMemzero(apa,cn*(sizeof(MatScalar)+2*sizeof(int)));CHKERRQ(ierr);
595eb9c0419SKris Buschelman 
596d20bfe6fSHong Zhang   apj      = (int *)(apa + cn);
597d20bfe6fSHong Zhang   apjdense = apj + cn;
598d20bfe6fSHong Zhang 
599d20bfe6fSHong Zhang   /* Clear old values in C */
600d20bfe6fSHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
601d20bfe6fSHong Zhang 
602d20bfe6fSHong Zhang   for (i=0;i<am;i++) {
603d20bfe6fSHong Zhang     /* Form sparse row of A*P */
604d20bfe6fSHong Zhang     anzi  = ai[i+1] - ai[i];
605d20bfe6fSHong Zhang     apnzj = 0;
606d20bfe6fSHong Zhang     for (j=0;j<anzi;j++) {
607d20bfe6fSHong Zhang       prow = *aj++;
608d20bfe6fSHong Zhang       pnzj = pi[prow+1] - pi[prow];
609d20bfe6fSHong Zhang       pjj  = pj + pi[prow];
610d20bfe6fSHong Zhang       paj  = pa + pi[prow];
611d20bfe6fSHong Zhang       for (k=0;k<pnzj;k++) {
612d20bfe6fSHong Zhang         if (!apjdense[pjj[k]]) {
613d20bfe6fSHong Zhang           apjdense[pjj[k]] = -1;
614d20bfe6fSHong Zhang           apj[apnzj++]     = pjj[k];
615d20bfe6fSHong Zhang         }
616d20bfe6fSHong Zhang         apa[pjj[k]] += (*aa)*paj[k];
617d20bfe6fSHong Zhang       }
618d20bfe6fSHong Zhang       flops += 2*pnzj;
619d20bfe6fSHong Zhang       aa++;
620d20bfe6fSHong Zhang     }
621d20bfe6fSHong Zhang 
622d20bfe6fSHong Zhang     /* Sort the j index array for quick sparse axpy. */
623d20bfe6fSHong Zhang     ierr = PetscSortInt(apnzj,apj);CHKERRQ(ierr);
624d20bfe6fSHong Zhang 
625d20bfe6fSHong Zhang     /* Compute P^T*A*P using outer product (P^T)[:,j]*(A*P)[j,:]. */
626d20bfe6fSHong Zhang     pnzi = pi[i+1] - pi[i];
627d20bfe6fSHong Zhang     for (j=0;j<pnzi;j++) {
628d20bfe6fSHong Zhang       nextap = 0;
629d20bfe6fSHong Zhang       crow   = *pJ++;
630d20bfe6fSHong Zhang       cjj    = cj + ci[crow];
631d20bfe6fSHong Zhang       caj    = ca + ci[crow];
632d20bfe6fSHong Zhang       /* Perform sparse axpy operation.  Note cjj includes apj. */
633d20bfe6fSHong Zhang       for (k=0;nextap<apnzj;k++) {
634d20bfe6fSHong Zhang         if (cjj[k]==apj[nextap]) {
635d20bfe6fSHong Zhang           caj[k] += (*pA)*apa[apj[nextap++]];
636d20bfe6fSHong Zhang         }
637d20bfe6fSHong Zhang       }
638d20bfe6fSHong Zhang       flops += 2*apnzj;
639d20bfe6fSHong Zhang       pA++;
640d20bfe6fSHong Zhang     }
641d20bfe6fSHong Zhang 
642d20bfe6fSHong Zhang     /* Zero the current row info for A*P */
643d20bfe6fSHong Zhang     for (j=0;j<apnzj;j++) {
644d20bfe6fSHong Zhang       apa[apj[j]]      = 0.;
645d20bfe6fSHong Zhang       apjdense[apj[j]] = 0;
646d20bfe6fSHong Zhang     }
647d20bfe6fSHong Zhang   }
648d20bfe6fSHong Zhang 
649d20bfe6fSHong Zhang   /* Assemble the final matrix and clean up */
650d20bfe6fSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
651d20bfe6fSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
652d20bfe6fSHong Zhang   ierr = PetscFree(apa);CHKERRQ(ierr);
653d20bfe6fSHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
654d20bfe6fSHong Zhang 
655eb9c0419SKris Buschelman   PetscFunctionReturn(0);
656eb9c0419SKris Buschelman }
6570e36024fSHong Zhang 
6580e36024fSHong Zhang /* Compute C = P[rstart:rend,:]^T * A * P of seqaij matrices - used by MatPtAP_MPIAIJ_MPIAIJ() */
6590e36024fSHong Zhang 
6600e36024fSHong Zhang #undef __FUNCT__
6610e36024fSHong Zhang #define __FUNCT__ "MatPtAPNumeric_SeqAIJ_SeqAIJ_ReducedPt"
6620e36024fSHong Zhang PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_ReducedPt(Mat A,Mat P,int prstart,int prend,Mat C)
6630e36024fSHong Zhang {
6640e36024fSHong Zhang   PetscErrorCode ierr;
6650e36024fSHong Zhang   int            flops=0;
6660e36024fSHong Zhang   Mat_SeqAIJ     *a  = (Mat_SeqAIJ *) A->data;
6670e36024fSHong Zhang   Mat_SeqAIJ     *p  = (Mat_SeqAIJ *) P->data;
6680e36024fSHong Zhang   Mat_SeqAIJ     *c  = (Mat_SeqAIJ *) C->data;
6690e36024fSHong Zhang   int            *ai=a->i,*aj=a->j,*apj,*apjdense,*pi=p->i,*pj=p->j,*pJ=p->j,*pjj;
6700e36024fSHong Zhang   int            *ci=c->i,*cj=c->j,*cjj;
6710e36024fSHong Zhang   int            am=A->M,cn=C->N,cm=C->M;
6720e36024fSHong Zhang   int            i,j,k,anzi,pnzi,apnzj,nextap,pnzj,prow,crow;
6730e36024fSHong Zhang   MatScalar      *aa=a->a,*apa,*pa=p->a,*pA=p->a,*paj,*ca=c->a,*caj;
6740e36024fSHong Zhang 
6750e36024fSHong Zhang   PetscFunctionBegin;
6760e36024fSHong Zhang   pA=p->a+pi[prstart];
6770e36024fSHong Zhang   /* Allocate temporary array for storage of one row of A*P */
6780e36024fSHong Zhang   ierr = PetscMalloc(cn*(sizeof(MatScalar)+2*sizeof(int)),&apa);CHKERRQ(ierr);
6790e36024fSHong Zhang   ierr = PetscMemzero(apa,cn*(sizeof(MatScalar)+2*sizeof(int)));CHKERRQ(ierr);
6800e36024fSHong Zhang 
6810e36024fSHong Zhang   apj      = (int *)(apa + cn);
6820e36024fSHong Zhang   apjdense = apj + cn;
6830e36024fSHong Zhang 
6840e36024fSHong Zhang   /* Clear old values in C */
6850e36024fSHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
6860e36024fSHong Zhang 
6870e36024fSHong Zhang   for (i=0;i<am;i++) {
6880e36024fSHong Zhang     /* Form sparse row of A*P */
6890e36024fSHong Zhang     anzi  = ai[i+1] - ai[i];
6900e36024fSHong Zhang     apnzj = 0;
6910e36024fSHong Zhang     for (j=0;j<anzi;j++) {
6920e36024fSHong Zhang       prow = *aj++;
6930e36024fSHong Zhang       pnzj = pi[prow+1] - pi[prow];
6940e36024fSHong Zhang       pjj  = pj + pi[prow];
6950e36024fSHong Zhang       paj  = pa + pi[prow];
6960e36024fSHong Zhang       for (k=0;k<pnzj;k++) {
6970e36024fSHong Zhang         if (!apjdense[pjj[k]]) {
6980e36024fSHong Zhang           apjdense[pjj[k]] = -1;
6990e36024fSHong Zhang           apj[apnzj++]     = pjj[k];
7000e36024fSHong Zhang         }
7010e36024fSHong Zhang         apa[pjj[k]] += (*aa)*paj[k];
7020e36024fSHong Zhang       }
7030e36024fSHong Zhang       flops += 2*pnzj;
7040e36024fSHong Zhang       aa++;
7050e36024fSHong Zhang     }
7060e36024fSHong Zhang 
7070e36024fSHong Zhang     /* Sort the j index array for quick sparse axpy. */
7080e36024fSHong Zhang     ierr = PetscSortInt(apnzj,apj);CHKERRQ(ierr);
7090e36024fSHong Zhang 
7100e36024fSHong Zhang     /* Compute P[[prstart:prend,:]^T*A*P using outer product (P^T)[:,j+prstart]*(A*P)[j,:]. */
7110e36024fSHong Zhang     pnzi = pi[i+1+prstart] - pi[i+prstart];
7120e36024fSHong Zhang     for (j=0;j<pnzi;j++) {
7130e36024fSHong Zhang       nextap = 0;
7140e36024fSHong Zhang       crow   = *pJ++;
7150e36024fSHong Zhang       cjj    = cj + ci[crow];
7160e36024fSHong Zhang       caj    = ca + ci[crow];
7170e36024fSHong Zhang       /* Perform sparse axpy operation.  Note cjj includes apj. */
7180e36024fSHong Zhang       for (k=0;nextap<apnzj;k++) {
7190e36024fSHong Zhang         if (cjj[k]==apj[nextap]) {
7200e36024fSHong Zhang           caj[k] += (*pA)*apa[apj[nextap++]];
7210e36024fSHong Zhang         }
7220e36024fSHong Zhang       }
7230e36024fSHong Zhang       flops += 2*apnzj;
7240e36024fSHong Zhang       pA++;
7250e36024fSHong Zhang     }
7260e36024fSHong Zhang 
7270e36024fSHong Zhang     /* Zero the current row info for A*P */
7280e36024fSHong Zhang     for (j=0;j<apnzj;j++) {
7290e36024fSHong Zhang       apa[apj[j]]      = 0.;
7300e36024fSHong Zhang       apjdense[apj[j]] = 0;
7310e36024fSHong Zhang     }
7320e36024fSHong Zhang   }
7330e36024fSHong Zhang 
7340e36024fSHong Zhang   /* Assemble the final matrix and clean up */
7350e36024fSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7360e36024fSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7370e36024fSHong Zhang   ierr = PetscFree(apa);CHKERRQ(ierr);
7380e36024fSHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
7390e36024fSHong Zhang 
7400e36024fSHong Zhang   PetscFunctionReturn(0);
7410e36024fSHong Zhang }
7420e36024fSHong Zhang 
7430e36024fSHong Zhang #undef __FUNCT__
7440e36024fSHong Zhang #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqAIJ_ReducedPt"
7450e36024fSHong Zhang PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_ReducedPt(Mat A,Mat P,PetscReal fill,int prstart,int prend,Mat *C) {
7460e36024fSHong Zhang   PetscErrorCode ierr;
7470e36024fSHong Zhang   FreeSpaceList  free_space=PETSC_NULL,current_space=PETSC_NULL;
7480e36024fSHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*p=(Mat_SeqAIJ*)P->data,*c;
7490e36024fSHong Zhang   int            *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj;
7500e36024fSHong Zhang   int            *ci,*cj,*denserow,*sparserow,*ptadenserow,*ptasparserow,*ptaj;
7510e36024fSHong Zhang   int            an=A->N,am=A->M,pn=P->N,pm=P->M;
7520e36024fSHong Zhang   int            i,j,k,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi;
7530e36024fSHong Zhang   MatScalar      *ca;
7540e36024fSHong Zhang   Mat *psub,P_sub;
7550e36024fSHong Zhang   IS  isrow,iscol;
7560e36024fSHong Zhang   int m = prend - prstart;
757*0b89d903Svictorle 
758*0b89d903Svictorle   PetscFunctionBegin;
759*0b89d903Svictorle   /* Get ij structure of P[rstart:rend,:]^T */
7600e36024fSHong Zhang   ierr = ISCreateStride(PETSC_COMM_SELF,m,prstart,1,&isrow);CHKERRQ(ierr);
7610e36024fSHong Zhang   ierr = ISCreateStride(PETSC_COMM_SELF,P->n,0,1,&iscol);CHKERRQ(ierr);
7620e36024fSHong Zhang   ierr = MatGetSubMatrices(P,1,&isrow,&iscol,MAT_INITIAL_MATRIX,&psub);CHKERRQ(ierr);
7630e36024fSHong Zhang   ierr = ISDestroy(isrow);CHKERRQ(ierr);
7640e36024fSHong Zhang   ierr = ISDestroy(iscol);CHKERRQ(ierr);
7650e36024fSHong Zhang   P_sub = psub[0];
7660e36024fSHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(P_sub,&pti,&ptj);CHKERRQ(ierr);
7670e36024fSHong Zhang   ierr = MatDestroyMatrices(1,&psub);CHKERRQ(ierr);
7680e36024fSHong Zhang   ptJ=ptj;
7690e36024fSHong Zhang 
7700e36024fSHong Zhang   /* Allocate ci array, arrays for fill computation and */
7710e36024fSHong Zhang   /* free space for accumulating nonzero column info */
7720e36024fSHong Zhang   ierr = PetscMalloc((pn+1)*sizeof(int),&ci);CHKERRQ(ierr);
7730e36024fSHong Zhang   ci[0] = 0;
7740e36024fSHong Zhang 
7750e36024fSHong Zhang   ierr = PetscMalloc((2*pn+2*an+1)*sizeof(int),&ptadenserow);CHKERRQ(ierr);
7760e36024fSHong Zhang   ierr = PetscMemzero(ptadenserow,(2*pn+2*an+1)*sizeof(int));CHKERRQ(ierr);
7770e36024fSHong Zhang   ptasparserow = ptadenserow  + an;
7780e36024fSHong Zhang   denserow     = ptasparserow + an;
7790e36024fSHong Zhang   sparserow    = denserow     + pn;
7800e36024fSHong Zhang 
7810e36024fSHong Zhang   /* Set initial free space to be nnz(A) scaled by aspect ratio of P. */
7820e36024fSHong Zhang   /* This should be reasonable if sparsity of PtAP is similar to that of A. */
7830e36024fSHong Zhang   ierr          = GetMoreSpace((ai[am]/pm)*pn,&free_space);
7840e36024fSHong Zhang   current_space = free_space;
7850e36024fSHong Zhang 
7860e36024fSHong Zhang   /* Determine symbolic info for each row of C: */
7870e36024fSHong Zhang   for (i=0;i<pn;i++) {
7880e36024fSHong Zhang     ptnzi  = pti[i+1] - pti[i];
7890e36024fSHong Zhang     ptanzi = 0;
7900e36024fSHong Zhang     /* Determine symbolic row of PtA_reduced: */
7910e36024fSHong Zhang     for (j=0;j<ptnzi;j++) {
7920e36024fSHong Zhang       arow = *ptJ++;
7930e36024fSHong Zhang       anzj = ai[arow+1] - ai[arow];
7940e36024fSHong Zhang       ajj  = aj + ai[arow];
7950e36024fSHong Zhang       for (k=0;k<anzj;k++) {
7960e36024fSHong Zhang         if (!ptadenserow[ajj[k]]) {
7970e36024fSHong Zhang           ptadenserow[ajj[k]]    = -1;
7980e36024fSHong Zhang           ptasparserow[ptanzi++] = ajj[k];
7990e36024fSHong Zhang         }
8000e36024fSHong Zhang       }
8010e36024fSHong Zhang     }
8020e36024fSHong Zhang       /* Using symbolic info for row of PtA, determine symbolic info for row of C: */
8030e36024fSHong Zhang     ptaj = ptasparserow;
8040e36024fSHong Zhang     cnzi   = 0;
8050e36024fSHong Zhang     for (j=0;j<ptanzi;j++) {
8060e36024fSHong Zhang       prow = *ptaj++;
8070e36024fSHong Zhang       pnzj = pi[prow+1] - pi[prow];
8080e36024fSHong Zhang       pjj  = pj + pi[prow];
8090e36024fSHong Zhang       for (k=0;k<pnzj;k++) {
8100e36024fSHong Zhang         if (!denserow[pjj[k]]) {
8110e36024fSHong Zhang             denserow[pjj[k]]  = -1;
8120e36024fSHong Zhang             sparserow[cnzi++] = pjj[k];
8130e36024fSHong Zhang         }
8140e36024fSHong Zhang       }
8150e36024fSHong Zhang     }
8160e36024fSHong Zhang 
8170e36024fSHong Zhang     /* sort sparserow */
8180e36024fSHong Zhang     ierr = PetscSortInt(cnzi,sparserow);CHKERRQ(ierr);
8190e36024fSHong Zhang 
8200e36024fSHong Zhang     /* If free space is not available, make more free space */
8210e36024fSHong Zhang     /* Double the amount of total space in the list */
8220e36024fSHong Zhang     if (current_space->local_remaining<cnzi) {
8230e36024fSHong Zhang       ierr = GetMoreSpace(current_space->total_array_size,&current_space);CHKERRQ(ierr);
8240e36024fSHong Zhang     }
8250e36024fSHong Zhang 
8260e36024fSHong Zhang     /* Copy data into free space, and zero out denserows */
8270e36024fSHong Zhang     ierr = PetscMemcpy(current_space->array,sparserow,cnzi*sizeof(int));CHKERRQ(ierr);
8280e36024fSHong Zhang     current_space->array           += cnzi;
8290e36024fSHong Zhang     current_space->local_used      += cnzi;
8300e36024fSHong Zhang     current_space->local_remaining -= cnzi;
8310e36024fSHong Zhang 
8320e36024fSHong Zhang     for (j=0;j<ptanzi;j++) {
8330e36024fSHong Zhang       ptadenserow[ptasparserow[j]] = 0;
8340e36024fSHong Zhang     }
8350e36024fSHong Zhang     for (j=0;j<cnzi;j++) {
8360e36024fSHong Zhang       denserow[sparserow[j]] = 0;
8370e36024fSHong Zhang     }
8380e36024fSHong Zhang       /* Aside: Perhaps we should save the pta info for the numerical factorization. */
8390e36024fSHong Zhang       /*        For now, we will recompute what is needed. */
8400e36024fSHong Zhang     ci[i+1] = ci[i] + cnzi;
8410e36024fSHong Zhang   }
8420e36024fSHong Zhang   /* nnz is now stored in ci[ptm], column indices are in the list of free space */
8430e36024fSHong Zhang   /* Allocate space for cj, initialize cj, and */
8440e36024fSHong Zhang   /* destroy list of free space and other temporary array(s) */
8450e36024fSHong Zhang   ierr = PetscMalloc((ci[pn]+1)*sizeof(int),&cj);CHKERRQ(ierr);
8460e36024fSHong Zhang   ierr = MakeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
8470e36024fSHong Zhang   ierr = PetscFree(ptadenserow);CHKERRQ(ierr);
8480e36024fSHong Zhang 
8490e36024fSHong Zhang   /* Allocate space for ca */
8500e36024fSHong Zhang   ierr = PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
8510e36024fSHong Zhang   ierr = PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));CHKERRQ(ierr);
8520e36024fSHong Zhang 
8530e36024fSHong Zhang   /* put together the new matrix */
8540e36024fSHong Zhang   ierr = MatCreateSeqAIJWithArrays(A->comm,pn,pn,ci,cj,ca,C);CHKERRQ(ierr);
8550e36024fSHong Zhang 
8560e36024fSHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
8570e36024fSHong Zhang   /* Since these are PETSc arrays, change flags to free them as necessary. */
8580e36024fSHong Zhang   c = (Mat_SeqAIJ *)((*C)->data);
8590e36024fSHong Zhang   c->freedata = PETSC_TRUE;
8600e36024fSHong Zhang   c->nonew    = 0;
8610e36024fSHong Zhang 
8620e36024fSHong Zhang   /* Clean up. */
8630e36024fSHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
8640e36024fSHong Zhang 
8650e36024fSHong Zhang   PetscFunctionReturn(0);
8660e36024fSHong Zhang }
8670e36024fSHong Zhang 
8680e36024fSHong Zhang #undef __FUNCT__
8690e36024fSHong Zhang #define __FUNCT__ "MatPtAP_SeqAIJ_SeqAIJ_ReducedPt"
8700e36024fSHong Zhang PetscErrorCode MatPtAP_SeqAIJ_SeqAIJ_ReducedPt(Mat A,Mat P,MatReuse scall,PetscReal fill,int prstart,int prend,Mat *C)
8710e36024fSHong Zhang {
8720e36024fSHong Zhang   PetscErrorCode ierr;
8730e36024fSHong Zhang   PetscFunctionBegin;
8740e36024fSHong Zhang   if (A->m != prend-prstart) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",A->m,prend-prstart);
8750e36024fSHong 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);
8760e36024fSHong Zhang   if (scall == MAT_INITIAL_MATRIX){
8770e36024fSHong Zhang     ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ_ReducedPt(A,P,fill,prstart,prend,C);CHKERRQ(ierr);
8780e36024fSHong Zhang   }
8790e36024fSHong Zhang 
8800e36024fSHong Zhang   ierr = MatPtAPNumeric_SeqAIJ_SeqAIJ_ReducedPt(A,P,prstart,prend,*C);CHKERRQ(ierr);
8810e36024fSHong Zhang 
8820e36024fSHong Zhang   PetscFunctionReturn(0);
8830e36024fSHong Zhang }
8840e36024fSHong Zhang 
885