xref: /petsc/src/mat/impls/aij/seq/matptap.c (revision dadf0e6bf3fbf8619c99f30f0124ab4e08bd7289)
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);
15253565b12SHong Zhang 
15353565b12SHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
15453565b12SHong Zhang   /* Since these are PETSc arrays, change flags to free them as necessary. */
15553565b12SHong Zhang   c = (Mat_SeqAIJ *)((*C)->data);
15653565b12SHong Zhang   c->free_a  = PETSC_TRUE;
15753565b12SHong Zhang   c->free_ij = PETSC_TRUE;
15853565b12SHong Zhang   c->nonew   = 0;
15953565b12SHong Zhang   A->ops->ptapnumeric = MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy2; /* should use *C->ops until PtAP insterface is updated to double dispatch as MatMatMult() */
16053565b12SHong Zhang 
16153565b12SHong Zhang   /* Clean up. */
16253565b12SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
16353565b12SHong Zhang #if defined(PETSC_USE_INFO)
16453565b12SHong Zhang   if (ci[pn] != 0) {
16553565b12SHong Zhang     PetscReal afill = ((PetscReal)ci[pn])/ai[am];
16653565b12SHong Zhang     if (afill < 1.0) afill = 1.0;
16753565b12SHong Zhang     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
16853565b12SHong Zhang     ierr = PetscInfo1((*C),"Use MatPtAP(A,P,MatReuse,%G,&C) for best performance.\n",afill);CHKERRQ(ierr);
16953565b12SHong Zhang   } else {
17053565b12SHong Zhang     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
17153565b12SHong Zhang   }
172f79958ebSHong Zhang #endif
17353565b12SHong Zhang   PetscFunctionReturn(0);
17453565b12SHong Zhang }
17553565b12SHong Zhang 
17653565b12SHong Zhang #undef __FUNCT__
17753565b12SHong Zhang #define __FUNCT__ "MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy2"
17853565b12SHong Zhang PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy2(Mat A,Mat P,Mat C)
17953565b12SHong Zhang {
18053565b12SHong Zhang   PetscErrorCode ierr;
18153565b12SHong Zhang   Mat_SeqAIJ     *a  = (Mat_SeqAIJ *) A->data;
18253565b12SHong Zhang   Mat_SeqAIJ     *p  = (Mat_SeqAIJ *) P->data;
18353565b12SHong Zhang   Mat_SeqAIJ     *c  = (Mat_SeqAIJ *) C->data;
18453565b12SHong Zhang   PetscInt       *ai=a->i,*aj=a->j,*apj,*apjdense,*pi=p->i,*pj=p->j,*pJ=p->j,*pjj;
18553565b12SHong Zhang   PetscInt       *ci=c->i,*cj=c->j,*cjj;
18653565b12SHong Zhang   PetscInt       am=A->rmap->N,cn=C->cmap->N,cm=C->rmap->N;
18753565b12SHong Zhang   PetscInt       i,j,k,anzi,pnzi,apnzj,nextap,pnzj,prow,crow;
18853565b12SHong Zhang   MatScalar      *aa=a->a,*apa,*pa=p->a,*pA=p->a,*paj,*ca=c->a,*caj;
18953565b12SHong Zhang 
19053565b12SHong Zhang   PetscFunctionBegin;
19153565b12SHong Zhang   /* Allocate temporary array for storage of one row of A*P */
192db1ebd87SHong Zhang   ierr = PetscMalloc(cn*(sizeof(MatScalar)+sizeof(PetscInt))+c->rmax*sizeof(PetscInt),&apa);CHKERRQ(ierr);
19353565b12SHong Zhang 
194db1ebd87SHong Zhang   apjdense = (PetscInt *)(apa + cn);
195db1ebd87SHong Zhang   apj      = apjdense + cn;
196db1ebd87SHong Zhang   ierr = PetscMemzero(apa,cn*(sizeof(MatScalar)+sizeof(PetscInt)));CHKERRQ(ierr);
19753565b12SHong Zhang 
19853565b12SHong Zhang   /* Clear old values in C */
19953565b12SHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
20053565b12SHong Zhang 
20153565b12SHong Zhang   for (i=0; i<am; i++) {
20253565b12SHong Zhang     /* Form sparse row of A*P */
20353565b12SHong Zhang     anzi  = ai[i+1] - ai[i];
20453565b12SHong Zhang     apnzj = 0;
20553565b12SHong Zhang     for (j=0; j<anzi; j++) {
20653565b12SHong Zhang       prow = *aj++;
20753565b12SHong Zhang       pnzj = pi[prow+1] - pi[prow];
20853565b12SHong Zhang       pjj  = pj + pi[prow];
20953565b12SHong Zhang       paj  = pa + pi[prow];
21053565b12SHong Zhang       for (k=0;k<pnzj;k++) {
21153565b12SHong Zhang         if (!apjdense[pjj[k]]) {
21253565b12SHong Zhang           apjdense[pjj[k]] = -1;
21353565b12SHong Zhang           apj[apnzj++]     = pjj[k];
21453565b12SHong Zhang         }
21553565b12SHong Zhang         apa[pjj[k]] += (*aa)*paj[k];
21653565b12SHong Zhang       }
21753565b12SHong Zhang       ierr = PetscLogFlops(2.0*pnzj);CHKERRQ(ierr);
21853565b12SHong Zhang       aa++;
21953565b12SHong Zhang     }
22053565b12SHong Zhang 
22153565b12SHong Zhang     /* Sort the j index array for quick sparse axpy. */
22253565b12SHong Zhang     /* Note: a array does not need sorting as it is in dense storage locations. */
22353565b12SHong Zhang     ierr = PetscSortInt(apnzj,apj);CHKERRQ(ierr);
22453565b12SHong Zhang 
22553565b12SHong Zhang     /* Compute P^T*A*P using outer product (P^T)[:,j]*(A*P)[j,:]. */
22653565b12SHong Zhang     pnzi = pi[i+1] - pi[i];
22753565b12SHong Zhang     for (j=0; j<pnzi; j++) {
22853565b12SHong Zhang       nextap = 0;
22953565b12SHong Zhang       crow   = *pJ++;
23053565b12SHong Zhang       cjj    = cj + ci[crow];
23153565b12SHong Zhang       caj    = ca + ci[crow];
23253565b12SHong Zhang       /* Perform sparse axpy operation.  Note cjj includes apj. */
23353565b12SHong Zhang       for (k=0;nextap<apnzj;k++) {
23453565b12SHong Zhang #if defined(PETSC_USE_DEBUG)
23553565b12SHong Zhang         if (k >= ci[crow+1] - ci[crow]) {
23653565b12SHong Zhang           SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"k too large k %d, crow %d",k,crow);
23753565b12SHong Zhang         }
23853565b12SHong Zhang #endif
23953565b12SHong Zhang         if (cjj[k]==apj[nextap]) {
24053565b12SHong Zhang           caj[k] += (*pA)*apa[apj[nextap++]];
24153565b12SHong Zhang         }
24253565b12SHong Zhang       }
24353565b12SHong Zhang       ierr = PetscLogFlops(2.0*apnzj);CHKERRQ(ierr);
24453565b12SHong Zhang       pA++;
24553565b12SHong Zhang     }
24653565b12SHong Zhang 
24753565b12SHong Zhang     /* Zero the current row info for A*P */
24853565b12SHong Zhang     for (j=0; j<apnzj; j++) {
24953565b12SHong Zhang       apa[apj[j]]      = 0.;
25053565b12SHong Zhang       apjdense[apj[j]] = 0;
25153565b12SHong Zhang     }
25253565b12SHong Zhang   }
25353565b12SHong Zhang 
25453565b12SHong Zhang   /* Assemble the final matrix and clean up */
25553565b12SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
25653565b12SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
25753565b12SHong Zhang   ierr = PetscFree(apa);CHKERRQ(ierr);
25853565b12SHong Zhang   PetscFunctionReturn(0);
25953565b12SHong Zhang }
26053565b12SHong Zhang 
26153565b12SHong Zhang #undef __FUNCT__
26253565b12SHong Zhang #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqAIJ"
26353565b12SHong Zhang PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
26453565b12SHong Zhang {
26553565b12SHong Zhang   PetscErrorCode     ierr;
26653565b12SHong Zhang   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data,*p = (Mat_SeqAIJ*)P->data,*c;
2679096154bSHong Zhang   PetscInt           *pti,*ptj,*ai=a->i,*aj=a->j,*pi=p->i,*pj=p->j,*api,*apj;
26853565b12SHong Zhang   PetscInt           *ci,*cj,ndouble_ap,ndouble_ptap;
26953565b12SHong Zhang   PetscInt           an=A->cmap->N,am=A->rmap->N,pn=P->cmap->N;
27053565b12SHong Zhang   MatScalar          *ca;
27153565b12SHong Zhang   Mat_PtAP           *ptap;
27253565b12SHong Zhang   PetscInt           sparse_axpy=0;
27377584fe6SHong Zhang   PetscLogDouble     t0,tf,time_Trans=0.0,time_GetSymbolic1=0.0,time_GetSymbolic2=0.0;
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);
28453565b12SHong Zhang     PetscFunctionReturn(0);
28553565b12SHong Zhang   }
28653565b12SHong Zhang 
28753565b12SHong Zhang   /* Get ij structure of Pt = P^T */
28877584fe6SHong Zhang   ierr = PetscGetTime(&t0);CHKERRQ(ierr);
28953565b12SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
29077584fe6SHong Zhang   ierr = PetscGetTime(&tf);CHKERRQ(ierr);
29177584fe6SHong Zhang   time_Trans += tf - t0;
29253565b12SHong Zhang 
29353565b12SHong Zhang   /* Get structure of AP = A*P */
29477584fe6SHong Zhang   ierr = PetscGetTime(&t0);CHKERRQ(ierr);
29553565b12SHong Zhang   ierr = MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ(am,ai,aj,an,pn,pi,pj,fill,&api,&apj,&ndouble_ap);CHKERRQ(ierr);
29677584fe6SHong Zhang   ierr = PetscGetTime(&tf);CHKERRQ(ierr);
29777584fe6SHong Zhang   time_GetSymbolic1 += tf - t0;
29853565b12SHong Zhang 
29953565b12SHong Zhang   /* Get structure of C = Pt*AP */
30077584fe6SHong Zhang   ierr = PetscGetTime(&t0);CHKERRQ(ierr);
30153565b12SHong Zhang   ierr = MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ(pn,pti,ptj,am,pn,api,apj,fill,&ci,&cj,&ndouble_ptap);CHKERRQ(ierr);
30277584fe6SHong Zhang   ierr = PetscGetTime(&tf);CHKERRQ(ierr);
30377584fe6SHong Zhang   time_GetSymbolic2 += tf - t0;
304f79958ebSHong Zhang 
305d20bfe6fSHong Zhang   /* Allocate space for ca */
306d20bfe6fSHong Zhang   ierr = PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
307d20bfe6fSHong Zhang   ierr = PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));CHKERRQ(ierr);
308d20bfe6fSHong Zhang 
309d20bfe6fSHong Zhang   /* put together the new matrix */
3107adad957SLisandro Dalcin   ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,pn,pn,ci,cj,ca,C);CHKERRQ(ierr);
311d20bfe6fSHong Zhang 
312d20bfe6fSHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
313d20bfe6fSHong Zhang   /* Since these are PETSc arrays, change flags to free them as necessary. */
3148d0a38e4SHong Zhang   c          = (Mat_SeqAIJ *)(*C)->data;
315e6b907acSBarry Smith   c->free_a  = PETSC_TRUE;
316e6b907acSBarry Smith   c->free_ij = PETSC_TRUE;
317d20bfe6fSHong Zhang   c->nonew   = 0;
318d20bfe6fSHong Zhang 
31953565b12SHong Zhang   /* create a supporting struct for reuse by MatPtAPNumeric() */
3208d0a38e4SHong Zhang   ierr = PetscNew(Mat_PtAP,&ptap);CHKERRQ(ierr);
32153565b12SHong Zhang   c->ptap            = ptap;
3228d0a38e4SHong Zhang   ptap->destroy      = (*C)->ops->destroy;
3238d0a38e4SHong Zhang   (*C)->ops->destroy = MatDestroy_SeqAIJ_PtAP;
3248d0a38e4SHong Zhang 
3258d0a38e4SHong Zhang   /* Allocate temporary array for storage of one row of A*P */
326db1ebd87SHong Zhang   ierr = PetscMalloc((pn+1)*sizeof(PetscScalar),&ptap->apa);CHKERRQ(ierr);
327db1ebd87SHong Zhang   ierr = PetscMemzero(ptap->apa,(pn+1)*sizeof(PetscScalar));CHKERRQ(ierr);
32853565b12SHong Zhang   if (sparse_axpy == 1){
32953565b12SHong Zhang     A->ops->ptapnumeric = MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy;
33053565b12SHong Zhang   } else {
33153565b12SHong Zhang     A->ops->ptapnumeric = MatPtAPNumeric_SeqAIJ_SeqAIJ;
33253565b12SHong Zhang   }
333f79958ebSHong Zhang   ptap->api = api;
334f79958ebSHong Zhang   ptap->apj = apj;
3358d0a38e4SHong Zhang 
336d20bfe6fSHong Zhang   /* Clean up. */
337d20bfe6fSHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
338f2b054eeSHong Zhang #if defined(PETSC_USE_INFO)
339f2b054eeSHong Zhang   if (ci[pn] != 0) {
340f79958ebSHong Zhang     PetscReal apfill,ptapfill;
341f79958ebSHong Zhang     apfill = ((PetscReal)api[am])/(ai[am]+pi[an]);
342f79958ebSHong Zhang     if (apfill < 1.0) apfill = 1.0;
343f79958ebSHong Zhang     ierr = PetscInfo3((*C),"A*P: Reallocs %D; Fill ratio: given %G needed %G.\n",ndouble_ap,fill,apfill);CHKERRQ(ierr);
344f79958ebSHong Zhang     ptapfill = ((PetscReal)ci[pn])/(pi[an]+api[am]);
345f79958ebSHong Zhang     if (ptapfill < 1.0) ptapfill = 1.0;
346f79958ebSHong Zhang     ierr = PetscInfo3((*C),"Pt*AP: Reallocs %D; Fill ratio: given %G needed %G.\n",ndouble_ptap,fill,ptapfill);CHKERRQ(ierr);
347f79958ebSHong Zhang 
348f79958ebSHong Zhang     ierr = PetscInfo1((*C),"Use MatPtAP(A,P,MatReuse,%G,&C) for best performance.\n",PetscMax(apfill,ptapfill));CHKERRQ(ierr);
349f79958ebSHong Zhang     ierr = PetscInfo4((*C),"nonzeros: A %D, P %D, A*P %D, C=PtAP %D\n",ai[am],pi[an],api[am],ci[pn]);CHKERRQ(ierr);
350f2b054eeSHong Zhang   } else {
351f2b054eeSHong Zhang     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
352f2b054eeSHong Zhang   }
353f2b054eeSHong Zhang #endif
354b8f276daSHong Zhang   /*
35577584fe6SHong Zhang   printf("MatPtAPSymbolic_SeqAIJ_SeqAIJ time %g + %g + %g\n",time_Trans,time_GetSymbolic1,time_GetSymbolic2);
356b8f276daSHong Zhang    */
357eb9c0419SKris Buschelman   PetscFunctionReturn(0);
358eb9c0419SKris Buschelman }
359eb9c0419SKris Buschelman 
360f2535c29SHong Zhang /* #define PROFILE_MatPtAPNumeric */
361eb9c0419SKris Buschelman #undef __FUNCT__
3629af31e4aSHong Zhang #define __FUNCT__ "MatPtAPNumeric_SeqAIJ_SeqAIJ"
363dfbe8321SBarry Smith PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat A,Mat P,Mat C)
364dfbe8321SBarry Smith {
365dfbe8321SBarry Smith   PetscErrorCode    ierr;
366d20bfe6fSHong Zhang   Mat_SeqAIJ        *a  = (Mat_SeqAIJ *) A->data;
367f79958ebSHong Zhang   Mat_SeqAIJ        *p  = (Mat_SeqAIJ *) P->data;
368d20bfe6fSHong Zhang   Mat_SeqAIJ        *c  = (Mat_SeqAIJ *) C->data;
369b8f276daSHong Zhang   const PetscInt    *ai=a->i,*aj=a->j,*pi=p->i,*pj=p->j,*ci=c->i,*cj=c->j;
370b8f276daSHong Zhang   const PetscScalar *aa=a->a,*pa=p->a,*pval;
371b8f276daSHong Zhang   const PetscInt    *apj,*pcol,*cjj;
372b8f276daSHong Zhang   const PetscInt    am=A->rmap->N,cm=C->rmap->N;
373b8f276daSHong Zhang   PetscInt          i,j,k,anz,apnz,pnz,prow,crow,cnz;
374b8f276daSHong Zhang   PetscScalar       *apa,*ca=c->a,*caj,pvalj;
37553565b12SHong Zhang   Mat_PtAP          *ptap = c->ptap;
376f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric)
377f2535c29SHong Zhang   PetscLogDouble    t0,tf,time_Cseq0=0.0,time_Cseq1=0.0;
378f2535c29SHong Zhang   PetscInt          flops0=0,flops1=0;
379f2535c29SHong Zhang #endif
380eb9c0419SKris Buschelman 
381eb9c0419SKris Buschelman   PetscFunctionBegin;
3828d0a38e4SHong Zhang   /* Get temporary array for storage of one row of A*P */
3838d0a38e4SHong Zhang   apa = ptap->apa;
384d20bfe6fSHong Zhang 
385d20bfe6fSHong Zhang   /* Clear old values in C */
386d20bfe6fSHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
387d20bfe6fSHong Zhang 
388d20bfe6fSHong Zhang   for (i=0;i<am;i++) {
38978844af4SHong Zhang     /* Form sparse row of AP[i,:] = A[i,:]*P */
390f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric)
391f2535c29SHong Zhang     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
392f2535c29SHong Zhang #endif
39378844af4SHong Zhang     anz  = ai[i+1] - ai[i];
39478844af4SHong Zhang     apnz = 0;
39578844af4SHong Zhang     for (j=0; j<anz; j++) {
39678844af4SHong Zhang       prow = aj[j];
39778844af4SHong Zhang       pnz  = pi[prow+1] - pi[prow];
39878844af4SHong Zhang       pcol = pj + pi[prow];
39978844af4SHong Zhang       pval = pa + pi[prow];
40078844af4SHong Zhang       for (k=0; k<pnz; k++) {
40178844af4SHong Zhang         apa[pcol[k]] += aa[j]*pval[k];
402d20bfe6fSHong Zhang       }
40378844af4SHong Zhang       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
404f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric)
405f2535c29SHong Zhang       flops0 += 2.0*pnz;
406f2535c29SHong Zhang #endif
40778844af4SHong Zhang     }
40878844af4SHong Zhang     aj += anz; aa += anz;
409f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric)
410f2535c29SHong Zhang     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
411f2535c29SHong Zhang     time_Cseq0 += tf - t0;
412f2535c29SHong Zhang #endif
41378844af4SHong Zhang 
41478844af4SHong Zhang     /* Compute P^T*A*P using outer product P[i,:]^T*AP[i,:]. */
415f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric)
416f2535c29SHong Zhang     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
417f2535c29SHong Zhang #endif
418f79958ebSHong Zhang     apj  = ptap->apj + ptap->api[i];
419f79958ebSHong Zhang     apnz = ptap->api[i+1] - ptap->api[i];
42078844af4SHong Zhang     pnz  = pi[i+1] - pi[i];
42178844af4SHong Zhang     pcol = pj + pi[i];
42278844af4SHong Zhang     pval = pa + pi[i];
4238d0a38e4SHong Zhang 
42453565b12SHong Zhang     /* Perform dense axpy */
42553565b12SHong Zhang     for (j=0; j<pnz; j++) {
42653565b12SHong Zhang       crow  = pcol[j];
42753565b12SHong Zhang       cjj   = cj + ci[crow];
42853565b12SHong Zhang       caj   = ca + ci[crow];
429f2535c29SHong Zhang       pvalj = pval[j];
43053565b12SHong Zhang       cnz   = ci[crow+1] - ci[crow];
43153565b12SHong Zhang       for (k=0; k<cnz; k++){
432f2535c29SHong Zhang         caj[k] += pvalj*apa[cjj[k]];
43353565b12SHong Zhang       }
43453565b12SHong Zhang       ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
435f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric)
436f2535c29SHong Zhang       flops1 += 2.0*cnz;
437f2535c29SHong Zhang #endif
43853565b12SHong Zhang     }
439f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric)
440f2535c29SHong Zhang     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
441f2535c29SHong Zhang     time_Cseq1 += tf - t0;
442f2535c29SHong Zhang #endif
44353565b12SHong Zhang 
44453565b12SHong Zhang     /* Zero the current row info for A*P */
4451d633602SHong Zhang     for (j=0; j<apnz; j++) apa[apj[j]] = 0.0;
44653565b12SHong Zhang   }
44753565b12SHong Zhang 
44853565b12SHong Zhang   /* Assemble the final matrix and clean up */
44953565b12SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
45053565b12SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
451f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric)
452f2535c29SHong Zhang   printf("PtAPNumeric_SeqAIJ time %g + %g, flops %d %d\n",time_Cseq0,time_Cseq1,flops0,flops1);
453f2535c29SHong Zhang #endif
45453565b12SHong Zhang   PetscFunctionReturn(0);
45553565b12SHong Zhang }
45653565b12SHong Zhang 
45753565b12SHong Zhang #undef __FUNCT__
45853565b12SHong Zhang #define __FUNCT__ "MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy"
45953565b12SHong Zhang PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat A,Mat P,Mat C)
46053565b12SHong Zhang {
46153565b12SHong Zhang   PetscErrorCode    ierr;
46253565b12SHong Zhang   Mat_SeqAIJ        *a  = (Mat_SeqAIJ *) A->data;
46353565b12SHong Zhang   Mat_SeqAIJ        *p  = (Mat_SeqAIJ *) P->data;
46453565b12SHong Zhang   Mat_SeqAIJ        *c  = (Mat_SeqAIJ *) C->data;
465b8f276daSHong Zhang   const PetscInt    *ai=a->i,*aj=a->j,*pi=p->i,*pj=p->j,*ci=c->i,*cj=c->j;
466b8f276daSHong Zhang   const PetscScalar *aa=a->a,*pa=p->a,*pval;
467b8f276daSHong Zhang   const PetscInt    *apj,*pcol,*cjj;
468b8f276daSHong Zhang   const PetscInt    am=A->rmap->N,cm=C->rmap->N;
46953565b12SHong Zhang   PetscInt          i,j,k,anz,apnz,pnz,prow,crow,apcol,nextap;
470b8f276daSHong Zhang   PetscScalar       *apa,*ca=c->a,*caj,pvalj;
47153565b12SHong Zhang   Mat_PtAP          *ptap = c->ptap;
472f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric)
473f2535c29SHong Zhang   PetscLogDouble    t0,tf,time_Cseq0=0.0,time_Cseq1=0.0;
474f2535c29SHong Zhang   PetscInt          flops0=0,flops1=0;
475f2535c29SHong Zhang #endif
47653565b12SHong Zhang 
47753565b12SHong Zhang   PetscFunctionBegin;
47853565b12SHong Zhang   /* Get temporary array for storage of one row of A*P */
47953565b12SHong Zhang   apa = ptap->apa;
48053565b12SHong Zhang 
48153565b12SHong Zhang   /* Clear old values in C */
48253565b12SHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
48353565b12SHong Zhang 
48453565b12SHong Zhang   for (i=0;i<am;i++) {
48553565b12SHong Zhang     /* Form sparse row of AP[i,:] = A[i,:]*P */
486f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric)
487f2535c29SHong Zhang     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
488f2535c29SHong Zhang #endif
48953565b12SHong Zhang     anz  = ai[i+1] - ai[i];
49053565b12SHong Zhang     apnz = 0;
49153565b12SHong Zhang     for (j=0; j<anz; j++) {
49253565b12SHong Zhang       prow = aj[j];
49353565b12SHong Zhang       pnz  = pi[prow+1] - pi[prow];
49453565b12SHong Zhang       pcol = pj + pi[prow];
49553565b12SHong Zhang       pval = pa + pi[prow];
49653565b12SHong Zhang       for (k=0; k<pnz; k++) {
49753565b12SHong Zhang         apa[pcol[k]] += aa[j]*pval[k];
49853565b12SHong Zhang       }
49953565b12SHong Zhang       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
500f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric)
501f2535c29SHong Zhang       flops0 += 2.0*pnz;
502f2535c29SHong Zhang #endif
50353565b12SHong Zhang     }
50453565b12SHong Zhang     aj += anz; aa += anz;
505f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric)
506f2535c29SHong Zhang     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
507f2535c29SHong Zhang     time_Cseq0 += tf - t0;
508f2535c29SHong Zhang #endif
50953565b12SHong Zhang 
51053565b12SHong Zhang     /* Compute P^T*A*P using outer product P[i,:]^T*AP[i,:]. */
511f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric)
512f2535c29SHong Zhang     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
513f2535c29SHong Zhang #endif
51453565b12SHong Zhang     apj  = ptap->apj + ptap->api[i];
51553565b12SHong Zhang     apnz = ptap->api[i+1] - ptap->api[i];
51653565b12SHong Zhang     pnz  = pi[i+1] - pi[i];
51753565b12SHong Zhang     pcol = pj + pi[i];
51853565b12SHong Zhang     pval = pa + pi[i];
51953565b12SHong Zhang 
52053565b12SHong Zhang     /* Perform sparse axpy */
52178844af4SHong Zhang     for (j=0; j<pnz; j++) {
52278844af4SHong Zhang       crow   = pcol[j];
523d20bfe6fSHong Zhang       cjj    = cj + ci[crow];
524d20bfe6fSHong Zhang       caj    = ca + ci[crow];
525f2535c29SHong Zhang       pvalj  = pval[j];
52678844af4SHong Zhang       nextap = 0;
527*b5b29cf7SHong Zhang       apcol  = apj[nextap];
52878844af4SHong Zhang       for (k=0; nextap<apnz; k++) {
529*b5b29cf7SHong Zhang         if (cjj[k] == apcol) {
530f2535c29SHong Zhang           caj[k] += pvalj*apa[apcol];
531b8f276daSHong Zhang           nextap++;
532b8f276daSHong Zhang           apcol   = apj[nextap]; /* apj is size of AP->i[m]+1, an extra space is available in case nextap=apnz */
533d20bfe6fSHong Zhang         }
534d20bfe6fSHong Zhang       }
53578844af4SHong Zhang       ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
536f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric)
537f2535c29SHong Zhang       flops1 += 2.0*apnz;
538f2535c29SHong Zhang #endif
539f79958ebSHong Zhang     }
540f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric)
541f2535c29SHong Zhang     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
542f2535c29SHong Zhang     time_Cseq1 += tf - t0;
543f2535c29SHong Zhang #endif
544d20bfe6fSHong Zhang 
545d20bfe6fSHong Zhang     /* Zero the current row info for A*P */
54678844af4SHong Zhang     for (j=0; j<apnz; j++) {
547*b5b29cf7SHong Zhang       apcol      = apj[j];
548*b5b29cf7SHong Zhang       apa[apcol] = 0.;
549d20bfe6fSHong Zhang     }
550d20bfe6fSHong Zhang   }
551d20bfe6fSHong Zhang 
552d20bfe6fSHong Zhang   /* Assemble the final matrix and clean up */
553d20bfe6fSHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
554d20bfe6fSHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
555f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric)
556f2535c29SHong Zhang   printf("MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy time %g + %g, flops %d %d\n",time_Cseq0,time_Cseq1,flops0,flops1);
557f2535c29SHong Zhang #endif
558eb9c0419SKris Buschelman   PetscFunctionReturn(0);
559eb9c0419SKris Buschelman }
560