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 = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 2065e8a0caSHong Zhang ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ(A,P,fill,C);CHKERRQ(ierr); 2165e8a0caSHong Zhang ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 227ba1a0bfSKris Buschelman } 2365e8a0caSHong Zhang ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 24ae32dd66SHong Zhang ierr = (*(*C)->ops->ptapnumeric)(A,P,*C);CHKERRQ(ierr); 2565e8a0caSHong Zhang ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 267ba1a0bfSKris Buschelman PetscFunctionReturn(0); 277ba1a0bfSKris Buschelman } 287ba1a0bfSKris Buschelman 297ba1a0bfSKris Buschelman #undef __FUNCT__ 3053565b12SHong Zhang #define __FUNCT__ "MatDestroy_SeqAIJ_PtAP" 3153565b12SHong Zhang PetscErrorCode MatDestroy_SeqAIJ_PtAP(Mat A) 328d0a38e4SHong Zhang { 338d0a38e4SHong Zhang PetscErrorCode ierr; 3453565b12SHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 3553565b12SHong Zhang Mat_PtAP *ptap = a->ptap; 368d0a38e4SHong Zhang 378d0a38e4SHong Zhang PetscFunctionBegin; 388d0a38e4SHong Zhang ierr = PetscFree(ptap->apa);CHKERRQ(ierr); 39f79958ebSHong Zhang ierr = PetscFree(ptap->api);CHKERRQ(ierr); 40f79958ebSHong Zhang ierr = PetscFree(ptap->apj);CHKERRQ(ierr); 4153565b12SHong Zhang ierr = (ptap->destroy)(A);CHKERRQ(ierr); 428d0a38e4SHong Zhang ierr = PetscFree(ptap);CHKERRQ(ierr); 438d0a38e4SHong Zhang PetscFunctionReturn(0); 448d0a38e4SHong Zhang } 458d0a38e4SHong Zhang 468d0a38e4SHong Zhang #undef __FUNCT__ 47ae32dd66SHong Zhang #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy" 48ae32dd66SHong Zhang PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(Mat A,Mat P,PetscReal fill,Mat *C) 498d0a38e4SHong Zhang { 508d0a38e4SHong Zhang PetscErrorCode ierr; 5153565b12SHong Zhang PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 52d20bfe6fSHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*p = (Mat_SeqAIJ*)P->data,*c; 5353565b12SHong Zhang PetscInt *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj; 5453565b12SHong Zhang PetscInt *ci,*cj,*ptadenserow,*ptasparserow,*ptaj,nspacedouble=0; 55d0f46423SBarry Smith PetscInt an=A->cmap->N,am=A->rmap->N,pn=P->cmap->N; 5653565b12SHong Zhang PetscInt i,j,k,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi,nlnk,*lnk; 57d20bfe6fSHong Zhang MatScalar *ca; 5853565b12SHong Zhang PetscBT lnkbt; 59eb9c0419SKris Buschelman 60eb9c0419SKris Buschelman PetscFunctionBegin; 6153565b12SHong Zhang /* Get ij structure of P^T */ 62eb9c0419SKris Buschelman ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr); 63d20bfe6fSHong Zhang ptJ=ptj; 64eb9c0419SKris Buschelman 65d20bfe6fSHong Zhang /* Allocate ci array, arrays for fill computation and */ 66d20bfe6fSHong Zhang /* free space for accumulating nonzero column info */ 6755a3bba9SHong Zhang ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 68d20bfe6fSHong Zhang ci[0] = 0; 69eb9c0419SKris Buschelman 7055a3bba9SHong Zhang ierr = PetscMalloc((2*an+1)*sizeof(PetscInt),&ptadenserow);CHKERRQ(ierr); 7155a3bba9SHong Zhang ierr = PetscMemzero(ptadenserow,(2*an+1)*sizeof(PetscInt));CHKERRQ(ierr); 72d20bfe6fSHong Zhang ptasparserow = ptadenserow + an; 7355a3bba9SHong Zhang 7455a3bba9SHong Zhang /* create and initialize a linked list */ 7555a3bba9SHong Zhang nlnk = pn+1; 7655a3bba9SHong Zhang ierr = PetscLLCreate(pn,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr); 77eb9c0419SKris Buschelman 78f2b054eeSHong Zhang /* Set initial free space to be fill*nnz(A). */ 79d20bfe6fSHong Zhang /* This should be reasonable if sparsity of PtAP is similar to that of A. */ 80*a2ea699eSBarry Smith ierr = PetscFreeSpaceGet((PetscInt)(fill*ai[am]),&free_space);CHKERRQ(ierr); 81d20bfe6fSHong Zhang current_space = free_space; 82d20bfe6fSHong Zhang 83d20bfe6fSHong Zhang /* Determine symbolic info for each row of C: */ 84d20bfe6fSHong Zhang for (i=0;i<pn;i++) { 85d20bfe6fSHong Zhang ptnzi = pti[i+1] - pti[i]; 86d20bfe6fSHong Zhang ptanzi = 0; 87d20bfe6fSHong Zhang /* Determine symbolic row of PtA: */ 88d20bfe6fSHong Zhang for (j=0;j<ptnzi;j++) { 89d20bfe6fSHong Zhang arow = *ptJ++; 90d20bfe6fSHong Zhang anzj = ai[arow+1] - ai[arow]; 91d20bfe6fSHong Zhang ajj = aj + ai[arow]; 92d20bfe6fSHong Zhang for (k=0;k<anzj;k++) { 93d20bfe6fSHong Zhang if (!ptadenserow[ajj[k]]) { 94d20bfe6fSHong Zhang ptadenserow[ajj[k]] = -1; 95d20bfe6fSHong Zhang ptasparserow[ptanzi++] = ajj[k]; 96d20bfe6fSHong Zhang } 97d20bfe6fSHong Zhang } 98d20bfe6fSHong Zhang } 99d20bfe6fSHong Zhang /* Using symbolic info for row of PtA, determine symbolic info for row of C: */ 100d20bfe6fSHong Zhang ptaj = ptasparserow; 101d20bfe6fSHong Zhang cnzi = 0; 102d20bfe6fSHong Zhang for (j=0;j<ptanzi;j++) { 103d20bfe6fSHong Zhang prow = *ptaj++; 104d20bfe6fSHong Zhang pnzj = pi[prow+1] - pi[prow]; 105d20bfe6fSHong Zhang pjj = pj + pi[prow]; 10655a3bba9SHong Zhang /* add non-zero cols of P into the sorted linked list lnk */ 107dadf0e6bSHong Zhang ierr = PetscLLAddSorted(pnzj,pjj,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr); 10855a3bba9SHong Zhang cnzi += nlnk; 109d20bfe6fSHong Zhang } 110d20bfe6fSHong Zhang 111d20bfe6fSHong Zhang /* If free space is not available, make more free space */ 112d20bfe6fSHong Zhang /* Double the amount of total space in the list */ 113d20bfe6fSHong Zhang if (current_space->local_remaining<cnzi) { 1144238b7adSHong Zhang ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 115f2b054eeSHong Zhang nspacedouble++; 116d20bfe6fSHong Zhang } 117d20bfe6fSHong Zhang 118d20bfe6fSHong Zhang /* Copy data into free space, and zero out denserows */ 11955a3bba9SHong Zhang ierr = PetscLLClean(pn,pn,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 120d20bfe6fSHong Zhang current_space->array += cnzi; 121d20bfe6fSHong Zhang current_space->local_used += cnzi; 122d20bfe6fSHong Zhang current_space->local_remaining -= cnzi; 123d20bfe6fSHong Zhang 124d20bfe6fSHong Zhang for (j=0;j<ptanzi;j++) { 125d20bfe6fSHong Zhang ptadenserow[ptasparserow[j]] = 0; 126d20bfe6fSHong Zhang } 127d20bfe6fSHong Zhang /* Aside: Perhaps we should save the pta info for the numerical factorization. */ 128d20bfe6fSHong Zhang /* For now, we will recompute what is needed. */ 129d20bfe6fSHong Zhang ci[i+1] = ci[i] + cnzi; 130d20bfe6fSHong Zhang } 131d20bfe6fSHong Zhang /* nnz is now stored in ci[ptm], column indices are in the list of free space */ 132d20bfe6fSHong Zhang /* Allocate space for cj, initialize cj, and */ 133d20bfe6fSHong Zhang /* destroy list of free space and other temporary array(s) */ 13455a3bba9SHong Zhang ierr = PetscMalloc((ci[pn]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr); 135a1a86e44SBarry Smith ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 136d20bfe6fSHong Zhang ierr = PetscFree(ptadenserow);CHKERRQ(ierr); 13755a3bba9SHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 138d20bfe6fSHong Zhang 13953565b12SHong Zhang /* Allocate space for ca */ 14053565b12SHong Zhang ierr = PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 14153565b12SHong Zhang ierr = PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));CHKERRQ(ierr); 14253565b12SHong Zhang 14353565b12SHong Zhang /* put together the new matrix */ 14453565b12SHong Zhang ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,pn,pn,ci,cj,ca,C);CHKERRQ(ierr); 145a2f3521dSMark F. Adams (*C)->rmap->bs = P->cmap->bs; 146a2f3521dSMark F. Adams (*C)->cmap->bs = P->cmap->bs; 147ae32dd66SHong Zhang 14853565b12SHong Zhang /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 14953565b12SHong Zhang /* Since these are PETSc arrays, change flags to free them as necessary. */ 15053565b12SHong Zhang c = (Mat_SeqAIJ *)((*C)->data); 15153565b12SHong Zhang c->free_a = PETSC_TRUE; 15253565b12SHong Zhang c->free_ij = PETSC_TRUE; 15353565b12SHong Zhang c->nonew = 0; 154f4a1e0b0SHong Zhang (*C)->ops->ptapnumeric = MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy; 15553565b12SHong Zhang 15653565b12SHong Zhang /* Clean up. */ 15753565b12SHong Zhang ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr); 15853565b12SHong Zhang #if defined(PETSC_USE_INFO) 15953565b12SHong Zhang if (ci[pn] != 0) { 16053565b12SHong Zhang PetscReal afill = ((PetscReal)ci[pn])/ai[am]; 16153565b12SHong Zhang if (afill < 1.0) afill = 1.0; 16253565b12SHong Zhang ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr); 16353565b12SHong Zhang ierr = PetscInfo1((*C),"Use MatPtAP(A,P,MatReuse,%G,&C) for best performance.\n",afill);CHKERRQ(ierr); 16453565b12SHong Zhang } else { 16553565b12SHong Zhang ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr); 16653565b12SHong Zhang } 167f79958ebSHong Zhang #endif 16853565b12SHong Zhang PetscFunctionReturn(0); 16953565b12SHong Zhang } 17053565b12SHong Zhang 17153565b12SHong Zhang #undef __FUNCT__ 172ae32dd66SHong Zhang #define __FUNCT__ "MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy" 173ae32dd66SHong Zhang PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat A,Mat P,Mat C) 17453565b12SHong Zhang { 17553565b12SHong Zhang PetscErrorCode ierr; 17653565b12SHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 17753565b12SHong Zhang Mat_SeqAIJ *p = (Mat_SeqAIJ *) P->data; 17853565b12SHong Zhang Mat_SeqAIJ *c = (Mat_SeqAIJ *) C->data; 17953565b12SHong Zhang PetscInt *ai=a->i,*aj=a->j,*apj,*apjdense,*pi=p->i,*pj=p->j,*pJ=p->j,*pjj; 18053565b12SHong Zhang PetscInt *ci=c->i,*cj=c->j,*cjj; 18153565b12SHong Zhang PetscInt am=A->rmap->N,cn=C->cmap->N,cm=C->rmap->N; 18253565b12SHong Zhang PetscInt i,j,k,anzi,pnzi,apnzj,nextap,pnzj,prow,crow; 18353565b12SHong Zhang MatScalar *aa=a->a,*apa,*pa=p->a,*pA=p->a,*paj,*ca=c->a,*caj; 18453565b12SHong Zhang 18553565b12SHong Zhang PetscFunctionBegin; 186ae32dd66SHong Zhang /* Allocate temporary array for storage of one row of A*P (cn: non-scalable) */ 187db1ebd87SHong Zhang ierr = PetscMalloc(cn*(sizeof(MatScalar)+sizeof(PetscInt))+c->rmax*sizeof(PetscInt),&apa);CHKERRQ(ierr); 18853565b12SHong Zhang 189db1ebd87SHong Zhang apjdense = (PetscInt *)(apa + cn); 190db1ebd87SHong Zhang apj = apjdense + cn; 191db1ebd87SHong Zhang ierr = PetscMemzero(apa,cn*(sizeof(MatScalar)+sizeof(PetscInt)));CHKERRQ(ierr); 19253565b12SHong Zhang 19353565b12SHong Zhang /* Clear old values in C */ 19453565b12SHong Zhang ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 19553565b12SHong Zhang 19653565b12SHong Zhang for (i=0; i<am; i++) { 19753565b12SHong Zhang /* Form sparse row of A*P */ 19853565b12SHong Zhang anzi = ai[i+1] - ai[i]; 19953565b12SHong Zhang apnzj = 0; 20053565b12SHong Zhang for (j=0; j<anzi; j++) { 20153565b12SHong Zhang prow = *aj++; 20253565b12SHong Zhang pnzj = pi[prow+1] - pi[prow]; 20353565b12SHong Zhang pjj = pj + pi[prow]; 20453565b12SHong Zhang paj = pa + pi[prow]; 20553565b12SHong Zhang for (k=0;k<pnzj;k++) { 20653565b12SHong Zhang if (!apjdense[pjj[k]]) { 20753565b12SHong Zhang apjdense[pjj[k]] = -1; 20853565b12SHong Zhang apj[apnzj++] = pjj[k]; 20953565b12SHong Zhang } 21053565b12SHong Zhang apa[pjj[k]] += (*aa)*paj[k]; 21153565b12SHong Zhang } 21253565b12SHong Zhang ierr = PetscLogFlops(2.0*pnzj);CHKERRQ(ierr); 21353565b12SHong Zhang aa++; 21453565b12SHong Zhang } 21553565b12SHong Zhang 21653565b12SHong Zhang /* Sort the j index array for quick sparse axpy. */ 21753565b12SHong Zhang /* Note: a array does not need sorting as it is in dense storage locations. */ 21853565b12SHong Zhang ierr = PetscSortInt(apnzj,apj);CHKERRQ(ierr); 21953565b12SHong Zhang 22053565b12SHong Zhang /* Compute P^T*A*P using outer product (P^T)[:,j]*(A*P)[j,:]. */ 22153565b12SHong Zhang pnzi = pi[i+1] - pi[i]; 22253565b12SHong Zhang for (j=0; j<pnzi; j++) { 22353565b12SHong Zhang nextap = 0; 22453565b12SHong Zhang crow = *pJ++; 22553565b12SHong Zhang cjj = cj + ci[crow]; 22653565b12SHong Zhang caj = ca + ci[crow]; 22753565b12SHong Zhang /* Perform sparse axpy operation. Note cjj includes apj. */ 22853565b12SHong Zhang for (k=0;nextap<apnzj;k++) { 22953565b12SHong Zhang #if defined(PETSC_USE_DEBUG) 23053565b12SHong Zhang if (k >= ci[crow+1] - ci[crow]) { 23153565b12SHong Zhang SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"k too large k %d, crow %d",k,crow); 23253565b12SHong Zhang } 23353565b12SHong Zhang #endif 23453565b12SHong Zhang if (cjj[k]==apj[nextap]) { 23553565b12SHong Zhang caj[k] += (*pA)*apa[apj[nextap++]]; 23653565b12SHong Zhang } 23753565b12SHong Zhang } 23853565b12SHong Zhang ierr = PetscLogFlops(2.0*apnzj);CHKERRQ(ierr); 23953565b12SHong Zhang pA++; 24053565b12SHong Zhang } 24153565b12SHong Zhang 24253565b12SHong Zhang /* Zero the current row info for A*P */ 24353565b12SHong Zhang for (j=0; j<apnzj; j++) { 24453565b12SHong Zhang apa[apj[j]] = 0.; 24553565b12SHong Zhang apjdense[apj[j]] = 0; 24653565b12SHong Zhang } 24753565b12SHong Zhang } 24853565b12SHong Zhang 24953565b12SHong Zhang /* Assemble the final matrix and clean up */ 25053565b12SHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 25153565b12SHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 252ae32dd66SHong Zhang 25353565b12SHong Zhang ierr = PetscFree(apa);CHKERRQ(ierr); 25453565b12SHong Zhang PetscFunctionReturn(0); 25553565b12SHong Zhang } 25653565b12SHong Zhang 25753565b12SHong Zhang #undef __FUNCT__ 25853565b12SHong Zhang #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqAIJ" 25953565b12SHong Zhang PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat P,PetscReal fill,Mat *C) 26053565b12SHong Zhang { 26153565b12SHong Zhang PetscErrorCode ierr; 262afff960fSHong Zhang Mat_SeqAIJ *ap,*c; 263ae32dd66SHong Zhang PetscInt *api,*apj,*ci,pn=P->cmap->N; 26453565b12SHong Zhang MatScalar *ca; 26553565b12SHong Zhang Mat_PtAP *ptap; 266971236abSHong Zhang Mat Pt,AP; 267ae32dd66SHong Zhang PetscBool sparse_axpy=PETSC_TRUE; 268b8f276daSHong Zhang 26953565b12SHong Zhang PetscFunctionBegin; 270ae32dd66SHong Zhang ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr); 27153565b12SHong Zhang /* flag 'sparse_axpy' determines which implementations to be used: 272ae32dd66SHong Zhang 0: do dense axpy in MatPtAPNumeric() - fastest, but requires storage of struct A*P; 273ae32dd66SHong Zhang 1: do two sparse axpy in MatPtAPNumeric() - slowest, does not store structure of A*P. */ 274ae32dd66SHong Zhang ierr = PetscOptionsBool("-matptap_scalable","Use sparse axpy but slower MatPtAPNumeric()","",sparse_axpy,&sparse_axpy,PETSC_NULL);CHKERRQ(ierr); 275ae32dd66SHong Zhang ierr = PetscOptionsEnd();CHKERRQ(ierr); 276ae32dd66SHong Zhang if (sparse_axpy){ 277ae32dd66SHong Zhang ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(A,P,fill,C);CHKERRQ(ierr); 27853565b12SHong Zhang PetscFunctionReturn(0); 27953565b12SHong Zhang } 28053565b12SHong Zhang 281c0e3927dSHong Zhang /* Get symbolic Pt = P^T */ 282c0e3927dSHong Zhang ierr = MatTransposeSymbolic_SeqAIJ(P,&Pt);CHKERRQ(ierr); 28353565b12SHong Zhang 284c0e3927dSHong Zhang /* Get symbolic AP = A*P */ 285c0e3927dSHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,P,fill,&AP);CHKERRQ(ierr); 286a2f3521dSMark F. Adams 287c0e3927dSHong Zhang ap = (Mat_SeqAIJ*)AP->data; 288c0e3927dSHong Zhang api = ap->i; 289c0e3927dSHong Zhang apj = ap->j; 290c0e3927dSHong Zhang ap->free_ij = PETSC_FALSE; /* api and apj are kept in struct ptap, cannot be destroyed with AP */ 29153565b12SHong Zhang 292c0e3927dSHong Zhang /* Get C = Pt*AP */ 293c0e3927dSHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(Pt,AP,fill,C);CHKERRQ(ierr); 294a2f3521dSMark F. Adams 295c0e3927dSHong Zhang c = (Mat_SeqAIJ*)(*C)->data; 296c0e3927dSHong Zhang ci = c->i; 297d20bfe6fSHong Zhang ierr = PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 298d20bfe6fSHong Zhang ierr = PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));CHKERRQ(ierr); 299c0e3927dSHong Zhang c->a = ca; 300e6b907acSBarry Smith c->free_a = PETSC_TRUE; 301d20bfe6fSHong Zhang 302c0e3927dSHong Zhang /* Create a supporting struct for reuse by MatPtAPNumeric() */ 3038d0a38e4SHong Zhang ierr = PetscNew(Mat_PtAP,&ptap);CHKERRQ(ierr); 30453565b12SHong Zhang c->ptap = ptap; 3058d0a38e4SHong Zhang ptap->destroy = (*C)->ops->destroy; 3068d0a38e4SHong Zhang (*C)->ops->destroy = MatDestroy_SeqAIJ_PtAP; 3078d0a38e4SHong Zhang 3088d0a38e4SHong Zhang /* Allocate temporary array for storage of one row of A*P */ 309db1ebd87SHong Zhang ierr = PetscMalloc((pn+1)*sizeof(PetscScalar),&ptap->apa);CHKERRQ(ierr); 310db1ebd87SHong Zhang ierr = PetscMemzero(ptap->apa,(pn+1)*sizeof(PetscScalar));CHKERRQ(ierr); 311ae32dd66SHong Zhang (*C)->ops->ptapnumeric = MatPtAPNumeric_SeqAIJ_SeqAIJ; 312ae32dd66SHong Zhang 313f79958ebSHong Zhang ptap->api = api; 314f79958ebSHong Zhang ptap->apj = apj; 3158d0a38e4SHong Zhang 316d20bfe6fSHong Zhang /* Clean up. */ 317971236abSHong Zhang ierr = MatDestroy(&Pt);CHKERRQ(ierr); 318971236abSHong Zhang ierr = MatDestroy(&AP);CHKERRQ(ierr); 319ae32dd66SHong Zhang #if defined(PETSC_USE_INFO) 320ae32dd66SHong Zhang ierr = PetscInfo2((*C),"given fill %G, use scalable %d\n",fill,sparse_axpy);CHKERRQ(ierr); 321ae32dd66SHong Zhang #endif 322eb9c0419SKris Buschelman PetscFunctionReturn(0); 323eb9c0419SKris Buschelman } 324eb9c0419SKris Buschelman 325f2535c29SHong Zhang /* #define PROFILE_MatPtAPNumeric */ 326eb9c0419SKris Buschelman #undef __FUNCT__ 3279af31e4aSHong Zhang #define __FUNCT__ "MatPtAPNumeric_SeqAIJ_SeqAIJ" 328dfbe8321SBarry Smith PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat A,Mat P,Mat C) 329dfbe8321SBarry Smith { 330dfbe8321SBarry Smith PetscErrorCode ierr; 331d20bfe6fSHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 332f79958ebSHong Zhang Mat_SeqAIJ *p = (Mat_SeqAIJ *) P->data; 333d20bfe6fSHong Zhang Mat_SeqAIJ *c = (Mat_SeqAIJ *) C->data; 334b8f276daSHong Zhang const PetscInt *ai=a->i,*aj=a->j,*pi=p->i,*pj=p->j,*ci=c->i,*cj=c->j; 335b8f276daSHong Zhang const PetscScalar *aa=a->a,*pa=p->a,*pval; 336b8f276daSHong Zhang const PetscInt *apj,*pcol,*cjj; 337b8f276daSHong Zhang const PetscInt am=A->rmap->N,cm=C->rmap->N; 338b8f276daSHong Zhang PetscInt i,j,k,anz,apnz,pnz,prow,crow,cnz; 339b8f276daSHong Zhang PetscScalar *apa,*ca=c->a,*caj,pvalj; 34053565b12SHong Zhang Mat_PtAP *ptap = c->ptap; 341f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric) 342f2535c29SHong Zhang PetscLogDouble t0,tf,time_Cseq0=0.0,time_Cseq1=0.0; 343f2535c29SHong Zhang PetscInt flops0=0,flops1=0; 344f2535c29SHong Zhang #endif 345eb9c0419SKris Buschelman 346eb9c0419SKris Buschelman PetscFunctionBegin; 3478d0a38e4SHong Zhang /* Get temporary array for storage of one row of A*P */ 3488d0a38e4SHong Zhang apa = ptap->apa; 349d20bfe6fSHong Zhang 350d20bfe6fSHong Zhang /* Clear old values in C */ 351d20bfe6fSHong Zhang ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 352d20bfe6fSHong Zhang 353d20bfe6fSHong Zhang for (i=0;i<am;i++) { 35478844af4SHong Zhang /* Form sparse row of AP[i,:] = A[i,:]*P */ 355f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric) 356f2535c29SHong Zhang ierr = PetscGetTime(&t0);CHKERRQ(ierr); 357f2535c29SHong Zhang #endif 35878844af4SHong Zhang anz = ai[i+1] - ai[i]; 35978844af4SHong Zhang apnz = 0; 36078844af4SHong Zhang for (j=0; j<anz; j++) { 36178844af4SHong Zhang prow = aj[j]; 36278844af4SHong Zhang pnz = pi[prow+1] - pi[prow]; 36378844af4SHong Zhang pcol = pj + pi[prow]; 36478844af4SHong Zhang pval = pa + pi[prow]; 36578844af4SHong Zhang for (k=0; k<pnz; k++) { 36678844af4SHong Zhang apa[pcol[k]] += aa[j]*pval[k]; 367d20bfe6fSHong Zhang } 36878844af4SHong Zhang ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 369f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric) 370f2535c29SHong Zhang flops0 += 2.0*pnz; 371f2535c29SHong Zhang #endif 37278844af4SHong Zhang } 37378844af4SHong Zhang aj += anz; aa += anz; 374f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric) 375f2535c29SHong Zhang ierr = PetscGetTime(&tf);CHKERRQ(ierr); 376f2535c29SHong Zhang time_Cseq0 += tf - t0; 377f2535c29SHong Zhang #endif 37878844af4SHong Zhang 37978844af4SHong Zhang /* Compute P^T*A*P using outer product P[i,:]^T*AP[i,:]. */ 380f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric) 381f2535c29SHong Zhang ierr = PetscGetTime(&t0);CHKERRQ(ierr); 382f2535c29SHong Zhang #endif 383f79958ebSHong Zhang apj = ptap->apj + ptap->api[i]; 384f79958ebSHong Zhang apnz = ptap->api[i+1] - ptap->api[i]; 38578844af4SHong Zhang pnz = pi[i+1] - pi[i]; 38678844af4SHong Zhang pcol = pj + pi[i]; 38778844af4SHong Zhang pval = pa + pi[i]; 3888d0a38e4SHong Zhang 38953565b12SHong Zhang /* Perform dense axpy */ 39053565b12SHong Zhang for (j=0; j<pnz; j++) { 39153565b12SHong Zhang crow = pcol[j]; 39253565b12SHong Zhang cjj = cj + ci[crow]; 39353565b12SHong Zhang caj = ca + ci[crow]; 394f2535c29SHong Zhang pvalj = pval[j]; 39553565b12SHong Zhang cnz = ci[crow+1] - ci[crow]; 39653565b12SHong Zhang for (k=0; k<cnz; k++){ 397f2535c29SHong Zhang caj[k] += pvalj*apa[cjj[k]]; 39853565b12SHong Zhang } 39953565b12SHong Zhang ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 400f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric) 401f2535c29SHong Zhang flops1 += 2.0*cnz; 402f2535c29SHong Zhang #endif 40353565b12SHong Zhang } 404f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric) 405f2535c29SHong Zhang ierr = PetscGetTime(&tf);CHKERRQ(ierr); 406f2535c29SHong Zhang time_Cseq1 += tf - t0; 407f2535c29SHong Zhang #endif 40853565b12SHong Zhang 40953565b12SHong Zhang /* Zero the current row info for A*P */ 4101d633602SHong Zhang for (j=0; j<apnz; j++) apa[apj[j]] = 0.0; 41153565b12SHong Zhang } 41253565b12SHong Zhang 41353565b12SHong Zhang /* Assemble the final matrix and clean up */ 41453565b12SHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 41553565b12SHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 416f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric) 417f2535c29SHong Zhang printf("PtAPNumeric_SeqAIJ time %g + %g, flops %d %d\n",time_Cseq0,time_Cseq1,flops0,flops1); 418f2535c29SHong Zhang #endif 41953565b12SHong Zhang PetscFunctionReturn(0); 42053565b12SHong Zhang } 421