/* Defines projective product routines where A is a AIJ matrix C = P^T * A * P */ #include "src/mat/impls/aij/seq/aij.h" /*I "petscmat.h" I*/ #include "src/mat/utils/freespace.h" #include "src/mat/impls/aij/mpi/mpiaij.h" EXTERN PetscErrorCode RegisterMatMatMultRoutines_Private(Mat); static PetscEvent MAT_PtAPSymbolic = 0; static PetscEvent MAT_PtAPNumeric = 0; #undef __FUNCT__ #define __FUNCT__ "MatPtAP" /*@ MatPtAP - Creates the matrix projection C = P^T * A * P Collective on Mat Input Parameters: + A - the matrix . P - the projection matrix . scall - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX - fill - expected fill as ratio of nnz(C)/(nnz(A) + nnz(P)) Output Parameters: . C - the product matrix Notes: C will be created and must be destroyed by the user with MatDestroy(). This routine is currently only implemented for pairs of SeqAIJ matrices and classes which inherit from SeqAIJ. C will be of type MATSEQAIJ. Level: intermediate .seealso: MatPtAPSymbolic(),MatPtAPNumeric(),MatMatMult() @*/ PetscErrorCode MatPtAP(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C) { PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(A,MAT_COOKIE,1); PetscValidType(A,1); MatPreallocated(A); if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); PetscValidHeaderSpecific(P,MAT_COOKIE,2); PetscValidType(P,2); MatPreallocated(P); if (!P->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); if (P->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); PetscValidPointer(C,3); if (P->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->M,A->N); if (fill <=0.0) SETERRQ1(PETSC_ERR_ARG_SIZ,"fill=%g must be > 0.0",fill); ierr = PetscLogEventBegin(MAT_PtAP,A,P,0,0);CHKERRQ(ierr); ierr = (*A->ops->matptap)(A,P,scall,fill,C);CHKERRQ(ierr); ierr = PetscLogEventEnd(MAT_PtAP,A,P,0,0);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatPtAP_SeqAIJ_SeqAIJ" PetscErrorCode MatPtAP_SeqAIJ_SeqAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C) { PetscErrorCode ierr; PetscFunctionBegin; if (scall == MAT_INITIAL_MATRIX){ ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ(A,P,fill,C);CHKERRQ(ierr); } ierr = MatPtAPNumeric_SeqAIJ_SeqAIJ(A,P,*C);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatPtAPSymbolic" /* MatPtAPSymbolic - Creates the (i,j) structure of the matrix projection C = P^T * A * P Collective on Mat Input Parameters: + A - the matrix - P - the projection matrix Output Parameters: . C - the (i,j) structure of the product matrix Notes: C will be created and must be destroyed by the user with MatDestroy(). This routine is currently only implemented for pairs of SeqAIJ matrices and classes which inherit from SeqAIJ. C will be of type MATSEQAIJ. The product is computed using this (i,j) structure by calling MatPtAPNumeric(). Level: intermediate .seealso: MatPtAP(),MatPtAPNumeric(),MatMatMultSymbolic() */ PetscErrorCode MatPtAPSymbolic(Mat A,Mat P,PetscReal fill,Mat *C) { PetscErrorCode ierr; char funct[80]; PetscErrorCode (*f)(Mat,Mat,Mat*); PetscFunctionBegin; PetscValidHeaderSpecific(A,MAT_COOKIE,1); PetscValidType(A,1); MatPreallocated(A); if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); PetscValidHeaderSpecific(P,MAT_COOKIE,2); PetscValidType(P,2); MatPreallocated(P); if (!P->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); if (P->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); PetscValidPointer(C,3); if (P->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->M,A->N); if (A->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix 'A' must be square, %d != %d",A->M,A->N); /* Currently only _seqaij_seqaij is implemented, so just query for it. */ /* When other implementations exist, attack the multiple dispatch problem. */ ierr = PetscStrcpy(funct,"MatPtAPSymbolic_seqaij_seqaij");CHKERRQ(ierr); ierr = PetscObjectQueryFunction((PetscObject)P,funct,(PetscVoidFunction)&f);CHKERRQ(ierr); if (!f) SETERRQ1(PETSC_ERR_SUP,"MatPtAPSymbolic is not supported for P of type %s",P->type_name); ierr = PetscObjectQueryFunction((PetscObject)A,funct,(PetscVoidFunction)&f);CHKERRQ(ierr); if (!f) SETERRQ1(PETSC_ERR_SUP,"MatPtAPSymbolic is not supported for A of type %s",A->type_name); ierr = (*f)(A,P,C);CHKERRQ(ierr); PetscFunctionReturn(0); } typedef struct { Mat symAP; } Mat_PtAPstruct; #undef __FUNCT__ #define __FUNCT__ "MatDestroy_SeqAIJ_PtAP" PetscErrorCode MatDestroy_SeqAIJ_PtAP(Mat A) { PetscErrorCode ierr; Mat_PtAPstruct *ptap=(Mat_PtAPstruct*)A->spptr; PetscFunctionBegin; ierr = MatDestroy(ptap->symAP);CHKERRQ(ierr); ierr = PetscFree(ptap);CHKERRQ(ierr); ierr = MatDestroy(A);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqAIJ" PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat P,PetscReal fill,Mat *C) { PetscErrorCode ierr; int *pti,*ptj; Mat Pt,AP; Mat_PtAPstruct *ptap; PetscFunctionBegin; /* create symbolic Pt */ ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr); ierr = MatCreateSeqAIJWithArrays(P->comm,P->N,P->M,pti,ptj,PETSC_NULL,&Pt);CHKERRQ(ierr); /* get symbolic AP=A*P and C=Pt*AP */ ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,P,fill,&AP);CHKERRQ(ierr); ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(Pt,AP,fill,C);CHKERRQ(ierr); /* clean up */ ierr = MatRestoreSymbolicTranspose_SeqAIJ(Pt,&pti,&ptj);CHKERRQ(ierr); ierr = MatDestroy(Pt);CHKERRQ(ierr); /* save symbolic AP - to be used by MatPtAPNumeric_SeqAIJ_SeqAIJ() */ ierr = PetscNew(Mat_PtAPstruct,&ptap);CHKERRQ(ierr); ptap->symAP = AP; (*C)->spptr = (void*)ptap; (*C)->ops->destroy = MatDestroy_SeqAIJ_PtAP; PetscFunctionReturn(0); } #include "src/mat/impls/maij/maij.h" EXTERN_C_BEGIN #undef __FUNCT__ #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqMAIJ" PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqMAIJ(Mat A,Mat PP,Mat *C) { /* This routine requires testing -- I don't think it works. */ PetscErrorCode ierr; FreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; Mat_SeqMAIJ *pp=(Mat_SeqMAIJ*)PP->data; Mat P=pp->AIJ; Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*p=(Mat_SeqAIJ*)P->data,*c; int *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj; int *ci,*cj,*denserow,*sparserow,*ptadenserow,*ptasparserow,*ptaj; int an=A->N,am=A->M,pn=P->N,pm=P->M,ppdof=pp->dof; int i,j,k,dof,pdof,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi; MatScalar *ca; PetscFunctionBegin; /* Start timer */ ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,PP,0,0);CHKERRQ(ierr); /* Get ij structure of P^T */ ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr); /* Allocate ci array, arrays for fill computation and */ /* free space for accumulating nonzero column info */ ierr = PetscMalloc((pn+1)*sizeof(int),&ci);CHKERRQ(ierr); ci[0] = 0; ierr = PetscMalloc((2*pn+2*an+1)*sizeof(int),&ptadenserow);CHKERRQ(ierr); ierr = PetscMemzero(ptadenserow,(2*pn+2*an+1)*sizeof(int));CHKERRQ(ierr); ptasparserow = ptadenserow + an; denserow = ptasparserow + an; sparserow = denserow + pn; /* Set initial free space to be nnz(A) scaled by aspect ratio of P. */ /* This should be reasonable if sparsity of PtAP is similar to that of A. */ ierr = GetMoreSpace((ai[am]/pm)*pn,&free_space); current_space = free_space; /* Determine symbolic info for each row of C: */ for (i=0;ilocal_remainingtotal_array_size,¤t_space);CHKERRQ(ierr); } /* Copy data into free space, and zero out denserows */ ierr = PetscMemcpy(current_space->array,sparserow,cnzi*sizeof(int));CHKERRQ(ierr); current_space->array += cnzi; current_space->local_used += cnzi; current_space->local_remaining -= cnzi; for (j=0;jcomm,pn,pn,ci,cj,ca,C);CHKERRQ(ierr); /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ /* Since these are PETSc arrays, change flags to free them as necessary. */ c = (Mat_SeqAIJ *)((*C)->data); c->freedata = PETSC_TRUE; c->nonew = 0; /* Clean up. */ ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr); ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,PP,0,0);CHKERRQ(ierr); PetscFunctionReturn(0); } EXTERN_C_END #undef __FUNCT__ #define __FUNCT__ "MatPtAPNumeric" /* MatPtAPNumeric - Computes the matrix projection C = P^T * A * P Collective on Mat Input Parameters: + A - the matrix - P - the projection matrix Output Parameters: . C - the product matrix Notes: C must have been created by calling MatPtAPSymbolic and must be destroyed by the user using MatDeatroy(). This routine is currently only implemented for pairs of SeqAIJ matrices and classes which inherit from SeqAIJ. C will be of type MATSEQAIJ. Level: intermediate .seealso: MatPtAP(),MatPtAPSymbolic(),MatMatMultNumeric() */ PetscErrorCode MatPtAPNumeric(Mat A,Mat P,Mat C) { PetscErrorCode ierr; char funct[80]; PetscErrorCode (*f)(Mat,Mat,Mat); PetscFunctionBegin; PetscValidHeaderSpecific(A,MAT_COOKIE,1); PetscValidType(A,1); MatPreallocated(A); if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); PetscValidHeaderSpecific(P,MAT_COOKIE,2); PetscValidType(P,2); MatPreallocated(P); if (!P->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); if (P->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); PetscValidHeaderSpecific(C,MAT_COOKIE,3); PetscValidType(C,3); MatPreallocated(C); if (!C->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); if (C->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); if (P->N!=C->M) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->N,C->M); if (P->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->M,A->N); if (A->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix 'A' must be square, %d != %d",A->M,A->N); if (P->N!=C->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",P->N,C->N); /* Currently only _seqaij_seqaij is implemented, so just query for it. */ /* When other implementations exist, attack the multiple dispatch problem. */ ierr = PetscStrcpy(funct,"MatPtAPNumeric_seqaij_seqaij");CHKERRQ(ierr); ierr = PetscObjectQueryFunction((PetscObject)P,funct,(PetscVoidFunction)&f);CHKERRQ(ierr); if (!f) SETERRQ1(PETSC_ERR_SUP,"MatPtAPNumeric is not supported for P of type %s",P->type_name); ierr = PetscObjectQueryFunction((PetscObject)A,funct,(PetscVoidFunction)&f);CHKERRQ(ierr); if (!f) SETERRQ1(PETSC_ERR_SUP,"MatPtAPNumeric is not supported for A of type %s",A->type_name); ierr = (*f)(A,P,C);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatPtAPNumeric_SeqAIJ_SeqAIJ" PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat A,Mat P,Mat C) { PetscErrorCode ierr; int flops=0; Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; Mat_SeqAIJ *p = (Mat_SeqAIJ *) P->data; Mat_SeqAIJ *c = (Mat_SeqAIJ *) C->data; int *ai=a->i,*aj=a->j,*pi=p->i,*pj=p->j,*pJ=p->j,*pjj; int *ci=c->i,*cj=c->j,*cjj; int am=A->M,cn=C->N,cm=C->M; int i,j,k,anzi,pnzi,apnzj,nextap,pnzj,prow,crow; MatScalar *aa=a->a,*apa,*pa=p->a,*pA=p->a,*paj,*ca=c->a,*caj; Mat_PtAPstruct *ptap=(Mat_PtAPstruct*)C->spptr; Mat_SeqAIJ *ap = (Mat_SeqAIJ *)(ptap->symAP)->data; int *api=ap->i,*apj=ap->j,apj_nextap; PetscFunctionBegin; /* Allocate temporary array for storage of one row of A*P */ ierr = PetscMalloc(cn*sizeof(MatScalar),&apa);CHKERRQ(ierr); ierr = PetscMemzero(apa,cn*sizeof(MatScalar));CHKERRQ(ierr); /* Clear old values in C */ ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); for (i=0;ij + api[i]; apnzj = api[i+1] - api[i]; pnzi = pi[i+1] - pi[i]; for (j=0;j