/* Defines projective product routines where A is a SeqAIJ 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 "petscbt.h" #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); } #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 #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;i