xref: /petsc/src/mat/impls/aij/seq/matptap.c (revision 0390132c8f796715fb0f5ee590874d59a39ac360)
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"
955a3bba9SHong Zhang #include "petscbt.h"
10eb9c0419SKris Buschelman 
11eb9c0419SKris Buschelman #undef __FUNCT__
129af31e4aSHong Zhang #define __FUNCT__ "MatPtAP"
134d3841fdSKris Buschelman /*@
149af31e4aSHong Zhang    MatPtAP - Creates the matrix projection C = P^T * A * P
154d3841fdSKris Buschelman 
164d3841fdSKris Buschelman    Collective on Mat
174d3841fdSKris Buschelman 
184d3841fdSKris Buschelman    Input Parameters:
194d3841fdSKris Buschelman +  A - the matrix
20f747e1aeSHong Zhang .  P - the projection matrix
21f747e1aeSHong Zhang .  scall - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX
22f747e1aeSHong Zhang -  fill - expected fill as ratio of nnz(C)/(nnz(A) + nnz(P))
234d3841fdSKris Buschelman 
244d3841fdSKris Buschelman    Output Parameters:
254d3841fdSKris Buschelman .  C - the product matrix
264d3841fdSKris Buschelman 
274d3841fdSKris Buschelman    Notes:
284d3841fdSKris Buschelman    C will be created and must be destroyed by the user with MatDestroy().
294d3841fdSKris Buschelman 
304d3841fdSKris Buschelman    This routine is currently only implemented for pairs of SeqAIJ matrices and classes
314d3841fdSKris Buschelman    which inherit from SeqAIJ.  C will be of type MATSEQAIJ.
324d3841fdSKris Buschelman 
334d3841fdSKris Buschelman    Level: intermediate
344d3841fdSKris Buschelman 
359af31e4aSHong Zhang .seealso: MatPtAPSymbolic(),MatPtAPNumeric(),MatMatMult()
364d3841fdSKris Buschelman @*/
3755a3bba9SHong Zhang PetscErrorCode MatPtAP(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
3855a3bba9SHong Zhang {
39dfbe8321SBarry Smith   PetscErrorCode ierr;
40534c1384SKris Buschelman   PetscErrorCode (*fA)(Mat,Mat,MatReuse,PetscReal,Mat *);
41534c1384SKris Buschelman   PetscErrorCode (*fP)(Mat,Mat,MatReuse,PetscReal,Mat *);
42eb9c0419SKris Buschelman 
43eb9c0419SKris Buschelman   PetscFunctionBegin;
449af31e4aSHong Zhang   PetscValidHeaderSpecific(A,MAT_COOKIE,1);
459af31e4aSHong Zhang   PetscValidType(A,1);
469af31e4aSHong Zhang   MatPreallocated(A);
479af31e4aSHong Zhang   if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
489af31e4aSHong Zhang   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
499af31e4aSHong Zhang   PetscValidHeaderSpecific(P,MAT_COOKIE,2);
509af31e4aSHong Zhang   PetscValidType(P,2);
519af31e4aSHong Zhang   MatPreallocated(P);
529af31e4aSHong Zhang   if (!P->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
539af31e4aSHong Zhang   if (P->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
549af31e4aSHong Zhang   PetscValidPointer(C,3);
55534c1384SKris Buschelman 
569af31e4aSHong Zhang   if (P->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->M,A->N);
57eb9c0419SKris Buschelman 
589af31e4aSHong Zhang   if (fill <=0.0) SETERRQ1(PETSC_ERR_ARG_SIZ,"fill=%g must be > 0.0",fill);
59eb9c0419SKris Buschelman 
60534c1384SKris Buschelman   /* For now, we do not dispatch based on the type of A and P */
61534c1384SKris Buschelman   /* When implementations like _SeqAIJ_MAIJ exist, attack the multiple dispatch problem. */
62534c1384SKris Buschelman   fA = A->ops->ptap;
63534c1384SKris Buschelman   if (!fA) SETERRQ1(PETSC_ERR_SUP,"MatPtAP not supported for A of type %s",A->type_name);
64534c1384SKris Buschelman   fP = P->ops->ptap;
65534c1384SKris Buschelman   if (!fP) SETERRQ1(PETSC_ERR_SUP,"MatPtAP not supported for P of type %s",P->type_name);
66534c1384SKris 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);
67534c1384SKris Buschelman 
689af31e4aSHong Zhang   ierr = PetscLogEventBegin(MAT_PtAP,A,P,0,0);CHKERRQ(ierr);
69534c1384SKris Buschelman   ierr = (*fA)(A,P,scall,fill,C);CHKERRQ(ierr);
709af31e4aSHong Zhang   ierr = PetscLogEventEnd(MAT_PtAP,A,P,0,0);CHKERRQ(ierr);
71eb9c0419SKris Buschelman   PetscFunctionReturn(0);
72eb9c0419SKris Buschelman }
73eb9c0419SKris Buschelman 
74*0390132cSHong Zhang EXTERN PetscErrorCode MatPtAPReducedPt_MPIAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,PetscInt,PetscInt,Mat*);
75*0390132cSHong Zhang EXTERN PetscErrorCode MatPtAPReducedPtSymbolic_MPIAIJ_SeqAIJ(Mat,Mat,PetscReal,PetscInt,PetscInt,Mat*);
76*0390132cSHong Zhang EXTERN PetscErrorCode MatPtAPReducedPtNumeric_MPIAIJ_SeqAIJ(Mat,Mat,PetscInt,PetscInt,Mat);
7722e3da61SHong Zhang 
78eb9c0419SKris Buschelman #undef __FUNCT__
79ff134f7aSHong Zhang #define __FUNCT__ "MatPtAP_MPIAIJ_MPIAIJ"
80ff134f7aSHong Zhang PetscErrorCode MatPtAP_MPIAIJ_MPIAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
81ff134f7aSHong Zhang {
82ff134f7aSHong Zhang   PetscErrorCode    ierr;
83b90dcfe3SHong Zhang 
84b90dcfe3SHong Zhang   PetscFunctionBegin;
85b90dcfe3SHong Zhang   if (scall == MAT_INITIAL_MATRIX){
864c627768SHong Zhang     ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
87b90dcfe3SHong Zhang     ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ(A,P,fill,C);CHKERRQ(ierr);/* numeric product is computed as well */
884c627768SHong Zhang     ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
89b90dcfe3SHong Zhang   } else if (scall == MAT_REUSE_MATRIX){
904c627768SHong Zhang     ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
91b90dcfe3SHong Zhang     ierr = MatPtAPNumeric_MPIAIJ_MPIAIJ(A,P,*C);CHKERRQ(ierr);
924c627768SHong Zhang     ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
93b90dcfe3SHong Zhang   } else {
94b90dcfe3SHong Zhang     SETERRQ1(PETSC_ERR_ARG_WRONG,"Invalid MatReuse %d",scall);
95b90dcfe3SHong Zhang   }
96b90dcfe3SHong Zhang   PetscFunctionReturn(0);
97b90dcfe3SHong Zhang }
98b90dcfe3SHong Zhang 
99b90dcfe3SHong Zhang #undef __FUNCT__
100b90dcfe3SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ"
101b90dcfe3SHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
102b90dcfe3SHong Zhang {
103b90dcfe3SHong Zhang   PetscErrorCode    ierr;
104*0390132cSHong Zhang   Mat               P_seq,C_seq;
105b1d57f15SBarry Smith   PetscInt          prstart,prend,m=P->m;
106ff134f7aSHong Zhang 
107ff134f7aSHong Zhang   PetscFunctionBegin;
10825616d81SHong Zhang   /* get P_seq = submatrix of P by taking rows of P that equal to nonzero col of A */
10922e3da61SHong Zhang   ierr = MatGetBrowsOfAcols(A,P,MAT_INITIAL_MATRIX,PETSC_NULL,PETSC_NULL,&prstart,&P_seq);CHKERRQ(ierr);
1100e36024fSHong Zhang 
11125616d81SHong Zhang   /* compute C_seq = P_loc^T * A_loc * P_seq */
112ff134f7aSHong Zhang   prend  = prstart + m;
113*0390132cSHong Zhang   ierr = MatPtAPReducedPt_MPIAIJ_SeqAIJ(A,P_seq,MAT_INITIAL_MATRIX,fill,prstart,prend,&C_seq);CHKERRQ(ierr);
114*0390132cSHong Zhang 
11525616d81SHong Zhang   ierr = MatDestroy(P_seq);CHKERRQ(ierr);
116b90dcfe3SHong Zhang 
117b90dcfe3SHong Zhang   /* add C_seq into mpi C */
118b90dcfe3SHong Zhang   ierr = MatMerge_SeqsToMPI(A->comm,C_seq,P->n,P->n,MAT_INITIAL_MATRIX,C);CHKERRQ(ierr);
119b90dcfe3SHong Zhang 
120ff134f7aSHong Zhang   PetscFunctionReturn(0);
121ff134f7aSHong Zhang }
122ff134f7aSHong Zhang 
123ff134f7aSHong Zhang #undef __FUNCT__
124ff134f7aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ"
125b90dcfe3SHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
126ff134f7aSHong Zhang {
127b90dcfe3SHong Zhang   PetscErrorCode       ierr;
128*0390132cSHong Zhang   Mat                  P_seq,C_seq;
129b1d57f15SBarry Smith   PetscInt             prstart,prend,m=P->m;
130671beff6SHong Zhang   Mat_Merge_SeqsToMPI  *merge;
131671beff6SHong Zhang   PetscObjectContainer container;
132ff134f7aSHong Zhang 
133ff134f7aSHong Zhang   PetscFunctionBegin;
134671beff6SHong Zhang   ierr = PetscObjectQuery((PetscObject)C,"MatMergeSeqsToMPI",(PetscObject *)&container);CHKERRQ(ierr);
135671beff6SHong Zhang   if (container) {
1367f79fd58SMatthew Knepley     ierr  = PetscObjectContainerGetPointer(container,(void **)&merge);CHKERRQ(ierr);
1377f79fd58SMatthew Knepley   } else {
1387f79fd58SMatthew Knepley     SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Matrix C does not posses and object container");
139671beff6SHong Zhang   }
140671beff6SHong Zhang 
141b90dcfe3SHong Zhang   /* get P_seq = submatrix of P by taking rows of P that equal to nonzero col of A */
142*0390132cSHong Zhang   ierr = MatGetBrowsOfAcols(A,P,MAT_INITIAL_MATRIX,PETSC_NULL,PETSC_NULL,&prstart,&P_seq);CHKERRQ(ierr);
143ff134f7aSHong Zhang 
144b90dcfe3SHong Zhang   /* compute C_seq = P_loc^T * A_loc * P_seq */
145b90dcfe3SHong Zhang   prend = prstart + m;
146b90dcfe3SHong Zhang   C_seq = merge->C_seq;
147*0390132cSHong Zhang   ierr = MatPtAPReducedPt_MPIAIJ_SeqAIJ(A,P_seq,MAT_REUSE_MATRIX,1.0,prstart,prend,&C_seq);CHKERRQ(ierr);
148b90dcfe3SHong Zhang   ierr = MatDestroy(P_seq);CHKERRQ(ierr);
149b90dcfe3SHong Zhang 
150b90dcfe3SHong Zhang   /* add C_seq into mpi C */
151b90dcfe3SHong Zhang   ierr = MatMerge_SeqsToMPI(A->comm,C_seq,P->n,P->n,MAT_REUSE_MATRIX,&C);CHKERRQ(ierr);
152b90dcfe3SHong Zhang 
153ff134f7aSHong Zhang   PetscFunctionReturn(0);
154ff134f7aSHong Zhang }
155ff134f7aSHong Zhang 
156ff134f7aSHong Zhang #undef __FUNCT__
1579af31e4aSHong Zhang #define __FUNCT__ "MatPtAP_SeqAIJ_SeqAIJ"
158dfbe8321SBarry Smith PetscErrorCode MatPtAP_SeqAIJ_SeqAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
1599af31e4aSHong Zhang {
160dfbe8321SBarry Smith   PetscErrorCode ierr;
161b1d57f15SBarry Smith 
1629af31e4aSHong Zhang   PetscFunctionBegin;
1639af31e4aSHong Zhang   if (scall == MAT_INITIAL_MATRIX){
164d20bfe6fSHong Zhang     ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
1659af31e4aSHong Zhang     ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ(A,P,fill,C);CHKERRQ(ierr);
166d20bfe6fSHong Zhang     ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
1679af31e4aSHong Zhang   }
168d20bfe6fSHong Zhang   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
1699af31e4aSHong Zhang   ierr = MatPtAPNumeric_SeqAIJ_SeqAIJ(A,P,*C);CHKERRQ(ierr);
170d20bfe6fSHong Zhang   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
1719af31e4aSHong Zhang   PetscFunctionReturn(0);
1729af31e4aSHong Zhang }
1739af31e4aSHong Zhang 
1749af31e4aSHong Zhang #undef __FUNCT__
1759af31e4aSHong Zhang #define __FUNCT__ "MatPtAPSymbolic"
1766849ba73SBarry Smith /*
1779af31e4aSHong Zhang    MatPtAPSymbolic - Creates the (i,j) structure of the matrix projection C = P^T * A * P
1784d3841fdSKris Buschelman 
1794d3841fdSKris Buschelman    Collective on Mat
1804d3841fdSKris Buschelman 
1814d3841fdSKris Buschelman    Input Parameters:
1824d3841fdSKris Buschelman +  A - the matrix
1834d3841fdSKris Buschelman -  P - the projection matrix
1844d3841fdSKris Buschelman 
1854d3841fdSKris Buschelman    Output Parameters:
1864d3841fdSKris Buschelman .  C - the (i,j) structure of the product matrix
1874d3841fdSKris Buschelman 
1884d3841fdSKris Buschelman    Notes:
1894d3841fdSKris Buschelman    C will be created and must be destroyed by the user with MatDestroy().
1904d3841fdSKris Buschelman 
1914d3841fdSKris Buschelman    This routine is currently only implemented for pairs of SeqAIJ matrices and classes
1924d3841fdSKris Buschelman    which inherit from SeqAIJ.  C will be of type MATSEQAIJ.  The product is computed using
1939af31e4aSHong Zhang    this (i,j) structure by calling MatPtAPNumeric().
1944d3841fdSKris Buschelman 
1954d3841fdSKris Buschelman    Level: intermediate
1964d3841fdSKris Buschelman 
1979af31e4aSHong Zhang .seealso: MatPtAP(),MatPtAPNumeric(),MatMatMultSymbolic()
1986849ba73SBarry Smith */
19955a3bba9SHong Zhang PetscErrorCode MatPtAPSymbolic(Mat A,Mat P,PetscReal fill,Mat *C)
20055a3bba9SHong Zhang {
201dfbe8321SBarry Smith   PetscErrorCode ierr;
202534c1384SKris Buschelman   PetscErrorCode (*fA)(Mat,Mat,PetscReal,Mat*);
203534c1384SKris Buschelman   PetscErrorCode (*fP)(Mat,Mat,PetscReal,Mat*);
204eb9c0419SKris Buschelman 
205eb9c0419SKris Buschelman   PetscFunctionBegin;
206eb9c0419SKris Buschelman 
2074482741eSBarry Smith   PetscValidHeaderSpecific(A,MAT_COOKIE,1);
208c9780b6fSBarry Smith   PetscValidType(A,1);
209eb9c0419SKris Buschelman   MatPreallocated(A);
210eb9c0419SKris Buschelman   if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
211eb9c0419SKris Buschelman   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
212eb9c0419SKris Buschelman 
2134482741eSBarry Smith   PetscValidHeaderSpecific(P,MAT_COOKIE,2);
214c9780b6fSBarry Smith   PetscValidType(P,2);
215eb9c0419SKris Buschelman   MatPreallocated(P);
216eb9c0419SKris Buschelman   if (!P->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
217eb9c0419SKris Buschelman   if (P->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
218eb9c0419SKris Buschelman 
2194482741eSBarry Smith   PetscValidPointer(C,3);
2204482741eSBarry Smith 
221eb9c0419SKris Buschelman   if (P->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->M,A->N);
222eb9c0419SKris Buschelman   if (A->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix 'A' must be square, %d != %d",A->M,A->N);
223eb9c0419SKris Buschelman 
224534c1384SKris Buschelman   /* For now, we do not dispatch based on the type of A and P */
225534c1384SKris Buschelman   /* When implementations like _SeqAIJ_MAIJ exist, attack the multiple dispatch problem. */
226534c1384SKris Buschelman   fA = A->ops->ptapsymbolic;
227534c1384SKris Buschelman   if (!fA) SETERRQ1(PETSC_ERR_SUP,"MatPtAPSymbolic not supported for A of type %s",A->type_name);
228534c1384SKris Buschelman   fP = P->ops->ptapsymbolic;
229534c1384SKris Buschelman   if (!fP) SETERRQ1(PETSC_ERR_SUP,"MatPtAPSymbolic not supported for P of type %s",P->type_name);
230534c1384SKris 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);
2314d3841fdSKris Buschelman 
232534c1384SKris Buschelman   ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
233534c1384SKris Buschelman   ierr = (*fA)(A,P,fill,C);CHKERRQ(ierr);
234534c1384SKris Buschelman   ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
235eb9c0419SKris Buschelman 
236eb9c0419SKris Buschelman   PetscFunctionReturn(0);
237eb9c0419SKris Buschelman }
238eb9c0419SKris Buschelman 
239f747e1aeSHong Zhang typedef struct {
240f747e1aeSHong Zhang   Mat    symAP;
241f747e1aeSHong Zhang } Mat_PtAPstruct;
242f747e1aeSHong Zhang 
24378a80504SBarry Smith EXTERN PetscErrorCode MatDestroy_SeqAIJ(Mat);
24478a80504SBarry Smith 
245f747e1aeSHong Zhang #undef __FUNCT__
246f747e1aeSHong Zhang #define __FUNCT__ "MatDestroy_SeqAIJ_PtAP"
247f4a850bbSBarry Smith PetscErrorCode MatDestroy_SeqAIJ_PtAP(Mat A)
248f747e1aeSHong Zhang {
249f4a850bbSBarry Smith   PetscErrorCode    ierr;
250f747e1aeSHong Zhang   Mat_PtAPstruct    *ptap=(Mat_PtAPstruct*)A->spptr;
251f747e1aeSHong Zhang 
252f747e1aeSHong Zhang   PetscFunctionBegin;
253f747e1aeSHong Zhang   ierr = MatDestroy(ptap->symAP);CHKERRQ(ierr);
254f747e1aeSHong Zhang   ierr = PetscFree(ptap);CHKERRQ(ierr);
25578a80504SBarry Smith   ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr);
256f747e1aeSHong Zhang   PetscFunctionReturn(0);
257f747e1aeSHong Zhang }
258f747e1aeSHong Zhang 
259eb9c0419SKris Buschelman #undef __FUNCT__
2609af31e4aSHong Zhang #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqAIJ"
26155a3bba9SHong Zhang PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
26255a3bba9SHong Zhang {
263dfbe8321SBarry Smith   PetscErrorCode ierr;
264d20bfe6fSHong Zhang   FreeSpaceList  free_space=PETSC_NULL,current_space=PETSC_NULL;
265d20bfe6fSHong Zhang   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data,*p = (Mat_SeqAIJ*)P->data,*c;
26655a3bba9SHong Zhang   PetscInt       *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj;
267b8374ebeSBarry Smith   PetscInt       *ci,*cj,*ptadenserow,*ptasparserow,*ptaj;
26855a3bba9SHong Zhang   PetscInt       an=A->N,am=A->M,pn=P->N,pm=P->M;
269b8374ebeSBarry Smith   PetscInt       i,j,k,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi,nlnk,*lnk;
270d20bfe6fSHong Zhang   MatScalar      *ca;
27155a3bba9SHong Zhang   PetscBT        lnkbt;
272eb9c0419SKris Buschelman 
273eb9c0419SKris Buschelman   PetscFunctionBegin;
274d20bfe6fSHong Zhang   /* Get ij structure of P^T */
275eb9c0419SKris Buschelman   ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
276d20bfe6fSHong Zhang   ptJ=ptj;
277eb9c0419SKris Buschelman 
278d20bfe6fSHong Zhang   /* Allocate ci array, arrays for fill computation and */
279d20bfe6fSHong Zhang   /* free space for accumulating nonzero column info */
28055a3bba9SHong Zhang   ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
281d20bfe6fSHong Zhang   ci[0] = 0;
282eb9c0419SKris Buschelman 
28355a3bba9SHong Zhang   ierr = PetscMalloc((2*an+1)*sizeof(PetscInt),&ptadenserow);CHKERRQ(ierr);
28455a3bba9SHong Zhang   ierr = PetscMemzero(ptadenserow,(2*an+1)*sizeof(PetscInt));CHKERRQ(ierr);
285d20bfe6fSHong Zhang   ptasparserow = ptadenserow  + an;
28655a3bba9SHong Zhang 
28755a3bba9SHong Zhang   /* create and initialize a linked list */
28855a3bba9SHong Zhang   nlnk = pn+1;
28955a3bba9SHong Zhang   ierr = PetscLLCreate(pn,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr);
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];
31955a3bba9SHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
32055a3bba9SHong Zhang       ierr = PetscLLAdd(pnzj,pjj,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr);
32155a3bba9SHong Zhang       cnzi += nlnk;
322d20bfe6fSHong Zhang     }
323d20bfe6fSHong Zhang 
324d20bfe6fSHong Zhang     /* If free space is not available, make more free space */
325d20bfe6fSHong Zhang     /* Double the amount of total space in the list */
326d20bfe6fSHong Zhang     if (current_space->local_remaining<cnzi) {
327d20bfe6fSHong Zhang       ierr = GetMoreSpace(current_space->total_array_size,&current_space);CHKERRQ(ierr);
328d20bfe6fSHong Zhang     }
329d20bfe6fSHong Zhang 
330d20bfe6fSHong Zhang     /* Copy data into free space, and zero out denserows */
33155a3bba9SHong Zhang     ierr = PetscLLClean(pn,pn,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
332d20bfe6fSHong Zhang     current_space->array           += cnzi;
333d20bfe6fSHong Zhang     current_space->local_used      += cnzi;
334d20bfe6fSHong Zhang     current_space->local_remaining -= cnzi;
335d20bfe6fSHong Zhang 
336d20bfe6fSHong Zhang     for (j=0;j<ptanzi;j++) {
337d20bfe6fSHong Zhang       ptadenserow[ptasparserow[j]] = 0;
338d20bfe6fSHong Zhang     }
339d20bfe6fSHong Zhang     /* Aside: Perhaps we should save the pta info for the numerical factorization. */
340d20bfe6fSHong Zhang     /*        For now, we will recompute what is needed. */
341d20bfe6fSHong Zhang     ci[i+1] = ci[i] + cnzi;
342d20bfe6fSHong Zhang   }
343d20bfe6fSHong Zhang   /* nnz is now stored in ci[ptm], column indices are in the list of free space */
344d20bfe6fSHong Zhang   /* Allocate space for cj, initialize cj, and */
345d20bfe6fSHong Zhang   /* destroy list of free space and other temporary array(s) */
34655a3bba9SHong Zhang   ierr = PetscMalloc((ci[pn]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
347d20bfe6fSHong Zhang   ierr = MakeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
348d20bfe6fSHong Zhang   ierr = PetscFree(ptadenserow);CHKERRQ(ierr);
34955a3bba9SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
350d20bfe6fSHong Zhang 
351d20bfe6fSHong Zhang   /* Allocate space for ca */
352d20bfe6fSHong Zhang   ierr = PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
353d20bfe6fSHong Zhang   ierr = PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));CHKERRQ(ierr);
354d20bfe6fSHong Zhang 
355d20bfe6fSHong Zhang   /* put together the new matrix */
356d20bfe6fSHong Zhang   ierr = MatCreateSeqAIJWithArrays(A->comm,pn,pn,ci,cj,ca,C);CHKERRQ(ierr);
357d20bfe6fSHong Zhang 
358d20bfe6fSHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
359d20bfe6fSHong Zhang   /* Since these are PETSc arrays, change flags to free them as necessary. */
360d20bfe6fSHong Zhang   c = (Mat_SeqAIJ *)((*C)->data);
361d20bfe6fSHong Zhang   c->freedata = PETSC_TRUE;
362d20bfe6fSHong Zhang   c->nonew    = 0;
363d20bfe6fSHong Zhang 
364d20bfe6fSHong Zhang   /* Clean up. */
365d20bfe6fSHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
366eb9c0419SKris Buschelman 
367eb9c0419SKris Buschelman   PetscFunctionReturn(0);
368eb9c0419SKris Buschelman }
369eb9c0419SKris Buschelman 
3703985e5eaSKris Buschelman #include "src/mat/impls/maij/maij.h"
3713985e5eaSKris Buschelman EXTERN_C_BEGIN
3723985e5eaSKris Buschelman #undef __FUNCT__
3739af31e4aSHong Zhang #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqMAIJ"
37455a3bba9SHong Zhang PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqMAIJ(Mat A,Mat PP,Mat *C)
37555a3bba9SHong Zhang {
3765c66b693SKris Buschelman   /* This routine requires testing -- I don't think it works. */
377dfbe8321SBarry Smith   PetscErrorCode ierr;
3783985e5eaSKris Buschelman   FreeSpaceList  free_space=PETSC_NULL,current_space=PETSC_NULL;
3793985e5eaSKris Buschelman   Mat_SeqMAIJ    *pp=(Mat_SeqMAIJ*)PP->data;
3803985e5eaSKris Buschelman   Mat            P=pp->AIJ;
3813985e5eaSKris Buschelman   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*p=(Mat_SeqAIJ*)P->data,*c;
382b1d57f15SBarry Smith   PetscInt       *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj;
383b1d57f15SBarry Smith   PetscInt       *ci,*cj,*denserow,*sparserow,*ptadenserow,*ptasparserow,*ptaj;
384b1d57f15SBarry Smith   PetscInt       an=A->N,am=A->M,pn=P->N,pm=P->M,ppdof=pp->dof;
385b1d57f15SBarry Smith   PetscInt       i,j,k,dof,pdof,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi;
3863985e5eaSKris Buschelman   MatScalar      *ca;
3873985e5eaSKris Buschelman 
3883985e5eaSKris Buschelman   PetscFunctionBegin;
3893985e5eaSKris Buschelman   /* Start timer */
3909af31e4aSHong Zhang   ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,PP,0,0);CHKERRQ(ierr);
3913985e5eaSKris Buschelman 
3923985e5eaSKris Buschelman   /* Get ij structure of P^T */
3933985e5eaSKris Buschelman   ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
3943985e5eaSKris Buschelman 
3953985e5eaSKris Buschelman   /* Allocate ci array, arrays for fill computation and */
3963985e5eaSKris Buschelman   /* free space for accumulating nonzero column info */
397b1d57f15SBarry Smith   ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
3983985e5eaSKris Buschelman   ci[0] = 0;
3993985e5eaSKris Buschelman 
400b1d57f15SBarry Smith   ierr = PetscMalloc((2*pn+2*an+1)*sizeof(PetscInt),&ptadenserow);CHKERRQ(ierr);
401b1d57f15SBarry Smith   ierr = PetscMemzero(ptadenserow,(2*pn+2*an+1)*sizeof(PetscInt));CHKERRQ(ierr);
4023985e5eaSKris Buschelman   ptasparserow = ptadenserow  + an;
4033985e5eaSKris Buschelman   denserow     = ptasparserow + an;
4043985e5eaSKris Buschelman   sparserow    = denserow     + pn;
4053985e5eaSKris Buschelman 
4063985e5eaSKris Buschelman   /* Set initial free space to be nnz(A) scaled by aspect ratio of P. */
4073985e5eaSKris Buschelman   /* This should be reasonable if sparsity of PtAP is similar to that of A. */
4083985e5eaSKris Buschelman   ierr          = GetMoreSpace((ai[am]/pm)*pn,&free_space);
4093985e5eaSKris Buschelman   current_space = free_space;
4103985e5eaSKris Buschelman 
4113985e5eaSKris Buschelman   /* Determine symbolic info for each row of C: */
4123985e5eaSKris Buschelman   for (i=0;i<pn/ppdof;i++) {
4133985e5eaSKris Buschelman     ptnzi  = pti[i+1] - pti[i];
4143985e5eaSKris Buschelman     ptanzi = 0;
4153985e5eaSKris Buschelman     ptJ    = ptj + pti[i];
4163985e5eaSKris Buschelman     for (dof=0;dof<ppdof;dof++) {
4173985e5eaSKris Buschelman     /* Determine symbolic row of PtA: */
4183985e5eaSKris Buschelman       for (j=0;j<ptnzi;j++) {
4193985e5eaSKris Buschelman         arow = ptJ[j] + dof;
4203985e5eaSKris Buschelman         anzj = ai[arow+1] - ai[arow];
4213985e5eaSKris Buschelman         ajj  = aj + ai[arow];
4223985e5eaSKris Buschelman         for (k=0;k<anzj;k++) {
4233985e5eaSKris Buschelman           if (!ptadenserow[ajj[k]]) {
4243985e5eaSKris Buschelman             ptadenserow[ajj[k]]    = -1;
4253985e5eaSKris Buschelman             ptasparserow[ptanzi++] = ajj[k];
4263985e5eaSKris Buschelman           }
4273985e5eaSKris Buschelman         }
4283985e5eaSKris Buschelman       }
4293985e5eaSKris Buschelman       /* Using symbolic info for row of PtA, determine symbolic info for row of C: */
4303985e5eaSKris Buschelman       ptaj = ptasparserow;
4313985e5eaSKris Buschelman       cnzi   = 0;
4323985e5eaSKris Buschelman       for (j=0;j<ptanzi;j++) {
433fe05a634SKris Buschelman         pdof = *ptaj%dof;
4343985e5eaSKris Buschelman         prow = (*ptaj++)/dof;
4353985e5eaSKris Buschelman         pnzj = pi[prow+1] - pi[prow];
4363985e5eaSKris Buschelman         pjj  = pj + pi[prow];
4373985e5eaSKris Buschelman         for (k=0;k<pnzj;k++) {
438fe05a634SKris Buschelman           if (!denserow[pjj[k]+pdof]) {
439fe05a634SKris Buschelman             denserow[pjj[k]+pdof] = -1;
440fe05a634SKris Buschelman             sparserow[cnzi++]     = pjj[k]+pdof;
4413985e5eaSKris Buschelman           }
4423985e5eaSKris Buschelman         }
4433985e5eaSKris Buschelman       }
4443985e5eaSKris Buschelman 
4453985e5eaSKris Buschelman       /* sort sparserow */
4463985e5eaSKris Buschelman       ierr = PetscSortInt(cnzi,sparserow);CHKERRQ(ierr);
4473985e5eaSKris Buschelman 
4483985e5eaSKris Buschelman       /* If free space is not available, make more free space */
4493985e5eaSKris Buschelman       /* Double the amount of total space in the list */
4503985e5eaSKris Buschelman       if (current_space->local_remaining<cnzi) {
4513985e5eaSKris Buschelman         ierr = GetMoreSpace(current_space->total_array_size,&current_space);CHKERRQ(ierr);
4523985e5eaSKris Buschelman       }
4533985e5eaSKris Buschelman 
4543985e5eaSKris Buschelman       /* Copy data into free space, and zero out denserows */
455b1d57f15SBarry Smith       ierr = PetscMemcpy(current_space->array,sparserow,cnzi*sizeof(PetscInt));CHKERRQ(ierr);
4563985e5eaSKris Buschelman       current_space->array           += cnzi;
4573985e5eaSKris Buschelman       current_space->local_used      += cnzi;
4583985e5eaSKris Buschelman       current_space->local_remaining -= cnzi;
4593985e5eaSKris Buschelman 
4603985e5eaSKris Buschelman       for (j=0;j<ptanzi;j++) {
4613985e5eaSKris Buschelman         ptadenserow[ptasparserow[j]] = 0;
4623985e5eaSKris Buschelman       }
4633985e5eaSKris Buschelman       for (j=0;j<cnzi;j++) {
4643985e5eaSKris Buschelman         denserow[sparserow[j]] = 0;
4653985e5eaSKris Buschelman       }
4663985e5eaSKris Buschelman       /* Aside: Perhaps we should save the pta info for the numerical factorization. */
4673985e5eaSKris Buschelman       /*        For now, we will recompute what is needed. */
4683985e5eaSKris Buschelman       ci[i+1+dof] = ci[i+dof] + cnzi;
4693985e5eaSKris Buschelman     }
4703985e5eaSKris Buschelman   }
4713985e5eaSKris Buschelman   /* nnz is now stored in ci[ptm], column indices are in the list of free space */
4723985e5eaSKris Buschelman   /* Allocate space for cj, initialize cj, and */
4733985e5eaSKris Buschelman   /* destroy list of free space and other temporary array(s) */
474b1d57f15SBarry Smith   ierr = PetscMalloc((ci[pn]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
4753985e5eaSKris Buschelman   ierr = MakeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
4763985e5eaSKris Buschelman   ierr = PetscFree(ptadenserow);CHKERRQ(ierr);
4773985e5eaSKris Buschelman 
4783985e5eaSKris Buschelman   /* Allocate space for ca */
4793985e5eaSKris Buschelman   ierr = PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
4803985e5eaSKris Buschelman   ierr = PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));CHKERRQ(ierr);
4813985e5eaSKris Buschelman 
4823985e5eaSKris Buschelman   /* put together the new matrix */
4833985e5eaSKris Buschelman   ierr = MatCreateSeqAIJWithArrays(A->comm,pn,pn,ci,cj,ca,C);CHKERRQ(ierr);
4843985e5eaSKris Buschelman 
4853985e5eaSKris Buschelman   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
4863985e5eaSKris Buschelman   /* Since these are PETSc arrays, change flags to free them as necessary. */
4873985e5eaSKris Buschelman   c = (Mat_SeqAIJ *)((*C)->data);
4883985e5eaSKris Buschelman   c->freedata = PETSC_TRUE;
4893985e5eaSKris Buschelman   c->nonew    = 0;
4903985e5eaSKris Buschelman 
4913985e5eaSKris Buschelman   /* Clean up. */
4923985e5eaSKris Buschelman   ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
4933985e5eaSKris Buschelman 
4949af31e4aSHong Zhang   ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,PP,0,0);CHKERRQ(ierr);
4953985e5eaSKris Buschelman   PetscFunctionReturn(0);
4963985e5eaSKris Buschelman }
4973985e5eaSKris Buschelman EXTERN_C_END
4983985e5eaSKris Buschelman 
499eb9c0419SKris Buschelman #undef __FUNCT__
5009af31e4aSHong Zhang #define __FUNCT__ "MatPtAPNumeric"
5016849ba73SBarry Smith /*
5029af31e4aSHong Zhang    MatPtAPNumeric - Computes the matrix projection C = P^T * A * P
5034d3841fdSKris Buschelman 
5044d3841fdSKris Buschelman    Collective on Mat
5054d3841fdSKris Buschelman 
5064d3841fdSKris Buschelman    Input Parameters:
5074d3841fdSKris Buschelman +  A - the matrix
5084d3841fdSKris Buschelman -  P - the projection matrix
5094d3841fdSKris Buschelman 
5104d3841fdSKris Buschelman    Output Parameters:
5114d3841fdSKris Buschelman .  C - the product matrix
5124d3841fdSKris Buschelman 
5134d3841fdSKris Buschelman    Notes:
5149af31e4aSHong Zhang    C must have been created by calling MatPtAPSymbolic and must be destroyed by
5154d3841fdSKris Buschelman    the user using MatDeatroy().
5164d3841fdSKris Buschelman 
517170ef064SHong Zhang    This routine is currently only implemented for pairs of AIJ matrices and classes
518170ef064SHong Zhang    which inherit from AIJ.  C will be of type MATAIJ.
5194d3841fdSKris Buschelman 
5204d3841fdSKris Buschelman    Level: intermediate
5214d3841fdSKris Buschelman 
5229af31e4aSHong Zhang .seealso: MatPtAP(),MatPtAPSymbolic(),MatMatMultNumeric()
5236849ba73SBarry Smith */
52455a3bba9SHong Zhang PetscErrorCode MatPtAPNumeric(Mat A,Mat P,Mat C)
52555a3bba9SHong Zhang {
526dfbe8321SBarry Smith   PetscErrorCode ierr;
527534c1384SKris Buschelman   PetscErrorCode (*fA)(Mat,Mat,Mat);
528534c1384SKris Buschelman   PetscErrorCode (*fP)(Mat,Mat,Mat);
529eb9c0419SKris Buschelman 
530eb9c0419SKris Buschelman   PetscFunctionBegin;
531eb9c0419SKris Buschelman 
5324482741eSBarry Smith   PetscValidHeaderSpecific(A,MAT_COOKIE,1);
533c9780b6fSBarry Smith   PetscValidType(A,1);
534eb9c0419SKris Buschelman   MatPreallocated(A);
535eb9c0419SKris Buschelman   if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
536eb9c0419SKris Buschelman   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
537eb9c0419SKris Buschelman 
5384482741eSBarry Smith   PetscValidHeaderSpecific(P,MAT_COOKIE,2);
539c9780b6fSBarry Smith   PetscValidType(P,2);
540eb9c0419SKris Buschelman   MatPreallocated(P);
541eb9c0419SKris Buschelman   if (!P->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
542eb9c0419SKris Buschelman   if (P->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
543eb9c0419SKris Buschelman 
5444482741eSBarry Smith   PetscValidHeaderSpecific(C,MAT_COOKIE,3);
545c9780b6fSBarry Smith   PetscValidType(C,3);
546eb9c0419SKris Buschelman   MatPreallocated(C);
547eb9c0419SKris Buschelman   if (!C->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
548eb9c0419SKris Buschelman   if (C->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
549eb9c0419SKris Buschelman 
550eb9c0419SKris Buschelman   if (P->N!=C->M) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->N,C->M);
551eb9c0419SKris Buschelman   if (P->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->M,A->N);
552eb9c0419SKris Buschelman   if (A->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix 'A' must be square, %d != %d",A->M,A->N);
553eb9c0419SKris Buschelman   if (P->N!=C->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->N,C->N);
554eb9c0419SKris Buschelman 
555534c1384SKris Buschelman   /* For now, we do not dispatch based on the type of A and P */
556534c1384SKris Buschelman   /* When implementations like _SeqAIJ_MAIJ exist, attack the multiple dispatch problem. */
557534c1384SKris Buschelman   fA = A->ops->ptapnumeric;
558534c1384SKris Buschelman   if (!fA) SETERRQ1(PETSC_ERR_SUP,"MatPtAPNumeric not supported for A of type %s",A->type_name);
559534c1384SKris Buschelman   fP = P->ops->ptapnumeric;
560534c1384SKris Buschelman   if (!fP) SETERRQ1(PETSC_ERR_SUP,"MatPtAPNumeric not supported for P of type %s",P->type_name);
561534c1384SKris Buschelman   if (fP!=fA) SETERRQ2(PETSC_ERR_ARG_INCOMP,"MatPtAPNumeric requires A, %s, to be compatible with P, %s",A->type_name,P->type_name);
5624d3841fdSKris Buschelman 
563534c1384SKris Buschelman   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
564534c1384SKris Buschelman   ierr = (*fA)(A,P,C);CHKERRQ(ierr);
565534c1384SKris Buschelman   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
566eb9c0419SKris Buschelman 
567eb9c0419SKris Buschelman   PetscFunctionReturn(0);
568eb9c0419SKris Buschelman }
569eb9c0419SKris Buschelman 
570eb9c0419SKris Buschelman #undef __FUNCT__
5719af31e4aSHong Zhang #define __FUNCT__ "MatPtAPNumeric_SeqAIJ_SeqAIJ"
572dfbe8321SBarry Smith PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat A,Mat P,Mat C)
573dfbe8321SBarry Smith {
574dfbe8321SBarry Smith   PetscErrorCode ierr;
575b1d57f15SBarry Smith   PetscInt       flops=0;
576d20bfe6fSHong Zhang   Mat_SeqAIJ     *a  = (Mat_SeqAIJ *) A->data;
577d20bfe6fSHong Zhang   Mat_SeqAIJ     *p  = (Mat_SeqAIJ *) P->data;
578d20bfe6fSHong Zhang   Mat_SeqAIJ     *c  = (Mat_SeqAIJ *) C->data;
579b1d57f15SBarry Smith   PetscInt       *ai=a->i,*aj=a->j,*apj,*apjdense,*pi=p->i,*pj=p->j,*pJ=p->j,*pjj;
580b1d57f15SBarry Smith   PetscInt       *ci=c->i,*cj=c->j,*cjj;
581b1d57f15SBarry Smith   PetscInt       am=A->M,cn=C->N,cm=C->M;
582b1d57f15SBarry Smith   PetscInt       i,j,k,anzi,pnzi,apnzj,nextap,pnzj,prow,crow;
583d20bfe6fSHong Zhang   MatScalar      *aa=a->a,*apa,*pa=p->a,*pA=p->a,*paj,*ca=c->a,*caj;
584eb9c0419SKris Buschelman 
585eb9c0419SKris Buschelman   PetscFunctionBegin;
586d20bfe6fSHong Zhang   /* Allocate temporary array for storage of one row of A*P */
587b1d57f15SBarry Smith   ierr = PetscMalloc(cn*(sizeof(MatScalar)+2*sizeof(PetscInt)),&apa);CHKERRQ(ierr);
588b1d57f15SBarry Smith   ierr = PetscMemzero(apa,cn*(sizeof(MatScalar)+2*sizeof(PetscInt)));CHKERRQ(ierr);
589eb9c0419SKris Buschelman 
590b1d57f15SBarry Smith   apj      = (PetscInt *)(apa + cn);
591d20bfe6fSHong Zhang   apjdense = apj + cn;
592d20bfe6fSHong Zhang 
593d20bfe6fSHong Zhang   /* Clear old values in C */
594d20bfe6fSHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
595d20bfe6fSHong Zhang 
596d20bfe6fSHong Zhang   for (i=0;i<am;i++) {
597d20bfe6fSHong Zhang     /* Form sparse row of A*P */
598d20bfe6fSHong Zhang     anzi  = ai[i+1] - ai[i];
599d20bfe6fSHong Zhang     apnzj = 0;
600d20bfe6fSHong Zhang     for (j=0;j<anzi;j++) {
601d20bfe6fSHong Zhang       prow = *aj++;
602d20bfe6fSHong Zhang       pnzj = pi[prow+1] - pi[prow];
603d20bfe6fSHong Zhang       pjj  = pj + pi[prow];
604d20bfe6fSHong Zhang       paj  = pa + pi[prow];
605d20bfe6fSHong Zhang       for (k=0;k<pnzj;k++) {
606d20bfe6fSHong Zhang         if (!apjdense[pjj[k]]) {
607d20bfe6fSHong Zhang           apjdense[pjj[k]] = -1;
608d20bfe6fSHong Zhang           apj[apnzj++]     = pjj[k];
609d20bfe6fSHong Zhang         }
610d20bfe6fSHong Zhang         apa[pjj[k]] += (*aa)*paj[k];
611d20bfe6fSHong Zhang       }
612d20bfe6fSHong Zhang       flops += 2*pnzj;
613d20bfe6fSHong Zhang       aa++;
614d20bfe6fSHong Zhang     }
615d20bfe6fSHong Zhang 
616d20bfe6fSHong Zhang     /* Sort the j index array for quick sparse axpy. */
617d20bfe6fSHong Zhang     ierr = PetscSortInt(apnzj,apj);CHKERRQ(ierr);
618d20bfe6fSHong Zhang 
619d20bfe6fSHong Zhang     /* Compute P^T*A*P using outer product (P^T)[:,j]*(A*P)[j,:]. */
620d20bfe6fSHong Zhang     pnzi = pi[i+1] - pi[i];
621d20bfe6fSHong Zhang     for (j=0;j<pnzi;j++) {
622d20bfe6fSHong Zhang       nextap = 0;
623d20bfe6fSHong Zhang       crow   = *pJ++;
624d20bfe6fSHong Zhang       cjj    = cj + ci[crow];
625d20bfe6fSHong Zhang       caj    = ca + ci[crow];
626d20bfe6fSHong Zhang       /* Perform sparse axpy operation.  Note cjj includes apj. */
627d20bfe6fSHong Zhang       for (k=0;nextap<apnzj;k++) {
628d20bfe6fSHong Zhang         if (cjj[k]==apj[nextap]) {
629d20bfe6fSHong Zhang           caj[k] += (*pA)*apa[apj[nextap++]];
630d20bfe6fSHong Zhang         }
631d20bfe6fSHong Zhang       }
632d20bfe6fSHong Zhang       flops += 2*apnzj;
633d20bfe6fSHong Zhang       pA++;
634d20bfe6fSHong Zhang     }
635d20bfe6fSHong Zhang 
636d20bfe6fSHong Zhang     /* Zero the current row info for A*P */
637d20bfe6fSHong Zhang     for (j=0;j<apnzj;j++) {
638d20bfe6fSHong Zhang       apa[apj[j]]      = 0.;
639d20bfe6fSHong Zhang       apjdense[apj[j]] = 0;
640d20bfe6fSHong Zhang     }
641d20bfe6fSHong Zhang   }
642d20bfe6fSHong Zhang 
643d20bfe6fSHong Zhang   /* Assemble the final matrix and clean up */
644d20bfe6fSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
645d20bfe6fSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
646d20bfe6fSHong Zhang   ierr = PetscFree(apa);CHKERRQ(ierr);
647d20bfe6fSHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
648d20bfe6fSHong Zhang 
649eb9c0419SKris Buschelman   PetscFunctionReturn(0);
650eb9c0419SKris Buschelman }
6510e36024fSHong Zhang 
652*0390132cSHong Zhang /* Compute local C = P[rstart:rend,:]^T * A * P - used by MatPtAP_MPIAIJ_MPIAIJ() */
6530e36024fSHong Zhang 
6540e36024fSHong Zhang #undef __FUNCT__
6550e36024fSHong Zhang #define __FUNCT__ "MatPtAPNumeric_SeqAIJ_SeqAIJ_ReducedPt"
6567f79fd58SMatthew Knepley /*@C
657e9b43d0fSSatish Balay    MatPtAPNumeric_SeqAIJ_SeqAIJ_ReducedPt - Creates C = P[rstart:rend,:]^T * A * P of seqaij matrices,
658b90dcfe3SHong Zhang                 used by MatPtAP_MPIAIJ_MPIAIJ()
659b90dcfe3SHong Zhang 
660b90dcfe3SHong Zhang    Collective on Mat
661b90dcfe3SHong Zhang 
662b90dcfe3SHong Zhang    Input Parameters:
663*0390132cSHong Zhang +  A - the matrix in mpiaij format
664b90dcfe3SHong Zhang .  P - the projection matrix in seqaij format
665b90dcfe3SHong Zhang .  scall - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX
666b90dcfe3SHong Zhang .  fill - expected fill (not being used when scall==MAT_REUSE_MATRIX)
667b90dcfe3SHong Zhang .  prstart, prend - the starting and ending-1-th row of the matrix P to be used for transpose
668b90dcfe3SHong Zhang 
669b90dcfe3SHong Zhang    Output Parameters:
670b90dcfe3SHong Zhang .  C - the product matrix in seqaij format
671b90dcfe3SHong Zhang 
672b90dcfe3SHong Zhang    Level: developer
673b90dcfe3SHong Zhang 
674b90dcfe3SHong Zhang .seealso: MatPtAPSymbolic(),MatPtAPNumeric(),MatMatMult()
675b90dcfe3SHong Zhang @*/
6760e36024fSHong Zhang 
6770e36024fSHong Zhang #undef __FUNCT__
678*0390132cSHong Zhang #define __FUNCT__ "MatPtAPReducedPtNumeric_MPIAIJ_SeqAIJ"
679*0390132cSHong Zhang PetscErrorCode MatPtAPReducedPtNumeric_MPIAIJ_SeqAIJ(Mat A,Mat P,PetscInt prstart,PetscInt prend,Mat C)
68022e3da61SHong Zhang {
68122e3da61SHong Zhang   PetscErrorCode ierr;
68222e3da61SHong Zhang   PetscInt       flops=0;
68322e3da61SHong Zhang   Mat_MPIAIJ     *a=(Mat_MPIAIJ*)A->data;
68422e3da61SHong Zhang   Mat_SeqAIJ     *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
68522e3da61SHong Zhang   Mat_SeqAIJ     *p  = (Mat_SeqAIJ *) P->data;
68622e3da61SHong Zhang   Mat_SeqAIJ     *c  = (Mat_SeqAIJ *) C->data;
68722e3da61SHong Zhang   PetscInt       *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,*apj,*apjdense,cstart=a->cstart,cend=a->cend,col;
68822e3da61SHong Zhang   PetscInt       *pi=p->i,*pj=p->j,*pJ=p->j+pi[prstart],*pjj;
68922e3da61SHong Zhang   PetscInt       *ci=c->i,*cj=c->j,*cjj;
69022e3da61SHong Zhang   PetscInt       am=A->m,cn=C->N,cm=C->M;
69122e3da61SHong Zhang   PetscInt       i,j,k,anzi,pnzi,apnzj,nextap,pnzj,prow,crow;
69222e3da61SHong Zhang   MatScalar      *ada=ad->a,*aoa=ao->a,*apa,*pa=p->a,*pA=p->a,*paj,*ca=c->a,*caj;
69322e3da61SHong Zhang 
69422e3da61SHong Zhang   PetscFunctionBegin;
69522e3da61SHong Zhang   pA=p->a+pi[prstart];
69622e3da61SHong Zhang   /* Allocate temporary array for storage of one row of A*P */
69722e3da61SHong Zhang   ierr = PetscMalloc(cn*(sizeof(MatScalar)+2*sizeof(PetscInt)),&apa);CHKERRQ(ierr);
69822e3da61SHong Zhang   ierr = PetscMemzero(apa,cn*(sizeof(MatScalar)+2*sizeof(PetscInt)));CHKERRQ(ierr);
69922e3da61SHong Zhang 
70022e3da61SHong Zhang   apj      = (PetscInt *)(apa + cn);
70122e3da61SHong Zhang   apjdense = apj + cn;
70222e3da61SHong Zhang 
70322e3da61SHong Zhang   /* Clear old values in C */
70422e3da61SHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
70522e3da61SHong Zhang 
70622e3da61SHong Zhang   for (i=0;i<am;i++) {
70722e3da61SHong Zhang     /* Form sparse row of A*P */
70822e3da61SHong Zhang     /* diagonal portion of A */
70922e3da61SHong Zhang     anzi  = adi[i+1] - adi[i];
71022e3da61SHong Zhang     apnzj = 0;
71122e3da61SHong Zhang     for (j=0;j<anzi;j++) {
71222e3da61SHong Zhang       prow = *adj + prstart;
71322e3da61SHong Zhang       adj++;
71422e3da61SHong Zhang       pnzj = pi[prow+1] - pi[prow];
71522e3da61SHong Zhang       pjj  = pj + pi[prow];
71622e3da61SHong Zhang       paj  = pa + pi[prow];
71722e3da61SHong Zhang       for (k=0;k<pnzj;k++) {
71822e3da61SHong Zhang         if (!apjdense[pjj[k]]) {
71922e3da61SHong Zhang           apjdense[pjj[k]] = -1;
72022e3da61SHong Zhang           apj[apnzj++]     = pjj[k];
72122e3da61SHong Zhang         }
72222e3da61SHong Zhang         apa[pjj[k]] += (*ada)*paj[k];
72322e3da61SHong Zhang       }
72422e3da61SHong Zhang       flops += 2*pnzj;
72522e3da61SHong Zhang       ada++;
72622e3da61SHong Zhang     }
72722e3da61SHong Zhang     /* off-diagonal portion of A */
72822e3da61SHong Zhang     anzi  = aoi[i+1] - aoi[i];
72922e3da61SHong Zhang     for (j=0;j<anzi;j++) {
73022e3da61SHong Zhang       col = a->garray[*aoj];
73122e3da61SHong Zhang       if (col < cstart){
73222e3da61SHong Zhang         prow = *aoj;
73322e3da61SHong Zhang       } else if (col >= cend){
73422e3da61SHong Zhang         prow = *aoj + prend-prstart;
73522e3da61SHong Zhang       } else {
73622e3da61SHong Zhang         SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Off-diagonal portion of A has wrong column map");
73722e3da61SHong Zhang       }
73822e3da61SHong Zhang       aoj++;
73922e3da61SHong Zhang       pnzj = pi[prow+1] - pi[prow];
74022e3da61SHong Zhang       pjj  = pj + pi[prow];
74122e3da61SHong Zhang       paj  = pa + pi[prow];
74222e3da61SHong Zhang       for (k=0;k<pnzj;k++) {
74322e3da61SHong Zhang         if (!apjdense[pjj[k]]) {
74422e3da61SHong Zhang           apjdense[pjj[k]] = -1;
74522e3da61SHong Zhang           apj[apnzj++]     = pjj[k];
74622e3da61SHong Zhang         }
74722e3da61SHong Zhang         apa[pjj[k]] += (*aoa)*paj[k];
74822e3da61SHong Zhang       }
74922e3da61SHong Zhang       flops += 2*pnzj;
75022e3da61SHong Zhang       aoa++;
75122e3da61SHong Zhang     }
75222e3da61SHong Zhang 
75322e3da61SHong Zhang     /* Sort the j index array for quick sparse axpy. */
75422e3da61SHong Zhang     ierr = PetscSortInt(apnzj,apj);CHKERRQ(ierr);
75522e3da61SHong Zhang 
75622e3da61SHong Zhang     /* Compute P[[prstart:prend,:]^T*A*P using outer product (P^T)[:,j+prstart]*(A*P)[j,:]. */
75722e3da61SHong Zhang     pnzi = pi[i+1+prstart] - pi[i+prstart];
75822e3da61SHong Zhang     for (j=0;j<pnzi;j++) {
75922e3da61SHong Zhang       nextap = 0;
76022e3da61SHong Zhang       crow   = *pJ++;
76122e3da61SHong Zhang       cjj    = cj + ci[crow];
76222e3da61SHong Zhang       caj    = ca + ci[crow];
76322e3da61SHong Zhang       /* Perform sparse axpy operation.  Note cjj includes apj. */
76422e3da61SHong Zhang       for (k=0;nextap<apnzj;k++) {
76522e3da61SHong Zhang         if (cjj[k]==apj[nextap]) {
76622e3da61SHong Zhang           caj[k] += (*pA)*apa[apj[nextap++]];
76722e3da61SHong Zhang         }
76822e3da61SHong Zhang       }
76922e3da61SHong Zhang       flops += 2*apnzj;
77022e3da61SHong Zhang       pA++;
77122e3da61SHong Zhang     }
77222e3da61SHong Zhang 
77322e3da61SHong Zhang     /* Zero the current row info for A*P */
77422e3da61SHong Zhang     for (j=0;j<apnzj;j++) {
77522e3da61SHong Zhang       apa[apj[j]]      = 0.;
77622e3da61SHong Zhang       apjdense[apj[j]] = 0;
77722e3da61SHong Zhang     }
77822e3da61SHong Zhang   }
77922e3da61SHong Zhang 
78022e3da61SHong Zhang   /* Assemble the final matrix and clean up */
78122e3da61SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
78222e3da61SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
78322e3da61SHong Zhang   ierr = PetscFree(apa);CHKERRQ(ierr);
78422e3da61SHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
78522e3da61SHong Zhang 
78622e3da61SHong Zhang   PetscFunctionReturn(0);
78722e3da61SHong Zhang }
78822e3da61SHong Zhang 
78922e3da61SHong Zhang #undef __FUNCT__
790*0390132cSHong Zhang #define __FUNCT__ "MatPtAPReducedPtSymbolic_MPIAIJ_SeqAIJ"
791*0390132cSHong Zhang PetscErrorCode MatPtAPReducedPtSymbolic_MPIAIJ_SeqAIJ(Mat A,Mat P,PetscReal fill,PetscInt prstart,PetscInt prend,Mat *C)
79222e3da61SHong Zhang {
79322e3da61SHong Zhang   PetscErrorCode ierr;
79422e3da61SHong Zhang   FreeSpaceList  free_space=PETSC_NULL,current_space=PETSC_NULL;
79522e3da61SHong Zhang   Mat_MPIAIJ     *a=(Mat_MPIAIJ*)A->data;
79622e3da61SHong Zhang   Mat_SeqAIJ     *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
79722e3da61SHong Zhang   Mat_SeqAIJ     *p=(Mat_SeqAIJ*)P->data,*c;
79822e3da61SHong Zhang   PetscInt       *pti,*ptj,*ptJ,*ajj,*pi=p->i,*pj=p->j,*pjj;
79922e3da61SHong Zhang   PetscInt       *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,nnz,cstart=a->cstart,cend=a->cend,col;
80022e3da61SHong Zhang   PetscInt       *ci,*cj,*ptadenserow,*ptasparserow,*ptaj;
80122e3da61SHong Zhang   PetscInt       an=A->N,am=A->m,pn=P->N,pm=P->M;
80222e3da61SHong Zhang   PetscInt       i,j,k,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi;
80322e3da61SHong Zhang   PetscInt       m=prend-prstart,nlnk,*lnk;
80422e3da61SHong Zhang   MatScalar      *ca;
80522e3da61SHong Zhang   PetscBT        lnkbt;
80622e3da61SHong Zhang 
80722e3da61SHong Zhang   PetscFunctionBegin;
80822e3da61SHong Zhang 
80922e3da61SHong Zhang   /* Get ij structure of P[rstart:rend,:]^T */
810e462e02cSHong Zhang   ierr = MatGetSymbolicTransposeReduced_SeqAIJ(P,prstart,prend,&pti,&ptj);CHKERRQ(ierr);
81122e3da61SHong Zhang   ptJ=ptj;
81222e3da61SHong Zhang 
81322e3da61SHong Zhang   /* Allocate ci array, arrays for fill computation and */
81422e3da61SHong Zhang   /* free space for accumulating nonzero column info */
81522e3da61SHong Zhang   ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
81622e3da61SHong Zhang   ci[0] = 0;
81722e3da61SHong Zhang 
81822e3da61SHong Zhang   ierr = PetscMalloc((2*an+1)*sizeof(PetscInt),&ptadenserow);CHKERRQ(ierr);
81922e3da61SHong Zhang   ierr = PetscMemzero(ptadenserow,(2*an+1)*sizeof(PetscInt));CHKERRQ(ierr);
82022e3da61SHong Zhang   ptasparserow = ptadenserow  + an;
82122e3da61SHong Zhang 
82222e3da61SHong Zhang   /* create and initialize a linked list */
82322e3da61SHong Zhang   nlnk = pn+1;
82422e3da61SHong Zhang   ierr = PetscLLCreate(pn,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr);
82522e3da61SHong Zhang 
82622e3da61SHong Zhang   /* Set initial free space to be nnz(A) scaled by aspect ratio of P. */
82722e3da61SHong Zhang   /* This should be reasonable if sparsity of PtAP is similar to that of A. */
82822e3da61SHong Zhang   nnz           = adi[am] + aoi[am];
82922e3da61SHong Zhang   ierr          = GetMoreSpace((PetscInt)(fill*nnz/pm)*pn,&free_space);
83022e3da61SHong Zhang   current_space = free_space;
83122e3da61SHong Zhang 
83222e3da61SHong Zhang   /* Determine symbolic info for each row of C: */
83322e3da61SHong Zhang   for (i=0;i<pn;i++) {
83422e3da61SHong Zhang     ptnzi  = pti[i+1] - pti[i];
83522e3da61SHong Zhang     ptanzi = 0;
83622e3da61SHong Zhang     /* Determine symbolic row of PtA_reduced: */
83722e3da61SHong Zhang     for (j=0;j<ptnzi;j++) {
83822e3da61SHong Zhang       arow = *ptJ++;
83922e3da61SHong Zhang       /* diagonal portion of A */
84022e3da61SHong Zhang       anzj = adi[arow+1] - adi[arow];
84122e3da61SHong Zhang       ajj  = adj + adi[arow];
84222e3da61SHong Zhang       for (k=0;k<anzj;k++) {
84322e3da61SHong Zhang         col = ajj[k]+prstart;
84422e3da61SHong Zhang         if (!ptadenserow[col]) {
84522e3da61SHong Zhang           ptadenserow[col]    = -1;
84622e3da61SHong Zhang           ptasparserow[ptanzi++] = col;
84722e3da61SHong Zhang         }
84822e3da61SHong Zhang       }
84922e3da61SHong Zhang       /* off-diagonal portion of A */
85022e3da61SHong Zhang       anzj = aoi[arow+1] - aoi[arow];
85122e3da61SHong Zhang       ajj  = aoj + aoi[arow];
85222e3da61SHong Zhang       for (k=0;k<anzj;k++) {
85322e3da61SHong Zhang         col = a->garray[ajj[k]];  /* global col */
85422e3da61SHong Zhang         if (col < cstart){
85522e3da61SHong Zhang           col = ajj[k];
85622e3da61SHong Zhang         } else if (col >= cend){
85722e3da61SHong Zhang           col = ajj[k] + m;
85822e3da61SHong Zhang         } else {
85922e3da61SHong Zhang           SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Off-diagonal portion of A has wrong column map");
86022e3da61SHong Zhang         }
86122e3da61SHong Zhang         if (!ptadenserow[col]) {
86222e3da61SHong Zhang           ptadenserow[col]    = -1;
86322e3da61SHong Zhang           ptasparserow[ptanzi++] = col;
86422e3da61SHong Zhang         }
86522e3da61SHong Zhang       }
86622e3da61SHong Zhang     }
86722e3da61SHong Zhang     /* Using symbolic info for row of PtA, determine symbolic info for row of C: */
86822e3da61SHong Zhang     ptaj = ptasparserow;
86922e3da61SHong Zhang     cnzi   = 0;
87022e3da61SHong Zhang     for (j=0;j<ptanzi;j++) {
87122e3da61SHong Zhang       prow = *ptaj++;
87222e3da61SHong Zhang       pnzj = pi[prow+1] - pi[prow];
87322e3da61SHong Zhang       pjj  = pj + pi[prow];
87422e3da61SHong Zhang       /* if (rank==1) printf(" prow: %d, pnzj: %d\n",prow,pnzj); */
87522e3da61SHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
87622e3da61SHong Zhang       ierr = PetscLLAdd(pnzj,pjj,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr);
87722e3da61SHong Zhang       cnzi += nlnk;
87822e3da61SHong Zhang     }
87922e3da61SHong Zhang 
88022e3da61SHong Zhang     /* If free space is not available, make more free space */
88122e3da61SHong Zhang     /* Double the amount of total space in the list */
88222e3da61SHong Zhang     if (current_space->local_remaining<cnzi) {
88322e3da61SHong Zhang       ierr = GetMoreSpace(current_space->total_array_size,&current_space);CHKERRQ(ierr);
88422e3da61SHong Zhang     }
88522e3da61SHong Zhang 
88622e3da61SHong Zhang     /* Copy data into free space, and zero out denserows */
88722e3da61SHong Zhang     ierr = PetscLLClean(pn,pn,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
88822e3da61SHong Zhang     current_space->array           += cnzi;
88922e3da61SHong Zhang     current_space->local_used      += cnzi;
89022e3da61SHong Zhang     current_space->local_remaining -= cnzi;
89122e3da61SHong Zhang 
89222e3da61SHong Zhang     for (j=0;j<ptanzi;j++) {
89322e3da61SHong Zhang       ptadenserow[ptasparserow[j]] = 0;
89422e3da61SHong Zhang     }
89522e3da61SHong Zhang 
89622e3da61SHong Zhang     /* Aside: Perhaps we should save the pta info for the numerical factorization. */
89722e3da61SHong Zhang     /*        For now, we will recompute what is needed. */
89822e3da61SHong Zhang     ci[i+1] = ci[i] + cnzi;
89922e3da61SHong Zhang   }
90022e3da61SHong Zhang   /* nnz is now stored in ci[ptm], column indices are in the list of free space */
90122e3da61SHong Zhang   /* Allocate space for cj, initialize cj, and */
90222e3da61SHong Zhang   /* destroy list of free space and other temporary array(s) */
90322e3da61SHong Zhang   ierr = PetscMalloc((ci[pn]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
90422e3da61SHong Zhang   ierr = MakeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
90522e3da61SHong Zhang   ierr = PetscFree(ptadenserow);CHKERRQ(ierr);
90622e3da61SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
90722e3da61SHong Zhang 
90822e3da61SHong Zhang   /* Allocate space for ca */
90922e3da61SHong Zhang   ierr = PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
91022e3da61SHong Zhang   ierr = PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));CHKERRQ(ierr);
91122e3da61SHong Zhang 
91222e3da61SHong Zhang   /* put together the new matrix */
91322e3da61SHong Zhang   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,pn,pn,ci,cj,ca,C);CHKERRQ(ierr);
91422e3da61SHong Zhang 
91522e3da61SHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
91622e3da61SHong Zhang   /* Since these are PETSc arrays, change flags to free them as necessary. */
91722e3da61SHong Zhang   c = (Mat_SeqAIJ *)((*C)->data);
91822e3da61SHong Zhang   c->freedata = PETSC_TRUE;
91922e3da61SHong Zhang   c->nonew    = 0;
92022e3da61SHong Zhang 
92122e3da61SHong Zhang   /* Clean up. */
92222e3da61SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
92322e3da61SHong Zhang 
92422e3da61SHong Zhang   PetscFunctionReturn(0);
92522e3da61SHong Zhang }
92622e3da61SHong Zhang 
927*0390132cSHong Zhang static PetscEvent logkey_matptapreducedpt = 0;
92822e3da61SHong Zhang #undef __FUNCT__
929*0390132cSHong Zhang #define __FUNCT__ "MatPtAPReducedPt_MPIAIJ_SeqAIJ"
930*0390132cSHong Zhang PetscErrorCode MatPtAPReducedPt_MPIAIJ_SeqAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,PetscInt prstart,PetscInt prend,Mat *C)
93122e3da61SHong Zhang {
93222e3da61SHong Zhang   PetscErrorCode ierr;
93322e3da61SHong Zhang 
93422e3da61SHong Zhang   PetscFunctionBegin;
93522e3da61SHong Zhang   if (A->m != prend-prstart) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",A->m,prend-prstart);
93622e3da61SHong 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);
937*0390132cSHong Zhang   if (!logkey_matptapreducedpt) {
938*0390132cSHong Zhang     ierr = PetscLogEventRegister(&logkey_matptapreducedpt,"MatPtAP_ReducedPt",MAT_COOKIE);
939e462e02cSHong Zhang   }
940*0390132cSHong Zhang   PetscLogEventBegin(logkey_matptapreducedpt,A,P,0,0);CHKERRQ(ierr);
94122e3da61SHong Zhang 
94222e3da61SHong Zhang   if (scall == MAT_INITIAL_MATRIX){
943*0390132cSHong Zhang     ierr = MatPtAPReducedPtSymbolic_MPIAIJ_SeqAIJ(A,P,fill,prstart,prend,C);CHKERRQ(ierr);
94422e3da61SHong Zhang   }
945*0390132cSHong Zhang   ierr = MatPtAPReducedPtNumeric_MPIAIJ_SeqAIJ(A,P,prstart,prend,*C);CHKERRQ(ierr);
946*0390132cSHong Zhang   ierr = PetscLogEventEnd(logkey_matptapreducedpt,A,P,0,0);CHKERRQ(ierr);
94722e3da61SHong Zhang 
94822e3da61SHong Zhang   PetscFunctionReturn(0);
94922e3da61SHong Zhang }
95022e3da61SHong Zhang 
951