1 2 /* 3 Defines matrix-matrix product routines for pairs of MPIAIJ matrices 4 C = A^T * B 5 The routines are slightly modified from MatTransposeMatMultxxx_SeqAIJ_SeqDense(). 6 */ 7 #include <../src/mat/impls/aij/seq/aij.h> /*I "petscmat.h" I*/ 8 #include <../src/mat/impls/aij/mpi/mpiaij.h> 9 #include <../src/mat/impls/dense/mpi/mpidense.h> 10 11 PetscErrorCode MatDestroy_MPIDense_MatTransMatMult(Mat A) 12 { 13 PetscErrorCode ierr; 14 Mat_MPIDense *a = (Mat_MPIDense*)A->data; 15 Mat_MatTransMatMult *atb = a->atb; 16 17 PetscFunctionBegin; 18 ierr = MatDestroy(&atb->mA);CHKERRQ(ierr); 19 ierr = VecDestroy(&atb->bt);CHKERRQ(ierr); 20 ierr = VecDestroy(&atb->ct);CHKERRQ(ierr); 21 ierr = (atb->destroy)(A);CHKERRQ(ierr); 22 ierr = PetscFree(atb);CHKERRQ(ierr); 23 PetscFunctionReturn(0); 24 } 25 26 PetscErrorCode MatTransposeMatMult_MPIAIJ_MPIDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 27 { 28 PetscErrorCode ierr; 29 30 PetscFunctionBegin; 31 if (scall == MAT_INITIAL_MATRIX) { 32 ierr = PetscLogEventBegin(MAT_TransposeMatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 33 ierr = MatTransposeMatMultSymbolic_MPIAIJ_MPIDense(A,B,fill,C);CHKERRQ(ierr); 34 ierr = PetscLogEventEnd(MAT_TransposeMatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 35 } 36 ierr = PetscLogEventBegin(MAT_TransposeMatMultNumeric,A,B,0,0);CHKERRQ(ierr); 37 ierr = MatTransposeMatMultNumeric_MPIAIJ_MPIDense(A,B,*C);CHKERRQ(ierr); 38 ierr = PetscLogEventEnd(MAT_TransposeMatMultNumeric,A,B,0,0);CHKERRQ(ierr); 39 PetscFunctionReturn(0); 40 } 41 42 PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C) 43 { 44 PetscErrorCode ierr; 45 PetscInt m=A->rmap->n,n=A->cmap->n,BN=B->cmap->N; 46 Mat_MatTransMatMult *atb; 47 Mat Cdense; 48 Vec bt,ct; 49 Mat_MPIDense *c; 50 51 PetscFunctionBegin; 52 ierr = PetscNew(&atb);CHKERRQ(ierr); 53 54 /* create output dense matrix C = A^T*B */ 55 ierr = MatCreate(PetscObjectComm((PetscObject)A),&Cdense);CHKERRQ(ierr); 56 ierr = MatSetSizes(Cdense,n,PETSC_DECIDE,PETSC_DECIDE,BN);CHKERRQ(ierr); 57 ierr = MatSetType(Cdense,MATMPIDENSE);CHKERRQ(ierr); 58 ierr = MatMPIDenseSetPreallocation(Cdense,NULL);CHKERRQ(ierr); 59 60 /* create vectors bt and ct to hold locally transposed arrays of B and C */ 61 ierr = VecCreate(PetscObjectComm((PetscObject)A),&bt);CHKERRQ(ierr); 62 ierr = VecSetSizes(bt,m*BN,PETSC_DECIDE);CHKERRQ(ierr); 63 ierr = VecSetType(bt,VECSTANDARD);CHKERRQ(ierr); 64 ierr = VecCreate(PetscObjectComm((PetscObject)A),&ct);CHKERRQ(ierr); 65 ierr = VecSetSizes(ct,n*BN,PETSC_DECIDE);CHKERRQ(ierr); 66 ierr = VecSetType(ct,VECSTANDARD);CHKERRQ(ierr); 67 atb->bt = bt; 68 atb->ct = ct; 69 70 *C = Cdense; 71 c = (Mat_MPIDense*)Cdense->data; 72 c->atb = atb; 73 atb->destroy = Cdense->ops->destroy; 74 Cdense->ops->destroy = MatDestroy_MPIDense_MatTransMatMult; 75 PetscFunctionReturn(0); 76 } 77 78 PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIDense(Mat A,Mat B,Mat C) 79 { 80 PetscErrorCode ierr; 81 PetscInt i,j,k,m=A->rmap->n,n=A->cmap->n,BN=B->cmap->N; 82 PetscScalar *Barray,*Carray,*btarray,*ctarray; 83 Mat_MPIDense *c=(Mat_MPIDense*)C->data; 84 Mat_MatTransMatMult *atb=c->atb; 85 Vec bt=atb->bt,ct=atb->ct; 86 87 PetscFunctionBegin; 88 /* create MAIJ matrix mA from A -- should be done in symbolic phase */ 89 ierr = MatDestroy(&atb->mA);CHKERRQ(ierr); 90 ierr = MatCreateMAIJ(A,BN,&atb->mA);CHKERRQ(ierr); 91 92 /* transpose local arry of B, then copy it to vector bt */ 93 ierr = MatDenseGetArray(B,&Barray);CHKERRQ(ierr); 94 ierr = VecGetArray(bt,&btarray);CHKERRQ(ierr); 95 96 k=0; 97 for (j=0; j<BN; j++) { 98 for (i=0; i<m; i++) btarray[i*BN + j] = Barray[k++]; 99 } 100 ierr = VecRestoreArray(bt,&btarray);CHKERRQ(ierr); 101 ierr = MatDenseRestoreArray(B,&Barray);CHKERRQ(ierr); 102 103 /* compute ct = mA^T * cb */ 104 ierr = MatMultTranspose(atb->mA,bt,ct);CHKERRQ(ierr); 105 106 /* transpose local arry of ct to matrix C */ 107 ierr = MatDenseGetArray(C,&Carray);CHKERRQ(ierr); 108 ierr = VecGetArray(ct,&ctarray);CHKERRQ(ierr); 109 k = 0; 110 for (j=0; j<BN; j++) { 111 for (i=0; i<n; i++) Carray[k++] = ctarray[i*BN + j]; 112 } 113 ierr = VecRestoreArray(ct,&ctarray);CHKERRQ(ierr); 114 ierr = MatDenseRestoreArray(C,&Carray);CHKERRQ(ierr); 115 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 116 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 117 PetscFunctionReturn(0); 118 } 119