xref: /petsc/src/mat/impls/aij/seq/matptap.c (revision cf9c20a2d58da010f7c4712defbcdf61cc8f72b5)
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)
134222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatPtAPSymbolic_AIJ_AIJ_wHYPRE(Mat,Mat,PetscReal,Mat);
146f231fbdSstefano_zampini #endif
156f231fbdSstefano_zampini 
164222ddf1SHong Zhang PetscErrorCode MatProductSymbolic_PtAP_SeqAIJ_SeqAIJ(Mat C)
177ba1a0bfSKris Buschelman {
187ba1a0bfSKris Buschelman   PetscErrorCode      ierr;
194222ddf1SHong Zhang   Mat_Product         *product = C->product;
204222ddf1SHong Zhang   Mat                 A=product->A,P=product->B;
214222ddf1SHong Zhang   MatProductAlgorithm alg=product->alg;
224222ddf1SHong Zhang   PetscReal           fill=product->fill;
234222ddf1SHong Zhang   PetscBool           flg;
248fa4b5a6SHong Zhang   Mat                 Pt;
258fa4b5a6SHong Zhang   Mat_MatTransMatMult *atb;
268fa4b5a6SHong Zhang   Mat_SeqAIJ          *c;
277ba1a0bfSKris Buschelman 
287ba1a0bfSKris Buschelman   PetscFunctionBegin;
294222ddf1SHong Zhang   /* "scalable" */
304222ddf1SHong Zhang   ierr = PetscStrcmp(alg,"scalable",&flg);CHKERRQ(ierr);
314222ddf1SHong Zhang   if (flg) {
324222ddf1SHong Zhang     ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(A,P,fill,C);CHKERRQ(ierr);
334222ddf1SHong Zhang     C->ops->productnumeric = MatProductNumeric_PtAP;
344222ddf1SHong Zhang     PetscFunctionReturn(0);
354222ddf1SHong Zhang   }
364222ddf1SHong Zhang 
374222ddf1SHong Zhang   /* "rap" */
384222ddf1SHong Zhang   ierr = PetscStrcmp(alg,"rap",&flg);CHKERRQ(ierr);
394222ddf1SHong Zhang   if (flg) { /* Set default algorithm */
408fa4b5a6SHong Zhang     ierr = PetscNew(&atb);CHKERRQ(ierr);
418fa4b5a6SHong Zhang     ierr = MatTranspose_SeqAIJ(P,MAT_INITIAL_MATRIX,&Pt);CHKERRQ(ierr);
424222ddf1SHong Zhang     ierr = MatMatMatMultSymbolic_SeqAIJ_SeqAIJ_SeqAIJ(Pt,A,P,fill,C);CHKERRQ(ierr);
438fa4b5a6SHong Zhang 
444222ddf1SHong Zhang     c                      = (Mat_SeqAIJ*)C->data;
458fa4b5a6SHong Zhang     c->atb                 = atb;
468fa4b5a6SHong Zhang     atb->At                = Pt;
474222ddf1SHong Zhang     atb->destroy           = C->ops->destroy;
484222ddf1SHong Zhang     C->ops->destroy        = MatDestroy_SeqAIJ_MatTransMatMult;
494222ddf1SHong Zhang     C->ops->ptapnumeric    = MatPtAPNumeric_SeqAIJ_SeqAIJ;
504222ddf1SHong Zhang     C->ops->productnumeric = MatProductNumeric_PtAP;
518fa4b5a6SHong Zhang     PetscFunctionReturn(0);
524222ddf1SHong Zhang   }
534222ddf1SHong Zhang 
544222ddf1SHong Zhang   /* hypre */
556f231fbdSstefano_zampini #if defined(PETSC_HAVE_HYPRE)
564222ddf1SHong Zhang   ierr = PetscStrcmp(alg,"hypre",&flg);CHKERRQ(ierr);
574222ddf1SHong Zhang   if (flg) {
586f231fbdSstefano_zampini     ierr = MatPtAPSymbolic_AIJ_AIJ_wHYPRE(A,P,fill,C);CHKERRQ(ierr);
594222ddf1SHong Zhang     PetscFunctionReturn(0);
604222ddf1SHong Zhang   }
616f231fbdSstefano_zampini #endif
624222ddf1SHong Zhang 
634222ddf1SHong Zhang   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatProductType is not supported");
647ba1a0bfSKris Buschelman   PetscFunctionReturn(0);
657ba1a0bfSKris Buschelman }
667ba1a0bfSKris Buschelman 
6753565b12SHong Zhang PetscErrorCode MatDestroy_SeqAIJ_PtAP(Mat A)
688d0a38e4SHong Zhang {
698d0a38e4SHong Zhang   PetscErrorCode ierr;
7053565b12SHong Zhang   Mat_SeqAIJ     *a  = (Mat_SeqAIJ*)A->data;
713cdca5ebSHong Zhang   Mat_AP         *ap = a->ap;
728d0a38e4SHong Zhang 
738d0a38e4SHong Zhang   PetscFunctionBegin;
743cdca5ebSHong Zhang   ierr = PetscFree(ap->apa);CHKERRQ(ierr);
753cdca5ebSHong Zhang   ierr = PetscFree(ap->api);CHKERRQ(ierr);
763cdca5ebSHong Zhang   ierr = PetscFree(ap->apj);CHKERRQ(ierr);
773cdca5ebSHong Zhang   ierr = (ap->destroy)(A);CHKERRQ(ierr);
783cdca5ebSHong Zhang   ierr = PetscFree(ap);CHKERRQ(ierr);
798d0a38e4SHong Zhang   PetscFunctionReturn(0);
808d0a38e4SHong Zhang }
818d0a38e4SHong Zhang 
824222ddf1SHong Zhang PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(Mat A,Mat P,PetscReal fill,Mat C)
838d0a38e4SHong Zhang {
848d0a38e4SHong Zhang   PetscErrorCode     ierr;
850298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
86d20bfe6fSHong Zhang   Mat_SeqAIJ         *a        = (Mat_SeqAIJ*)A->data,*p = (Mat_SeqAIJ*)P->data,*c;
8753565b12SHong Zhang   PetscInt           *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj;
8853565b12SHong Zhang   PetscInt           *ci,*cj,*ptadenserow,*ptasparserow,*ptaj,nspacedouble=0;
897d69fa63SHong Zhang   PetscInt           an=A->cmap->N,am=A->rmap->N,pn=P->cmap->N,pm=P->rmap->N;
9053565b12SHong Zhang   PetscInt           i,j,k,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi,nlnk,*lnk;
91d20bfe6fSHong Zhang   MatScalar          *ca;
9253565b12SHong Zhang   PetscBT            lnkbt;
937d69fa63SHong Zhang   PetscReal          afill;
94eb9c0419SKris Buschelman 
95eb9c0419SKris Buschelman   PetscFunctionBegin;
9653565b12SHong Zhang   /* Get ij structure of P^T */
97eb9c0419SKris Buschelman   ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
98d20bfe6fSHong Zhang   ptJ  = ptj;
99eb9c0419SKris Buschelman 
100d20bfe6fSHong Zhang   /* Allocate ci array, arrays for fill computation and */
101d20bfe6fSHong Zhang   /* free space for accumulating nonzero column info */
102854ce69bSBarry Smith   ierr  = PetscMalloc1(pn+1,&ci);CHKERRQ(ierr);
103d20bfe6fSHong Zhang   ci[0] = 0;
104eb9c0419SKris Buschelman 
1051795a4d1SJed Brown   ierr         = PetscCalloc1(2*an+1,&ptadenserow);CHKERRQ(ierr);
106d20bfe6fSHong Zhang   ptasparserow = ptadenserow  + an;
10755a3bba9SHong Zhang 
10855a3bba9SHong Zhang   /* create and initialize a linked list */
10955a3bba9SHong Zhang   nlnk = pn+1;
11055a3bba9SHong Zhang   ierr = PetscLLCreate(pn,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr);
111eb9c0419SKris Buschelman 
1127d69fa63SHong Zhang   /* Set initial free space to be fill*(nnz(A)+ nnz(P)) */
113f91af8c7SBarry Smith   ierr          = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],pi[pm])),&free_space);CHKERRQ(ierr);
114d20bfe6fSHong Zhang   current_space = free_space;
115d20bfe6fSHong Zhang 
116d20bfe6fSHong Zhang   /* Determine symbolic info for each row of C: */
117d20bfe6fSHong Zhang   for (i=0; i<pn; i++) {
118d20bfe6fSHong Zhang     ptnzi  = pti[i+1] - pti[i];
119d20bfe6fSHong Zhang     ptanzi = 0;
120d20bfe6fSHong Zhang     /* Determine symbolic row of PtA: */
121d20bfe6fSHong Zhang     for (j=0; j<ptnzi; j++) {
122d20bfe6fSHong Zhang       arow = *ptJ++;
123d20bfe6fSHong Zhang       anzj = ai[arow+1] - ai[arow];
124d20bfe6fSHong Zhang       ajj  = aj + ai[arow];
125d20bfe6fSHong Zhang       for (k=0; k<anzj; k++) {
126d20bfe6fSHong Zhang         if (!ptadenserow[ajj[k]]) {
127d20bfe6fSHong Zhang           ptadenserow[ajj[k]]    = -1;
128d20bfe6fSHong Zhang           ptasparserow[ptanzi++] = ajj[k];
129d20bfe6fSHong Zhang         }
130d20bfe6fSHong Zhang       }
131d20bfe6fSHong Zhang     }
132d20bfe6fSHong Zhang     /* Using symbolic info for row of PtA, determine symbolic info for row of C: */
133d20bfe6fSHong Zhang     ptaj = ptasparserow;
134d20bfe6fSHong Zhang     cnzi = 0;
135d20bfe6fSHong Zhang     for (j=0; j<ptanzi; j++) {
136d20bfe6fSHong Zhang       prow = *ptaj++;
137d20bfe6fSHong Zhang       pnzj = pi[prow+1] - pi[prow];
138d20bfe6fSHong Zhang       pjj  = pj + pi[prow];
13955a3bba9SHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
140dadf0e6bSHong Zhang       ierr  = PetscLLAddSorted(pnzj,pjj,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr);
14155a3bba9SHong Zhang       cnzi += nlnk;
142d20bfe6fSHong Zhang     }
143d20bfe6fSHong Zhang 
144d20bfe6fSHong Zhang     /* If free space is not available, make more free space */
145d20bfe6fSHong Zhang     /* Double the amount of total space in the list */
146d20bfe6fSHong Zhang     if (current_space->local_remaining<cnzi) {
147f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(cnzi,current_space->total_array_size),&current_space);CHKERRQ(ierr);
148f2b054eeSHong Zhang       nspacedouble++;
149d20bfe6fSHong Zhang     }
150d20bfe6fSHong Zhang 
151d20bfe6fSHong Zhang     /* Copy data into free space, and zero out denserows */
15255a3bba9SHong Zhang     ierr = PetscLLClean(pn,pn,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
1532205254eSKarl Rupp 
154d20bfe6fSHong Zhang     current_space->array           += cnzi;
155d20bfe6fSHong Zhang     current_space->local_used      += cnzi;
156d20bfe6fSHong Zhang     current_space->local_remaining -= cnzi;
157d20bfe6fSHong Zhang 
1582205254eSKarl Rupp     for (j=0; j<ptanzi; j++) ptadenserow[ptasparserow[j]] = 0;
1592205254eSKarl Rupp 
160d20bfe6fSHong Zhang     /* Aside: Perhaps we should save the pta info for the numerical factorization. */
161d20bfe6fSHong Zhang     /*        For now, we will recompute what is needed. */
162d20bfe6fSHong Zhang     ci[i+1] = ci[i] + cnzi;
163d20bfe6fSHong Zhang   }
164d20bfe6fSHong Zhang   /* nnz is now stored in ci[ptm], column indices are in the list of free space */
165d20bfe6fSHong Zhang   /* Allocate space for cj, initialize cj, and */
166d20bfe6fSHong Zhang   /* destroy list of free space and other temporary array(s) */
167854ce69bSBarry Smith   ierr = PetscMalloc1(ci[pn]+1,&cj);CHKERRQ(ierr);
168a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
169d20bfe6fSHong Zhang   ierr = PetscFree(ptadenserow);CHKERRQ(ierr);
17055a3bba9SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
171d20bfe6fSHong Zhang 
172854ce69bSBarry Smith   ierr = PetscCalloc1(ci[pn]+1,&ca);CHKERRQ(ierr);
17353565b12SHong Zhang 
17453565b12SHong Zhang   /* put together the new matrix */
1754222ddf1SHong Zhang   ierr = MatSetSeqAIJWithArrays_private(PetscObjectComm((PetscObject)A),pn,pn,ci,cj,ca,C);CHKERRQ(ierr);
1764222ddf1SHong Zhang   ierr = MatSetBlockSizes(C,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr);
1774222ddf1SHong Zhang   ierr = MatSetType(C,((PetscObject)A)->type_name);CHKERRQ(ierr);
178ae32dd66SHong Zhang 
17953565b12SHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
18053565b12SHong Zhang   /* Since these are PETSc arrays, change flags to free them as necessary. */
1814222ddf1SHong Zhang   c          = (Mat_SeqAIJ*)((C)->data);
18253565b12SHong Zhang   c->free_a  = PETSC_TRUE;
18353565b12SHong Zhang   c->free_ij = PETSC_TRUE;
18453565b12SHong Zhang   c->nonew   = 0;
1854222ddf1SHong Zhang   C->ops->ptapnumeric    = MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy;
18653565b12SHong Zhang 
1877d69fa63SHong Zhang   /* set MatInfo */
1887d69fa63SHong Zhang   afill = (PetscReal)ci[pn]/(ai[am]+pi[pm] + 1.e-5);
1897d69fa63SHong Zhang   if (afill < 1.0) afill = 1.0;
1907d69fa63SHong Zhang   c->maxnz                     = ci[pn];
1917d69fa63SHong Zhang   c->nz                        = ci[pn];
1924222ddf1SHong Zhang   C->info.mallocs           = nspacedouble;
1934222ddf1SHong Zhang   C->info.fill_ratio_given  = fill;
1944222ddf1SHong Zhang   C->info.fill_ratio_needed = afill;
1957d69fa63SHong Zhang 
19653565b12SHong Zhang   /* Clean up. */
19753565b12SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
19853565b12SHong Zhang #if defined(PETSC_USE_INFO)
19953565b12SHong Zhang   if (ci[pn] != 0) {
2004222ddf1SHong Zhang     ierr = PetscInfo3(C,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
2014222ddf1SHong Zhang     ierr = PetscInfo1(C,"Use MatPtAP(A,P,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr);
20253565b12SHong Zhang   } else {
2034222ddf1SHong Zhang     ierr = PetscInfo(C,"Empty matrix product\n");CHKERRQ(ierr);
20453565b12SHong Zhang   }
205f79958ebSHong Zhang #endif
20653565b12SHong Zhang   PetscFunctionReturn(0);
20753565b12SHong Zhang }
20853565b12SHong Zhang 
209ae32dd66SHong Zhang PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat A,Mat P,Mat C)
21053565b12SHong Zhang {
21153565b12SHong Zhang   PetscErrorCode ierr;
21253565b12SHong Zhang   Mat_SeqAIJ     *a = (Mat_SeqAIJ*) A->data;
21353565b12SHong Zhang   Mat_SeqAIJ     *p = (Mat_SeqAIJ*) P->data;
21453565b12SHong Zhang   Mat_SeqAIJ     *c = (Mat_SeqAIJ*) C->data;
21553565b12SHong Zhang   PetscInt       *ai=a->i,*aj=a->j,*apj,*apjdense,*pi=p->i,*pj=p->j,*pJ=p->j,*pjj;
21653565b12SHong Zhang   PetscInt       *ci=c->i,*cj=c->j,*cjj;
21753565b12SHong Zhang   PetscInt       am =A->rmap->N,cn=C->cmap->N,cm=C->rmap->N;
21853565b12SHong Zhang   PetscInt       i,j,k,anzi,pnzi,apnzj,nextap,pnzj,prow,crow;
21953565b12SHong Zhang   MatScalar      *aa=a->a,*apa,*pa=p->a,*pA=p->a,*paj,*ca=c->a,*caj;
22053565b12SHong Zhang 
22153565b12SHong Zhang   PetscFunctionBegin;
222ae32dd66SHong Zhang   /* Allocate temporary array for storage of one row of A*P (cn: non-scalable) */
2238f8f2f0dSBarry Smith   ierr = PetscCalloc2(cn,&apa,cn,&apjdense);CHKERRQ(ierr);
224857a15f1SBarry Smith   ierr = PetscMalloc1(cn,&apj);CHKERRQ(ierr);
22553565b12SHong Zhang 
22653565b12SHong Zhang   /* Clear old values in C */
227580bdb30SBarry Smith   ierr = PetscArrayzero(ca,ci[cm]);CHKERRQ(ierr);
22853565b12SHong Zhang 
22953565b12SHong Zhang   for (i=0; i<am; i++) {
23053565b12SHong Zhang     /* Form sparse row of A*P */
23153565b12SHong Zhang     anzi  = ai[i+1] - ai[i];
23253565b12SHong Zhang     apnzj = 0;
23353565b12SHong Zhang     for (j=0; j<anzi; j++) {
23453565b12SHong Zhang       prow = *aj++;
23553565b12SHong Zhang       pnzj = pi[prow+1] - pi[prow];
23653565b12SHong Zhang       pjj  = pj + pi[prow];
23753565b12SHong Zhang       paj  = pa + pi[prow];
23853565b12SHong Zhang       for (k=0; k<pnzj; k++) {
23953565b12SHong Zhang         if (!apjdense[pjj[k]]) {
24053565b12SHong Zhang           apjdense[pjj[k]] = -1;
24153565b12SHong Zhang           apj[apnzj++]     = pjj[k];
24253565b12SHong Zhang         }
24353565b12SHong Zhang         apa[pjj[k]] += (*aa)*paj[k];
24453565b12SHong Zhang       }
24553565b12SHong Zhang       ierr = PetscLogFlops(2.0*pnzj);CHKERRQ(ierr);
24653565b12SHong Zhang       aa++;
24753565b12SHong Zhang     }
24853565b12SHong Zhang 
24953565b12SHong Zhang     /* Sort the j index array for quick sparse axpy. */
25053565b12SHong Zhang     /* Note: a array does not need sorting as it is in dense storage locations. */
25153565b12SHong Zhang     ierr = PetscSortInt(apnzj,apj);CHKERRQ(ierr);
25253565b12SHong Zhang 
25353565b12SHong Zhang     /* Compute P^T*A*P using outer product (P^T)[:,j]*(A*P)[j,:]. */
25453565b12SHong Zhang     pnzi = pi[i+1] - pi[i];
25553565b12SHong Zhang     for (j=0; j<pnzi; j++) {
25653565b12SHong Zhang       nextap = 0;
25753565b12SHong Zhang       crow   = *pJ++;
25853565b12SHong Zhang       cjj    = cj + ci[crow];
25953565b12SHong Zhang       caj    = ca + ci[crow];
26053565b12SHong Zhang       /* Perform sparse axpy operation.  Note cjj includes apj. */
26153565b12SHong Zhang       for (k=0; nextap<apnzj; k++) {
262*cf9c20a2SJed Brown         if (PetscUnlikelyDebug(k >= ci[crow+1] - ci[crow])) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"k too large k %d, crow %d",k,crow);
26353565b12SHong Zhang         if (cjj[k]==apj[nextap]) {
26453565b12SHong Zhang           caj[k] += (*pA)*apa[apj[nextap++]];
26553565b12SHong Zhang         }
26653565b12SHong Zhang       }
26753565b12SHong Zhang       ierr = PetscLogFlops(2.0*apnzj);CHKERRQ(ierr);
26853565b12SHong Zhang       pA++;
26953565b12SHong Zhang     }
27053565b12SHong Zhang 
27153565b12SHong Zhang     /* Zero the current row info for A*P */
27253565b12SHong Zhang     for (j=0; j<apnzj; j++) {
27353565b12SHong Zhang       apa[apj[j]]      = 0.;
27453565b12SHong Zhang       apjdense[apj[j]] = 0;
27553565b12SHong Zhang     }
27653565b12SHong Zhang   }
27753565b12SHong Zhang 
27853565b12SHong Zhang   /* Assemble the final matrix and clean up */
27953565b12SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
28053565b12SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
281ae32dd66SHong Zhang 
282857a15f1SBarry Smith   ierr = PetscFree2(apa,apjdense);CHKERRQ(ierr);
283857a15f1SBarry Smith   ierr = PetscFree(apj);CHKERRQ(ierr);
28453565b12SHong Zhang   PetscFunctionReturn(0);
28553565b12SHong Zhang }
28653565b12SHong Zhang 
287dfbe8321SBarry Smith PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat A,Mat P,Mat C)
288dfbe8321SBarry Smith {
289dfbe8321SBarry Smith   PetscErrorCode      ierr;
290d20bfe6fSHong Zhang   Mat_SeqAIJ          *c = (Mat_SeqAIJ*)C->data;
2918fa4b5a6SHong Zhang   Mat_MatTransMatMult *atb = c->atb;
2928fa4b5a6SHong Zhang   Mat                 Pt = atb->At;
293eb9c0419SKris Buschelman 
294eb9c0419SKris Buschelman   PetscFunctionBegin;
2958fa4b5a6SHong Zhang   ierr = MatTranspose_SeqAIJ(P,MAT_REUSE_MATRIX,&Pt);CHKERRQ(ierr);
2968fa4b5a6SHong Zhang   ierr = (C->ops->matmatmultnumeric)(Pt,A,P,C);CHKERRQ(ierr);
29753565b12SHong Zhang   PetscFunctionReturn(0);
29853565b12SHong Zhang }
299