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" 9*55a3bba9SHong Zhang #include "petscbt.h" 10eb9c0419SKris Buschelman 11eb9c0419SKris Buschelman #undef __FUNCT__ 129af31e4aSHong Zhang #define __FUNCT__ "MatPtAP" 134d3841fdSKris Buschelman /*@ 149af31e4aSHong Zhang MatPtAP - Creates the matrix projection C = P^T * A * P 154d3841fdSKris Buschelman 164d3841fdSKris Buschelman Collective on Mat 174d3841fdSKris Buschelman 184d3841fdSKris Buschelman Input Parameters: 194d3841fdSKris Buschelman + A - the matrix 20f747e1aeSHong Zhang . P - the projection matrix 21f747e1aeSHong Zhang . scall - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX 22f747e1aeSHong Zhang - fill - expected fill as ratio of nnz(C)/(nnz(A) + nnz(P)) 234d3841fdSKris Buschelman 244d3841fdSKris Buschelman Output Parameters: 254d3841fdSKris Buschelman . C - the product matrix 264d3841fdSKris Buschelman 274d3841fdSKris Buschelman Notes: 284d3841fdSKris Buschelman C will be created and must be destroyed by the user with MatDestroy(). 294d3841fdSKris Buschelman 304d3841fdSKris Buschelman This routine is currently only implemented for pairs of SeqAIJ matrices and classes 314d3841fdSKris Buschelman which inherit from SeqAIJ. C will be of type MATSEQAIJ. 324d3841fdSKris Buschelman 334d3841fdSKris Buschelman Level: intermediate 344d3841fdSKris Buschelman 359af31e4aSHong Zhang .seealso: MatPtAPSymbolic(),MatPtAPNumeric(),MatMatMult() 364d3841fdSKris Buschelman @*/ 37*55a3bba9SHong Zhang PetscErrorCode MatPtAP(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C) 38*55a3bba9SHong Zhang { 39dfbe8321SBarry Smith PetscErrorCode ierr; 40534c1384SKris Buschelman PetscErrorCode (*fA)(Mat,Mat,MatReuse,PetscReal,Mat *); 41534c1384SKris Buschelman PetscErrorCode (*fP)(Mat,Mat,MatReuse,PetscReal,Mat *); 42eb9c0419SKris Buschelman 43eb9c0419SKris Buschelman PetscFunctionBegin; 449af31e4aSHong Zhang PetscValidHeaderSpecific(A,MAT_COOKIE,1); 459af31e4aSHong Zhang PetscValidType(A,1); 469af31e4aSHong Zhang MatPreallocated(A); 479af31e4aSHong Zhang if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 489af31e4aSHong Zhang if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 499af31e4aSHong Zhang PetscValidHeaderSpecific(P,MAT_COOKIE,2); 509af31e4aSHong Zhang PetscValidType(P,2); 519af31e4aSHong Zhang MatPreallocated(P); 529af31e4aSHong Zhang if (!P->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 539af31e4aSHong Zhang if (P->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 549af31e4aSHong Zhang PetscValidPointer(C,3); 55534c1384SKris Buschelman 569af31e4aSHong Zhang if (P->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->M,A->N); 57eb9c0419SKris Buschelman 589af31e4aSHong Zhang if (fill <=0.0) SETERRQ1(PETSC_ERR_ARG_SIZ,"fill=%g must be > 0.0",fill); 59eb9c0419SKris Buschelman 60534c1384SKris Buschelman /* For now, we do not dispatch based on the type of A and P */ 61534c1384SKris Buschelman /* When implementations like _SeqAIJ_MAIJ exist, attack the multiple dispatch problem. */ 62534c1384SKris Buschelman fA = A->ops->ptap; 63534c1384SKris Buschelman if (!fA) SETERRQ1(PETSC_ERR_SUP,"MatPtAP not supported for A of type %s",A->type_name); 64534c1384SKris Buschelman fP = P->ops->ptap; 65534c1384SKris Buschelman if (!fP) SETERRQ1(PETSC_ERR_SUP,"MatPtAP not supported for P of type %s",P->type_name); 66534c1384SKris 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); 67534c1384SKris Buschelman 689af31e4aSHong Zhang ierr = PetscLogEventBegin(MAT_PtAP,A,P,0,0);CHKERRQ(ierr); 69534c1384SKris Buschelman ierr = (*fA)(A,P,scall,fill,C);CHKERRQ(ierr); 709af31e4aSHong Zhang ierr = PetscLogEventEnd(MAT_PtAP,A,P,0,0);CHKERRQ(ierr); 71eb9c0419SKris Buschelman PetscFunctionReturn(0); 72eb9c0419SKris Buschelman } 73eb9c0419SKris Buschelman 740e36024fSHong Zhang EXTERN PetscErrorCode MatPtAP_SeqAIJ_SeqAIJ_ReducedPt(Mat,Mat,MatReuse,PetscReal,int,int,Mat*); 750e36024fSHong Zhang EXTERN PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_ReducedPt(Mat,Mat,PetscReal,int,int,Mat*); 760e36024fSHong Zhang EXTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_ReducedPt(Mat,Mat,int,int,Mat); 77b90dcfe3SHong Zhang 78eb9c0419SKris Buschelman #undef __FUNCT__ 79ff134f7aSHong Zhang #define __FUNCT__ "MatPtAP_MPIAIJ_MPIAIJ" 80ff134f7aSHong Zhang PetscErrorCode MatPtAP_MPIAIJ_MPIAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C) 81ff134f7aSHong Zhang { 82ff134f7aSHong Zhang PetscErrorCode ierr; 83b90dcfe3SHong Zhang 84b90dcfe3SHong Zhang PetscFunctionBegin; 85b90dcfe3SHong Zhang if (scall == MAT_INITIAL_MATRIX){ 864c627768SHong Zhang ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 87b90dcfe3SHong Zhang ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ(A,P,fill,C);CHKERRQ(ierr);/* numeric product is computed as well */ 884c627768SHong Zhang ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 89b90dcfe3SHong Zhang } else if (scall == MAT_REUSE_MATRIX){ 904c627768SHong Zhang ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 91b90dcfe3SHong Zhang ierr = MatPtAPNumeric_MPIAIJ_MPIAIJ(A,P,*C);CHKERRQ(ierr); 924c627768SHong Zhang ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 93b90dcfe3SHong Zhang } else { 94b90dcfe3SHong Zhang SETERRQ1(PETSC_ERR_ARG_WRONG,"Invalid MatReuse %d",scall); 95b90dcfe3SHong Zhang } 96b90dcfe3SHong Zhang PetscFunctionReturn(0); 97b90dcfe3SHong Zhang } 98b90dcfe3SHong Zhang 99b90dcfe3SHong Zhang #undef __FUNCT__ 100b90dcfe3SHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ" 101b90dcfe3SHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C) 102b90dcfe3SHong Zhang { 103b90dcfe3SHong Zhang PetscErrorCode ierr; 10425616d81SHong Zhang Mat P_seq,A_loc,C_seq; 1050e36024fSHong Zhang int prstart,prend,m=P->m; 10625616d81SHong Zhang IS isrowp,iscolp; 107ff134f7aSHong Zhang 108ff134f7aSHong Zhang PetscFunctionBegin; 10925616d81SHong Zhang /* get P_seq = submatrix of P by taking rows of P that equal to nonzero col of A */ 110b90dcfe3SHong Zhang ierr = MatGetBrowsOfAcols(A,P,MAT_INITIAL_MATRIX,&isrowp,&iscolp,&prstart,&P_seq);CHKERRQ(ierr); 11125616d81SHong Zhang ierr = ISDestroy(iscolp);CHKERRQ(ierr); 112ff134f7aSHong Zhang 11325616d81SHong Zhang /* get A_loc = submatrix of A by taking all local rows of A */ 114b90dcfe3SHong Zhang ierr = MatGetLocalMat(A,MAT_INITIAL_MATRIX,PETSC_NULL,&isrowp,&A_loc);CHKERRQ(ierr); 11525616d81SHong Zhang ierr = ISDestroy(isrowp);CHKERRQ(ierr); 1160e36024fSHong Zhang 11725616d81SHong Zhang /* compute C_seq = P_loc^T * A_loc * P_seq */ 118ff134f7aSHong Zhang prend = prstart + m; 119b90dcfe3SHong Zhang ierr = MatPtAP_SeqAIJ_SeqAIJ_ReducedPt(A_loc,P_seq,MAT_INITIAL_MATRIX,fill,prstart,prend,&C_seq);CHKERRQ(ierr); 12025616d81SHong Zhang ierr = MatDestroy(P_seq);CHKERRQ(ierr); 12125616d81SHong Zhang ierr = MatDestroy(A_loc);CHKERRQ(ierr); 122b90dcfe3SHong Zhang 123b90dcfe3SHong Zhang /* add C_seq into mpi C */ 124b90dcfe3SHong Zhang ierr = MatMerge_SeqsToMPI(A->comm,C_seq,P->n,P->n,MAT_INITIAL_MATRIX,C);CHKERRQ(ierr); 125b90dcfe3SHong Zhang 126ff134f7aSHong Zhang PetscFunctionReturn(0); 127ff134f7aSHong Zhang } 128ff134f7aSHong Zhang 129ff134f7aSHong Zhang #undef __FUNCT__ 130ff134f7aSHong Zhang #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ" 131b90dcfe3SHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C) 132ff134f7aSHong Zhang { 133b90dcfe3SHong Zhang PetscErrorCode ierr; 134b90dcfe3SHong Zhang Mat P_seq,A_loc,C_seq; 135b90dcfe3SHong Zhang int prstart,prend,m=P->m; 136b90dcfe3SHong Zhang IS isrowp,iscolp; 137671beff6SHong Zhang Mat_Merge_SeqsToMPI *merge; 138671beff6SHong Zhang PetscObjectContainer container; 139ff134f7aSHong Zhang 140ff134f7aSHong Zhang PetscFunctionBegin; 141671beff6SHong Zhang ierr = PetscObjectQuery((PetscObject)C,"MatMergeSeqsToMPI",(PetscObject *)&container);CHKERRQ(ierr); 142671beff6SHong Zhang if (container) { 1437f79fd58SMatthew Knepley ierr = PetscObjectContainerGetPointer(container,(void **)&merge);CHKERRQ(ierr); 1447f79fd58SMatthew Knepley } else { 1457f79fd58SMatthew Knepley SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Matrix C does not posses and object container"); 146671beff6SHong Zhang } 147671beff6SHong Zhang 148b90dcfe3SHong Zhang /* get P_seq = submatrix of P by taking rows of P that equal to nonzero col of A */ 149b90dcfe3SHong Zhang ierr = MatGetBrowsOfAcols(A,P,MAT_INITIAL_MATRIX,&isrowp,&iscolp,&prstart,&P_seq);CHKERRQ(ierr); 150b90dcfe3SHong Zhang ierr = ISDestroy(iscolp);CHKERRQ(ierr); 151ff134f7aSHong Zhang 152b90dcfe3SHong Zhang /* get A_loc = submatrix of A by taking all local rows of A */ 153b90dcfe3SHong Zhang ierr = MatGetLocalMat(A,MAT_INITIAL_MATRIX,PETSC_NULL,&isrowp,&A_loc);CHKERRQ(ierr); 154b90dcfe3SHong Zhang ierr = ISDestroy(isrowp);CHKERRQ(ierr); 155ff134f7aSHong Zhang 156b90dcfe3SHong Zhang /* compute C_seq = P_loc^T * A_loc * P_seq */ 157b90dcfe3SHong Zhang prend = prstart + m; 158b90dcfe3SHong Zhang C_seq = merge->C_seq; 159b90dcfe3SHong Zhang ierr = MatPtAP_SeqAIJ_SeqAIJ_ReducedPt(A_loc,P_seq,MAT_REUSE_MATRIX,1.0,prstart,prend,&C_seq);CHKERRQ(ierr); 160b90dcfe3SHong Zhang ierr = MatDestroy(P_seq);CHKERRQ(ierr); 161b90dcfe3SHong Zhang ierr = MatDestroy(A_loc);CHKERRQ(ierr); 162b90dcfe3SHong Zhang 163b90dcfe3SHong Zhang /* add C_seq into mpi C */ 164b90dcfe3SHong Zhang ierr = MatMerge_SeqsToMPI(A->comm,C_seq,P->n,P->n,MAT_REUSE_MATRIX,&C);CHKERRQ(ierr); 165b90dcfe3SHong Zhang 166ff134f7aSHong Zhang PetscFunctionReturn(0); 167ff134f7aSHong Zhang } 168ff134f7aSHong Zhang 169ff134f7aSHong Zhang #undef __FUNCT__ 1709af31e4aSHong Zhang #define __FUNCT__ "MatPtAP_SeqAIJ_SeqAIJ" 171dfbe8321SBarry Smith PetscErrorCode MatPtAP_SeqAIJ_SeqAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C) 1729af31e4aSHong Zhang { 173dfbe8321SBarry Smith PetscErrorCode ierr; 1749af31e4aSHong Zhang PetscFunctionBegin; 1759af31e4aSHong Zhang if (scall == MAT_INITIAL_MATRIX){ 176d20bfe6fSHong Zhang ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 1779af31e4aSHong Zhang ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ(A,P,fill,C);CHKERRQ(ierr); 178d20bfe6fSHong Zhang ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 1799af31e4aSHong Zhang } 180d20bfe6fSHong Zhang ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 1819af31e4aSHong Zhang ierr = MatPtAPNumeric_SeqAIJ_SeqAIJ(A,P,*C);CHKERRQ(ierr); 182d20bfe6fSHong Zhang ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 1839af31e4aSHong Zhang PetscFunctionReturn(0); 1849af31e4aSHong Zhang } 1859af31e4aSHong Zhang 1869af31e4aSHong Zhang #undef __FUNCT__ 1879af31e4aSHong Zhang #define __FUNCT__ "MatPtAPSymbolic" 1886849ba73SBarry Smith /* 1899af31e4aSHong Zhang MatPtAPSymbolic - Creates the (i,j) structure of the matrix projection C = P^T * A * P 1904d3841fdSKris Buschelman 1914d3841fdSKris Buschelman Collective on Mat 1924d3841fdSKris Buschelman 1934d3841fdSKris Buschelman Input Parameters: 1944d3841fdSKris Buschelman + A - the matrix 1954d3841fdSKris Buschelman - P - the projection matrix 1964d3841fdSKris Buschelman 1974d3841fdSKris Buschelman Output Parameters: 1984d3841fdSKris Buschelman . C - the (i,j) structure of the product matrix 1994d3841fdSKris Buschelman 2004d3841fdSKris Buschelman Notes: 2014d3841fdSKris Buschelman C will be created and must be destroyed by the user with MatDestroy(). 2024d3841fdSKris Buschelman 2034d3841fdSKris Buschelman This routine is currently only implemented for pairs of SeqAIJ matrices and classes 2044d3841fdSKris Buschelman which inherit from SeqAIJ. C will be of type MATSEQAIJ. The product is computed using 2059af31e4aSHong Zhang this (i,j) structure by calling MatPtAPNumeric(). 2064d3841fdSKris Buschelman 2074d3841fdSKris Buschelman Level: intermediate 2084d3841fdSKris Buschelman 2099af31e4aSHong Zhang .seealso: MatPtAP(),MatPtAPNumeric(),MatMatMultSymbolic() 2106849ba73SBarry Smith */ 211*55a3bba9SHong Zhang PetscErrorCode MatPtAPSymbolic(Mat A,Mat P,PetscReal fill,Mat *C) 212*55a3bba9SHong Zhang { 213dfbe8321SBarry Smith PetscErrorCode ierr; 214534c1384SKris Buschelman PetscErrorCode (*fA)(Mat,Mat,PetscReal,Mat*); 215534c1384SKris Buschelman PetscErrorCode (*fP)(Mat,Mat,PetscReal,Mat*); 216eb9c0419SKris Buschelman 217eb9c0419SKris Buschelman PetscFunctionBegin; 218eb9c0419SKris Buschelman 2194482741eSBarry Smith PetscValidHeaderSpecific(A,MAT_COOKIE,1); 220c9780b6fSBarry Smith PetscValidType(A,1); 221eb9c0419SKris Buschelman MatPreallocated(A); 222eb9c0419SKris Buschelman if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 223eb9c0419SKris Buschelman if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 224eb9c0419SKris Buschelman 2254482741eSBarry Smith PetscValidHeaderSpecific(P,MAT_COOKIE,2); 226c9780b6fSBarry Smith PetscValidType(P,2); 227eb9c0419SKris Buschelman MatPreallocated(P); 228eb9c0419SKris Buschelman if (!P->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 229eb9c0419SKris Buschelman if (P->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 230eb9c0419SKris Buschelman 2314482741eSBarry Smith PetscValidPointer(C,3); 2324482741eSBarry Smith 233eb9c0419SKris Buschelman if (P->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->M,A->N); 234eb9c0419SKris Buschelman if (A->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix 'A' must be square, %d != %d",A->M,A->N); 235eb9c0419SKris Buschelman 236534c1384SKris Buschelman /* For now, we do not dispatch based on the type of A and P */ 237534c1384SKris Buschelman /* When implementations like _SeqAIJ_MAIJ exist, attack the multiple dispatch problem. */ 238534c1384SKris Buschelman fA = A->ops->ptapsymbolic; 239534c1384SKris Buschelman if (!fA) SETERRQ1(PETSC_ERR_SUP,"MatPtAPSymbolic not supported for A of type %s",A->type_name); 240534c1384SKris Buschelman fP = P->ops->ptapsymbolic; 241534c1384SKris Buschelman if (!fP) SETERRQ1(PETSC_ERR_SUP,"MatPtAPSymbolic not supported for P of type %s",P->type_name); 242534c1384SKris 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); 2434d3841fdSKris Buschelman 244534c1384SKris Buschelman ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 245534c1384SKris Buschelman ierr = (*fA)(A,P,fill,C);CHKERRQ(ierr); 246534c1384SKris Buschelman ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 247eb9c0419SKris Buschelman 248eb9c0419SKris Buschelman PetscFunctionReturn(0); 249eb9c0419SKris Buschelman } 250eb9c0419SKris Buschelman 251f747e1aeSHong Zhang typedef struct { 252f747e1aeSHong Zhang Mat symAP; 253f747e1aeSHong Zhang } Mat_PtAPstruct; 254f747e1aeSHong Zhang 25578a80504SBarry Smith EXTERN PetscErrorCode MatDestroy_SeqAIJ(Mat); 25678a80504SBarry Smith 257f747e1aeSHong Zhang #undef __FUNCT__ 258f747e1aeSHong Zhang #define __FUNCT__ "MatDestroy_SeqAIJ_PtAP" 259f4a850bbSBarry Smith PetscErrorCode MatDestroy_SeqAIJ_PtAP(Mat A) 260f747e1aeSHong Zhang { 261f4a850bbSBarry Smith PetscErrorCode ierr; 262f747e1aeSHong Zhang Mat_PtAPstruct *ptap=(Mat_PtAPstruct*)A->spptr; 263f747e1aeSHong Zhang 264f747e1aeSHong Zhang PetscFunctionBegin; 265f747e1aeSHong Zhang ierr = MatDestroy(ptap->symAP);CHKERRQ(ierr); 266f747e1aeSHong Zhang ierr = PetscFree(ptap);CHKERRQ(ierr); 26778a80504SBarry Smith ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr); 268f747e1aeSHong Zhang PetscFunctionReturn(0); 269f747e1aeSHong Zhang } 270f747e1aeSHong Zhang 271eb9c0419SKris Buschelman #undef __FUNCT__ 2729af31e4aSHong Zhang #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqAIJ" 273*55a3bba9SHong Zhang PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat P,PetscReal fill,Mat *C) 274*55a3bba9SHong Zhang { 275dfbe8321SBarry Smith PetscErrorCode ierr; 276d20bfe6fSHong Zhang FreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 277d20bfe6fSHong Zhang Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*p=(Mat_SeqAIJ*)P->data,*c; 278*55a3bba9SHong Zhang PetscInt *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj; 279*55a3bba9SHong Zhang PetscInt *ci,*cj,*ptadenserow,*ptasparserow,*ptaj,nlnk,*lnk; 280*55a3bba9SHong Zhang PetscInt an=A->N,am=A->M,pn=P->N,pm=P->M; 281*55a3bba9SHong Zhang PetscInt i,j,k,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi; 282d20bfe6fSHong Zhang MatScalar *ca; 283*55a3bba9SHong Zhang PetscBT lnkbt; 284eb9c0419SKris Buschelman 285eb9c0419SKris Buschelman PetscFunctionBegin; 286d20bfe6fSHong Zhang /* Get ij structure of P^T */ 287eb9c0419SKris Buschelman ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr); 288d20bfe6fSHong Zhang ptJ=ptj; 289eb9c0419SKris Buschelman 290d20bfe6fSHong Zhang /* Allocate ci array, arrays for fill computation and */ 291d20bfe6fSHong Zhang /* free space for accumulating nonzero column info */ 292*55a3bba9SHong Zhang ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 293d20bfe6fSHong Zhang ci[0] = 0; 294eb9c0419SKris Buschelman 295*55a3bba9SHong Zhang ierr = PetscMalloc((2*an+1)*sizeof(PetscInt),&ptadenserow);CHKERRQ(ierr); 296*55a3bba9SHong Zhang ierr = PetscMemzero(ptadenserow,(2*an+1)*sizeof(PetscInt));CHKERRQ(ierr); 297d20bfe6fSHong Zhang ptasparserow = ptadenserow + an; 298*55a3bba9SHong Zhang 299*55a3bba9SHong Zhang /* create and initialize a linked list */ 300*55a3bba9SHong Zhang nlnk = pn+1; 301*55a3bba9SHong Zhang ierr = PetscLLCreate(pn,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr); 302eb9c0419SKris Buschelman 303d20bfe6fSHong Zhang /* Set initial free space to be nnz(A) scaled by aspect ratio of P. */ 304d20bfe6fSHong Zhang /* This should be reasonable if sparsity of PtAP is similar to that of A. */ 305d20bfe6fSHong Zhang ierr = GetMoreSpace((ai[am]/pm)*pn,&free_space); 306d20bfe6fSHong Zhang current_space = free_space; 307d20bfe6fSHong Zhang 308d20bfe6fSHong Zhang /* Determine symbolic info for each row of C: */ 309d20bfe6fSHong Zhang for (i=0;i<pn;i++) { 310d20bfe6fSHong Zhang ptnzi = pti[i+1] - pti[i]; 311d20bfe6fSHong Zhang ptanzi = 0; 312d20bfe6fSHong Zhang /* Determine symbolic row of PtA: */ 313d20bfe6fSHong Zhang for (j=0;j<ptnzi;j++) { 314d20bfe6fSHong Zhang arow = *ptJ++; 315d20bfe6fSHong Zhang anzj = ai[arow+1] - ai[arow]; 316d20bfe6fSHong Zhang ajj = aj + ai[arow]; 317d20bfe6fSHong Zhang for (k=0;k<anzj;k++) { 318d20bfe6fSHong Zhang if (!ptadenserow[ajj[k]]) { 319d20bfe6fSHong Zhang ptadenserow[ajj[k]] = -1; 320d20bfe6fSHong Zhang ptasparserow[ptanzi++] = ajj[k]; 321d20bfe6fSHong Zhang } 322d20bfe6fSHong Zhang } 323d20bfe6fSHong Zhang } 324d20bfe6fSHong Zhang /* Using symbolic info for row of PtA, determine symbolic info for row of C: */ 325d20bfe6fSHong Zhang ptaj = ptasparserow; 326d20bfe6fSHong Zhang cnzi = 0; 327d20bfe6fSHong Zhang for (j=0;j<ptanzi;j++) { 328d20bfe6fSHong Zhang prow = *ptaj++; 329d20bfe6fSHong Zhang pnzj = pi[prow+1] - pi[prow]; 330d20bfe6fSHong Zhang pjj = pj + pi[prow]; 331*55a3bba9SHong Zhang /* add non-zero cols of P into the sorted linked list lnk */ 332*55a3bba9SHong Zhang ierr = PetscLLAdd(pnzj,pjj,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr); 333*55a3bba9SHong Zhang cnzi += nlnk; 334d20bfe6fSHong Zhang } 335d20bfe6fSHong Zhang 336d20bfe6fSHong Zhang /* If free space is not available, make more free space */ 337d20bfe6fSHong Zhang /* Double the amount of total space in the list */ 338d20bfe6fSHong Zhang if (current_space->local_remaining<cnzi) { 339d20bfe6fSHong Zhang ierr = GetMoreSpace(current_space->total_array_size,¤t_space);CHKERRQ(ierr); 340d20bfe6fSHong Zhang } 341d20bfe6fSHong Zhang 342d20bfe6fSHong Zhang /* Copy data into free space, and zero out denserows */ 343*55a3bba9SHong Zhang ierr = PetscLLClean(pn,pn,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 344d20bfe6fSHong Zhang current_space->array += cnzi; 345d20bfe6fSHong Zhang current_space->local_used += cnzi; 346d20bfe6fSHong Zhang current_space->local_remaining -= cnzi; 347d20bfe6fSHong Zhang 348d20bfe6fSHong Zhang for (j=0;j<ptanzi;j++) { 349d20bfe6fSHong Zhang ptadenserow[ptasparserow[j]] = 0; 350d20bfe6fSHong Zhang } 351d20bfe6fSHong Zhang /* Aside: Perhaps we should save the pta info for the numerical factorization. */ 352d20bfe6fSHong Zhang /* For now, we will recompute what is needed. */ 353d20bfe6fSHong Zhang ci[i+1] = ci[i] + cnzi; 354d20bfe6fSHong Zhang } 355d20bfe6fSHong Zhang /* nnz is now stored in ci[ptm], column indices are in the list of free space */ 356d20bfe6fSHong Zhang /* Allocate space for cj, initialize cj, and */ 357d20bfe6fSHong Zhang /* destroy list of free space and other temporary array(s) */ 358*55a3bba9SHong Zhang ierr = PetscMalloc((ci[pn]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr); 359d20bfe6fSHong Zhang ierr = MakeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 360d20bfe6fSHong Zhang ierr = PetscFree(ptadenserow);CHKERRQ(ierr); 361*55a3bba9SHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 362d20bfe6fSHong Zhang 363d20bfe6fSHong Zhang /* Allocate space for ca */ 364d20bfe6fSHong Zhang ierr = PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 365d20bfe6fSHong Zhang ierr = PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));CHKERRQ(ierr); 366d20bfe6fSHong Zhang 367d20bfe6fSHong Zhang /* put together the new matrix */ 368d20bfe6fSHong Zhang ierr = MatCreateSeqAIJWithArrays(A->comm,pn,pn,ci,cj,ca,C);CHKERRQ(ierr); 369d20bfe6fSHong Zhang 370d20bfe6fSHong Zhang /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 371d20bfe6fSHong Zhang /* Since these are PETSc arrays, change flags to free them as necessary. */ 372d20bfe6fSHong Zhang c = (Mat_SeqAIJ *)((*C)->data); 373d20bfe6fSHong Zhang c->freedata = PETSC_TRUE; 374d20bfe6fSHong Zhang c->nonew = 0; 375d20bfe6fSHong Zhang 376d20bfe6fSHong Zhang /* Clean up. */ 377d20bfe6fSHong Zhang ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr); 378eb9c0419SKris Buschelman 379eb9c0419SKris Buschelman PetscFunctionReturn(0); 380eb9c0419SKris Buschelman } 381eb9c0419SKris Buschelman 3823985e5eaSKris Buschelman #include "src/mat/impls/maij/maij.h" 3833985e5eaSKris Buschelman EXTERN_C_BEGIN 3843985e5eaSKris Buschelman #undef __FUNCT__ 3859af31e4aSHong Zhang #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqMAIJ" 386*55a3bba9SHong Zhang PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqMAIJ(Mat A,Mat PP,Mat *C) 387*55a3bba9SHong Zhang { 3885c66b693SKris Buschelman /* This routine requires testing -- I don't think it works. */ 389dfbe8321SBarry Smith PetscErrorCode ierr; 3903985e5eaSKris Buschelman FreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 3913985e5eaSKris Buschelman Mat_SeqMAIJ *pp=(Mat_SeqMAIJ*)PP->data; 3923985e5eaSKris Buschelman Mat P=pp->AIJ; 3933985e5eaSKris Buschelman Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*p=(Mat_SeqAIJ*)P->data,*c; 3943985e5eaSKris Buschelman int *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj; 3953985e5eaSKris Buschelman int *ci,*cj,*denserow,*sparserow,*ptadenserow,*ptasparserow,*ptaj; 3963985e5eaSKris Buschelman int an=A->N,am=A->M,pn=P->N,pm=P->M,ppdof=pp->dof; 397fe05a634SKris Buschelman int i,j,k,dof,pdof,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi; 3983985e5eaSKris Buschelman MatScalar *ca; 3993985e5eaSKris Buschelman 4003985e5eaSKris Buschelman PetscFunctionBegin; 4013985e5eaSKris Buschelman /* Start timer */ 4029af31e4aSHong Zhang ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,PP,0,0);CHKERRQ(ierr); 4033985e5eaSKris Buschelman 4043985e5eaSKris Buschelman /* Get ij structure of P^T */ 4053985e5eaSKris Buschelman ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr); 4063985e5eaSKris Buschelman 4073985e5eaSKris Buschelman /* Allocate ci array, arrays for fill computation and */ 4083985e5eaSKris Buschelman /* free space for accumulating nonzero column info */ 4093985e5eaSKris Buschelman ierr = PetscMalloc((pn+1)*sizeof(int),&ci);CHKERRQ(ierr); 4103985e5eaSKris Buschelman ci[0] = 0; 4113985e5eaSKris Buschelman 4123985e5eaSKris Buschelman ierr = PetscMalloc((2*pn+2*an+1)*sizeof(int),&ptadenserow);CHKERRQ(ierr); 4133985e5eaSKris Buschelman ierr = PetscMemzero(ptadenserow,(2*pn+2*an+1)*sizeof(int));CHKERRQ(ierr); 4143985e5eaSKris Buschelman ptasparserow = ptadenserow + an; 4153985e5eaSKris Buschelman denserow = ptasparserow + an; 4163985e5eaSKris Buschelman sparserow = denserow + pn; 4173985e5eaSKris Buschelman 4183985e5eaSKris Buschelman /* Set initial free space to be nnz(A) scaled by aspect ratio of P. */ 4193985e5eaSKris Buschelman /* This should be reasonable if sparsity of PtAP is similar to that of A. */ 4203985e5eaSKris Buschelman ierr = GetMoreSpace((ai[am]/pm)*pn,&free_space); 4213985e5eaSKris Buschelman current_space = free_space; 4223985e5eaSKris Buschelman 4233985e5eaSKris Buschelman /* Determine symbolic info for each row of C: */ 4243985e5eaSKris Buschelman for (i=0;i<pn/ppdof;i++) { 4253985e5eaSKris Buschelman ptnzi = pti[i+1] - pti[i]; 4263985e5eaSKris Buschelman ptanzi = 0; 4273985e5eaSKris Buschelman ptJ = ptj + pti[i]; 4283985e5eaSKris Buschelman for (dof=0;dof<ppdof;dof++) { 4293985e5eaSKris Buschelman /* Determine symbolic row of PtA: */ 4303985e5eaSKris Buschelman for (j=0;j<ptnzi;j++) { 4313985e5eaSKris Buschelman arow = ptJ[j] + dof; 4323985e5eaSKris Buschelman anzj = ai[arow+1] - ai[arow]; 4333985e5eaSKris Buschelman ajj = aj + ai[arow]; 4343985e5eaSKris Buschelman for (k=0;k<anzj;k++) { 4353985e5eaSKris Buschelman if (!ptadenserow[ajj[k]]) { 4363985e5eaSKris Buschelman ptadenserow[ajj[k]] = -1; 4373985e5eaSKris Buschelman ptasparserow[ptanzi++] = ajj[k]; 4383985e5eaSKris Buschelman } 4393985e5eaSKris Buschelman } 4403985e5eaSKris Buschelman } 4413985e5eaSKris Buschelman /* Using symbolic info for row of PtA, determine symbolic info for row of C: */ 4423985e5eaSKris Buschelman ptaj = ptasparserow; 4433985e5eaSKris Buschelman cnzi = 0; 4443985e5eaSKris Buschelman for (j=0;j<ptanzi;j++) { 445fe05a634SKris Buschelman pdof = *ptaj%dof; 4463985e5eaSKris Buschelman prow = (*ptaj++)/dof; 4473985e5eaSKris Buschelman pnzj = pi[prow+1] - pi[prow]; 4483985e5eaSKris Buschelman pjj = pj + pi[prow]; 4493985e5eaSKris Buschelman for (k=0;k<pnzj;k++) { 450fe05a634SKris Buschelman if (!denserow[pjj[k]+pdof]) { 451fe05a634SKris Buschelman denserow[pjj[k]+pdof] = -1; 452fe05a634SKris Buschelman sparserow[cnzi++] = pjj[k]+pdof; 4533985e5eaSKris Buschelman } 4543985e5eaSKris Buschelman } 4553985e5eaSKris Buschelman } 4563985e5eaSKris Buschelman 4573985e5eaSKris Buschelman /* sort sparserow */ 4583985e5eaSKris Buschelman ierr = PetscSortInt(cnzi,sparserow);CHKERRQ(ierr); 4593985e5eaSKris Buschelman 4603985e5eaSKris Buschelman /* If free space is not available, make more free space */ 4613985e5eaSKris Buschelman /* Double the amount of total space in the list */ 4623985e5eaSKris Buschelman if (current_space->local_remaining<cnzi) { 4633985e5eaSKris Buschelman ierr = GetMoreSpace(current_space->total_array_size,¤t_space);CHKERRQ(ierr); 4643985e5eaSKris Buschelman } 4653985e5eaSKris Buschelman 4663985e5eaSKris Buschelman /* Copy data into free space, and zero out denserows */ 4673985e5eaSKris Buschelman ierr = PetscMemcpy(current_space->array,sparserow,cnzi*sizeof(int));CHKERRQ(ierr); 4683985e5eaSKris Buschelman current_space->array += cnzi; 4693985e5eaSKris Buschelman current_space->local_used += cnzi; 4703985e5eaSKris Buschelman current_space->local_remaining -= cnzi; 4713985e5eaSKris Buschelman 4723985e5eaSKris Buschelman for (j=0;j<ptanzi;j++) { 4733985e5eaSKris Buschelman ptadenserow[ptasparserow[j]] = 0; 4743985e5eaSKris Buschelman } 4753985e5eaSKris Buschelman for (j=0;j<cnzi;j++) { 4763985e5eaSKris Buschelman denserow[sparserow[j]] = 0; 4773985e5eaSKris Buschelman } 4783985e5eaSKris Buschelman /* Aside: Perhaps we should save the pta info for the numerical factorization. */ 4793985e5eaSKris Buschelman /* For now, we will recompute what is needed. */ 4803985e5eaSKris Buschelman ci[i+1+dof] = ci[i+dof] + cnzi; 4813985e5eaSKris Buschelman } 4823985e5eaSKris Buschelman } 4833985e5eaSKris Buschelman /* nnz is now stored in ci[ptm], column indices are in the list of free space */ 4843985e5eaSKris Buschelman /* Allocate space for cj, initialize cj, and */ 4853985e5eaSKris Buschelman /* destroy list of free space and other temporary array(s) */ 4863985e5eaSKris Buschelman ierr = PetscMalloc((ci[pn]+1)*sizeof(int),&cj);CHKERRQ(ierr); 4873985e5eaSKris Buschelman ierr = MakeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 4883985e5eaSKris Buschelman ierr = PetscFree(ptadenserow);CHKERRQ(ierr); 4893985e5eaSKris Buschelman 4903985e5eaSKris Buschelman /* Allocate space for ca */ 4913985e5eaSKris Buschelman ierr = PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 4923985e5eaSKris Buschelman ierr = PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));CHKERRQ(ierr); 4933985e5eaSKris Buschelman 4943985e5eaSKris Buschelman /* put together the new matrix */ 4953985e5eaSKris Buschelman ierr = MatCreateSeqAIJWithArrays(A->comm,pn,pn,ci,cj,ca,C);CHKERRQ(ierr); 4963985e5eaSKris Buschelman 4973985e5eaSKris Buschelman /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 4983985e5eaSKris Buschelman /* Since these are PETSc arrays, change flags to free them as necessary. */ 4993985e5eaSKris Buschelman c = (Mat_SeqAIJ *)((*C)->data); 5003985e5eaSKris Buschelman c->freedata = PETSC_TRUE; 5013985e5eaSKris Buschelman c->nonew = 0; 5023985e5eaSKris Buschelman 5033985e5eaSKris Buschelman /* Clean up. */ 5043985e5eaSKris Buschelman ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr); 5053985e5eaSKris Buschelman 5069af31e4aSHong Zhang ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,PP,0,0);CHKERRQ(ierr); 5073985e5eaSKris Buschelman PetscFunctionReturn(0); 5083985e5eaSKris Buschelman } 5093985e5eaSKris Buschelman EXTERN_C_END 5103985e5eaSKris Buschelman 511eb9c0419SKris Buschelman #undef __FUNCT__ 5129af31e4aSHong Zhang #define __FUNCT__ "MatPtAPNumeric" 5136849ba73SBarry Smith /* 5149af31e4aSHong Zhang MatPtAPNumeric - Computes the matrix projection C = P^T * A * P 5154d3841fdSKris Buschelman 5164d3841fdSKris Buschelman Collective on Mat 5174d3841fdSKris Buschelman 5184d3841fdSKris Buschelman Input Parameters: 5194d3841fdSKris Buschelman + A - the matrix 5204d3841fdSKris Buschelman - P - the projection matrix 5214d3841fdSKris Buschelman 5224d3841fdSKris Buschelman Output Parameters: 5234d3841fdSKris Buschelman . C - the product matrix 5244d3841fdSKris Buschelman 5254d3841fdSKris Buschelman Notes: 5269af31e4aSHong Zhang C must have been created by calling MatPtAPSymbolic and must be destroyed by 5274d3841fdSKris Buschelman the user using MatDeatroy(). 5284d3841fdSKris Buschelman 529170ef064SHong Zhang This routine is currently only implemented for pairs of AIJ matrices and classes 530170ef064SHong Zhang which inherit from AIJ. C will be of type MATAIJ. 5314d3841fdSKris Buschelman 5324d3841fdSKris Buschelman Level: intermediate 5334d3841fdSKris Buschelman 5349af31e4aSHong Zhang .seealso: MatPtAP(),MatPtAPSymbolic(),MatMatMultNumeric() 5356849ba73SBarry Smith */ 536*55a3bba9SHong Zhang PetscErrorCode MatPtAPNumeric(Mat A,Mat P,Mat C) 537*55a3bba9SHong Zhang { 538dfbe8321SBarry Smith PetscErrorCode ierr; 539534c1384SKris Buschelman PetscErrorCode (*fA)(Mat,Mat,Mat); 540534c1384SKris Buschelman PetscErrorCode (*fP)(Mat,Mat,Mat); 541eb9c0419SKris Buschelman 542eb9c0419SKris Buschelman PetscFunctionBegin; 543eb9c0419SKris Buschelman 5444482741eSBarry Smith PetscValidHeaderSpecific(A,MAT_COOKIE,1); 545c9780b6fSBarry Smith PetscValidType(A,1); 546eb9c0419SKris Buschelman MatPreallocated(A); 547eb9c0419SKris Buschelman if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 548eb9c0419SKris Buschelman if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 549eb9c0419SKris Buschelman 5504482741eSBarry Smith PetscValidHeaderSpecific(P,MAT_COOKIE,2); 551c9780b6fSBarry Smith PetscValidType(P,2); 552eb9c0419SKris Buschelman MatPreallocated(P); 553eb9c0419SKris Buschelman if (!P->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 554eb9c0419SKris Buschelman if (P->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 555eb9c0419SKris Buschelman 5564482741eSBarry Smith PetscValidHeaderSpecific(C,MAT_COOKIE,3); 557c9780b6fSBarry Smith PetscValidType(C,3); 558eb9c0419SKris Buschelman MatPreallocated(C); 559eb9c0419SKris Buschelman if (!C->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 560eb9c0419SKris Buschelman if (C->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 561eb9c0419SKris Buschelman 562eb9c0419SKris Buschelman if (P->N!=C->M) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->N,C->M); 563eb9c0419SKris Buschelman if (P->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->M,A->N); 564eb9c0419SKris Buschelman if (A->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix 'A' must be square, %d != %d",A->M,A->N); 565eb9c0419SKris Buschelman if (P->N!=C->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->N,C->N); 566eb9c0419SKris Buschelman 567534c1384SKris Buschelman /* For now, we do not dispatch based on the type of A and P */ 568534c1384SKris Buschelman /* When implementations like _SeqAIJ_MAIJ exist, attack the multiple dispatch problem. */ 569534c1384SKris Buschelman fA = A->ops->ptapnumeric; 570534c1384SKris Buschelman if (!fA) SETERRQ1(PETSC_ERR_SUP,"MatPtAPNumeric not supported for A of type %s",A->type_name); 571534c1384SKris Buschelman fP = P->ops->ptapnumeric; 572534c1384SKris Buschelman if (!fP) SETERRQ1(PETSC_ERR_SUP,"MatPtAPNumeric not supported for P of type %s",P->type_name); 573534c1384SKris 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); 5744d3841fdSKris Buschelman 575534c1384SKris Buschelman ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 576534c1384SKris Buschelman ierr = (*fA)(A,P,C);CHKERRQ(ierr); 577534c1384SKris Buschelman ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 578eb9c0419SKris Buschelman 579eb9c0419SKris Buschelman PetscFunctionReturn(0); 580eb9c0419SKris Buschelman } 581eb9c0419SKris Buschelman 582eb9c0419SKris Buschelman #undef __FUNCT__ 5839af31e4aSHong Zhang #define __FUNCT__ "MatPtAPNumeric_SeqAIJ_SeqAIJ" 584dfbe8321SBarry Smith PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat A,Mat P,Mat C) 585dfbe8321SBarry Smith { 586dfbe8321SBarry Smith PetscErrorCode ierr; 587d20bfe6fSHong Zhang int flops=0; 588d20bfe6fSHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 589d20bfe6fSHong Zhang Mat_SeqAIJ *p = (Mat_SeqAIJ *) P->data; 590d20bfe6fSHong Zhang Mat_SeqAIJ *c = (Mat_SeqAIJ *) C->data; 591d20bfe6fSHong Zhang int *ai=a->i,*aj=a->j,*apj,*apjdense,*pi=p->i,*pj=p->j,*pJ=p->j,*pjj; 592d20bfe6fSHong Zhang int *ci=c->i,*cj=c->j,*cjj; 593d20bfe6fSHong Zhang int am=A->M,cn=C->N,cm=C->M; 594d20bfe6fSHong Zhang int i,j,k,anzi,pnzi,apnzj,nextap,pnzj,prow,crow; 595d20bfe6fSHong Zhang MatScalar *aa=a->a,*apa,*pa=p->a,*pA=p->a,*paj,*ca=c->a,*caj; 596eb9c0419SKris Buschelman 597eb9c0419SKris Buschelman PetscFunctionBegin; 598d20bfe6fSHong Zhang /* Allocate temporary array for storage of one row of A*P */ 599d20bfe6fSHong Zhang ierr = PetscMalloc(cn*(sizeof(MatScalar)+2*sizeof(int)),&apa);CHKERRQ(ierr); 600d20bfe6fSHong Zhang ierr = PetscMemzero(apa,cn*(sizeof(MatScalar)+2*sizeof(int)));CHKERRQ(ierr); 601eb9c0419SKris Buschelman 602d20bfe6fSHong Zhang apj = (int *)(apa + cn); 603d20bfe6fSHong Zhang apjdense = apj + cn; 604d20bfe6fSHong Zhang 605d20bfe6fSHong Zhang /* Clear old values in C */ 606d20bfe6fSHong Zhang ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 607d20bfe6fSHong Zhang 608d20bfe6fSHong Zhang for (i=0;i<am;i++) { 609d20bfe6fSHong Zhang /* Form sparse row of A*P */ 610d20bfe6fSHong Zhang anzi = ai[i+1] - ai[i]; 611d20bfe6fSHong Zhang apnzj = 0; 612d20bfe6fSHong Zhang for (j=0;j<anzi;j++) { 613d20bfe6fSHong Zhang prow = *aj++; 614d20bfe6fSHong Zhang pnzj = pi[prow+1] - pi[prow]; 615d20bfe6fSHong Zhang pjj = pj + pi[prow]; 616d20bfe6fSHong Zhang paj = pa + pi[prow]; 617d20bfe6fSHong Zhang for (k=0;k<pnzj;k++) { 618d20bfe6fSHong Zhang if (!apjdense[pjj[k]]) { 619d20bfe6fSHong Zhang apjdense[pjj[k]] = -1; 620d20bfe6fSHong Zhang apj[apnzj++] = pjj[k]; 621d20bfe6fSHong Zhang } 622d20bfe6fSHong Zhang apa[pjj[k]] += (*aa)*paj[k]; 623d20bfe6fSHong Zhang } 624d20bfe6fSHong Zhang flops += 2*pnzj; 625d20bfe6fSHong Zhang aa++; 626d20bfe6fSHong Zhang } 627d20bfe6fSHong Zhang 628d20bfe6fSHong Zhang /* Sort the j index array for quick sparse axpy. */ 629d20bfe6fSHong Zhang ierr = PetscSortInt(apnzj,apj);CHKERRQ(ierr); 630d20bfe6fSHong Zhang 631d20bfe6fSHong Zhang /* Compute P^T*A*P using outer product (P^T)[:,j]*(A*P)[j,:]. */ 632d20bfe6fSHong Zhang pnzi = pi[i+1] - pi[i]; 633d20bfe6fSHong Zhang for (j=0;j<pnzi;j++) { 634d20bfe6fSHong Zhang nextap = 0; 635d20bfe6fSHong Zhang crow = *pJ++; 636d20bfe6fSHong Zhang cjj = cj + ci[crow]; 637d20bfe6fSHong Zhang caj = ca + ci[crow]; 638d20bfe6fSHong Zhang /* Perform sparse axpy operation. Note cjj includes apj. */ 639d20bfe6fSHong Zhang for (k=0;nextap<apnzj;k++) { 640d20bfe6fSHong Zhang if (cjj[k]==apj[nextap]) { 641d20bfe6fSHong Zhang caj[k] += (*pA)*apa[apj[nextap++]]; 642d20bfe6fSHong Zhang } 643d20bfe6fSHong Zhang } 644d20bfe6fSHong Zhang flops += 2*apnzj; 645d20bfe6fSHong Zhang pA++; 646d20bfe6fSHong Zhang } 647d20bfe6fSHong Zhang 648d20bfe6fSHong Zhang /* Zero the current row info for A*P */ 649d20bfe6fSHong Zhang for (j=0;j<apnzj;j++) { 650d20bfe6fSHong Zhang apa[apj[j]] = 0.; 651d20bfe6fSHong Zhang apjdense[apj[j]] = 0; 652d20bfe6fSHong Zhang } 653d20bfe6fSHong Zhang } 654d20bfe6fSHong Zhang 655d20bfe6fSHong Zhang /* Assemble the final matrix and clean up */ 656d20bfe6fSHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 657d20bfe6fSHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 658d20bfe6fSHong Zhang ierr = PetscFree(apa);CHKERRQ(ierr); 659d20bfe6fSHong Zhang ierr = PetscLogFlops(flops);CHKERRQ(ierr); 660d20bfe6fSHong Zhang 661eb9c0419SKris Buschelman PetscFunctionReturn(0); 662eb9c0419SKris Buschelman } 6630e36024fSHong Zhang 6640e36024fSHong Zhang /* Compute C = P[rstart:rend,:]^T * A * P of seqaij matrices - used by MatPtAP_MPIAIJ_MPIAIJ() */ 6650e36024fSHong Zhang 6660e36024fSHong Zhang #undef __FUNCT__ 6670e36024fSHong Zhang #define __FUNCT__ "MatPtAPNumeric_SeqAIJ_SeqAIJ_ReducedPt" 6687f79fd58SMatthew Knepley /*@C 669e9b43d0fSSatish Balay MatPtAPNumeric_SeqAIJ_SeqAIJ_ReducedPt - Creates C = P[rstart:rend,:]^T * A * P of seqaij matrices, 670b90dcfe3SHong Zhang used by MatPtAP_MPIAIJ_MPIAIJ() 671b90dcfe3SHong Zhang 672b90dcfe3SHong Zhang Collective on Mat 673b90dcfe3SHong Zhang 674b90dcfe3SHong Zhang Input Parameters: 675b90dcfe3SHong Zhang + A - the matrix in seqaij format 676b90dcfe3SHong Zhang . P - the projection matrix in seqaij format 677b90dcfe3SHong Zhang . scall - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX 678b90dcfe3SHong Zhang . fill - expected fill (not being used when scall==MAT_REUSE_MATRIX) 679b90dcfe3SHong Zhang . prstart, prend - the starting and ending-1-th row of the matrix P to be used for transpose 680b90dcfe3SHong Zhang 681b90dcfe3SHong Zhang Output Parameters: 682b90dcfe3SHong Zhang . C - the product matrix in seqaij format 683b90dcfe3SHong Zhang 684b90dcfe3SHong Zhang Level: developer 685b90dcfe3SHong Zhang 686b90dcfe3SHong Zhang .seealso: MatPtAPSymbolic(),MatPtAPNumeric(),MatMatMult() 687b90dcfe3SHong Zhang @*/ 6880e36024fSHong Zhang PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_ReducedPt(Mat A,Mat P,int prstart,int prend,Mat C) 6890e36024fSHong Zhang { 6900e36024fSHong Zhang PetscErrorCode ierr; 6910e36024fSHong Zhang int flops=0; 6920e36024fSHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 6930e36024fSHong Zhang Mat_SeqAIJ *p = (Mat_SeqAIJ *) P->data; 6940e36024fSHong Zhang Mat_SeqAIJ *c = (Mat_SeqAIJ *) C->data; 69520d4747cSHong Zhang int *ai=a->i,*aj=a->j,*apj,*apjdense; 69620d4747cSHong Zhang int *pi=p->i,*pj=p->j,*pJ=p->j+pi[prstart],*pjj; 6970e36024fSHong Zhang int *ci=c->i,*cj=c->j,*cjj; 6980e36024fSHong Zhang int am=A->M,cn=C->N,cm=C->M; 6990e36024fSHong Zhang int i,j,k,anzi,pnzi,apnzj,nextap,pnzj,prow,crow; 7000e36024fSHong Zhang MatScalar *aa=a->a,*apa,*pa=p->a,*pA=p->a,*paj,*ca=c->a,*caj; 7010e36024fSHong Zhang 7020e36024fSHong Zhang PetscFunctionBegin; 7030e36024fSHong Zhang pA=p->a+pi[prstart]; 7040e36024fSHong Zhang /* Allocate temporary array for storage of one row of A*P */ 7050e36024fSHong Zhang ierr = PetscMalloc(cn*(sizeof(MatScalar)+2*sizeof(int)),&apa);CHKERRQ(ierr); 7060e36024fSHong Zhang ierr = PetscMemzero(apa,cn*(sizeof(MatScalar)+2*sizeof(int)));CHKERRQ(ierr); 7070e36024fSHong Zhang 7080e36024fSHong Zhang apj = (int *)(apa + cn); 7090e36024fSHong Zhang apjdense = apj + cn; 7100e36024fSHong Zhang 7110e36024fSHong Zhang /* Clear old values in C */ 7120e36024fSHong Zhang ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 7130e36024fSHong Zhang 7140e36024fSHong Zhang for (i=0;i<am;i++) { 7150e36024fSHong Zhang /* Form sparse row of A*P */ 7160e36024fSHong Zhang anzi = ai[i+1] - ai[i]; 7170e36024fSHong Zhang apnzj = 0; 7180e36024fSHong Zhang for (j=0;j<anzi;j++) { 7190e36024fSHong Zhang prow = *aj++; 7200e36024fSHong Zhang pnzj = pi[prow+1] - pi[prow]; 7210e36024fSHong Zhang pjj = pj + pi[prow]; 7220e36024fSHong Zhang paj = pa + pi[prow]; 7230e36024fSHong Zhang for (k=0;k<pnzj;k++) { 7240e36024fSHong Zhang if (!apjdense[pjj[k]]) { 7250e36024fSHong Zhang apjdense[pjj[k]] = -1; 7260e36024fSHong Zhang apj[apnzj++] = pjj[k]; 7270e36024fSHong Zhang } 7280e36024fSHong Zhang apa[pjj[k]] += (*aa)*paj[k]; 7290e36024fSHong Zhang } 7300e36024fSHong Zhang flops += 2*pnzj; 7310e36024fSHong Zhang aa++; 7320e36024fSHong Zhang } 7330e36024fSHong Zhang 7340e36024fSHong Zhang /* Sort the j index array for quick sparse axpy. */ 7350e36024fSHong Zhang ierr = PetscSortInt(apnzj,apj);CHKERRQ(ierr); 7360e36024fSHong Zhang 7370e36024fSHong Zhang /* Compute P[[prstart:prend,:]^T*A*P using outer product (P^T)[:,j+prstart]*(A*P)[j,:]. */ 7380e36024fSHong Zhang pnzi = pi[i+1+prstart] - pi[i+prstart]; 7390e36024fSHong Zhang for (j=0;j<pnzi;j++) { 7400e36024fSHong Zhang nextap = 0; 7410e36024fSHong Zhang crow = *pJ++; 7420e36024fSHong Zhang cjj = cj + ci[crow]; 7430e36024fSHong Zhang caj = ca + ci[crow]; 7440e36024fSHong Zhang /* Perform sparse axpy operation. Note cjj includes apj. */ 7450e36024fSHong Zhang for (k=0;nextap<apnzj;k++) { 7460e36024fSHong Zhang if (cjj[k]==apj[nextap]) { 7470e36024fSHong Zhang caj[k] += (*pA)*apa[apj[nextap++]]; 7480e36024fSHong Zhang } 7490e36024fSHong Zhang } 7500e36024fSHong Zhang flops += 2*apnzj; 7510e36024fSHong Zhang pA++; 7520e36024fSHong Zhang } 7530e36024fSHong Zhang 7540e36024fSHong Zhang /* Zero the current row info for A*P */ 7550e36024fSHong Zhang for (j=0;j<apnzj;j++) { 7560e36024fSHong Zhang apa[apj[j]] = 0.; 7570e36024fSHong Zhang apjdense[apj[j]] = 0; 7580e36024fSHong Zhang } 7590e36024fSHong Zhang } 7600e36024fSHong Zhang 7610e36024fSHong Zhang /* Assemble the final matrix and clean up */ 7620e36024fSHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 7630e36024fSHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 7640e36024fSHong Zhang ierr = PetscFree(apa);CHKERRQ(ierr); 7650e36024fSHong Zhang ierr = PetscLogFlops(flops);CHKERRQ(ierr); 7660e36024fSHong Zhang 7670e36024fSHong Zhang PetscFunctionReturn(0); 7680e36024fSHong Zhang } 7690e36024fSHong Zhang 7700e36024fSHong Zhang #undef __FUNCT__ 7710e36024fSHong Zhang #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqAIJ_ReducedPt" 772*55a3bba9SHong Zhang PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_ReducedPt(Mat A,Mat P,PetscReal fill,int prstart,int prend,Mat *C) 773*55a3bba9SHong Zhang { 7740e36024fSHong Zhang PetscErrorCode ierr; 7750e36024fSHong Zhang FreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 7760e36024fSHong Zhang Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*p=(Mat_SeqAIJ*)P->data,*c; 777*55a3bba9SHong Zhang PetscInt *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj; 778*55a3bba9SHong Zhang PetscInt *ci,*cj,*ptadenserow,*ptasparserow,*ptaj,nlnk,*lnk; 779*55a3bba9SHong Zhang PetscInt an=A->N,am=A->M,pn=P->N,pm=P->M; 780*55a3bba9SHong Zhang PetscInt i,j,k,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi,m=prend - prstart; 7810e36024fSHong Zhang MatScalar *ca; 7820e36024fSHong Zhang Mat *psub,P_sub; 7830e36024fSHong Zhang IS isrow,iscol; 784*55a3bba9SHong Zhang PetscBT lnkbt; 7850b89d903Svictorle 7860b89d903Svictorle PetscFunctionBegin; 7870b89d903Svictorle /* Get ij structure of P[rstart:rend,:]^T */ 7880e36024fSHong Zhang ierr = ISCreateStride(PETSC_COMM_SELF,m,prstart,1,&isrow);CHKERRQ(ierr); 7890e36024fSHong Zhang ierr = ISCreateStride(PETSC_COMM_SELF,P->n,0,1,&iscol);CHKERRQ(ierr); 7900e36024fSHong Zhang ierr = MatGetSubMatrices(P,1,&isrow,&iscol,MAT_INITIAL_MATRIX,&psub);CHKERRQ(ierr); 7910e36024fSHong Zhang ierr = ISDestroy(isrow);CHKERRQ(ierr); 7920e36024fSHong Zhang ierr = ISDestroy(iscol);CHKERRQ(ierr); 7930e36024fSHong Zhang P_sub = psub[0]; 7940e36024fSHong Zhang ierr = MatGetSymbolicTranspose_SeqAIJ(P_sub,&pti,&ptj);CHKERRQ(ierr); 7950e36024fSHong Zhang ierr = MatDestroyMatrices(1,&psub);CHKERRQ(ierr); 7960e36024fSHong Zhang ptJ=ptj; 7970e36024fSHong Zhang 7980e36024fSHong Zhang /* Allocate ci array, arrays for fill computation and */ 7990e36024fSHong Zhang /* free space for accumulating nonzero column info */ 800*55a3bba9SHong Zhang ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 8010e36024fSHong Zhang ci[0] = 0; 8020e36024fSHong Zhang 803*55a3bba9SHong Zhang ierr = PetscMalloc((2*an+1)*sizeof(PetscInt),&ptadenserow);CHKERRQ(ierr); 804*55a3bba9SHong Zhang ierr = PetscMemzero(ptadenserow,(2*an+1)*sizeof(PetscInt));CHKERRQ(ierr); 8050e36024fSHong Zhang ptasparserow = ptadenserow + an; 806*55a3bba9SHong Zhang 807*55a3bba9SHong Zhang /* create and initialize a linked list */ 808*55a3bba9SHong Zhang nlnk = pn+1; 809*55a3bba9SHong Zhang ierr = PetscLLCreate(pn,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr); 8100e36024fSHong Zhang 8110e36024fSHong Zhang /* Set initial free space to be nnz(A) scaled by aspect ratio of P. */ 8120e36024fSHong Zhang /* This should be reasonable if sparsity of PtAP is similar to that of A. */ 813*55a3bba9SHong Zhang ierr = GetMoreSpace((PetscInt)(fill*ai[am]/pm)*pn,&free_space); 8140e36024fSHong Zhang current_space = free_space; 8150e36024fSHong Zhang 8160e36024fSHong Zhang /* Determine symbolic info for each row of C: */ 8170e36024fSHong Zhang for (i=0;i<pn;i++) { 8180e36024fSHong Zhang ptnzi = pti[i+1] - pti[i]; 8190e36024fSHong Zhang ptanzi = 0; 8200e36024fSHong Zhang /* Determine symbolic row of PtA_reduced: */ 8210e36024fSHong Zhang for (j=0;j<ptnzi;j++) { 8220e36024fSHong Zhang arow = *ptJ++; 8230e36024fSHong Zhang anzj = ai[arow+1] - ai[arow]; 8240e36024fSHong Zhang ajj = aj + ai[arow]; 8250e36024fSHong Zhang for (k=0;k<anzj;k++) { 8260e36024fSHong Zhang if (!ptadenserow[ajj[k]]) { 8270e36024fSHong Zhang ptadenserow[ajj[k]] = -1; 8280e36024fSHong Zhang ptasparserow[ptanzi++] = ajj[k]; 8290e36024fSHong Zhang } 8300e36024fSHong Zhang } 8310e36024fSHong Zhang } 8320e36024fSHong Zhang /* Using symbolic info for row of PtA, determine symbolic info for row of C: */ 8330e36024fSHong Zhang ptaj = ptasparserow; 8340e36024fSHong Zhang cnzi = 0; 8350e36024fSHong Zhang for (j=0;j<ptanzi;j++) { 8360e36024fSHong Zhang prow = *ptaj++; 8370e36024fSHong Zhang pnzj = pi[prow+1] - pi[prow]; 8380e36024fSHong Zhang pjj = pj + pi[prow]; 839*55a3bba9SHong Zhang /* add non-zero cols of P into the sorted linked list lnk */ 840*55a3bba9SHong Zhang ierr = PetscLLAdd(pnzj,pjj,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr); 841*55a3bba9SHong Zhang cnzi += nlnk; 8420e36024fSHong Zhang } 8430e36024fSHong Zhang 8440e36024fSHong Zhang /* If free space is not available, make more free space */ 8450e36024fSHong Zhang /* Double the amount of total space in the list */ 8460e36024fSHong Zhang if (current_space->local_remaining<cnzi) { 8470e36024fSHong Zhang ierr = GetMoreSpace(current_space->total_array_size,¤t_space);CHKERRQ(ierr); 8480e36024fSHong Zhang } 8490e36024fSHong Zhang 8500e36024fSHong Zhang /* Copy data into free space, and zero out denserows */ 851*55a3bba9SHong Zhang ierr = PetscLLClean(pn,pn,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 8520e36024fSHong Zhang current_space->array += cnzi; 8530e36024fSHong Zhang current_space->local_used += cnzi; 8540e36024fSHong Zhang current_space->local_remaining -= cnzi; 8550e36024fSHong Zhang 8560e36024fSHong Zhang for (j=0;j<ptanzi;j++) { 8570e36024fSHong Zhang ptadenserow[ptasparserow[j]] = 0; 8580e36024fSHong Zhang } 859*55a3bba9SHong Zhang 8600e36024fSHong Zhang /* Aside: Perhaps we should save the pta info for the numerical factorization. */ 8610e36024fSHong Zhang /* For now, we will recompute what is needed. */ 8620e36024fSHong Zhang ci[i+1] = ci[i] + cnzi; 8630e36024fSHong Zhang } 8640e36024fSHong Zhang /* nnz is now stored in ci[ptm], column indices are in the list of free space */ 8650e36024fSHong Zhang /* Allocate space for cj, initialize cj, and */ 8660e36024fSHong Zhang /* destroy list of free space and other temporary array(s) */ 867*55a3bba9SHong Zhang ierr = PetscMalloc((ci[pn]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr); 8680e36024fSHong Zhang ierr = MakeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 8690e36024fSHong Zhang ierr = PetscFree(ptadenserow);CHKERRQ(ierr); 870*55a3bba9SHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 8710e36024fSHong Zhang 8720e36024fSHong Zhang /* Allocate space for ca */ 8730e36024fSHong Zhang ierr = PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 8740e36024fSHong Zhang ierr = PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));CHKERRQ(ierr); 8750e36024fSHong Zhang 8760e36024fSHong Zhang /* put together the new matrix */ 8770e36024fSHong Zhang ierr = MatCreateSeqAIJWithArrays(A->comm,pn,pn,ci,cj,ca,C);CHKERRQ(ierr); 8780e36024fSHong Zhang 8790e36024fSHong Zhang /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 8800e36024fSHong Zhang /* Since these are PETSc arrays, change flags to free them as necessary. */ 8810e36024fSHong Zhang c = (Mat_SeqAIJ *)((*C)->data); 8820e36024fSHong Zhang c->freedata = PETSC_TRUE; 8830e36024fSHong Zhang c->nonew = 0; 8840e36024fSHong Zhang 8850e36024fSHong Zhang /* Clean up. */ 8860e36024fSHong Zhang ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr); 8870e36024fSHong Zhang 8880e36024fSHong Zhang PetscFunctionReturn(0); 8890e36024fSHong Zhang } 8900e36024fSHong Zhang 8910e36024fSHong Zhang #undef __FUNCT__ 8920e36024fSHong Zhang #define __FUNCT__ "MatPtAP_SeqAIJ_SeqAIJ_ReducedPt" 893*55a3bba9SHong Zhang PetscErrorCode MatPtAP_SeqAIJ_SeqAIJ_ReducedPt(Mat A,Mat P,MatReuse scall,PetscReal fill,PetscInt prstart,PetscInt prend,Mat *C) 8940e36024fSHong Zhang { 8950e36024fSHong Zhang PetscErrorCode ierr; 8960e36024fSHong Zhang PetscFunctionBegin; 8970e36024fSHong Zhang if (A->m != prend-prstart) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",A->m,prend-prstart); 8980e36024fSHong Zhang if (prend-prstart > P->m) SETERRQ2(PETSC_ERR_ARG_SIZ," prend-prstart %d cannot be larger than P->m %d",prend-prstart,P->m); 8990e36024fSHong Zhang if (scall == MAT_INITIAL_MATRIX){ 9000e36024fSHong Zhang ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ_ReducedPt(A,P,fill,prstart,prend,C);CHKERRQ(ierr); 9010e36024fSHong Zhang } 9020e36024fSHong Zhang 9030e36024fSHong Zhang ierr = MatPtAPNumeric_SeqAIJ_SeqAIJ_ReducedPt(A,P,prstart,prend,*C);CHKERRQ(ierr); 9040e36024fSHong Zhang 9050e36024fSHong Zhang PetscFunctionReturn(0); 9060e36024fSHong Zhang } 9070e36024fSHong Zhang 908