xref: /petsc/src/mat/impls/aij/seq/matptap.c (revision 8f8f2f0d3bbfe99bf6fe84d4337dd6d7e3a5b04f)
1be1d678aSKris Buschelman 
2eb9c0419SKris Buschelman /*
342c57489SHong Zhang   Defines projective product routines where A is a SeqAIJ matrix
4eb9c0419SKris Buschelman           C = P^T * A * P
5eb9c0419SKris Buschelman */
6eb9c0419SKris Buschelman 
7c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h>   /*I "petscmat.h" I*/
8c6db04a5SJed Brown #include <../src/mat/utils/freespace.h>
9c6db04a5SJed Brown #include <petscbt.h>
108563dfccSBarry Smith #include <petsctime.h>
11eb9c0419SKris Buschelman 
126f231fbdSstefano_zampini #if defined(PETSC_HAVE_HYPRE)
136f231fbdSstefano_zampini PETSC_INTERN PetscErrorCode MatPtAPSymbolic_AIJ_AIJ_wHYPRE(Mat,Mat,PetscReal,Mat*);
146f231fbdSstefano_zampini #endif
156f231fbdSstefano_zampini 
16150d2497SBarry Smith PETSC_INTERN PetscErrorCode MatPtAP_SeqAIJ_SeqAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
177ba1a0bfSKris Buschelman {
187ba1a0bfSKris Buschelman   PetscErrorCode      ierr;
196f231fbdSstefano_zampini #if !defined(PETSC_HAVE_HYPRE)
208fa4b5a6SHong Zhang   const char          *algTypes[2] = {"scalable","rap"};
216f231fbdSstefano_zampini   PetscInt            nalg = 2;
226f231fbdSstefano_zampini #else
238fa4b5a6SHong Zhang   const char          *algTypes[3] = {"scalable","rap","hypre"};
246f231fbdSstefano_zampini   PetscInt            nalg = 3;
256f231fbdSstefano_zampini #endif
268fa4b5a6SHong Zhang   PetscInt            alg = 1; /* set default algorithm */
278fa4b5a6SHong Zhang   Mat                 Pt;
288fa4b5a6SHong Zhang   Mat_MatTransMatMult *atb;
298fa4b5a6SHong Zhang   Mat_SeqAIJ          *c;
307ba1a0bfSKris Buschelman 
317ba1a0bfSKris Buschelman   PetscFunctionBegin;
3265e8a0caSHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
334a1b09b7SHong Zhang     /*
344a1b09b7SHong Zhang      Alg 'scalable' determines which implementations to be used:
358fa4b5a6SHong Zhang        "rap":      Pt = P^T and C = Pt*A*P
368fa4b5a6SHong Zhang        "scalable": do outer product and two sparse axpy in MatPtAPNumeric() - might slow, does not store structure of A*P.
376f231fbdSstefano_zampini        "hypre":    use boomerAMGBuildCoarseOperator.
384a1b09b7SHong Zhang      */
39715a5346SHong Zhang     ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MatPtAP","Mat");CHKERRQ(ierr);
406f231fbdSstefano_zampini     ierr = PetscOptionsEList("-matptap_via","Algorithmic approach","MatPtAP",algTypes,nalg,algTypes[0],&alg,NULL);CHKERRQ(ierr);
414a1b09b7SHong Zhang     ierr = PetscOptionsEnd();CHKERRQ(ierr);
424a1b09b7SHong Zhang     switch (alg) {
434a1b09b7SHong Zhang     case 1:
448fa4b5a6SHong Zhang       ierr = PetscNew(&atb);CHKERRQ(ierr);
458fa4b5a6SHong Zhang       ierr = MatTranspose_SeqAIJ(P,MAT_INITIAL_MATRIX,&Pt);CHKERRQ(ierr);
468fa4b5a6SHong Zhang       ierr = MatMatMatMult(Pt,A,P,MAT_INITIAL_MATRIX,fill,C);CHKERRQ(ierr);
478fa4b5a6SHong Zhang 
488fa4b5a6SHong Zhang       c                      = (Mat_SeqAIJ*)(*C)->data;
498fa4b5a6SHong Zhang       c->atb                 = atb;
508fa4b5a6SHong Zhang       atb->At                = Pt;
518fa4b5a6SHong Zhang       atb->destroy           = (*C)->ops->destroy;
528fa4b5a6SHong Zhang       (*C)->ops->destroy     = MatDestroy_SeqAIJ_MatTransMatMult;
538fa4b5a6SHong Zhang       (*C)->ops->ptapnumeric = MatPtAPNumeric_SeqAIJ_SeqAIJ;
548fa4b5a6SHong Zhang       PetscFunctionReturn(0);
554a1b09b7SHong Zhang       break;
566f231fbdSstefano_zampini #if defined(PETSC_HAVE_HYPRE)
576f231fbdSstefano_zampini     case 2:
580abfe76eSFande Kong       ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
596f231fbdSstefano_zampini       ierr = MatPtAPSymbolic_AIJ_AIJ_wHYPRE(A,P,fill,C);CHKERRQ(ierr);
600abfe76eSFande Kong       ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
616f231fbdSstefano_zampini       break;
626f231fbdSstefano_zampini #endif
634a1b09b7SHong Zhang     default:
640abfe76eSFande Kong       ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
654a1b09b7SHong Zhang       ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(A,P,fill,C);CHKERRQ(ierr);
660abfe76eSFande Kong       ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
674a1b09b7SHong Zhang       break;
684a1b09b7SHong Zhang     }
697ba1a0bfSKris Buschelman   }
700abfe76eSFande Kong   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
71ae32dd66SHong Zhang   ierr = (*(*C)->ops->ptapnumeric)(A,P,*C);CHKERRQ(ierr);
720abfe76eSFande Kong   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
737ba1a0bfSKris Buschelman   PetscFunctionReturn(0);
747ba1a0bfSKris Buschelman }
757ba1a0bfSKris Buschelman 
7653565b12SHong Zhang PetscErrorCode MatDestroy_SeqAIJ_PtAP(Mat A)
778d0a38e4SHong Zhang {
788d0a38e4SHong Zhang   PetscErrorCode ierr;
7953565b12SHong Zhang   Mat_SeqAIJ     *a  = (Mat_SeqAIJ*)A->data;
803cdca5ebSHong Zhang   Mat_AP         *ap = a->ap;
818d0a38e4SHong Zhang 
828d0a38e4SHong Zhang   PetscFunctionBegin;
833cdca5ebSHong Zhang   ierr = PetscFree(ap->apa);CHKERRQ(ierr);
843cdca5ebSHong Zhang   ierr = PetscFree(ap->api);CHKERRQ(ierr);
853cdca5ebSHong Zhang   ierr = PetscFree(ap->apj);CHKERRQ(ierr);
863cdca5ebSHong Zhang   ierr = (ap->destroy)(A);CHKERRQ(ierr);
873cdca5ebSHong Zhang   ierr = PetscFree(ap);CHKERRQ(ierr);
888d0a38e4SHong Zhang   PetscFunctionReturn(0);
898d0a38e4SHong Zhang }
908d0a38e4SHong Zhang 
91ae32dd66SHong Zhang PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(Mat A,Mat P,PetscReal fill,Mat *C)
928d0a38e4SHong Zhang {
938d0a38e4SHong Zhang   PetscErrorCode     ierr;
940298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
95d20bfe6fSHong Zhang   Mat_SeqAIJ         *a        = (Mat_SeqAIJ*)A->data,*p = (Mat_SeqAIJ*)P->data,*c;
9653565b12SHong Zhang   PetscInt           *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj;
9753565b12SHong Zhang   PetscInt           *ci,*cj,*ptadenserow,*ptasparserow,*ptaj,nspacedouble=0;
987d69fa63SHong Zhang   PetscInt           an=A->cmap->N,am=A->rmap->N,pn=P->cmap->N,pm=P->rmap->N;
9953565b12SHong Zhang   PetscInt           i,j,k,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi,nlnk,*lnk;
100d20bfe6fSHong Zhang   MatScalar          *ca;
10153565b12SHong Zhang   PetscBT            lnkbt;
1027d69fa63SHong Zhang   PetscReal          afill;
103eb9c0419SKris Buschelman 
104eb9c0419SKris Buschelman   PetscFunctionBegin;
10553565b12SHong Zhang   /* Get ij structure of P^T */
106eb9c0419SKris Buschelman   ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
107d20bfe6fSHong Zhang   ptJ  = ptj;
108eb9c0419SKris Buschelman 
109d20bfe6fSHong Zhang   /* Allocate ci array, arrays for fill computation and */
110d20bfe6fSHong Zhang   /* free space for accumulating nonzero column info */
111854ce69bSBarry Smith   ierr  = PetscMalloc1(pn+1,&ci);CHKERRQ(ierr);
112d20bfe6fSHong Zhang   ci[0] = 0;
113eb9c0419SKris Buschelman 
1141795a4d1SJed Brown   ierr         = PetscCalloc1(2*an+1,&ptadenserow);CHKERRQ(ierr);
115d20bfe6fSHong Zhang   ptasparserow = ptadenserow  + an;
11655a3bba9SHong Zhang 
11755a3bba9SHong Zhang   /* create and initialize a linked list */
11855a3bba9SHong Zhang   nlnk = pn+1;
11955a3bba9SHong Zhang   ierr = PetscLLCreate(pn,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr);
120eb9c0419SKris Buschelman 
1217d69fa63SHong Zhang   /* Set initial free space to be fill*(nnz(A)+ nnz(P)) */
122f91af8c7SBarry Smith   ierr          = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],pi[pm])),&free_space);CHKERRQ(ierr);
123d20bfe6fSHong Zhang   current_space = free_space;
124d20bfe6fSHong Zhang 
125d20bfe6fSHong Zhang   /* Determine symbolic info for each row of C: */
126d20bfe6fSHong Zhang   for (i=0; i<pn; i++) {
127d20bfe6fSHong Zhang     ptnzi  = pti[i+1] - pti[i];
128d20bfe6fSHong Zhang     ptanzi = 0;
129d20bfe6fSHong Zhang     /* Determine symbolic row of PtA: */
130d20bfe6fSHong Zhang     for (j=0; j<ptnzi; j++) {
131d20bfe6fSHong Zhang       arow = *ptJ++;
132d20bfe6fSHong Zhang       anzj = ai[arow+1] - ai[arow];
133d20bfe6fSHong Zhang       ajj  = aj + ai[arow];
134d20bfe6fSHong Zhang       for (k=0; k<anzj; k++) {
135d20bfe6fSHong Zhang         if (!ptadenserow[ajj[k]]) {
136d20bfe6fSHong Zhang           ptadenserow[ajj[k]]    = -1;
137d20bfe6fSHong Zhang           ptasparserow[ptanzi++] = ajj[k];
138d20bfe6fSHong Zhang         }
139d20bfe6fSHong Zhang       }
140d20bfe6fSHong Zhang     }
141d20bfe6fSHong Zhang     /* Using symbolic info for row of PtA, determine symbolic info for row of C: */
142d20bfe6fSHong Zhang     ptaj = ptasparserow;
143d20bfe6fSHong Zhang     cnzi = 0;
144d20bfe6fSHong Zhang     for (j=0; j<ptanzi; j++) {
145d20bfe6fSHong Zhang       prow = *ptaj++;
146d20bfe6fSHong Zhang       pnzj = pi[prow+1] - pi[prow];
147d20bfe6fSHong Zhang       pjj  = pj + pi[prow];
14855a3bba9SHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
149dadf0e6bSHong Zhang       ierr  = PetscLLAddSorted(pnzj,pjj,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr);
15055a3bba9SHong Zhang       cnzi += nlnk;
151d20bfe6fSHong Zhang     }
152d20bfe6fSHong Zhang 
153d20bfe6fSHong Zhang     /* If free space is not available, make more free space */
154d20bfe6fSHong Zhang     /* Double the amount of total space in the list */
155d20bfe6fSHong Zhang     if (current_space->local_remaining<cnzi) {
156f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(cnzi,current_space->total_array_size),&current_space);CHKERRQ(ierr);
157f2b054eeSHong Zhang       nspacedouble++;
158d20bfe6fSHong Zhang     }
159d20bfe6fSHong Zhang 
160d20bfe6fSHong Zhang     /* Copy data into free space, and zero out denserows */
16155a3bba9SHong Zhang     ierr = PetscLLClean(pn,pn,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
1622205254eSKarl Rupp 
163d20bfe6fSHong Zhang     current_space->array           += cnzi;
164d20bfe6fSHong Zhang     current_space->local_used      += cnzi;
165d20bfe6fSHong Zhang     current_space->local_remaining -= cnzi;
166d20bfe6fSHong Zhang 
1672205254eSKarl Rupp     for (j=0; j<ptanzi; j++) ptadenserow[ptasparserow[j]] = 0;
1682205254eSKarl Rupp 
169d20bfe6fSHong Zhang     /* Aside: Perhaps we should save the pta info for the numerical factorization. */
170d20bfe6fSHong Zhang     /*        For now, we will recompute what is needed. */
171d20bfe6fSHong Zhang     ci[i+1] = ci[i] + cnzi;
172d20bfe6fSHong Zhang   }
173d20bfe6fSHong Zhang   /* nnz is now stored in ci[ptm], column indices are in the list of free space */
174d20bfe6fSHong Zhang   /* Allocate space for cj, initialize cj, and */
175d20bfe6fSHong Zhang   /* destroy list of free space and other temporary array(s) */
176854ce69bSBarry Smith   ierr = PetscMalloc1(ci[pn]+1,&cj);CHKERRQ(ierr);
177a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
178d20bfe6fSHong Zhang   ierr = PetscFree(ptadenserow);CHKERRQ(ierr);
17955a3bba9SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
180d20bfe6fSHong Zhang 
181854ce69bSBarry Smith   ierr = PetscCalloc1(ci[pn]+1,&ca);CHKERRQ(ierr);
18253565b12SHong Zhang 
18353565b12SHong Zhang   /* put together the new matrix */
184ce94432eSBarry Smith   ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),pn,pn,ci,cj,ca,C);CHKERRQ(ierr);
18533d57670SJed Brown   ierr = MatSetBlockSizes(*C,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr);
186ae32dd66SHong Zhang 
18753565b12SHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
18853565b12SHong Zhang   /* Since these are PETSc arrays, change flags to free them as necessary. */
18953565b12SHong Zhang   c          = (Mat_SeqAIJ*)((*C)->data);
19053565b12SHong Zhang   c->free_a  = PETSC_TRUE;
19153565b12SHong Zhang   c->free_ij = PETSC_TRUE;
19253565b12SHong Zhang   c->nonew   = 0;
193f4a1e0b0SHong Zhang   (*C)->ops->ptapnumeric = MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy;
19453565b12SHong Zhang 
1957d69fa63SHong Zhang   /* set MatInfo */
1967d69fa63SHong Zhang   afill = (PetscReal)ci[pn]/(ai[am]+pi[pm] + 1.e-5);
1977d69fa63SHong Zhang   if (afill < 1.0) afill = 1.0;
1987d69fa63SHong Zhang   c->maxnz                     = ci[pn];
1997d69fa63SHong Zhang   c->nz                        = ci[pn];
2007d69fa63SHong Zhang   (*C)->info.mallocs           = nspacedouble;
2017d69fa63SHong Zhang   (*C)->info.fill_ratio_given  = fill;
2027d69fa63SHong Zhang   (*C)->info.fill_ratio_needed = afill;
2037d69fa63SHong Zhang 
20453565b12SHong Zhang   /* Clean up. */
20553565b12SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
20653565b12SHong Zhang #if defined(PETSC_USE_INFO)
20753565b12SHong Zhang   if (ci[pn] != 0) {
20857622a8eSBarry Smith     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
20957622a8eSBarry Smith     ierr = PetscInfo1((*C),"Use MatPtAP(A,P,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr);
21053565b12SHong Zhang   } else {
21153565b12SHong Zhang     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
21253565b12SHong Zhang   }
213f79958ebSHong Zhang #endif
21453565b12SHong Zhang   PetscFunctionReturn(0);
21553565b12SHong Zhang }
21653565b12SHong Zhang 
217ae32dd66SHong Zhang PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat A,Mat P,Mat C)
21853565b12SHong Zhang {
21953565b12SHong Zhang   PetscErrorCode ierr;
22053565b12SHong Zhang   Mat_SeqAIJ     *a = (Mat_SeqAIJ*) A->data;
22153565b12SHong Zhang   Mat_SeqAIJ     *p = (Mat_SeqAIJ*) P->data;
22253565b12SHong Zhang   Mat_SeqAIJ     *c = (Mat_SeqAIJ*) C->data;
22353565b12SHong Zhang   PetscInt       *ai=a->i,*aj=a->j,*apj,*apjdense,*pi=p->i,*pj=p->j,*pJ=p->j,*pjj;
22453565b12SHong Zhang   PetscInt       *ci=c->i,*cj=c->j,*cjj;
22553565b12SHong Zhang   PetscInt       am =A->rmap->N,cn=C->cmap->N,cm=C->rmap->N;
22653565b12SHong Zhang   PetscInt       i,j,k,anzi,pnzi,apnzj,nextap,pnzj,prow,crow;
22753565b12SHong Zhang   MatScalar      *aa=a->a,*apa,*pa=p->a,*pA=p->a,*paj,*ca=c->a,*caj;
22853565b12SHong Zhang 
22953565b12SHong Zhang   PetscFunctionBegin;
230ae32dd66SHong Zhang   /* Allocate temporary array for storage of one row of A*P (cn: non-scalable) */
231*8f8f2f0dSBarry Smith   ierr = PetscCalloc2(cn,&apa,cn,&apjdense);CHKERRQ(ierr);
232857a15f1SBarry Smith   ierr = PetscMalloc1(cn,&apj);CHKERRQ(ierr);
23353565b12SHong Zhang 
23453565b12SHong Zhang   /* Clear old values in C */
235580bdb30SBarry Smith   ierr = PetscArrayzero(ca,ci[cm]);CHKERRQ(ierr);
23653565b12SHong Zhang 
23753565b12SHong Zhang   for (i=0; i<am; i++) {
23853565b12SHong Zhang     /* Form sparse row of A*P */
23953565b12SHong Zhang     anzi  = ai[i+1] - ai[i];
24053565b12SHong Zhang     apnzj = 0;
24153565b12SHong Zhang     for (j=0; j<anzi; j++) {
24253565b12SHong Zhang       prow = *aj++;
24353565b12SHong Zhang       pnzj = pi[prow+1] - pi[prow];
24453565b12SHong Zhang       pjj  = pj + pi[prow];
24553565b12SHong Zhang       paj  = pa + pi[prow];
24653565b12SHong Zhang       for (k=0; k<pnzj; k++) {
24753565b12SHong Zhang         if (!apjdense[pjj[k]]) {
24853565b12SHong Zhang           apjdense[pjj[k]] = -1;
24953565b12SHong Zhang           apj[apnzj++]     = pjj[k];
25053565b12SHong Zhang         }
25153565b12SHong Zhang         apa[pjj[k]] += (*aa)*paj[k];
25253565b12SHong Zhang       }
25353565b12SHong Zhang       ierr = PetscLogFlops(2.0*pnzj);CHKERRQ(ierr);
25453565b12SHong Zhang       aa++;
25553565b12SHong Zhang     }
25653565b12SHong Zhang 
25753565b12SHong Zhang     /* Sort the j index array for quick sparse axpy. */
25853565b12SHong Zhang     /* Note: a array does not need sorting as it is in dense storage locations. */
25953565b12SHong Zhang     ierr = PetscSortInt(apnzj,apj);CHKERRQ(ierr);
26053565b12SHong Zhang 
26153565b12SHong Zhang     /* Compute P^T*A*P using outer product (P^T)[:,j]*(A*P)[j,:]. */
26253565b12SHong Zhang     pnzi = pi[i+1] - pi[i];
26353565b12SHong Zhang     for (j=0; j<pnzi; j++) {
26453565b12SHong Zhang       nextap = 0;
26553565b12SHong Zhang       crow   = *pJ++;
26653565b12SHong Zhang       cjj    = cj + ci[crow];
26753565b12SHong Zhang       caj    = ca + ci[crow];
26853565b12SHong Zhang       /* Perform sparse axpy operation.  Note cjj includes apj. */
26953565b12SHong Zhang       for (k=0; nextap<apnzj; k++) {
27053565b12SHong Zhang #if defined(PETSC_USE_DEBUG)
271f23aa3ddSBarry Smith         if (k >= ci[crow+1] - ci[crow]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"k too large k %d, crow %d",k,crow);
27253565b12SHong Zhang #endif
27353565b12SHong Zhang         if (cjj[k]==apj[nextap]) {
27453565b12SHong Zhang           caj[k] += (*pA)*apa[apj[nextap++]];
27553565b12SHong Zhang         }
27653565b12SHong Zhang       }
27753565b12SHong Zhang       ierr = PetscLogFlops(2.0*apnzj);CHKERRQ(ierr);
27853565b12SHong Zhang       pA++;
27953565b12SHong Zhang     }
28053565b12SHong Zhang 
28153565b12SHong Zhang     /* Zero the current row info for A*P */
28253565b12SHong Zhang     for (j=0; j<apnzj; j++) {
28353565b12SHong Zhang       apa[apj[j]]      = 0.;
28453565b12SHong Zhang       apjdense[apj[j]] = 0;
28553565b12SHong Zhang     }
28653565b12SHong Zhang   }
28753565b12SHong Zhang 
28853565b12SHong Zhang   /* Assemble the final matrix and clean up */
28953565b12SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
29053565b12SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
291ae32dd66SHong Zhang 
292857a15f1SBarry Smith   ierr = PetscFree2(apa,apjdense);CHKERRQ(ierr);
293857a15f1SBarry Smith   ierr = PetscFree(apj);CHKERRQ(ierr);
29453565b12SHong Zhang   PetscFunctionReturn(0);
29553565b12SHong Zhang }
29653565b12SHong Zhang 
297dfbe8321SBarry Smith PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat A,Mat P,Mat C)
298dfbe8321SBarry Smith {
299dfbe8321SBarry Smith   PetscErrorCode      ierr;
300d20bfe6fSHong Zhang   Mat_SeqAIJ          *c = (Mat_SeqAIJ*)C->data;
3018fa4b5a6SHong Zhang   Mat_MatTransMatMult *atb = c->atb;
3028fa4b5a6SHong Zhang   Mat                 Pt = atb->At;
303eb9c0419SKris Buschelman 
304eb9c0419SKris Buschelman   PetscFunctionBegin;
3058fa4b5a6SHong Zhang   ierr = MatTranspose_SeqAIJ(P,MAT_REUSE_MATRIX,&Pt);CHKERRQ(ierr);
3068fa4b5a6SHong Zhang   ierr = (C->ops->matmatmultnumeric)(Pt,A,P,C);CHKERRQ(ierr);
30753565b12SHong Zhang   PetscFunctionReturn(0);
30853565b12SHong Zhang }
309