Actual source code: mpimattransposematmult.c
petsc-3.6.1 2015-08-06
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>
13: PetscErrorCode MatDestroy_MPIDense_MatTransMatMult(Mat A)
14: {
15: PetscErrorCode ierr;
16: Mat_MPIDense *a = (Mat_MPIDense*)A->data;
17: Mat_MatTransMatMult *atb = a->atb;
20: MatDestroy(&atb->mA);
21: VecDestroy(&atb->bt);
22: VecDestroy(&atb->ct);
23: (atb->destroy)(A);
24: PetscFree(atb);
25: return(0);
26: }
30: PetscErrorCode MatTransposeMatMult_MPIAIJ_MPIDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
31: {
35: if (scall == MAT_INITIAL_MATRIX) {
36: PetscLogEventBegin(MAT_TransposeMatMultSymbolic,A,B,0,0);
37: MatTransposeMatMultSymbolic_MPIAIJ_MPIDense(A,B,fill,C);
38: PetscLogEventEnd(MAT_TransposeMatMultSymbolic,A,B,0,0);
39: }
40: PetscLogEventBegin(MAT_TransposeMatMultNumeric,A,B,0,0);
41: MatTransposeMatMultNumeric_MPIAIJ_MPIDense(A,B,*C);
42: PetscLogEventEnd(MAT_TransposeMatMultNumeric,A,B,0,0);
43: return(0);
44: }
48: PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C)
49: {
50: PetscErrorCode ierr;
51: PetscInt m=A->rmap->n,n=A->cmap->n,BN=B->cmap->N;
52: Mat_MatTransMatMult *atb;
53: Mat Cdense;
54: Vec bt,ct;
55: Mat_MPIDense *c;
56:
58: PetscNew(&atb);
60: /* create output dense matrix C = A^T*B */
61: MatCreate(PetscObjectComm((PetscObject)A),&Cdense);
62: MatSetSizes(Cdense,n,PETSC_DECIDE,PETSC_DECIDE,BN);
63: MatSetType(Cdense,MATMPIDENSE);
64: MatMPIDenseSetPreallocation(Cdense,NULL);
66: /* create vectors bt and ct to hold locally transposed arrays of B and C */
67: VecCreate(PetscObjectComm((PetscObject)A),&bt);
68: VecSetSizes(bt,m*BN,PETSC_DECIDE);
69: VecSetType(bt,VECSTANDARD);
70: VecCreate(PetscObjectComm((PetscObject)A),&ct);
71: VecSetSizes(ct,n*BN,PETSC_DECIDE);
72: VecSetType(ct,VECSTANDARD);
73: atb->bt = bt;
74: atb->ct = ct;
76: *C = Cdense;
77: c = (Mat_MPIDense*)Cdense->data;
78: c->atb = atb;
79: atb->destroy = Cdense->ops->destroy;
80: Cdense->ops->destroy = MatDestroy_MPIDense_MatTransMatMult;
81: return(0);
82: }
86: PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIDense(Mat A,Mat B,Mat C)
87: {
88: PetscErrorCode ierr;
89: PetscInt i,j,k,m=A->rmap->n,n=A->cmap->n,BN=B->cmap->N;
90: PetscScalar *Barray,*Carray,*btarray,*ctarray;
91: Mat_MPIDense *c=(Mat_MPIDense*)C->data;
92: Mat_MatTransMatMult *atb=c->atb;
93: Vec bt=atb->bt,ct=atb->ct;
96: /* create MAIJ matrix mA from A -- should be done in symbolic phase */
97: MatDestroy(&atb->mA);
98: MatCreateMAIJ(A,BN,&atb->mA);
100: /* transpose local arry of B, then copy it to vector bt */
101: MatDenseGetArray(B,&Barray);
102: VecGetArray(bt,&btarray);
104: k=0;
105: for (j=0; j<BN; j++) {
106: for (i=0; i<m; i++) btarray[i*BN + j] = Barray[k++];
107: }
108: VecRestoreArray(bt,&btarray);
109: MatDenseRestoreArray(B,&Barray);
111: /* compute ct = mA^T * cb */
112: MatMultTranspose(atb->mA,bt,ct);
114: /* transpose local arry of ct to matrix C */
115: MatDenseGetArray(C,&Carray);
116: VecGetArray(ct,&ctarray);
117: k = 0;
118: for (j=0; j<BN; j++) {
119: for (i=0; i<n; i++) Carray[k++] = ctarray[i*BN + j];
120: }
121: VecRestoreArray(ct,&ctarray);
122: MatDenseRestoreArray(C,&Carray);
123: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
124: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
125: return(0);
126: }