xref: /petsc/src/mat/impls/aij/seq/matptap.c (revision 8fa4b5a6cf9ec79f104c6bc6f70641dcee59642e)
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)
20*8fa4b5a6SHong Zhang   const char          *algTypes[2] = {"scalable","rap"};
216f231fbdSstefano_zampini   PetscInt            nalg = 2;
226f231fbdSstefano_zampini #else
23*8fa4b5a6SHong Zhang   const char          *algTypes[3] = {"scalable","rap","hypre"};
246f231fbdSstefano_zampini   PetscInt            nalg = 3;
256f231fbdSstefano_zampini #endif
26*8fa4b5a6SHong Zhang   PetscInt            alg = 1; /* set default algorithm */
27*8fa4b5a6SHong Zhang   Mat                 Pt;
28*8fa4b5a6SHong Zhang   Mat_MatTransMatMult *atb;
29*8fa4b5a6SHong 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:
35*8fa4b5a6SHong Zhang        "rap":      Pt = P^T and C = Pt*A*P
36*8fa4b5a6SHong 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:
45*8fa4b5a6SHong Zhang       ierr = PetscNew(&atb);CHKERRQ(ierr);
46*8fa4b5a6SHong Zhang       ierr = MatTranspose_SeqAIJ(P,MAT_INITIAL_MATRIX,&Pt);CHKERRQ(ierr);
47*8fa4b5a6SHong Zhang       ierr = MatMatMatMult(Pt,A,P,MAT_INITIAL_MATRIX,fill,C);CHKERRQ(ierr);
48*8fa4b5a6SHong Zhang 
49*8fa4b5a6SHong Zhang       c                      = (Mat_SeqAIJ*)(*C)->data;
50*8fa4b5a6SHong Zhang       c->atb                 = atb;
51*8fa4b5a6SHong Zhang       atb->At                = Pt;
52*8fa4b5a6SHong Zhang       atb->destroy           = (*C)->ops->destroy;
53*8fa4b5a6SHong Zhang       (*C)->ops->destroy     = MatDestroy_SeqAIJ_MatTransMatMult;
54*8fa4b5a6SHong Zhang       (*C)->ops->ptapnumeric = MatPtAPNumeric_SeqAIJ_SeqAIJ;
55*8fa4b5a6SHong Zhang       PetscFunctionReturn(0);
564a1b09b7SHong Zhang       break;
576f231fbdSstefano_zampini #if defined(PETSC_HAVE_HYPRE)
586f231fbdSstefano_zampini     case 2:
596f231fbdSstefano_zampini       ierr = MatPtAPSymbolic_AIJ_AIJ_wHYPRE(A,P,fill,C);CHKERRQ(ierr);
606f231fbdSstefano_zampini       break;
616f231fbdSstefano_zampini #endif
624a1b09b7SHong Zhang     default:
634a1b09b7SHong Zhang       ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(A,P,fill,C);CHKERRQ(ierr);
644a1b09b7SHong Zhang       break;
654a1b09b7SHong Zhang     }
667ba1a0bfSKris Buschelman   }
67ae32dd66SHong Zhang   ierr = (*(*C)->ops->ptapnumeric)(A,P,*C);CHKERRQ(ierr);
687ba1a0bfSKris Buschelman   PetscFunctionReturn(0);
697ba1a0bfSKris Buschelman }
707ba1a0bfSKris Buschelman 
7153565b12SHong Zhang PetscErrorCode MatDestroy_SeqAIJ_PtAP(Mat A)
728d0a38e4SHong Zhang {
738d0a38e4SHong Zhang   PetscErrorCode ierr;
7453565b12SHong Zhang   Mat_SeqAIJ     *a    = (Mat_SeqAIJ*)A->data;
7553565b12SHong Zhang   Mat_PtAP       *ptap = a->ptap;
768d0a38e4SHong Zhang 
778d0a38e4SHong Zhang   PetscFunctionBegin;
788d0a38e4SHong Zhang   ierr = PetscFree(ptap->apa);CHKERRQ(ierr);
79f79958ebSHong Zhang   ierr = PetscFree(ptap->api);CHKERRQ(ierr);
80f79958ebSHong Zhang   ierr = PetscFree(ptap->apj);CHKERRQ(ierr);
8153565b12SHong Zhang   ierr = (ptap->destroy)(A);CHKERRQ(ierr);
828d0a38e4SHong Zhang   ierr = PetscFree(ptap);CHKERRQ(ierr);
838d0a38e4SHong Zhang   PetscFunctionReturn(0);
848d0a38e4SHong Zhang }
858d0a38e4SHong Zhang 
86ae32dd66SHong Zhang PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(Mat A,Mat P,PetscReal fill,Mat *C)
878d0a38e4SHong Zhang {
888d0a38e4SHong Zhang   PetscErrorCode     ierr;
890298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
90d20bfe6fSHong Zhang   Mat_SeqAIJ         *a        = (Mat_SeqAIJ*)A->data,*p = (Mat_SeqAIJ*)P->data,*c;
9153565b12SHong Zhang   PetscInt           *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj;
9253565b12SHong Zhang   PetscInt           *ci,*cj,*ptadenserow,*ptasparserow,*ptaj,nspacedouble=0;
937d69fa63SHong Zhang   PetscInt           an=A->cmap->N,am=A->rmap->N,pn=P->cmap->N,pm=P->rmap->N;
9453565b12SHong Zhang   PetscInt           i,j,k,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi,nlnk,*lnk;
95d20bfe6fSHong Zhang   MatScalar          *ca;
9653565b12SHong Zhang   PetscBT            lnkbt;
977d69fa63SHong Zhang   PetscReal          afill;
98eb9c0419SKris Buschelman 
99eb9c0419SKris Buschelman   PetscFunctionBegin;
10053565b12SHong Zhang   /* Get ij structure of P^T */
101eb9c0419SKris Buschelman   ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
102d20bfe6fSHong Zhang   ptJ  = ptj;
103eb9c0419SKris Buschelman 
104d20bfe6fSHong Zhang   /* Allocate ci array, arrays for fill computation and */
105d20bfe6fSHong Zhang   /* free space for accumulating nonzero column info */
106854ce69bSBarry Smith   ierr  = PetscMalloc1(pn+1,&ci);CHKERRQ(ierr);
107d20bfe6fSHong Zhang   ci[0] = 0;
108eb9c0419SKris Buschelman 
1091795a4d1SJed Brown   ierr         = PetscCalloc1(2*an+1,&ptadenserow);CHKERRQ(ierr);
110d20bfe6fSHong Zhang   ptasparserow = ptadenserow  + an;
11155a3bba9SHong Zhang 
11255a3bba9SHong Zhang   /* create and initialize a linked list */
11355a3bba9SHong Zhang   nlnk = pn+1;
11455a3bba9SHong Zhang   ierr = PetscLLCreate(pn,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr);
115eb9c0419SKris Buschelman 
1167d69fa63SHong Zhang   /* Set initial free space to be fill*(nnz(A)+ nnz(P)) */
117f91af8c7SBarry Smith   ierr          = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],pi[pm])),&free_space);CHKERRQ(ierr);
118d20bfe6fSHong Zhang   current_space = free_space;
119d20bfe6fSHong Zhang 
120d20bfe6fSHong Zhang   /* Determine symbolic info for each row of C: */
121d20bfe6fSHong Zhang   for (i=0; i<pn; i++) {
122d20bfe6fSHong Zhang     ptnzi  = pti[i+1] - pti[i];
123d20bfe6fSHong Zhang     ptanzi = 0;
124d20bfe6fSHong Zhang     /* Determine symbolic row of PtA: */
125d20bfe6fSHong Zhang     for (j=0; j<ptnzi; j++) {
126d20bfe6fSHong Zhang       arow = *ptJ++;
127d20bfe6fSHong Zhang       anzj = ai[arow+1] - ai[arow];
128d20bfe6fSHong Zhang       ajj  = aj + ai[arow];
129d20bfe6fSHong Zhang       for (k=0; k<anzj; k++) {
130d20bfe6fSHong Zhang         if (!ptadenserow[ajj[k]]) {
131d20bfe6fSHong Zhang           ptadenserow[ajj[k]]    = -1;
132d20bfe6fSHong Zhang           ptasparserow[ptanzi++] = ajj[k];
133d20bfe6fSHong Zhang         }
134d20bfe6fSHong Zhang       }
135d20bfe6fSHong Zhang     }
136d20bfe6fSHong Zhang     /* Using symbolic info for row of PtA, determine symbolic info for row of C: */
137d20bfe6fSHong Zhang     ptaj = ptasparserow;
138d20bfe6fSHong Zhang     cnzi = 0;
139d20bfe6fSHong Zhang     for (j=0; j<ptanzi; j++) {
140d20bfe6fSHong Zhang       prow = *ptaj++;
141d20bfe6fSHong Zhang       pnzj = pi[prow+1] - pi[prow];
142d20bfe6fSHong Zhang       pjj  = pj + pi[prow];
14355a3bba9SHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
144dadf0e6bSHong Zhang       ierr  = PetscLLAddSorted(pnzj,pjj,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr);
14555a3bba9SHong Zhang       cnzi += nlnk;
146d20bfe6fSHong Zhang     }
147d20bfe6fSHong Zhang 
148d20bfe6fSHong Zhang     /* If free space is not available, make more free space */
149d20bfe6fSHong Zhang     /* Double the amount of total space in the list */
150d20bfe6fSHong Zhang     if (current_space->local_remaining<cnzi) {
151f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(cnzi,current_space->total_array_size),&current_space);CHKERRQ(ierr);
152f2b054eeSHong Zhang       nspacedouble++;
153d20bfe6fSHong Zhang     }
154d20bfe6fSHong Zhang 
155d20bfe6fSHong Zhang     /* Copy data into free space, and zero out denserows */
15655a3bba9SHong Zhang     ierr = PetscLLClean(pn,pn,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
1572205254eSKarl Rupp 
158d20bfe6fSHong Zhang     current_space->array           += cnzi;
159d20bfe6fSHong Zhang     current_space->local_used      += cnzi;
160d20bfe6fSHong Zhang     current_space->local_remaining -= cnzi;
161d20bfe6fSHong Zhang 
1622205254eSKarl Rupp     for (j=0; j<ptanzi; j++) ptadenserow[ptasparserow[j]] = 0;
1632205254eSKarl Rupp 
164d20bfe6fSHong Zhang     /* Aside: Perhaps we should save the pta info for the numerical factorization. */
165d20bfe6fSHong Zhang     /*        For now, we will recompute what is needed. */
166d20bfe6fSHong Zhang     ci[i+1] = ci[i] + cnzi;
167d20bfe6fSHong Zhang   }
168d20bfe6fSHong Zhang   /* nnz is now stored in ci[ptm], column indices are in the list of free space */
169d20bfe6fSHong Zhang   /* Allocate space for cj, initialize cj, and */
170d20bfe6fSHong Zhang   /* destroy list of free space and other temporary array(s) */
171854ce69bSBarry Smith   ierr = PetscMalloc1(ci[pn]+1,&cj);CHKERRQ(ierr);
172a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
173d20bfe6fSHong Zhang   ierr = PetscFree(ptadenserow);CHKERRQ(ierr);
17455a3bba9SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
175d20bfe6fSHong Zhang 
176854ce69bSBarry Smith   ierr = PetscCalloc1(ci[pn]+1,&ca);CHKERRQ(ierr);
17753565b12SHong Zhang 
17853565b12SHong Zhang   /* put together the new matrix */
179ce94432eSBarry Smith   ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),pn,pn,ci,cj,ca,C);CHKERRQ(ierr);
18033d57670SJed Brown   ierr = MatSetBlockSizes(*C,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr);
181ae32dd66SHong Zhang 
18253565b12SHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
18353565b12SHong Zhang   /* Since these are PETSc arrays, change flags to free them as necessary. */
18453565b12SHong Zhang   c          = (Mat_SeqAIJ*)((*C)->data);
18553565b12SHong Zhang   c->free_a  = PETSC_TRUE;
18653565b12SHong Zhang   c->free_ij = PETSC_TRUE;
18753565b12SHong Zhang   c->nonew   = 0;
188f4a1e0b0SHong Zhang   (*C)->ops->ptapnumeric = MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy;
18953565b12SHong Zhang 
1907d69fa63SHong Zhang   /* set MatInfo */
1917d69fa63SHong Zhang   afill = (PetscReal)ci[pn]/(ai[am]+pi[pm] + 1.e-5);
1927d69fa63SHong Zhang   if (afill < 1.0) afill = 1.0;
1937d69fa63SHong Zhang   c->maxnz                     = ci[pn];
1947d69fa63SHong Zhang   c->nz                        = ci[pn];
1957d69fa63SHong Zhang   (*C)->info.mallocs           = nspacedouble;
1967d69fa63SHong Zhang   (*C)->info.fill_ratio_given  = fill;
1977d69fa63SHong Zhang   (*C)->info.fill_ratio_needed = afill;
1987d69fa63SHong Zhang 
19953565b12SHong Zhang   /* Clean up. */
20053565b12SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
20153565b12SHong Zhang #if defined(PETSC_USE_INFO)
20253565b12SHong Zhang   if (ci[pn] != 0) {
20357622a8eSBarry Smith     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
20457622a8eSBarry Smith     ierr = PetscInfo1((*C),"Use MatPtAP(A,P,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr);
20553565b12SHong Zhang   } else {
20653565b12SHong Zhang     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
20753565b12SHong Zhang   }
208f79958ebSHong Zhang #endif
20953565b12SHong Zhang   PetscFunctionReturn(0);
21053565b12SHong Zhang }
21153565b12SHong Zhang 
212ae32dd66SHong Zhang PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat A,Mat P,Mat C)
21353565b12SHong Zhang {
21453565b12SHong Zhang   PetscErrorCode ierr;
21553565b12SHong Zhang   Mat_SeqAIJ     *a = (Mat_SeqAIJ*) A->data;
21653565b12SHong Zhang   Mat_SeqAIJ     *p = (Mat_SeqAIJ*) P->data;
21753565b12SHong Zhang   Mat_SeqAIJ     *c = (Mat_SeqAIJ*) C->data;
21853565b12SHong Zhang   PetscInt       *ai=a->i,*aj=a->j,*apj,*apjdense,*pi=p->i,*pj=p->j,*pJ=p->j,*pjj;
21953565b12SHong Zhang   PetscInt       *ci=c->i,*cj=c->j,*cjj;
22053565b12SHong Zhang   PetscInt       am =A->rmap->N,cn=C->cmap->N,cm=C->rmap->N;
22153565b12SHong Zhang   PetscInt       i,j,k,anzi,pnzi,apnzj,nextap,pnzj,prow,crow;
22253565b12SHong Zhang   MatScalar      *aa=a->a,*apa,*pa=p->a,*pA=p->a,*paj,*ca=c->a,*caj;
22353565b12SHong Zhang 
22453565b12SHong Zhang   PetscFunctionBegin;
225ae32dd66SHong Zhang   /* Allocate temporary array for storage of one row of A*P (cn: non-scalable) */
2264bcca364SHong Zhang   ierr = PetscMalloc3(cn,&apa,cn,&apjdense,cn,&apj);CHKERRQ(ierr);
227eb9baa12SBarry Smith   ierr = PetscMemzero(apa,cn*sizeof(MatScalar));CHKERRQ(ierr);
228eb9baa12SBarry Smith   ierr = PetscMemzero(apjdense,cn*sizeof(PetscInt));CHKERRQ(ierr);
22953565b12SHong Zhang 
23053565b12SHong Zhang   /* Clear old values in C */
23153565b12SHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
23253565b12SHong Zhang 
23353565b12SHong Zhang   for (i=0; i<am; i++) {
23453565b12SHong Zhang     /* Form sparse row of A*P */
23553565b12SHong Zhang     anzi  = ai[i+1] - ai[i];
23653565b12SHong Zhang     apnzj = 0;
23753565b12SHong Zhang     for (j=0; j<anzi; j++) {
23853565b12SHong Zhang       prow = *aj++;
23953565b12SHong Zhang       pnzj = pi[prow+1] - pi[prow];
24053565b12SHong Zhang       pjj  = pj + pi[prow];
24153565b12SHong Zhang       paj  = pa + pi[prow];
24253565b12SHong Zhang       for (k=0; k<pnzj; k++) {
24353565b12SHong Zhang         if (!apjdense[pjj[k]]) {
24453565b12SHong Zhang           apjdense[pjj[k]] = -1;
24553565b12SHong Zhang           apj[apnzj++]     = pjj[k];
24653565b12SHong Zhang         }
24753565b12SHong Zhang         apa[pjj[k]] += (*aa)*paj[k];
24853565b12SHong Zhang       }
24953565b12SHong Zhang       ierr = PetscLogFlops(2.0*pnzj);CHKERRQ(ierr);
25053565b12SHong Zhang       aa++;
25153565b12SHong Zhang     }
25253565b12SHong Zhang 
25353565b12SHong Zhang     /* Sort the j index array for quick sparse axpy. */
25453565b12SHong Zhang     /* Note: a array does not need sorting as it is in dense storage locations. */
25553565b12SHong Zhang     ierr = PetscSortInt(apnzj,apj);CHKERRQ(ierr);
25653565b12SHong Zhang 
25753565b12SHong Zhang     /* Compute P^T*A*P using outer product (P^T)[:,j]*(A*P)[j,:]. */
25853565b12SHong Zhang     pnzi = pi[i+1] - pi[i];
25953565b12SHong Zhang     for (j=0; j<pnzi; j++) {
26053565b12SHong Zhang       nextap = 0;
26153565b12SHong Zhang       crow   = *pJ++;
26253565b12SHong Zhang       cjj    = cj + ci[crow];
26353565b12SHong Zhang       caj    = ca + ci[crow];
26453565b12SHong Zhang       /* Perform sparse axpy operation.  Note cjj includes apj. */
26553565b12SHong Zhang       for (k=0; nextap<apnzj; k++) {
26653565b12SHong Zhang #if defined(PETSC_USE_DEBUG)
267f23aa3ddSBarry Smith         if (k >= ci[crow+1] - ci[crow]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"k too large k %d, crow %d",k,crow);
26853565b12SHong Zhang #endif
26953565b12SHong Zhang         if (cjj[k]==apj[nextap]) {
27053565b12SHong Zhang           caj[k] += (*pA)*apa[apj[nextap++]];
27153565b12SHong Zhang         }
27253565b12SHong Zhang       }
27353565b12SHong Zhang       ierr = PetscLogFlops(2.0*apnzj);CHKERRQ(ierr);
27453565b12SHong Zhang       pA++;
27553565b12SHong Zhang     }
27653565b12SHong Zhang 
27753565b12SHong Zhang     /* Zero the current row info for A*P */
27853565b12SHong Zhang     for (j=0; j<apnzj; j++) {
27953565b12SHong Zhang       apa[apj[j]]      = 0.;
28053565b12SHong Zhang       apjdense[apj[j]] = 0;
28153565b12SHong Zhang     }
28253565b12SHong Zhang   }
28353565b12SHong Zhang 
28453565b12SHong Zhang   /* Assemble the final matrix and clean up */
28553565b12SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
28653565b12SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
287ae32dd66SHong Zhang 
288eb9baa12SBarry Smith   ierr = PetscFree3(apa,apjdense,apj);CHKERRQ(ierr);
28953565b12SHong Zhang   PetscFunctionReturn(0);
29053565b12SHong Zhang }
29153565b12SHong Zhang 
292dfbe8321SBarry Smith PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat A,Mat P,Mat C)
293dfbe8321SBarry Smith {
294dfbe8321SBarry Smith   PetscErrorCode      ierr;
295d20bfe6fSHong Zhang   Mat_SeqAIJ          *c = (Mat_SeqAIJ*)C->data;
296*8fa4b5a6SHong Zhang   Mat_MatTransMatMult *atb = c->atb;
297*8fa4b5a6SHong Zhang   Mat                 Pt = atb->At;
298eb9c0419SKris Buschelman 
299eb9c0419SKris Buschelman   PetscFunctionBegin;
300*8fa4b5a6SHong Zhang   ierr = MatTranspose_SeqAIJ(P,MAT_REUSE_MATRIX,&Pt);CHKERRQ(ierr);
301*8fa4b5a6SHong Zhang   ierr = (C->ops->matmatmultnumeric)(Pt,A,P,C);CHKERRQ(ierr);
30253565b12SHong Zhang   PetscFunctionReturn(0);
30353565b12SHong Zhang }
304