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