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 116718818eSStefano Zampini PetscErrorCode MatDestroy_MPIDense_MatTransMatMult(void *data) 128949adfdSHong Zhang { 136718818eSStefano Zampini Mat_MatTransMatMult *atb = (Mat_MatTransMatMult*)data; 148949adfdSHong Zhang 158949adfdSHong Zhang PetscFunctionBegin; 165f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&atb->mA)); 175f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&atb->bt)); 185f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&atb->ct)); 195f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(atb)); 208949adfdSHong Zhang PetscFunctionReturn(0); 218949adfdSHong Zhang } 228949adfdSHong Zhang 236718818eSStefano Zampini static PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIDense(Mat,Mat,Mat); 246718818eSStefano Zampini 256718818eSStefano Zampini PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIDense(Mat A,Mat B,PetscReal fill,Mat C) 268949adfdSHong Zhang { 278949adfdSHong Zhang Mat_MatTransMatMult *atb; 286718818eSStefano Zampini PetscBool cisdense; 298949adfdSHong Zhang 308949adfdSHong Zhang PetscFunctionBegin; 316718818eSStefano Zampini MatCheckProduct(C,4); 32*28b400f6SJacob Faibussowitsch PetscCheck(!C->product->data,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Extra product struct not empty"); 338949adfdSHong Zhang 348949adfdSHong Zhang /* create output dense matrix C = A^T*B */ 355f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(C,A->cmap->n,B->cmap->n,A->cmap->N,B->cmap->N)); 365f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompareAny((PetscObject)C,&cisdense,MATMPIDENSE,MATMPIDENSECUDA,"")); 376718818eSStefano Zampini if (!cisdense) { 385f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(C,((PetscObject)B)->type_name)); 396718818eSStefano Zampini } 405f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(C)); 418949adfdSHong Zhang 426718818eSStefano Zampini /* create additional data structure for the product */ 435f80ce2aSJacob Faibussowitsch CHKERRQ(PetscNew(&atb)); 446718818eSStefano Zampini if (B->cmap->N) { 455f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateMAIJ(A,B->cmap->N,&atb->mA)); 46445ca090SPierre Jolivet if (!atb->mA->assembled) { 475f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(atb->mA,MAT_FINAL_ASSEMBLY)); 485f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(atb->mA,MAT_FINAL_ASSEMBLY)); 49445ca090SPierre Jolivet } 505f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(atb->mA,&atb->ct,&atb->bt)); 516718818eSStefano Zampini } 526718818eSStefano Zampini C->product->data = atb; 536718818eSStefano Zampini C->product->destroy = MatDestroy_MPIDense_MatTransMatMult; 548949adfdSHong Zhang 554222ddf1SHong Zhang C->ops->transposematmultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIDense; 568949adfdSHong Zhang PetscFunctionReturn(0); 578949adfdSHong Zhang } 588949adfdSHong Zhang 596718818eSStefano Zampini static PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIDense(Mat A,Mat B,Mat C) 608949adfdSHong Zhang { 611683a169SBarry Smith const PetscScalar *Barray,*ctarray; 621683a169SBarry Smith PetscScalar *Carray,*btarray; 63b45e3bf4SStefano Zampini PetscInt i,j,m=A->rmap->n,n=A->cmap->n,ldb,BN=B->cmap->N,ldc; 646718818eSStefano Zampini Mat_MatTransMatMult *atb; 656718818eSStefano Zampini Vec bt,ct; 668949adfdSHong Zhang 678949adfdSHong Zhang PetscFunctionBegin; 686718818eSStefano Zampini MatCheckProduct(C,3); 696718818eSStefano Zampini atb = (Mat_MatTransMatMult *)C->product->data; 702c71b3e2SJacob Faibussowitsch PetscCheckFalse(!atb,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing product struct"); 716718818eSStefano Zampini if (!BN) { 725f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 735f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 746718818eSStefano Zampini PetscFunctionReturn(0); 756718818eSStefano Zampini } 766718818eSStefano Zampini bt = atb->bt; 776718818eSStefano Zampini ct = atb->ct; 788949adfdSHong Zhang 79b45e3bf4SStefano Zampini /* transpose local array of B, then copy it to vector bt */ 805f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseGetArrayRead(B,&Barray)); 815f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseGetLDA(B,&ldb)); 825f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(bt,&btarray)); 83b45e3bf4SStefano Zampini for (j=0; j<BN; j++) 84b45e3bf4SStefano Zampini for (i=0; i<m; i++) 85b45e3bf4SStefano Zampini btarray[i*BN + j] = Barray[j*ldb + i]; 865f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(bt,&btarray)); 875f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseRestoreArrayRead(B,&Barray)); 888949adfdSHong Zhang 898949adfdSHong Zhang /* compute ct = mA^T * cb */ 905f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTranspose(atb->mA,bt,ct)); 918949adfdSHong Zhang 92905b3b74SHong Zhang /* transpose local array of ct to matrix C */ 935f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseGetArray(C,&Carray)); 945f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseGetLDA(C,&ldc)); 955f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(ct,&ctarray)); 96b45e3bf4SStefano Zampini for (j=0; j<BN; j++) 97b45e3bf4SStefano Zampini for (i=0; i<n; i++) 98b45e3bf4SStefano Zampini Carray[j*ldc + i] = ctarray[i*BN + j]; 995f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(ct,&ctarray)); 1005f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseRestoreArray(C,&Carray)); 1015f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 1025f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 1038949adfdSHong Zhang PetscFunctionReturn(0); 1048949adfdSHong Zhang } 105