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 118949adfdSHong Zhang #undef __FUNCT__ 128949adfdSHong Zhang #define __FUNCT__ "MatDestroy_MPIDense_MatTransMatMult" 138949adfdSHong Zhang PetscErrorCode MatDestroy_MPIDense_MatTransMatMult(Mat A) 148949adfdSHong Zhang { 158949adfdSHong Zhang PetscErrorCode ierr; 168949adfdSHong Zhang Mat_MPIDense *a = (Mat_MPIDense*)A->data; 178949adfdSHong Zhang Mat_MatTransMatMult *atb = a->atb; 188949adfdSHong Zhang 198949adfdSHong Zhang PetscFunctionBegin; 208949adfdSHong Zhang ierr = MatDestroy(&atb->mA);CHKERRQ(ierr); 218949adfdSHong Zhang ierr = VecDestroy(&atb->bt);CHKERRQ(ierr); 228949adfdSHong Zhang ierr = VecDestroy(&atb->ct);CHKERRQ(ierr); 238949adfdSHong Zhang ierr = (atb->destroy)(A);CHKERRQ(ierr); 248949adfdSHong Zhang ierr = PetscFree(atb);CHKERRQ(ierr); 258949adfdSHong Zhang PetscFunctionReturn(0); 268949adfdSHong Zhang } 278949adfdSHong Zhang 288949adfdSHong Zhang #undef __FUNCT__ 298949adfdSHong Zhang #define __FUNCT__ "MatTransposeMatMult_MPIAIJ_MPIDense" 308949adfdSHong Zhang PetscErrorCode MatTransposeMatMult_MPIAIJ_MPIDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 318949adfdSHong Zhang { 328949adfdSHong Zhang PetscErrorCode ierr; 338949adfdSHong Zhang 348949adfdSHong Zhang PetscFunctionBegin; 358949adfdSHong Zhang if (scall == MAT_INITIAL_MATRIX) { 368949adfdSHong Zhang ierr = PetscLogEventBegin(MAT_TransposeMatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 378949adfdSHong Zhang ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIDense(A,B,fill,C);CHKERRQ(ierr); 388949adfdSHong Zhang ierr = PetscLogEventEnd(MAT_TransposeMatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 398949adfdSHong Zhang } 408949adfdSHong Zhang ierr = PetscLogEventBegin(MAT_TransposeMatMultNumeric,A,B,0,0);CHKERRQ(ierr); 418949adfdSHong Zhang ierr = MatTransposeMatMultNumeric_MPIAIJ_MPIDense(A,B,*C);CHKERRQ(ierr); 428949adfdSHong Zhang ierr = PetscLogEventEnd(MAT_TransposeMatMultNumeric,A,B,0,0);CHKERRQ(ierr); 438949adfdSHong Zhang PetscFunctionReturn(0); 448949adfdSHong Zhang } 458949adfdSHong Zhang 468949adfdSHong Zhang #undef __FUNCT__ 478949adfdSHong Zhang #define __FUNCT__ "MatTransposeMatMultSymbolic_MPIAIJ_MPIDense" 488949adfdSHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C) 498949adfdSHong Zhang { 508949adfdSHong Zhang PetscErrorCode ierr; 518949adfdSHong Zhang PetscInt m=A->rmap->n,n=A->cmap->n,BN=B->cmap->N; 528949adfdSHong Zhang Mat_MatTransMatMult *atb; 538949adfdSHong Zhang Mat Cdense; 548949adfdSHong Zhang Vec bt,ct; 558949adfdSHong Zhang Mat_MPIDense *c; 568949adfdSHong Zhang 578949adfdSHong Zhang PetscFunctionBegin; 58*b00a9115SJed Brown ierr = PetscNew(&atb);CHKERRQ(ierr); 598949adfdSHong Zhang 608949adfdSHong Zhang /* create output dense matrix C = A^T*B */ 618949adfdSHong Zhang ierr = MatCreate(PetscObjectComm((PetscObject)A),&Cdense);CHKERRQ(ierr); 628949adfdSHong Zhang ierr = MatSetSizes(Cdense,n,PETSC_DECIDE,PETSC_DECIDE,BN);CHKERRQ(ierr); 638949adfdSHong Zhang ierr = MatSetType(Cdense,MATMPIDENSE);CHKERRQ(ierr); 648949adfdSHong Zhang ierr = MatMPIDenseSetPreallocation(Cdense,NULL);CHKERRQ(ierr); 658949adfdSHong Zhang 668949adfdSHong Zhang /* create vectors bt and ct to hold locally transposed arrays of B and C */ 678949adfdSHong Zhang ierr = VecCreate(PetscObjectComm((PetscObject)A),&bt);CHKERRQ(ierr); 688949adfdSHong Zhang ierr = VecSetSizes(bt,m*BN,PETSC_DECIDE);CHKERRQ(ierr); 69c0dedaeaSBarry Smith ierr = VecSetType(bt,VECSTANDARD);CHKERRQ(ierr); 708949adfdSHong Zhang ierr = VecCreate(PetscObjectComm((PetscObject)A),&ct);CHKERRQ(ierr); 718949adfdSHong Zhang ierr = VecSetSizes(ct,n*BN,PETSC_DECIDE);CHKERRQ(ierr); 72c0dedaeaSBarry Smith ierr = VecSetType(ct,VECSTANDARD);CHKERRQ(ierr); 738949adfdSHong Zhang atb->bt = bt; 748949adfdSHong Zhang atb->ct = ct; 758949adfdSHong Zhang 768949adfdSHong Zhang *C = Cdense; 778949adfdSHong Zhang c = (Mat_MPIDense*)Cdense->data; 788949adfdSHong Zhang c->atb = atb; 798949adfdSHong Zhang atb->destroy = Cdense->ops->destroy; 808949adfdSHong Zhang Cdense->ops->destroy = MatDestroy_MPIDense_MatTransMatMult; 818949adfdSHong Zhang PetscFunctionReturn(0); 828949adfdSHong Zhang } 838949adfdSHong Zhang 848949adfdSHong Zhang #undef __FUNCT__ 858949adfdSHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_MPIAIJ_MPIDense" 868949adfdSHong Zhang PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIDense(Mat A,Mat B,Mat C) 878949adfdSHong Zhang { 888949adfdSHong Zhang PetscErrorCode ierr; 898949adfdSHong Zhang PetscInt i,j,k,m=A->rmap->n,n=A->cmap->n,BN=B->cmap->N; 908949adfdSHong Zhang PetscScalar *Barray,*Carray,*btarray,*ctarray; 918949adfdSHong Zhang Mat_MPIDense *c=(Mat_MPIDense*)C->data; 928949adfdSHong Zhang Mat_MatTransMatMult *atb=c->atb; 938949adfdSHong Zhang Vec bt=atb->bt,ct=atb->ct; 948949adfdSHong Zhang 958949adfdSHong Zhang PetscFunctionBegin; 968949adfdSHong Zhang /* create MAIJ matrix mA from A -- should be done in symbolic phase */ 978949adfdSHong Zhang ierr = MatDestroy(&atb->mA);CHKERRQ(ierr); 988949adfdSHong Zhang ierr = MatCreateMAIJ(A,BN,&atb->mA);CHKERRQ(ierr); 998949adfdSHong Zhang 1008949adfdSHong Zhang /* transpose local arry of B, then copy it to vector bt */ 1018949adfdSHong Zhang ierr = MatDenseGetArray(B,&Barray);CHKERRQ(ierr); 1028949adfdSHong Zhang ierr = VecGetArray(bt,&btarray);CHKERRQ(ierr); 1038949adfdSHong Zhang 1048949adfdSHong Zhang k=0; 1058949adfdSHong Zhang for (j=0; j<BN; j++) { 1068949adfdSHong Zhang for (i=0; i<m; i++) btarray[i*BN + j] = Barray[k++]; 1078949adfdSHong Zhang } 1088949adfdSHong Zhang ierr = VecRestoreArray(bt,&btarray);CHKERRQ(ierr); 1098949adfdSHong Zhang ierr = MatDenseRestoreArray(B,&Barray);CHKERRQ(ierr); 1108949adfdSHong Zhang 1118949adfdSHong Zhang /* compute ct = mA^T * cb */ 1128949adfdSHong Zhang ierr = MatMultTranspose(atb->mA,bt,ct);CHKERRQ(ierr); 1138949adfdSHong Zhang 1148949adfdSHong Zhang /* transpose local arry of ct to matrix C */ 1158949adfdSHong Zhang ierr = MatDenseGetArray(C,&Carray);CHKERRQ(ierr); 1168949adfdSHong Zhang ierr = VecGetArray(ct,&ctarray);CHKERRQ(ierr); 1178949adfdSHong Zhang k = 0; 1188949adfdSHong Zhang for (j=0; j<BN; j++) { 1198949adfdSHong Zhang for (i=0; i<n; i++) Carray[k++] = ctarray[i*BN + j]; 1208949adfdSHong Zhang } 1218949adfdSHong Zhang ierr = VecRestoreArray(ct,&ctarray);CHKERRQ(ierr); 1228949adfdSHong Zhang ierr = MatDenseRestoreArray(C,&Carray);CHKERRQ(ierr); 1238949adfdSHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1248949adfdSHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1258949adfdSHong Zhang PetscFunctionReturn(0); 1268949adfdSHong Zhang } 127