xref: /petsc/src/mat/impls/aij/seq/matptap.c (revision 68eacb73b84ae7f3fd7363217d47f23a8f967155)
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)
204a1b09b7SHong Zhang   const char     *algTypes[2] = {"scalable","nonscalable"};
216f231fbdSstefano_zampini   PetscInt       nalg = 2;
226f231fbdSstefano_zampini #else
236f231fbdSstefano_zampini   const char     *algTypes[3] = {"scalable","nonscalable","hypre"};
246f231fbdSstefano_zampini   PetscInt       nalg = 3;
256f231fbdSstefano_zampini #endif
264a1b09b7SHong Zhang   PetscInt       alg = 0; /* set default algorithm */
277ba1a0bfSKris Buschelman 
287ba1a0bfSKris Buschelman   PetscFunctionBegin;
2965e8a0caSHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
304a1b09b7SHong Zhang     /*
314a1b09b7SHong Zhang      Alg 'scalable' determines which implementations to be used:
324a1b09b7SHong Zhang        "nonscalable": do dense axpy in MatPtAPNumeric() - fastest, but requires storage of struct A*P;
334a1b09b7SHong Zhang        "scalable":    do two sparse axpy in MatPtAPNumeric() - might slow, does not store structure of A*P.
346f231fbdSstefano_zampini        "hypre":    use boomerAMGBuildCoarseOperator.
354a1b09b7SHong Zhang      */
364a1b09b7SHong Zhang     ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
37*68eacb73SHong Zhang     PetscOptionsObject->alreadyprinted = PETSC_FALSE; /* a hack to ensure the option shows in '-help' */
386f231fbdSstefano_zampini     ierr = PetscOptionsEList("-matptap_via","Algorithmic approach","MatPtAP",algTypes,nalg,algTypes[0],&alg,NULL);CHKERRQ(ierr);
394a1b09b7SHong Zhang     ierr = PetscOptionsEnd();CHKERRQ(ierr);
403ff4c91cSHong Zhang     ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
414a1b09b7SHong Zhang     switch (alg) {
424a1b09b7SHong Zhang     case 1:
434a1b09b7SHong Zhang       ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ_DenseAxpy(A,P,fill,C);CHKERRQ(ierr);
444a1b09b7SHong Zhang       break;
456f231fbdSstefano_zampini #if defined(PETSC_HAVE_HYPRE)
466f231fbdSstefano_zampini     case 2:
476f231fbdSstefano_zampini       ierr = MatPtAPSymbolic_AIJ_AIJ_wHYPRE(A,P,fill,C);CHKERRQ(ierr);
486f231fbdSstefano_zampini       break;
496f231fbdSstefano_zampini #endif
504a1b09b7SHong Zhang     default:
514a1b09b7SHong Zhang       ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(A,P,fill,C);CHKERRQ(ierr);
524a1b09b7SHong Zhang       break;
534a1b09b7SHong Zhang     }
543ff4c91cSHong Zhang     ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
557ba1a0bfSKris Buschelman   }
563ff4c91cSHong Zhang   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
57ae32dd66SHong Zhang   ierr = (*(*C)->ops->ptapnumeric)(A,P,*C);CHKERRQ(ierr);
583ff4c91cSHong Zhang   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
597ba1a0bfSKris Buschelman   PetscFunctionReturn(0);
607ba1a0bfSKris Buschelman }
617ba1a0bfSKris Buschelman 
6253565b12SHong Zhang PetscErrorCode MatDestroy_SeqAIJ_PtAP(Mat A)
638d0a38e4SHong Zhang {
648d0a38e4SHong Zhang   PetscErrorCode ierr;
6553565b12SHong Zhang   Mat_SeqAIJ     *a    = (Mat_SeqAIJ*)A->data;
6653565b12SHong Zhang   Mat_PtAP       *ptap = a->ptap;
678d0a38e4SHong Zhang 
688d0a38e4SHong Zhang   PetscFunctionBegin;
698d0a38e4SHong Zhang   ierr = PetscFree(ptap->apa);CHKERRQ(ierr);
70f79958ebSHong Zhang   ierr = PetscFree(ptap->api);CHKERRQ(ierr);
71f79958ebSHong Zhang   ierr = PetscFree(ptap->apj);CHKERRQ(ierr);
7253565b12SHong Zhang   ierr = (ptap->destroy)(A);CHKERRQ(ierr);
738d0a38e4SHong Zhang   ierr = PetscFree(ptap);CHKERRQ(ierr);
748d0a38e4SHong Zhang   PetscFunctionReturn(0);
758d0a38e4SHong Zhang }
768d0a38e4SHong Zhang 
77ae32dd66SHong Zhang PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(Mat A,Mat P,PetscReal fill,Mat *C)
788d0a38e4SHong Zhang {
798d0a38e4SHong Zhang   PetscErrorCode     ierr;
800298fd71SBarry Smith   PetscFreeSpaceList free_space=NULL,current_space=NULL;
81d20bfe6fSHong Zhang   Mat_SeqAIJ         *a        = (Mat_SeqAIJ*)A->data,*p = (Mat_SeqAIJ*)P->data,*c;
8253565b12SHong Zhang   PetscInt           *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj;
8353565b12SHong Zhang   PetscInt           *ci,*cj,*ptadenserow,*ptasparserow,*ptaj,nspacedouble=0;
847d69fa63SHong Zhang   PetscInt           an=A->cmap->N,am=A->rmap->N,pn=P->cmap->N,pm=P->rmap->N;
8553565b12SHong Zhang   PetscInt           i,j,k,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi,nlnk,*lnk;
86d20bfe6fSHong Zhang   MatScalar          *ca;
8753565b12SHong Zhang   PetscBT            lnkbt;
887d69fa63SHong Zhang   PetscReal          afill;
89eb9c0419SKris Buschelman 
90eb9c0419SKris Buschelman   PetscFunctionBegin;
9153565b12SHong Zhang   /* Get ij structure of P^T */
92eb9c0419SKris Buschelman   ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
93d20bfe6fSHong Zhang   ptJ  = ptj;
94eb9c0419SKris Buschelman 
95d20bfe6fSHong Zhang   /* Allocate ci array, arrays for fill computation and */
96d20bfe6fSHong Zhang   /* free space for accumulating nonzero column info */
97854ce69bSBarry Smith   ierr  = PetscMalloc1(pn+1,&ci);CHKERRQ(ierr);
98d20bfe6fSHong Zhang   ci[0] = 0;
99eb9c0419SKris Buschelman 
1001795a4d1SJed Brown   ierr         = PetscCalloc1(2*an+1,&ptadenserow);CHKERRQ(ierr);
101d20bfe6fSHong Zhang   ptasparserow = ptadenserow  + an;
10255a3bba9SHong Zhang 
10355a3bba9SHong Zhang   /* create and initialize a linked list */
10455a3bba9SHong Zhang   nlnk = pn+1;
10555a3bba9SHong Zhang   ierr = PetscLLCreate(pn,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr);
106eb9c0419SKris Buschelman 
1077d69fa63SHong Zhang   /* Set initial free space to be fill*(nnz(A)+ nnz(P)) */
108f91af8c7SBarry Smith   ierr          = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],pi[pm])),&free_space);CHKERRQ(ierr);
109d20bfe6fSHong Zhang   current_space = free_space;
110d20bfe6fSHong Zhang 
111d20bfe6fSHong Zhang   /* Determine symbolic info for each row of C: */
112d20bfe6fSHong Zhang   for (i=0; i<pn; i++) {
113d20bfe6fSHong Zhang     ptnzi  = pti[i+1] - pti[i];
114d20bfe6fSHong Zhang     ptanzi = 0;
115d20bfe6fSHong Zhang     /* Determine symbolic row of PtA: */
116d20bfe6fSHong Zhang     for (j=0; j<ptnzi; j++) {
117d20bfe6fSHong Zhang       arow = *ptJ++;
118d20bfe6fSHong Zhang       anzj = ai[arow+1] - ai[arow];
119d20bfe6fSHong Zhang       ajj  = aj + ai[arow];
120d20bfe6fSHong Zhang       for (k=0; k<anzj; k++) {
121d20bfe6fSHong Zhang         if (!ptadenserow[ajj[k]]) {
122d20bfe6fSHong Zhang           ptadenserow[ajj[k]]    = -1;
123d20bfe6fSHong Zhang           ptasparserow[ptanzi++] = ajj[k];
124d20bfe6fSHong Zhang         }
125d20bfe6fSHong Zhang       }
126d20bfe6fSHong Zhang     }
127d20bfe6fSHong Zhang     /* Using symbolic info for row of PtA, determine symbolic info for row of C: */
128d20bfe6fSHong Zhang     ptaj = ptasparserow;
129d20bfe6fSHong Zhang     cnzi = 0;
130d20bfe6fSHong Zhang     for (j=0; j<ptanzi; j++) {
131d20bfe6fSHong Zhang       prow = *ptaj++;
132d20bfe6fSHong Zhang       pnzj = pi[prow+1] - pi[prow];
133d20bfe6fSHong Zhang       pjj  = pj + pi[prow];
13455a3bba9SHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
135dadf0e6bSHong Zhang       ierr  = PetscLLAddSorted(pnzj,pjj,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr);
13655a3bba9SHong Zhang       cnzi += nlnk;
137d20bfe6fSHong Zhang     }
138d20bfe6fSHong Zhang 
139d20bfe6fSHong Zhang     /* If free space is not available, make more free space */
140d20bfe6fSHong Zhang     /* Double the amount of total space in the list */
141d20bfe6fSHong Zhang     if (current_space->local_remaining<cnzi) {
142f91af8c7SBarry Smith       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(cnzi,current_space->total_array_size),&current_space);CHKERRQ(ierr);
143f2b054eeSHong Zhang       nspacedouble++;
144d20bfe6fSHong Zhang     }
145d20bfe6fSHong Zhang 
146d20bfe6fSHong Zhang     /* Copy data into free space, and zero out denserows */
14755a3bba9SHong Zhang     ierr = PetscLLClean(pn,pn,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
1482205254eSKarl Rupp 
149d20bfe6fSHong Zhang     current_space->array           += cnzi;
150d20bfe6fSHong Zhang     current_space->local_used      += cnzi;
151d20bfe6fSHong Zhang     current_space->local_remaining -= cnzi;
152d20bfe6fSHong Zhang 
1532205254eSKarl Rupp     for (j=0; j<ptanzi; j++) ptadenserow[ptasparserow[j]] = 0;
1542205254eSKarl Rupp 
155d20bfe6fSHong Zhang     /* Aside: Perhaps we should save the pta info for the numerical factorization. */
156d20bfe6fSHong Zhang     /*        For now, we will recompute what is needed. */
157d20bfe6fSHong Zhang     ci[i+1] = ci[i] + cnzi;
158d20bfe6fSHong Zhang   }
159d20bfe6fSHong Zhang   /* nnz is now stored in ci[ptm], column indices are in the list of free space */
160d20bfe6fSHong Zhang   /* Allocate space for cj, initialize cj, and */
161d20bfe6fSHong Zhang   /* destroy list of free space and other temporary array(s) */
162854ce69bSBarry Smith   ierr = PetscMalloc1(ci[pn]+1,&cj);CHKERRQ(ierr);
163a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
164d20bfe6fSHong Zhang   ierr = PetscFree(ptadenserow);CHKERRQ(ierr);
16555a3bba9SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
166d20bfe6fSHong Zhang 
167854ce69bSBarry Smith   ierr = PetscCalloc1(ci[pn]+1,&ca);CHKERRQ(ierr);
16853565b12SHong Zhang 
16953565b12SHong Zhang   /* put together the new matrix */
170ce94432eSBarry Smith   ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),pn,pn,ci,cj,ca,C);CHKERRQ(ierr);
17133d57670SJed Brown   ierr = MatSetBlockSizes(*C,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr);
172ae32dd66SHong Zhang 
17353565b12SHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
17453565b12SHong Zhang   /* Since these are PETSc arrays, change flags to free them as necessary. */
17553565b12SHong Zhang   c          = (Mat_SeqAIJ*)((*C)->data);
17653565b12SHong Zhang   c->free_a  = PETSC_TRUE;
17753565b12SHong Zhang   c->free_ij = PETSC_TRUE;
17853565b12SHong Zhang   c->nonew   = 0;
179f4a1e0b0SHong Zhang   (*C)->ops->ptapnumeric = MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy;
18053565b12SHong Zhang 
1817d69fa63SHong Zhang   /* set MatInfo */
1827d69fa63SHong Zhang   afill = (PetscReal)ci[pn]/(ai[am]+pi[pm] + 1.e-5);
1837d69fa63SHong Zhang   if (afill < 1.0) afill = 1.0;
1847d69fa63SHong Zhang   c->maxnz                     = ci[pn];
1857d69fa63SHong Zhang   c->nz                        = ci[pn];
1867d69fa63SHong Zhang   (*C)->info.mallocs           = nspacedouble;
1877d69fa63SHong Zhang   (*C)->info.fill_ratio_given  = fill;
1887d69fa63SHong Zhang   (*C)->info.fill_ratio_needed = afill;
1897d69fa63SHong Zhang 
19053565b12SHong Zhang   /* Clean up. */
19153565b12SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
19253565b12SHong Zhang #if defined(PETSC_USE_INFO)
19353565b12SHong Zhang   if (ci[pn] != 0) {
19457622a8eSBarry Smith     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
19557622a8eSBarry Smith     ierr = PetscInfo1((*C),"Use MatPtAP(A,P,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr);
19653565b12SHong Zhang   } else {
19753565b12SHong Zhang     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
19853565b12SHong Zhang   }
199f79958ebSHong Zhang #endif
20053565b12SHong Zhang   PetscFunctionReturn(0);
20153565b12SHong Zhang }
20253565b12SHong Zhang 
203ae32dd66SHong Zhang PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat A,Mat P,Mat C)
20453565b12SHong Zhang {
20553565b12SHong Zhang   PetscErrorCode ierr;
20653565b12SHong Zhang   Mat_SeqAIJ     *a = (Mat_SeqAIJ*) A->data;
20753565b12SHong Zhang   Mat_SeqAIJ     *p = (Mat_SeqAIJ*) P->data;
20853565b12SHong Zhang   Mat_SeqAIJ     *c = (Mat_SeqAIJ*) C->data;
20953565b12SHong Zhang   PetscInt       *ai=a->i,*aj=a->j,*apj,*apjdense,*pi=p->i,*pj=p->j,*pJ=p->j,*pjj;
21053565b12SHong Zhang   PetscInt       *ci=c->i,*cj=c->j,*cjj;
21153565b12SHong Zhang   PetscInt       am =A->rmap->N,cn=C->cmap->N,cm=C->rmap->N;
21253565b12SHong Zhang   PetscInt       i,j,k,anzi,pnzi,apnzj,nextap,pnzj,prow,crow;
21353565b12SHong Zhang   MatScalar      *aa=a->a,*apa,*pa=p->a,*pA=p->a,*paj,*ca=c->a,*caj;
21453565b12SHong Zhang 
21553565b12SHong Zhang   PetscFunctionBegin;
216ae32dd66SHong Zhang   /* Allocate temporary array for storage of one row of A*P (cn: non-scalable) */
2174bcca364SHong Zhang   ierr = PetscMalloc3(cn,&apa,cn,&apjdense,cn,&apj);CHKERRQ(ierr);
218eb9baa12SBarry Smith   ierr = PetscMemzero(apa,cn*sizeof(MatScalar));CHKERRQ(ierr);
219eb9baa12SBarry Smith   ierr = PetscMemzero(apjdense,cn*sizeof(PetscInt));CHKERRQ(ierr);
22053565b12SHong Zhang 
22153565b12SHong Zhang   /* Clear old values in C */
22253565b12SHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
22353565b12SHong Zhang 
22453565b12SHong Zhang   for (i=0; i<am; i++) {
22553565b12SHong Zhang     /* Form sparse row of A*P */
22653565b12SHong Zhang     anzi  = ai[i+1] - ai[i];
22753565b12SHong Zhang     apnzj = 0;
22853565b12SHong Zhang     for (j=0; j<anzi; j++) {
22953565b12SHong Zhang       prow = *aj++;
23053565b12SHong Zhang       pnzj = pi[prow+1] - pi[prow];
23153565b12SHong Zhang       pjj  = pj + pi[prow];
23253565b12SHong Zhang       paj  = pa + pi[prow];
23353565b12SHong Zhang       for (k=0; k<pnzj; k++) {
23453565b12SHong Zhang         if (!apjdense[pjj[k]]) {
23553565b12SHong Zhang           apjdense[pjj[k]] = -1;
23653565b12SHong Zhang           apj[apnzj++]     = pjj[k];
23753565b12SHong Zhang         }
23853565b12SHong Zhang         apa[pjj[k]] += (*aa)*paj[k];
23953565b12SHong Zhang       }
24053565b12SHong Zhang       ierr = PetscLogFlops(2.0*pnzj);CHKERRQ(ierr);
24153565b12SHong Zhang       aa++;
24253565b12SHong Zhang     }
24353565b12SHong Zhang 
24453565b12SHong Zhang     /* Sort the j index array for quick sparse axpy. */
24553565b12SHong Zhang     /* Note: a array does not need sorting as it is in dense storage locations. */
24653565b12SHong Zhang     ierr = PetscSortInt(apnzj,apj);CHKERRQ(ierr);
24753565b12SHong Zhang 
24853565b12SHong Zhang     /* Compute P^T*A*P using outer product (P^T)[:,j]*(A*P)[j,:]. */
24953565b12SHong Zhang     pnzi = pi[i+1] - pi[i];
25053565b12SHong Zhang     for (j=0; j<pnzi; j++) {
25153565b12SHong Zhang       nextap = 0;
25253565b12SHong Zhang       crow   = *pJ++;
25353565b12SHong Zhang       cjj    = cj + ci[crow];
25453565b12SHong Zhang       caj    = ca + ci[crow];
25553565b12SHong Zhang       /* Perform sparse axpy operation.  Note cjj includes apj. */
25653565b12SHong Zhang       for (k=0; nextap<apnzj; k++) {
25753565b12SHong Zhang #if defined(PETSC_USE_DEBUG)
258f23aa3ddSBarry Smith         if (k >= ci[crow+1] - ci[crow]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"k too large k %d, crow %d",k,crow);
25953565b12SHong Zhang #endif
26053565b12SHong Zhang         if (cjj[k]==apj[nextap]) {
26153565b12SHong Zhang           caj[k] += (*pA)*apa[apj[nextap++]];
26253565b12SHong Zhang         }
26353565b12SHong Zhang       }
26453565b12SHong Zhang       ierr = PetscLogFlops(2.0*apnzj);CHKERRQ(ierr);
26553565b12SHong Zhang       pA++;
26653565b12SHong Zhang     }
26753565b12SHong Zhang 
26853565b12SHong Zhang     /* Zero the current row info for A*P */
26953565b12SHong Zhang     for (j=0; j<apnzj; j++) {
27053565b12SHong Zhang       apa[apj[j]]      = 0.;
27153565b12SHong Zhang       apjdense[apj[j]] = 0;
27253565b12SHong Zhang     }
27353565b12SHong Zhang   }
27453565b12SHong Zhang 
27553565b12SHong Zhang   /* Assemble the final matrix and clean up */
27653565b12SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
27753565b12SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
278ae32dd66SHong Zhang 
279eb9baa12SBarry Smith   ierr = PetscFree3(apa,apjdense,apj);CHKERRQ(ierr);
28053565b12SHong Zhang   PetscFunctionReturn(0);
28153565b12SHong Zhang }
28253565b12SHong Zhang 
2834a1b09b7SHong Zhang PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_DenseAxpy(Mat A,Mat P,PetscReal fill,Mat *C)
28453565b12SHong Zhang {
28553565b12SHong Zhang   PetscErrorCode ierr;
286afff960fSHong Zhang   Mat_SeqAIJ     *ap,*c;
287ae32dd66SHong Zhang   PetscInt       *api,*apj,*ci,pn=P->cmap->N;
28853565b12SHong Zhang   MatScalar      *ca;
28953565b12SHong Zhang   Mat_PtAP       *ptap;
290971236abSHong Zhang   Mat            Pt,AP;
291b8f276daSHong Zhang 
29253565b12SHong Zhang   PetscFunctionBegin;
293c0e3927dSHong Zhang   /* Get symbolic Pt = P^T */
294c0e3927dSHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(P,&Pt);CHKERRQ(ierr);
29553565b12SHong Zhang 
296c0e3927dSHong Zhang   /* Get symbolic AP = A*P */
297c0e3927dSHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,P,fill,&AP);CHKERRQ(ierr);
298a2f3521dSMark F. Adams 
299c0e3927dSHong Zhang   ap          = (Mat_SeqAIJ*)AP->data;
300c0e3927dSHong Zhang   api         = ap->i;
301c0e3927dSHong Zhang   apj         = ap->j;
302c0e3927dSHong Zhang   ap->free_ij = PETSC_FALSE; /* api and apj are kept in struct ptap, cannot be destroyed with AP */
30353565b12SHong Zhang 
304c0e3927dSHong Zhang   /* Get C = Pt*AP */
305c0e3927dSHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(Pt,AP,fill,C);CHKERRQ(ierr);
306a2f3521dSMark F. Adams 
307c0e3927dSHong Zhang   c         = (Mat_SeqAIJ*)(*C)->data;
308c0e3927dSHong Zhang   ci        = c->i;
3091795a4d1SJed Brown   ierr      = PetscCalloc1(ci[pn]+1,&ca);CHKERRQ(ierr);
310c0e3927dSHong Zhang   c->a      = ca;
311e6b907acSBarry Smith   c->free_a = PETSC_TRUE;
312d20bfe6fSHong Zhang 
313c0e3927dSHong Zhang   /* Create a supporting struct for reuse by MatPtAPNumeric() */
314b00a9115SJed Brown   ierr = PetscNew(&ptap);CHKERRQ(ierr);
3152205254eSKarl Rupp 
31653565b12SHong Zhang   c->ptap            = ptap;
3178d0a38e4SHong Zhang   ptap->destroy      = (*C)->ops->destroy;
3188d0a38e4SHong Zhang   (*C)->ops->destroy = MatDestroy_SeqAIJ_PtAP;
3198d0a38e4SHong Zhang 
3208d0a38e4SHong Zhang   /* Allocate temporary array for storage of one row of A*P */
3211795a4d1SJed Brown   ierr = PetscCalloc1(pn+1,&ptap->apa);CHKERRQ(ierr);
3222205254eSKarl Rupp 
323ae32dd66SHong Zhang   (*C)->ops->ptapnumeric = MatPtAPNumeric_SeqAIJ_SeqAIJ;
324ae32dd66SHong Zhang 
325f79958ebSHong Zhang   ptap->api = api;
326f79958ebSHong Zhang   ptap->apj = apj;
3278d0a38e4SHong Zhang 
328d20bfe6fSHong Zhang   /* Clean up. */
329971236abSHong Zhang   ierr = MatDestroy(&Pt);CHKERRQ(ierr);
330971236abSHong Zhang   ierr = MatDestroy(&AP);CHKERRQ(ierr);
331ae32dd66SHong Zhang #if defined(PETSC_USE_INFO)
33257622a8eSBarry Smith   ierr = PetscInfo1((*C),"given fill %g\n",(double)fill);CHKERRQ(ierr);
333ae32dd66SHong Zhang #endif
334eb9c0419SKris Buschelman   PetscFunctionReturn(0);
335eb9c0419SKris Buschelman }
336eb9c0419SKris Buschelman 
337f2535c29SHong Zhang /* #define PROFILE_MatPtAPNumeric */
338dfbe8321SBarry Smith PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat A,Mat P,Mat C)
339dfbe8321SBarry Smith {
340dfbe8321SBarry Smith   PetscErrorCode    ierr;
341d20bfe6fSHong Zhang   Mat_SeqAIJ        *a = (Mat_SeqAIJ*) A->data;
342f79958ebSHong Zhang   Mat_SeqAIJ        *p = (Mat_SeqAIJ*) P->data;
343d20bfe6fSHong Zhang   Mat_SeqAIJ        *c = (Mat_SeqAIJ*) C->data;
344b8f276daSHong Zhang   const PetscInt    *ai=a->i,*aj=a->j,*pi=p->i,*pj=p->j,*ci=c->i,*cj=c->j;
345b8f276daSHong Zhang   const PetscScalar *aa=a->a,*pa=p->a,*pval;
346b8f276daSHong Zhang   const PetscInt    *apj,*pcol,*cjj;
347b8f276daSHong Zhang   const PetscInt    am=A->rmap->N,cm=C->rmap->N;
348b8f276daSHong Zhang   PetscInt          i,j,k,anz,apnz,pnz,prow,crow,cnz;
349b8f276daSHong Zhang   PetscScalar       *apa,*ca=c->a,*caj,pvalj;
35053565b12SHong Zhang   Mat_PtAP          *ptap = c->ptap;
351f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric)
352f2535c29SHong Zhang   PetscLogDouble t0,tf,time_Cseq0=0.0,time_Cseq1=0.0;
353f2535c29SHong Zhang   PetscInt       flops0=0,flops1=0;
354f2535c29SHong Zhang #endif
355eb9c0419SKris Buschelman 
356eb9c0419SKris Buschelman   PetscFunctionBegin;
3578d0a38e4SHong Zhang   /* Get temporary array for storage of one row of A*P */
3588d0a38e4SHong Zhang   apa = ptap->apa;
359d20bfe6fSHong Zhang 
360d20bfe6fSHong Zhang   /* Clear old values in C */
361d20bfe6fSHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
362d20bfe6fSHong Zhang 
363d20bfe6fSHong Zhang   for (i=0; i<am; i++) {
36478844af4SHong Zhang     /* Form sparse row of AP[i,:] = A[i,:]*P */
365f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric)
3668563dfccSBarry Smith     ierr = PetscTime(&t0);CHKERRQ(ierr);
367f2535c29SHong Zhang #endif
36878844af4SHong Zhang     anz  = ai[i+1] - ai[i];
36978844af4SHong Zhang     for (j=0; j<anz; j++) {
37078844af4SHong Zhang       prow = aj[j];
37178844af4SHong Zhang       pnz  = pi[prow+1] - pi[prow];
37278844af4SHong Zhang       pcol = pj + pi[prow];
37378844af4SHong Zhang       pval = pa + pi[prow];
37478844af4SHong Zhang       for (k=0; k<pnz; k++) {
37578844af4SHong Zhang         apa[pcol[k]] += aa[j]*pval[k];
376d20bfe6fSHong Zhang       }
37778844af4SHong Zhang       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
378f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric)
379f2535c29SHong Zhang       flops0 += 2.0*pnz;
380f2535c29SHong Zhang #endif
38178844af4SHong Zhang     }
38278844af4SHong Zhang     aj += anz; aa += anz;
383f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric)
3848563dfccSBarry Smith     ierr = PetscTime(&tf);CHKERRQ(ierr);
3852205254eSKarl Rupp 
386f2535c29SHong Zhang     time_Cseq0 += tf - t0;
387f2535c29SHong Zhang #endif
38878844af4SHong Zhang 
38978844af4SHong Zhang     /* Compute P^T*A*P using outer product P[i,:]^T*AP[i,:]. */
390f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric)
3918563dfccSBarry Smith     ierr = PetscTime(&t0);CHKERRQ(ierr);
392f2535c29SHong Zhang #endif
393f79958ebSHong Zhang     apj  = ptap->apj + ptap->api[i];
394f79958ebSHong Zhang     apnz = ptap->api[i+1] - ptap->api[i];
39578844af4SHong Zhang     pnz  = pi[i+1] - pi[i];
39678844af4SHong Zhang     pcol = pj + pi[i];
39778844af4SHong Zhang     pval = pa + pi[i];
3988d0a38e4SHong Zhang 
39953565b12SHong Zhang     /* Perform dense axpy */
40053565b12SHong Zhang     for (j=0; j<pnz; j++) {
40153565b12SHong Zhang       crow  = pcol[j];
40253565b12SHong Zhang       cjj   = cj + ci[crow];
40353565b12SHong Zhang       caj   = ca + ci[crow];
404f2535c29SHong Zhang       pvalj = pval[j];
40553565b12SHong Zhang       cnz   = ci[crow+1] - ci[crow];
4062205254eSKarl Rupp       for (k=0; k<cnz; k++) caj[k] += pvalj*apa[cjj[k]];
40753565b12SHong Zhang       ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
408f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric)
409f2535c29SHong Zhang       flops1 += 2.0*cnz;
410f2535c29SHong Zhang #endif
41153565b12SHong Zhang     }
412f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric)
4138563dfccSBarry Smith     ierr        = PetscTime(&tf);CHKERRQ(ierr);
414f2535c29SHong Zhang     time_Cseq1 += tf - t0;
415f2535c29SHong Zhang #endif
41653565b12SHong Zhang 
41753565b12SHong Zhang     /* Zero the current row info for A*P */
4181d633602SHong Zhang     for (j=0; j<apnz; j++) apa[apj[j]] = 0.0;
41953565b12SHong Zhang   }
42053565b12SHong Zhang 
42153565b12SHong Zhang   /* Assemble the final matrix and clean up */
42253565b12SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
42353565b12SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
424f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric)
425f2535c29SHong Zhang   printf("PtAPNumeric_SeqAIJ time %g + %g, flops %d %d\n",time_Cseq0,time_Cseq1,flops0,flops1);
426f2535c29SHong Zhang #endif
42753565b12SHong Zhang   PetscFunctionReturn(0);
42853565b12SHong Zhang }
429