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> 10*8563dfccSBarry Smith #include <petsctime.h> 11eb9c0419SKris Buschelman 12ff134f7aSHong Zhang #undef __FUNCT__ 1365e8a0caSHong Zhang #define __FUNCT__ "MatPtAP_SeqAIJ_SeqAIJ" 1465e8a0caSHong Zhang PetscErrorCode MatPtAP_SeqAIJ_SeqAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C) 157ba1a0bfSKris Buschelman { 167ba1a0bfSKris Buschelman PetscErrorCode ierr; 177ba1a0bfSKris Buschelman 187ba1a0bfSKris Buschelman PetscFunctionBegin; 1965e8a0caSHong Zhang if (scall == MAT_INITIAL_MATRIX) { 2065e8a0caSHong Zhang ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ(A,P,fill,C);CHKERRQ(ierr); 217ba1a0bfSKris Buschelman } 22ae32dd66SHong Zhang ierr = (*(*C)->ops->ptapnumeric)(A,P,*C);CHKERRQ(ierr); 237ba1a0bfSKris Buschelman PetscFunctionReturn(0); 247ba1a0bfSKris Buschelman } 257ba1a0bfSKris Buschelman 267ba1a0bfSKris Buschelman #undef __FUNCT__ 2753565b12SHong Zhang #define __FUNCT__ "MatDestroy_SeqAIJ_PtAP" 2853565b12SHong Zhang PetscErrorCode MatDestroy_SeqAIJ_PtAP(Mat A) 298d0a38e4SHong Zhang { 308d0a38e4SHong Zhang PetscErrorCode ierr; 3153565b12SHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 3253565b12SHong Zhang Mat_PtAP *ptap = a->ptap; 338d0a38e4SHong Zhang 348d0a38e4SHong Zhang PetscFunctionBegin; 358d0a38e4SHong Zhang ierr = PetscFree(ptap->apa);CHKERRQ(ierr); 36f79958ebSHong Zhang ierr = PetscFree(ptap->api);CHKERRQ(ierr); 37f79958ebSHong Zhang ierr = PetscFree(ptap->apj);CHKERRQ(ierr); 3853565b12SHong Zhang ierr = (ptap->destroy)(A);CHKERRQ(ierr); 398d0a38e4SHong Zhang ierr = PetscFree(ptap);CHKERRQ(ierr); 408d0a38e4SHong Zhang PetscFunctionReturn(0); 418d0a38e4SHong Zhang } 428d0a38e4SHong Zhang 438d0a38e4SHong Zhang #undef __FUNCT__ 44ae32dd66SHong Zhang #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy" 45ae32dd66SHong Zhang PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(Mat A,Mat P,PetscReal fill,Mat *C) 468d0a38e4SHong Zhang { 478d0a38e4SHong Zhang PetscErrorCode ierr; 480298fd71SBarry Smith PetscFreeSpaceList free_space=NULL,current_space=NULL; 49d20bfe6fSHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*p = (Mat_SeqAIJ*)P->data,*c; 5053565b12SHong Zhang PetscInt *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj; 5153565b12SHong Zhang PetscInt *ci,*cj,*ptadenserow,*ptasparserow,*ptaj,nspacedouble=0; 527d69fa63SHong Zhang PetscInt an=A->cmap->N,am=A->rmap->N,pn=P->cmap->N,pm=P->rmap->N; 5353565b12SHong Zhang PetscInt i,j,k,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi,nlnk,*lnk; 54d20bfe6fSHong Zhang MatScalar *ca; 5553565b12SHong Zhang PetscBT lnkbt; 567d69fa63SHong Zhang PetscReal afill; 57eb9c0419SKris Buschelman 58eb9c0419SKris Buschelman PetscFunctionBegin; 5953565b12SHong Zhang /* Get ij structure of P^T */ 60eb9c0419SKris Buschelman ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr); 61d20bfe6fSHong Zhang ptJ = ptj; 62eb9c0419SKris Buschelman 63d20bfe6fSHong Zhang /* Allocate ci array, arrays for fill computation and */ 64d20bfe6fSHong Zhang /* free space for accumulating nonzero column info */ 6555a3bba9SHong Zhang ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 66d20bfe6fSHong Zhang ci[0] = 0; 67eb9c0419SKris Buschelman 6855a3bba9SHong Zhang ierr = PetscMalloc((2*an+1)*sizeof(PetscInt),&ptadenserow);CHKERRQ(ierr); 6955a3bba9SHong Zhang ierr = PetscMemzero(ptadenserow,(2*an+1)*sizeof(PetscInt));CHKERRQ(ierr); 70d20bfe6fSHong Zhang ptasparserow = ptadenserow + an; 7155a3bba9SHong Zhang 7255a3bba9SHong Zhang /* create and initialize a linked list */ 7355a3bba9SHong Zhang nlnk = pn+1; 7455a3bba9SHong Zhang ierr = PetscLLCreate(pn,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr); 75eb9c0419SKris Buschelman 767d69fa63SHong Zhang /* Set initial free space to be fill*(nnz(A)+ nnz(P)) */ 777d69fa63SHong Zhang ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+pi[pm])),&free_space);CHKERRQ(ierr); 78d20bfe6fSHong Zhang current_space = free_space; 79d20bfe6fSHong Zhang 80d20bfe6fSHong Zhang /* Determine symbolic info for each row of C: */ 81d20bfe6fSHong Zhang for (i=0; i<pn; i++) { 82d20bfe6fSHong Zhang ptnzi = pti[i+1] - pti[i]; 83d20bfe6fSHong Zhang ptanzi = 0; 84d20bfe6fSHong Zhang /* Determine symbolic row of PtA: */ 85d20bfe6fSHong Zhang for (j=0; j<ptnzi; j++) { 86d20bfe6fSHong Zhang arow = *ptJ++; 87d20bfe6fSHong Zhang anzj = ai[arow+1] - ai[arow]; 88d20bfe6fSHong Zhang ajj = aj + ai[arow]; 89d20bfe6fSHong Zhang for (k=0; k<anzj; k++) { 90d20bfe6fSHong Zhang if (!ptadenserow[ajj[k]]) { 91d20bfe6fSHong Zhang ptadenserow[ajj[k]] = -1; 92d20bfe6fSHong Zhang ptasparserow[ptanzi++] = ajj[k]; 93d20bfe6fSHong Zhang } 94d20bfe6fSHong Zhang } 95d20bfe6fSHong Zhang } 96d20bfe6fSHong Zhang /* Using symbolic info for row of PtA, determine symbolic info for row of C: */ 97d20bfe6fSHong Zhang ptaj = ptasparserow; 98d20bfe6fSHong Zhang cnzi = 0; 99d20bfe6fSHong Zhang for (j=0; j<ptanzi; j++) { 100d20bfe6fSHong Zhang prow = *ptaj++; 101d20bfe6fSHong Zhang pnzj = pi[prow+1] - pi[prow]; 102d20bfe6fSHong Zhang pjj = pj + pi[prow]; 10355a3bba9SHong Zhang /* add non-zero cols of P into the sorted linked list lnk */ 104dadf0e6bSHong Zhang ierr = PetscLLAddSorted(pnzj,pjj,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr); 10555a3bba9SHong Zhang cnzi += nlnk; 106d20bfe6fSHong Zhang } 107d20bfe6fSHong Zhang 108d20bfe6fSHong Zhang /* If free space is not available, make more free space */ 109d20bfe6fSHong Zhang /* Double the amount of total space in the list */ 110d20bfe6fSHong Zhang if (current_space->local_remaining<cnzi) { 1114238b7adSHong Zhang ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 112f2b054eeSHong Zhang nspacedouble++; 113d20bfe6fSHong Zhang } 114d20bfe6fSHong Zhang 115d20bfe6fSHong Zhang /* Copy data into free space, and zero out denserows */ 11655a3bba9SHong Zhang ierr = PetscLLClean(pn,pn,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 1172205254eSKarl Rupp 118d20bfe6fSHong Zhang current_space->array += cnzi; 119d20bfe6fSHong Zhang current_space->local_used += cnzi; 120d20bfe6fSHong Zhang current_space->local_remaining -= cnzi; 121d20bfe6fSHong Zhang 1222205254eSKarl Rupp for (j=0; j<ptanzi; j++) ptadenserow[ptasparserow[j]] = 0; 1232205254eSKarl Rupp 124d20bfe6fSHong Zhang /* Aside: Perhaps we should save the pta info for the numerical factorization. */ 125d20bfe6fSHong Zhang /* For now, we will recompute what is needed. */ 126d20bfe6fSHong Zhang ci[i+1] = ci[i] + cnzi; 127d20bfe6fSHong Zhang } 128d20bfe6fSHong Zhang /* nnz is now stored in ci[ptm], column indices are in the list of free space */ 129d20bfe6fSHong Zhang /* Allocate space for cj, initialize cj, and */ 130d20bfe6fSHong Zhang /* destroy list of free space and other temporary array(s) */ 13155a3bba9SHong Zhang ierr = PetscMalloc((ci[pn]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr); 132a1a86e44SBarry Smith ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 133d20bfe6fSHong Zhang ierr = PetscFree(ptadenserow);CHKERRQ(ierr); 13455a3bba9SHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 135d20bfe6fSHong Zhang 13653565b12SHong Zhang /* Allocate space for ca */ 13753565b12SHong Zhang ierr = PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 13853565b12SHong Zhang ierr = PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));CHKERRQ(ierr); 13953565b12SHong Zhang 14053565b12SHong Zhang /* put together the new matrix */ 141ce94432eSBarry Smith ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),pn,pn,ci,cj,ca,C);CHKERRQ(ierr); 1422205254eSKarl Rupp 143a2f3521dSMark F. Adams (*C)->rmap->bs = P->cmap->bs; 144a2f3521dSMark F. Adams (*C)->cmap->bs = P->cmap->bs; 145ae32dd66SHong Zhang 14653565b12SHong Zhang /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 14753565b12SHong Zhang /* Since these are PETSc arrays, change flags to free them as necessary. */ 14853565b12SHong Zhang c = (Mat_SeqAIJ*)((*C)->data); 14953565b12SHong Zhang c->free_a = PETSC_TRUE; 15053565b12SHong Zhang c->free_ij = PETSC_TRUE; 15153565b12SHong Zhang c->nonew = 0; 152f4a1e0b0SHong Zhang (*C)->ops->ptapnumeric = MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy; 15353565b12SHong Zhang 1547d69fa63SHong Zhang /* set MatInfo */ 1557d69fa63SHong Zhang afill = (PetscReal)ci[pn]/(ai[am]+pi[pm] + 1.e-5); 1567d69fa63SHong Zhang if (afill < 1.0) afill = 1.0; 1577d69fa63SHong Zhang c->maxnz = ci[pn]; 1587d69fa63SHong Zhang c->nz = ci[pn]; 1597d69fa63SHong Zhang (*C)->info.mallocs = nspacedouble; 1607d69fa63SHong Zhang (*C)->info.fill_ratio_given = fill; 1617d69fa63SHong Zhang (*C)->info.fill_ratio_needed = afill; 1627d69fa63SHong Zhang 16353565b12SHong Zhang /* Clean up. */ 16453565b12SHong Zhang ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr); 16553565b12SHong Zhang #if defined(PETSC_USE_INFO) 16653565b12SHong Zhang if (ci[pn] != 0) { 16753565b12SHong Zhang 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__ 177ae32dd66SHong Zhang #define __FUNCT__ "MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy" 178ae32dd66SHong Zhang PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy(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; 191ae32dd66SHong Zhang /* Allocate temporary array for storage of one row of A*P (cn: non-scalable) */ 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) 235f23aa3ddSBarry Smith if (k >= ci[crow+1] - ci[crow]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"k too large k %d, crow %d",k,crow); 23653565b12SHong Zhang #endif 23753565b12SHong Zhang if (cjj[k]==apj[nextap]) { 23853565b12SHong Zhang caj[k] += (*pA)*apa[apj[nextap++]]; 23953565b12SHong Zhang } 24053565b12SHong Zhang } 24153565b12SHong Zhang ierr = PetscLogFlops(2.0*apnzj);CHKERRQ(ierr); 24253565b12SHong Zhang pA++; 24353565b12SHong Zhang } 24453565b12SHong Zhang 24553565b12SHong Zhang /* Zero the current row info for A*P */ 24653565b12SHong Zhang for (j=0; j<apnzj; j++) { 24753565b12SHong Zhang apa[apj[j]] = 0.; 24853565b12SHong Zhang apjdense[apj[j]] = 0; 24953565b12SHong Zhang } 25053565b12SHong Zhang } 25153565b12SHong Zhang 25253565b12SHong Zhang /* Assemble the final matrix and clean up */ 25353565b12SHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 25453565b12SHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 255ae32dd66SHong Zhang 25653565b12SHong Zhang ierr = PetscFree(apa);CHKERRQ(ierr); 25753565b12SHong Zhang PetscFunctionReturn(0); 25853565b12SHong Zhang } 25953565b12SHong Zhang 26053565b12SHong Zhang #undef __FUNCT__ 26153565b12SHong Zhang #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqAIJ" 26253565b12SHong Zhang PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat P,PetscReal fill,Mat *C) 26353565b12SHong Zhang { 26453565b12SHong Zhang PetscErrorCode ierr; 265afff960fSHong Zhang Mat_SeqAIJ *ap,*c; 266ae32dd66SHong Zhang PetscInt *api,*apj,*ci,pn=P->cmap->N; 26753565b12SHong Zhang MatScalar *ca; 26853565b12SHong Zhang Mat_PtAP *ptap; 269971236abSHong Zhang Mat Pt,AP; 270ae32dd66SHong Zhang PetscBool sparse_axpy=PETSC_TRUE; 271b8f276daSHong Zhang 27253565b12SHong Zhang PetscFunctionBegin; 273ae32dd66SHong Zhang ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr); 27453565b12SHong Zhang /* flag 'sparse_axpy' determines which implementations to be used: 275ae32dd66SHong Zhang 0: do dense axpy in MatPtAPNumeric() - fastest, but requires storage of struct A*P; 276ae32dd66SHong Zhang 1: do two sparse axpy in MatPtAPNumeric() - slowest, does not store structure of A*P. */ 2770298fd71SBarry Smith ierr = PetscOptionsBool("-matptap_scalable","Use sparse axpy but slower MatPtAPNumeric()","",sparse_axpy,&sparse_axpy,NULL);CHKERRQ(ierr); 278ae32dd66SHong Zhang ierr = PetscOptionsEnd();CHKERRQ(ierr); 279ae32dd66SHong Zhang if (sparse_axpy) { 280ae32dd66SHong Zhang ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(A,P,fill,C);CHKERRQ(ierr); 28153565b12SHong Zhang PetscFunctionReturn(0); 28253565b12SHong Zhang } 28353565b12SHong Zhang 284c0e3927dSHong Zhang /* Get symbolic Pt = P^T */ 285c0e3927dSHong Zhang ierr = MatTransposeSymbolic_SeqAIJ(P,&Pt);CHKERRQ(ierr); 28653565b12SHong Zhang 287c0e3927dSHong Zhang /* Get symbolic AP = A*P */ 288c0e3927dSHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,P,fill,&AP);CHKERRQ(ierr); 289a2f3521dSMark F. Adams 290c0e3927dSHong Zhang ap = (Mat_SeqAIJ*)AP->data; 291c0e3927dSHong Zhang api = ap->i; 292c0e3927dSHong Zhang apj = ap->j; 293c0e3927dSHong Zhang ap->free_ij = PETSC_FALSE; /* api and apj are kept in struct ptap, cannot be destroyed with AP */ 29453565b12SHong Zhang 295c0e3927dSHong Zhang /* Get C = Pt*AP */ 296c0e3927dSHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(Pt,AP,fill,C);CHKERRQ(ierr); 297a2f3521dSMark F. Adams 298c0e3927dSHong Zhang c = (Mat_SeqAIJ*)(*C)->data; 299c0e3927dSHong Zhang ci = c->i; 300d20bfe6fSHong Zhang ierr = PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 301d20bfe6fSHong Zhang ierr = PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));CHKERRQ(ierr); 302c0e3927dSHong Zhang c->a = ca; 303e6b907acSBarry Smith c->free_a = PETSC_TRUE; 304d20bfe6fSHong Zhang 305c0e3927dSHong Zhang /* Create a supporting struct for reuse by MatPtAPNumeric() */ 3068d0a38e4SHong Zhang ierr = PetscNew(Mat_PtAP,&ptap);CHKERRQ(ierr); 3072205254eSKarl Rupp 30853565b12SHong Zhang c->ptap = ptap; 3098d0a38e4SHong Zhang ptap->destroy = (*C)->ops->destroy; 3108d0a38e4SHong Zhang (*C)->ops->destroy = MatDestroy_SeqAIJ_PtAP; 3118d0a38e4SHong Zhang 3128d0a38e4SHong Zhang /* Allocate temporary array for storage of one row of A*P */ 313db1ebd87SHong Zhang ierr = PetscMalloc((pn+1)*sizeof(PetscScalar),&ptap->apa);CHKERRQ(ierr); 314db1ebd87SHong Zhang ierr = PetscMemzero(ptap->apa,(pn+1)*sizeof(PetscScalar));CHKERRQ(ierr); 3152205254eSKarl Rupp 316ae32dd66SHong Zhang (*C)->ops->ptapnumeric = MatPtAPNumeric_SeqAIJ_SeqAIJ; 317ae32dd66SHong Zhang 318f79958ebSHong Zhang ptap->api = api; 319f79958ebSHong Zhang ptap->apj = apj; 3208d0a38e4SHong Zhang 321d20bfe6fSHong Zhang /* Clean up. */ 322971236abSHong Zhang ierr = MatDestroy(&Pt);CHKERRQ(ierr); 323971236abSHong Zhang ierr = MatDestroy(&AP);CHKERRQ(ierr); 324ae32dd66SHong Zhang #if defined(PETSC_USE_INFO) 325ae32dd66SHong Zhang ierr = PetscInfo2((*C),"given fill %G, use scalable %d\n",fill,sparse_axpy);CHKERRQ(ierr); 326ae32dd66SHong Zhang #endif 327eb9c0419SKris Buschelman PetscFunctionReturn(0); 328eb9c0419SKris Buschelman } 329eb9c0419SKris Buschelman 330f2535c29SHong Zhang /* #define PROFILE_MatPtAPNumeric */ 331eb9c0419SKris Buschelman #undef __FUNCT__ 3329af31e4aSHong Zhang #define __FUNCT__ "MatPtAPNumeric_SeqAIJ_SeqAIJ" 333dfbe8321SBarry Smith PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat A,Mat P,Mat C) 334dfbe8321SBarry Smith { 335dfbe8321SBarry Smith PetscErrorCode ierr; 336d20bfe6fSHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ*) A->data; 337f79958ebSHong Zhang Mat_SeqAIJ *p = (Mat_SeqAIJ*) P->data; 338d20bfe6fSHong Zhang Mat_SeqAIJ *c = (Mat_SeqAIJ*) C->data; 339b8f276daSHong Zhang const PetscInt *ai=a->i,*aj=a->j,*pi=p->i,*pj=p->j,*ci=c->i,*cj=c->j; 340b8f276daSHong Zhang const PetscScalar *aa=a->a,*pa=p->a,*pval; 341b8f276daSHong Zhang const PetscInt *apj,*pcol,*cjj; 342b8f276daSHong Zhang const PetscInt am=A->rmap->N,cm=C->rmap->N; 343b8f276daSHong Zhang PetscInt i,j,k,anz,apnz,pnz,prow,crow,cnz; 344b8f276daSHong Zhang PetscScalar *apa,*ca=c->a,*caj,pvalj; 34553565b12SHong Zhang Mat_PtAP *ptap = c->ptap; 346f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric) 347f2535c29SHong Zhang PetscLogDouble t0,tf,time_Cseq0=0.0,time_Cseq1=0.0; 348f2535c29SHong Zhang PetscInt flops0=0,flops1=0; 349f2535c29SHong Zhang #endif 350eb9c0419SKris Buschelman 351eb9c0419SKris Buschelman PetscFunctionBegin; 3528d0a38e4SHong Zhang /* Get temporary array for storage of one row of A*P */ 3538d0a38e4SHong Zhang apa = ptap->apa; 354d20bfe6fSHong Zhang 355d20bfe6fSHong Zhang /* Clear old values in C */ 356d20bfe6fSHong Zhang ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 357d20bfe6fSHong Zhang 358d20bfe6fSHong Zhang for (i=0; i<am; i++) { 35978844af4SHong Zhang /* Form sparse row of AP[i,:] = A[i,:]*P */ 360f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric) 361*8563dfccSBarry Smith ierr = PetscTime(&t0);CHKERRQ(ierr); 362f2535c29SHong Zhang #endif 36378844af4SHong Zhang anz = ai[i+1] - ai[i]; 36478844af4SHong Zhang apnz = 0; 36578844af4SHong Zhang for (j=0; j<anz; j++) { 36678844af4SHong Zhang prow = aj[j]; 36778844af4SHong Zhang pnz = pi[prow+1] - pi[prow]; 36878844af4SHong Zhang pcol = pj + pi[prow]; 36978844af4SHong Zhang pval = pa + pi[prow]; 37078844af4SHong Zhang for (k=0; k<pnz; k++) { 37178844af4SHong Zhang apa[pcol[k]] += aa[j]*pval[k]; 372d20bfe6fSHong Zhang } 37378844af4SHong Zhang ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 374f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric) 375f2535c29SHong Zhang flops0 += 2.0*pnz; 376f2535c29SHong Zhang #endif 37778844af4SHong Zhang } 37878844af4SHong Zhang aj += anz; aa += anz; 379f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric) 380*8563dfccSBarry Smith ierr = PetscTime(&tf);CHKERRQ(ierr); 3812205254eSKarl Rupp 382f2535c29SHong Zhang time_Cseq0 += tf - t0; 383f2535c29SHong Zhang #endif 38478844af4SHong Zhang 38578844af4SHong Zhang /* Compute P^T*A*P using outer product P[i,:]^T*AP[i,:]. */ 386f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric) 387*8563dfccSBarry Smith ierr = PetscTime(&t0);CHKERRQ(ierr); 388f2535c29SHong Zhang #endif 389f79958ebSHong Zhang apj = ptap->apj + ptap->api[i]; 390f79958ebSHong Zhang apnz = ptap->api[i+1] - ptap->api[i]; 39178844af4SHong Zhang pnz = pi[i+1] - pi[i]; 39278844af4SHong Zhang pcol = pj + pi[i]; 39378844af4SHong Zhang pval = pa + pi[i]; 3948d0a38e4SHong Zhang 39553565b12SHong Zhang /* Perform dense axpy */ 39653565b12SHong Zhang for (j=0; j<pnz; j++) { 39753565b12SHong Zhang crow = pcol[j]; 39853565b12SHong Zhang cjj = cj + ci[crow]; 39953565b12SHong Zhang caj = ca + ci[crow]; 400f2535c29SHong Zhang pvalj = pval[j]; 40153565b12SHong Zhang cnz = ci[crow+1] - ci[crow]; 4022205254eSKarl Rupp for (k=0; k<cnz; k++) caj[k] += pvalj*apa[cjj[k]]; 40353565b12SHong Zhang ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 404f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric) 405f2535c29SHong Zhang flops1 += 2.0*cnz; 406f2535c29SHong Zhang #endif 40753565b12SHong Zhang } 408f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric) 409*8563dfccSBarry Smith ierr = PetscTime(&tf);CHKERRQ(ierr); 410f2535c29SHong Zhang time_Cseq1 += tf - t0; 411f2535c29SHong Zhang #endif 41253565b12SHong Zhang 41353565b12SHong Zhang /* Zero the current row info for A*P */ 4141d633602SHong Zhang for (j=0; j<apnz; j++) apa[apj[j]] = 0.0; 41553565b12SHong Zhang } 41653565b12SHong Zhang 41753565b12SHong Zhang /* Assemble the final matrix and clean up */ 41853565b12SHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 41953565b12SHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 420f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric) 421f2535c29SHong Zhang printf("PtAPNumeric_SeqAIJ time %g + %g, flops %d %d\n",time_Cseq0,time_Cseq1,flops0,flops1); 422f2535c29SHong Zhang #endif 42353565b12SHong Zhang PetscFunctionReturn(0); 42453565b12SHong Zhang } 425