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