xref: /petsc/src/mat/impls/aij/seq/matptap.c (revision ae32dd66e07cee58953b871fdd996e27c5d9a395)
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>
10eb9c0419SKris Buschelman 
11ff134f7aSHong Zhang #undef __FUNCT__
1265e8a0caSHong Zhang #define __FUNCT__ "MatPtAP_SeqAIJ_SeqAIJ"
1365e8a0caSHong Zhang PetscErrorCode MatPtAP_SeqAIJ_SeqAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
147ba1a0bfSKris Buschelman {
157ba1a0bfSKris Buschelman   PetscErrorCode ierr;
167ba1a0bfSKris Buschelman 
177ba1a0bfSKris Buschelman   PetscFunctionBegin;
1865e8a0caSHong Zhang   if (scall == MAT_INITIAL_MATRIX){
1965e8a0caSHong Zhang     ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
2065e8a0caSHong Zhang     ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ(A,P,fill,C);CHKERRQ(ierr);
2165e8a0caSHong Zhang     ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
227ba1a0bfSKris Buschelman   }
2365e8a0caSHong Zhang   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
24*ae32dd66SHong Zhang   //ierr = (*(A)->ops->ptapnumeric)(A,P,*C);CHKERRQ(ierr); /* should use matrix *C, not A! */
25*ae32dd66SHong Zhang   ierr = (*(*C)->ops->ptapnumeric)(A,P,*C);CHKERRQ(ierr);
2665e8a0caSHong Zhang   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
277ba1a0bfSKris Buschelman   PetscFunctionReturn(0);
287ba1a0bfSKris Buschelman }
297ba1a0bfSKris Buschelman 
307ba1a0bfSKris Buschelman #undef __FUNCT__
3153565b12SHong Zhang #define __FUNCT__ "MatDestroy_SeqAIJ_PtAP"
3253565b12SHong Zhang PetscErrorCode MatDestroy_SeqAIJ_PtAP(Mat A)
338d0a38e4SHong Zhang {
348d0a38e4SHong Zhang   PetscErrorCode ierr;
3553565b12SHong Zhang   Mat_SeqAIJ     *a = (Mat_SeqAIJ *)A->data;
3653565b12SHong Zhang   Mat_PtAP       *ptap = a->ptap;
378d0a38e4SHong Zhang 
388d0a38e4SHong Zhang   PetscFunctionBegin;
3953565b12SHong Zhang   /* free ptap, then A */
408d0a38e4SHong Zhang   ierr = PetscFree(ptap->apa);CHKERRQ(ierr);
41f79958ebSHong Zhang   ierr = PetscFree(ptap->api);CHKERRQ(ierr);
42f79958ebSHong Zhang   ierr = PetscFree(ptap->apj);CHKERRQ(ierr);
4353565b12SHong Zhang   ierr = (ptap->destroy)(A);CHKERRQ(ierr);
448d0a38e4SHong Zhang   ierr = PetscFree(ptap);CHKERRQ(ierr);
458d0a38e4SHong Zhang   PetscFunctionReturn(0);
468d0a38e4SHong Zhang }
478d0a38e4SHong Zhang 
488d0a38e4SHong Zhang #undef __FUNCT__
49*ae32dd66SHong Zhang #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy"
50*ae32dd66SHong Zhang PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(Mat A,Mat P,PetscReal fill,Mat *C)
518d0a38e4SHong Zhang {
528d0a38e4SHong Zhang   PetscErrorCode     ierr;
5353565b12SHong Zhang   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
54d20bfe6fSHong Zhang   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data,*p = (Mat_SeqAIJ*)P->data,*c;
5553565b12SHong Zhang   PetscInt           *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj;
5653565b12SHong Zhang   PetscInt           *ci,*cj,*ptadenserow,*ptasparserow,*ptaj,nspacedouble=0;
57d0f46423SBarry Smith   PetscInt           an=A->cmap->N,am=A->rmap->N,pn=P->cmap->N;
5853565b12SHong Zhang   PetscInt           i,j,k,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi,nlnk,*lnk;
59d20bfe6fSHong Zhang   MatScalar          *ca;
6053565b12SHong Zhang   PetscBT            lnkbt;
61eb9c0419SKris Buschelman 
62eb9c0419SKris Buschelman   PetscFunctionBegin;
6353565b12SHong Zhang   /* Get ij structure of P^T */
64eb9c0419SKris Buschelman   ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
65d20bfe6fSHong Zhang   ptJ=ptj;
66eb9c0419SKris Buschelman 
67d20bfe6fSHong Zhang   /* Allocate ci array, arrays for fill computation and */
68d20bfe6fSHong Zhang   /* free space for accumulating nonzero column info */
6955a3bba9SHong Zhang   ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
70d20bfe6fSHong Zhang   ci[0] = 0;
71eb9c0419SKris Buschelman 
7255a3bba9SHong Zhang   ierr = PetscMalloc((2*an+1)*sizeof(PetscInt),&ptadenserow);CHKERRQ(ierr);
7355a3bba9SHong Zhang   ierr = PetscMemzero(ptadenserow,(2*an+1)*sizeof(PetscInt));CHKERRQ(ierr);
74d20bfe6fSHong Zhang   ptasparserow = ptadenserow  + an;
7555a3bba9SHong Zhang 
7655a3bba9SHong Zhang   /* create and initialize a linked list */
7755a3bba9SHong Zhang   nlnk = pn+1;
7855a3bba9SHong Zhang   ierr = PetscLLCreate(pn,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr);
79eb9c0419SKris Buschelman 
80f2b054eeSHong Zhang   /* Set initial free space to be fill*nnz(A). */
81d20bfe6fSHong Zhang   /* This should be reasonable if sparsity of PtAP is similar to that of A. */
82f2b054eeSHong Zhang   ierr          = PetscFreeSpaceGet((PetscInt)(fill*ai[am]),&free_space);
83d20bfe6fSHong Zhang   current_space = free_space;
84d20bfe6fSHong Zhang 
85d20bfe6fSHong Zhang   /* Determine symbolic info for each row of C: */
86d20bfe6fSHong Zhang   for (i=0;i<pn;i++) {
87d20bfe6fSHong Zhang     ptnzi  = pti[i+1] - pti[i];
88d20bfe6fSHong Zhang     ptanzi = 0;
89d20bfe6fSHong Zhang     /* Determine symbolic row of PtA: */
90d20bfe6fSHong Zhang     for (j=0;j<ptnzi;j++) {
91d20bfe6fSHong Zhang       arow = *ptJ++;
92d20bfe6fSHong Zhang       anzj = ai[arow+1] - ai[arow];
93d20bfe6fSHong Zhang       ajj  = aj + ai[arow];
94d20bfe6fSHong Zhang       for (k=0;k<anzj;k++) {
95d20bfe6fSHong Zhang         if (!ptadenserow[ajj[k]]) {
96d20bfe6fSHong Zhang           ptadenserow[ajj[k]]    = -1;
97d20bfe6fSHong Zhang           ptasparserow[ptanzi++] = ajj[k];
98d20bfe6fSHong Zhang         }
99d20bfe6fSHong Zhang       }
100d20bfe6fSHong Zhang     }
101d20bfe6fSHong Zhang     /* Using symbolic info for row of PtA, determine symbolic info for row of C: */
102d20bfe6fSHong Zhang     ptaj = ptasparserow;
103d20bfe6fSHong Zhang     cnzi   = 0;
104d20bfe6fSHong Zhang     for (j=0;j<ptanzi;j++) {
105d20bfe6fSHong Zhang       prow = *ptaj++;
106d20bfe6fSHong Zhang       pnzj = pi[prow+1] - pi[prow];
107d20bfe6fSHong Zhang       pjj  = pj + pi[prow];
10855a3bba9SHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
109dadf0e6bSHong Zhang       ierr = PetscLLAddSorted(pnzj,pjj,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr);
11055a3bba9SHong Zhang       cnzi += nlnk;
111d20bfe6fSHong Zhang     }
112d20bfe6fSHong Zhang 
113d20bfe6fSHong Zhang     /* If free space is not available, make more free space */
114d20bfe6fSHong Zhang     /* Double the amount of total space in the list */
115d20bfe6fSHong Zhang     if (current_space->local_remaining<cnzi) {
1164238b7adSHong Zhang       ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,&current_space);CHKERRQ(ierr);
117f2b054eeSHong Zhang       nspacedouble++;
118d20bfe6fSHong Zhang     }
119d20bfe6fSHong Zhang 
120d20bfe6fSHong Zhang     /* Copy data into free space, and zero out denserows */
12155a3bba9SHong Zhang     ierr = PetscLLClean(pn,pn,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
122d20bfe6fSHong Zhang     current_space->array           += cnzi;
123d20bfe6fSHong Zhang     current_space->local_used      += cnzi;
124d20bfe6fSHong Zhang     current_space->local_remaining -= cnzi;
125d20bfe6fSHong Zhang 
126d20bfe6fSHong Zhang     for (j=0;j<ptanzi;j++) {
127d20bfe6fSHong Zhang       ptadenserow[ptasparserow[j]] = 0;
128d20bfe6fSHong Zhang     }
129d20bfe6fSHong Zhang     /* Aside: Perhaps we should save the pta info for the numerical factorization. */
130d20bfe6fSHong Zhang     /*        For now, we will recompute what is needed. */
131d20bfe6fSHong Zhang     ci[i+1] = ci[i] + cnzi;
132d20bfe6fSHong Zhang   }
133d20bfe6fSHong Zhang   /* nnz is now stored in ci[ptm], column indices are in the list of free space */
134d20bfe6fSHong Zhang   /* Allocate space for cj, initialize cj, and */
135d20bfe6fSHong Zhang   /* destroy list of free space and other temporary array(s) */
13655a3bba9SHong Zhang   ierr = PetscMalloc((ci[pn]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
137a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
138d20bfe6fSHong Zhang   ierr = PetscFree(ptadenserow);CHKERRQ(ierr);
13955a3bba9SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
140d20bfe6fSHong Zhang 
14153565b12SHong Zhang   /* Allocate space for ca */
14253565b12SHong Zhang   ierr = PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
14353565b12SHong Zhang   ierr = PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));CHKERRQ(ierr);
14453565b12SHong Zhang 
14553565b12SHong Zhang   /* put together the new matrix */
14653565b12SHong Zhang   ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,pn,pn,ci,cj,ca,C);CHKERRQ(ierr);
147a2f3521dSMark F. Adams   (*C)->rmap->bs = P->cmap->bs;
148a2f3521dSMark F. Adams   (*C)->cmap->bs = P->cmap->bs;
149*ae32dd66SHong Zhang 
15053565b12SHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
15153565b12SHong Zhang   /* Since these are PETSc arrays, change flags to free them as necessary. */
15253565b12SHong Zhang   c = (Mat_SeqAIJ *)((*C)->data);
15353565b12SHong Zhang   c->free_a  = PETSC_TRUE;
15453565b12SHong Zhang   c->free_ij = PETSC_TRUE;
15553565b12SHong Zhang   c->nonew   = 0;
156*ae32dd66SHong Zhang   (*C)->ops->ptapnumeric = MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy; /* should use *C->ops until PtAP insterface is updated to double dispatch as MatMatMult() */
15753565b12SHong Zhang 
15853565b12SHong Zhang   /* Clean up. */
15953565b12SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
16053565b12SHong Zhang #if defined(PETSC_USE_INFO)
16153565b12SHong Zhang   if (ci[pn] != 0) {
16253565b12SHong Zhang     PetscReal afill = ((PetscReal)ci[pn])/ai[am];
16353565b12SHong Zhang     if (afill < 1.0) afill = 1.0;
16453565b12SHong Zhang     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
16553565b12SHong Zhang     ierr = PetscInfo1((*C),"Use MatPtAP(A,P,MatReuse,%G,&C) for best performance.\n",afill);CHKERRQ(ierr);
16653565b12SHong Zhang   } else {
16753565b12SHong Zhang     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
16853565b12SHong Zhang   }
169f79958ebSHong Zhang #endif
17053565b12SHong Zhang   PetscFunctionReturn(0);
17153565b12SHong Zhang }
17253565b12SHong Zhang 
17353565b12SHong Zhang #undef __FUNCT__
174*ae32dd66SHong Zhang #define __FUNCT__ "MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy"
175*ae32dd66SHong Zhang PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat A,Mat P,Mat C)
17653565b12SHong Zhang {
17753565b12SHong Zhang   PetscErrorCode ierr;
17853565b12SHong Zhang   Mat_SeqAIJ     *a  = (Mat_SeqAIJ *) A->data;
17953565b12SHong Zhang   Mat_SeqAIJ     *p  = (Mat_SeqAIJ *) P->data;
18053565b12SHong Zhang   Mat_SeqAIJ     *c  = (Mat_SeqAIJ *) C->data;
18153565b12SHong Zhang   PetscInt       *ai=a->i,*aj=a->j,*apj,*apjdense,*pi=p->i,*pj=p->j,*pJ=p->j,*pjj;
18253565b12SHong Zhang   PetscInt       *ci=c->i,*cj=c->j,*cjj;
18353565b12SHong Zhang   PetscInt       am=A->rmap->N,cn=C->cmap->N,cm=C->rmap->N;
18453565b12SHong Zhang   PetscInt       i,j,k,anzi,pnzi,apnzj,nextap,pnzj,prow,crow;
18553565b12SHong Zhang   MatScalar      *aa=a->a,*apa,*pa=p->a,*pA=p->a,*paj,*ca=c->a,*caj;
18653565b12SHong Zhang 
18753565b12SHong Zhang   PetscFunctionBegin;
188*ae32dd66SHong Zhang   /* Allocate temporary array for storage of one row of A*P (cn: non-scalable) */
189db1ebd87SHong Zhang   ierr = PetscMalloc(cn*(sizeof(MatScalar)+sizeof(PetscInt))+c->rmax*sizeof(PetscInt),&apa);CHKERRQ(ierr);
19053565b12SHong Zhang 
191db1ebd87SHong Zhang   apjdense = (PetscInt *)(apa + cn);
192db1ebd87SHong Zhang   apj      = apjdense + cn;
193db1ebd87SHong Zhang   ierr = PetscMemzero(apa,cn*(sizeof(MatScalar)+sizeof(PetscInt)));CHKERRQ(ierr);
19453565b12SHong Zhang 
19553565b12SHong Zhang   /* Clear old values in C */
19653565b12SHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
19753565b12SHong Zhang 
19853565b12SHong Zhang   for (i=0; i<am; i++) {
19953565b12SHong Zhang     /* Form sparse row of A*P */
20053565b12SHong Zhang     anzi  = ai[i+1] - ai[i];
20153565b12SHong Zhang     apnzj = 0;
20253565b12SHong Zhang     for (j=0; j<anzi; j++) {
20353565b12SHong Zhang       prow = *aj++;
20453565b12SHong Zhang       pnzj = pi[prow+1] - pi[prow];
20553565b12SHong Zhang       pjj  = pj + pi[prow];
20653565b12SHong Zhang       paj  = pa + pi[prow];
20753565b12SHong Zhang       for (k=0;k<pnzj;k++) {
20853565b12SHong Zhang         if (!apjdense[pjj[k]]) {
20953565b12SHong Zhang           apjdense[pjj[k]] = -1;
21053565b12SHong Zhang           apj[apnzj++]     = pjj[k];
21153565b12SHong Zhang         }
21253565b12SHong Zhang         apa[pjj[k]] += (*aa)*paj[k];
21353565b12SHong Zhang       }
21453565b12SHong Zhang       ierr = PetscLogFlops(2.0*pnzj);CHKERRQ(ierr);
21553565b12SHong Zhang       aa++;
21653565b12SHong Zhang     }
21753565b12SHong Zhang 
21853565b12SHong Zhang     /* Sort the j index array for quick sparse axpy. */
21953565b12SHong Zhang     /* Note: a array does not need sorting as it is in dense storage locations. */
22053565b12SHong Zhang     ierr = PetscSortInt(apnzj,apj);CHKERRQ(ierr);
22153565b12SHong Zhang 
22253565b12SHong Zhang     /* Compute P^T*A*P using outer product (P^T)[:,j]*(A*P)[j,:]. */
22353565b12SHong Zhang     pnzi = pi[i+1] - pi[i];
22453565b12SHong Zhang     for (j=0; j<pnzi; j++) {
22553565b12SHong Zhang       nextap = 0;
22653565b12SHong Zhang       crow   = *pJ++;
22753565b12SHong Zhang       cjj    = cj + ci[crow];
22853565b12SHong Zhang       caj    = ca + ci[crow];
22953565b12SHong Zhang       /* Perform sparse axpy operation.  Note cjj includes apj. */
23053565b12SHong Zhang       for (k=0;nextap<apnzj;k++) {
23153565b12SHong Zhang #if defined(PETSC_USE_DEBUG)
23253565b12SHong Zhang         if (k >= ci[crow+1] - ci[crow]) {
23353565b12SHong Zhang           SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"k too large k %d, crow %d",k,crow);
23453565b12SHong Zhang         }
23553565b12SHong Zhang #endif
23653565b12SHong Zhang         if (cjj[k]==apj[nextap]) {
23753565b12SHong Zhang           caj[k] += (*pA)*apa[apj[nextap++]];
23853565b12SHong Zhang         }
23953565b12SHong Zhang       }
24053565b12SHong Zhang       ierr = PetscLogFlops(2.0*apnzj);CHKERRQ(ierr);
24153565b12SHong Zhang       pA++;
24253565b12SHong Zhang     }
24353565b12SHong Zhang 
24453565b12SHong Zhang     /* Zero the current row info for A*P */
24553565b12SHong Zhang     for (j=0; j<apnzj; j++) {
24653565b12SHong Zhang       apa[apj[j]]      = 0.;
24753565b12SHong Zhang       apjdense[apj[j]] = 0;
24853565b12SHong Zhang     }
24953565b12SHong Zhang   }
25053565b12SHong Zhang 
25153565b12SHong Zhang   /* Assemble the final matrix and clean up */
25253565b12SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
25353565b12SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
254*ae32dd66SHong Zhang 
25553565b12SHong Zhang   ierr = PetscFree(apa);CHKERRQ(ierr);
25653565b12SHong Zhang   PetscFunctionReturn(0);
25753565b12SHong Zhang }
25853565b12SHong Zhang 
25953565b12SHong Zhang #undef __FUNCT__
26053565b12SHong Zhang #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqAIJ"
26153565b12SHong Zhang PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
26253565b12SHong Zhang {
26353565b12SHong Zhang   PetscErrorCode     ierr;
264afff960fSHong Zhang   Mat_SeqAIJ         *ap,*c;
265*ae32dd66SHong Zhang   PetscInt           *api,*apj,*ci,pn=P->cmap->N;
26653565b12SHong Zhang   MatScalar          *ca;
26753565b12SHong Zhang   Mat_PtAP           *ptap;
268971236abSHong Zhang   Mat                Pt,AP;
269*ae32dd66SHong Zhang   PetscBool          sparse_axpy=PETSC_TRUE;
270b8f276daSHong Zhang 
27153565b12SHong Zhang   PetscFunctionBegin;
272*ae32dd66SHong Zhang   ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
27353565b12SHong Zhang   /* flag 'sparse_axpy' determines which implementations to be used:
274*ae32dd66SHong Zhang        0: do dense axpy in MatPtAPNumeric() - fastest, but requires storage of struct A*P;
275*ae32dd66SHong Zhang        1: do two sparse axpy in MatPtAPNumeric() - slowest, does not store structure of A*P. */
276*ae32dd66SHong Zhang   ierr = PetscOptionsBool("-matptap_scalable","Use sparse axpy but slower MatPtAPNumeric()","",sparse_axpy,&sparse_axpy,PETSC_NULL);CHKERRQ(ierr);
277*ae32dd66SHong Zhang   ierr = PetscOptionsEnd();CHKERRQ(ierr);
278*ae32dd66SHong Zhang   if (sparse_axpy){
279*ae32dd66SHong Zhang     ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(A,P,fill,C);CHKERRQ(ierr);
28053565b12SHong Zhang     PetscFunctionReturn(0);
28153565b12SHong Zhang   }
28253565b12SHong Zhang 
283c0e3927dSHong Zhang   /* Get symbolic Pt = P^T */
284c0e3927dSHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(P,&Pt);CHKERRQ(ierr);
28553565b12SHong Zhang 
286c0e3927dSHong Zhang   /* Get symbolic AP = A*P */
287c0e3927dSHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,P,fill,&AP);CHKERRQ(ierr);
288a2f3521dSMark F. Adams 
289c0e3927dSHong Zhang   ap          = (Mat_SeqAIJ*)AP->data;
290c0e3927dSHong Zhang   api         = ap->i;
291c0e3927dSHong Zhang   apj         = ap->j;
292c0e3927dSHong Zhang   ap->free_ij = PETSC_FALSE; /* api and apj are kept in struct ptap, cannot be destroyed with AP */
29353565b12SHong Zhang 
294c0e3927dSHong Zhang   /* Get C = Pt*AP */
295c0e3927dSHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(Pt,AP,fill,C);CHKERRQ(ierr);
296a2f3521dSMark F. Adams 
297c0e3927dSHong Zhang   c  = (Mat_SeqAIJ*)(*C)->data;
298c0e3927dSHong Zhang   ci = c->i;
299d20bfe6fSHong Zhang   ierr = PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
300d20bfe6fSHong Zhang   ierr = PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));CHKERRQ(ierr);
301c0e3927dSHong Zhang   c->a       = ca;
302e6b907acSBarry Smith   c->free_a  = PETSC_TRUE;
303d20bfe6fSHong Zhang 
304c0e3927dSHong Zhang   /* Create a supporting struct for reuse by MatPtAPNumeric() */
3058d0a38e4SHong Zhang   ierr = PetscNew(Mat_PtAP,&ptap);CHKERRQ(ierr);
30653565b12SHong Zhang   c->ptap            = ptap;
3078d0a38e4SHong Zhang   ptap->destroy      = (*C)->ops->destroy;
3088d0a38e4SHong Zhang   (*C)->ops->destroy = MatDestroy_SeqAIJ_PtAP;
3098d0a38e4SHong Zhang 
3108d0a38e4SHong Zhang   /* Allocate temporary array for storage of one row of A*P */
311db1ebd87SHong Zhang   ierr = PetscMalloc((pn+1)*sizeof(PetscScalar),&ptap->apa);CHKERRQ(ierr);
312db1ebd87SHong Zhang   ierr = PetscMemzero(ptap->apa,(pn+1)*sizeof(PetscScalar));CHKERRQ(ierr);
313*ae32dd66SHong Zhang   (*C)->ops->ptapnumeric = MatPtAPNumeric_SeqAIJ_SeqAIJ;
314*ae32dd66SHong Zhang 
315f79958ebSHong Zhang   ptap->api = api;
316f79958ebSHong Zhang   ptap->apj = apj;
3178d0a38e4SHong Zhang 
318d20bfe6fSHong Zhang   /* Clean up. */
319971236abSHong Zhang   ierr = MatDestroy(&Pt);CHKERRQ(ierr);
320971236abSHong Zhang   ierr = MatDestroy(&AP);CHKERRQ(ierr);
321*ae32dd66SHong Zhang #if defined(PETSC_USE_INFO)
322*ae32dd66SHong Zhang   ierr = PetscInfo2((*C),"given fill %G, use scalable %d\n",fill,sparse_axpy);CHKERRQ(ierr);
323*ae32dd66SHong Zhang #endif
324eb9c0419SKris Buschelman   PetscFunctionReturn(0);
325eb9c0419SKris Buschelman }
326eb9c0419SKris Buschelman 
327f2535c29SHong Zhang /* #define PROFILE_MatPtAPNumeric */
328eb9c0419SKris Buschelman #undef __FUNCT__
3299af31e4aSHong Zhang #define __FUNCT__ "MatPtAPNumeric_SeqAIJ_SeqAIJ"
330dfbe8321SBarry Smith PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat A,Mat P,Mat C)
331dfbe8321SBarry Smith {
332dfbe8321SBarry Smith   PetscErrorCode    ierr;
333d20bfe6fSHong Zhang   Mat_SeqAIJ        *a  = (Mat_SeqAIJ *) A->data;
334f79958ebSHong Zhang   Mat_SeqAIJ        *p  = (Mat_SeqAIJ *) P->data;
335d20bfe6fSHong Zhang   Mat_SeqAIJ        *c  = (Mat_SeqAIJ *) C->data;
336b8f276daSHong Zhang   const PetscInt    *ai=a->i,*aj=a->j,*pi=p->i,*pj=p->j,*ci=c->i,*cj=c->j;
337b8f276daSHong Zhang   const PetscScalar *aa=a->a,*pa=p->a,*pval;
338b8f276daSHong Zhang   const PetscInt    *apj,*pcol,*cjj;
339b8f276daSHong Zhang   const PetscInt    am=A->rmap->N,cm=C->rmap->N;
340b8f276daSHong Zhang   PetscInt          i,j,k,anz,apnz,pnz,prow,crow,cnz;
341b8f276daSHong Zhang   PetscScalar       *apa,*ca=c->a,*caj,pvalj;
34253565b12SHong Zhang   Mat_PtAP          *ptap = c->ptap;
343f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric)
344f2535c29SHong Zhang   PetscLogDouble    t0,tf,time_Cseq0=0.0,time_Cseq1=0.0;
345f2535c29SHong Zhang   PetscInt          flops0=0,flops1=0;
346f2535c29SHong Zhang #endif
347eb9c0419SKris Buschelman 
348eb9c0419SKris Buschelman   PetscFunctionBegin;
3498d0a38e4SHong Zhang   /* Get temporary array for storage of one row of A*P */
3508d0a38e4SHong Zhang   apa = ptap->apa;
351d20bfe6fSHong Zhang 
352d20bfe6fSHong Zhang   /* Clear old values in C */
353d20bfe6fSHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
354d20bfe6fSHong Zhang 
355d20bfe6fSHong Zhang   for (i=0;i<am;i++) {
35678844af4SHong Zhang     /* Form sparse row of AP[i,:] = A[i,:]*P */
357f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric)
358f2535c29SHong Zhang     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
359f2535c29SHong Zhang #endif
36078844af4SHong Zhang     anz  = ai[i+1] - ai[i];
36178844af4SHong Zhang     apnz = 0;
36278844af4SHong Zhang     for (j=0; j<anz; j++) {
36378844af4SHong Zhang       prow = aj[j];
36478844af4SHong Zhang       pnz  = pi[prow+1] - pi[prow];
36578844af4SHong Zhang       pcol = pj + pi[prow];
36678844af4SHong Zhang       pval = pa + pi[prow];
36778844af4SHong Zhang       for (k=0; k<pnz; k++) {
36878844af4SHong Zhang         apa[pcol[k]] += aa[j]*pval[k];
369d20bfe6fSHong Zhang       }
37078844af4SHong Zhang       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
371f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric)
372f2535c29SHong Zhang       flops0 += 2.0*pnz;
373f2535c29SHong Zhang #endif
37478844af4SHong Zhang     }
37578844af4SHong Zhang     aj += anz; aa += anz;
376f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric)
377f2535c29SHong Zhang     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
378f2535c29SHong Zhang     time_Cseq0 += tf - t0;
379f2535c29SHong Zhang #endif
38078844af4SHong Zhang 
38178844af4SHong Zhang     /* Compute P^T*A*P using outer product P[i,:]^T*AP[i,:]. */
382f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric)
383f2535c29SHong Zhang     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
384f2535c29SHong Zhang #endif
385f79958ebSHong Zhang     apj  = ptap->apj + ptap->api[i];
386f79958ebSHong Zhang     apnz = ptap->api[i+1] - ptap->api[i];
38778844af4SHong Zhang     pnz  = pi[i+1] - pi[i];
38878844af4SHong Zhang     pcol = pj + pi[i];
38978844af4SHong Zhang     pval = pa + pi[i];
3908d0a38e4SHong Zhang 
39153565b12SHong Zhang     /* Perform dense axpy */
39253565b12SHong Zhang     for (j=0; j<pnz; j++) {
39353565b12SHong Zhang       crow  = pcol[j];
39453565b12SHong Zhang       cjj   = cj + ci[crow];
39553565b12SHong Zhang       caj   = ca + ci[crow];
396f2535c29SHong Zhang       pvalj = pval[j];
39753565b12SHong Zhang       cnz   = ci[crow+1] - ci[crow];
39853565b12SHong Zhang       for (k=0; k<cnz; k++){
399f2535c29SHong Zhang         caj[k] += pvalj*apa[cjj[k]];
40053565b12SHong Zhang       }
40153565b12SHong Zhang       ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
402f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric)
403f2535c29SHong Zhang       flops1 += 2.0*cnz;
404f2535c29SHong Zhang #endif
40553565b12SHong Zhang     }
406f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric)
407f2535c29SHong Zhang     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
408f2535c29SHong Zhang     time_Cseq1 += tf - t0;
409f2535c29SHong Zhang #endif
41053565b12SHong Zhang 
41153565b12SHong Zhang     /* Zero the current row info for A*P */
4121d633602SHong Zhang     for (j=0; j<apnz; j++) apa[apj[j]] = 0.0;
41353565b12SHong Zhang   }
41453565b12SHong Zhang 
41553565b12SHong Zhang   /* Assemble the final matrix and clean up */
41653565b12SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
41753565b12SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
418f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric)
419f2535c29SHong Zhang   printf("PtAPNumeric_SeqAIJ time %g + %g, flops %d %d\n",time_Cseq0,time_Cseq1,flops0,flops1);
420f2535c29SHong Zhang #endif
42153565b12SHong Zhang   PetscFunctionReturn(0);
42253565b12SHong Zhang }
423