/* Defines matrix-matrix product routines C = A^T * B */ #include <../src/mat/impls/aij/seq/aij.h> /*I "petscmat.h" I*/ #include <../src/mat/impls/dense/seq/dense.h> #undef __FUNCT__ #define __FUNCT__ "MatDestroy_SeqDense_MatTransMatMult" PetscErrorCode MatDestroy_SeqDense_MatTransMatMult(Mat A) { PetscErrorCode ierr; Mat_SeqDense *a = (Mat_SeqDense*)A->data; Mat_MatTransMatMult *atb = a->atb; PetscFunctionBegin; ierr = MatDestroy(&atb->mA);CHKERRQ(ierr); ierr = VecDestroy(&atb->bt);CHKERRQ(ierr); ierr = VecDestroy(&atb->ct);CHKERRQ(ierr); ierr = (atb->destroy)(A);CHKERRQ(ierr); ierr = PetscFree(atb);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatTransposeMatMult_SeqAIJ_SeqDense" PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) { PetscErrorCode ierr; PetscFunctionBegin; if (scall == MAT_INITIAL_MATRIX) { ierr = PetscLogEventBegin(MAT_TransposeMatMultSymbolic,A,B,0,0);CHKERRQ(ierr); ierr = MatTransposeMatMultSymbolic_SeqAIJ_SeqDense(A,B,fill,C);CHKERRQ(ierr); ierr = PetscLogEventEnd(MAT_TransposeMatMultSymbolic,A,B,0,0);CHKERRQ(ierr); } ierr = PetscLogEventBegin(MAT_TransposeMatMultNumeric,A,B,0,0);CHKERRQ(ierr); ierr = MatTransposeMatMultNumeric_SeqAIJ_SeqDense(A,B,*C);CHKERRQ(ierr); ierr = PetscLogEventEnd(MAT_TransposeMatMultNumeric,A,B,0,0);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatTransposeMatMultSymbolic_SeqAIJ_SeqDense" PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C) { PetscErrorCode ierr; PetscInt m=A->rmap->n,n=A->cmap->n,BN=B->cmap->N; Mat_MatTransMatMult *atb; Mat Cdense; Vec bt,ct; Mat_SeqDense *c; PetscFunctionBegin; ierr = PetscNew(&atb);CHKERRQ(ierr); /* create output dense matrix C = A^T*B */ ierr = MatCreate(PETSC_COMM_SELF,&Cdense);CHKERRQ(ierr); ierr = MatSetSizes(Cdense,n,BN,n,BN);CHKERRQ(ierr); ierr = MatSetType(Cdense,MATSEQDENSE);CHKERRQ(ierr); ierr = MatSeqDenseSetPreallocation(Cdense,NULL);CHKERRQ(ierr); /* create vectors bt and ct to hold locally transposed arrays of B and C */ ierr = VecCreate(PETSC_COMM_SELF,&bt);CHKERRQ(ierr); ierr = VecSetSizes(bt,m*BN,m*BN);CHKERRQ(ierr); ierr = VecSetType(bt,VECSTANDARD);CHKERRQ(ierr); ierr = VecCreate(PETSC_COMM_SELF,&ct);CHKERRQ(ierr); ierr = VecSetSizes(ct,n*BN,n*BN);CHKERRQ(ierr); ierr = VecSetType(ct,VECSTANDARD);CHKERRQ(ierr); atb->bt = bt; atb->ct = ct; *C = Cdense; c = (Mat_SeqDense*)Cdense->data; c->atb = atb; atb->destroy = Cdense->ops->destroy; Cdense->ops->destroy = MatDestroy_SeqDense_MatTransMatMult; PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatTransposeMatMultNumeric_SeqAIJ_SeqDense" PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C) { PetscErrorCode ierr; PetscInt i,j,k,m=A->rmap->n,n=A->cmap->n,BN=B->cmap->N; PetscScalar *Barray,*Carray,*btarray,*ctarray; Mat_SeqDense *c=(Mat_SeqDense*)C->data; Mat_MatTransMatMult *atb=c->atb; Vec bt=atb->bt,ct=atb->ct; PetscFunctionBegin; /* create MAIJ matrix mA from A -- should be done in symbolic phase */ ierr = MatDestroy(&atb->mA);CHKERRQ(ierr); ierr = MatCreateMAIJ(A,BN,&atb->mA);CHKERRQ(ierr); /* transpose local arry of B, then copy it to vector bt */ ierr = MatDenseGetArray(B,&Barray);CHKERRQ(ierr); ierr = VecGetArray(bt,&btarray);CHKERRQ(ierr); k=0; for (j=0; jmA,bt,ct);CHKERRQ(ierr); /* transpose local arry of ct to matrix C */ ierr = MatDenseGetArray(C,&Carray);CHKERRQ(ierr); ierr = VecGetArray(ct,&ctarray);CHKERRQ(ierr); k = 0; for (j=0; j