xref: /petsc/src/mat/impls/aij/seq/matptap.c (revision 2205254efee3a00a594e5e2a3a70f74dcb40bc03)
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 = MatPtAPSymbolic_SeqAIJ_SeqAIJ(A,P,fill,C);CHKERRQ(ierr);
207ba1a0bfSKris Buschelman   }
21ae32dd66SHong Zhang   ierr = (*(*C)->ops->ptapnumeric)(A,P,*C);CHKERRQ(ierr);
227ba1a0bfSKris Buschelman   PetscFunctionReturn(0);
237ba1a0bfSKris Buschelman }
247ba1a0bfSKris Buschelman 
257ba1a0bfSKris Buschelman #undef __FUNCT__
2653565b12SHong Zhang #define __FUNCT__ "MatDestroy_SeqAIJ_PtAP"
2753565b12SHong Zhang PetscErrorCode MatDestroy_SeqAIJ_PtAP(Mat A)
288d0a38e4SHong Zhang {
298d0a38e4SHong Zhang   PetscErrorCode ierr;
3053565b12SHong Zhang   Mat_SeqAIJ     *a    = (Mat_SeqAIJ*)A->data;
3153565b12SHong Zhang   Mat_PtAP       *ptap = a->ptap;
328d0a38e4SHong Zhang 
338d0a38e4SHong Zhang   PetscFunctionBegin;
348d0a38e4SHong Zhang   ierr = PetscFree(ptap->apa);CHKERRQ(ierr);
35f79958ebSHong Zhang   ierr = PetscFree(ptap->api);CHKERRQ(ierr);
36f79958ebSHong Zhang   ierr = PetscFree(ptap->apj);CHKERRQ(ierr);
3753565b12SHong Zhang   ierr = (ptap->destroy)(A);CHKERRQ(ierr);
388d0a38e4SHong Zhang   ierr = PetscFree(ptap);CHKERRQ(ierr);
398d0a38e4SHong Zhang   PetscFunctionReturn(0);
408d0a38e4SHong Zhang }
418d0a38e4SHong Zhang 
428d0a38e4SHong Zhang #undef __FUNCT__
43ae32dd66SHong Zhang #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy"
44ae32dd66SHong Zhang PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(Mat A,Mat P,PetscReal fill,Mat *C)
458d0a38e4SHong Zhang {
468d0a38e4SHong Zhang   PetscErrorCode     ierr;
4753565b12SHong Zhang   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
48d20bfe6fSHong Zhang   Mat_SeqAIJ         *a        = (Mat_SeqAIJ*)A->data,*p = (Mat_SeqAIJ*)P->data,*c;
4953565b12SHong Zhang   PetscInt           *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj;
5053565b12SHong Zhang   PetscInt           *ci,*cj,*ptadenserow,*ptasparserow,*ptaj,nspacedouble=0;
51d0f46423SBarry Smith   PetscInt           an=A->cmap->N,am=A->rmap->N,pn=P->cmap->N;
5253565b12SHong Zhang   PetscInt           i,j,k,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi,nlnk,*lnk;
53d20bfe6fSHong Zhang   MatScalar          *ca;
5453565b12SHong Zhang   PetscBT            lnkbt;
55eb9c0419SKris Buschelman 
56eb9c0419SKris Buschelman   PetscFunctionBegin;
5753565b12SHong Zhang   /* Get ij structure of P^T */
58eb9c0419SKris Buschelman   ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
59d20bfe6fSHong Zhang   ptJ  = ptj;
60eb9c0419SKris Buschelman 
61d20bfe6fSHong Zhang   /* Allocate ci array, arrays for fill computation and */
62d20bfe6fSHong Zhang   /* free space for accumulating nonzero column info */
6355a3bba9SHong Zhang   ierr  = PetscMalloc((pn+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
64d20bfe6fSHong Zhang   ci[0] = 0;
65eb9c0419SKris Buschelman 
6655a3bba9SHong Zhang   ierr         = PetscMalloc((2*an+1)*sizeof(PetscInt),&ptadenserow);CHKERRQ(ierr);
6755a3bba9SHong Zhang   ierr         = PetscMemzero(ptadenserow,(2*an+1)*sizeof(PetscInt));CHKERRQ(ierr);
68d20bfe6fSHong Zhang   ptasparserow = ptadenserow  + an;
6955a3bba9SHong Zhang 
7055a3bba9SHong Zhang   /* create and initialize a linked list */
7155a3bba9SHong Zhang   nlnk = pn+1;
7255a3bba9SHong Zhang   ierr = PetscLLCreate(pn,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr);
73eb9c0419SKris Buschelman 
74f2b054eeSHong Zhang   /* Set initial free space to be fill*nnz(A). */
75d20bfe6fSHong Zhang   /* This should be reasonable if sparsity of PtAP is similar to that of A. */
76a2ea699eSBarry Smith   ierr          = PetscFreeSpaceGet((PetscInt)(fill*ai[am]),&free_space);CHKERRQ(ierr);
77d20bfe6fSHong Zhang   current_space = free_space;
78d20bfe6fSHong Zhang 
79d20bfe6fSHong Zhang   /* Determine symbolic info for each row of C: */
80d20bfe6fSHong Zhang   for (i=0;i<pn;i++) {
81d20bfe6fSHong Zhang     ptnzi  = pti[i+1] - pti[i];
82d20bfe6fSHong Zhang     ptanzi = 0;
83d20bfe6fSHong Zhang     /* Determine symbolic row of PtA: */
84d20bfe6fSHong Zhang     for (j=0;j<ptnzi;j++) {
85d20bfe6fSHong Zhang       arow = *ptJ++;
86d20bfe6fSHong Zhang       anzj = ai[arow+1] - ai[arow];
87d20bfe6fSHong Zhang       ajj  = aj + ai[arow];
88d20bfe6fSHong Zhang       for (k=0;k<anzj;k++) {
89d20bfe6fSHong Zhang         if (!ptadenserow[ajj[k]]) {
90d20bfe6fSHong Zhang           ptadenserow[ajj[k]]    = -1;
91d20bfe6fSHong Zhang           ptasparserow[ptanzi++] = ajj[k];
92d20bfe6fSHong Zhang         }
93d20bfe6fSHong Zhang       }
94d20bfe6fSHong Zhang     }
95d20bfe6fSHong Zhang     /* Using symbolic info for row of PtA, determine symbolic info for row of C: */
96d20bfe6fSHong Zhang     ptaj = ptasparserow;
97d20bfe6fSHong Zhang     cnzi = 0;
98d20bfe6fSHong Zhang     for (j=0; j<ptanzi; j++) {
99d20bfe6fSHong Zhang       prow = *ptaj++;
100d20bfe6fSHong Zhang       pnzj = pi[prow+1] - pi[prow];
101d20bfe6fSHong Zhang       pjj  = pj + pi[prow];
10255a3bba9SHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
103dadf0e6bSHong Zhang       ierr  = PetscLLAddSorted(pnzj,pjj,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr);
10455a3bba9SHong Zhang       cnzi += nlnk;
105d20bfe6fSHong Zhang     }
106d20bfe6fSHong Zhang 
107d20bfe6fSHong Zhang     /* If free space is not available, make more free space */
108d20bfe6fSHong Zhang     /* Double the amount of total space in the list */
109d20bfe6fSHong Zhang     if (current_space->local_remaining<cnzi) {
1104238b7adSHong Zhang       ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,&current_space);CHKERRQ(ierr);
111f2b054eeSHong Zhang       nspacedouble++;
112d20bfe6fSHong Zhang     }
113d20bfe6fSHong Zhang 
114d20bfe6fSHong Zhang     /* Copy data into free space, and zero out denserows */
11555a3bba9SHong Zhang     ierr = PetscLLClean(pn,pn,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
116*2205254eSKarl Rupp 
117d20bfe6fSHong Zhang     current_space->array           += cnzi;
118d20bfe6fSHong Zhang     current_space->local_used      += cnzi;
119d20bfe6fSHong Zhang     current_space->local_remaining -= cnzi;
120d20bfe6fSHong Zhang 
121*2205254eSKarl Rupp     for (j=0;j<ptanzi;j++) ptadenserow[ptasparserow[j]] = 0;
122*2205254eSKarl Rupp 
123d20bfe6fSHong Zhang     /* Aside: Perhaps we should save the pta info for the numerical factorization. */
124d20bfe6fSHong Zhang     /*        For now, we will recompute what is needed. */
125d20bfe6fSHong Zhang     ci[i+1] = ci[i] + cnzi;
126d20bfe6fSHong Zhang   }
127d20bfe6fSHong Zhang   /* nnz is now stored in ci[ptm], column indices are in the list of free space */
128d20bfe6fSHong Zhang   /* Allocate space for cj, initialize cj, and */
129d20bfe6fSHong Zhang   /* destroy list of free space and other temporary array(s) */
13055a3bba9SHong Zhang   ierr = PetscMalloc((ci[pn]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
131a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
132d20bfe6fSHong Zhang   ierr = PetscFree(ptadenserow);CHKERRQ(ierr);
13355a3bba9SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
134d20bfe6fSHong Zhang 
13553565b12SHong Zhang   /* Allocate space for ca */
13653565b12SHong Zhang   ierr = PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
13753565b12SHong Zhang   ierr = PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));CHKERRQ(ierr);
13853565b12SHong Zhang 
13953565b12SHong Zhang   /* put together the new matrix */
14053565b12SHong Zhang   ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,pn,pn,ci,cj,ca,C);CHKERRQ(ierr);
141*2205254eSKarl Rupp 
142a2f3521dSMark F. Adams   (*C)->rmap->bs = P->cmap->bs;
143a2f3521dSMark F. Adams   (*C)->cmap->bs = P->cmap->bs;
144ae32dd66SHong Zhang 
14553565b12SHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
14653565b12SHong Zhang   /* Since these are PETSc arrays, change flags to free them as necessary. */
14753565b12SHong Zhang   c          = (Mat_SeqAIJ*)((*C)->data);
14853565b12SHong Zhang   c->free_a  = PETSC_TRUE;
14953565b12SHong Zhang   c->free_ij = PETSC_TRUE;
15053565b12SHong Zhang   c->nonew   = 0;
151*2205254eSKarl Rupp 
152f4a1e0b0SHong Zhang   (*C)->ops->ptapnumeric = MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy;
15353565b12SHong Zhang 
15453565b12SHong Zhang   /* Clean up. */
15553565b12SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
15653565b12SHong Zhang #if defined(PETSC_USE_INFO)
15753565b12SHong Zhang   if (ci[pn] != 0) {
15853565b12SHong Zhang     PetscReal afill = ((PetscReal)ci[pn])/ai[am];
15953565b12SHong Zhang     if (afill < 1.0) afill = 1.0;
16053565b12SHong Zhang     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
16153565b12SHong Zhang     ierr = PetscInfo1((*C),"Use MatPtAP(A,P,MatReuse,%G,&C) for best performance.\n",afill);CHKERRQ(ierr);
16253565b12SHong Zhang   } else {
16353565b12SHong Zhang     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
16453565b12SHong Zhang   }
165f79958ebSHong Zhang #endif
16653565b12SHong Zhang   PetscFunctionReturn(0);
16753565b12SHong Zhang }
16853565b12SHong Zhang 
16953565b12SHong Zhang #undef __FUNCT__
170ae32dd66SHong Zhang #define __FUNCT__ "MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy"
171ae32dd66SHong Zhang PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat A,Mat P,Mat C)
17253565b12SHong Zhang {
17353565b12SHong Zhang   PetscErrorCode ierr;
17453565b12SHong Zhang   Mat_SeqAIJ     *a = (Mat_SeqAIJ*) A->data;
17553565b12SHong Zhang   Mat_SeqAIJ     *p = (Mat_SeqAIJ*) P->data;
17653565b12SHong Zhang   Mat_SeqAIJ     *c = (Mat_SeqAIJ*) C->data;
17753565b12SHong Zhang   PetscInt       *ai=a->i,*aj=a->j,*apj,*apjdense,*pi=p->i,*pj=p->j,*pJ=p->j,*pjj;
17853565b12SHong Zhang   PetscInt       *ci=c->i,*cj=c->j,*cjj;
17953565b12SHong Zhang   PetscInt       am =A->rmap->N,cn=C->cmap->N,cm=C->rmap->N;
18053565b12SHong Zhang   PetscInt       i,j,k,anzi,pnzi,apnzj,nextap,pnzj,prow,crow;
18153565b12SHong Zhang   MatScalar      *aa=a->a,*apa,*pa=p->a,*pA=p->a,*paj,*ca=c->a,*caj;
18253565b12SHong Zhang 
18353565b12SHong Zhang   PetscFunctionBegin;
184ae32dd66SHong Zhang   /* Allocate temporary array for storage of one row of A*P (cn: non-scalable) */
185db1ebd87SHong Zhang   ierr = PetscMalloc(cn*(sizeof(MatScalar)+sizeof(PetscInt))+c->rmax*sizeof(PetscInt),&apa);CHKERRQ(ierr);
18653565b12SHong Zhang 
187db1ebd87SHong Zhang   apjdense = (PetscInt*)(apa + cn);
188db1ebd87SHong Zhang   apj      = apjdense + cn;
189db1ebd87SHong Zhang   ierr     = PetscMemzero(apa,cn*(sizeof(MatScalar)+sizeof(PetscInt)));CHKERRQ(ierr);
19053565b12SHong Zhang 
19153565b12SHong Zhang   /* Clear old values in C */
19253565b12SHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
19353565b12SHong Zhang 
19453565b12SHong Zhang   for (i=0; i<am; i++) {
19553565b12SHong Zhang     /* Form sparse row of A*P */
19653565b12SHong Zhang     anzi  = ai[i+1] - ai[i];
19753565b12SHong Zhang     apnzj = 0;
19853565b12SHong Zhang     for (j=0; j<anzi; j++) {
19953565b12SHong Zhang       prow = *aj++;
20053565b12SHong Zhang       pnzj = pi[prow+1] - pi[prow];
20153565b12SHong Zhang       pjj  = pj + pi[prow];
20253565b12SHong Zhang       paj  = pa + pi[prow];
20353565b12SHong Zhang       for (k=0; k<pnzj; k++) {
20453565b12SHong Zhang         if (!apjdense[pjj[k]]) {
20553565b12SHong Zhang           apjdense[pjj[k]] = -1;
20653565b12SHong Zhang           apj[apnzj++]     = pjj[k];
20753565b12SHong Zhang         }
20853565b12SHong Zhang         apa[pjj[k]] += (*aa)*paj[k];
20953565b12SHong Zhang       }
21053565b12SHong Zhang       ierr = PetscLogFlops(2.0*pnzj);CHKERRQ(ierr);
21153565b12SHong Zhang       aa++;
21253565b12SHong Zhang     }
21353565b12SHong Zhang 
21453565b12SHong Zhang     /* Sort the j index array for quick sparse axpy. */
21553565b12SHong Zhang     /* Note: a array does not need sorting as it is in dense storage locations. */
21653565b12SHong Zhang     ierr = PetscSortInt(apnzj,apj);CHKERRQ(ierr);
21753565b12SHong Zhang 
21853565b12SHong Zhang     /* Compute P^T*A*P using outer product (P^T)[:,j]*(A*P)[j,:]. */
21953565b12SHong Zhang     pnzi = pi[i+1] - pi[i];
22053565b12SHong Zhang     for (j=0; j<pnzi; j++) {
22153565b12SHong Zhang       nextap = 0;
22253565b12SHong Zhang       crow   = *pJ++;
22353565b12SHong Zhang       cjj    = cj + ci[crow];
22453565b12SHong Zhang       caj    = ca + ci[crow];
22553565b12SHong Zhang       /* Perform sparse axpy operation.  Note cjj includes apj. */
22653565b12SHong Zhang       for (k=0; nextap<apnzj; k++) {
22753565b12SHong Zhang #if defined(PETSC_USE_DEBUG)
228f23aa3ddSBarry Smith         if (k >= ci[crow+1] - ci[crow]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"k too large k %d, crow %d",k,crow);
22953565b12SHong Zhang #endif
23053565b12SHong Zhang         if (cjj[k]==apj[nextap]) {
23153565b12SHong Zhang           caj[k] += (*pA)*apa[apj[nextap++]];
23253565b12SHong Zhang         }
23353565b12SHong Zhang       }
23453565b12SHong Zhang       ierr = PetscLogFlops(2.0*apnzj);CHKERRQ(ierr);
23553565b12SHong Zhang       pA++;
23653565b12SHong Zhang     }
23753565b12SHong Zhang 
23853565b12SHong Zhang     /* Zero the current row info for A*P */
23953565b12SHong Zhang     for (j=0; j<apnzj; j++) {
24053565b12SHong Zhang       apa[apj[j]]      = 0.;
24153565b12SHong Zhang       apjdense[apj[j]] = 0;
24253565b12SHong Zhang     }
24353565b12SHong Zhang   }
24453565b12SHong Zhang 
24553565b12SHong Zhang   /* Assemble the final matrix and clean up */
24653565b12SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
24753565b12SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
248ae32dd66SHong Zhang 
24953565b12SHong Zhang   ierr = PetscFree(apa);CHKERRQ(ierr);
25053565b12SHong Zhang   PetscFunctionReturn(0);
25153565b12SHong Zhang }
25253565b12SHong Zhang 
25353565b12SHong Zhang #undef __FUNCT__
25453565b12SHong Zhang #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqAIJ"
25553565b12SHong Zhang PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
25653565b12SHong Zhang {
25753565b12SHong Zhang   PetscErrorCode ierr;
258afff960fSHong Zhang   Mat_SeqAIJ     *ap,*c;
259ae32dd66SHong Zhang   PetscInt       *api,*apj,*ci,pn=P->cmap->N;
26053565b12SHong Zhang   MatScalar      *ca;
26153565b12SHong Zhang   Mat_PtAP       *ptap;
262971236abSHong Zhang   Mat            Pt,AP;
263ae32dd66SHong Zhang   PetscBool      sparse_axpy=PETSC_TRUE;
264b8f276daSHong Zhang 
26553565b12SHong Zhang   PetscFunctionBegin;
266ae32dd66SHong Zhang   ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
26753565b12SHong Zhang   /* flag 'sparse_axpy' determines which implementations to be used:
268ae32dd66SHong Zhang        0: do dense axpy in MatPtAPNumeric() - fastest, but requires storage of struct A*P;
269ae32dd66SHong Zhang        1: do two sparse axpy in MatPtAPNumeric() - slowest, does not store structure of A*P. */
270ae32dd66SHong Zhang   ierr = PetscOptionsBool("-matptap_scalable","Use sparse axpy but slower MatPtAPNumeric()","",sparse_axpy,&sparse_axpy,PETSC_NULL);CHKERRQ(ierr);
271ae32dd66SHong Zhang   ierr = PetscOptionsEnd();CHKERRQ(ierr);
272ae32dd66SHong Zhang   if (sparse_axpy) {
273ae32dd66SHong Zhang     ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(A,P,fill,C);CHKERRQ(ierr);
27453565b12SHong Zhang     PetscFunctionReturn(0);
27553565b12SHong Zhang   }
27653565b12SHong Zhang 
277c0e3927dSHong Zhang   /* Get symbolic Pt = P^T */
278c0e3927dSHong Zhang   ierr = MatTransposeSymbolic_SeqAIJ(P,&Pt);CHKERRQ(ierr);
27953565b12SHong Zhang 
280c0e3927dSHong Zhang   /* Get symbolic AP = A*P */
281c0e3927dSHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,P,fill,&AP);CHKERRQ(ierr);
282a2f3521dSMark F. Adams 
283c0e3927dSHong Zhang   ap          = (Mat_SeqAIJ*)AP->data;
284c0e3927dSHong Zhang   api         = ap->i;
285c0e3927dSHong Zhang   apj         = ap->j;
286c0e3927dSHong Zhang   ap->free_ij = PETSC_FALSE; /* api and apj are kept in struct ptap, cannot be destroyed with AP */
28753565b12SHong Zhang 
288c0e3927dSHong Zhang   /* Get C = Pt*AP */
289c0e3927dSHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(Pt,AP,fill,C);CHKERRQ(ierr);
290a2f3521dSMark F. Adams 
291c0e3927dSHong Zhang   c         = (Mat_SeqAIJ*)(*C)->data;
292c0e3927dSHong Zhang   ci        = c->i;
293d20bfe6fSHong Zhang   ierr      = PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
294d20bfe6fSHong Zhang   ierr      = PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));CHKERRQ(ierr);
295c0e3927dSHong Zhang   c->a      = ca;
296e6b907acSBarry Smith   c->free_a = PETSC_TRUE;
297d20bfe6fSHong Zhang 
298c0e3927dSHong Zhang   /* Create a supporting struct for reuse by MatPtAPNumeric() */
2998d0a38e4SHong Zhang   ierr = PetscNew(Mat_PtAP,&ptap);CHKERRQ(ierr);
300*2205254eSKarl Rupp 
30153565b12SHong Zhang   c->ptap            = ptap;
3028d0a38e4SHong Zhang   ptap->destroy      = (*C)->ops->destroy;
3038d0a38e4SHong Zhang   (*C)->ops->destroy = MatDestroy_SeqAIJ_PtAP;
3048d0a38e4SHong Zhang 
3058d0a38e4SHong Zhang   /* Allocate temporary array for storage of one row of A*P */
306db1ebd87SHong Zhang   ierr = PetscMalloc((pn+1)*sizeof(PetscScalar),&ptap->apa);CHKERRQ(ierr);
307db1ebd87SHong Zhang   ierr = PetscMemzero(ptap->apa,(pn+1)*sizeof(PetscScalar));CHKERRQ(ierr);
308*2205254eSKarl Rupp 
309ae32dd66SHong Zhang   (*C)->ops->ptapnumeric = MatPtAPNumeric_SeqAIJ_SeqAIJ;
310ae32dd66SHong Zhang 
311f79958ebSHong Zhang   ptap->api = api;
312f79958ebSHong Zhang   ptap->apj = apj;
3138d0a38e4SHong Zhang 
314d20bfe6fSHong Zhang   /* Clean up. */
315971236abSHong Zhang   ierr = MatDestroy(&Pt);CHKERRQ(ierr);
316971236abSHong Zhang   ierr = MatDestroy(&AP);CHKERRQ(ierr);
317ae32dd66SHong Zhang #if defined(PETSC_USE_INFO)
318ae32dd66SHong Zhang   ierr = PetscInfo2((*C),"given fill %G, use scalable %d\n",fill,sparse_axpy);CHKERRQ(ierr);
319ae32dd66SHong Zhang #endif
320eb9c0419SKris Buschelman   PetscFunctionReturn(0);
321eb9c0419SKris Buschelman }
322eb9c0419SKris Buschelman 
323f2535c29SHong Zhang /* #define PROFILE_MatPtAPNumeric */
324eb9c0419SKris Buschelman #undef __FUNCT__
3259af31e4aSHong Zhang #define __FUNCT__ "MatPtAPNumeric_SeqAIJ_SeqAIJ"
326dfbe8321SBarry Smith PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat A,Mat P,Mat C)
327dfbe8321SBarry Smith {
328dfbe8321SBarry Smith   PetscErrorCode    ierr;
329d20bfe6fSHong Zhang   Mat_SeqAIJ        *a = (Mat_SeqAIJ*) A->data;
330f79958ebSHong Zhang   Mat_SeqAIJ        *p = (Mat_SeqAIJ*) P->data;
331d20bfe6fSHong Zhang   Mat_SeqAIJ        *c = (Mat_SeqAIJ*) C->data;
332b8f276daSHong Zhang   const PetscInt    *ai=a->i,*aj=a->j,*pi=p->i,*pj=p->j,*ci=c->i,*cj=c->j;
333b8f276daSHong Zhang   const PetscScalar *aa=a->a,*pa=p->a,*pval;
334b8f276daSHong Zhang   const PetscInt    *apj,*pcol,*cjj;
335b8f276daSHong Zhang   const PetscInt    am=A->rmap->N,cm=C->rmap->N;
336b8f276daSHong Zhang   PetscInt          i,j,k,anz,apnz,pnz,prow,crow,cnz;
337b8f276daSHong Zhang   PetscScalar       *apa,*ca=c->a,*caj,pvalj;
33853565b12SHong Zhang   Mat_PtAP          *ptap = c->ptap;
339f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric)
340f2535c29SHong Zhang   PetscLogDouble t0,tf,time_Cseq0=0.0,time_Cseq1=0.0;
341f2535c29SHong Zhang   PetscInt       flops0=0,flops1=0;
342f2535c29SHong Zhang #endif
343eb9c0419SKris Buschelman 
344eb9c0419SKris Buschelman   PetscFunctionBegin;
3458d0a38e4SHong Zhang   /* Get temporary array for storage of one row of A*P */
3468d0a38e4SHong Zhang   apa = ptap->apa;
347d20bfe6fSHong Zhang 
348d20bfe6fSHong Zhang   /* Clear old values in C */
349d20bfe6fSHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
350d20bfe6fSHong Zhang 
351d20bfe6fSHong Zhang   for (i=0; i<am; i++) {
35278844af4SHong Zhang     /* Form sparse row of AP[i,:] = A[i,:]*P */
353f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric)
354f2535c29SHong Zhang     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
355f2535c29SHong Zhang #endif
35678844af4SHong Zhang     anz  = ai[i+1] - ai[i];
35778844af4SHong Zhang     apnz = 0;
35878844af4SHong Zhang     for (j=0; j<anz; j++) {
35978844af4SHong Zhang       prow = aj[j];
36078844af4SHong Zhang       pnz  = pi[prow+1] - pi[prow];
36178844af4SHong Zhang       pcol = pj + pi[prow];
36278844af4SHong Zhang       pval = pa + pi[prow];
36378844af4SHong Zhang       for (k=0; k<pnz; k++) {
36478844af4SHong Zhang         apa[pcol[k]] += aa[j]*pval[k];
365d20bfe6fSHong Zhang       }
36678844af4SHong Zhang       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
367f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric)
368f2535c29SHong Zhang       flops0 += 2.0*pnz;
369f2535c29SHong Zhang #endif
37078844af4SHong Zhang     }
37178844af4SHong Zhang     aj += anz; aa += anz;
372f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric)
373f2535c29SHong Zhang     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
374*2205254eSKarl Rupp 
375f2535c29SHong Zhang     time_Cseq0 += tf - t0;
376f2535c29SHong Zhang #endif
37778844af4SHong Zhang 
37878844af4SHong Zhang     /* Compute P^T*A*P using outer product P[i,:]^T*AP[i,:]. */
379f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric)
380f2535c29SHong Zhang     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
381f2535c29SHong Zhang #endif
382f79958ebSHong Zhang     apj  = ptap->apj + ptap->api[i];
383f79958ebSHong Zhang     apnz = ptap->api[i+1] - ptap->api[i];
38478844af4SHong Zhang     pnz  = pi[i+1] - pi[i];
38578844af4SHong Zhang     pcol = pj + pi[i];
38678844af4SHong Zhang     pval = pa + pi[i];
3878d0a38e4SHong Zhang 
38853565b12SHong Zhang     /* Perform dense axpy */
38953565b12SHong Zhang     for (j=0; j<pnz; j++) {
39053565b12SHong Zhang       crow  = pcol[j];
39153565b12SHong Zhang       cjj   = cj + ci[crow];
39253565b12SHong Zhang       caj   = ca + ci[crow];
393f2535c29SHong Zhang       pvalj = pval[j];
39453565b12SHong Zhang       cnz   = ci[crow+1] - ci[crow];
395*2205254eSKarl Rupp       for (k=0; k<cnz; k++) caj[k] += pvalj*apa[cjj[k]];
39653565b12SHong Zhang       ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
397f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric)
398f2535c29SHong Zhang       flops1 += 2.0*cnz;
399f2535c29SHong Zhang #endif
40053565b12SHong Zhang     }
401f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric)
402f2535c29SHong Zhang     ierr        = PetscGetTime(&tf);CHKERRQ(ierr);
403f2535c29SHong Zhang     time_Cseq1 += tf - t0;
404f2535c29SHong Zhang #endif
40553565b12SHong Zhang 
40653565b12SHong Zhang     /* Zero the current row info for A*P */
4071d633602SHong Zhang     for (j=0; j<apnz; j++) apa[apj[j]] = 0.0;
40853565b12SHong Zhang   }
40953565b12SHong Zhang 
41053565b12SHong Zhang   /* Assemble the final matrix and clean up */
41153565b12SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
41253565b12SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
413f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric)
414f2535c29SHong Zhang   printf("PtAPNumeric_SeqAIJ time %g + %g, flops %d %d\n",time_Cseq0,time_Cseq1,flops0,flops1);
415f2535c29SHong Zhang #endif
41653565b12SHong Zhang   PetscFunctionReturn(0);
41753565b12SHong Zhang }
418