1eb9c0419SKris Buschelman /* 29af31e4aSHong Zhang Defines projective product routines where A is a AIJ 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" 89af31e4aSHong Zhang #include "src/mat/impls/aij/mpi/mpiaij.h" 9eb9c0419SKris Buschelman 10eb9c0419SKris Buschelman #undef __FUNCT__ 119af31e4aSHong Zhang #define __FUNCT__ "MatPtAP" 124d3841fdSKris Buschelman /*@ 139af31e4aSHong Zhang MatPtAP - Creates the matrix projection C = P^T * A * P 144d3841fdSKris Buschelman 154d3841fdSKris Buschelman Collective on Mat 164d3841fdSKris Buschelman 174d3841fdSKris Buschelman Input Parameters: 184d3841fdSKris Buschelman + A - the matrix 19f747e1aeSHong Zhang . P - the projection matrix 20f747e1aeSHong Zhang . scall - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX 21f747e1aeSHong Zhang - fill - expected fill as ratio of nnz(C)/(nnz(A) + nnz(P)) 224d3841fdSKris Buschelman 234d3841fdSKris Buschelman Output Parameters: 244d3841fdSKris Buschelman . C - the product matrix 254d3841fdSKris Buschelman 264d3841fdSKris Buschelman Notes: 274d3841fdSKris Buschelman C will be created and must be destroyed by the user with MatDestroy(). 284d3841fdSKris Buschelman 294d3841fdSKris Buschelman This routine is currently only implemented for pairs of SeqAIJ matrices and classes 304d3841fdSKris Buschelman which inherit from SeqAIJ. C will be of type MATSEQAIJ. 314d3841fdSKris Buschelman 324d3841fdSKris Buschelman Level: intermediate 334d3841fdSKris Buschelman 349af31e4aSHong Zhang .seealso: MatPtAPSymbolic(),MatPtAPNumeric(),MatMatMult() 354d3841fdSKris Buschelman @*/ 36dfbe8321SBarry Smith PetscErrorCode MatPtAP(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C) { 37dfbe8321SBarry Smith PetscErrorCode ierr; 38534c1384SKris Buschelman PetscErrorCode (*fA)(Mat,Mat,MatReuse,PetscReal,Mat *); 39534c1384SKris Buschelman PetscErrorCode (*fP)(Mat,Mat,MatReuse,PetscReal,Mat *); 40eb9c0419SKris Buschelman 41eb9c0419SKris Buschelman PetscFunctionBegin; 429af31e4aSHong Zhang PetscValidHeaderSpecific(A,MAT_COOKIE,1); 439af31e4aSHong Zhang PetscValidType(A,1); 449af31e4aSHong Zhang MatPreallocated(A); 459af31e4aSHong Zhang if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 469af31e4aSHong Zhang if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 479af31e4aSHong Zhang PetscValidHeaderSpecific(P,MAT_COOKIE,2); 489af31e4aSHong Zhang PetscValidType(P,2); 499af31e4aSHong Zhang MatPreallocated(P); 509af31e4aSHong Zhang if (!P->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 519af31e4aSHong Zhang if (P->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 529af31e4aSHong Zhang PetscValidPointer(C,3); 53534c1384SKris Buschelman 549af31e4aSHong Zhang if (P->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->M,A->N); 55eb9c0419SKris Buschelman 569af31e4aSHong Zhang if (fill <=0.0) SETERRQ1(PETSC_ERR_ARG_SIZ,"fill=%g must be > 0.0",fill); 57eb9c0419SKris Buschelman 58534c1384SKris Buschelman /* For now, we do not dispatch based on the type of A and P */ 59534c1384SKris Buschelman /* When implementations like _SeqAIJ_MAIJ exist, attack the multiple dispatch problem. */ 60534c1384SKris Buschelman fA = A->ops->ptap; 61534c1384SKris Buschelman if (!fA) SETERRQ1(PETSC_ERR_SUP,"MatPtAP not supported for A of type %s",A->type_name); 62534c1384SKris Buschelman fP = P->ops->ptap; 63534c1384SKris Buschelman if (!fP) SETERRQ1(PETSC_ERR_SUP,"MatPtAP not supported for P of type %s",P->type_name); 64534c1384SKris Buschelman if (fP!=fA) SETERRQ2(PETSC_ERR_ARG_INCOMP,"MatPtAP requires A, %s, to be compatible with P, %s",A->type_name,P->type_name); 65534c1384SKris Buschelman 669af31e4aSHong Zhang ierr = PetscLogEventBegin(MAT_PtAP,A,P,0,0);CHKERRQ(ierr); 67534c1384SKris Buschelman ierr = (*fA)(A,P,scall,fill,C);CHKERRQ(ierr); 689af31e4aSHong Zhang ierr = PetscLogEventEnd(MAT_PtAP,A,P,0,0);CHKERRQ(ierr); 69eb9c0419SKris Buschelman PetscFunctionReturn(0); 70eb9c0419SKris Buschelman } 71eb9c0419SKris Buschelman 72eb9c0419SKris Buschelman #undef __FUNCT__ 73*ff134f7aSHong Zhang #define __FUNCT__ "MatPtAP_MPIAIJ_MPIAIJ" 74*ff134f7aSHong Zhang PetscErrorCode MatPtAP_MPIAIJ_MPIAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C) 75*ff134f7aSHong Zhang { 76*ff134f7aSHong Zhang PetscErrorCode ierr; 77*ff134f7aSHong Zhang Mat C_mpi,AP_seq,P_seq,P_subseq,*psubseq; 78*ff134f7aSHong Zhang Mat_MPIAIJ *p = (Mat_MPIAIJ*)P->data; 79*ff134f7aSHong Zhang Mat_MatMatMultMPI *mult; 80*ff134f7aSHong Zhang int i,prow,prstart,prend,m=P->m,pncols; 81*ff134f7aSHong Zhang const int *pcols; 82*ff134f7aSHong Zhang const PetscScalar *pvals; 83*ff134f7aSHong Zhang int rank; 84*ff134f7aSHong Zhang 85*ff134f7aSHong Zhang PetscFunctionBegin; 86*ff134f7aSHong Zhang ierr = MPI_Comm_rank(A->comm,&rank);CHKERRQ(ierr); 87*ff134f7aSHong Zhang 88*ff134f7aSHong Zhang ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ(A,P,fill,&C_mpi);CHKERRQ(ierr); 89*ff134f7aSHong Zhang mult = (Mat_MatMatMultMPI*)C_mpi->spptr; 90*ff134f7aSHong Zhang P_seq = mult->bseq[0]; 91*ff134f7aSHong Zhang AP_seq = mult->C_seq; 92*ff134f7aSHong Zhang prstart = mult->brstart; 93*ff134f7aSHong Zhang prend = prstart + m; 94*ff134f7aSHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," [%d] prstart: %d, prend: %d, dim of P_seq: %d, %d\n",rank,prstart,prend,P_seq->m,P_seq->n); 95*ff134f7aSHong Zhang 96*ff134f7aSHong Zhang /* get seq matrix P_subseq by taking local rows of P */ 97*ff134f7aSHong Zhang IS isrow,iscol; 98*ff134f7aSHong Zhang ierr = ISCreateStride(PETSC_COMM_SELF,m,prstart,1,&isrow);CHKERRQ(ierr); 99*ff134f7aSHong Zhang ierr = ISCreateStride(PETSC_COMM_SELF,P_seq->n,0,1,&iscol);CHKERRQ(ierr); 100*ff134f7aSHong Zhang ierr = MatGetSubMatrices(P_seq,1,&isrow,&iscol,MAT_INITIAL_MATRIX,&psubseq);CHKERRQ(ierr); 101*ff134f7aSHong Zhang P_subseq = psubseq[0]; 102*ff134f7aSHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," [%d] dim of P_subseq: %d, %d\n",rank,P_subseq->m,P_subseq->n); 103*ff134f7aSHong Zhang 104*ff134f7aSHong Zhang /* Compute P_subseq^T*C_seq using outer product (P_loc^T)[:,i]*C_seq[i,:]. */ 105*ff134f7aSHong Zhang for (i=0;i<m;i++) { 106*ff134f7aSHong Zhang prow = prstart + i; 107*ff134f7aSHong Zhang ierr = MatGetRow(P_seq,prow,&pncols,&pcols,&pvals);CHKERRQ(ierr); 108*ff134f7aSHong Zhang ierr = MatRestoreRow(P_seq,prow,&pncols,&pcols,&pvals);CHKERRQ(ierr); 109*ff134f7aSHong Zhang } 110*ff134f7aSHong Zhang 111*ff134f7aSHong Zhang *C = C_mpi; /* to be removed! */ 112*ff134f7aSHong Zhang PetscFunctionReturn(0); 113*ff134f7aSHong Zhang } 114*ff134f7aSHong Zhang 115*ff134f7aSHong Zhang #undef __FUNCT__ 116*ff134f7aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ" 117*ff134f7aSHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat B,PetscReal fill,Mat *C) 118*ff134f7aSHong Zhang { 119*ff134f7aSHong Zhang Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data; 120*ff134f7aSHong Zhang PetscErrorCode ierr; 121*ff134f7aSHong Zhang 122*ff134f7aSHong Zhang PetscFunctionBegin; 123*ff134f7aSHong Zhang PetscFunctionReturn(0); 124*ff134f7aSHong Zhang } 125*ff134f7aSHong Zhang 126*ff134f7aSHong Zhang #undef __FUNCT__ 127*ff134f7aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ" 128*ff134f7aSHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat B,Mat C) 129*ff134f7aSHong Zhang { 130*ff134f7aSHong Zhang Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data; 131*ff134f7aSHong Zhang PetscErrorCode ierr; 132*ff134f7aSHong Zhang 133*ff134f7aSHong Zhang PetscFunctionBegin; 134*ff134f7aSHong Zhang PetscFunctionReturn(0); 135*ff134f7aSHong Zhang } 136*ff134f7aSHong Zhang 137*ff134f7aSHong Zhang #undef __FUNCT__ 1389af31e4aSHong Zhang #define __FUNCT__ "MatPtAP_SeqAIJ_SeqAIJ" 139dfbe8321SBarry Smith PetscErrorCode MatPtAP_SeqAIJ_SeqAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C) 1409af31e4aSHong Zhang { 141dfbe8321SBarry Smith PetscErrorCode ierr; 1429af31e4aSHong Zhang PetscFunctionBegin; 1439af31e4aSHong Zhang if (scall == MAT_INITIAL_MATRIX){ 1449af31e4aSHong Zhang ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ(A,P,fill,C);CHKERRQ(ierr); 1459af31e4aSHong Zhang } 1469af31e4aSHong Zhang ierr = MatPtAPNumeric_SeqAIJ_SeqAIJ(A,P,*C);CHKERRQ(ierr); 1479af31e4aSHong Zhang PetscFunctionReturn(0); 1489af31e4aSHong Zhang } 1499af31e4aSHong Zhang 1509af31e4aSHong Zhang #undef __FUNCT__ 1519af31e4aSHong Zhang #define __FUNCT__ "MatPtAPSymbolic" 1526849ba73SBarry Smith /* 1539af31e4aSHong Zhang MatPtAPSymbolic - Creates the (i,j) structure of the matrix projection C = P^T * A * P 1544d3841fdSKris Buschelman 1554d3841fdSKris Buschelman Collective on Mat 1564d3841fdSKris Buschelman 1574d3841fdSKris Buschelman Input Parameters: 1584d3841fdSKris Buschelman + A - the matrix 1594d3841fdSKris Buschelman - P - the projection matrix 1604d3841fdSKris Buschelman 1614d3841fdSKris Buschelman Output Parameters: 1624d3841fdSKris Buschelman . C - the (i,j) structure of the product matrix 1634d3841fdSKris Buschelman 1644d3841fdSKris Buschelman Notes: 1654d3841fdSKris Buschelman C will be created and must be destroyed by the user with MatDestroy(). 1664d3841fdSKris Buschelman 1674d3841fdSKris Buschelman This routine is currently only implemented for pairs of SeqAIJ matrices and classes 1684d3841fdSKris Buschelman which inherit from SeqAIJ. C will be of type MATSEQAIJ. The product is computed using 1699af31e4aSHong Zhang this (i,j) structure by calling MatPtAPNumeric(). 1704d3841fdSKris Buschelman 1714d3841fdSKris Buschelman Level: intermediate 1724d3841fdSKris Buschelman 1739af31e4aSHong Zhang .seealso: MatPtAP(),MatPtAPNumeric(),MatMatMultSymbolic() 1746849ba73SBarry Smith */ 175dfbe8321SBarry Smith PetscErrorCode MatPtAPSymbolic(Mat A,Mat P,PetscReal fill,Mat *C) { 176dfbe8321SBarry Smith PetscErrorCode ierr; 177534c1384SKris Buschelman PetscErrorCode (*fA)(Mat,Mat,PetscReal,Mat*); 178534c1384SKris Buschelman PetscErrorCode (*fP)(Mat,Mat,PetscReal,Mat*); 179eb9c0419SKris Buschelman 180eb9c0419SKris Buschelman PetscFunctionBegin; 181eb9c0419SKris Buschelman 1824482741eSBarry Smith PetscValidHeaderSpecific(A,MAT_COOKIE,1); 183c9780b6fSBarry Smith PetscValidType(A,1); 184eb9c0419SKris Buschelman MatPreallocated(A); 185eb9c0419SKris Buschelman if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 186eb9c0419SKris Buschelman if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 187eb9c0419SKris Buschelman 1884482741eSBarry Smith PetscValidHeaderSpecific(P,MAT_COOKIE,2); 189c9780b6fSBarry Smith PetscValidType(P,2); 190eb9c0419SKris Buschelman MatPreallocated(P); 191eb9c0419SKris Buschelman if (!P->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 192eb9c0419SKris Buschelman if (P->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 193eb9c0419SKris Buschelman 1944482741eSBarry Smith PetscValidPointer(C,3); 1954482741eSBarry Smith 196eb9c0419SKris Buschelman if (P->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->M,A->N); 197eb9c0419SKris Buschelman if (A->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix 'A' must be square, %d != %d",A->M,A->N); 198eb9c0419SKris Buschelman 199534c1384SKris Buschelman /* For now, we do not dispatch based on the type of A and P */ 200534c1384SKris Buschelman /* When implementations like _SeqAIJ_MAIJ exist, attack the multiple dispatch problem. */ 201534c1384SKris Buschelman fA = A->ops->ptapsymbolic; 202534c1384SKris Buschelman if (!fA) SETERRQ1(PETSC_ERR_SUP,"MatPtAPSymbolic not supported for A of type %s",A->type_name); 203534c1384SKris Buschelman fP = P->ops->ptapsymbolic; 204534c1384SKris Buschelman if (!fP) SETERRQ1(PETSC_ERR_SUP,"MatPtAPSymbolic not supported for P of type %s",P->type_name); 205534c1384SKris Buschelman if (fP!=fA) SETERRQ2(PETSC_ERR_ARG_INCOMP,"MatPtAPSymbolic requires A, %s, to be compatible with P, %s",A->type_name,P->type_name); 2064d3841fdSKris Buschelman 207534c1384SKris Buschelman ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 208534c1384SKris Buschelman ierr = (*fA)(A,P,fill,C);CHKERRQ(ierr); 209534c1384SKris Buschelman ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 210eb9c0419SKris Buschelman 211eb9c0419SKris Buschelman PetscFunctionReturn(0); 212eb9c0419SKris Buschelman } 213eb9c0419SKris Buschelman 214f747e1aeSHong Zhang typedef struct { 215f747e1aeSHong Zhang Mat symAP; 216f747e1aeSHong Zhang } Mat_PtAPstruct; 217f747e1aeSHong Zhang 21878a80504SBarry Smith EXTERN PetscErrorCode MatDestroy_SeqAIJ(Mat); 21978a80504SBarry Smith 220f747e1aeSHong Zhang #undef __FUNCT__ 221f747e1aeSHong Zhang #define __FUNCT__ "MatDestroy_SeqAIJ_PtAP" 222f4a850bbSBarry Smith PetscErrorCode MatDestroy_SeqAIJ_PtAP(Mat A) 223f747e1aeSHong Zhang { 224f4a850bbSBarry Smith PetscErrorCode ierr; 225f747e1aeSHong Zhang Mat_PtAPstruct *ptap=(Mat_PtAPstruct*)A->spptr; 226f747e1aeSHong Zhang 227f747e1aeSHong Zhang PetscFunctionBegin; 228f747e1aeSHong Zhang ierr = MatDestroy(ptap->symAP);CHKERRQ(ierr); 229f747e1aeSHong Zhang ierr = PetscFree(ptap);CHKERRQ(ierr); 23078a80504SBarry Smith ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr); 231f747e1aeSHong Zhang PetscFunctionReturn(0); 232f747e1aeSHong Zhang } 233f747e1aeSHong Zhang 234eb9c0419SKris Buschelman #undef __FUNCT__ 2359af31e4aSHong Zhang #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqAIJ" 236dfbe8321SBarry Smith PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat P,PetscReal fill,Mat *C) { 237dfbe8321SBarry Smith PetscErrorCode ierr; 238f4a850bbSBarry Smith int *pti,*ptj; 239f747e1aeSHong Zhang Mat Pt,AP; 240f747e1aeSHong Zhang Mat_PtAPstruct *ptap; 241eb9c0419SKris Buschelman 242eb9c0419SKris Buschelman PetscFunctionBegin; 243f747e1aeSHong Zhang /* create symbolic Pt */ 244eb9c0419SKris Buschelman ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr); 245f747e1aeSHong Zhang ierr = MatCreateSeqAIJWithArrays(P->comm,P->N,P->M,pti,ptj,PETSC_NULL,&Pt);CHKERRQ(ierr); 246eb9c0419SKris Buschelman 247f747e1aeSHong Zhang /* get symbolic AP=A*P and C=Pt*AP */ 248f747e1aeSHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,P,fill,&AP);CHKERRQ(ierr); 249f747e1aeSHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(Pt,AP,fill,C);CHKERRQ(ierr); 250eb9c0419SKris Buschelman 251f747e1aeSHong Zhang /* clean up */ 252f747e1aeSHong Zhang ierr = MatRestoreSymbolicTranspose_SeqAIJ(Pt,&pti,&ptj);CHKERRQ(ierr); 253f747e1aeSHong Zhang ierr = MatDestroy(Pt);CHKERRQ(ierr); 254eb9c0419SKris Buschelman 255f747e1aeSHong Zhang /* save symbolic AP - to be used by MatPtAPNumeric_SeqAIJ_SeqAIJ() */ 256f747e1aeSHong Zhang ierr = PetscNew(Mat_PtAPstruct,&ptap);CHKERRQ(ierr); 257f747e1aeSHong Zhang ptap->symAP = AP; 258f747e1aeSHong Zhang (*C)->spptr = (void*)ptap; 259f747e1aeSHong Zhang (*C)->ops->destroy = MatDestroy_SeqAIJ_PtAP; 260eb9c0419SKris Buschelman 261eb9c0419SKris Buschelman PetscFunctionReturn(0); 262eb9c0419SKris Buschelman } 263eb9c0419SKris Buschelman 2643985e5eaSKris Buschelman #include "src/mat/impls/maij/maij.h" 2653985e5eaSKris Buschelman EXTERN_C_BEGIN 2663985e5eaSKris Buschelman #undef __FUNCT__ 2679af31e4aSHong Zhang #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqMAIJ" 268dfbe8321SBarry Smith PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqMAIJ(Mat A,Mat PP,Mat *C) { 2695c66b693SKris Buschelman /* This routine requires testing -- I don't think it works. */ 270dfbe8321SBarry Smith PetscErrorCode ierr; 2713985e5eaSKris Buschelman FreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 2723985e5eaSKris Buschelman Mat_SeqMAIJ *pp=(Mat_SeqMAIJ*)PP->data; 2733985e5eaSKris Buschelman Mat P=pp->AIJ; 2743985e5eaSKris Buschelman Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*p=(Mat_SeqAIJ*)P->data,*c; 2753985e5eaSKris Buschelman int *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj; 2763985e5eaSKris Buschelman int *ci,*cj,*denserow,*sparserow,*ptadenserow,*ptasparserow,*ptaj; 2773985e5eaSKris Buschelman int an=A->N,am=A->M,pn=P->N,pm=P->M,ppdof=pp->dof; 278fe05a634SKris Buschelman int i,j,k,dof,pdof,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi; 2793985e5eaSKris Buschelman MatScalar *ca; 2803985e5eaSKris Buschelman 2813985e5eaSKris Buschelman PetscFunctionBegin; 2823985e5eaSKris Buschelman /* Start timer */ 2839af31e4aSHong Zhang ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,PP,0,0);CHKERRQ(ierr); 2843985e5eaSKris Buschelman 2853985e5eaSKris Buschelman /* Get ij structure of P^T */ 2863985e5eaSKris Buschelman ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr); 2873985e5eaSKris Buschelman 2883985e5eaSKris Buschelman /* Allocate ci array, arrays for fill computation and */ 2893985e5eaSKris Buschelman /* free space for accumulating nonzero column info */ 2903985e5eaSKris Buschelman ierr = PetscMalloc((pn+1)*sizeof(int),&ci);CHKERRQ(ierr); 2913985e5eaSKris Buschelman ci[0] = 0; 2923985e5eaSKris Buschelman 2933985e5eaSKris Buschelman ierr = PetscMalloc((2*pn+2*an+1)*sizeof(int),&ptadenserow);CHKERRQ(ierr); 2943985e5eaSKris Buschelman ierr = PetscMemzero(ptadenserow,(2*pn+2*an+1)*sizeof(int));CHKERRQ(ierr); 2953985e5eaSKris Buschelman ptasparserow = ptadenserow + an; 2963985e5eaSKris Buschelman denserow = ptasparserow + an; 2973985e5eaSKris Buschelman sparserow = denserow + pn; 2983985e5eaSKris Buschelman 2993985e5eaSKris Buschelman /* Set initial free space to be nnz(A) scaled by aspect ratio of P. */ 3003985e5eaSKris Buschelman /* This should be reasonable if sparsity of PtAP is similar to that of A. */ 3013985e5eaSKris Buschelman ierr = GetMoreSpace((ai[am]/pm)*pn,&free_space); 3023985e5eaSKris Buschelman current_space = free_space; 3033985e5eaSKris Buschelman 3043985e5eaSKris Buschelman /* Determine symbolic info for each row of C: */ 3053985e5eaSKris Buschelman for (i=0;i<pn/ppdof;i++) { 3063985e5eaSKris Buschelman ptnzi = pti[i+1] - pti[i]; 3073985e5eaSKris Buschelman ptanzi = 0; 3083985e5eaSKris Buschelman ptJ = ptj + pti[i]; 3093985e5eaSKris Buschelman for (dof=0;dof<ppdof;dof++) { 3103985e5eaSKris Buschelman /* Determine symbolic row of PtA: */ 3113985e5eaSKris Buschelman for (j=0;j<ptnzi;j++) { 3123985e5eaSKris Buschelman arow = ptJ[j] + dof; 3133985e5eaSKris Buschelman anzj = ai[arow+1] - ai[arow]; 3143985e5eaSKris Buschelman ajj = aj + ai[arow]; 3153985e5eaSKris Buschelman for (k=0;k<anzj;k++) { 3163985e5eaSKris Buschelman if (!ptadenserow[ajj[k]]) { 3173985e5eaSKris Buschelman ptadenserow[ajj[k]] = -1; 3183985e5eaSKris Buschelman ptasparserow[ptanzi++] = ajj[k]; 3193985e5eaSKris Buschelman } 3203985e5eaSKris Buschelman } 3213985e5eaSKris Buschelman } 3223985e5eaSKris Buschelman /* Using symbolic info for row of PtA, determine symbolic info for row of C: */ 3233985e5eaSKris Buschelman ptaj = ptasparserow; 3243985e5eaSKris Buschelman cnzi = 0; 3253985e5eaSKris Buschelman for (j=0;j<ptanzi;j++) { 326fe05a634SKris Buschelman pdof = *ptaj%dof; 3273985e5eaSKris Buschelman prow = (*ptaj++)/dof; 3283985e5eaSKris Buschelman pnzj = pi[prow+1] - pi[prow]; 3293985e5eaSKris Buschelman pjj = pj + pi[prow]; 3303985e5eaSKris Buschelman for (k=0;k<pnzj;k++) { 331fe05a634SKris Buschelman if (!denserow[pjj[k]+pdof]) { 332fe05a634SKris Buschelman denserow[pjj[k]+pdof] = -1; 333fe05a634SKris Buschelman sparserow[cnzi++] = pjj[k]+pdof; 3343985e5eaSKris Buschelman } 3353985e5eaSKris Buschelman } 3363985e5eaSKris Buschelman } 3373985e5eaSKris Buschelman 3383985e5eaSKris Buschelman /* sort sparserow */ 3393985e5eaSKris Buschelman ierr = PetscSortInt(cnzi,sparserow);CHKERRQ(ierr); 3403985e5eaSKris Buschelman 3413985e5eaSKris Buschelman /* If free space is not available, make more free space */ 3423985e5eaSKris Buschelman /* Double the amount of total space in the list */ 3433985e5eaSKris Buschelman if (current_space->local_remaining<cnzi) { 3443985e5eaSKris Buschelman ierr = GetMoreSpace(current_space->total_array_size,¤t_space);CHKERRQ(ierr); 3453985e5eaSKris Buschelman } 3463985e5eaSKris Buschelman 3473985e5eaSKris Buschelman /* Copy data into free space, and zero out denserows */ 3483985e5eaSKris Buschelman ierr = PetscMemcpy(current_space->array,sparserow,cnzi*sizeof(int));CHKERRQ(ierr); 3493985e5eaSKris Buschelman current_space->array += cnzi; 3503985e5eaSKris Buschelman current_space->local_used += cnzi; 3513985e5eaSKris Buschelman current_space->local_remaining -= cnzi; 3523985e5eaSKris Buschelman 3533985e5eaSKris Buschelman for (j=0;j<ptanzi;j++) { 3543985e5eaSKris Buschelman ptadenserow[ptasparserow[j]] = 0; 3553985e5eaSKris Buschelman } 3563985e5eaSKris Buschelman for (j=0;j<cnzi;j++) { 3573985e5eaSKris Buschelman denserow[sparserow[j]] = 0; 3583985e5eaSKris Buschelman } 3593985e5eaSKris Buschelman /* Aside: Perhaps we should save the pta info for the numerical factorization. */ 3603985e5eaSKris Buschelman /* For now, we will recompute what is needed. */ 3613985e5eaSKris Buschelman ci[i+1+dof] = ci[i+dof] + cnzi; 3623985e5eaSKris Buschelman } 3633985e5eaSKris Buschelman } 3643985e5eaSKris Buschelman /* nnz is now stored in ci[ptm], column indices are in the list of free space */ 3653985e5eaSKris Buschelman /* Allocate space for cj, initialize cj, and */ 3663985e5eaSKris Buschelman /* destroy list of free space and other temporary array(s) */ 3673985e5eaSKris Buschelman ierr = PetscMalloc((ci[pn]+1)*sizeof(int),&cj);CHKERRQ(ierr); 3683985e5eaSKris Buschelman ierr = MakeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 3693985e5eaSKris Buschelman ierr = PetscFree(ptadenserow);CHKERRQ(ierr); 3703985e5eaSKris Buschelman 3713985e5eaSKris Buschelman /* Allocate space for ca */ 3723985e5eaSKris Buschelman ierr = PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 3733985e5eaSKris Buschelman ierr = PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));CHKERRQ(ierr); 3743985e5eaSKris Buschelman 3753985e5eaSKris Buschelman /* put together the new matrix */ 3763985e5eaSKris Buschelman ierr = MatCreateSeqAIJWithArrays(A->comm,pn,pn,ci,cj,ca,C);CHKERRQ(ierr); 3773985e5eaSKris Buschelman 3783985e5eaSKris Buschelman /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 3793985e5eaSKris Buschelman /* Since these are PETSc arrays, change flags to free them as necessary. */ 3803985e5eaSKris Buschelman c = (Mat_SeqAIJ *)((*C)->data); 3813985e5eaSKris Buschelman c->freedata = PETSC_TRUE; 3823985e5eaSKris Buschelman c->nonew = 0; 3833985e5eaSKris Buschelman 3843985e5eaSKris Buschelman /* Clean up. */ 3853985e5eaSKris Buschelman ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr); 3863985e5eaSKris Buschelman 3879af31e4aSHong Zhang ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,PP,0,0);CHKERRQ(ierr); 3883985e5eaSKris Buschelman PetscFunctionReturn(0); 3893985e5eaSKris Buschelman } 3903985e5eaSKris Buschelman EXTERN_C_END 3913985e5eaSKris Buschelman 392eb9c0419SKris Buschelman #undef __FUNCT__ 3939af31e4aSHong Zhang #define __FUNCT__ "MatPtAPNumeric" 3946849ba73SBarry Smith /* 3959af31e4aSHong Zhang MatPtAPNumeric - Computes the matrix projection C = P^T * A * P 3964d3841fdSKris Buschelman 3974d3841fdSKris Buschelman Collective on Mat 3984d3841fdSKris Buschelman 3994d3841fdSKris Buschelman Input Parameters: 4004d3841fdSKris Buschelman + A - the matrix 4014d3841fdSKris Buschelman - P - the projection matrix 4024d3841fdSKris Buschelman 4034d3841fdSKris Buschelman Output Parameters: 4044d3841fdSKris Buschelman . C - the product matrix 4054d3841fdSKris Buschelman 4064d3841fdSKris Buschelman Notes: 4079af31e4aSHong Zhang C must have been created by calling MatPtAPSymbolic and must be destroyed by 4084d3841fdSKris Buschelman the user using MatDeatroy(). 4094d3841fdSKris Buschelman 4104d3841fdSKris Buschelman This routine is currently only implemented for pairs of SeqAIJ matrices and classes 4114d3841fdSKris Buschelman which inherit from SeqAIJ. C will be of type MATSEQAIJ. 4124d3841fdSKris Buschelman 4134d3841fdSKris Buschelman Level: intermediate 4144d3841fdSKris Buschelman 4159af31e4aSHong Zhang .seealso: MatPtAP(),MatPtAPSymbolic(),MatMatMultNumeric() 4166849ba73SBarry Smith */ 417dfbe8321SBarry Smith PetscErrorCode MatPtAPNumeric(Mat A,Mat P,Mat C) { 418dfbe8321SBarry Smith PetscErrorCode ierr; 419534c1384SKris Buschelman PetscErrorCode (*fA)(Mat,Mat,Mat); 420534c1384SKris Buschelman PetscErrorCode (*fP)(Mat,Mat,Mat); 421eb9c0419SKris Buschelman 422eb9c0419SKris Buschelman PetscFunctionBegin; 423eb9c0419SKris Buschelman 4244482741eSBarry Smith PetscValidHeaderSpecific(A,MAT_COOKIE,1); 425c9780b6fSBarry Smith PetscValidType(A,1); 426eb9c0419SKris Buschelman MatPreallocated(A); 427eb9c0419SKris Buschelman if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 428eb9c0419SKris Buschelman if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 429eb9c0419SKris Buschelman 4304482741eSBarry Smith PetscValidHeaderSpecific(P,MAT_COOKIE,2); 431c9780b6fSBarry Smith PetscValidType(P,2); 432eb9c0419SKris Buschelman MatPreallocated(P); 433eb9c0419SKris Buschelman if (!P->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 434eb9c0419SKris Buschelman if (P->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 435eb9c0419SKris Buschelman 4364482741eSBarry Smith PetscValidHeaderSpecific(C,MAT_COOKIE,3); 437c9780b6fSBarry Smith PetscValidType(C,3); 438eb9c0419SKris Buschelman MatPreallocated(C); 439eb9c0419SKris Buschelman if (!C->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 440eb9c0419SKris Buschelman if (C->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 441eb9c0419SKris Buschelman 442eb9c0419SKris Buschelman if (P->N!=C->M) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->N,C->M); 443eb9c0419SKris Buschelman if (P->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->M,A->N); 444eb9c0419SKris Buschelman if (A->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix 'A' must be square, %d != %d",A->M,A->N); 445eb9c0419SKris Buschelman if (P->N!=C->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->N,C->N); 446eb9c0419SKris Buschelman 447534c1384SKris Buschelman /* For now, we do not dispatch based on the type of A and P */ 448534c1384SKris Buschelman /* When implementations like _SeqAIJ_MAIJ exist, attack the multiple dispatch problem. */ 449534c1384SKris Buschelman fA = A->ops->ptapnumeric; 450534c1384SKris Buschelman if (!fA) SETERRQ1(PETSC_ERR_SUP,"MatPtAPNumeric not supported for A of type %s",A->type_name); 451534c1384SKris Buschelman fP = P->ops->ptapnumeric; 452534c1384SKris Buschelman if (!fP) SETERRQ1(PETSC_ERR_SUP,"MatPtAPNumeric not supported for P of type %s",P->type_name); 453534c1384SKris Buschelman if (fP!=fA) SETERRQ2(PETSC_ERR_ARG_INCOMP,"MatPtAPNumeric requires A, %s, to be compatible with P, %s",A->type_name,P->type_name); 4544d3841fdSKris Buschelman 455534c1384SKris Buschelman ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 456534c1384SKris Buschelman ierr = (*fA)(A,P,C);CHKERRQ(ierr); 457534c1384SKris Buschelman ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 458eb9c0419SKris Buschelman 459eb9c0419SKris Buschelman PetscFunctionReturn(0); 460eb9c0419SKris Buschelman } 461eb9c0419SKris Buschelman 462eb9c0419SKris Buschelman #undef __FUNCT__ 4639af31e4aSHong Zhang #define __FUNCT__ "MatPtAPNumeric_SeqAIJ_SeqAIJ" 464dfbe8321SBarry Smith PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat A,Mat P,Mat C) 465dfbe8321SBarry Smith { 466dfbe8321SBarry Smith PetscErrorCode ierr; 467dfbe8321SBarry Smith int flops=0; 468eb9c0419SKris Buschelman Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 469eb9c0419SKris Buschelman Mat_SeqAIJ *p = (Mat_SeqAIJ *) P->data; 470eb9c0419SKris Buschelman Mat_SeqAIJ *c = (Mat_SeqAIJ *) C->data; 471f4a850bbSBarry Smith int *ai=a->i,*aj=a->j,*pi=p->i,*pj=p->j,*pJ=p->j,*pjj; 472eb9c0419SKris Buschelman int *ci=c->i,*cj=c->j,*cjj; 473eb9c0419SKris Buschelman int am=A->M,cn=C->N,cm=C->M; 474eb9c0419SKris Buschelman int i,j,k,anzi,pnzi,apnzj,nextap,pnzj,prow,crow; 475eb9c0419SKris Buschelman MatScalar *aa=a->a,*apa,*pa=p->a,*pA=p->a,*paj,*ca=c->a,*caj; 476f747e1aeSHong Zhang Mat_PtAPstruct *ptap=(Mat_PtAPstruct*)C->spptr; 477f747e1aeSHong Zhang Mat_SeqAIJ *ap = (Mat_SeqAIJ *)(ptap->symAP)->data; 478f747e1aeSHong Zhang int *api=ap->i,*apj=ap->j,apj_nextap; 479eb9c0419SKris Buschelman 480eb9c0419SKris Buschelman PetscFunctionBegin; 481eb9c0419SKris Buschelman /* Allocate temporary array for storage of one row of A*P */ 482f747e1aeSHong Zhang ierr = PetscMalloc(cn*sizeof(MatScalar),&apa);CHKERRQ(ierr); 483f747e1aeSHong Zhang ierr = PetscMemzero(apa,cn*sizeof(MatScalar));CHKERRQ(ierr); 484eb9c0419SKris Buschelman 485eb9c0419SKris Buschelman /* Clear old values in C */ 486eb9c0419SKris Buschelman ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 487eb9c0419SKris Buschelman 488eb9c0419SKris Buschelman for (i=0;i<am;i++) { 489f747e1aeSHong Zhang /* Get sparse values of A*P[i,:] */ 490eb9c0419SKris Buschelman anzi = ai[i+1] - ai[i]; 491eb9c0419SKris Buschelman apnzj = 0; 492eb9c0419SKris Buschelman for (j=0;j<anzi;j++) { 493eb9c0419SKris Buschelman prow = *aj++; 494eb9c0419SKris Buschelman pnzj = pi[prow+1] - pi[prow]; 495eb9c0419SKris Buschelman pjj = pj + pi[prow]; 496eb9c0419SKris Buschelman paj = pa + pi[prow]; 497eb9c0419SKris Buschelman for (k=0;k<pnzj;k++) { 498eb9c0419SKris Buschelman apa[pjj[k]] += (*aa)*paj[k]; 499eb9c0419SKris Buschelman } 500eb9c0419SKris Buschelman flops += 2*pnzj; 501eb9c0419SKris Buschelman aa++; 502eb9c0419SKris Buschelman } 503eb9c0419SKris Buschelman 504eb9c0419SKris Buschelman /* Compute P^T*A*P using outer product (P^T)[:,j]*(A*P)[j,:]. */ 505f747e1aeSHong Zhang apj = ap->j + api[i]; 506f747e1aeSHong Zhang apnzj = api[i+1] - api[i]; 507eb9c0419SKris Buschelman pnzi = pi[i+1] - pi[i]; 508eb9c0419SKris Buschelman for (j=0;j<pnzi;j++) { 509eb9c0419SKris Buschelman nextap = 0; 510eb9c0419SKris Buschelman crow = *pJ++; 511eb9c0419SKris Buschelman cjj = cj + ci[crow]; 512eb9c0419SKris Buschelman caj = ca + ci[crow]; 513eb9c0419SKris Buschelman /* Perform sparse axpy operation. Note cjj includes apj. */ 514eb9c0419SKris Buschelman for (k=0; nextap<apnzj; k++) { 515f747e1aeSHong Zhang apj_nextap = *(apj+nextap); 516f747e1aeSHong Zhang if (cjj[k]==apj_nextap) { 517f747e1aeSHong Zhang caj[k] += (*pA)*apa[apj_nextap]; 518f747e1aeSHong Zhang nextap++; 519eb9c0419SKris Buschelman } 520eb9c0419SKris Buschelman } 521eb9c0419SKris Buschelman flops += 2*apnzj; 522eb9c0419SKris Buschelman pA++; 523eb9c0419SKris Buschelman } 524eb9c0419SKris Buschelman 525f747e1aeSHong Zhang /* Zero the current row values for A*P */ 526f747e1aeSHong Zhang for (j=0;j<apnzj;j++) apa[apj[j]] = 0.0; 527eb9c0419SKris Buschelman } 528eb9c0419SKris Buschelman 529eb9c0419SKris Buschelman /* Assemble the final matrix and clean up */ 530eb9c0419SKris Buschelman ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 531eb9c0419SKris Buschelman ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 532eb9c0419SKris Buschelman ierr = PetscFree(apa);CHKERRQ(ierr); 533eb9c0419SKris Buschelman ierr = PetscLogFlops(flops);CHKERRQ(ierr); 534eb9c0419SKris Buschelman 535eb9c0419SKris Buschelman PetscFunctionReturn(0); 536eb9c0419SKris Buschelman } 537