1eb9c0419SKris Buschelman /* 2*42c57489SHong Zhang Defines projective product routines where A is a SeqAIJ 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" 855a3bba9SHong Zhang #include "petscbt.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 @*/ 3655a3bba9SHong Zhang PetscErrorCode MatPtAP(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C) 3755a3bba9SHong Zhang { 38dfbe8321SBarry Smith PetscErrorCode ierr; 39534c1384SKris Buschelman PetscErrorCode (*fA)(Mat,Mat,MatReuse,PetscReal,Mat *); 40534c1384SKris Buschelman PetscErrorCode (*fP)(Mat,Mat,MatReuse,PetscReal,Mat *); 41eb9c0419SKris Buschelman 42eb9c0419SKris Buschelman PetscFunctionBegin; 439af31e4aSHong Zhang PetscValidHeaderSpecific(A,MAT_COOKIE,1); 449af31e4aSHong Zhang PetscValidType(A,1); 459af31e4aSHong Zhang MatPreallocated(A); 469af31e4aSHong Zhang if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 479af31e4aSHong Zhang if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 489af31e4aSHong Zhang PetscValidHeaderSpecific(P,MAT_COOKIE,2); 499af31e4aSHong Zhang PetscValidType(P,2); 509af31e4aSHong Zhang MatPreallocated(P); 519af31e4aSHong Zhang if (!P->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 529af31e4aSHong Zhang if (P->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 539af31e4aSHong Zhang PetscValidPointer(C,3); 54534c1384SKris Buschelman 5577431f27SBarry Smith if (P->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %D != %D",P->M,A->N); 56eb9c0419SKris Buschelman 579af31e4aSHong Zhang if (fill <=0.0) SETERRQ1(PETSC_ERR_ARG_SIZ,"fill=%g must be > 0.0",fill); 58eb9c0419SKris Buschelman 59534c1384SKris Buschelman /* For now, we do not dispatch based on the type of A and P */ 60534c1384SKris Buschelman /* When implementations like _SeqAIJ_MAIJ exist, attack the multiple dispatch problem. */ 61534c1384SKris Buschelman fA = A->ops->ptap; 62534c1384SKris Buschelman if (!fA) SETERRQ1(PETSC_ERR_SUP,"MatPtAP not supported for A of type %s",A->type_name); 63534c1384SKris Buschelman fP = P->ops->ptap; 64534c1384SKris Buschelman if (!fP) SETERRQ1(PETSC_ERR_SUP,"MatPtAP not supported for P of type %s",P->type_name); 65534c1384SKris 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); 66534c1384SKris Buschelman 679af31e4aSHong Zhang ierr = PetscLogEventBegin(MAT_PtAP,A,P,0,0);CHKERRQ(ierr); 68534c1384SKris Buschelman ierr = (*fA)(A,P,scall,fill,C);CHKERRQ(ierr); 699af31e4aSHong Zhang ierr = PetscLogEventEnd(MAT_PtAP,A,P,0,0);CHKERRQ(ierr); 70eb9c0419SKris Buschelman PetscFunctionReturn(0); 71eb9c0419SKris Buschelman } 72eb9c0419SKris Buschelman 73ff134f7aSHong Zhang #undef __FUNCT__ 749af31e4aSHong Zhang #define __FUNCT__ "MatPtAP_SeqAIJ_SeqAIJ" 75dfbe8321SBarry Smith PetscErrorCode MatPtAP_SeqAIJ_SeqAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C) 769af31e4aSHong Zhang { 77dfbe8321SBarry Smith PetscErrorCode ierr; 78b1d57f15SBarry Smith 799af31e4aSHong Zhang PetscFunctionBegin; 809af31e4aSHong Zhang if (scall == MAT_INITIAL_MATRIX){ 81d20bfe6fSHong Zhang ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 829af31e4aSHong Zhang ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ(A,P,fill,C);CHKERRQ(ierr); 83d20bfe6fSHong Zhang ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 849af31e4aSHong Zhang } 85d20bfe6fSHong Zhang ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 869af31e4aSHong Zhang ierr = MatPtAPNumeric_SeqAIJ_SeqAIJ(A,P,*C);CHKERRQ(ierr); 87d20bfe6fSHong Zhang ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 889af31e4aSHong Zhang PetscFunctionReturn(0); 899af31e4aSHong Zhang } 909af31e4aSHong Zhang 916849ba73SBarry Smith /* 929af31e4aSHong Zhang MatPtAPSymbolic - Creates the (i,j) structure of the matrix projection C = P^T * A * P 934d3841fdSKris Buschelman 944d3841fdSKris Buschelman Collective on Mat 954d3841fdSKris Buschelman 964d3841fdSKris Buschelman Input Parameters: 974d3841fdSKris Buschelman + A - the matrix 984d3841fdSKris Buschelman - P - the projection matrix 994d3841fdSKris Buschelman 1004d3841fdSKris Buschelman Output Parameters: 1014d3841fdSKris Buschelman . C - the (i,j) structure of the product matrix 1024d3841fdSKris Buschelman 1034d3841fdSKris Buschelman Notes: 1044d3841fdSKris Buschelman C will be created and must be destroyed by the user with MatDestroy(). 1054d3841fdSKris Buschelman 1064d3841fdSKris Buschelman This routine is currently only implemented for pairs of SeqAIJ matrices and classes 1074d3841fdSKris Buschelman which inherit from SeqAIJ. C will be of type MATSEQAIJ. The product is computed using 1089af31e4aSHong Zhang this (i,j) structure by calling MatPtAPNumeric(). 1094d3841fdSKris Buschelman 1104d3841fdSKris Buschelman Level: intermediate 1114d3841fdSKris Buschelman 1129af31e4aSHong Zhang .seealso: MatPtAP(),MatPtAPNumeric(),MatMatMultSymbolic() 1136849ba73SBarry Smith */ 114e240928fSHong Zhang #undef __FUNCT__ 115e240928fSHong Zhang #define __FUNCT__ "MatPtAPSymbolic" 11655a3bba9SHong Zhang PetscErrorCode MatPtAPSymbolic(Mat A,Mat P,PetscReal fill,Mat *C) 11755a3bba9SHong Zhang { 118dfbe8321SBarry Smith PetscErrorCode ierr; 119534c1384SKris Buschelman PetscErrorCode (*fA)(Mat,Mat,PetscReal,Mat*); 120534c1384SKris Buschelman PetscErrorCode (*fP)(Mat,Mat,PetscReal,Mat*); 121eb9c0419SKris Buschelman 122eb9c0419SKris Buschelman PetscFunctionBegin; 123eb9c0419SKris Buschelman 1244482741eSBarry Smith PetscValidHeaderSpecific(A,MAT_COOKIE,1); 125c9780b6fSBarry Smith PetscValidType(A,1); 126eb9c0419SKris Buschelman MatPreallocated(A); 127eb9c0419SKris Buschelman if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 128eb9c0419SKris Buschelman if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 129eb9c0419SKris Buschelman 1304482741eSBarry Smith PetscValidHeaderSpecific(P,MAT_COOKIE,2); 131c9780b6fSBarry Smith PetscValidType(P,2); 132eb9c0419SKris Buschelman MatPreallocated(P); 133eb9c0419SKris Buschelman if (!P->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 134eb9c0419SKris Buschelman if (P->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 135eb9c0419SKris Buschelman 1364482741eSBarry Smith PetscValidPointer(C,3); 1374482741eSBarry Smith 13877431f27SBarry Smith if (P->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %D != %D",P->M,A->N); 13977431f27SBarry Smith if (A->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix 'A' must be square, %D != %D",A->M,A->N); 140eb9c0419SKris Buschelman 141534c1384SKris Buschelman /* For now, we do not dispatch based on the type of A and P */ 142534c1384SKris Buschelman /* When implementations like _SeqAIJ_MAIJ exist, attack the multiple dispatch problem. */ 143534c1384SKris Buschelman fA = A->ops->ptapsymbolic; 144534c1384SKris Buschelman if (!fA) SETERRQ1(PETSC_ERR_SUP,"MatPtAPSymbolic not supported for A of type %s",A->type_name); 145534c1384SKris Buschelman fP = P->ops->ptapsymbolic; 146534c1384SKris Buschelman if (!fP) SETERRQ1(PETSC_ERR_SUP,"MatPtAPSymbolic not supported for P of type %s",P->type_name); 147534c1384SKris 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); 1484d3841fdSKris Buschelman 149534c1384SKris Buschelman ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 150534c1384SKris Buschelman ierr = (*fA)(A,P,fill,C);CHKERRQ(ierr); 151534c1384SKris Buschelman ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 152eb9c0419SKris Buschelman 153eb9c0419SKris Buschelman PetscFunctionReturn(0); 154eb9c0419SKris Buschelman } 155f747e1aeSHong Zhang 156eb9c0419SKris Buschelman #undef __FUNCT__ 1579af31e4aSHong Zhang #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqAIJ" 15855a3bba9SHong Zhang PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat P,PetscReal fill,Mat *C) 15955a3bba9SHong Zhang { 160dfbe8321SBarry Smith PetscErrorCode ierr; 161d20bfe6fSHong Zhang FreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 162d20bfe6fSHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*p = (Mat_SeqAIJ*)P->data,*c; 16355a3bba9SHong Zhang PetscInt *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj; 164b8374ebeSBarry Smith PetscInt *ci,*cj,*ptadenserow,*ptasparserow,*ptaj; 16555a3bba9SHong Zhang PetscInt an=A->N,am=A->M,pn=P->N,pm=P->M; 166b8374ebeSBarry Smith PetscInt i,j,k,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi,nlnk,*lnk; 167d20bfe6fSHong Zhang MatScalar *ca; 16855a3bba9SHong Zhang PetscBT lnkbt; 169eb9c0419SKris Buschelman 170eb9c0419SKris Buschelman PetscFunctionBegin; 171d20bfe6fSHong Zhang /* Get ij structure of P^T */ 172eb9c0419SKris Buschelman ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr); 173d20bfe6fSHong Zhang ptJ=ptj; 174eb9c0419SKris Buschelman 175d20bfe6fSHong Zhang /* Allocate ci array, arrays for fill computation and */ 176d20bfe6fSHong Zhang /* free space for accumulating nonzero column info */ 17755a3bba9SHong Zhang ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 178d20bfe6fSHong Zhang ci[0] = 0; 179eb9c0419SKris Buschelman 18055a3bba9SHong Zhang ierr = PetscMalloc((2*an+1)*sizeof(PetscInt),&ptadenserow);CHKERRQ(ierr); 18155a3bba9SHong Zhang ierr = PetscMemzero(ptadenserow,(2*an+1)*sizeof(PetscInt));CHKERRQ(ierr); 182d20bfe6fSHong Zhang ptasparserow = ptadenserow + an; 18355a3bba9SHong Zhang 18455a3bba9SHong Zhang /* create and initialize a linked list */ 18555a3bba9SHong Zhang nlnk = pn+1; 18655a3bba9SHong Zhang ierr = PetscLLCreate(pn,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr); 187eb9c0419SKris Buschelman 188d20bfe6fSHong Zhang /* Set initial free space to be nnz(A) scaled by aspect ratio of P. */ 189d20bfe6fSHong Zhang /* This should be reasonable if sparsity of PtAP is similar to that of A. */ 190d20bfe6fSHong Zhang ierr = GetMoreSpace((ai[am]/pm)*pn,&free_space); 191d20bfe6fSHong Zhang current_space = free_space; 192d20bfe6fSHong Zhang 193d20bfe6fSHong Zhang /* Determine symbolic info for each row of C: */ 194d20bfe6fSHong Zhang for (i=0;i<pn;i++) { 195d20bfe6fSHong Zhang ptnzi = pti[i+1] - pti[i]; 196d20bfe6fSHong Zhang ptanzi = 0; 197d20bfe6fSHong Zhang /* Determine symbolic row of PtA: */ 198d20bfe6fSHong Zhang for (j=0;j<ptnzi;j++) { 199d20bfe6fSHong Zhang arow = *ptJ++; 200d20bfe6fSHong Zhang anzj = ai[arow+1] - ai[arow]; 201d20bfe6fSHong Zhang ajj = aj + ai[arow]; 202d20bfe6fSHong Zhang for (k=0;k<anzj;k++) { 203d20bfe6fSHong Zhang if (!ptadenserow[ajj[k]]) { 204d20bfe6fSHong Zhang ptadenserow[ajj[k]] = -1; 205d20bfe6fSHong Zhang ptasparserow[ptanzi++] = ajj[k]; 206d20bfe6fSHong Zhang } 207d20bfe6fSHong Zhang } 208d20bfe6fSHong Zhang } 209d20bfe6fSHong Zhang /* Using symbolic info for row of PtA, determine symbolic info for row of C: */ 210d20bfe6fSHong Zhang ptaj = ptasparserow; 211d20bfe6fSHong Zhang cnzi = 0; 212d20bfe6fSHong Zhang for (j=0;j<ptanzi;j++) { 213d20bfe6fSHong Zhang prow = *ptaj++; 214d20bfe6fSHong Zhang pnzj = pi[prow+1] - pi[prow]; 215d20bfe6fSHong Zhang pjj = pj + pi[prow]; 21655a3bba9SHong Zhang /* add non-zero cols of P into the sorted linked list lnk */ 21755a3bba9SHong Zhang ierr = PetscLLAdd(pnzj,pjj,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr); 21855a3bba9SHong Zhang cnzi += nlnk; 219d20bfe6fSHong Zhang } 220d20bfe6fSHong Zhang 221d20bfe6fSHong Zhang /* If free space is not available, make more free space */ 222d20bfe6fSHong Zhang /* Double the amount of total space in the list */ 223d20bfe6fSHong Zhang if (current_space->local_remaining<cnzi) { 224d20bfe6fSHong Zhang ierr = GetMoreSpace(current_space->total_array_size,¤t_space);CHKERRQ(ierr); 225d20bfe6fSHong Zhang } 226d20bfe6fSHong Zhang 227d20bfe6fSHong Zhang /* Copy data into free space, and zero out denserows */ 22855a3bba9SHong Zhang ierr = PetscLLClean(pn,pn,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 229d20bfe6fSHong Zhang current_space->array += cnzi; 230d20bfe6fSHong Zhang current_space->local_used += cnzi; 231d20bfe6fSHong Zhang current_space->local_remaining -= cnzi; 232d20bfe6fSHong Zhang 233d20bfe6fSHong Zhang for (j=0;j<ptanzi;j++) { 234d20bfe6fSHong Zhang ptadenserow[ptasparserow[j]] = 0; 235d20bfe6fSHong Zhang } 236d20bfe6fSHong Zhang /* Aside: Perhaps we should save the pta info for the numerical factorization. */ 237d20bfe6fSHong Zhang /* For now, we will recompute what is needed. */ 238d20bfe6fSHong Zhang ci[i+1] = ci[i] + cnzi; 239d20bfe6fSHong Zhang } 240d20bfe6fSHong Zhang /* nnz is now stored in ci[ptm], column indices are in the list of free space */ 241d20bfe6fSHong Zhang /* Allocate space for cj, initialize cj, and */ 242d20bfe6fSHong Zhang /* destroy list of free space and other temporary array(s) */ 24355a3bba9SHong Zhang ierr = PetscMalloc((ci[pn]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr); 244d20bfe6fSHong Zhang ierr = MakeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 245d20bfe6fSHong Zhang ierr = PetscFree(ptadenserow);CHKERRQ(ierr); 24655a3bba9SHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 247d20bfe6fSHong Zhang 248d20bfe6fSHong Zhang /* Allocate space for ca */ 249d20bfe6fSHong Zhang ierr = PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 250d20bfe6fSHong Zhang ierr = PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));CHKERRQ(ierr); 251d20bfe6fSHong Zhang 252d20bfe6fSHong Zhang /* put together the new matrix */ 253d20bfe6fSHong Zhang ierr = MatCreateSeqAIJWithArrays(A->comm,pn,pn,ci,cj,ca,C);CHKERRQ(ierr); 254d20bfe6fSHong Zhang 255d20bfe6fSHong Zhang /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 256d20bfe6fSHong Zhang /* Since these are PETSc arrays, change flags to free them as necessary. */ 257d20bfe6fSHong Zhang c = (Mat_SeqAIJ *)((*C)->data); 258d20bfe6fSHong Zhang c->freedata = PETSC_TRUE; 259d20bfe6fSHong Zhang c->nonew = 0; 260d20bfe6fSHong Zhang 261d20bfe6fSHong Zhang /* Clean up. */ 262d20bfe6fSHong Zhang ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr); 263eb9c0419SKris Buschelman 264eb9c0419SKris Buschelman PetscFunctionReturn(0); 265eb9c0419SKris Buschelman } 266eb9c0419SKris Buschelman 2673985e5eaSKris Buschelman #include "src/mat/impls/maij/maij.h" 2683985e5eaSKris Buschelman EXTERN_C_BEGIN 2693985e5eaSKris Buschelman #undef __FUNCT__ 2709af31e4aSHong Zhang #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqMAIJ" 27155a3bba9SHong Zhang PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqMAIJ(Mat A,Mat PP,Mat *C) 27255a3bba9SHong Zhang { 2735c66b693SKris Buschelman /* This routine requires testing -- I don't think it works. */ 274dfbe8321SBarry Smith PetscErrorCode ierr; 2753985e5eaSKris Buschelman FreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 2763985e5eaSKris Buschelman Mat_SeqMAIJ *pp=(Mat_SeqMAIJ*)PP->data; 2773985e5eaSKris Buschelman Mat P=pp->AIJ; 2783985e5eaSKris Buschelman Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*p=(Mat_SeqAIJ*)P->data,*c; 279b1d57f15SBarry Smith PetscInt *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj; 280b1d57f15SBarry Smith PetscInt *ci,*cj,*denserow,*sparserow,*ptadenserow,*ptasparserow,*ptaj; 281b1d57f15SBarry Smith PetscInt an=A->N,am=A->M,pn=P->N,pm=P->M,ppdof=pp->dof; 282b1d57f15SBarry Smith PetscInt i,j,k,dof,pdof,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi; 2833985e5eaSKris Buschelman MatScalar *ca; 2843985e5eaSKris Buschelman 2853985e5eaSKris Buschelman PetscFunctionBegin; 2863985e5eaSKris Buschelman /* Start timer */ 2879af31e4aSHong Zhang ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,PP,0,0);CHKERRQ(ierr); 2883985e5eaSKris Buschelman 2893985e5eaSKris Buschelman /* Get ij structure of P^T */ 2903985e5eaSKris Buschelman ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr); 2913985e5eaSKris Buschelman 2923985e5eaSKris Buschelman /* Allocate ci array, arrays for fill computation and */ 2933985e5eaSKris Buschelman /* free space for accumulating nonzero column info */ 294b1d57f15SBarry Smith ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 2953985e5eaSKris Buschelman ci[0] = 0; 2963985e5eaSKris Buschelman 297b1d57f15SBarry Smith ierr = PetscMalloc((2*pn+2*an+1)*sizeof(PetscInt),&ptadenserow);CHKERRQ(ierr); 298b1d57f15SBarry Smith ierr = PetscMemzero(ptadenserow,(2*pn+2*an+1)*sizeof(PetscInt));CHKERRQ(ierr); 2993985e5eaSKris Buschelman ptasparserow = ptadenserow + an; 3003985e5eaSKris Buschelman denserow = ptasparserow + an; 3013985e5eaSKris Buschelman sparserow = denserow + pn; 3023985e5eaSKris Buschelman 3033985e5eaSKris Buschelman /* Set initial free space to be nnz(A) scaled by aspect ratio of P. */ 3043985e5eaSKris Buschelman /* This should be reasonable if sparsity of PtAP is similar to that of A. */ 3053985e5eaSKris Buschelman ierr = GetMoreSpace((ai[am]/pm)*pn,&free_space); 3063985e5eaSKris Buschelman current_space = free_space; 3073985e5eaSKris Buschelman 3083985e5eaSKris Buschelman /* Determine symbolic info for each row of C: */ 3093985e5eaSKris Buschelman for (i=0;i<pn/ppdof;i++) { 3103985e5eaSKris Buschelman ptnzi = pti[i+1] - pti[i]; 3113985e5eaSKris Buschelman ptanzi = 0; 3123985e5eaSKris Buschelman ptJ = ptj + pti[i]; 3133985e5eaSKris Buschelman for (dof=0;dof<ppdof;dof++) { 3143985e5eaSKris Buschelman /* Determine symbolic row of PtA: */ 3153985e5eaSKris Buschelman for (j=0;j<ptnzi;j++) { 3163985e5eaSKris Buschelman arow = ptJ[j] + dof; 3173985e5eaSKris Buschelman anzj = ai[arow+1] - ai[arow]; 3183985e5eaSKris Buschelman ajj = aj + ai[arow]; 3193985e5eaSKris Buschelman for (k=0;k<anzj;k++) { 3203985e5eaSKris Buschelman if (!ptadenserow[ajj[k]]) { 3213985e5eaSKris Buschelman ptadenserow[ajj[k]] = -1; 3223985e5eaSKris Buschelman ptasparserow[ptanzi++] = ajj[k]; 3233985e5eaSKris Buschelman } 3243985e5eaSKris Buschelman } 3253985e5eaSKris Buschelman } 3263985e5eaSKris Buschelman /* Using symbolic info for row of PtA, determine symbolic info for row of C: */ 3273985e5eaSKris Buschelman ptaj = ptasparserow; 3283985e5eaSKris Buschelman cnzi = 0; 3293985e5eaSKris Buschelman for (j=0;j<ptanzi;j++) { 330fe05a634SKris Buschelman pdof = *ptaj%dof; 3313985e5eaSKris Buschelman prow = (*ptaj++)/dof; 3323985e5eaSKris Buschelman pnzj = pi[prow+1] - pi[prow]; 3333985e5eaSKris Buschelman pjj = pj + pi[prow]; 3343985e5eaSKris Buschelman for (k=0;k<pnzj;k++) { 335fe05a634SKris Buschelman if (!denserow[pjj[k]+pdof]) { 336fe05a634SKris Buschelman denserow[pjj[k]+pdof] = -1; 337fe05a634SKris Buschelman sparserow[cnzi++] = pjj[k]+pdof; 3383985e5eaSKris Buschelman } 3393985e5eaSKris Buschelman } 3403985e5eaSKris Buschelman } 3413985e5eaSKris Buschelman 3423985e5eaSKris Buschelman /* sort sparserow */ 3433985e5eaSKris Buschelman ierr = PetscSortInt(cnzi,sparserow);CHKERRQ(ierr); 3443985e5eaSKris Buschelman 3453985e5eaSKris Buschelman /* If free space is not available, make more free space */ 3463985e5eaSKris Buschelman /* Double the amount of total space in the list */ 3473985e5eaSKris Buschelman if (current_space->local_remaining<cnzi) { 3483985e5eaSKris Buschelman ierr = GetMoreSpace(current_space->total_array_size,¤t_space);CHKERRQ(ierr); 3493985e5eaSKris Buschelman } 3503985e5eaSKris Buschelman 3513985e5eaSKris Buschelman /* Copy data into free space, and zero out denserows */ 352b1d57f15SBarry Smith ierr = PetscMemcpy(current_space->array,sparserow,cnzi*sizeof(PetscInt));CHKERRQ(ierr); 3533985e5eaSKris Buschelman current_space->array += cnzi; 3543985e5eaSKris Buschelman current_space->local_used += cnzi; 3553985e5eaSKris Buschelman current_space->local_remaining -= cnzi; 3563985e5eaSKris Buschelman 3573985e5eaSKris Buschelman for (j=0;j<ptanzi;j++) { 3583985e5eaSKris Buschelman ptadenserow[ptasparserow[j]] = 0; 3593985e5eaSKris Buschelman } 3603985e5eaSKris Buschelman for (j=0;j<cnzi;j++) { 3613985e5eaSKris Buschelman denserow[sparserow[j]] = 0; 3623985e5eaSKris Buschelman } 3633985e5eaSKris Buschelman /* Aside: Perhaps we should save the pta info for the numerical factorization. */ 3643985e5eaSKris Buschelman /* For now, we will recompute what is needed. */ 3653985e5eaSKris Buschelman ci[i+1+dof] = ci[i+dof] + cnzi; 3663985e5eaSKris Buschelman } 3673985e5eaSKris Buschelman } 3683985e5eaSKris Buschelman /* nnz is now stored in ci[ptm], column indices are in the list of free space */ 3693985e5eaSKris Buschelman /* Allocate space for cj, initialize cj, and */ 3703985e5eaSKris Buschelman /* destroy list of free space and other temporary array(s) */ 371b1d57f15SBarry Smith ierr = PetscMalloc((ci[pn]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr); 3723985e5eaSKris Buschelman ierr = MakeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 3733985e5eaSKris Buschelman ierr = PetscFree(ptadenserow);CHKERRQ(ierr); 3743985e5eaSKris Buschelman 3753985e5eaSKris Buschelman /* Allocate space for ca */ 3763985e5eaSKris Buschelman ierr = PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 3773985e5eaSKris Buschelman ierr = PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));CHKERRQ(ierr); 3783985e5eaSKris Buschelman 3793985e5eaSKris Buschelman /* put together the new matrix */ 3803985e5eaSKris Buschelman ierr = MatCreateSeqAIJWithArrays(A->comm,pn,pn,ci,cj,ca,C);CHKERRQ(ierr); 3813985e5eaSKris Buschelman 3823985e5eaSKris Buschelman /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 3833985e5eaSKris Buschelman /* Since these are PETSc arrays, change flags to free them as necessary. */ 3843985e5eaSKris Buschelman c = (Mat_SeqAIJ *)((*C)->data); 3853985e5eaSKris Buschelman c->freedata = PETSC_TRUE; 3863985e5eaSKris Buschelman c->nonew = 0; 3873985e5eaSKris Buschelman 3883985e5eaSKris Buschelman /* Clean up. */ 3893985e5eaSKris Buschelman ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr); 3903985e5eaSKris Buschelman 3919af31e4aSHong Zhang ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,PP,0,0);CHKERRQ(ierr); 3923985e5eaSKris Buschelman PetscFunctionReturn(0); 3933985e5eaSKris Buschelman } 3943985e5eaSKris Buschelman EXTERN_C_END 3953985e5eaSKris Buschelman 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 412170ef064SHong Zhang This routine is currently only implemented for pairs of AIJ matrices and classes 413170ef064SHong 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 */ 419e240928fSHong Zhang #undef __FUNCT__ 420e240928fSHong Zhang #define __FUNCT__ "MatPtAPNumeric" 42155a3bba9SHong Zhang PetscErrorCode MatPtAPNumeric(Mat A,Mat P,Mat C) 42255a3bba9SHong Zhang { 423dfbe8321SBarry Smith PetscErrorCode ierr; 424534c1384SKris Buschelman PetscErrorCode (*fA)(Mat,Mat,Mat); 425534c1384SKris Buschelman PetscErrorCode (*fP)(Mat,Mat,Mat); 426eb9c0419SKris Buschelman 427eb9c0419SKris Buschelman PetscFunctionBegin; 428eb9c0419SKris Buschelman 4294482741eSBarry Smith PetscValidHeaderSpecific(A,MAT_COOKIE,1); 430c9780b6fSBarry Smith PetscValidType(A,1); 431eb9c0419SKris Buschelman MatPreallocated(A); 432eb9c0419SKris Buschelman if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 433eb9c0419SKris Buschelman if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 434eb9c0419SKris Buschelman 4354482741eSBarry Smith PetscValidHeaderSpecific(P,MAT_COOKIE,2); 436c9780b6fSBarry Smith PetscValidType(P,2); 437eb9c0419SKris Buschelman MatPreallocated(P); 438eb9c0419SKris Buschelman if (!P->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 439eb9c0419SKris Buschelman if (P->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 440eb9c0419SKris Buschelman 4414482741eSBarry Smith PetscValidHeaderSpecific(C,MAT_COOKIE,3); 442c9780b6fSBarry Smith PetscValidType(C,3); 443eb9c0419SKris Buschelman MatPreallocated(C); 444eb9c0419SKris Buschelman if (!C->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 445eb9c0419SKris Buschelman if (C->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 446eb9c0419SKris Buschelman 44777431f27SBarry Smith if (P->N!=C->M) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %D != %D",P->N,C->M); 44877431f27SBarry Smith if (P->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %D != %D",P->M,A->N); 44977431f27SBarry Smith if (A->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix 'A' must be square, %D != %D",A->M,A->N); 45077431f27SBarry Smith if (P->N!=C->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %D != %D",P->N,C->N); 451eb9c0419SKris Buschelman 452534c1384SKris Buschelman /* For now, we do not dispatch based on the type of A and P */ 453534c1384SKris Buschelman /* When implementations like _SeqAIJ_MAIJ exist, attack the multiple dispatch problem. */ 454534c1384SKris Buschelman fA = A->ops->ptapnumeric; 455534c1384SKris Buschelman if (!fA) SETERRQ1(PETSC_ERR_SUP,"MatPtAPNumeric not supported for A of type %s",A->type_name); 456534c1384SKris Buschelman fP = P->ops->ptapnumeric; 457534c1384SKris Buschelman if (!fP) SETERRQ1(PETSC_ERR_SUP,"MatPtAPNumeric not supported for P of type %s",P->type_name); 458534c1384SKris 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); 4594d3841fdSKris Buschelman 460534c1384SKris Buschelman ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 461534c1384SKris Buschelman ierr = (*fA)(A,P,C);CHKERRQ(ierr); 462534c1384SKris Buschelman ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 463eb9c0419SKris Buschelman 464eb9c0419SKris Buschelman PetscFunctionReturn(0); 465eb9c0419SKris Buschelman } 466eb9c0419SKris Buschelman 467eb9c0419SKris Buschelman #undef __FUNCT__ 4689af31e4aSHong Zhang #define __FUNCT__ "MatPtAPNumeric_SeqAIJ_SeqAIJ" 469dfbe8321SBarry Smith PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat A,Mat P,Mat C) 470dfbe8321SBarry Smith { 471dfbe8321SBarry Smith PetscErrorCode ierr; 472b1d57f15SBarry Smith PetscInt flops=0; 473d20bfe6fSHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 474d20bfe6fSHong Zhang Mat_SeqAIJ *p = (Mat_SeqAIJ *) P->data; 475d20bfe6fSHong Zhang Mat_SeqAIJ *c = (Mat_SeqAIJ *) C->data; 476b1d57f15SBarry Smith PetscInt *ai=a->i,*aj=a->j,*apj,*apjdense,*pi=p->i,*pj=p->j,*pJ=p->j,*pjj; 477b1d57f15SBarry Smith PetscInt *ci=c->i,*cj=c->j,*cjj; 478b1d57f15SBarry Smith PetscInt am=A->M,cn=C->N,cm=C->M; 479b1d57f15SBarry Smith PetscInt i,j,k,anzi,pnzi,apnzj,nextap,pnzj,prow,crow; 480d20bfe6fSHong Zhang MatScalar *aa=a->a,*apa,*pa=p->a,*pA=p->a,*paj,*ca=c->a,*caj; 481eb9c0419SKris Buschelman 482eb9c0419SKris Buschelman PetscFunctionBegin; 483d20bfe6fSHong Zhang /* Allocate temporary array for storage of one row of A*P */ 484b1d57f15SBarry Smith ierr = PetscMalloc(cn*(sizeof(MatScalar)+2*sizeof(PetscInt)),&apa);CHKERRQ(ierr); 485b1d57f15SBarry Smith ierr = PetscMemzero(apa,cn*(sizeof(MatScalar)+2*sizeof(PetscInt)));CHKERRQ(ierr); 486eb9c0419SKris Buschelman 487b1d57f15SBarry Smith apj = (PetscInt *)(apa + cn); 488d20bfe6fSHong Zhang apjdense = apj + cn; 489d20bfe6fSHong Zhang 490d20bfe6fSHong Zhang /* Clear old values in C */ 491d20bfe6fSHong Zhang ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 492d20bfe6fSHong Zhang 493d20bfe6fSHong Zhang for (i=0;i<am;i++) { 494d20bfe6fSHong Zhang /* Form sparse row of A*P */ 495d20bfe6fSHong Zhang anzi = ai[i+1] - ai[i]; 496d20bfe6fSHong Zhang apnzj = 0; 497d20bfe6fSHong Zhang for (j=0;j<anzi;j++) { 498d20bfe6fSHong Zhang prow = *aj++; 499d20bfe6fSHong Zhang pnzj = pi[prow+1] - pi[prow]; 500d20bfe6fSHong Zhang pjj = pj + pi[prow]; 501d20bfe6fSHong Zhang paj = pa + pi[prow]; 502d20bfe6fSHong Zhang for (k=0;k<pnzj;k++) { 503d20bfe6fSHong Zhang if (!apjdense[pjj[k]]) { 504d20bfe6fSHong Zhang apjdense[pjj[k]] = -1; 505d20bfe6fSHong Zhang apj[apnzj++] = pjj[k]; 506d20bfe6fSHong Zhang } 507d20bfe6fSHong Zhang apa[pjj[k]] += (*aa)*paj[k]; 508d20bfe6fSHong Zhang } 509d20bfe6fSHong Zhang flops += 2*pnzj; 510d20bfe6fSHong Zhang aa++; 511d20bfe6fSHong Zhang } 512d20bfe6fSHong Zhang 513d20bfe6fSHong Zhang /* Sort the j index array for quick sparse axpy. */ 514d20bfe6fSHong Zhang ierr = PetscSortInt(apnzj,apj);CHKERRQ(ierr); 515d20bfe6fSHong Zhang 516d20bfe6fSHong Zhang /* Compute P^T*A*P using outer product (P^T)[:,j]*(A*P)[j,:]. */ 517d20bfe6fSHong Zhang pnzi = pi[i+1] - pi[i]; 518d20bfe6fSHong Zhang for (j=0;j<pnzi;j++) { 519d20bfe6fSHong Zhang nextap = 0; 520d20bfe6fSHong Zhang crow = *pJ++; 521d20bfe6fSHong Zhang cjj = cj + ci[crow]; 522d20bfe6fSHong Zhang caj = ca + ci[crow]; 523d20bfe6fSHong Zhang /* Perform sparse axpy operation. Note cjj includes apj. */ 524d20bfe6fSHong Zhang for (k=0;nextap<apnzj;k++) { 525d20bfe6fSHong Zhang if (cjj[k]==apj[nextap]) { 526d20bfe6fSHong Zhang caj[k] += (*pA)*apa[apj[nextap++]]; 527d20bfe6fSHong Zhang } 528d20bfe6fSHong Zhang } 529d20bfe6fSHong Zhang flops += 2*apnzj; 530d20bfe6fSHong Zhang pA++; 531d20bfe6fSHong Zhang } 532d20bfe6fSHong Zhang 533d20bfe6fSHong Zhang /* Zero the current row info for A*P */ 534d20bfe6fSHong Zhang for (j=0;j<apnzj;j++) { 535d20bfe6fSHong Zhang apa[apj[j]] = 0.; 536d20bfe6fSHong Zhang apjdense[apj[j]] = 0; 537d20bfe6fSHong Zhang } 538d20bfe6fSHong Zhang } 539d20bfe6fSHong Zhang 540d20bfe6fSHong Zhang /* Assemble the final matrix and clean up */ 541d20bfe6fSHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 542d20bfe6fSHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 543d20bfe6fSHong Zhang ierr = PetscFree(apa);CHKERRQ(ierr); 544d20bfe6fSHong Zhang ierr = PetscLogFlops(flops);CHKERRQ(ierr); 545d20bfe6fSHong Zhang 546eb9c0419SKris Buschelman PetscFunctionReturn(0); 547eb9c0419SKris Buschelman } 548