xref: /petsc/src/mat/impls/aij/seq/matptap.c (revision a2f3521de2a9d7869700f4c17e26c23fcfeaa6f6)
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__
127ba1a0bfSKris Buschelman #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ"
137ba1a0bfSKris Buschelman PetscErrorCode MatPtAPSymbolic_SeqAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
147ba1a0bfSKris Buschelman {
157ba1a0bfSKris Buschelman   PetscErrorCode ierr;
167ba1a0bfSKris Buschelman 
177ba1a0bfSKris Buschelman   PetscFunctionBegin;
1865e19b50SBarry Smith   if (!P->ops->ptapsymbolic_seqaij) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_SUP,"Not implemented for A=%s and P=%s",((PetscObject)A)->type_name,((PetscObject)P)->type_name);
197ba1a0bfSKris Buschelman   ierr = (*P->ops->ptapsymbolic_seqaij)(A,P,fill,C);CHKERRQ(ierr);
207ba1a0bfSKris Buschelman   PetscFunctionReturn(0);
217ba1a0bfSKris Buschelman }
227ba1a0bfSKris Buschelman 
237ba1a0bfSKris Buschelman #undef __FUNCT__
247ba1a0bfSKris Buschelman #define __FUNCT__ "MatPtAPNumeric_SeqAIJ"
257ba1a0bfSKris Buschelman PetscErrorCode MatPtAPNumeric_SeqAIJ(Mat A,Mat P,Mat C)
267ba1a0bfSKris Buschelman {
277ba1a0bfSKris Buschelman   PetscErrorCode ierr;
287ba1a0bfSKris Buschelman 
297ba1a0bfSKris Buschelman   PetscFunctionBegin;
3065e19b50SBarry Smith   if (!P->ops->ptapnumeric_seqaij) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_SUP,"Not implemented for A=%s and P=%s",((PetscObject)A)->type_name,((PetscObject)P)->type_name);
317ba1a0bfSKris Buschelman   ierr = (*P->ops->ptapnumeric_seqaij)(A,P,C);CHKERRQ(ierr);
327ba1a0bfSKris Buschelman   PetscFunctionReturn(0);
337ba1a0bfSKris Buschelman }
347ba1a0bfSKris Buschelman 
357ba1a0bfSKris Buschelman #undef __FUNCT__
3653565b12SHong Zhang #define __FUNCT__ "MatDestroy_SeqAIJ_PtAP"
3753565b12SHong Zhang PetscErrorCode MatDestroy_SeqAIJ_PtAP(Mat A)
388d0a38e4SHong Zhang {
398d0a38e4SHong Zhang   PetscErrorCode ierr;
4053565b12SHong Zhang   Mat_SeqAIJ     *a = (Mat_SeqAIJ *)A->data;
4153565b12SHong Zhang   Mat_PtAP       *ptap = a->ptap;
428d0a38e4SHong Zhang 
438d0a38e4SHong Zhang   PetscFunctionBegin;
4453565b12SHong Zhang   /* free ptap, then A */
458d0a38e4SHong Zhang   ierr = PetscFree(ptap->apa);CHKERRQ(ierr);
46f79958ebSHong Zhang   ierr = PetscFree(ptap->api);CHKERRQ(ierr);
47f79958ebSHong Zhang   ierr = PetscFree(ptap->apj);CHKERRQ(ierr);
4853565b12SHong Zhang   ierr = (ptap->destroy)(A);CHKERRQ(ierr);
498d0a38e4SHong Zhang   ierr = PetscFree(ptap);CHKERRQ(ierr);
508d0a38e4SHong Zhang   PetscFunctionReturn(0);
518d0a38e4SHong Zhang }
528d0a38e4SHong Zhang 
538d0a38e4SHong Zhang #undef __FUNCT__
5453565b12SHong Zhang #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy2"
5553565b12SHong Zhang PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy2(Mat A,Mat P,PetscReal fill,Mat *C)
568d0a38e4SHong Zhang {
578d0a38e4SHong Zhang   PetscErrorCode     ierr;
5853565b12SHong Zhang   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
59d20bfe6fSHong Zhang   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data,*p = (Mat_SeqAIJ*)P->data,*c;
6053565b12SHong Zhang   PetscInt           *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj;
6153565b12SHong Zhang   PetscInt           *ci,*cj,*ptadenserow,*ptasparserow,*ptaj,nspacedouble=0;
62d0f46423SBarry Smith   PetscInt           an=A->cmap->N,am=A->rmap->N,pn=P->cmap->N;
6353565b12SHong Zhang   PetscInt           i,j,k,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi,nlnk,*lnk;
64d20bfe6fSHong Zhang   MatScalar          *ca;
6553565b12SHong Zhang   PetscBT            lnkbt;
66eb9c0419SKris Buschelman 
67eb9c0419SKris Buschelman   PetscFunctionBegin;
6853565b12SHong Zhang   /* Get ij structure of P^T */
69eb9c0419SKris Buschelman   ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
70d20bfe6fSHong Zhang   ptJ=ptj;
71eb9c0419SKris Buschelman 
72d20bfe6fSHong Zhang   /* Allocate ci array, arrays for fill computation and */
73d20bfe6fSHong Zhang   /* free space for accumulating nonzero column info */
7455a3bba9SHong Zhang   ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
75d20bfe6fSHong Zhang   ci[0] = 0;
76eb9c0419SKris Buschelman 
7755a3bba9SHong Zhang   ierr = PetscMalloc((2*an+1)*sizeof(PetscInt),&ptadenserow);CHKERRQ(ierr);
7855a3bba9SHong Zhang   ierr = PetscMemzero(ptadenserow,(2*an+1)*sizeof(PetscInt));CHKERRQ(ierr);
79d20bfe6fSHong Zhang   ptasparserow = ptadenserow  + an;
8055a3bba9SHong Zhang 
8155a3bba9SHong Zhang   /* create and initialize a linked list */
8255a3bba9SHong Zhang   nlnk = pn+1;
8355a3bba9SHong Zhang   ierr = PetscLLCreate(pn,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr);
84eb9c0419SKris Buschelman 
85f2b054eeSHong Zhang   /* Set initial free space to be fill*nnz(A). */
86d20bfe6fSHong Zhang   /* This should be reasonable if sparsity of PtAP is similar to that of A. */
87f2b054eeSHong Zhang   ierr          = PetscFreeSpaceGet((PetscInt)(fill*ai[am]),&free_space);
88d20bfe6fSHong Zhang   current_space = free_space;
89d20bfe6fSHong Zhang 
90d20bfe6fSHong Zhang   /* Determine symbolic info for each row of C: */
91d20bfe6fSHong Zhang   for (i=0;i<pn;i++) {
92d20bfe6fSHong Zhang     ptnzi  = pti[i+1] - pti[i];
93d20bfe6fSHong Zhang     ptanzi = 0;
94d20bfe6fSHong Zhang     /* Determine symbolic row of PtA: */
95d20bfe6fSHong Zhang     for (j=0;j<ptnzi;j++) {
96d20bfe6fSHong Zhang       arow = *ptJ++;
97d20bfe6fSHong Zhang       anzj = ai[arow+1] - ai[arow];
98d20bfe6fSHong Zhang       ajj  = aj + ai[arow];
99d20bfe6fSHong Zhang       for (k=0;k<anzj;k++) {
100d20bfe6fSHong Zhang         if (!ptadenserow[ajj[k]]) {
101d20bfe6fSHong Zhang           ptadenserow[ajj[k]]    = -1;
102d20bfe6fSHong Zhang           ptasparserow[ptanzi++] = ajj[k];
103d20bfe6fSHong Zhang         }
104d20bfe6fSHong Zhang       }
105d20bfe6fSHong Zhang     }
106d20bfe6fSHong Zhang     /* Using symbolic info for row of PtA, determine symbolic info for row of C: */
107d20bfe6fSHong Zhang     ptaj = ptasparserow;
108d20bfe6fSHong Zhang     cnzi   = 0;
109d20bfe6fSHong Zhang     for (j=0;j<ptanzi;j++) {
110d20bfe6fSHong Zhang       prow = *ptaj++;
111d20bfe6fSHong Zhang       pnzj = pi[prow+1] - pi[prow];
112d20bfe6fSHong Zhang       pjj  = pj + pi[prow];
11355a3bba9SHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
114dadf0e6bSHong Zhang       ierr = PetscLLAddSorted(pnzj,pjj,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr);
11555a3bba9SHong Zhang       cnzi += nlnk;
116d20bfe6fSHong Zhang     }
117d20bfe6fSHong Zhang 
118d20bfe6fSHong Zhang     /* If free space is not available, make more free space */
119d20bfe6fSHong Zhang     /* Double the amount of total space in the list */
120d20bfe6fSHong Zhang     if (current_space->local_remaining<cnzi) {
1214238b7adSHong Zhang       ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,&current_space);CHKERRQ(ierr);
122f2b054eeSHong Zhang       nspacedouble++;
123d20bfe6fSHong Zhang     }
124d20bfe6fSHong Zhang 
125d20bfe6fSHong Zhang     /* Copy data into free space, and zero out denserows */
12655a3bba9SHong Zhang     ierr = PetscLLClean(pn,pn,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
127d20bfe6fSHong Zhang     current_space->array           += cnzi;
128d20bfe6fSHong Zhang     current_space->local_used      += cnzi;
129d20bfe6fSHong Zhang     current_space->local_remaining -= cnzi;
130d20bfe6fSHong Zhang 
131d20bfe6fSHong Zhang     for (j=0;j<ptanzi;j++) {
132d20bfe6fSHong Zhang       ptadenserow[ptasparserow[j]] = 0;
133d20bfe6fSHong Zhang     }
134d20bfe6fSHong Zhang     /* Aside: Perhaps we should save the pta info for the numerical factorization. */
135d20bfe6fSHong Zhang     /*        For now, we will recompute what is needed. */
136d20bfe6fSHong Zhang     ci[i+1] = ci[i] + cnzi;
137d20bfe6fSHong Zhang   }
138d20bfe6fSHong Zhang   /* nnz is now stored in ci[ptm], column indices are in the list of free space */
139d20bfe6fSHong Zhang   /* Allocate space for cj, initialize cj, and */
140d20bfe6fSHong Zhang   /* destroy list of free space and other temporary array(s) */
14155a3bba9SHong Zhang   ierr = PetscMalloc((ci[pn]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
142a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
143d20bfe6fSHong Zhang   ierr = PetscFree(ptadenserow);CHKERRQ(ierr);
14455a3bba9SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
145d20bfe6fSHong Zhang 
14653565b12SHong Zhang   /* Allocate space for ca */
14753565b12SHong Zhang   ierr = PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
14853565b12SHong Zhang   ierr = PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));CHKERRQ(ierr);
14953565b12SHong Zhang 
15053565b12SHong Zhang   /* put together the new matrix */
15153565b12SHong Zhang   ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,pn,pn,ci,cj,ca,C);CHKERRQ(ierr);
152*a2f3521dSMark F. Adams   (*C)->rmap->bs = P->cmap->bs;
153*a2f3521dSMark F. Adams   (*C)->cmap->bs = P->cmap->bs;
154*a2f3521dSMark F. Adams PetscPrintf(PETSC_COMM_SELF,"************%s C.bs=%d,%d\n",__FUNCT__,(*C)->rmap->bs,(*C)->cmap->bs);
15553565b12SHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
15653565b12SHong Zhang   /* Since these are PETSc arrays, change flags to free them as necessary. */
15753565b12SHong Zhang   c = (Mat_SeqAIJ *)((*C)->data);
15853565b12SHong Zhang   c->free_a  = PETSC_TRUE;
15953565b12SHong Zhang   c->free_ij = PETSC_TRUE;
16053565b12SHong Zhang   c->nonew   = 0;
16153565b12SHong Zhang   A->ops->ptapnumeric = MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy2; /* should use *C->ops until PtAP insterface is updated to double dispatch as MatMatMult() */
16253565b12SHong Zhang 
16353565b12SHong Zhang   /* Clean up. */
16453565b12SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
16553565b12SHong Zhang #if defined(PETSC_USE_INFO)
16653565b12SHong Zhang   if (ci[pn] != 0) {
16753565b12SHong Zhang     PetscReal afill = ((PetscReal)ci[pn])/ai[am];
16853565b12SHong Zhang     if (afill < 1.0) afill = 1.0;
16953565b12SHong Zhang     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
17053565b12SHong Zhang     ierr = PetscInfo1((*C),"Use MatPtAP(A,P,MatReuse,%G,&C) for best performance.\n",afill);CHKERRQ(ierr);
17153565b12SHong Zhang   } else {
17253565b12SHong Zhang     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
17353565b12SHong Zhang   }
174f79958ebSHong Zhang #endif
17553565b12SHong Zhang   PetscFunctionReturn(0);
17653565b12SHong Zhang }
17753565b12SHong Zhang 
17853565b12SHong Zhang #undef __FUNCT__
17953565b12SHong Zhang #define __FUNCT__ "MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy2"
18053565b12SHong Zhang PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy2(Mat A,Mat P,Mat C)
18153565b12SHong Zhang {
18253565b12SHong Zhang   PetscErrorCode ierr;
18353565b12SHong Zhang   Mat_SeqAIJ     *a  = (Mat_SeqAIJ *) A->data;
18453565b12SHong Zhang   Mat_SeqAIJ     *p  = (Mat_SeqAIJ *) P->data;
18553565b12SHong Zhang   Mat_SeqAIJ     *c  = (Mat_SeqAIJ *) C->data;
18653565b12SHong Zhang   PetscInt       *ai=a->i,*aj=a->j,*apj,*apjdense,*pi=p->i,*pj=p->j,*pJ=p->j,*pjj;
18753565b12SHong Zhang   PetscInt       *ci=c->i,*cj=c->j,*cjj;
18853565b12SHong Zhang   PetscInt       am=A->rmap->N,cn=C->cmap->N,cm=C->rmap->N;
18953565b12SHong Zhang   PetscInt       i,j,k,anzi,pnzi,apnzj,nextap,pnzj,prow,crow;
19053565b12SHong Zhang   MatScalar      *aa=a->a,*apa,*pa=p->a,*pA=p->a,*paj,*ca=c->a,*caj;
19153565b12SHong Zhang 
19253565b12SHong Zhang   PetscFunctionBegin;
19353565b12SHong Zhang   /* Allocate temporary array for storage of one row of A*P */
194db1ebd87SHong Zhang   ierr = PetscMalloc(cn*(sizeof(MatScalar)+sizeof(PetscInt))+c->rmax*sizeof(PetscInt),&apa);CHKERRQ(ierr);
19553565b12SHong Zhang 
196db1ebd87SHong Zhang   apjdense = (PetscInt *)(apa + cn);
197db1ebd87SHong Zhang   apj      = apjdense + cn;
198db1ebd87SHong Zhang   ierr = PetscMemzero(apa,cn*(sizeof(MatScalar)+sizeof(PetscInt)));CHKERRQ(ierr);
19953565b12SHong Zhang 
20053565b12SHong Zhang   /* Clear old values in C */
20153565b12SHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
20253565b12SHong Zhang 
20353565b12SHong Zhang   for (i=0; i<am; i++) {
20453565b12SHong Zhang     /* Form sparse row of A*P */
20553565b12SHong Zhang     anzi  = ai[i+1] - ai[i];
20653565b12SHong Zhang     apnzj = 0;
20753565b12SHong Zhang     for (j=0; j<anzi; j++) {
20853565b12SHong Zhang       prow = *aj++;
20953565b12SHong Zhang       pnzj = pi[prow+1] - pi[prow];
21053565b12SHong Zhang       pjj  = pj + pi[prow];
21153565b12SHong Zhang       paj  = pa + pi[prow];
21253565b12SHong Zhang       for (k=0;k<pnzj;k++) {
21353565b12SHong Zhang         if (!apjdense[pjj[k]]) {
21453565b12SHong Zhang           apjdense[pjj[k]] = -1;
21553565b12SHong Zhang           apj[apnzj++]     = pjj[k];
21653565b12SHong Zhang         }
21753565b12SHong Zhang         apa[pjj[k]] += (*aa)*paj[k];
21853565b12SHong Zhang       }
21953565b12SHong Zhang       ierr = PetscLogFlops(2.0*pnzj);CHKERRQ(ierr);
22053565b12SHong Zhang       aa++;
22153565b12SHong Zhang     }
22253565b12SHong Zhang 
22353565b12SHong Zhang     /* Sort the j index array for quick sparse axpy. */
22453565b12SHong Zhang     /* Note: a array does not need sorting as it is in dense storage locations. */
22553565b12SHong Zhang     ierr = PetscSortInt(apnzj,apj);CHKERRQ(ierr);
22653565b12SHong Zhang 
22753565b12SHong Zhang     /* Compute P^T*A*P using outer product (P^T)[:,j]*(A*P)[j,:]. */
22853565b12SHong Zhang     pnzi = pi[i+1] - pi[i];
22953565b12SHong Zhang     for (j=0; j<pnzi; j++) {
23053565b12SHong Zhang       nextap = 0;
23153565b12SHong Zhang       crow   = *pJ++;
23253565b12SHong Zhang       cjj    = cj + ci[crow];
23353565b12SHong Zhang       caj    = ca + ci[crow];
23453565b12SHong Zhang       /* Perform sparse axpy operation.  Note cjj includes apj. */
23553565b12SHong Zhang       for (k=0;nextap<apnzj;k++) {
23653565b12SHong Zhang #if defined(PETSC_USE_DEBUG)
23753565b12SHong Zhang         if (k >= ci[crow+1] - ci[crow]) {
23853565b12SHong Zhang           SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"k too large k %d, crow %d",k,crow);
23953565b12SHong Zhang         }
24053565b12SHong Zhang #endif
24153565b12SHong Zhang         if (cjj[k]==apj[nextap]) {
24253565b12SHong Zhang           caj[k] += (*pA)*apa[apj[nextap++]];
24353565b12SHong Zhang         }
24453565b12SHong Zhang       }
24553565b12SHong Zhang       ierr = PetscLogFlops(2.0*apnzj);CHKERRQ(ierr);
24653565b12SHong Zhang       pA++;
24753565b12SHong Zhang     }
24853565b12SHong Zhang 
24953565b12SHong Zhang     /* Zero the current row info for A*P */
25053565b12SHong Zhang     for (j=0; j<apnzj; j++) {
25153565b12SHong Zhang       apa[apj[j]]      = 0.;
25253565b12SHong Zhang       apjdense[apj[j]] = 0;
25353565b12SHong Zhang     }
25453565b12SHong Zhang   }
25553565b12SHong Zhang 
25653565b12SHong Zhang   /* Assemble the final matrix and clean up */
25753565b12SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
25853565b12SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
259*a2f3521dSMark F. Adams PetscPrintf(PETSC_COMM_SELF,"************%s C.bs=%d,%d\n",__FUNCT__,C->rmap->bs,C->cmap->bs);
26053565b12SHong Zhang   ierr = PetscFree(apa);CHKERRQ(ierr);
26153565b12SHong Zhang   PetscFunctionReturn(0);
26253565b12SHong Zhang }
26353565b12SHong Zhang 
26453565b12SHong Zhang #undef __FUNCT__
26553565b12SHong Zhang #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqAIJ"
26653565b12SHong Zhang PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
26753565b12SHong Zhang {
26853565b12SHong Zhang   PetscErrorCode     ierr;
269afff960fSHong Zhang   Mat_SeqAIJ         *ap,*c;
27081b6fd65SJed Brown   PetscInt           *api,*apj,*ci,pn=P->cmap->N,sparse_axpy=0;
27153565b12SHong Zhang   MatScalar          *ca;
27253565b12SHong Zhang   Mat_PtAP           *ptap;
273971236abSHong Zhang   Mat                Pt,AP;
274b8f276daSHong Zhang 
27553565b12SHong Zhang   PetscFunctionBegin;
27653565b12SHong Zhang   /* flag 'sparse_axpy' determines which implementations to be used:
277db1ebd87SHong Zhang        0: do dense axpy in MatPtAPNumeric() - fastest, but requires storage of struct A*P; (default)
278db1ebd87SHong Zhang        1: do one sparse axpy - uses same memory as sparse_axpy=0 and might execute less flops
279db1ebd87SHong Zhang           (apnz vs. cnz in the outerproduct), slower than case '0' when cnz is not too large than apnz;
280db1ebd87SHong Zhang        2: do two sparse axpy in MatPtAPNumeric() - slowest, does not store structure of A*P. */
28153565b12SHong Zhang   ierr = PetscOptionsGetInt(PETSC_NULL,"-matptap_sparseaxpy",&sparse_axpy,PETSC_NULL);CHKERRQ(ierr);
28253565b12SHong Zhang   if (sparse_axpy == 2){
28353565b12SHong Zhang     ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy2(A,P,fill,C);CHKERRQ(ierr);
284*a2f3521dSMark F. Adams 
28553565b12SHong Zhang     PetscFunctionReturn(0);
28653565b12SHong Zhang   }
28753565b12SHong Zhang 
288c0e3927dSHong Zhang   /* Get symbolic Pt = P^T */
289c0e3927dSHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(P,&Pt);CHKERRQ(ierr);
29053565b12SHong Zhang 
291c0e3927dSHong Zhang   /* Get symbolic AP = A*P */
292c0e3927dSHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,P,fill,&AP);CHKERRQ(ierr);
293*a2f3521dSMark F. Adams 
294c0e3927dSHong Zhang   ap          = (Mat_SeqAIJ*)AP->data;
295c0e3927dSHong Zhang   api         = ap->i;
296c0e3927dSHong Zhang   apj         = ap->j;
297c0e3927dSHong Zhang   ap->free_ij = PETSC_FALSE; /* api and apj are kept in struct ptap, cannot be destroyed with AP */
29853565b12SHong Zhang 
299c0e3927dSHong Zhang   /* Get C = Pt*AP */
300c0e3927dSHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(Pt,AP,fill,C);CHKERRQ(ierr);
301*a2f3521dSMark F. Adams 
302c0e3927dSHong Zhang   c  = (Mat_SeqAIJ*)(*C)->data;
303c0e3927dSHong Zhang   ci = c->i;
304d20bfe6fSHong Zhang   ierr = PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
305d20bfe6fSHong Zhang   ierr = PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));CHKERRQ(ierr);
306c0e3927dSHong Zhang   c->a       = ca;
307e6b907acSBarry Smith   c->free_a  = PETSC_TRUE;
308d20bfe6fSHong Zhang 
309c0e3927dSHong Zhang   /* Create a supporting struct for reuse by MatPtAPNumeric() */
3108d0a38e4SHong Zhang   ierr = PetscNew(Mat_PtAP,&ptap);CHKERRQ(ierr);
31153565b12SHong Zhang   c->ptap            = ptap;
3128d0a38e4SHong Zhang   ptap->destroy      = (*C)->ops->destroy;
3138d0a38e4SHong Zhang   (*C)->ops->destroy = MatDestroy_SeqAIJ_PtAP;
3148d0a38e4SHong Zhang 
3158d0a38e4SHong Zhang   /* Allocate temporary array for storage of one row of A*P */
316db1ebd87SHong Zhang   ierr = PetscMalloc((pn+1)*sizeof(PetscScalar),&ptap->apa);CHKERRQ(ierr);
317db1ebd87SHong Zhang   ierr = PetscMemzero(ptap->apa,(pn+1)*sizeof(PetscScalar));CHKERRQ(ierr);
31853565b12SHong Zhang   if (sparse_axpy == 1){
31953565b12SHong Zhang     A->ops->ptapnumeric = MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy;
32053565b12SHong Zhang   } else {
32153565b12SHong Zhang     A->ops->ptapnumeric = MatPtAPNumeric_SeqAIJ_SeqAIJ;
32253565b12SHong Zhang   }
323f79958ebSHong Zhang   ptap->api = api;
324f79958ebSHong Zhang   ptap->apj = apj;
3258d0a38e4SHong Zhang 
326d20bfe6fSHong Zhang   /* Clean up. */
327971236abSHong Zhang   ierr = MatDestroy(&Pt);CHKERRQ(ierr);
328971236abSHong Zhang   ierr = MatDestroy(&AP);CHKERRQ(ierr);
329eb9c0419SKris Buschelman   PetscFunctionReturn(0);
330eb9c0419SKris Buschelman }
331eb9c0419SKris Buschelman 
332f2535c29SHong Zhang /* #define PROFILE_MatPtAPNumeric */
333eb9c0419SKris Buschelman #undef __FUNCT__
3349af31e4aSHong Zhang #define __FUNCT__ "MatPtAPNumeric_SeqAIJ_SeqAIJ"
335dfbe8321SBarry Smith PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat A,Mat P,Mat C)
336dfbe8321SBarry Smith {
337dfbe8321SBarry Smith   PetscErrorCode    ierr;
338d20bfe6fSHong Zhang   Mat_SeqAIJ        *a  = (Mat_SeqAIJ *) A->data;
339f79958ebSHong Zhang   Mat_SeqAIJ        *p  = (Mat_SeqAIJ *) P->data;
340d20bfe6fSHong Zhang   Mat_SeqAIJ        *c  = (Mat_SeqAIJ *) C->data;
341b8f276daSHong Zhang   const PetscInt    *ai=a->i,*aj=a->j,*pi=p->i,*pj=p->j,*ci=c->i,*cj=c->j;
342b8f276daSHong Zhang   const PetscScalar *aa=a->a,*pa=p->a,*pval;
343b8f276daSHong Zhang   const PetscInt    *apj,*pcol,*cjj;
344b8f276daSHong Zhang   const PetscInt    am=A->rmap->N,cm=C->rmap->N;
345b8f276daSHong Zhang   PetscInt          i,j,k,anz,apnz,pnz,prow,crow,cnz;
346b8f276daSHong Zhang   PetscScalar       *apa,*ca=c->a,*caj,pvalj;
34753565b12SHong Zhang   Mat_PtAP          *ptap = c->ptap;
348f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric)
349f2535c29SHong Zhang   PetscLogDouble    t0,tf,time_Cseq0=0.0,time_Cseq1=0.0;
350f2535c29SHong Zhang   PetscInt          flops0=0,flops1=0;
351f2535c29SHong Zhang #endif
352eb9c0419SKris Buschelman 
353eb9c0419SKris Buschelman   PetscFunctionBegin;
3548d0a38e4SHong Zhang   /* Get temporary array for storage of one row of A*P */
3558d0a38e4SHong Zhang   apa = ptap->apa;
356d20bfe6fSHong Zhang 
357d20bfe6fSHong Zhang   /* Clear old values in C */
358d20bfe6fSHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
359d20bfe6fSHong Zhang 
360d20bfe6fSHong Zhang   for (i=0;i<am;i++) {
36178844af4SHong Zhang     /* Form sparse row of AP[i,:] = A[i,:]*P */
362f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric)
363f2535c29SHong Zhang     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
364f2535c29SHong Zhang #endif
36578844af4SHong Zhang     anz  = ai[i+1] - ai[i];
36678844af4SHong Zhang     apnz = 0;
36778844af4SHong Zhang     for (j=0; j<anz; j++) {
36878844af4SHong Zhang       prow = aj[j];
36978844af4SHong Zhang       pnz  = pi[prow+1] - pi[prow];
37078844af4SHong Zhang       pcol = pj + pi[prow];
37178844af4SHong Zhang       pval = pa + pi[prow];
37278844af4SHong Zhang       for (k=0; k<pnz; k++) {
37378844af4SHong Zhang         apa[pcol[k]] += aa[j]*pval[k];
374d20bfe6fSHong Zhang       }
37578844af4SHong Zhang       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
376f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric)
377f2535c29SHong Zhang       flops0 += 2.0*pnz;
378f2535c29SHong Zhang #endif
37978844af4SHong Zhang     }
38078844af4SHong Zhang     aj += anz; aa += anz;
381f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric)
382f2535c29SHong Zhang     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
383f2535c29SHong Zhang     time_Cseq0 += tf - t0;
384f2535c29SHong Zhang #endif
38578844af4SHong Zhang 
38678844af4SHong Zhang     /* Compute P^T*A*P using outer product P[i,:]^T*AP[i,:]. */
387f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric)
388f2535c29SHong Zhang     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
389f2535c29SHong Zhang #endif
390f79958ebSHong Zhang     apj  = ptap->apj + ptap->api[i];
391f79958ebSHong Zhang     apnz = ptap->api[i+1] - ptap->api[i];
39278844af4SHong Zhang     pnz  = pi[i+1] - pi[i];
39378844af4SHong Zhang     pcol = pj + pi[i];
39478844af4SHong Zhang     pval = pa + pi[i];
3958d0a38e4SHong Zhang 
39653565b12SHong Zhang     /* Perform dense axpy */
39753565b12SHong Zhang     for (j=0; j<pnz; j++) {
39853565b12SHong Zhang       crow  = pcol[j];
39953565b12SHong Zhang       cjj   = cj + ci[crow];
40053565b12SHong Zhang       caj   = ca + ci[crow];
401f2535c29SHong Zhang       pvalj = pval[j];
40253565b12SHong Zhang       cnz   = ci[crow+1] - ci[crow];
40353565b12SHong Zhang       for (k=0; k<cnz; k++){
404f2535c29SHong Zhang         caj[k] += pvalj*apa[cjj[k]];
40553565b12SHong Zhang       }
40653565b12SHong Zhang       ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
407f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric)
408f2535c29SHong Zhang       flops1 += 2.0*cnz;
409f2535c29SHong Zhang #endif
41053565b12SHong Zhang     }
411f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric)
412f2535c29SHong Zhang     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
413f2535c29SHong Zhang     time_Cseq1 += tf - t0;
414f2535c29SHong Zhang #endif
41553565b12SHong Zhang 
41653565b12SHong Zhang     /* Zero the current row info for A*P */
4171d633602SHong Zhang     for (j=0; j<apnz; j++) apa[apj[j]] = 0.0;
41853565b12SHong Zhang   }
41953565b12SHong Zhang 
42053565b12SHong Zhang   /* Assemble the final matrix and clean up */
42153565b12SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
42253565b12SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
423f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric)
424f2535c29SHong Zhang   printf("PtAPNumeric_SeqAIJ time %g + %g, flops %d %d\n",time_Cseq0,time_Cseq1,flops0,flops1);
425f2535c29SHong Zhang #endif
42653565b12SHong Zhang   PetscFunctionReturn(0);
42753565b12SHong Zhang }
42853565b12SHong Zhang 
42953565b12SHong Zhang #undef __FUNCT__
43053565b12SHong Zhang #define __FUNCT__ "MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy"
43153565b12SHong Zhang PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat A,Mat P,Mat C)
43253565b12SHong Zhang {
43353565b12SHong Zhang   PetscErrorCode    ierr;
43453565b12SHong Zhang   Mat_SeqAIJ        *a  = (Mat_SeqAIJ *) A->data;
43553565b12SHong Zhang   Mat_SeqAIJ        *p  = (Mat_SeqAIJ *) P->data;
43653565b12SHong Zhang   Mat_SeqAIJ        *c  = (Mat_SeqAIJ *) C->data;
437b8f276daSHong Zhang   const PetscInt    *ai=a->i,*aj=a->j,*pi=p->i,*pj=p->j,*ci=c->i,*cj=c->j;
438b8f276daSHong Zhang   const PetscScalar *aa=a->a,*pa=p->a,*pval;
439b8f276daSHong Zhang   const PetscInt    *apj,*pcol,*cjj;
440b8f276daSHong Zhang   const PetscInt    am=A->rmap->N,cm=C->rmap->N;
44153565b12SHong Zhang   PetscInt          i,j,k,anz,apnz,pnz,prow,crow,apcol,nextap;
442b8f276daSHong Zhang   PetscScalar       *apa,*ca=c->a,*caj,pvalj;
44353565b12SHong Zhang   Mat_PtAP          *ptap = c->ptap;
444f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric)
445f2535c29SHong Zhang   PetscLogDouble    t0,tf,time_Cseq0=0.0,time_Cseq1=0.0;
446f2535c29SHong Zhang   PetscInt          flops0=0,flops1=0;
447f2535c29SHong Zhang #endif
44853565b12SHong Zhang 
44953565b12SHong Zhang   PetscFunctionBegin;
45053565b12SHong Zhang   /* Get temporary array for storage of one row of A*P */
45153565b12SHong Zhang   apa = ptap->apa;
45253565b12SHong Zhang 
45353565b12SHong Zhang   /* Clear old values in C */
45453565b12SHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
45553565b12SHong Zhang 
45653565b12SHong Zhang   for (i=0;i<am;i++) {
45753565b12SHong Zhang     /* Form sparse row of AP[i,:] = A[i,:]*P */
458f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric)
459f2535c29SHong Zhang     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
460f2535c29SHong Zhang #endif
46153565b12SHong Zhang     anz  = ai[i+1] - ai[i];
46253565b12SHong Zhang     apnz = 0;
46353565b12SHong Zhang     for (j=0; j<anz; j++) {
46453565b12SHong Zhang       prow = aj[j];
46553565b12SHong Zhang       pnz  = pi[prow+1] - pi[prow];
46653565b12SHong Zhang       pcol = pj + pi[prow];
46753565b12SHong Zhang       pval = pa + pi[prow];
46853565b12SHong Zhang       for (k=0; k<pnz; k++) {
46953565b12SHong Zhang         apa[pcol[k]] += aa[j]*pval[k];
47053565b12SHong Zhang       }
47153565b12SHong Zhang       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
472f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric)
473f2535c29SHong Zhang       flops0 += 2.0*pnz;
474f2535c29SHong Zhang #endif
47553565b12SHong Zhang     }
47653565b12SHong Zhang     aj += anz; aa += anz;
477f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric)
478f2535c29SHong Zhang     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
479f2535c29SHong Zhang     time_Cseq0 += tf - t0;
480f2535c29SHong Zhang #endif
48153565b12SHong Zhang 
48253565b12SHong Zhang     /* Compute P^T*A*P using outer product P[i,:]^T*AP[i,:]. */
483f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric)
484f2535c29SHong Zhang     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
485f2535c29SHong Zhang #endif
48653565b12SHong Zhang     apj  = ptap->apj + ptap->api[i];
48753565b12SHong Zhang     apnz = ptap->api[i+1] - ptap->api[i];
48853565b12SHong Zhang     pnz  = pi[i+1] - pi[i];
48953565b12SHong Zhang     pcol = pj + pi[i];
49053565b12SHong Zhang     pval = pa + pi[i];
49153565b12SHong Zhang 
49253565b12SHong Zhang     /* Perform sparse axpy */
49378844af4SHong Zhang     for (j=0; j<pnz; j++) {
49478844af4SHong Zhang       crow   = pcol[j];
495d20bfe6fSHong Zhang       cjj    = cj + ci[crow];
496d20bfe6fSHong Zhang       caj    = ca + ci[crow];
497f2535c29SHong Zhang       pvalj  = pval[j];
4980d04baf8SBarry Smith       nextap = 1;
499b5b29cf7SHong Zhang       apcol  = apj[nextap];
50078844af4SHong Zhang       for (k=0; nextap<apnz; k++) {
501b5b29cf7SHong Zhang         if (cjj[k] == apcol) {
502f2535c29SHong Zhang           caj[k] += pvalj*apa[apcol];
5030d04baf8SBarry Smith           apcol   = apj[nextap++];
504d20bfe6fSHong Zhang         }
505d20bfe6fSHong Zhang       }
50678844af4SHong Zhang       ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
507f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric)
508f2535c29SHong Zhang       flops1 += 2.0*apnz;
509f2535c29SHong Zhang #endif
510f79958ebSHong Zhang     }
511f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric)
512f2535c29SHong Zhang     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
513f2535c29SHong Zhang     time_Cseq1 += tf - t0;
514f2535c29SHong Zhang #endif
515d20bfe6fSHong Zhang 
516d20bfe6fSHong Zhang     /* Zero the current row info for A*P */
51778844af4SHong Zhang     for (j=0; j<apnz; j++) {
518b5b29cf7SHong Zhang       apcol      = apj[j];
519b5b29cf7SHong Zhang       apa[apcol] = 0.;
520d20bfe6fSHong Zhang     }
521d20bfe6fSHong Zhang   }
522d20bfe6fSHong Zhang 
523d20bfe6fSHong Zhang   /* Assemble the final matrix and clean up */
524d20bfe6fSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
525d20bfe6fSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
526f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric)
527f2535c29SHong Zhang   printf("MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy time %g + %g, flops %d %d\n",time_Cseq0,time_Cseq1,flops0,flops1);
528f2535c29SHong Zhang #endif
529eb9c0419SKris Buschelman   PetscFunctionReturn(0);
530eb9c0419SKris Buschelman }
531