1eb9c0419SKris Buschelman /* 242c57489SHong Zhang Defines projective product routines where A is a SeqAIJ matrix 3eb9c0419SKris Buschelman C = P^T * A * P 4eb9c0419SKris Buschelman */ 5eb9c0419SKris Buschelman 6231952e2SKris Buschelman #include "src/mat/impls/aij/seq/aij.h" /*I "petscmat.h" I*/ 7eb9c0419SKris Buschelman #include "src/mat/utils/freespace.h" 855a3bba9SHong Zhang #include "petscbt.h" 9eb9c0419SKris Buschelman 10eb9c0419SKris Buschelman 11ff134f7aSHong Zhang #undef __FUNCT__ 129af31e4aSHong Zhang #define __FUNCT__ "MatPtAP_SeqAIJ_SeqAIJ" 13dfbe8321SBarry Smith PetscErrorCode MatPtAP_SeqAIJ_SeqAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C) 149af31e4aSHong Zhang { 15dfbe8321SBarry Smith PetscErrorCode ierr; 16b1d57f15SBarry Smith 179af31e4aSHong Zhang PetscFunctionBegin; 189af31e4aSHong Zhang if (scall == MAT_INITIAL_MATRIX){ 19d20bfe6fSHong Zhang ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 209af31e4aSHong Zhang ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ(A,P,fill,C);CHKERRQ(ierr); 21d20bfe6fSHong Zhang ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 229af31e4aSHong Zhang } 23d20bfe6fSHong Zhang ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 249af31e4aSHong Zhang ierr = MatPtAPNumeric_SeqAIJ_SeqAIJ(A,P,*C);CHKERRQ(ierr); 25d20bfe6fSHong Zhang ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 269af31e4aSHong Zhang PetscFunctionReturn(0); 279af31e4aSHong Zhang } 289af31e4aSHong Zhang 29eb9c0419SKris Buschelman #undef __FUNCT__ 309af31e4aSHong Zhang #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqAIJ" 3155a3bba9SHong Zhang PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat P,PetscReal fill,Mat *C) 3255a3bba9SHong Zhang { 33dfbe8321SBarry Smith PetscErrorCode ierr; 34d20bfe6fSHong Zhang FreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 35d20bfe6fSHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*p = (Mat_SeqAIJ*)P->data,*c; 3655a3bba9SHong Zhang PetscInt *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj; 37b8374ebeSBarry Smith PetscInt *ci,*cj,*ptadenserow,*ptasparserow,*ptaj; 3855a3bba9SHong Zhang PetscInt an=A->N,am=A->M,pn=P->N,pm=P->M; 39b8374ebeSBarry Smith PetscInt i,j,k,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi,nlnk,*lnk; 40d20bfe6fSHong Zhang MatScalar *ca; 4155a3bba9SHong Zhang PetscBT lnkbt; 42eb9c0419SKris Buschelman 43eb9c0419SKris Buschelman PetscFunctionBegin; 44d20bfe6fSHong Zhang /* Get ij structure of P^T */ 45eb9c0419SKris Buschelman ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr); 46d20bfe6fSHong Zhang ptJ=ptj; 47eb9c0419SKris Buschelman 48d20bfe6fSHong Zhang /* Allocate ci array, arrays for fill computation and */ 49d20bfe6fSHong Zhang /* free space for accumulating nonzero column info */ 5055a3bba9SHong Zhang ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 51d20bfe6fSHong Zhang ci[0] = 0; 52eb9c0419SKris Buschelman 5355a3bba9SHong Zhang ierr = PetscMalloc((2*an+1)*sizeof(PetscInt),&ptadenserow);CHKERRQ(ierr); 5455a3bba9SHong Zhang ierr = PetscMemzero(ptadenserow,(2*an+1)*sizeof(PetscInt));CHKERRQ(ierr); 55d20bfe6fSHong Zhang ptasparserow = ptadenserow + an; 5655a3bba9SHong Zhang 5755a3bba9SHong Zhang /* create and initialize a linked list */ 5855a3bba9SHong Zhang nlnk = pn+1; 5955a3bba9SHong Zhang ierr = PetscLLCreate(pn,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr); 60eb9c0419SKris Buschelman 61d20bfe6fSHong Zhang /* Set initial free space to be nnz(A) scaled by aspect ratio of P. */ 62d20bfe6fSHong Zhang /* This should be reasonable if sparsity of PtAP is similar to that of A. */ 63d20bfe6fSHong Zhang ierr = GetMoreSpace((ai[am]/pm)*pn,&free_space); 64d20bfe6fSHong Zhang current_space = free_space; 65d20bfe6fSHong Zhang 66d20bfe6fSHong Zhang /* Determine symbolic info for each row of C: */ 67d20bfe6fSHong Zhang for (i=0;i<pn;i++) { 68d20bfe6fSHong Zhang ptnzi = pti[i+1] - pti[i]; 69d20bfe6fSHong Zhang ptanzi = 0; 70d20bfe6fSHong Zhang /* Determine symbolic row of PtA: */ 71d20bfe6fSHong Zhang for (j=0;j<ptnzi;j++) { 72d20bfe6fSHong Zhang arow = *ptJ++; 73d20bfe6fSHong Zhang anzj = ai[arow+1] - ai[arow]; 74d20bfe6fSHong Zhang ajj = aj + ai[arow]; 75d20bfe6fSHong Zhang for (k=0;k<anzj;k++) { 76d20bfe6fSHong Zhang if (!ptadenserow[ajj[k]]) { 77d20bfe6fSHong Zhang ptadenserow[ajj[k]] = -1; 78d20bfe6fSHong Zhang ptasparserow[ptanzi++] = ajj[k]; 79d20bfe6fSHong Zhang } 80d20bfe6fSHong Zhang } 81d20bfe6fSHong Zhang } 82d20bfe6fSHong Zhang /* Using symbolic info for row of PtA, determine symbolic info for row of C: */ 83d20bfe6fSHong Zhang ptaj = ptasparserow; 84d20bfe6fSHong Zhang cnzi = 0; 85d20bfe6fSHong Zhang for (j=0;j<ptanzi;j++) { 86d20bfe6fSHong Zhang prow = *ptaj++; 87d20bfe6fSHong Zhang pnzj = pi[prow+1] - pi[prow]; 88d20bfe6fSHong Zhang pjj = pj + pi[prow]; 8955a3bba9SHong Zhang /* add non-zero cols of P into the sorted linked list lnk */ 9055a3bba9SHong Zhang ierr = PetscLLAdd(pnzj,pjj,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr); 9155a3bba9SHong Zhang cnzi += nlnk; 92d20bfe6fSHong Zhang } 93d20bfe6fSHong Zhang 94d20bfe6fSHong Zhang /* If free space is not available, make more free space */ 95d20bfe6fSHong Zhang /* Double the amount of total space in the list */ 96d20bfe6fSHong Zhang if (current_space->local_remaining<cnzi) { 97d20bfe6fSHong Zhang ierr = GetMoreSpace(current_space->total_array_size,¤t_space);CHKERRQ(ierr); 98d20bfe6fSHong Zhang } 99d20bfe6fSHong Zhang 100d20bfe6fSHong Zhang /* Copy data into free space, and zero out denserows */ 10155a3bba9SHong Zhang ierr = PetscLLClean(pn,pn,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 102d20bfe6fSHong Zhang current_space->array += cnzi; 103d20bfe6fSHong Zhang current_space->local_used += cnzi; 104d20bfe6fSHong Zhang current_space->local_remaining -= cnzi; 105d20bfe6fSHong Zhang 106d20bfe6fSHong Zhang for (j=0;j<ptanzi;j++) { 107d20bfe6fSHong Zhang ptadenserow[ptasparserow[j]] = 0; 108d20bfe6fSHong Zhang } 109d20bfe6fSHong Zhang /* Aside: Perhaps we should save the pta info for the numerical factorization. */ 110d20bfe6fSHong Zhang /* For now, we will recompute what is needed. */ 111d20bfe6fSHong Zhang ci[i+1] = ci[i] + cnzi; 112d20bfe6fSHong Zhang } 113d20bfe6fSHong Zhang /* nnz is now stored in ci[ptm], column indices are in the list of free space */ 114d20bfe6fSHong Zhang /* Allocate space for cj, initialize cj, and */ 115d20bfe6fSHong Zhang /* destroy list of free space and other temporary array(s) */ 11655a3bba9SHong Zhang ierr = PetscMalloc((ci[pn]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr); 117d20bfe6fSHong Zhang ierr = MakeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 118d20bfe6fSHong Zhang ierr = PetscFree(ptadenserow);CHKERRQ(ierr); 11955a3bba9SHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 120d20bfe6fSHong Zhang 121d20bfe6fSHong Zhang /* Allocate space for ca */ 122d20bfe6fSHong Zhang ierr = PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 123d20bfe6fSHong Zhang ierr = PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));CHKERRQ(ierr); 124d20bfe6fSHong Zhang 125d20bfe6fSHong Zhang /* put together the new matrix */ 126d20bfe6fSHong Zhang ierr = MatCreateSeqAIJWithArrays(A->comm,pn,pn,ci,cj,ca,C);CHKERRQ(ierr); 127d20bfe6fSHong Zhang 128d20bfe6fSHong Zhang /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 129d20bfe6fSHong Zhang /* Since these are PETSc arrays, change flags to free them as necessary. */ 130d20bfe6fSHong Zhang c = (Mat_SeqAIJ *)((*C)->data); 131d20bfe6fSHong Zhang c->freedata = PETSC_TRUE; 132d20bfe6fSHong Zhang c->nonew = 0; 133d20bfe6fSHong Zhang 134d20bfe6fSHong Zhang /* Clean up. */ 135d20bfe6fSHong Zhang ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr); 136eb9c0419SKris Buschelman 137eb9c0419SKris Buschelman PetscFunctionReturn(0); 138eb9c0419SKris Buschelman } 139eb9c0419SKris Buschelman 1403985e5eaSKris Buschelman #include "src/mat/impls/maij/maij.h" 1413985e5eaSKris Buschelman #undef __FUNCT__ 1429af31e4aSHong Zhang #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqMAIJ" 14355a3bba9SHong Zhang PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqMAIJ(Mat A,Mat PP,Mat *C) 14455a3bba9SHong Zhang { 145*7e692031SKris Buschelman /* This routine requires testing -- but it's getting better. */ 146dfbe8321SBarry Smith PetscErrorCode ierr; 1473985e5eaSKris Buschelman FreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 1483985e5eaSKris Buschelman Mat_SeqMAIJ *pp=(Mat_SeqMAIJ*)PP->data; 1493985e5eaSKris Buschelman Mat P=pp->AIJ; 1503985e5eaSKris Buschelman Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*p=(Mat_SeqAIJ*)P->data,*c; 151b1d57f15SBarry Smith PetscInt *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj; 152*7e692031SKris Buschelman PetscInt *ci,*cj,*ptadenserow,*ptasparserow,*denserow,*sparserow,*ptaj; 153*7e692031SKris Buschelman PetscInt an=A->N,am=A->M,pn=P->N,pm=P->M,ppdof=pp->dof,cn; 154*7e692031SKris Buschelman PetscInt i,j,k,dof,pshift,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi; 1553985e5eaSKris Buschelman MatScalar *ca; 1563985e5eaSKris Buschelman 1573985e5eaSKris Buschelman PetscFunctionBegin; 1583985e5eaSKris Buschelman /* Start timer */ 1599af31e4aSHong Zhang ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,PP,0,0);CHKERRQ(ierr); 1603985e5eaSKris Buschelman 1613985e5eaSKris Buschelman /* Get ij structure of P^T */ 1623985e5eaSKris Buschelman ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr); 1633985e5eaSKris Buschelman 164*7e692031SKris Buschelman cn = pn*ppdof; 1653985e5eaSKris Buschelman /* Allocate ci array, arrays for fill computation and */ 1663985e5eaSKris Buschelman /* free space for accumulating nonzero column info */ 167*7e692031SKris Buschelman ierr = PetscMalloc((cn+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 1683985e5eaSKris Buschelman ci[0] = 0; 1693985e5eaSKris Buschelman 170*7e692031SKris Buschelman /* Work arrays for rows of P^T*A */ 171*7e692031SKris Buschelman ierr = PetscMalloc((2*cn+2*an+1)*sizeof(PetscInt),&ptadenserow);CHKERRQ(ierr); 172*7e692031SKris Buschelman ierr = PetscMemzero(ptadenserow,(2*cn+2*an+1)*sizeof(PetscInt));CHKERRQ(ierr); 1733985e5eaSKris Buschelman ptasparserow = ptadenserow + an; 1743985e5eaSKris Buschelman denserow = ptasparserow + an; 175*7e692031SKris Buschelman sparserow = denserow + cn; 1763985e5eaSKris Buschelman 1773985e5eaSKris Buschelman /* Set initial free space to be nnz(A) scaled by aspect ratio of P. */ 1783985e5eaSKris Buschelman /* This should be reasonable if sparsity of PtAP is similar to that of A. */ 179*7e692031SKris Buschelman /* Note, aspect ratio of P is the same as the aspect ratio of SeqAIJ inside P */ 1803985e5eaSKris Buschelman ierr = GetMoreSpace((ai[am]/pm)*pn,&free_space); 1813985e5eaSKris Buschelman current_space = free_space; 1823985e5eaSKris Buschelman 1833985e5eaSKris Buschelman /* Determine symbolic info for each row of C: */ 184*7e692031SKris Buschelman for (i=0;i<pn;i++) { 1853985e5eaSKris Buschelman ptnzi = pti[i+1] - pti[i]; 1863985e5eaSKris Buschelman ptanzi = 0; 1873985e5eaSKris Buschelman ptJ = ptj + pti[i]; 1883985e5eaSKris Buschelman for (dof=0;dof<ppdof;dof++) { 1893985e5eaSKris Buschelman /* Determine symbolic row of PtA: */ 1903985e5eaSKris Buschelman for (j=0;j<ptnzi;j++) { 191*7e692031SKris Buschelman /* Expand ptJ[j] by block size and shift by dof to get the right row of A */ 192*7e692031SKris Buschelman arow = ptJ[j]*ppdof + dof; 193*7e692031SKris Buschelman /* Nonzeros of P^T*A will be in same locations as any element of A in that row */ 1943985e5eaSKris Buschelman anzj = ai[arow+1] - ai[arow]; 1953985e5eaSKris Buschelman ajj = aj + ai[arow]; 1963985e5eaSKris Buschelman for (k=0;k<anzj;k++) { 1973985e5eaSKris Buschelman if (!ptadenserow[ajj[k]]) { 1983985e5eaSKris Buschelman ptadenserow[ajj[k]] = -1; 1993985e5eaSKris Buschelman ptasparserow[ptanzi++] = ajj[k]; 2003985e5eaSKris Buschelman } 2013985e5eaSKris Buschelman } 2023985e5eaSKris Buschelman } 2033985e5eaSKris Buschelman /* Using symbolic info for row of PtA, determine symbolic info for row of C: */ 2043985e5eaSKris Buschelman ptaj = ptasparserow; 2053985e5eaSKris Buschelman cnzi = 0; 2063985e5eaSKris Buschelman for (j=0;j<ptanzi;j++) { 207*7e692031SKris Buschelman /* Get offset within block of P */ 208*7e692031SKris Buschelman pshift = *ptaj%ppdof; 209*7e692031SKris Buschelman /* Get block row of P */ 210*7e692031SKris Buschelman prow = (*ptaj++)/ppdof; 211*7e692031SKris Buschelman /* P has same number of nonzeros per row as the compressed form */ 2123985e5eaSKris Buschelman pnzj = pi[prow+1] - pi[prow]; 2133985e5eaSKris Buschelman pjj = pj + pi[prow]; 2143985e5eaSKris Buschelman for (k=0;k<pnzj;k++) { 215*7e692031SKris Buschelman /* Locations in C are shifted by the offset within the block */ 216*7e692031SKris Buschelman if (!denserow[pjj[k]+pshift]) { 217*7e692031SKris Buschelman denserow[pjj[k]+pshift] = -1; 218*7e692031SKris Buschelman sparserow[cnzi++] = pjj[k]+pshift; 2193985e5eaSKris Buschelman } 2203985e5eaSKris Buschelman } 2213985e5eaSKris Buschelman } 2223985e5eaSKris Buschelman 2233985e5eaSKris Buschelman /* sort sparserow */ 2243985e5eaSKris Buschelman ierr = PetscSortInt(cnzi,sparserow);CHKERRQ(ierr); 2253985e5eaSKris Buschelman 2263985e5eaSKris Buschelman /* If free space is not available, make more free space */ 2273985e5eaSKris Buschelman /* Double the amount of total space in the list */ 2283985e5eaSKris Buschelman if (current_space->local_remaining<cnzi) { 2293985e5eaSKris Buschelman ierr = GetMoreSpace(current_space->total_array_size,¤t_space);CHKERRQ(ierr); 2303985e5eaSKris Buschelman } 2313985e5eaSKris Buschelman 2323985e5eaSKris Buschelman /* Copy data into free space, and zero out denserows */ 233b1d57f15SBarry Smith ierr = PetscMemcpy(current_space->array,sparserow,cnzi*sizeof(PetscInt));CHKERRQ(ierr); 2343985e5eaSKris Buschelman current_space->array += cnzi; 2353985e5eaSKris Buschelman current_space->local_used += cnzi; 2363985e5eaSKris Buschelman current_space->local_remaining -= cnzi; 2373985e5eaSKris Buschelman 2383985e5eaSKris Buschelman for (j=0;j<ptanzi;j++) { 2393985e5eaSKris Buschelman ptadenserow[ptasparserow[j]] = 0; 2403985e5eaSKris Buschelman } 2413985e5eaSKris Buschelman for (j=0;j<cnzi;j++) { 2423985e5eaSKris Buschelman denserow[sparserow[j]] = 0; 2433985e5eaSKris Buschelman } 2443985e5eaSKris Buschelman /* Aside: Perhaps we should save the pta info for the numerical factorization. */ 2453985e5eaSKris Buschelman /* For now, we will recompute what is needed. */ 246*7e692031SKris Buschelman ci[i*ppdof+1+dof] = ci[i*ppdof+dof] + cnzi; 2473985e5eaSKris Buschelman } 2483985e5eaSKris Buschelman } 2493985e5eaSKris Buschelman /* nnz is now stored in ci[ptm], column indices are in the list of free space */ 2503985e5eaSKris Buschelman /* Allocate space for cj, initialize cj, and */ 2513985e5eaSKris Buschelman /* destroy list of free space and other temporary array(s) */ 252*7e692031SKris Buschelman ierr = PetscMalloc((ci[cn]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr); 2533985e5eaSKris Buschelman ierr = MakeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 2543985e5eaSKris Buschelman ierr = PetscFree(ptadenserow);CHKERRQ(ierr); 2553985e5eaSKris Buschelman 2563985e5eaSKris Buschelman /* Allocate space for ca */ 257*7e692031SKris Buschelman ierr = PetscMalloc((ci[cn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 258*7e692031SKris Buschelman ierr = PetscMemzero(ca,(ci[cn]+1)*sizeof(MatScalar));CHKERRQ(ierr); 2593985e5eaSKris Buschelman 2603985e5eaSKris Buschelman /* put together the new matrix */ 261*7e692031SKris Buschelman ierr = MatCreateSeqAIJWithArrays(A->comm,cn,cn,ci,cj,ca,C);CHKERRQ(ierr); 2623985e5eaSKris Buschelman 2633985e5eaSKris Buschelman /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 2643985e5eaSKris Buschelman /* Since these are PETSc arrays, change flags to free them as necessary. */ 2653985e5eaSKris Buschelman c = (Mat_SeqAIJ *)((*C)->data); 2663985e5eaSKris Buschelman c->freedata = PETSC_TRUE; 2673985e5eaSKris Buschelman c->nonew = 0; 2683985e5eaSKris Buschelman 2693985e5eaSKris Buschelman /* Clean up. */ 2703985e5eaSKris Buschelman ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr); 2713985e5eaSKris Buschelman 2729af31e4aSHong Zhang ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,PP,0,0);CHKERRQ(ierr); 2733985e5eaSKris Buschelman PetscFunctionReturn(0); 2743985e5eaSKris Buschelman } 275eb9c0419SKris Buschelman 276eb9c0419SKris Buschelman #undef __FUNCT__ 2779af31e4aSHong Zhang #define __FUNCT__ "MatPtAPNumeric_SeqAIJ_SeqAIJ" 278dfbe8321SBarry Smith PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat A,Mat P,Mat C) 279dfbe8321SBarry Smith { 280dfbe8321SBarry Smith PetscErrorCode ierr; 281b1d57f15SBarry Smith PetscInt flops=0; 282d20bfe6fSHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 283d20bfe6fSHong Zhang Mat_SeqAIJ *p = (Mat_SeqAIJ *) P->data; 284d20bfe6fSHong Zhang Mat_SeqAIJ *c = (Mat_SeqAIJ *) C->data; 285b1d57f15SBarry Smith PetscInt *ai=a->i,*aj=a->j,*apj,*apjdense,*pi=p->i,*pj=p->j,*pJ=p->j,*pjj; 286b1d57f15SBarry Smith PetscInt *ci=c->i,*cj=c->j,*cjj; 287b1d57f15SBarry Smith PetscInt am=A->M,cn=C->N,cm=C->M; 288b1d57f15SBarry Smith PetscInt i,j,k,anzi,pnzi,apnzj,nextap,pnzj,prow,crow; 289d20bfe6fSHong Zhang MatScalar *aa=a->a,*apa,*pa=p->a,*pA=p->a,*paj,*ca=c->a,*caj; 290eb9c0419SKris Buschelman 291eb9c0419SKris Buschelman PetscFunctionBegin; 292d20bfe6fSHong Zhang /* Allocate temporary array for storage of one row of A*P */ 293b1d57f15SBarry Smith ierr = PetscMalloc(cn*(sizeof(MatScalar)+2*sizeof(PetscInt)),&apa);CHKERRQ(ierr); 294b1d57f15SBarry Smith ierr = PetscMemzero(apa,cn*(sizeof(MatScalar)+2*sizeof(PetscInt)));CHKERRQ(ierr); 295eb9c0419SKris Buschelman 296b1d57f15SBarry Smith apj = (PetscInt *)(apa + cn); 297d20bfe6fSHong Zhang apjdense = apj + cn; 298d20bfe6fSHong Zhang 299d20bfe6fSHong Zhang /* Clear old values in C */ 300d20bfe6fSHong Zhang ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 301d20bfe6fSHong Zhang 302d20bfe6fSHong Zhang for (i=0;i<am;i++) { 303d20bfe6fSHong Zhang /* Form sparse row of A*P */ 304d20bfe6fSHong Zhang anzi = ai[i+1] - ai[i]; 305d20bfe6fSHong Zhang apnzj = 0; 306d20bfe6fSHong Zhang for (j=0;j<anzi;j++) { 307d20bfe6fSHong Zhang prow = *aj++; 308d20bfe6fSHong Zhang pnzj = pi[prow+1] - pi[prow]; 309d20bfe6fSHong Zhang pjj = pj + pi[prow]; 310d20bfe6fSHong Zhang paj = pa + pi[prow]; 311d20bfe6fSHong Zhang for (k=0;k<pnzj;k++) { 312d20bfe6fSHong Zhang if (!apjdense[pjj[k]]) { 313d20bfe6fSHong Zhang apjdense[pjj[k]] = -1; 314d20bfe6fSHong Zhang apj[apnzj++] = pjj[k]; 315d20bfe6fSHong Zhang } 316d20bfe6fSHong Zhang apa[pjj[k]] += (*aa)*paj[k]; 317d20bfe6fSHong Zhang } 318d20bfe6fSHong Zhang flops += 2*pnzj; 319d20bfe6fSHong Zhang aa++; 320d20bfe6fSHong Zhang } 321d20bfe6fSHong Zhang 322d20bfe6fSHong Zhang /* Sort the j index array for quick sparse axpy. */ 323d20bfe6fSHong Zhang ierr = PetscSortInt(apnzj,apj);CHKERRQ(ierr); 324d20bfe6fSHong Zhang 325d20bfe6fSHong Zhang /* Compute P^T*A*P using outer product (P^T)[:,j]*(A*P)[j,:]. */ 326d20bfe6fSHong Zhang pnzi = pi[i+1] - pi[i]; 327d20bfe6fSHong Zhang for (j=0;j<pnzi;j++) { 328d20bfe6fSHong Zhang nextap = 0; 329d20bfe6fSHong Zhang crow = *pJ++; 330d20bfe6fSHong Zhang cjj = cj + ci[crow]; 331d20bfe6fSHong Zhang caj = ca + ci[crow]; 332d20bfe6fSHong Zhang /* Perform sparse axpy operation. Note cjj includes apj. */ 333d20bfe6fSHong Zhang for (k=0;nextap<apnzj;k++) { 334d20bfe6fSHong Zhang if (cjj[k]==apj[nextap]) { 335d20bfe6fSHong Zhang caj[k] += (*pA)*apa[apj[nextap++]]; 336d20bfe6fSHong Zhang } 337d20bfe6fSHong Zhang } 338d20bfe6fSHong Zhang flops += 2*apnzj; 339d20bfe6fSHong Zhang pA++; 340d20bfe6fSHong Zhang } 341d20bfe6fSHong Zhang 342d20bfe6fSHong Zhang /* Zero the current row info for A*P */ 343d20bfe6fSHong Zhang for (j=0;j<apnzj;j++) { 344d20bfe6fSHong Zhang apa[apj[j]] = 0.; 345d20bfe6fSHong Zhang apjdense[apj[j]] = 0; 346d20bfe6fSHong Zhang } 347d20bfe6fSHong Zhang } 348d20bfe6fSHong Zhang 349d20bfe6fSHong Zhang /* Assemble the final matrix and clean up */ 350d20bfe6fSHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 351d20bfe6fSHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 352d20bfe6fSHong Zhang ierr = PetscFree(apa);CHKERRQ(ierr); 353d20bfe6fSHong Zhang ierr = PetscLogFlops(flops);CHKERRQ(ierr); 354d20bfe6fSHong Zhang 355eb9c0419SKris Buschelman PetscFunctionReturn(0); 356eb9c0419SKris Buschelman } 357