xref: /petsc/src/mat/impls/aij/seq/matptap.c (revision 65e4b4d46e77b7a054a0705b877751ea429d3d2b)
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;
257ba1a0bfSKris Buschelman 
267ba1a0bfSKris Buschelman   PetscFunctionBegin;
274222ddf1SHong Zhang   /* "scalable" */
284222ddf1SHong Zhang   ierr = PetscStrcmp(alg,"scalable",&flg);CHKERRQ(ierr);
294222ddf1SHong Zhang   if (flg) {
304222ddf1SHong Zhang     ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(A,P,fill,C);CHKERRQ(ierr);
314222ddf1SHong Zhang     C->ops->productnumeric = MatProductNumeric_PtAP;
324222ddf1SHong Zhang     PetscFunctionReturn(0);
334222ddf1SHong Zhang   }
344222ddf1SHong Zhang 
354222ddf1SHong Zhang   /* "rap" */
364222ddf1SHong Zhang   ierr = PetscStrcmp(alg,"rap",&flg);CHKERRQ(ierr);
376718818eSStefano Zampini   if (flg) {
386718818eSStefano Zampini     Mat_MatTransMatMult *atb;
396718818eSStefano Zampini 
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 
448fa4b5a6SHong Zhang     atb->At                = Pt;
456718818eSStefano Zampini     atb->data              = C->product->data;
466718818eSStefano Zampini     atb->destroy           = C->product->destroy;
476718818eSStefano Zampini     C->product->data       = atb;
486718818eSStefano Zampini     C->product->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 }
657ba1a0bfSKris Buschelman 
664222ddf1SHong Zhang PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(Mat A,Mat P,PetscReal fill,Mat C)
678d0a38e4SHong Zhang {
688d0a38e4SHong Zhang   PetscErrorCode     ierr;
690298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
70d20bfe6fSHong Zhang   Mat_SeqAIJ         *a        = (Mat_SeqAIJ*)A->data,*p = (Mat_SeqAIJ*)P->data,*c;
7153565b12SHong Zhang   PetscInt           *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj;
7253565b12SHong Zhang   PetscInt           *ci,*cj,*ptadenserow,*ptasparserow,*ptaj,nspacedouble=0;
737d69fa63SHong Zhang   PetscInt           an=A->cmap->N,am=A->rmap->N,pn=P->cmap->N,pm=P->rmap->N;
7453565b12SHong Zhang   PetscInt           i,j,k,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi,nlnk,*lnk;
75d20bfe6fSHong Zhang   MatScalar          *ca;
7653565b12SHong Zhang   PetscBT            lnkbt;
777d69fa63SHong Zhang   PetscReal          afill;
78eb9c0419SKris Buschelman 
79eb9c0419SKris Buschelman   PetscFunctionBegin;
8053565b12SHong Zhang   /* Get ij structure of P^T */
81eb9c0419SKris Buschelman   ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
82d20bfe6fSHong Zhang   ptJ  = ptj;
83eb9c0419SKris Buschelman 
84d20bfe6fSHong Zhang   /* Allocate ci array, arrays for fill computation and */
85d20bfe6fSHong Zhang   /* free space for accumulating nonzero column info */
86854ce69bSBarry Smith   ierr  = PetscMalloc1(pn+1,&ci);CHKERRQ(ierr);
87d20bfe6fSHong Zhang   ci[0] = 0;
88eb9c0419SKris Buschelman 
891795a4d1SJed Brown   ierr         = PetscCalloc1(2*an+1,&ptadenserow);CHKERRQ(ierr);
90d20bfe6fSHong Zhang   ptasparserow = ptadenserow  + an;
9155a3bba9SHong Zhang 
9255a3bba9SHong Zhang   /* create and initialize a linked list */
9355a3bba9SHong Zhang   nlnk = pn+1;
9455a3bba9SHong Zhang   ierr = PetscLLCreate(pn,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr);
95eb9c0419SKris Buschelman 
967d69fa63SHong Zhang   /* Set initial free space to be fill*(nnz(A)+ nnz(P)) */
97f91af8c7SBarry Smith   ierr          = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],pi[pm])),&free_space);CHKERRQ(ierr);
98d20bfe6fSHong Zhang   current_space = free_space;
99d20bfe6fSHong Zhang 
100d20bfe6fSHong Zhang   /* Determine symbolic info for each row of C: */
101d20bfe6fSHong Zhang   for (i=0; i<pn; i++) {
102d20bfe6fSHong Zhang     ptnzi  = pti[i+1] - pti[i];
103d20bfe6fSHong Zhang     ptanzi = 0;
104d20bfe6fSHong Zhang     /* Determine symbolic row of PtA: */
105d20bfe6fSHong Zhang     for (j=0; j<ptnzi; j++) {
106d20bfe6fSHong Zhang       arow = *ptJ++;
107d20bfe6fSHong Zhang       anzj = ai[arow+1] - ai[arow];
108d20bfe6fSHong Zhang       ajj  = aj + ai[arow];
109d20bfe6fSHong Zhang       for (k=0; k<anzj; k++) {
110d20bfe6fSHong Zhang         if (!ptadenserow[ajj[k]]) {
111d20bfe6fSHong Zhang           ptadenserow[ajj[k]]    = -1;
112d20bfe6fSHong Zhang           ptasparserow[ptanzi++] = ajj[k];
113d20bfe6fSHong Zhang         }
114d20bfe6fSHong Zhang       }
115d20bfe6fSHong Zhang     }
116d20bfe6fSHong Zhang     /* Using symbolic info for row of PtA, determine symbolic info for row of C: */
117d20bfe6fSHong Zhang     ptaj = ptasparserow;
118d20bfe6fSHong Zhang     cnzi = 0;
119d20bfe6fSHong Zhang     for (j=0; j<ptanzi; j++) {
120d20bfe6fSHong Zhang       prow = *ptaj++;
121d20bfe6fSHong Zhang       pnzj = pi[prow+1] - pi[prow];
122d20bfe6fSHong Zhang       pjj  = pj + pi[prow];
12355a3bba9SHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
124dadf0e6bSHong Zhang       ierr  = PetscLLAddSorted(pnzj,pjj,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr);
12555a3bba9SHong Zhang       cnzi += nlnk;
126d20bfe6fSHong Zhang     }
127d20bfe6fSHong Zhang 
128d20bfe6fSHong Zhang     /* If free space is not available, make more free space */
129d20bfe6fSHong Zhang     /* Double the amount of total space in the list */
130d20bfe6fSHong Zhang     if (current_space->local_remaining<cnzi) {
131f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(cnzi,current_space->total_array_size),&current_space);CHKERRQ(ierr);
132f2b054eeSHong Zhang       nspacedouble++;
133d20bfe6fSHong Zhang     }
134d20bfe6fSHong Zhang 
135d20bfe6fSHong Zhang     /* Copy data into free space, and zero out denserows */
13655a3bba9SHong Zhang     ierr = PetscLLClean(pn,pn,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
1372205254eSKarl Rupp 
138d20bfe6fSHong Zhang     current_space->array           += cnzi;
139d20bfe6fSHong Zhang     current_space->local_used      += cnzi;
140d20bfe6fSHong Zhang     current_space->local_remaining -= cnzi;
141d20bfe6fSHong Zhang 
1422205254eSKarl Rupp     for (j=0; j<ptanzi; j++) ptadenserow[ptasparserow[j]] = 0;
1432205254eSKarl Rupp 
144d20bfe6fSHong Zhang     /* Aside: Perhaps we should save the pta info for the numerical factorization. */
145d20bfe6fSHong Zhang     /*        For now, we will recompute what is needed. */
146d20bfe6fSHong Zhang     ci[i+1] = ci[i] + cnzi;
147d20bfe6fSHong Zhang   }
148d20bfe6fSHong Zhang   /* nnz is now stored in ci[ptm], column indices are in the list of free space */
149d20bfe6fSHong Zhang   /* Allocate space for cj, initialize cj, and */
150d20bfe6fSHong Zhang   /* destroy list of free space and other temporary array(s) */
151854ce69bSBarry Smith   ierr = PetscMalloc1(ci[pn]+1,&cj);CHKERRQ(ierr);
152a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
153d20bfe6fSHong Zhang   ierr = PetscFree(ptadenserow);CHKERRQ(ierr);
15455a3bba9SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
155d20bfe6fSHong Zhang 
156854ce69bSBarry Smith   ierr = PetscCalloc1(ci[pn]+1,&ca);CHKERRQ(ierr);
15753565b12SHong Zhang 
15853565b12SHong Zhang   /* put together the new matrix */
159e4e71118SStefano Zampini   ierr = MatSetSeqAIJWithArrays_private(PetscObjectComm((PetscObject)A),pn,pn,ci,cj,ca,((PetscObject)A)->type_name,C);CHKERRQ(ierr);
1604222ddf1SHong Zhang   ierr = MatSetBlockSizes(C,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr);
161ae32dd66SHong Zhang 
16253565b12SHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
16353565b12SHong Zhang   /* Since these are PETSc arrays, change flags to free them as necessary. */
1644222ddf1SHong Zhang   c          = (Mat_SeqAIJ*)((C)->data);
16553565b12SHong Zhang   c->free_a  = PETSC_TRUE;
16653565b12SHong Zhang   c->free_ij = PETSC_TRUE;
16753565b12SHong Zhang   c->nonew   = 0;
1686718818eSStefano Zampini 
1694222ddf1SHong Zhang   C->ops->ptapnumeric = MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy;
17053565b12SHong Zhang 
1717d69fa63SHong Zhang   /* set MatInfo */
1727d69fa63SHong Zhang   afill = (PetscReal)ci[pn]/(ai[am]+pi[pm] + 1.e-5);
1737d69fa63SHong Zhang   if (afill < 1.0) afill = 1.0;
1747d69fa63SHong Zhang   c->maxnz                     = ci[pn];
1757d69fa63SHong Zhang   c->nz                        = ci[pn];
1764222ddf1SHong Zhang   C->info.mallocs           = nspacedouble;
1774222ddf1SHong Zhang   C->info.fill_ratio_given  = fill;
1784222ddf1SHong Zhang   C->info.fill_ratio_needed = afill;
1797d69fa63SHong Zhang 
18053565b12SHong Zhang   /* Clean up. */
18153565b12SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
18253565b12SHong Zhang #if defined(PETSC_USE_INFO)
18353565b12SHong Zhang   if (ci[pn] != 0) {
1844222ddf1SHong Zhang     ierr = PetscInfo3(C,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
1854222ddf1SHong Zhang     ierr = PetscInfo1(C,"Use MatPtAP(A,P,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr);
18653565b12SHong Zhang   } else {
1874222ddf1SHong Zhang     ierr = PetscInfo(C,"Empty matrix product\n");CHKERRQ(ierr);
18853565b12SHong Zhang   }
189f79958ebSHong Zhang #endif
19053565b12SHong Zhang   PetscFunctionReturn(0);
19153565b12SHong Zhang }
19253565b12SHong Zhang 
193ae32dd66SHong Zhang PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat A,Mat P,Mat C)
19453565b12SHong Zhang {
19553565b12SHong Zhang   PetscErrorCode ierr;
19653565b12SHong Zhang   Mat_SeqAIJ     *a = (Mat_SeqAIJ*) A->data;
19753565b12SHong Zhang   Mat_SeqAIJ     *p = (Mat_SeqAIJ*) P->data;
19853565b12SHong Zhang   Mat_SeqAIJ     *c = (Mat_SeqAIJ*) C->data;
19953565b12SHong Zhang   PetscInt       *ai=a->i,*aj=a->j,*apj,*apjdense,*pi=p->i,*pj=p->j,*pJ=p->j,*pjj;
20053565b12SHong Zhang   PetscInt       *ci=c->i,*cj=c->j,*cjj;
20153565b12SHong Zhang   PetscInt       am =A->rmap->N,cn=C->cmap->N,cm=C->rmap->N;
20253565b12SHong Zhang   PetscInt       i,j,k,anzi,pnzi,apnzj,nextap,pnzj,prow,crow;
203*65e4b4d4SStefano Zampini   MatScalar      *aa,*apa,*pa,*pA,*paj,*ca,*caj;
20453565b12SHong Zhang 
20553565b12SHong Zhang   PetscFunctionBegin;
206ae32dd66SHong Zhang   /* Allocate temporary array for storage of one row of A*P (cn: non-scalable) */
2078f8f2f0dSBarry Smith   ierr = PetscCalloc2(cn,&apa,cn,&apjdense);CHKERRQ(ierr);
208857a15f1SBarry Smith   ierr = PetscMalloc1(cn,&apj);CHKERRQ(ierr);
209*65e4b4d4SStefano Zampini   /* trigger CPU copies if needed and flag CPU mask for C */
210*65e4b4d4SStefano Zampini #if defined(PETSC_HAVE_DEVICE)
211*65e4b4d4SStefano Zampini   {
212*65e4b4d4SStefano Zampini     const PetscScalar *dummy;
213*65e4b4d4SStefano Zampini     ierr = MatSeqAIJGetArrayRead(A,&dummy);CHKERRQ(ierr);
214*65e4b4d4SStefano Zampini     ierr = MatSeqAIJRestoreArrayRead(A,&dummy);CHKERRQ(ierr);
215*65e4b4d4SStefano Zampini     ierr = MatSeqAIJGetArrayRead(P,&dummy);CHKERRQ(ierr);
216*65e4b4d4SStefano Zampini     ierr = MatSeqAIJRestoreArrayRead(P,&dummy);CHKERRQ(ierr);
217*65e4b4d4SStefano Zampini     if (C->offloadmask != PETSC_OFFLOAD_UNALLOCATED) C->offloadmask = PETSC_OFFLOAD_CPU;
218*65e4b4d4SStefano Zampini   }
219*65e4b4d4SStefano Zampini #endif
220*65e4b4d4SStefano Zampini   aa = a->a;
221*65e4b4d4SStefano Zampini   pa = p->a;
222*65e4b4d4SStefano Zampini   pA = p->a;
223*65e4b4d4SStefano Zampini   ca = c->a;
22453565b12SHong Zhang 
22553565b12SHong Zhang   /* Clear old values in C */
226580bdb30SBarry Smith   ierr = PetscArrayzero(ca,ci[cm]);CHKERRQ(ierr);
22753565b12SHong Zhang 
22853565b12SHong Zhang   for (i=0; i<am; i++) {
22953565b12SHong Zhang     /* Form sparse row of A*P */
23053565b12SHong Zhang     anzi  = ai[i+1] - ai[i];
23153565b12SHong Zhang     apnzj = 0;
23253565b12SHong Zhang     for (j=0; j<anzi; j++) {
23353565b12SHong Zhang       prow = *aj++;
23453565b12SHong Zhang       pnzj = pi[prow+1] - pi[prow];
23553565b12SHong Zhang       pjj  = pj + pi[prow];
23653565b12SHong Zhang       paj  = pa + pi[prow];
23753565b12SHong Zhang       for (k=0; k<pnzj; k++) {
23853565b12SHong Zhang         if (!apjdense[pjj[k]]) {
23953565b12SHong Zhang           apjdense[pjj[k]] = -1;
24053565b12SHong Zhang           apj[apnzj++]     = pjj[k];
24153565b12SHong Zhang         }
24253565b12SHong Zhang         apa[pjj[k]] += (*aa)*paj[k];
24353565b12SHong Zhang       }
24453565b12SHong Zhang       ierr = PetscLogFlops(2.0*pnzj);CHKERRQ(ierr);
24553565b12SHong Zhang       aa++;
24653565b12SHong Zhang     }
24753565b12SHong Zhang 
24853565b12SHong Zhang     /* Sort the j index array for quick sparse axpy. */
24953565b12SHong Zhang     /* Note: a array does not need sorting as it is in dense storage locations. */
25053565b12SHong Zhang     ierr = PetscSortInt(apnzj,apj);CHKERRQ(ierr);
25153565b12SHong Zhang 
25253565b12SHong Zhang     /* Compute P^T*A*P using outer product (P^T)[:,j]*(A*P)[j,:]. */
25353565b12SHong Zhang     pnzi = pi[i+1] - pi[i];
25453565b12SHong Zhang     for (j=0; j<pnzi; j++) {
25553565b12SHong Zhang       nextap = 0;
25653565b12SHong Zhang       crow   = *pJ++;
25753565b12SHong Zhang       cjj    = cj + ci[crow];
25853565b12SHong Zhang       caj    = ca + ci[crow];
25953565b12SHong Zhang       /* Perform sparse axpy operation.  Note cjj includes apj. */
26053565b12SHong Zhang       for (k=0; nextap<apnzj; k++) {
261cf9c20a2SJed 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);
26253565b12SHong Zhang         if (cjj[k]==apj[nextap]) {
26353565b12SHong Zhang           caj[k] += (*pA)*apa[apj[nextap++]];
26453565b12SHong Zhang         }
26553565b12SHong Zhang       }
26653565b12SHong Zhang       ierr = PetscLogFlops(2.0*apnzj);CHKERRQ(ierr);
26753565b12SHong Zhang       pA++;
26853565b12SHong Zhang     }
26953565b12SHong Zhang 
27053565b12SHong Zhang     /* Zero the current row info for A*P */
27153565b12SHong Zhang     for (j=0; j<apnzj; j++) {
27253565b12SHong Zhang       apa[apj[j]]      = 0.;
27353565b12SHong Zhang       apjdense[apj[j]] = 0;
27453565b12SHong Zhang     }
27553565b12SHong Zhang   }
27653565b12SHong Zhang 
27753565b12SHong Zhang   /* Assemble the final matrix and clean up */
27853565b12SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
27953565b12SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
280ae32dd66SHong Zhang 
281857a15f1SBarry Smith   ierr = PetscFree2(apa,apjdense);CHKERRQ(ierr);
282857a15f1SBarry Smith   ierr = PetscFree(apj);CHKERRQ(ierr);
28353565b12SHong Zhang   PetscFunctionReturn(0);
28453565b12SHong Zhang }
28553565b12SHong Zhang 
286dfbe8321SBarry Smith PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat A,Mat P,Mat C)
287dfbe8321SBarry Smith {
288dfbe8321SBarry Smith   PetscErrorCode      ierr;
2896718818eSStefano Zampini   Mat_MatTransMatMult *atb;
290eb9c0419SKris Buschelman 
291eb9c0419SKris Buschelman   PetscFunctionBegin;
2926718818eSStefano Zampini   MatCheckProduct(C,3);
2936718818eSStefano Zampini   atb  = (Mat_MatTransMatMult*)C->product->data;
2946718818eSStefano Zampini   if (!atb) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_PLIB,"Missing data structure");
2956718818eSStefano Zampini   ierr = MatTranspose_SeqAIJ(P,MAT_REUSE_MATRIX,&atb->At);CHKERRQ(ierr);
2966718818eSStefano Zampini   if (!C->ops->matmultnumeric) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_PLIB,"Missing numeric operation");
2976718818eSStefano Zampini   /* when using rap, MatMatMatMultSymbolic used a different data */
2986718818eSStefano Zampini   if (atb->data) C->product->data = atb->data;
2996718818eSStefano Zampini   ierr = (*C->ops->matmatmultnumeric)(atb->At,A,P,C);CHKERRQ(ierr);
3006718818eSStefano Zampini   C->product->data = atb;
30153565b12SHong Zhang   PetscFunctionReturn(0);
30253565b12SHong Zhang }
303