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__ 73ff134f7aSHong Zhang #define __FUNCT__ "MatPtAP_MPIAIJ_MPIAIJ" 74ff134f7aSHong Zhang PetscErrorCode MatPtAP_MPIAIJ_MPIAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C) 75ff134f7aSHong Zhang { 76*170ef064SHong Zhang #ifdef TMP 77ff134f7aSHong Zhang PetscErrorCode ierr; 78ff134f7aSHong Zhang Mat C_mpi,AP_seq,P_seq,P_subseq,*psubseq; 79ff134f7aSHong Zhang Mat_MPIAIJ *p = (Mat_MPIAIJ*)P->data; 80ff134f7aSHong Zhang Mat_MatMatMultMPI *mult; 81ff134f7aSHong Zhang int i,prow,prstart,prend,m=P->m,pncols; 82ff134f7aSHong Zhang const int *pcols; 83ff134f7aSHong Zhang const PetscScalar *pvals; 84ff134f7aSHong Zhang int rank; 85ff134f7aSHong Zhang 86ff134f7aSHong Zhang PetscFunctionBegin; 87ff134f7aSHong Zhang ierr = MPI_Comm_rank(A->comm,&rank);CHKERRQ(ierr); 88ff134f7aSHong Zhang 89ff134f7aSHong Zhang ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ(A,P,fill,&C_mpi);CHKERRQ(ierr); 90ff134f7aSHong Zhang mult = (Mat_MatMatMultMPI*)C_mpi->spptr; 91ff134f7aSHong Zhang P_seq = mult->bseq[0]; 92ff134f7aSHong Zhang AP_seq = mult->C_seq; 93ff134f7aSHong Zhang prstart = mult->brstart; 94ff134f7aSHong Zhang prend = prstart + m; 95ff134f7aSHong 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); 96ff134f7aSHong Zhang 97ff134f7aSHong Zhang /* get seq matrix P_subseq by taking local rows of P */ 98ff134f7aSHong Zhang IS isrow,iscol; 99ff134f7aSHong Zhang ierr = ISCreateStride(PETSC_COMM_SELF,m,prstart,1,&isrow);CHKERRQ(ierr); 100ff134f7aSHong Zhang ierr = ISCreateStride(PETSC_COMM_SELF,P_seq->n,0,1,&iscol);CHKERRQ(ierr); 101ff134f7aSHong Zhang ierr = MatGetSubMatrices(P_seq,1,&isrow,&iscol,MAT_INITIAL_MATRIX,&psubseq);CHKERRQ(ierr); 102ff134f7aSHong Zhang P_subseq = psubseq[0]; 103ff134f7aSHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," [%d] dim of P_subseq: %d, %d\n",rank,P_subseq->m,P_subseq->n); 104ff134f7aSHong Zhang 105ff134f7aSHong Zhang /* Compute P_subseq^T*C_seq using outer product (P_loc^T)[:,i]*C_seq[i,:]. */ 106ff134f7aSHong Zhang for (i=0;i<m;i++) { 107ff134f7aSHong Zhang prow = prstart + i; 108ff134f7aSHong Zhang ierr = MatGetRow(P_seq,prow,&pncols,&pcols,&pvals);CHKERRQ(ierr); 109ff134f7aSHong Zhang ierr = MatRestoreRow(P_seq,prow,&pncols,&pcols,&pvals);CHKERRQ(ierr); 110ff134f7aSHong Zhang } 111ff134f7aSHong Zhang 112ff134f7aSHong Zhang *C = C_mpi; /* to be removed! */ 113*170ef064SHong Zhang #endif /* TMP */ 114ff134f7aSHong Zhang PetscFunctionReturn(0); 115ff134f7aSHong Zhang } 116ff134f7aSHong Zhang 117ff134f7aSHong Zhang #undef __FUNCT__ 118ff134f7aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ" 119ff134f7aSHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat B,PetscReal fill,Mat *C) 120ff134f7aSHong Zhang { 121ff134f7aSHong Zhang Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data; 122ff134f7aSHong Zhang PetscErrorCode ierr; 123ff134f7aSHong Zhang 124ff134f7aSHong Zhang PetscFunctionBegin; 125ff134f7aSHong Zhang PetscFunctionReturn(0); 126ff134f7aSHong Zhang } 127ff134f7aSHong Zhang 128ff134f7aSHong Zhang #undef __FUNCT__ 129ff134f7aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ" 130ff134f7aSHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat B,Mat C) 131ff134f7aSHong Zhang { 132ff134f7aSHong Zhang Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data; 133ff134f7aSHong Zhang PetscErrorCode ierr; 134ff134f7aSHong Zhang 135ff134f7aSHong Zhang PetscFunctionBegin; 136ff134f7aSHong Zhang PetscFunctionReturn(0); 137ff134f7aSHong Zhang } 138ff134f7aSHong Zhang 139ff134f7aSHong Zhang #undef __FUNCT__ 1409af31e4aSHong Zhang #define __FUNCT__ "MatPtAP_SeqAIJ_SeqAIJ" 141dfbe8321SBarry Smith PetscErrorCode MatPtAP_SeqAIJ_SeqAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C) 1429af31e4aSHong Zhang { 143dfbe8321SBarry Smith PetscErrorCode ierr; 1449af31e4aSHong Zhang PetscFunctionBegin; 1459af31e4aSHong Zhang if (scall == MAT_INITIAL_MATRIX){ 1469af31e4aSHong Zhang ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ(A,P,fill,C);CHKERRQ(ierr); 1479af31e4aSHong Zhang } 1489af31e4aSHong Zhang ierr = MatPtAPNumeric_SeqAIJ_SeqAIJ(A,P,*C);CHKERRQ(ierr); 1499af31e4aSHong Zhang PetscFunctionReturn(0); 1509af31e4aSHong Zhang } 1519af31e4aSHong Zhang 1529af31e4aSHong Zhang #undef __FUNCT__ 1539af31e4aSHong Zhang #define __FUNCT__ "MatPtAPSymbolic" 1546849ba73SBarry Smith /* 1559af31e4aSHong Zhang MatPtAPSymbolic - Creates the (i,j) structure of the matrix projection C = P^T * A * P 1564d3841fdSKris Buschelman 1574d3841fdSKris Buschelman Collective on Mat 1584d3841fdSKris Buschelman 1594d3841fdSKris Buschelman Input Parameters: 1604d3841fdSKris Buschelman + A - the matrix 1614d3841fdSKris Buschelman - P - the projection matrix 1624d3841fdSKris Buschelman 1634d3841fdSKris Buschelman Output Parameters: 1644d3841fdSKris Buschelman . C - the (i,j) structure of the product matrix 1654d3841fdSKris Buschelman 1664d3841fdSKris Buschelman Notes: 1674d3841fdSKris Buschelman C will be created and must be destroyed by the user with MatDestroy(). 1684d3841fdSKris Buschelman 1694d3841fdSKris Buschelman This routine is currently only implemented for pairs of SeqAIJ matrices and classes 1704d3841fdSKris Buschelman which inherit from SeqAIJ. C will be of type MATSEQAIJ. The product is computed using 1719af31e4aSHong Zhang this (i,j) structure by calling MatPtAPNumeric(). 1724d3841fdSKris Buschelman 1734d3841fdSKris Buschelman Level: intermediate 1744d3841fdSKris Buschelman 1759af31e4aSHong Zhang .seealso: MatPtAP(),MatPtAPNumeric(),MatMatMultSymbolic() 1766849ba73SBarry Smith */ 177dfbe8321SBarry Smith PetscErrorCode MatPtAPSymbolic(Mat A,Mat P,PetscReal fill,Mat *C) { 178dfbe8321SBarry Smith PetscErrorCode ierr; 179534c1384SKris Buschelman PetscErrorCode (*fA)(Mat,Mat,PetscReal,Mat*); 180534c1384SKris Buschelman PetscErrorCode (*fP)(Mat,Mat,PetscReal,Mat*); 181eb9c0419SKris Buschelman 182eb9c0419SKris Buschelman PetscFunctionBegin; 183eb9c0419SKris Buschelman 1844482741eSBarry Smith PetscValidHeaderSpecific(A,MAT_COOKIE,1); 185c9780b6fSBarry Smith PetscValidType(A,1); 186eb9c0419SKris Buschelman MatPreallocated(A); 187eb9c0419SKris Buschelman if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 188eb9c0419SKris Buschelman if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 189eb9c0419SKris Buschelman 1904482741eSBarry Smith PetscValidHeaderSpecific(P,MAT_COOKIE,2); 191c9780b6fSBarry Smith PetscValidType(P,2); 192eb9c0419SKris Buschelman MatPreallocated(P); 193eb9c0419SKris Buschelman if (!P->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 194eb9c0419SKris Buschelman if (P->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 195eb9c0419SKris Buschelman 1964482741eSBarry Smith PetscValidPointer(C,3); 1974482741eSBarry Smith 198eb9c0419SKris Buschelman if (P->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->M,A->N); 199eb9c0419SKris Buschelman if (A->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix 'A' must be square, %d != %d",A->M,A->N); 200eb9c0419SKris Buschelman 201534c1384SKris Buschelman /* For now, we do not dispatch based on the type of A and P */ 202534c1384SKris Buschelman /* When implementations like _SeqAIJ_MAIJ exist, attack the multiple dispatch problem. */ 203534c1384SKris Buschelman fA = A->ops->ptapsymbolic; 204534c1384SKris Buschelman if (!fA) SETERRQ1(PETSC_ERR_SUP,"MatPtAPSymbolic not supported for A of type %s",A->type_name); 205534c1384SKris Buschelman fP = P->ops->ptapsymbolic; 206534c1384SKris Buschelman if (!fP) SETERRQ1(PETSC_ERR_SUP,"MatPtAPSymbolic not supported for P of type %s",P->type_name); 207534c1384SKris 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); 2084d3841fdSKris Buschelman 209534c1384SKris Buschelman ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 210534c1384SKris Buschelman ierr = (*fA)(A,P,fill,C);CHKERRQ(ierr); 211534c1384SKris Buschelman ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 212eb9c0419SKris Buschelman 213eb9c0419SKris Buschelman PetscFunctionReturn(0); 214eb9c0419SKris Buschelman } 215eb9c0419SKris Buschelman 216f747e1aeSHong Zhang typedef struct { 217f747e1aeSHong Zhang Mat symAP; 218f747e1aeSHong Zhang } Mat_PtAPstruct; 219f747e1aeSHong Zhang 22078a80504SBarry Smith EXTERN PetscErrorCode MatDestroy_SeqAIJ(Mat); 22178a80504SBarry Smith 222f747e1aeSHong Zhang #undef __FUNCT__ 223f747e1aeSHong Zhang #define __FUNCT__ "MatDestroy_SeqAIJ_PtAP" 224f4a850bbSBarry Smith PetscErrorCode MatDestroy_SeqAIJ_PtAP(Mat A) 225f747e1aeSHong Zhang { 226f4a850bbSBarry Smith PetscErrorCode ierr; 227f747e1aeSHong Zhang Mat_PtAPstruct *ptap=(Mat_PtAPstruct*)A->spptr; 228f747e1aeSHong Zhang 229f747e1aeSHong Zhang PetscFunctionBegin; 230f747e1aeSHong Zhang ierr = MatDestroy(ptap->symAP);CHKERRQ(ierr); 231f747e1aeSHong Zhang ierr = PetscFree(ptap);CHKERRQ(ierr); 23278a80504SBarry Smith ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr); 233f747e1aeSHong Zhang PetscFunctionReturn(0); 234f747e1aeSHong Zhang } 235f747e1aeSHong Zhang 236eb9c0419SKris Buschelman #undef __FUNCT__ 2379af31e4aSHong Zhang #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqAIJ" 238dfbe8321SBarry Smith PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat P,PetscReal fill,Mat *C) { 239dfbe8321SBarry Smith PetscErrorCode ierr; 240f4a850bbSBarry Smith int *pti,*ptj; 241f747e1aeSHong Zhang Mat Pt,AP; 242f747e1aeSHong Zhang Mat_PtAPstruct *ptap; 243eb9c0419SKris Buschelman 244eb9c0419SKris Buschelman PetscFunctionBegin; 245f747e1aeSHong Zhang /* create symbolic Pt */ 246eb9c0419SKris Buschelman ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr); 247f747e1aeSHong Zhang ierr = MatCreateSeqAIJWithArrays(P->comm,P->N,P->M,pti,ptj,PETSC_NULL,&Pt);CHKERRQ(ierr); 248eb9c0419SKris Buschelman 249f747e1aeSHong Zhang /* get symbolic AP=A*P and C=Pt*AP */ 250f747e1aeSHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,P,fill,&AP);CHKERRQ(ierr); 251f747e1aeSHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(Pt,AP,fill,C);CHKERRQ(ierr); 252eb9c0419SKris Buschelman 253f747e1aeSHong Zhang /* clean up */ 254f747e1aeSHong Zhang ierr = MatRestoreSymbolicTranspose_SeqAIJ(Pt,&pti,&ptj);CHKERRQ(ierr); 255f747e1aeSHong Zhang ierr = MatDestroy(Pt);CHKERRQ(ierr); 256eb9c0419SKris Buschelman 257f747e1aeSHong Zhang /* save symbolic AP - to be used by MatPtAPNumeric_SeqAIJ_SeqAIJ() */ 258f747e1aeSHong Zhang ierr = PetscNew(Mat_PtAPstruct,&ptap);CHKERRQ(ierr); 259f747e1aeSHong Zhang ptap->symAP = AP; 260f747e1aeSHong Zhang (*C)->spptr = (void*)ptap; 261f747e1aeSHong Zhang (*C)->ops->destroy = MatDestroy_SeqAIJ_PtAP; 262eb9c0419SKris Buschelman 263eb9c0419SKris Buschelman PetscFunctionReturn(0); 264eb9c0419SKris Buschelman } 265eb9c0419SKris Buschelman 2663985e5eaSKris Buschelman #include "src/mat/impls/maij/maij.h" 2673985e5eaSKris Buschelman EXTERN_C_BEGIN 2683985e5eaSKris Buschelman #undef __FUNCT__ 2699af31e4aSHong Zhang #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqMAIJ" 270dfbe8321SBarry Smith PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqMAIJ(Mat A,Mat PP,Mat *C) { 2715c66b693SKris Buschelman /* This routine requires testing -- I don't think it works. */ 272dfbe8321SBarry Smith PetscErrorCode ierr; 2733985e5eaSKris Buschelman FreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 2743985e5eaSKris Buschelman Mat_SeqMAIJ *pp=(Mat_SeqMAIJ*)PP->data; 2753985e5eaSKris Buschelman Mat P=pp->AIJ; 2763985e5eaSKris Buschelman Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*p=(Mat_SeqAIJ*)P->data,*c; 2773985e5eaSKris Buschelman int *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj; 2783985e5eaSKris Buschelman int *ci,*cj,*denserow,*sparserow,*ptadenserow,*ptasparserow,*ptaj; 2793985e5eaSKris Buschelman int an=A->N,am=A->M,pn=P->N,pm=P->M,ppdof=pp->dof; 280fe05a634SKris Buschelman int i,j,k,dof,pdof,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi; 2813985e5eaSKris Buschelman MatScalar *ca; 2823985e5eaSKris Buschelman 2833985e5eaSKris Buschelman PetscFunctionBegin; 2843985e5eaSKris Buschelman /* Start timer */ 2859af31e4aSHong Zhang ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,PP,0,0);CHKERRQ(ierr); 2863985e5eaSKris Buschelman 2873985e5eaSKris Buschelman /* Get ij structure of P^T */ 2883985e5eaSKris Buschelman ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr); 2893985e5eaSKris Buschelman 2903985e5eaSKris Buschelman /* Allocate ci array, arrays for fill computation and */ 2913985e5eaSKris Buschelman /* free space for accumulating nonzero column info */ 2923985e5eaSKris Buschelman ierr = PetscMalloc((pn+1)*sizeof(int),&ci);CHKERRQ(ierr); 2933985e5eaSKris Buschelman ci[0] = 0; 2943985e5eaSKris Buschelman 2953985e5eaSKris Buschelman ierr = PetscMalloc((2*pn+2*an+1)*sizeof(int),&ptadenserow);CHKERRQ(ierr); 2963985e5eaSKris Buschelman ierr = PetscMemzero(ptadenserow,(2*pn+2*an+1)*sizeof(int));CHKERRQ(ierr); 2973985e5eaSKris Buschelman ptasparserow = ptadenserow + an; 2983985e5eaSKris Buschelman denserow = ptasparserow + an; 2993985e5eaSKris Buschelman sparserow = denserow + pn; 3003985e5eaSKris Buschelman 3013985e5eaSKris Buschelman /* Set initial free space to be nnz(A) scaled by aspect ratio of P. */ 3023985e5eaSKris Buschelman /* This should be reasonable if sparsity of PtAP is similar to that of A. */ 3033985e5eaSKris Buschelman ierr = GetMoreSpace((ai[am]/pm)*pn,&free_space); 3043985e5eaSKris Buschelman current_space = free_space; 3053985e5eaSKris Buschelman 3063985e5eaSKris Buschelman /* Determine symbolic info for each row of C: */ 3073985e5eaSKris Buschelman for (i=0;i<pn/ppdof;i++) { 3083985e5eaSKris Buschelman ptnzi = pti[i+1] - pti[i]; 3093985e5eaSKris Buschelman ptanzi = 0; 3103985e5eaSKris Buschelman ptJ = ptj + pti[i]; 3113985e5eaSKris Buschelman for (dof=0;dof<ppdof;dof++) { 3123985e5eaSKris Buschelman /* Determine symbolic row of PtA: */ 3133985e5eaSKris Buschelman for (j=0;j<ptnzi;j++) { 3143985e5eaSKris Buschelman arow = ptJ[j] + dof; 3153985e5eaSKris Buschelman anzj = ai[arow+1] - ai[arow]; 3163985e5eaSKris Buschelman ajj = aj + ai[arow]; 3173985e5eaSKris Buschelman for (k=0;k<anzj;k++) { 3183985e5eaSKris Buschelman if (!ptadenserow[ajj[k]]) { 3193985e5eaSKris Buschelman ptadenserow[ajj[k]] = -1; 3203985e5eaSKris Buschelman ptasparserow[ptanzi++] = ajj[k]; 3213985e5eaSKris Buschelman } 3223985e5eaSKris Buschelman } 3233985e5eaSKris Buschelman } 3243985e5eaSKris Buschelman /* Using symbolic info for row of PtA, determine symbolic info for row of C: */ 3253985e5eaSKris Buschelman ptaj = ptasparserow; 3263985e5eaSKris Buschelman cnzi = 0; 3273985e5eaSKris Buschelman for (j=0;j<ptanzi;j++) { 328fe05a634SKris Buschelman pdof = *ptaj%dof; 3293985e5eaSKris Buschelman prow = (*ptaj++)/dof; 3303985e5eaSKris Buschelman pnzj = pi[prow+1] - pi[prow]; 3313985e5eaSKris Buschelman pjj = pj + pi[prow]; 3323985e5eaSKris Buschelman for (k=0;k<pnzj;k++) { 333fe05a634SKris Buschelman if (!denserow[pjj[k]+pdof]) { 334fe05a634SKris Buschelman denserow[pjj[k]+pdof] = -1; 335fe05a634SKris Buschelman sparserow[cnzi++] = pjj[k]+pdof; 3363985e5eaSKris Buschelman } 3373985e5eaSKris Buschelman } 3383985e5eaSKris Buschelman } 3393985e5eaSKris Buschelman 3403985e5eaSKris Buschelman /* sort sparserow */ 3413985e5eaSKris Buschelman ierr = PetscSortInt(cnzi,sparserow);CHKERRQ(ierr); 3423985e5eaSKris Buschelman 3433985e5eaSKris Buschelman /* If free space is not available, make more free space */ 3443985e5eaSKris Buschelman /* Double the amount of total space in the list */ 3453985e5eaSKris Buschelman if (current_space->local_remaining<cnzi) { 3463985e5eaSKris Buschelman ierr = GetMoreSpace(current_space->total_array_size,¤t_space);CHKERRQ(ierr); 3473985e5eaSKris Buschelman } 3483985e5eaSKris Buschelman 3493985e5eaSKris Buschelman /* Copy data into free space, and zero out denserows */ 3503985e5eaSKris Buschelman ierr = PetscMemcpy(current_space->array,sparserow,cnzi*sizeof(int));CHKERRQ(ierr); 3513985e5eaSKris Buschelman current_space->array += cnzi; 3523985e5eaSKris Buschelman current_space->local_used += cnzi; 3533985e5eaSKris Buschelman current_space->local_remaining -= cnzi; 3543985e5eaSKris Buschelman 3553985e5eaSKris Buschelman for (j=0;j<ptanzi;j++) { 3563985e5eaSKris Buschelman ptadenserow[ptasparserow[j]] = 0; 3573985e5eaSKris Buschelman } 3583985e5eaSKris Buschelman for (j=0;j<cnzi;j++) { 3593985e5eaSKris Buschelman denserow[sparserow[j]] = 0; 3603985e5eaSKris Buschelman } 3613985e5eaSKris Buschelman /* Aside: Perhaps we should save the pta info for the numerical factorization. */ 3623985e5eaSKris Buschelman /* For now, we will recompute what is needed. */ 3633985e5eaSKris Buschelman ci[i+1+dof] = ci[i+dof] + cnzi; 3643985e5eaSKris Buschelman } 3653985e5eaSKris Buschelman } 3663985e5eaSKris Buschelman /* nnz is now stored in ci[ptm], column indices are in the list of free space */ 3673985e5eaSKris Buschelman /* Allocate space for cj, initialize cj, and */ 3683985e5eaSKris Buschelman /* destroy list of free space and other temporary array(s) */ 3693985e5eaSKris Buschelman ierr = PetscMalloc((ci[pn]+1)*sizeof(int),&cj);CHKERRQ(ierr); 3703985e5eaSKris Buschelman ierr = MakeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 3713985e5eaSKris Buschelman ierr = PetscFree(ptadenserow);CHKERRQ(ierr); 3723985e5eaSKris Buschelman 3733985e5eaSKris Buschelman /* Allocate space for ca */ 3743985e5eaSKris Buschelman ierr = PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 3753985e5eaSKris Buschelman ierr = PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));CHKERRQ(ierr); 3763985e5eaSKris Buschelman 3773985e5eaSKris Buschelman /* put together the new matrix */ 3783985e5eaSKris Buschelman ierr = MatCreateSeqAIJWithArrays(A->comm,pn,pn,ci,cj,ca,C);CHKERRQ(ierr); 3793985e5eaSKris Buschelman 3803985e5eaSKris Buschelman /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 3813985e5eaSKris Buschelman /* Since these are PETSc arrays, change flags to free them as necessary. */ 3823985e5eaSKris Buschelman c = (Mat_SeqAIJ *)((*C)->data); 3833985e5eaSKris Buschelman c->freedata = PETSC_TRUE; 3843985e5eaSKris Buschelman c->nonew = 0; 3853985e5eaSKris Buschelman 3863985e5eaSKris Buschelman /* Clean up. */ 3873985e5eaSKris Buschelman ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr); 3883985e5eaSKris Buschelman 3899af31e4aSHong Zhang ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,PP,0,0);CHKERRQ(ierr); 3903985e5eaSKris Buschelman PetscFunctionReturn(0); 3913985e5eaSKris Buschelman } 3923985e5eaSKris Buschelman EXTERN_C_END 3933985e5eaSKris Buschelman 394eb9c0419SKris Buschelman #undef __FUNCT__ 3959af31e4aSHong Zhang #define __FUNCT__ "MatPtAPNumeric" 3966849ba73SBarry Smith /* 3979af31e4aSHong Zhang MatPtAPNumeric - Computes the matrix projection C = P^T * A * P 3984d3841fdSKris Buschelman 3994d3841fdSKris Buschelman Collective on Mat 4004d3841fdSKris Buschelman 4014d3841fdSKris Buschelman Input Parameters: 4024d3841fdSKris Buschelman + A - the matrix 4034d3841fdSKris Buschelman - P - the projection matrix 4044d3841fdSKris Buschelman 4054d3841fdSKris Buschelman Output Parameters: 4064d3841fdSKris Buschelman . C - the product matrix 4074d3841fdSKris Buschelman 4084d3841fdSKris Buschelman Notes: 4099af31e4aSHong Zhang C must have been created by calling MatPtAPSymbolic and must be destroyed by 4104d3841fdSKris Buschelman the user using MatDeatroy(). 4114d3841fdSKris Buschelman 412*170ef064SHong Zhang This routine is currently only implemented for pairs of AIJ matrices and classes 413*170ef064SHong Zhang which inherit from AIJ. C will be of type MATAIJ. 4144d3841fdSKris Buschelman 4154d3841fdSKris Buschelman Level: intermediate 4164d3841fdSKris Buschelman 4179af31e4aSHong Zhang .seealso: MatPtAP(),MatPtAPSymbolic(),MatMatMultNumeric() 4186849ba73SBarry Smith */ 419dfbe8321SBarry Smith PetscErrorCode MatPtAPNumeric(Mat A,Mat P,Mat C) { 420dfbe8321SBarry Smith PetscErrorCode ierr; 421534c1384SKris Buschelman PetscErrorCode (*fA)(Mat,Mat,Mat); 422534c1384SKris Buschelman PetscErrorCode (*fP)(Mat,Mat,Mat); 423eb9c0419SKris Buschelman 424eb9c0419SKris Buschelman PetscFunctionBegin; 425eb9c0419SKris Buschelman 4264482741eSBarry Smith PetscValidHeaderSpecific(A,MAT_COOKIE,1); 427c9780b6fSBarry Smith PetscValidType(A,1); 428eb9c0419SKris Buschelman MatPreallocated(A); 429eb9c0419SKris Buschelman if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 430eb9c0419SKris Buschelman if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 431eb9c0419SKris Buschelman 4324482741eSBarry Smith PetscValidHeaderSpecific(P,MAT_COOKIE,2); 433c9780b6fSBarry Smith PetscValidType(P,2); 434eb9c0419SKris Buschelman MatPreallocated(P); 435eb9c0419SKris Buschelman if (!P->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 436eb9c0419SKris Buschelman if (P->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 437eb9c0419SKris Buschelman 4384482741eSBarry Smith PetscValidHeaderSpecific(C,MAT_COOKIE,3); 439c9780b6fSBarry Smith PetscValidType(C,3); 440eb9c0419SKris Buschelman MatPreallocated(C); 441eb9c0419SKris Buschelman if (!C->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 442eb9c0419SKris Buschelman if (C->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 443eb9c0419SKris Buschelman 444eb9c0419SKris Buschelman if (P->N!=C->M) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->N,C->M); 445eb9c0419SKris Buschelman if (P->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->M,A->N); 446eb9c0419SKris Buschelman if (A->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix 'A' must be square, %d != %d",A->M,A->N); 447eb9c0419SKris Buschelman if (P->N!=C->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->N,C->N); 448eb9c0419SKris Buschelman 449534c1384SKris Buschelman /* For now, we do not dispatch based on the type of A and P */ 450534c1384SKris Buschelman /* When implementations like _SeqAIJ_MAIJ exist, attack the multiple dispatch problem. */ 451534c1384SKris Buschelman fA = A->ops->ptapnumeric; 452534c1384SKris Buschelman if (!fA) SETERRQ1(PETSC_ERR_SUP,"MatPtAPNumeric not supported for A of type %s",A->type_name); 453534c1384SKris Buschelman fP = P->ops->ptapnumeric; 454534c1384SKris Buschelman if (!fP) SETERRQ1(PETSC_ERR_SUP,"MatPtAPNumeric not supported for P of type %s",P->type_name); 455534c1384SKris 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); 4564d3841fdSKris Buschelman 457534c1384SKris Buschelman ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 458534c1384SKris Buschelman ierr = (*fA)(A,P,C);CHKERRQ(ierr); 459534c1384SKris Buschelman ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 460eb9c0419SKris Buschelman 461eb9c0419SKris Buschelman PetscFunctionReturn(0); 462eb9c0419SKris Buschelman } 463eb9c0419SKris Buschelman 464eb9c0419SKris Buschelman #undef __FUNCT__ 4659af31e4aSHong Zhang #define __FUNCT__ "MatPtAPNumeric_SeqAIJ_SeqAIJ" 466dfbe8321SBarry Smith PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat A,Mat P,Mat C) 467dfbe8321SBarry Smith { 468dfbe8321SBarry Smith PetscErrorCode ierr; 469f747e1aeSHong Zhang Mat_PtAPstruct *ptap=(Mat_PtAPstruct*)C->spptr; 470*170ef064SHong Zhang Mat AP=ptap->symAP; 471eb9c0419SKris Buschelman 472eb9c0419SKris Buschelman PetscFunctionBegin; 473*170ef064SHong Zhang /* compute numeric AP = A*P */ 474*170ef064SHong Zhang ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ(A,P,AP);CHKERRQ(ierr); 475eb9c0419SKris Buschelman 476*170ef064SHong Zhang /* compute numeric P^T*AP */ 477*170ef064SHong Zhang ierr = MatMatMultTransposeNumeric_SeqAIJ_SeqAIJ(P,AP,C);CHKERRQ(ierr); 478eb9c0419SKris Buschelman PetscFunctionReturn(0); 479eb9c0419SKris Buschelman } 480