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> 108563dfccSBarry 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; 174a1b09b7SHong Zhang const char *algTypes[2] = {"scalable","nonscalable"}; 184a1b09b7SHong Zhang PetscInt alg=0; /* set default algorithm */ 197ba1a0bfSKris Buschelman 207ba1a0bfSKris Buschelman PetscFunctionBegin; 2165e8a0caSHong Zhang if (scall == MAT_INITIAL_MATRIX) { 224a1b09b7SHong Zhang /* 234a1b09b7SHong Zhang Alg 'scalable' determines which implementations to be used: 244a1b09b7SHong Zhang "nonscalable": do dense axpy in MatPtAPNumeric() - fastest, but requires storage of struct A*P; 254a1b09b7SHong Zhang "scalable": do two sparse axpy in MatPtAPNumeric() - might slow, does not store structure of A*P. 264a1b09b7SHong Zhang */ 274a1b09b7SHong Zhang ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr); 284a1b09b7SHong Zhang ierr = PetscOptionsEList("-matptap_via","Algorithmic approach","MatPtAP",algTypes,2,algTypes[0],&alg,NULL);CHKERRQ(ierr); 294a1b09b7SHong Zhang ierr = PetscOptionsEnd();CHKERRQ(ierr); 303ff4c91cSHong Zhang ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 314a1b09b7SHong Zhang switch (alg) { 324a1b09b7SHong Zhang case 1: 334a1b09b7SHong Zhang ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ_DenseAxpy(A,P,fill,C);CHKERRQ(ierr); 344a1b09b7SHong Zhang break; 354a1b09b7SHong Zhang default: 364a1b09b7SHong Zhang ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(A,P,fill,C);CHKERRQ(ierr); 374a1b09b7SHong Zhang break; 384a1b09b7SHong Zhang } 393ff4c91cSHong Zhang ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 407ba1a0bfSKris Buschelman } 413ff4c91cSHong Zhang ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 42ae32dd66SHong Zhang ierr = (*(*C)->ops->ptapnumeric)(A,P,*C);CHKERRQ(ierr); 433ff4c91cSHong Zhang ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 447ba1a0bfSKris Buschelman PetscFunctionReturn(0); 457ba1a0bfSKris Buschelman } 467ba1a0bfSKris Buschelman 477ba1a0bfSKris Buschelman #undef __FUNCT__ 4853565b12SHong Zhang #define __FUNCT__ "MatDestroy_SeqAIJ_PtAP" 4953565b12SHong Zhang PetscErrorCode MatDestroy_SeqAIJ_PtAP(Mat A) 508d0a38e4SHong Zhang { 518d0a38e4SHong Zhang PetscErrorCode ierr; 5253565b12SHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 5353565b12SHong Zhang Mat_PtAP *ptap = a->ptap; 548d0a38e4SHong Zhang 558d0a38e4SHong Zhang PetscFunctionBegin; 568d0a38e4SHong Zhang ierr = PetscFree(ptap->apa);CHKERRQ(ierr); 57f79958ebSHong Zhang ierr = PetscFree(ptap->api);CHKERRQ(ierr); 58f79958ebSHong Zhang ierr = PetscFree(ptap->apj);CHKERRQ(ierr); 5953565b12SHong Zhang ierr = (ptap->destroy)(A);CHKERRQ(ierr); 608d0a38e4SHong Zhang ierr = PetscFree(ptap);CHKERRQ(ierr); 618d0a38e4SHong Zhang PetscFunctionReturn(0); 628d0a38e4SHong Zhang } 638d0a38e4SHong Zhang 648d0a38e4SHong Zhang #undef __FUNCT__ 65ae32dd66SHong Zhang #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy" 66ae32dd66SHong Zhang PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(Mat A,Mat P,PetscReal fill,Mat *C) 678d0a38e4SHong Zhang { 688d0a38e4SHong Zhang PetscErrorCode ierr; 690298fd71SBarry Smith PetscFreeSpaceList free_space=NULL,current_space=NULL; 70d20bfe6fSHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*p = (Mat_SeqAIJ*)P->data,*c; 7153565b12SHong Zhang PetscInt *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj; 7253565b12SHong Zhang PetscInt *ci,*cj,*ptadenserow,*ptasparserow,*ptaj,nspacedouble=0; 737d69fa63SHong Zhang PetscInt an=A->cmap->N,am=A->rmap->N,pn=P->cmap->N,pm=P->rmap->N; 7453565b12SHong Zhang PetscInt i,j,k,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi,nlnk,*lnk; 75d20bfe6fSHong Zhang MatScalar *ca; 7653565b12SHong Zhang PetscBT lnkbt; 777d69fa63SHong Zhang PetscReal afill; 78eb9c0419SKris Buschelman 79eb9c0419SKris Buschelman PetscFunctionBegin; 8053565b12SHong Zhang /* Get ij structure of P^T */ 81eb9c0419SKris Buschelman ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr); 82d20bfe6fSHong Zhang ptJ = ptj; 83eb9c0419SKris Buschelman 84d20bfe6fSHong Zhang /* Allocate ci array, arrays for fill computation and */ 85d20bfe6fSHong Zhang /* free space for accumulating nonzero column info */ 86785e854fSJed Brown ierr = PetscMalloc1((pn+1),&ci);CHKERRQ(ierr); 87d20bfe6fSHong Zhang ci[0] = 0; 88eb9c0419SKris Buschelman 89*1795a4d1SJed Brown ierr = PetscCalloc1(2*an+1,&ptadenserow);CHKERRQ(ierr); 90d20bfe6fSHong Zhang ptasparserow = ptadenserow + an; 9155a3bba9SHong Zhang 9255a3bba9SHong Zhang /* create and initialize a linked list */ 9355a3bba9SHong Zhang nlnk = pn+1; 9455a3bba9SHong Zhang ierr = PetscLLCreate(pn,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr); 95eb9c0419SKris Buschelman 967d69fa63SHong Zhang /* Set initial free space to be fill*(nnz(A)+ nnz(P)) */ 977d69fa63SHong Zhang ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+pi[pm])),&free_space);CHKERRQ(ierr); 98d20bfe6fSHong Zhang current_space = free_space; 99d20bfe6fSHong Zhang 100d20bfe6fSHong Zhang /* Determine symbolic info for each row of C: */ 101d20bfe6fSHong Zhang for (i=0; i<pn; i++) { 102d20bfe6fSHong Zhang ptnzi = pti[i+1] - pti[i]; 103d20bfe6fSHong Zhang ptanzi = 0; 104d20bfe6fSHong Zhang /* Determine symbolic row of PtA: */ 105d20bfe6fSHong Zhang for (j=0; j<ptnzi; j++) { 106d20bfe6fSHong Zhang arow = *ptJ++; 107d20bfe6fSHong Zhang anzj = ai[arow+1] - ai[arow]; 108d20bfe6fSHong Zhang ajj = aj + ai[arow]; 109d20bfe6fSHong Zhang for (k=0; k<anzj; k++) { 110d20bfe6fSHong Zhang if (!ptadenserow[ajj[k]]) { 111d20bfe6fSHong Zhang ptadenserow[ajj[k]] = -1; 112d20bfe6fSHong Zhang ptasparserow[ptanzi++] = ajj[k]; 113d20bfe6fSHong Zhang } 114d20bfe6fSHong Zhang } 115d20bfe6fSHong Zhang } 116d20bfe6fSHong Zhang /* Using symbolic info for row of PtA, determine symbolic info for row of C: */ 117d20bfe6fSHong Zhang ptaj = ptasparserow; 118d20bfe6fSHong Zhang cnzi = 0; 119d20bfe6fSHong Zhang for (j=0; j<ptanzi; j++) { 120d20bfe6fSHong Zhang prow = *ptaj++; 121d20bfe6fSHong Zhang pnzj = pi[prow+1] - pi[prow]; 122d20bfe6fSHong Zhang pjj = pj + pi[prow]; 12355a3bba9SHong Zhang /* add non-zero cols of P into the sorted linked list lnk */ 124dadf0e6bSHong Zhang ierr = PetscLLAddSorted(pnzj,pjj,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr); 12555a3bba9SHong Zhang cnzi += nlnk; 126d20bfe6fSHong Zhang } 127d20bfe6fSHong Zhang 128d20bfe6fSHong Zhang /* If free space is not available, make more free space */ 129d20bfe6fSHong Zhang /* Double the amount of total space in the list */ 130d20bfe6fSHong Zhang if (current_space->local_remaining<cnzi) { 1314238b7adSHong Zhang ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 132f2b054eeSHong Zhang nspacedouble++; 133d20bfe6fSHong Zhang } 134d20bfe6fSHong Zhang 135d20bfe6fSHong Zhang /* Copy data into free space, and zero out denserows */ 13655a3bba9SHong Zhang ierr = PetscLLClean(pn,pn,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 1372205254eSKarl Rupp 138d20bfe6fSHong Zhang current_space->array += cnzi; 139d20bfe6fSHong Zhang current_space->local_used += cnzi; 140d20bfe6fSHong Zhang current_space->local_remaining -= cnzi; 141d20bfe6fSHong Zhang 1422205254eSKarl Rupp for (j=0; j<ptanzi; j++) ptadenserow[ptasparserow[j]] = 0; 1432205254eSKarl Rupp 144d20bfe6fSHong Zhang /* Aside: Perhaps we should save the pta info for the numerical factorization. */ 145d20bfe6fSHong Zhang /* For now, we will recompute what is needed. */ 146d20bfe6fSHong Zhang ci[i+1] = ci[i] + cnzi; 147d20bfe6fSHong Zhang } 148d20bfe6fSHong Zhang /* nnz is now stored in ci[ptm], column indices are in the list of free space */ 149d20bfe6fSHong Zhang /* Allocate space for cj, initialize cj, and */ 150d20bfe6fSHong Zhang /* destroy list of free space and other temporary array(s) */ 151785e854fSJed Brown ierr = PetscMalloc1((ci[pn]+1),&cj);CHKERRQ(ierr); 152a1a86e44SBarry Smith ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 153d20bfe6fSHong Zhang ierr = PetscFree(ptadenserow);CHKERRQ(ierr); 15455a3bba9SHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 155d20bfe6fSHong Zhang 156*1795a4d1SJed Brown ierr = PetscCalloc1((ci[pn]+1),&ca);CHKERRQ(ierr); 15753565b12SHong Zhang 15853565b12SHong Zhang /* put together the new matrix */ 159ce94432eSBarry Smith ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),pn,pn,ci,cj,ca,C);CHKERRQ(ierr); 1602205254eSKarl Rupp 161a2f3521dSMark F. Adams (*C)->rmap->bs = P->cmap->bs; 162a2f3521dSMark F. Adams (*C)->cmap->bs = P->cmap->bs; 163ae32dd66SHong Zhang 16453565b12SHong Zhang /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 16553565b12SHong Zhang /* Since these are PETSc arrays, change flags to free them as necessary. */ 16653565b12SHong Zhang c = (Mat_SeqAIJ*)((*C)->data); 16753565b12SHong Zhang c->free_a = PETSC_TRUE; 16853565b12SHong Zhang c->free_ij = PETSC_TRUE; 16953565b12SHong Zhang c->nonew = 0; 170f4a1e0b0SHong Zhang (*C)->ops->ptapnumeric = MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy; 17153565b12SHong Zhang 1727d69fa63SHong Zhang /* set MatInfo */ 1737d69fa63SHong Zhang afill = (PetscReal)ci[pn]/(ai[am]+pi[pm] + 1.e-5); 1747d69fa63SHong Zhang if (afill < 1.0) afill = 1.0; 1757d69fa63SHong Zhang c->maxnz = ci[pn]; 1767d69fa63SHong Zhang c->nz = ci[pn]; 1777d69fa63SHong Zhang (*C)->info.mallocs = nspacedouble; 1787d69fa63SHong Zhang (*C)->info.fill_ratio_given = fill; 1797d69fa63SHong Zhang (*C)->info.fill_ratio_needed = afill; 1807d69fa63SHong Zhang 18153565b12SHong Zhang /* Clean up. */ 18253565b12SHong Zhang ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr); 18353565b12SHong Zhang #if defined(PETSC_USE_INFO) 18453565b12SHong Zhang if (ci[pn] != 0) { 18553565b12SHong Zhang ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr); 18653565b12SHong Zhang ierr = PetscInfo1((*C),"Use MatPtAP(A,P,MatReuse,%G,&C) for best performance.\n",afill);CHKERRQ(ierr); 18753565b12SHong Zhang } else { 18853565b12SHong Zhang ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr); 18953565b12SHong Zhang } 190f79958ebSHong Zhang #endif 19153565b12SHong Zhang PetscFunctionReturn(0); 19253565b12SHong Zhang } 19353565b12SHong Zhang 19453565b12SHong Zhang #undef __FUNCT__ 195ae32dd66SHong Zhang #define __FUNCT__ "MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy" 196ae32dd66SHong Zhang PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat A,Mat P,Mat C) 19753565b12SHong Zhang { 19853565b12SHong Zhang PetscErrorCode ierr; 19953565b12SHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ*) A->data; 20053565b12SHong Zhang Mat_SeqAIJ *p = (Mat_SeqAIJ*) P->data; 20153565b12SHong Zhang Mat_SeqAIJ *c = (Mat_SeqAIJ*) C->data; 20253565b12SHong Zhang PetscInt *ai=a->i,*aj=a->j,*apj,*apjdense,*pi=p->i,*pj=p->j,*pJ=p->j,*pjj; 20353565b12SHong Zhang PetscInt *ci=c->i,*cj=c->j,*cjj; 20453565b12SHong Zhang PetscInt am =A->rmap->N,cn=C->cmap->N,cm=C->rmap->N; 20553565b12SHong Zhang PetscInt i,j,k,anzi,pnzi,apnzj,nextap,pnzj,prow,crow; 20653565b12SHong Zhang MatScalar *aa=a->a,*apa,*pa=p->a,*pA=p->a,*paj,*ca=c->a,*caj; 20753565b12SHong Zhang 20853565b12SHong Zhang PetscFunctionBegin; 209ae32dd66SHong Zhang /* Allocate temporary array for storage of one row of A*P (cn: non-scalable) */ 210dcca6d9dSJed Brown ierr = PetscMalloc3(cn,&apa,cn,&apjdense,c->rmax,&apj);CHKERRQ(ierr); 211eb9baa12SBarry Smith ierr = PetscMemzero(apa,cn*sizeof(MatScalar));CHKERRQ(ierr); 212eb9baa12SBarry Smith ierr = PetscMemzero(apjdense,cn*sizeof(PetscInt));CHKERRQ(ierr); 21353565b12SHong Zhang 21453565b12SHong Zhang /* Clear old values in C */ 21553565b12SHong Zhang ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 21653565b12SHong Zhang 21753565b12SHong Zhang for (i=0; i<am; i++) { 21853565b12SHong Zhang /* Form sparse row of A*P */ 21953565b12SHong Zhang anzi = ai[i+1] - ai[i]; 22053565b12SHong Zhang apnzj = 0; 22153565b12SHong Zhang for (j=0; j<anzi; j++) { 22253565b12SHong Zhang prow = *aj++; 22353565b12SHong Zhang pnzj = pi[prow+1] - pi[prow]; 22453565b12SHong Zhang pjj = pj + pi[prow]; 22553565b12SHong Zhang paj = pa + pi[prow]; 22653565b12SHong Zhang for (k=0; k<pnzj; k++) { 22753565b12SHong Zhang if (!apjdense[pjj[k]]) { 22853565b12SHong Zhang apjdense[pjj[k]] = -1; 22953565b12SHong Zhang apj[apnzj++] = pjj[k]; 23053565b12SHong Zhang } 23153565b12SHong Zhang apa[pjj[k]] += (*aa)*paj[k]; 23253565b12SHong Zhang } 23353565b12SHong Zhang ierr = PetscLogFlops(2.0*pnzj);CHKERRQ(ierr); 23453565b12SHong Zhang aa++; 23553565b12SHong Zhang } 23653565b12SHong Zhang 23753565b12SHong Zhang /* Sort the j index array for quick sparse axpy. */ 23853565b12SHong Zhang /* Note: a array does not need sorting as it is in dense storage locations. */ 23953565b12SHong Zhang ierr = PetscSortInt(apnzj,apj);CHKERRQ(ierr); 24053565b12SHong Zhang 24153565b12SHong Zhang /* Compute P^T*A*P using outer product (P^T)[:,j]*(A*P)[j,:]. */ 24253565b12SHong Zhang pnzi = pi[i+1] - pi[i]; 24353565b12SHong Zhang for (j=0; j<pnzi; j++) { 24453565b12SHong Zhang nextap = 0; 24553565b12SHong Zhang crow = *pJ++; 24653565b12SHong Zhang cjj = cj + ci[crow]; 24753565b12SHong Zhang caj = ca + ci[crow]; 24853565b12SHong Zhang /* Perform sparse axpy operation. Note cjj includes apj. */ 24953565b12SHong Zhang for (k=0; nextap<apnzj; k++) { 25053565b12SHong Zhang #if defined(PETSC_USE_DEBUG) 251f23aa3ddSBarry Smith if (k >= ci[crow+1] - ci[crow]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"k too large k %d, crow %d",k,crow); 25253565b12SHong Zhang #endif 25353565b12SHong Zhang if (cjj[k]==apj[nextap]) { 25453565b12SHong Zhang caj[k] += (*pA)*apa[apj[nextap++]]; 25553565b12SHong Zhang } 25653565b12SHong Zhang } 25753565b12SHong Zhang ierr = PetscLogFlops(2.0*apnzj);CHKERRQ(ierr); 25853565b12SHong Zhang pA++; 25953565b12SHong Zhang } 26053565b12SHong Zhang 26153565b12SHong Zhang /* Zero the current row info for A*P */ 26253565b12SHong Zhang for (j=0; j<apnzj; j++) { 26353565b12SHong Zhang apa[apj[j]] = 0.; 26453565b12SHong Zhang apjdense[apj[j]] = 0; 26553565b12SHong Zhang } 26653565b12SHong Zhang } 26753565b12SHong Zhang 26853565b12SHong Zhang /* Assemble the final matrix and clean up */ 26953565b12SHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 27053565b12SHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 271ae32dd66SHong Zhang 272eb9baa12SBarry Smith ierr = PetscFree3(apa,apjdense,apj);CHKERRQ(ierr); 27353565b12SHong Zhang PetscFunctionReturn(0); 27453565b12SHong Zhang } 27553565b12SHong Zhang 27653565b12SHong Zhang #undef __FUNCT__ 2774a1b09b7SHong Zhang #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqAIJ_DenseAxpy" 2784a1b09b7SHong Zhang PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_DenseAxpy(Mat A,Mat P,PetscReal fill,Mat *C) 27953565b12SHong Zhang { 28053565b12SHong Zhang PetscErrorCode ierr; 281afff960fSHong Zhang Mat_SeqAIJ *ap,*c; 282ae32dd66SHong Zhang PetscInt *api,*apj,*ci,pn=P->cmap->N; 28353565b12SHong Zhang MatScalar *ca; 28453565b12SHong Zhang Mat_PtAP *ptap; 285971236abSHong Zhang Mat Pt,AP; 286b8f276daSHong Zhang 28753565b12SHong Zhang PetscFunctionBegin; 288c0e3927dSHong Zhang /* Get symbolic Pt = P^T */ 289c0e3927dSHong Zhang ierr = MatTransposeSymbolic_SeqAIJ(P,&Pt);CHKERRQ(ierr); 29053565b12SHong Zhang 291c0e3927dSHong Zhang /* Get symbolic AP = A*P */ 292c0e3927dSHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,P,fill,&AP);CHKERRQ(ierr); 293a2f3521dSMark F. Adams 294c0e3927dSHong Zhang ap = (Mat_SeqAIJ*)AP->data; 295c0e3927dSHong Zhang api = ap->i; 296c0e3927dSHong Zhang apj = ap->j; 297c0e3927dSHong Zhang ap->free_ij = PETSC_FALSE; /* api and apj are kept in struct ptap, cannot be destroyed with AP */ 29853565b12SHong Zhang 299c0e3927dSHong Zhang /* Get C = Pt*AP */ 300c0e3927dSHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(Pt,AP,fill,C);CHKERRQ(ierr); 301a2f3521dSMark F. Adams 302c0e3927dSHong Zhang c = (Mat_SeqAIJ*)(*C)->data; 303c0e3927dSHong Zhang ci = c->i; 304*1795a4d1SJed Brown ierr = PetscCalloc1(ci[pn]+1,&ca);CHKERRQ(ierr); 305c0e3927dSHong Zhang c->a = ca; 306e6b907acSBarry Smith c->free_a = PETSC_TRUE; 307d20bfe6fSHong Zhang 308c0e3927dSHong Zhang /* Create a supporting struct for reuse by MatPtAPNumeric() */ 3098d0a38e4SHong Zhang ierr = PetscNew(Mat_PtAP,&ptap);CHKERRQ(ierr); 3102205254eSKarl Rupp 31153565b12SHong Zhang c->ptap = ptap; 3128d0a38e4SHong Zhang ptap->destroy = (*C)->ops->destroy; 3138d0a38e4SHong Zhang (*C)->ops->destroy = MatDestroy_SeqAIJ_PtAP; 3148d0a38e4SHong Zhang 3158d0a38e4SHong Zhang /* Allocate temporary array for storage of one row of A*P */ 316*1795a4d1SJed Brown ierr = PetscCalloc1(pn+1,&ptap->apa);CHKERRQ(ierr); 3172205254eSKarl Rupp 318ae32dd66SHong Zhang (*C)->ops->ptapnumeric = MatPtAPNumeric_SeqAIJ_SeqAIJ; 319ae32dd66SHong Zhang 320f79958ebSHong Zhang ptap->api = api; 321f79958ebSHong Zhang ptap->apj = apj; 3228d0a38e4SHong Zhang 323d20bfe6fSHong Zhang /* Clean up. */ 324971236abSHong Zhang ierr = MatDestroy(&Pt);CHKERRQ(ierr); 325971236abSHong Zhang ierr = MatDestroy(&AP);CHKERRQ(ierr); 326ae32dd66SHong Zhang #if defined(PETSC_USE_INFO) 3274a1b09b7SHong Zhang ierr = PetscInfo1((*C),"given fill %G\n",fill);CHKERRQ(ierr); 328ae32dd66SHong Zhang #endif 329eb9c0419SKris Buschelman PetscFunctionReturn(0); 330eb9c0419SKris Buschelman } 331eb9c0419SKris Buschelman 332f2535c29SHong Zhang /* #define PROFILE_MatPtAPNumeric */ 333eb9c0419SKris Buschelman #undef __FUNCT__ 3349af31e4aSHong Zhang #define __FUNCT__ "MatPtAPNumeric_SeqAIJ_SeqAIJ" 335dfbe8321SBarry Smith PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat A,Mat P,Mat C) 336dfbe8321SBarry Smith { 337dfbe8321SBarry Smith PetscErrorCode ierr; 338d20bfe6fSHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ*) A->data; 339f79958ebSHong Zhang Mat_SeqAIJ *p = (Mat_SeqAIJ*) P->data; 340d20bfe6fSHong Zhang Mat_SeqAIJ *c = (Mat_SeqAIJ*) C->data; 341b8f276daSHong Zhang const PetscInt *ai=a->i,*aj=a->j,*pi=p->i,*pj=p->j,*ci=c->i,*cj=c->j; 342b8f276daSHong Zhang const PetscScalar *aa=a->a,*pa=p->a,*pval; 343b8f276daSHong Zhang const PetscInt *apj,*pcol,*cjj; 344b8f276daSHong Zhang const PetscInt am=A->rmap->N,cm=C->rmap->N; 345b8f276daSHong Zhang PetscInt i,j,k,anz,apnz,pnz,prow,crow,cnz; 346b8f276daSHong Zhang PetscScalar *apa,*ca=c->a,*caj,pvalj; 34753565b12SHong Zhang Mat_PtAP *ptap = c->ptap; 348f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric) 349f2535c29SHong Zhang PetscLogDouble t0,tf,time_Cseq0=0.0,time_Cseq1=0.0; 350f2535c29SHong Zhang PetscInt flops0=0,flops1=0; 351f2535c29SHong Zhang #endif 352eb9c0419SKris Buschelman 353eb9c0419SKris Buschelman PetscFunctionBegin; 3548d0a38e4SHong Zhang /* Get temporary array for storage of one row of A*P */ 3558d0a38e4SHong Zhang apa = ptap->apa; 356d20bfe6fSHong Zhang 357d20bfe6fSHong Zhang /* Clear old values in C */ 358d20bfe6fSHong Zhang ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 359d20bfe6fSHong Zhang 360d20bfe6fSHong Zhang for (i=0; i<am; i++) { 36178844af4SHong Zhang /* Form sparse row of AP[i,:] = A[i,:]*P */ 362f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric) 3638563dfccSBarry Smith ierr = PetscTime(&t0);CHKERRQ(ierr); 364f2535c29SHong Zhang #endif 36578844af4SHong Zhang anz = ai[i+1] - ai[i]; 36678844af4SHong Zhang apnz = 0; 36778844af4SHong Zhang for (j=0; j<anz; j++) { 36878844af4SHong Zhang prow = aj[j]; 36978844af4SHong Zhang pnz = pi[prow+1] - pi[prow]; 37078844af4SHong Zhang pcol = pj + pi[prow]; 37178844af4SHong Zhang pval = pa + pi[prow]; 37278844af4SHong Zhang for (k=0; k<pnz; k++) { 37378844af4SHong Zhang apa[pcol[k]] += aa[j]*pval[k]; 374d20bfe6fSHong Zhang } 37578844af4SHong Zhang ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 376f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric) 377f2535c29SHong Zhang flops0 += 2.0*pnz; 378f2535c29SHong Zhang #endif 37978844af4SHong Zhang } 38078844af4SHong Zhang aj += anz; aa += anz; 381f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric) 3828563dfccSBarry Smith ierr = PetscTime(&tf);CHKERRQ(ierr); 3832205254eSKarl Rupp 384f2535c29SHong Zhang time_Cseq0 += tf - t0; 385f2535c29SHong Zhang #endif 38678844af4SHong Zhang 38778844af4SHong Zhang /* Compute P^T*A*P using outer product P[i,:]^T*AP[i,:]. */ 388f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric) 3898563dfccSBarry Smith ierr = PetscTime(&t0);CHKERRQ(ierr); 390f2535c29SHong Zhang #endif 391f79958ebSHong Zhang apj = ptap->apj + ptap->api[i]; 392f79958ebSHong Zhang apnz = ptap->api[i+1] - ptap->api[i]; 39378844af4SHong Zhang pnz = pi[i+1] - pi[i]; 39478844af4SHong Zhang pcol = pj + pi[i]; 39578844af4SHong Zhang pval = pa + pi[i]; 3968d0a38e4SHong Zhang 39753565b12SHong Zhang /* Perform dense axpy */ 39853565b12SHong Zhang for (j=0; j<pnz; j++) { 39953565b12SHong Zhang crow = pcol[j]; 40053565b12SHong Zhang cjj = cj + ci[crow]; 40153565b12SHong Zhang caj = ca + ci[crow]; 402f2535c29SHong Zhang pvalj = pval[j]; 40353565b12SHong Zhang cnz = ci[crow+1] - ci[crow]; 4042205254eSKarl Rupp for (k=0; k<cnz; k++) caj[k] += pvalj*apa[cjj[k]]; 40553565b12SHong Zhang ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 406f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric) 407f2535c29SHong Zhang flops1 += 2.0*cnz; 408f2535c29SHong Zhang #endif 40953565b12SHong Zhang } 410f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric) 4118563dfccSBarry Smith ierr = PetscTime(&tf);CHKERRQ(ierr); 412f2535c29SHong Zhang time_Cseq1 += tf - t0; 413f2535c29SHong Zhang #endif 41453565b12SHong Zhang 41553565b12SHong Zhang /* Zero the current row info for A*P */ 4161d633602SHong Zhang for (j=0; j<apnz; j++) apa[apj[j]] = 0.0; 41753565b12SHong Zhang } 41853565b12SHong Zhang 41953565b12SHong Zhang /* Assemble the final matrix and clean up */ 42053565b12SHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 42153565b12SHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 422f2535c29SHong Zhang #if defined(PROFILE_MatPtAPNumeric) 423f2535c29SHong Zhang printf("PtAPNumeric_SeqAIJ time %g + %g, flops %d %d\n",time_Cseq0,time_Cseq1,flops0,flops1); 424f2535c29SHong Zhang #endif 42553565b12SHong Zhang PetscFunctionReturn(0); 42653565b12SHong Zhang } 427