18949adfdSHong Zhang 28949adfdSHong Zhang /* 38949adfdSHong Zhang Defines matrix-matrix product routines for pairs of MPIAIJ matrices 48949adfdSHong Zhang C = A^T * B 58949adfdSHong Zhang The routines are slightly modified from MatTransposeMatMultxxx_SeqAIJ_SeqDense(). 68949adfdSHong Zhang */ 78949adfdSHong Zhang #include <../src/mat/impls/aij/seq/aij.h> /*I "petscmat.h" I*/ 88949adfdSHong Zhang #include <../src/mat/impls/aij/mpi/mpiaij.h> 98949adfdSHong Zhang #include <../src/mat/impls/dense/mpi/mpidense.h> 108949adfdSHong Zhang 11*6718818eSStefano Zampini PetscErrorCode MatDestroy_MPIDense_MatTransMatMult(void *data) 128949adfdSHong Zhang { 138949adfdSHong Zhang PetscErrorCode ierr; 14*6718818eSStefano Zampini Mat_MatTransMatMult *atb = (Mat_MatTransMatMult*)data; 158949adfdSHong Zhang 168949adfdSHong Zhang PetscFunctionBegin; 178949adfdSHong Zhang ierr = MatDestroy(&atb->mA);CHKERRQ(ierr); 188949adfdSHong Zhang ierr = VecDestroy(&atb->bt);CHKERRQ(ierr); 198949adfdSHong Zhang ierr = VecDestroy(&atb->ct);CHKERRQ(ierr); 208949adfdSHong Zhang ierr = PetscFree(atb);CHKERRQ(ierr); 218949adfdSHong Zhang PetscFunctionReturn(0); 228949adfdSHong Zhang } 238949adfdSHong Zhang 24*6718818eSStefano Zampini static PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIDense(Mat,Mat,Mat); 25*6718818eSStefano Zampini 26*6718818eSStefano Zampini PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIDense(Mat A,Mat B,PetscReal fill,Mat C) 278949adfdSHong Zhang { 288949adfdSHong Zhang PetscErrorCode ierr; 298949adfdSHong Zhang Mat_MatTransMatMult *atb; 30*6718818eSStefano Zampini PetscBool cisdense; 318949adfdSHong Zhang 328949adfdSHong Zhang PetscFunctionBegin; 33*6718818eSStefano Zampini MatCheckProduct(C,4); 34*6718818eSStefano Zampini if (C->product->data) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Extra product struct not empty"); 358949adfdSHong Zhang 368949adfdSHong Zhang /* create output dense matrix C = A^T*B */ 37*6718818eSStefano Zampini ierr = MatSetSizes(C,A->cmap->n,B->cmap->n,A->cmap->N,B->cmap->N);CHKERRQ(ierr); 38*6718818eSStefano Zampini ierr = PetscObjectTypeCompareAny((PetscObject)C,&cisdense,MATMPIDENSE,MATMPIDENSECUDA,"");CHKERRQ(ierr); 39*6718818eSStefano Zampini if (!cisdense) { 40*6718818eSStefano Zampini ierr = MatSetType(C,((PetscObject)B)->type_name);CHKERRQ(ierr); 41*6718818eSStefano Zampini } 424eeee6adSStefano Zampini ierr = MatSetUp(C);CHKERRQ(ierr); 438949adfdSHong Zhang 44*6718818eSStefano Zampini /* create additional data structure for the product */ 45*6718818eSStefano Zampini ierr = PetscNew(&atb);CHKERRQ(ierr); 46*6718818eSStefano Zampini if (B->cmap->N) { 47*6718818eSStefano Zampini ierr = MatCreateMAIJ(A,B->cmap->N,&atb->mA);CHKERRQ(ierr); 48*6718818eSStefano Zampini ierr = MatCreateVecs(atb->mA,&atb->ct,&atb->bt);CHKERRQ(ierr); 49*6718818eSStefano Zampini } 50*6718818eSStefano Zampini C->product->data = atb; 51*6718818eSStefano Zampini C->product->destroy = MatDestroy_MPIDense_MatTransMatMult; 528949adfdSHong Zhang 534222ddf1SHong Zhang C->ops->transposematmultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIDense; 548949adfdSHong Zhang PetscFunctionReturn(0); 558949adfdSHong Zhang } 568949adfdSHong Zhang 57*6718818eSStefano Zampini static PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIDense(Mat A,Mat B,Mat C) 588949adfdSHong Zhang { 598949adfdSHong Zhang PetscErrorCode ierr; 601683a169SBarry Smith const PetscScalar *Barray,*ctarray; 611683a169SBarry Smith PetscScalar *Carray,*btarray; 62905b3b74SHong Zhang Mat_MPIDense *b=(Mat_MPIDense*)B->data,*c=(Mat_MPIDense*)C->data; 63905b3b74SHong Zhang Mat_SeqDense *bseq=(Mat_SeqDense*)(b->A)->data,*cseq=(Mat_SeqDense*)(c->A)->data; 64905b3b74SHong Zhang PetscInt i,j,m=A->rmap->n,n=A->cmap->n,ldb=bseq->lda,BN=B->cmap->N,ldc=cseq->lda; 65*6718818eSStefano Zampini Mat_MatTransMatMult *atb; 66*6718818eSStefano Zampini Vec bt,ct; 678949adfdSHong Zhang 688949adfdSHong Zhang PetscFunctionBegin; 69*6718818eSStefano Zampini MatCheckProduct(C,3); 70*6718818eSStefano Zampini atb=(Mat_MatTransMatMult *)C->product->data; 71*6718818eSStefano Zampini if (!atb) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing product struct"); 72*6718818eSStefano Zampini if (!BN) { 73*6718818eSStefano Zampini ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 74*6718818eSStefano Zampini ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 75*6718818eSStefano Zampini PetscFunctionReturn(0); 76*6718818eSStefano Zampini } 77*6718818eSStefano Zampini bt = atb->bt; 78*6718818eSStefano Zampini ct = atb->ct; 798949adfdSHong Zhang /* transpose local arry of B, then copy it to vector bt */ 801683a169SBarry Smith ierr = MatDenseGetArrayRead(B,&Barray);CHKERRQ(ierr); 818949adfdSHong Zhang ierr = VecGetArray(bt,&btarray);CHKERRQ(ierr); 828949adfdSHong Zhang 838949adfdSHong Zhang for (j=0; j<BN; j++) { 84905b3b74SHong Zhang for (i=0; i<m; i++) btarray[i*BN + j] = Barray[j*ldb + i]; 858949adfdSHong Zhang } 868949adfdSHong Zhang ierr = VecRestoreArray(bt,&btarray);CHKERRQ(ierr); 871683a169SBarry Smith ierr = MatDenseRestoreArrayRead(B,&Barray);CHKERRQ(ierr); 888949adfdSHong Zhang 898949adfdSHong Zhang /* compute ct = mA^T * cb */ 908949adfdSHong Zhang ierr = MatMultTranspose(atb->mA,bt,ct);CHKERRQ(ierr); 918949adfdSHong Zhang 92905b3b74SHong Zhang /* transpose local array of ct to matrix C */ 938949adfdSHong Zhang ierr = MatDenseGetArray(C,&Carray);CHKERRQ(ierr); 941683a169SBarry Smith ierr = VecGetArrayRead(ct,&ctarray);CHKERRQ(ierr); 95905b3b74SHong Zhang 968949adfdSHong Zhang for (j=0; j<BN; j++) { 97905b3b74SHong Zhang for (i=0; i<n; i++) Carray[j*ldc + i] = ctarray[i*BN + j]; 988949adfdSHong Zhang } 991683a169SBarry Smith ierr = VecRestoreArrayRead(ct,&ctarray);CHKERRQ(ierr); 1008949adfdSHong Zhang ierr = MatDenseRestoreArray(C,&Carray);CHKERRQ(ierr); 1018949adfdSHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1028949adfdSHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1038949adfdSHong Zhang PetscFunctionReturn(0); 1048949adfdSHong Zhang } 105