Actual source code: mpimattransposematmult.c
petsc-3.13.6 2020-09-29
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>
8: #include <../src/mat/impls/aij/mpi/mpiaij.h>
9: #include <../src/mat/impls/dense/mpi/mpidense.h>
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;
18: MatDestroy(&atb->mA);
19: VecDestroy(&atb->bt);
20: VecDestroy(&atb->ct);
21: (atb->destroy)(A);
22: PetscFree(atb);
23: return(0);
24: }
26: PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIDense(Mat A,Mat B,PetscReal fill,Mat C)
27: {
28: PetscErrorCode ierr;
29: PetscInt m=A->rmap->n,n=A->cmap->n,BN=B->cmap->N;
30: Mat_MatTransMatMult *atb;
31: Vec bt,ct;
32: Mat_MPIDense *c;
35: PetscNew(&atb);
37: /* create output dense matrix C = A^T*B */
38: MatSetSizes(C,n,B->cmap->n,A->cmap->N,BN);
39: MatSetType(C,MATMPIDENSE);
40: MatSetUp(C);
42: /* create vectors bt and ct to hold locally transposed arrays of B and C */
43: VecCreate(PetscObjectComm((PetscObject)A),&bt);
44: VecSetSizes(bt,m*BN,PETSC_DECIDE);
45: VecSetType(bt,VECSTANDARD);
46: VecCreate(PetscObjectComm((PetscObject)A),&ct);
47: VecSetSizes(ct,n*BN,PETSC_DECIDE);
48: VecSetType(ct,VECSTANDARD);
49: atb->bt = bt;
50: atb->ct = ct;
52: c = (Mat_MPIDense*)C->data;
53: c->atb = atb;
54: atb->destroy = C->ops->destroy;
55: C->ops->destroy = MatDestroy_MPIDense_MatTransMatMult;
56: C->ops->transposematmultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIDense;
57: return(0);
58: }
60: PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIDense(Mat A,Mat B,Mat C)
61: {
62: PetscErrorCode ierr;
63: const PetscScalar *Barray,*ctarray;
64: PetscScalar *Carray,*btarray;
65: Mat_MPIDense *b=(Mat_MPIDense*)B->data,*c=(Mat_MPIDense*)C->data;
66: Mat_SeqDense *bseq=(Mat_SeqDense*)(b->A)->data,*cseq=(Mat_SeqDense*)(c->A)->data;
67: PetscInt i,j,m=A->rmap->n,n=A->cmap->n,ldb=bseq->lda,BN=B->cmap->N,ldc=cseq->lda;
68: Mat_MatTransMatMult *atb=c->atb;
69: Vec bt=atb->bt,ct=atb->ct;
72: /* create MAIJ matrix mA from A -- should be done in symbolic phase */
73: MatDestroy(&atb->mA);
74: MatCreateMAIJ(A,BN,&atb->mA);
76: /* transpose local arry of B, then copy it to vector bt */
77: MatDenseGetArrayRead(B,&Barray);
78: VecGetArray(bt,&btarray);
80: for (j=0; j<BN; j++) {
81: for (i=0; i<m; i++) btarray[i*BN + j] = Barray[j*ldb + i];
82: }
83: VecRestoreArray(bt,&btarray);
84: MatDenseRestoreArrayRead(B,&Barray);
86: /* compute ct = mA^T * cb */
87: MatMultTranspose(atb->mA,bt,ct);
89: /* transpose local array of ct to matrix C */
90: MatDenseGetArray(C,&Carray);
91: VecGetArrayRead(ct,&ctarray);
93: for (j=0; j<BN; j++) {
94: for (i=0; i<n; i++) Carray[j*ldc + i] = ctarray[i*BN + j];
95: }
96: VecRestoreArrayRead(ct,&ctarray);
97: MatDenseRestoreArray(C,&Carray);
98: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
99: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
100: return(0);
101: }