1be1d678aSKris Buschelman #define PETSCMAT_DLL 2be1d678aSKris Buschelman 3eb9c0419SKris Buschelman /* 442c57489SHong Zhang Defines projective product routines where A is a SeqAIJ matrix 5eb9c0419SKris Buschelman C = P^T * A * P 6eb9c0419SKris Buschelman */ 7eb9c0419SKris Buschelman 8231952e2SKris Buschelman #include "src/mat/impls/aij/seq/aij.h" /*I "petscmat.h" I*/ 9eb9c0419SKris Buschelman #include "src/mat/utils/freespace.h" 1055a3bba9SHong Zhang #include "petscbt.h" 11eb9c0419SKris Buschelman 12ff134f7aSHong Zhang #undef __FUNCT__ 137ba1a0bfSKris Buschelman #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ" 147ba1a0bfSKris Buschelman PetscErrorCode MatPtAPSymbolic_SeqAIJ(Mat A,Mat P,PetscReal fill,Mat *C) 157ba1a0bfSKris Buschelman { 167ba1a0bfSKris Buschelman PetscErrorCode ierr; 177ba1a0bfSKris Buschelman 187ba1a0bfSKris Buschelman PetscFunctionBegin; 197ba1a0bfSKris Buschelman if (!P->ops->ptapsymbolic_seqaij) { 207ba1a0bfSKris Buschelman SETERRQ2(PETSC_ERR_SUP,"Not implemented for A=%s and P=%s",A->type_name,P->type_name); 217ba1a0bfSKris Buschelman } 227ba1a0bfSKris Buschelman ierr = (*P->ops->ptapsymbolic_seqaij)(A,P,fill,C);CHKERRQ(ierr); 237ba1a0bfSKris Buschelman PetscFunctionReturn(0); 247ba1a0bfSKris Buschelman } 257ba1a0bfSKris Buschelman 267ba1a0bfSKris Buschelman #undef __FUNCT__ 277ba1a0bfSKris Buschelman #define __FUNCT__ "MatPtAPNumeric_SeqAIJ" 287ba1a0bfSKris Buschelman PetscErrorCode MatPtAPNumeric_SeqAIJ(Mat A,Mat P,Mat C) 297ba1a0bfSKris Buschelman { 307ba1a0bfSKris Buschelman PetscErrorCode ierr; 317ba1a0bfSKris Buschelman 327ba1a0bfSKris Buschelman PetscFunctionBegin; 337ba1a0bfSKris Buschelman if (!P->ops->ptapnumeric_seqaij) { 347ba1a0bfSKris Buschelman SETERRQ2(PETSC_ERR_SUP,"Not implemented for A=%s and P=%s",A->type_name,P->type_name); 357ba1a0bfSKris Buschelman } 367ba1a0bfSKris Buschelman ierr = (*P->ops->ptapnumeric_seqaij)(A,P,C);CHKERRQ(ierr); 377ba1a0bfSKris Buschelman PetscFunctionReturn(0); 387ba1a0bfSKris Buschelman } 397ba1a0bfSKris Buschelman 407ba1a0bfSKris Buschelman #undef __FUNCT__ 419af31e4aSHong Zhang #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqAIJ" 4255a3bba9SHong Zhang PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat P,PetscReal fill,Mat *C) 4355a3bba9SHong Zhang { 44dfbe8321SBarry Smith PetscErrorCode ierr; 45a1a86e44SBarry Smith PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 46d20bfe6fSHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*p = (Mat_SeqAIJ*)P->data,*c; 4755a3bba9SHong Zhang PetscInt *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj; 48*f2b054eeSHong Zhang PetscInt *ci,*cj,*ptadenserow,*ptasparserow,*ptaj,nspacedouble=0; 49899cda47SBarry Smith PetscInt an=A->cmap.N,am=A->rmap.N,pn=P->cmap.N,pm=P->rmap.N; 50b8374ebeSBarry Smith PetscInt i,j,k,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi,nlnk,*lnk; 51d20bfe6fSHong Zhang MatScalar *ca; 5255a3bba9SHong Zhang PetscBT lnkbt; 53eb9c0419SKris Buschelman 54eb9c0419SKris Buschelman PetscFunctionBegin; 55d20bfe6fSHong Zhang /* Get ij structure of P^T */ 56eb9c0419SKris Buschelman ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr); 57d20bfe6fSHong Zhang ptJ=ptj; 58eb9c0419SKris Buschelman 59d20bfe6fSHong Zhang /* Allocate ci array, arrays for fill computation and */ 60d20bfe6fSHong Zhang /* free space for accumulating nonzero column info */ 6155a3bba9SHong Zhang ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 62d20bfe6fSHong Zhang ci[0] = 0; 63eb9c0419SKris Buschelman 6455a3bba9SHong Zhang ierr = PetscMalloc((2*an+1)*sizeof(PetscInt),&ptadenserow);CHKERRQ(ierr); 6555a3bba9SHong Zhang ierr = PetscMemzero(ptadenserow,(2*an+1)*sizeof(PetscInt));CHKERRQ(ierr); 66d20bfe6fSHong Zhang ptasparserow = ptadenserow + an; 6755a3bba9SHong Zhang 6855a3bba9SHong Zhang /* create and initialize a linked list */ 6955a3bba9SHong Zhang nlnk = pn+1; 7055a3bba9SHong Zhang ierr = PetscLLCreate(pn,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr); 71eb9c0419SKris Buschelman 72*f2b054eeSHong Zhang /* Set initial free space to be fill*nnz(A). */ 73d20bfe6fSHong Zhang /* This should be reasonable if sparsity of PtAP is similar to that of A. */ 74*f2b054eeSHong Zhang ierr = PetscFreeSpaceGet((PetscInt)(fill*ai[am]),&free_space); 75d20bfe6fSHong Zhang current_space = free_space; 76d20bfe6fSHong Zhang 77d20bfe6fSHong Zhang /* Determine symbolic info for each row of C: */ 78d20bfe6fSHong Zhang for (i=0;i<pn;i++) { 79d20bfe6fSHong Zhang ptnzi = pti[i+1] - pti[i]; 80d20bfe6fSHong Zhang ptanzi = 0; 81d20bfe6fSHong Zhang /* Determine symbolic row of PtA: */ 82d20bfe6fSHong Zhang for (j=0;j<ptnzi;j++) { 83d20bfe6fSHong Zhang arow = *ptJ++; 84d20bfe6fSHong Zhang anzj = ai[arow+1] - ai[arow]; 85d20bfe6fSHong Zhang ajj = aj + ai[arow]; 86d20bfe6fSHong Zhang for (k=0;k<anzj;k++) { 87d20bfe6fSHong Zhang if (!ptadenserow[ajj[k]]) { 88d20bfe6fSHong Zhang ptadenserow[ajj[k]] = -1; 89d20bfe6fSHong Zhang ptasparserow[ptanzi++] = ajj[k]; 90d20bfe6fSHong Zhang } 91d20bfe6fSHong Zhang } 92d20bfe6fSHong Zhang } 93d20bfe6fSHong Zhang /* Using symbolic info for row of PtA, determine symbolic info for row of C: */ 94d20bfe6fSHong Zhang ptaj = ptasparserow; 95d20bfe6fSHong Zhang cnzi = 0; 96d20bfe6fSHong Zhang for (j=0;j<ptanzi;j++) { 97d20bfe6fSHong Zhang prow = *ptaj++; 98d20bfe6fSHong Zhang pnzj = pi[prow+1] - pi[prow]; 99d20bfe6fSHong Zhang pjj = pj + pi[prow]; 10055a3bba9SHong Zhang /* add non-zero cols of P into the sorted linked list lnk */ 10155a3bba9SHong Zhang ierr = PetscLLAdd(pnzj,pjj,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr); 10255a3bba9SHong Zhang cnzi += nlnk; 103d20bfe6fSHong Zhang } 104d20bfe6fSHong Zhang 105d20bfe6fSHong Zhang /* If free space is not available, make more free space */ 106d20bfe6fSHong Zhang /* Double the amount of total space in the list */ 107d20bfe6fSHong Zhang if (current_space->local_remaining<cnzi) { 108a1a86e44SBarry Smith ierr = PetscFreeSpaceGet(current_space->total_array_size,¤t_space);CHKERRQ(ierr); 109*f2b054eeSHong Zhang nspacedouble++; 110d20bfe6fSHong Zhang } 111d20bfe6fSHong Zhang 112d20bfe6fSHong Zhang /* Copy data into free space, and zero out denserows */ 11355a3bba9SHong Zhang ierr = PetscLLClean(pn,pn,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 114d20bfe6fSHong Zhang current_space->array += cnzi; 115d20bfe6fSHong Zhang current_space->local_used += cnzi; 116d20bfe6fSHong Zhang current_space->local_remaining -= cnzi; 117d20bfe6fSHong Zhang 118d20bfe6fSHong Zhang for (j=0;j<ptanzi;j++) { 119d20bfe6fSHong Zhang ptadenserow[ptasparserow[j]] = 0; 120d20bfe6fSHong Zhang } 121d20bfe6fSHong Zhang /* Aside: Perhaps we should save the pta info for the numerical factorization. */ 122d20bfe6fSHong Zhang /* For now, we will recompute what is needed. */ 123d20bfe6fSHong Zhang ci[i+1] = ci[i] + cnzi; 124d20bfe6fSHong Zhang } 125d20bfe6fSHong Zhang /* nnz is now stored in ci[ptm], column indices are in the list of free space */ 126d20bfe6fSHong Zhang /* Allocate space for cj, initialize cj, and */ 127d20bfe6fSHong Zhang /* destroy list of free space and other temporary array(s) */ 12855a3bba9SHong Zhang ierr = PetscMalloc((ci[pn]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr); 129a1a86e44SBarry Smith ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 130d20bfe6fSHong Zhang ierr = PetscFree(ptadenserow);CHKERRQ(ierr); 13155a3bba9SHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 132d20bfe6fSHong Zhang 133d20bfe6fSHong Zhang /* Allocate space for ca */ 134d20bfe6fSHong Zhang ierr = PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 135d20bfe6fSHong Zhang ierr = PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));CHKERRQ(ierr); 136d20bfe6fSHong Zhang 137d20bfe6fSHong Zhang /* put together the new matrix */ 138d20bfe6fSHong Zhang ierr = MatCreateSeqAIJWithArrays(A->comm,pn,pn,ci,cj,ca,C);CHKERRQ(ierr); 139d20bfe6fSHong Zhang 140d20bfe6fSHong Zhang /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 141d20bfe6fSHong Zhang /* Since these are PETSc arrays, change flags to free them as necessary. */ 142d20bfe6fSHong Zhang c = (Mat_SeqAIJ *)((*C)->data); 143e6b907acSBarry Smith c->free_a = PETSC_TRUE; 144e6b907acSBarry Smith c->free_ij = PETSC_TRUE; 145d20bfe6fSHong Zhang c->nonew = 0; 146d20bfe6fSHong Zhang 147d20bfe6fSHong Zhang /* Clean up. */ 148d20bfe6fSHong Zhang ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr); 149*f2b054eeSHong Zhang #if defined(PETSC_USE_INFO) 150*f2b054eeSHong Zhang if (ci[pn] != 0) { 151*f2b054eeSHong Zhang PetscReal afill = ((PetscReal)ci[pn])/ai[am]; 152*f2b054eeSHong Zhang if (afill < 1.0) afill = 1.0; 153*f2b054eeSHong Zhang ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr); 154*f2b054eeSHong Zhang ierr = PetscInfo1((*C),"Use MatPtAP(A,P,MatReuse,%G,&C) for best performance.\n",afill);CHKERRQ(ierr); 155*f2b054eeSHong Zhang } else { 156*f2b054eeSHong Zhang ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr); 157*f2b054eeSHong Zhang } 158*f2b054eeSHong Zhang #endif 159eb9c0419SKris Buschelman PetscFunctionReturn(0); 160eb9c0419SKris Buschelman } 161eb9c0419SKris Buschelman 162eb9c0419SKris Buschelman #undef __FUNCT__ 1639af31e4aSHong Zhang #define __FUNCT__ "MatPtAPNumeric_SeqAIJ_SeqAIJ" 164dfbe8321SBarry Smith PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat A,Mat P,Mat C) 165dfbe8321SBarry Smith { 166dfbe8321SBarry Smith PetscErrorCode ierr; 167b1d57f15SBarry Smith PetscInt flops=0; 168d20bfe6fSHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 169d20bfe6fSHong Zhang Mat_SeqAIJ *p = (Mat_SeqAIJ *) P->data; 170d20bfe6fSHong Zhang Mat_SeqAIJ *c = (Mat_SeqAIJ *) C->data; 171b1d57f15SBarry Smith PetscInt *ai=a->i,*aj=a->j,*apj,*apjdense,*pi=p->i,*pj=p->j,*pJ=p->j,*pjj; 172b1d57f15SBarry Smith PetscInt *ci=c->i,*cj=c->j,*cjj; 173899cda47SBarry Smith PetscInt am=A->rmap.N,cn=C->cmap.N,cm=C->rmap.N; 174b1d57f15SBarry Smith PetscInt i,j,k,anzi,pnzi,apnzj,nextap,pnzj,prow,crow; 175d20bfe6fSHong Zhang MatScalar *aa=a->a,*apa,*pa=p->a,*pA=p->a,*paj,*ca=c->a,*caj; 176eb9c0419SKris Buschelman 177eb9c0419SKris Buschelman PetscFunctionBegin; 178d20bfe6fSHong Zhang /* Allocate temporary array for storage of one row of A*P */ 179b1d57f15SBarry Smith ierr = PetscMalloc(cn*(sizeof(MatScalar)+2*sizeof(PetscInt)),&apa);CHKERRQ(ierr); 180b1d57f15SBarry Smith ierr = PetscMemzero(apa,cn*(sizeof(MatScalar)+2*sizeof(PetscInt)));CHKERRQ(ierr); 181eb9c0419SKris Buschelman 182b1d57f15SBarry Smith apj = (PetscInt *)(apa + cn); 183d20bfe6fSHong Zhang apjdense = apj + cn; 184d20bfe6fSHong Zhang 185d20bfe6fSHong Zhang /* Clear old values in C */ 186d20bfe6fSHong Zhang ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 187d20bfe6fSHong Zhang 188d20bfe6fSHong Zhang for (i=0;i<am;i++) { 189d20bfe6fSHong Zhang /* Form sparse row of A*P */ 190d20bfe6fSHong Zhang anzi = ai[i+1] - ai[i]; 191d20bfe6fSHong Zhang apnzj = 0; 192d20bfe6fSHong Zhang for (j=0;j<anzi;j++) { 193d20bfe6fSHong Zhang prow = *aj++; 194d20bfe6fSHong Zhang pnzj = pi[prow+1] - pi[prow]; 195d20bfe6fSHong Zhang pjj = pj + pi[prow]; 196d20bfe6fSHong Zhang paj = pa + pi[prow]; 197d20bfe6fSHong Zhang for (k=0;k<pnzj;k++) { 198d20bfe6fSHong Zhang if (!apjdense[pjj[k]]) { 199d20bfe6fSHong Zhang apjdense[pjj[k]] = -1; 200d20bfe6fSHong Zhang apj[apnzj++] = pjj[k]; 201d20bfe6fSHong Zhang } 202d20bfe6fSHong Zhang apa[pjj[k]] += (*aa)*paj[k]; 203d20bfe6fSHong Zhang } 204d20bfe6fSHong Zhang flops += 2*pnzj; 205d20bfe6fSHong Zhang aa++; 206d20bfe6fSHong Zhang } 207d20bfe6fSHong Zhang 208d20bfe6fSHong Zhang /* Sort the j index array for quick sparse axpy. */ 209d6509036SKris Buschelman /* Note: a array does not need sorting as it is in dense storage locations. */ 210d20bfe6fSHong Zhang ierr = PetscSortInt(apnzj,apj);CHKERRQ(ierr); 211d20bfe6fSHong Zhang 212d20bfe6fSHong Zhang /* Compute P^T*A*P using outer product (P^T)[:,j]*(A*P)[j,:]. */ 213d20bfe6fSHong Zhang pnzi = pi[i+1] - pi[i]; 214d20bfe6fSHong Zhang for (j=0;j<pnzi;j++) { 215d20bfe6fSHong Zhang nextap = 0; 216d20bfe6fSHong Zhang crow = *pJ++; 217d20bfe6fSHong Zhang cjj = cj + ci[crow]; 218d20bfe6fSHong Zhang caj = ca + ci[crow]; 219d20bfe6fSHong Zhang /* Perform sparse axpy operation. Note cjj includes apj. */ 220d20bfe6fSHong Zhang for (k=0;nextap<apnzj;k++) { 221f9f3b1caSBarry Smith #if defined(PETSC_USE_DEBUG) 222f9f3b1caSBarry Smith if (k >= ci[crow+1] - ci[crow]) { 223f9f3b1caSBarry Smith SETERRQ2(PETSC_ERR_PLIB,"k too large k %d, crow %d",k,crow); 224f9f3b1caSBarry Smith } 225f9f3b1caSBarry Smith #endif 226d20bfe6fSHong Zhang if (cjj[k]==apj[nextap]) { 227d20bfe6fSHong Zhang caj[k] += (*pA)*apa[apj[nextap++]]; 228d20bfe6fSHong Zhang } 229d20bfe6fSHong Zhang } 230d20bfe6fSHong Zhang flops += 2*apnzj; 231d20bfe6fSHong Zhang pA++; 232d20bfe6fSHong Zhang } 233d20bfe6fSHong Zhang 234d20bfe6fSHong Zhang /* Zero the current row info for A*P */ 235d20bfe6fSHong Zhang for (j=0;j<apnzj;j++) { 236d20bfe6fSHong Zhang apa[apj[j]] = 0.; 237d20bfe6fSHong Zhang apjdense[apj[j]] = 0; 238d20bfe6fSHong Zhang } 239d20bfe6fSHong Zhang } 240d20bfe6fSHong Zhang 241d20bfe6fSHong Zhang /* Assemble the final matrix and clean up */ 242d20bfe6fSHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 243d20bfe6fSHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 244d20bfe6fSHong Zhang ierr = PetscFree(apa);CHKERRQ(ierr); 245d20bfe6fSHong Zhang ierr = PetscLogFlops(flops);CHKERRQ(ierr); 246d20bfe6fSHong Zhang 247eb9c0419SKris Buschelman PetscFunctionReturn(0); 248eb9c0419SKris Buschelman } 249