/* 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" #include "petscbt.h" #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; PetscErrorCode (*fA)(Mat,Mat,MatReuse,PetscReal,Mat *); PetscErrorCode (*fP)(Mat,Mat,MatReuse,PetscReal,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 (fill <=0.0) SETERRQ1(PETSC_ERR_ARG_SIZ,"fill=%g must be > 0.0",fill); /* For now, we do not dispatch based on the type of A and P */ /* When implementations like _SeqAIJ_MAIJ exist, attack the multiple dispatch problem. */ fA = A->ops->ptap; if (!fA) SETERRQ1(PETSC_ERR_SUP,"MatPtAP not supported for A of type %s",A->type_name); fP = P->ops->ptap; if (!fP) SETERRQ1(PETSC_ERR_SUP,"MatPtAP not supported for P of type %s",P->type_name); if (fP!=fA) SETERRQ2(PETSC_ERR_ARG_INCOMP,"MatPtAP requires A, %s, to be compatible with P, %s",A->type_name,P->type_name); ierr = PetscLogEventBegin(MAT_PtAP,A,P,0,0);CHKERRQ(ierr); ierr = (*fA)(A,P,scall,fill,C);CHKERRQ(ierr); ierr = PetscLogEventEnd(MAT_PtAP,A,P,0,0);CHKERRQ(ierr); PetscFunctionReturn(0); } EXTERN PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_Local(Mat,Mat,Mat,PetscInt,PetscReal,Mat*); EXTERN PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_Local(Mat,Mat,Mat,PetscInt,Mat C); EXTERN PetscErrorCode MatDestroy_MPIAIJ(Mat); #undef __FUNCT__ #define __FUNCT__ "MatDestroy_MPIAIJ_PtAP" PetscErrorCode MatDestroy_MPIAIJ_PtAP(Mat A) { PetscErrorCode ierr; Mat_PtAPMPI *ptap; PetscObjectContainer container; PetscFunctionBegin; ierr = PetscObjectQuery((PetscObject)A,"MatPtAPMPI",(PetscObject *)&container);CHKERRQ(ierr); if (container) { ierr = PetscObjectContainerGetPointer(container,(void **)&ptap);CHKERRQ(ierr); ierr = MatDestroyMatrices(1,&ptap->ploc);CHKERRQ(ierr); ierr = MatDestroy(ptap->P_oth);CHKERRQ(ierr); ierr = PetscObjectContainerDestroy(container);CHKERRQ(ierr); ierr = PetscObjectCompose((PetscObject)A,"MatPtApMPI",0);CHKERRQ(ierr); } ierr = PetscFree(ptap);CHKERRQ(ierr); ierr = MatDestroy_MPIAIJ(A);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatPtAP_MPIAIJ_MPIAIJ" PetscErrorCode MatPtAP_MPIAIJ_MPIAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C) { PetscErrorCode ierr; Mat P_loc,P_oth,*ploc; PetscInt prstart,low; IS isrow,iscol; Mat_PtAPMPI *ptap; PetscObjectContainer container; PetscFunctionBegin; if (scall == MAT_INITIAL_MATRIX){ ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ ierr = MatGetBrowsOfAoCols(A,P,MAT_INITIAL_MATRIX,PETSC_NULL,PETSC_NULL,&prstart,&P_oth);CHKERRQ(ierr); /* get P_loc by taking all local rows of P */ ierr = MatGetOwnershipRange(P,&low,PETSC_NULL);CHKERRQ(ierr); ierr = ISCreateStride(PETSC_COMM_SELF,P->m,low,1,&isrow);CHKERRQ(ierr); ierr = ISCreateStride(PETSC_COMM_SELF,P->N,0,1,&iscol);CHKERRQ(ierr); ierr = MatGetSubMatrices(P,1,&isrow,&iscol,MAT_INITIAL_MATRIX,&ploc);CHKERRQ(ierr); P_loc = ploc[0]; ierr = ISDestroy(isrow);CHKERRQ(ierr); ierr = ISDestroy(iscol);CHKERRQ(ierr); /* attach the supporting struct to P for reuse */ ierr = PetscNew(Mat_PtAPMPI,&ptap);CHKERRQ(ierr); ptap->ploc = ploc; ptap->P_oth = P_oth; ptap->prstart = prstart; P->ops->destroy = MatDestroy_MPIAIJ_PtAP; ierr = PetscObjectContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr); ierr = PetscObjectContainerSetPointer(container,ptap);CHKERRQ(ierr); ierr = PetscObjectCompose((PetscObject)P,"MatPtAPMPI",(PetscObject)container);CHKERRQ(ierr); /* now, compute symbolic product */ ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ(A,P,fill,C);CHKERRQ(ierr);/* numeric product is computed as well */ ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); } else if (scall == MAT_REUSE_MATRIX){ ierr = PetscObjectQuery((PetscObject)P,"MatPtAPMPI",(PetscObject *)&container);CHKERRQ(ierr); if (container) { ierr = PetscObjectContainerGetPointer(container,(void **)&ptap);CHKERRQ(ierr); } else { SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Matrix P does not posses an object container"); } /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ /* improve! */ ierr = MatDestroy(ptap->P_oth);CHKERRQ(ierr); ierr = MatGetBrowsOfAoCols(A,P,MAT_INITIAL_MATRIX,PETSC_NULL,PETSC_NULL,&prstart,&ptap->P_oth);CHKERRQ(ierr); /* get P_loc by taking all local rows of P */ ierr = MatGetOwnershipRange(P,&low,PETSC_NULL);CHKERRQ(ierr); ierr = ISCreateStride(PETSC_COMM_SELF,P->m,low,1,&isrow);CHKERRQ(ierr); ierr = ISCreateStride(PETSC_COMM_SELF,P->N,0,1,&iscol);CHKERRQ(ierr); ierr = MatGetSubMatrices(P,1,&isrow,&iscol,MAT_REUSE_MATRIX,&ptap->ploc);CHKERRQ(ierr); ierr = ISDestroy(isrow);CHKERRQ(ierr); ierr = ISDestroy(iscol);CHKERRQ(ierr); } else { SETERRQ1(PETSC_ERR_ARG_WRONG,"Invalid MatReuse %d",scall); } /* now, compute numeric product */ ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); ierr = MatPtAPNumeric_MPIAIJ_MPIAIJ(A,P,*C);CHKERRQ(ierr); ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatPtAPSymbolic_MPIAIJ_MPIAIJ" PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C) { PetscErrorCode ierr; Mat C_seq; Mat_PtAPMPI *ptap; PetscObjectContainer container; PetscFunctionBegin; ierr = PetscObjectQuery((PetscObject)P,"MatPtAPMPI",(PetscObject *)&container);CHKERRQ(ierr); if (container) { ierr = PetscObjectContainerGetPointer(container,(void **)&ptap);CHKERRQ(ierr); } else { SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Matrix P does not posses an object container"); } /* compute C_seq = P_loc^T * A_loc * P_seq */ ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ_Local(A,ptap->ploc[0],ptap->P_oth,ptap->prstart,fill,&C_seq);CHKERRQ(ierr); /* add C_seq into mpi C */ ierr = MatMerge_SeqsToMPISymbolic(A->comm,C_seq,P->n,P->n,C);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatPtAPNumeric_MPIAIJ_MPIAIJ" PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C) { PetscErrorCode ierr; Mat_Merge_SeqsToMPI *merge; Mat_PtAPMPI *ptap; PetscObjectContainer cont_merge,cont_ptap; PetscFunctionBegin; ierr = PetscObjectQuery((PetscObject)C,"MatMergeSeqsToMPI",(PetscObject *)&cont_merge);CHKERRQ(ierr); if (cont_merge) { ierr = PetscObjectContainerGetPointer(cont_merge,(void **)&merge);CHKERRQ(ierr); } else { SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Matrix C does not posses an object container"); } ierr = PetscObjectQuery((PetscObject)P,"MatPtAPMPI",(PetscObject *)&cont_ptap);CHKERRQ(ierr); if (cont_ptap) { ierr = PetscObjectContainerGetPointer(cont_ptap,(void **)&ptap);CHKERRQ(ierr); } else { SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Matrix P does not posses an object container"); } ierr = MatPtAPNumeric_MPIAIJ_MPIAIJ_Local(A,ptap->ploc[0],ptap->P_oth,ptap->prstart,merge->C_seq);CHKERRQ(ierr); ierr = MatMerge_SeqsToMPINumeric(merge->C_seq,C);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 = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ(A,P,fill,C);CHKERRQ(ierr); ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); } ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); ierr = MatPtAPNumeric_SeqAIJ_SeqAIJ(A,P,*C);CHKERRQ(ierr); ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); PetscFunctionReturn(0); } /* 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() */ #undef __FUNCT__ #define __FUNCT__ "MatPtAPSymbolic" PetscErrorCode MatPtAPSymbolic(Mat A,Mat P,PetscReal fill,Mat *C) { PetscErrorCode ierr; PetscErrorCode (*fA)(Mat,Mat,PetscReal,Mat*); PetscErrorCode (*fP)(Mat,Mat,PetscReal,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); /* For now, we do not dispatch based on the type of A and P */ /* When implementations like _SeqAIJ_MAIJ exist, attack the multiple dispatch problem. */ fA = A->ops->ptapsymbolic; if (!fA) SETERRQ1(PETSC_ERR_SUP,"MatPtAPSymbolic not supported for A of type %s",A->type_name); fP = P->ops->ptapsymbolic; if (!fP) SETERRQ1(PETSC_ERR_SUP,"MatPtAPSymbolic not supported for P of type %s",P->type_name); if (fP!=fA) SETERRQ2(PETSC_ERR_ARG_INCOMP,"MatPtAPSymbolic requires A, %s, to be compatible with P, %s",A->type_name,P->type_name); ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); ierr = (*fA)(A,P,fill,C);CHKERRQ(ierr); ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); PetscFunctionReturn(0); } typedef struct { Mat symAP; } Mat_PtAPstruct; EXTERN PetscErrorCode MatDestroy_SeqAIJ(Mat); #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_SeqAIJ(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; FreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*p = (Mat_SeqAIJ*)P->data,*c; PetscInt *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj; PetscInt *ci,*cj,*ptadenserow,*ptasparserow,*ptaj; PetscInt an=A->N,am=A->M,pn=P->N,pm=P->M; PetscInt i,j,k,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi,nlnk,*lnk; MatScalar *ca; PetscBT lnkbt; PetscFunctionBegin; /* Get ij structure of P^T */ ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr); ptJ=ptj; /* Allocate ci array, arrays for fill computation and */ /* free space for accumulating nonzero column info */ ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr); ci[0] = 0; ierr = PetscMalloc((2*an+1)*sizeof(PetscInt),&ptadenserow);CHKERRQ(ierr); ierr = PetscMemzero(ptadenserow,(2*an+1)*sizeof(PetscInt));CHKERRQ(ierr); ptasparserow = ptadenserow + an; /* create and initialize a linked list */ nlnk = pn+1; ierr = PetscLLCreate(pn,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr); /* 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 = PetscLLClean(pn,pn,cnzi,lnk,current_space->array,lnkbt);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); 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; PetscInt *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj; PetscInt *ci,*cj,*denserow,*sparserow,*ptadenserow,*ptasparserow,*ptaj; PetscInt an=A->N,am=A->M,pn=P->N,pm=P->M,ppdof=pp->dof; PetscInt 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(PetscInt),&ci);CHKERRQ(ierr); ci[0] = 0; ierr = PetscMalloc((2*pn+2*an+1)*sizeof(PetscInt),&ptadenserow);CHKERRQ(ierr); ierr = PetscMemzero(ptadenserow,(2*pn+2*an+1)*sizeof(PetscInt));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(PetscInt));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 /* 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 AIJ matrices and classes which inherit from AIJ. C will be of type MATAIJ. Level: intermediate .seealso: MatPtAP(),MatPtAPSymbolic(),MatMatMultNumeric() */ #undef __FUNCT__ #define __FUNCT__ "MatPtAPNumeric" PetscErrorCode MatPtAPNumeric(Mat A,Mat P,Mat C) { PetscErrorCode ierr; PetscErrorCode (*fA)(Mat,Mat,Mat); PetscErrorCode (*fP)(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); /* For now, we do not dispatch based on the type of A and P */ /* When implementations like _SeqAIJ_MAIJ exist, attack the multiple dispatch problem. */ fA = A->ops->ptapnumeric; if (!fA) SETERRQ1(PETSC_ERR_SUP,"MatPtAPNumeric not supported for A of type %s",A->type_name); fP = P->ops->ptapnumeric; if (!fP) SETERRQ1(PETSC_ERR_SUP,"MatPtAPNumeric not supported for P of type %s",P->type_name); if (fP!=fA) SETERRQ2(PETSC_ERR_ARG_INCOMP,"MatPtAPNumeric requires A, %s, to be compatible with P, %s",A->type_name,P->type_name); ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); ierr = (*fA)(A,P,C);CHKERRQ(ierr); ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatPtAPNumeric_SeqAIJ_SeqAIJ" PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat A,Mat P,Mat C) { PetscErrorCode ierr; PetscInt flops=0; Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; Mat_SeqAIJ *p = (Mat_SeqAIJ *) P->data; Mat_SeqAIJ *c = (Mat_SeqAIJ *) C->data; PetscInt *ai=a->i,*aj=a->j,*apj,*apjdense,*pi=p->i,*pj=p->j,*pJ=p->j,*pjj; PetscInt *ci=c->i,*cj=c->j,*cjj; PetscInt am=A->M,cn=C->N,cm=C->M; PetscInt 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; PetscFunctionBegin; /* Allocate temporary array for storage of one row of A*P */ ierr = PetscMalloc(cn*(sizeof(MatScalar)+2*sizeof(PetscInt)),&apa);CHKERRQ(ierr); ierr = PetscMemzero(apa,cn*(sizeof(MatScalar)+2*sizeof(PetscInt)));CHKERRQ(ierr); apj = (PetscInt *)(apa + cn); apjdense = apj + cn; /* Clear old values in C */ ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); for (i=0;idata; Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; Mat_SeqAIJ *c = (Mat_SeqAIJ *) C->data; Mat_SeqAIJ *p_loc=(Mat_SeqAIJ*)P_loc->data,*p_oth=(Mat_SeqAIJ*)P_oth->data; PetscInt *pi_loc=p_loc->i,*pj_loc=p_loc->j,*pi_oth=p_oth->i,*pj_oth=p_oth->j; PetscInt *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,*apj,*apjdense,cstart=a->cstart,cend=a->cend,col; PetscInt *pJ=pj_loc,*pjj; PetscInt *ci=c->i,*cj=c->j,*cjj; PetscInt am=A->m,cn=C->N,cm=C->M; PetscInt i,j,k,anzi,pnzi,apnzj,nextap,pnzj,prow,crow; MatScalar *ada=ad->a,*aoa=ao->a,*apa,*paj,*ca=c->a,*caj; MatScalar *pa_loc=p_loc->a,*pA=pa_loc,*pa_oth=p_oth->a; PetscFunctionBegin; if (!logkey_matptapnumeric_local) { ierr = PetscLogEventRegister(&logkey_matptapnumeric_local,"MatPtAPNumeric_MPIAIJ_MPIAIJ_Local",MAT_COOKIE);CHKERRQ(ierr); } ierr = PetscLogEventBegin(logkey_matptapnumeric_local,A,0,0,0);CHKERRQ(ierr); /* Allocate temporary array for storage of one row of A*P */ ierr = PetscMalloc(cn*(sizeof(MatScalar)+2*sizeof(PetscInt)),&apa);CHKERRQ(ierr); ierr = PetscMemzero(apa,cn*(sizeof(MatScalar)+2*sizeof(PetscInt)));CHKERRQ(ierr); apj = (PetscInt *)(apa + cn); apjdense = apj + cn; /* Clear old values in C */ ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); for (i=0;igarray[*aoj]; if (col < cstart){ prow = *aoj; } else if (col >= cend){ prow = *aoj; } else { SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Off-diagonal portion of A has wrong column map"); } aoj++; pnzj = pi_oth[prow+1] - pi_oth[prow]; pjj = pj_oth + pi_oth[prow]; paj = pa_oth + pi_oth[prow]; for (k=0;kdata; Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*c; PetscInt *pti,*ptj,*ptJ,*ajj,*pjj; PetscInt *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,nnz,cstart=a->cstart,cend=a->cend,col; PetscInt *ci,*cj,*ptaj; PetscInt an=A->N,am=A->m,pN=P_loc->N; PetscInt i,j,k,ptnzi,arow,anzj,prow,pnzj,cnzi; PetscInt pm=P_loc->m,nlnk,*lnk; MatScalar *ca; PetscBT lnkbt; PetscInt prend,nprow_loc,nprow_oth; PetscInt *ptadenserow_loc,*ptadenserow_oth,*ptasparserow_loc,*ptasparserow_oth; Mat_SeqAIJ *p_loc=(Mat_SeqAIJ*)P_loc->data,*p_oth=(Mat_SeqAIJ*)P_oth->data; PetscInt *pi_loc=p_loc->i,*pj_loc=p_loc->j,*pi_oth=p_oth->i,*pj_oth=p_oth->j; PetscFunctionBegin; if (!logkey_matptapsymbolic_local) { ierr = PetscLogEventRegister(&logkey_matptapsymbolic_local,"MatPtAPSymbolic_MPIAIJ_MPIAIJ_Local",MAT_COOKIE); } ierr = PetscLogEventBegin(logkey_matptapsymbolic_local,A,P_loc,P_oth,0);CHKERRQ(ierr); prend = prstart + pm; /* get ij structure of P_loc^T */ ierr = MatGetSymbolicTranspose_SeqAIJ(P_loc,&pti,&ptj);CHKERRQ(ierr); ptJ=ptj; /* Allocate ci array, arrays for fill computation and */ /* free space for accumulating nonzero column info */ ierr = PetscMalloc((pN+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr); ci[0] = 0; ierr = PetscMalloc((4*an+1)*sizeof(PetscInt),&ptadenserow_loc);CHKERRQ(ierr); ierr = PetscMemzero(ptadenserow_loc,(4*an+1)*sizeof(PetscInt));CHKERRQ(ierr); ptasparserow_loc = ptadenserow_loc + an; ptadenserow_oth = ptasparserow_loc + an; ptasparserow_oth = ptadenserow_oth + an; /* create and initialize a linked list */ nlnk = pN+1; ierr = PetscLLCreate(pN,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr); /* Set initial free space to be nnz(A) scaled by aspect ratio of P (P->M=A->N). */ /* This should be reasonable if sparsity of PtAP is similar to that of A. */ nnz = adi[am] + aoi[am]; ierr = GetMoreSpace((PetscInt)(fill*nnz/A->N)*pN,&free_space); current_space = free_space; /* determine symbolic info for each row of C: */ for (i=0;igarray[ajj[k]]; /* global col */ if (col < cstart){ col = ajj[k]; } else if (col >= cend){ col = ajj[k] + pm; } else { SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Off-diagonal portion of A has wrong column map"); } if (!ptadenserow_oth[col]) { ptadenserow_oth[col] = -1; ptasparserow_oth[nprow_oth++] = col; } } } /* using symbolic info of local PtA, determine symbolic info for row of C: */ cnzi = 0; /* rows of P_loc */ ptaj = ptasparserow_loc; for (j=0; j= prend) prow -= pm; /* rm */ pnzj = pi_oth[prow+1] - pi_oth[prow]; pjj = pj_oth + pi_oth[prow]; /* add non-zero cols of P into the sorted linked list lnk */ ierr = PetscLLAdd(pnzj,pjj,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr); cnzi += nlnk; } /* If free space is not available, make more free space */ /* Double the amount of total space in the list */ if (current_space->local_remainingtotal_array_size,¤t_space);CHKERRQ(ierr); } /* Copy data into free space, and zero out denserows */ ierr = PetscLLClean(pN,pN,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr); current_space->array += cnzi; current_space->local_used += cnzi; current_space->local_remaining -= cnzi; for (j=0;jdata); c->freedata = PETSC_TRUE; c->nonew = 0; /* Clean up. */ ierr = MatRestoreSymbolicTranspose_SeqAIJ(P_loc,&pti,&ptj);CHKERRQ(ierr); ierr = PetscLogEventEnd(logkey_matptapsymbolic_local,A,P_loc,P_oth,0);CHKERRQ(ierr); PetscFunctionReturn(0); }