xref: /petsc/src/mat/impls/aij/seq/matptap.c (revision 0abfe76ec6b4e1cef922a23993ccad7b6ad6771d)
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      */
394a1b09b7SHong Zhang     ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
4068eacb73SHong Zhang     PetscOptionsObject->alreadyprinted = PETSC_FALSE; /* a hack to ensure the option shows in '-help' */
416f231fbdSstefano_zampini     ierr = PetscOptionsEList("-matptap_via","Algorithmic approach","MatPtAP",algTypes,nalg,algTypes[0],&alg,NULL);CHKERRQ(ierr);
424a1b09b7SHong Zhang     ierr = PetscOptionsEnd();CHKERRQ(ierr);
434a1b09b7SHong Zhang     switch (alg) {
444a1b09b7SHong Zhang     case 1:
458fa4b5a6SHong Zhang       ierr = PetscNew(&atb);CHKERRQ(ierr);
468fa4b5a6SHong Zhang       ierr = MatTranspose_SeqAIJ(P,MAT_INITIAL_MATRIX,&Pt);CHKERRQ(ierr);
478fa4b5a6SHong Zhang       ierr = MatMatMatMult(Pt,A,P,MAT_INITIAL_MATRIX,fill,C);CHKERRQ(ierr);
488fa4b5a6SHong Zhang 
498fa4b5a6SHong Zhang       c                      = (Mat_SeqAIJ*)(*C)->data;
508fa4b5a6SHong Zhang       c->atb                 = atb;
518fa4b5a6SHong Zhang       atb->At                = Pt;
528fa4b5a6SHong Zhang       atb->destroy           = (*C)->ops->destroy;
538fa4b5a6SHong Zhang       (*C)->ops->destroy     = MatDestroy_SeqAIJ_MatTransMatMult;
548fa4b5a6SHong Zhang       (*C)->ops->ptapnumeric = MatPtAPNumeric_SeqAIJ_SeqAIJ;
558fa4b5a6SHong Zhang       PetscFunctionReturn(0);
564a1b09b7SHong Zhang       break;
576f231fbdSstefano_zampini #if defined(PETSC_HAVE_HYPRE)
586f231fbdSstefano_zampini     case 2:
59*0abfe76eSFande Kong       ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
606f231fbdSstefano_zampini       ierr = MatPtAPSymbolic_AIJ_AIJ_wHYPRE(A,P,fill,C);CHKERRQ(ierr);
61*0abfe76eSFande Kong       ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
626f231fbdSstefano_zampini       break;
636f231fbdSstefano_zampini #endif
644a1b09b7SHong Zhang     default:
65*0abfe76eSFande Kong       ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
664a1b09b7SHong Zhang       ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(A,P,fill,C);CHKERRQ(ierr);
67*0abfe76eSFande Kong       ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
684a1b09b7SHong Zhang       break;
694a1b09b7SHong Zhang     }
707ba1a0bfSKris Buschelman   }
71*0abfe76eSFande Kong   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
72ae32dd66SHong Zhang   ierr = (*(*C)->ops->ptapnumeric)(A,P,*C);CHKERRQ(ierr);
73*0abfe76eSFande Kong   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
747ba1a0bfSKris Buschelman   PetscFunctionReturn(0);
757ba1a0bfSKris Buschelman }
767ba1a0bfSKris Buschelman 
7753565b12SHong Zhang PetscErrorCode MatDestroy_SeqAIJ_PtAP(Mat A)
788d0a38e4SHong Zhang {
798d0a38e4SHong Zhang   PetscErrorCode ierr;
8053565b12SHong Zhang   Mat_SeqAIJ     *a    = (Mat_SeqAIJ*)A->data;
8153565b12SHong Zhang   Mat_PtAP       *ptap = a->ptap;
828d0a38e4SHong Zhang 
838d0a38e4SHong Zhang   PetscFunctionBegin;
848d0a38e4SHong Zhang   ierr = PetscFree(ptap->apa);CHKERRQ(ierr);
85f79958ebSHong Zhang   ierr = PetscFree(ptap->api);CHKERRQ(ierr);
86f79958ebSHong Zhang   ierr = PetscFree(ptap->apj);CHKERRQ(ierr);
8753565b12SHong Zhang   ierr = (ptap->destroy)(A);CHKERRQ(ierr);
888d0a38e4SHong Zhang   ierr = PetscFree(ptap);CHKERRQ(ierr);
898d0a38e4SHong Zhang   PetscFunctionReturn(0);
908d0a38e4SHong Zhang }
918d0a38e4SHong Zhang 
92ae32dd66SHong Zhang PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(Mat A,Mat P,PetscReal fill,Mat *C)
938d0a38e4SHong Zhang {
948d0a38e4SHong Zhang   PetscErrorCode     ierr;
950298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
96d20bfe6fSHong Zhang   Mat_SeqAIJ         *a        = (Mat_SeqAIJ*)A->data,*p = (Mat_SeqAIJ*)P->data,*c;
9753565b12SHong Zhang   PetscInt           *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj;
9853565b12SHong Zhang   PetscInt           *ci,*cj,*ptadenserow,*ptasparserow,*ptaj,nspacedouble=0;
997d69fa63SHong Zhang   PetscInt           an=A->cmap->N,am=A->rmap->N,pn=P->cmap->N,pm=P->rmap->N;
10053565b12SHong Zhang   PetscInt           i,j,k,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi,nlnk,*lnk;
101d20bfe6fSHong Zhang   MatScalar          *ca;
10253565b12SHong Zhang   PetscBT            lnkbt;
1037d69fa63SHong Zhang   PetscReal          afill;
104eb9c0419SKris Buschelman 
105eb9c0419SKris Buschelman   PetscFunctionBegin;
10653565b12SHong Zhang   /* Get ij structure of P^T */
107eb9c0419SKris Buschelman   ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
108d20bfe6fSHong Zhang   ptJ  = ptj;
109eb9c0419SKris Buschelman 
110d20bfe6fSHong Zhang   /* Allocate ci array, arrays for fill computation and */
111d20bfe6fSHong Zhang   /* free space for accumulating nonzero column info */
112854ce69bSBarry Smith   ierr  = PetscMalloc1(pn+1,&ci);CHKERRQ(ierr);
113d20bfe6fSHong Zhang   ci[0] = 0;
114eb9c0419SKris Buschelman 
1151795a4d1SJed Brown   ierr         = PetscCalloc1(2*an+1,&ptadenserow);CHKERRQ(ierr);
116d20bfe6fSHong Zhang   ptasparserow = ptadenserow  + an;
11755a3bba9SHong Zhang 
11855a3bba9SHong Zhang   /* create and initialize a linked list */
11955a3bba9SHong Zhang   nlnk = pn+1;
12055a3bba9SHong Zhang   ierr = PetscLLCreate(pn,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr);
121eb9c0419SKris Buschelman 
1227d69fa63SHong Zhang   /* Set initial free space to be fill*(nnz(A)+ nnz(P)) */
123f91af8c7SBarry Smith   ierr          = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],pi[pm])),&free_space);CHKERRQ(ierr);
124d20bfe6fSHong Zhang   current_space = free_space;
125d20bfe6fSHong Zhang 
126d20bfe6fSHong Zhang   /* Determine symbolic info for each row of C: */
127d20bfe6fSHong Zhang   for (i=0; i<pn; i++) {
128d20bfe6fSHong Zhang     ptnzi  = pti[i+1] - pti[i];
129d20bfe6fSHong Zhang     ptanzi = 0;
130d20bfe6fSHong Zhang     /* Determine symbolic row of PtA: */
131d20bfe6fSHong Zhang     for (j=0; j<ptnzi; j++) {
132d20bfe6fSHong Zhang       arow = *ptJ++;
133d20bfe6fSHong Zhang       anzj = ai[arow+1] - ai[arow];
134d20bfe6fSHong Zhang       ajj  = aj + ai[arow];
135d20bfe6fSHong Zhang       for (k=0; k<anzj; k++) {
136d20bfe6fSHong Zhang         if (!ptadenserow[ajj[k]]) {
137d20bfe6fSHong Zhang           ptadenserow[ajj[k]]    = -1;
138d20bfe6fSHong Zhang           ptasparserow[ptanzi++] = ajj[k];
139d20bfe6fSHong Zhang         }
140d20bfe6fSHong Zhang       }
141d20bfe6fSHong Zhang     }
142d20bfe6fSHong Zhang     /* Using symbolic info for row of PtA, determine symbolic info for row of C: */
143d20bfe6fSHong Zhang     ptaj = ptasparserow;
144d20bfe6fSHong Zhang     cnzi = 0;
145d20bfe6fSHong Zhang     for (j=0; j<ptanzi; j++) {
146d20bfe6fSHong Zhang       prow = *ptaj++;
147d20bfe6fSHong Zhang       pnzj = pi[prow+1] - pi[prow];
148d20bfe6fSHong Zhang       pjj  = pj + pi[prow];
14955a3bba9SHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
150dadf0e6bSHong Zhang       ierr  = PetscLLAddSorted(pnzj,pjj,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr);
15155a3bba9SHong Zhang       cnzi += nlnk;
152d20bfe6fSHong Zhang     }
153d20bfe6fSHong Zhang 
154d20bfe6fSHong Zhang     /* If free space is not available, make more free space */
155d20bfe6fSHong Zhang     /* Double the amount of total space in the list */
156d20bfe6fSHong Zhang     if (current_space->local_remaining<cnzi) {
157f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(cnzi,current_space->total_array_size),&current_space);CHKERRQ(ierr);
158f2b054eeSHong Zhang       nspacedouble++;
159d20bfe6fSHong Zhang     }
160d20bfe6fSHong Zhang 
161d20bfe6fSHong Zhang     /* Copy data into free space, and zero out denserows */
16255a3bba9SHong Zhang     ierr = PetscLLClean(pn,pn,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
1632205254eSKarl Rupp 
164d20bfe6fSHong Zhang     current_space->array           += cnzi;
165d20bfe6fSHong Zhang     current_space->local_used      += cnzi;
166d20bfe6fSHong Zhang     current_space->local_remaining -= cnzi;
167d20bfe6fSHong Zhang 
1682205254eSKarl Rupp     for (j=0; j<ptanzi; j++) ptadenserow[ptasparserow[j]] = 0;
1692205254eSKarl Rupp 
170d20bfe6fSHong Zhang     /* Aside: Perhaps we should save the pta info for the numerical factorization. */
171d20bfe6fSHong Zhang     /*        For now, we will recompute what is needed. */
172d20bfe6fSHong Zhang     ci[i+1] = ci[i] + cnzi;
173d20bfe6fSHong Zhang   }
174d20bfe6fSHong Zhang   /* nnz is now stored in ci[ptm], column indices are in the list of free space */
175d20bfe6fSHong Zhang   /* Allocate space for cj, initialize cj, and */
176d20bfe6fSHong Zhang   /* destroy list of free space and other temporary array(s) */
177854ce69bSBarry Smith   ierr = PetscMalloc1(ci[pn]+1,&cj);CHKERRQ(ierr);
178a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
179d20bfe6fSHong Zhang   ierr = PetscFree(ptadenserow);CHKERRQ(ierr);
18055a3bba9SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
181d20bfe6fSHong Zhang 
182854ce69bSBarry Smith   ierr = PetscCalloc1(ci[pn]+1,&ca);CHKERRQ(ierr);
18353565b12SHong Zhang 
18453565b12SHong Zhang   /* put together the new matrix */
185ce94432eSBarry Smith   ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),pn,pn,ci,cj,ca,C);CHKERRQ(ierr);
18633d57670SJed Brown   ierr = MatSetBlockSizes(*C,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr);
187ae32dd66SHong Zhang 
18853565b12SHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
18953565b12SHong Zhang   /* Since these are PETSc arrays, change flags to free them as necessary. */
19053565b12SHong Zhang   c          = (Mat_SeqAIJ*)((*C)->data);
19153565b12SHong Zhang   c->free_a  = PETSC_TRUE;
19253565b12SHong Zhang   c->free_ij = PETSC_TRUE;
19353565b12SHong Zhang   c->nonew   = 0;
194f4a1e0b0SHong Zhang   (*C)->ops->ptapnumeric = MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy;
19553565b12SHong Zhang 
1967d69fa63SHong Zhang   /* set MatInfo */
1977d69fa63SHong Zhang   afill = (PetscReal)ci[pn]/(ai[am]+pi[pm] + 1.e-5);
1987d69fa63SHong Zhang   if (afill < 1.0) afill = 1.0;
1997d69fa63SHong Zhang   c->maxnz                     = ci[pn];
2007d69fa63SHong Zhang   c->nz                        = ci[pn];
2017d69fa63SHong Zhang   (*C)->info.mallocs           = nspacedouble;
2027d69fa63SHong Zhang   (*C)->info.fill_ratio_given  = fill;
2037d69fa63SHong Zhang   (*C)->info.fill_ratio_needed = afill;
2047d69fa63SHong Zhang 
20553565b12SHong Zhang   /* Clean up. */
20653565b12SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
20753565b12SHong Zhang #if defined(PETSC_USE_INFO)
20853565b12SHong Zhang   if (ci[pn] != 0) {
20957622a8eSBarry Smith     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
21057622a8eSBarry Smith     ierr = PetscInfo1((*C),"Use MatPtAP(A,P,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr);
21153565b12SHong Zhang   } else {
21253565b12SHong Zhang     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
21353565b12SHong Zhang   }
214f79958ebSHong Zhang #endif
21553565b12SHong Zhang   PetscFunctionReturn(0);
21653565b12SHong Zhang }
21753565b12SHong Zhang 
218ae32dd66SHong Zhang PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat A,Mat P,Mat C)
21953565b12SHong Zhang {
22053565b12SHong Zhang   PetscErrorCode ierr;
22153565b12SHong Zhang   Mat_SeqAIJ     *a = (Mat_SeqAIJ*) A->data;
22253565b12SHong Zhang   Mat_SeqAIJ     *p = (Mat_SeqAIJ*) P->data;
22353565b12SHong Zhang   Mat_SeqAIJ     *c = (Mat_SeqAIJ*) C->data;
22453565b12SHong Zhang   PetscInt       *ai=a->i,*aj=a->j,*apj,*apjdense,*pi=p->i,*pj=p->j,*pJ=p->j,*pjj;
22553565b12SHong Zhang   PetscInt       *ci=c->i,*cj=c->j,*cjj;
22653565b12SHong Zhang   PetscInt       am =A->rmap->N,cn=C->cmap->N,cm=C->rmap->N;
22753565b12SHong Zhang   PetscInt       i,j,k,anzi,pnzi,apnzj,nextap,pnzj,prow,crow;
22853565b12SHong Zhang   MatScalar      *aa=a->a,*apa,*pa=p->a,*pA=p->a,*paj,*ca=c->a,*caj;
22953565b12SHong Zhang 
23053565b12SHong Zhang   PetscFunctionBegin;
231ae32dd66SHong Zhang   /* Allocate temporary array for storage of one row of A*P (cn: non-scalable) */
2324bcca364SHong Zhang   ierr = PetscMalloc3(cn,&apa,cn,&apjdense,cn,&apj);CHKERRQ(ierr);
233eb9baa12SBarry Smith   ierr = PetscMemzero(apa,cn*sizeof(MatScalar));CHKERRQ(ierr);
234eb9baa12SBarry Smith   ierr = PetscMemzero(apjdense,cn*sizeof(PetscInt));CHKERRQ(ierr);
23553565b12SHong Zhang 
23653565b12SHong Zhang   /* Clear old values in C */
23753565b12SHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
23853565b12SHong Zhang 
23953565b12SHong Zhang   for (i=0; i<am; i++) {
24053565b12SHong Zhang     /* Form sparse row of A*P */
24153565b12SHong Zhang     anzi  = ai[i+1] - ai[i];
24253565b12SHong Zhang     apnzj = 0;
24353565b12SHong Zhang     for (j=0; j<anzi; j++) {
24453565b12SHong Zhang       prow = *aj++;
24553565b12SHong Zhang       pnzj = pi[prow+1] - pi[prow];
24653565b12SHong Zhang       pjj  = pj + pi[prow];
24753565b12SHong Zhang       paj  = pa + pi[prow];
24853565b12SHong Zhang       for (k=0; k<pnzj; k++) {
24953565b12SHong Zhang         if (!apjdense[pjj[k]]) {
25053565b12SHong Zhang           apjdense[pjj[k]] = -1;
25153565b12SHong Zhang           apj[apnzj++]     = pjj[k];
25253565b12SHong Zhang         }
25353565b12SHong Zhang         apa[pjj[k]] += (*aa)*paj[k];
25453565b12SHong Zhang       }
25553565b12SHong Zhang       ierr = PetscLogFlops(2.0*pnzj);CHKERRQ(ierr);
25653565b12SHong Zhang       aa++;
25753565b12SHong Zhang     }
25853565b12SHong Zhang 
25953565b12SHong Zhang     /* Sort the j index array for quick sparse axpy. */
26053565b12SHong Zhang     /* Note: a array does not need sorting as it is in dense storage locations. */
26153565b12SHong Zhang     ierr = PetscSortInt(apnzj,apj);CHKERRQ(ierr);
26253565b12SHong Zhang 
26353565b12SHong Zhang     /* Compute P^T*A*P using outer product (P^T)[:,j]*(A*P)[j,:]. */
26453565b12SHong Zhang     pnzi = pi[i+1] - pi[i];
26553565b12SHong Zhang     for (j=0; j<pnzi; j++) {
26653565b12SHong Zhang       nextap = 0;
26753565b12SHong Zhang       crow   = *pJ++;
26853565b12SHong Zhang       cjj    = cj + ci[crow];
26953565b12SHong Zhang       caj    = ca + ci[crow];
27053565b12SHong Zhang       /* Perform sparse axpy operation.  Note cjj includes apj. */
27153565b12SHong Zhang       for (k=0; nextap<apnzj; k++) {
27253565b12SHong Zhang #if defined(PETSC_USE_DEBUG)
273f23aa3ddSBarry Smith         if (k >= ci[crow+1] - ci[crow]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"k too large k %d, crow %d",k,crow);
27453565b12SHong Zhang #endif
27553565b12SHong Zhang         if (cjj[k]==apj[nextap]) {
27653565b12SHong Zhang           caj[k] += (*pA)*apa[apj[nextap++]];
27753565b12SHong Zhang         }
27853565b12SHong Zhang       }
27953565b12SHong Zhang       ierr = PetscLogFlops(2.0*apnzj);CHKERRQ(ierr);
28053565b12SHong Zhang       pA++;
28153565b12SHong Zhang     }
28253565b12SHong Zhang 
28353565b12SHong Zhang     /* Zero the current row info for A*P */
28453565b12SHong Zhang     for (j=0; j<apnzj; j++) {
28553565b12SHong Zhang       apa[apj[j]]      = 0.;
28653565b12SHong Zhang       apjdense[apj[j]] = 0;
28753565b12SHong Zhang     }
28853565b12SHong Zhang   }
28953565b12SHong Zhang 
29053565b12SHong Zhang   /* Assemble the final matrix and clean up */
29153565b12SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
29253565b12SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
293ae32dd66SHong Zhang 
294eb9baa12SBarry Smith   ierr = PetscFree3(apa,apjdense,apj);CHKERRQ(ierr);
29553565b12SHong Zhang   PetscFunctionReturn(0);
29653565b12SHong Zhang }
29753565b12SHong Zhang 
298dfbe8321SBarry Smith PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat A,Mat P,Mat C)
299dfbe8321SBarry Smith {
300dfbe8321SBarry Smith   PetscErrorCode      ierr;
301d20bfe6fSHong Zhang   Mat_SeqAIJ          *c = (Mat_SeqAIJ*)C->data;
3028fa4b5a6SHong Zhang   Mat_MatTransMatMult *atb = c->atb;
3038fa4b5a6SHong Zhang   Mat                 Pt = atb->At;
304eb9c0419SKris Buschelman 
305eb9c0419SKris Buschelman   PetscFunctionBegin;
3068fa4b5a6SHong Zhang   ierr = MatTranspose_SeqAIJ(P,MAT_REUSE_MATRIX,&Pt);CHKERRQ(ierr);
3078fa4b5a6SHong Zhang   ierr = (C->ops->matmatmultnumeric)(Pt,A,P,C);CHKERRQ(ierr);
30853565b12SHong Zhang   PetscFunctionReturn(0);
30953565b12SHong Zhang }
310