xref: /petsc/src/mat/impls/aij/mpi/mpimattransposematmult.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
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   Cdense->ops->transposematmultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIDense;
76   PetscFunctionReturn(0);
77 }
78 
79 PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIDense(Mat A,Mat B,Mat C)
80 {
81   PetscErrorCode      ierr;
82   const PetscScalar   *Barray,*ctarray;
83   PetscScalar         *Carray,*btarray;
84   Mat_MPIDense        *b=(Mat_MPIDense*)B->data,*c=(Mat_MPIDense*)C->data;
85   Mat_SeqDense        *bseq=(Mat_SeqDense*)(b->A)->data,*cseq=(Mat_SeqDense*)(c->A)->data;
86   PetscInt            i,j,m=A->rmap->n,n=A->cmap->n,ldb=bseq->lda,BN=B->cmap->N,ldc=cseq->lda;
87   Mat_MatTransMatMult *atb=c->atb;
88   Vec                 bt=atb->bt,ct=atb->ct;
89 
90   PetscFunctionBegin;
91   /* create MAIJ matrix mA from A -- should be done in symbolic phase */
92   ierr = MatDestroy(&atb->mA);CHKERRQ(ierr);
93   ierr = MatCreateMAIJ(A,BN,&atb->mA);CHKERRQ(ierr);
94 
95   /* transpose local arry of B, then copy it to vector bt */
96   ierr = MatDenseGetArrayRead(B,&Barray);CHKERRQ(ierr);
97   ierr = VecGetArray(bt,&btarray);CHKERRQ(ierr);
98 
99   for (j=0; j<BN; j++) {
100     for (i=0; i<m; i++) btarray[i*BN + j] = Barray[j*ldb + i];
101   }
102   ierr = VecRestoreArray(bt,&btarray);CHKERRQ(ierr);
103   ierr = MatDenseRestoreArrayRead(B,&Barray);CHKERRQ(ierr);
104 
105   /* compute ct = mA^T * cb */
106   ierr = MatMultTranspose(atb->mA,bt,ct);CHKERRQ(ierr);
107 
108   /* transpose local array of ct to matrix C */
109   ierr = MatDenseGetArray(C,&Carray);CHKERRQ(ierr);
110   ierr = VecGetArrayRead(ct,&ctarray);CHKERRQ(ierr);
111 
112   for (j=0; j<BN; j++) {
113     for (i=0; i<n; i++) Carray[j*ldc + i] = ctarray[i*BN + j];
114   }
115   ierr = VecRestoreArrayRead(ct,&ctarray);CHKERRQ(ierr);
116   ierr = MatDenseRestoreArray(C,&Carray);CHKERRQ(ierr);
117   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
118   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
119   PetscFunctionReturn(0);
120 }
121