Actual source code: mattransposematmult.c
petsc-3.9.4 2018-09-11
2: /*
3: Defines matrix-matrix product routines
4: C = A^T * B
5: */
7: #include <../src/mat/impls/aij/seq/aij.h>
8: #include <../src/mat/impls/dense/seq/dense.h>
10: PetscErrorCode MatDestroy_SeqDense_MatTransMatMult(Mat A)
11: {
12: PetscErrorCode ierr;
13: Mat_SeqDense *a = (Mat_SeqDense*)A->data;
14: Mat_MatTransMatMult *atb = a->atb;
17: MatDestroy(&atb->mA);
18: VecDestroy(&atb->bt);
19: VecDestroy(&atb->ct);
20: (atb->destroy)(A);
21: PetscFree(atb);
22: return(0);
23: }
25: PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
26: {
28:
30: if (scall == MAT_INITIAL_MATRIX) {
31: PetscLogEventBegin(MAT_TransposeMatMultSymbolic,A,B,0,0);
32: MatTransposeMatMultSymbolic_SeqAIJ_SeqDense(A,B,fill,C);
33: PetscLogEventEnd(MAT_TransposeMatMultSymbolic,A,B,0,0);
34: }
35: PetscLogEventBegin(MAT_TransposeMatMultNumeric,A,B,0,0);
36: MatTransposeMatMultNumeric_SeqAIJ_SeqDense(A,B,*C);
37: PetscLogEventEnd(MAT_TransposeMatMultNumeric,A,B,0,0);
38: return(0);
39: }
41: PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
42: {
43: PetscErrorCode ierr;
44: PetscInt m=A->rmap->n,n=A->cmap->n,BN=B->cmap->N;
45: Mat_MatTransMatMult *atb;
46: Mat Cdense;
47: Vec bt,ct;
48: Mat_SeqDense *c;
51: PetscNew(&atb);
53: /* create output dense matrix C = A^T*B */
54: MatCreate(PETSC_COMM_SELF,&Cdense);
55: MatSetSizes(Cdense,n,BN,n,BN);
56: MatSetType(Cdense,MATSEQDENSE);
57: MatSeqDenseSetPreallocation(Cdense,NULL);
59: /* create vectors bt and ct to hold locally transposed arrays of B and C */
60: VecCreate(PETSC_COMM_SELF,&bt);
61: VecSetSizes(bt,m*BN,m*BN);
62: VecSetType(bt,VECSTANDARD);
63: VecCreate(PETSC_COMM_SELF,&ct);
64: VecSetSizes(ct,n*BN,n*BN);
65: VecSetType(ct,VECSTANDARD);
66: atb->bt = bt;
67: atb->ct = ct;
69: *C = Cdense;
70: c = (Mat_SeqDense*)Cdense->data;
71: c->atb = atb;
72: atb->destroy = Cdense->ops->destroy;
73: Cdense->ops->destroy = MatDestroy_SeqDense_MatTransMatMult;
74: return(0);
75: }
77: PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
78: {
79: PetscErrorCode ierr;
80: PetscInt i,j,k,m=A->rmap->n,n=A->cmap->n,BN=B->cmap->N;
81: PetscScalar *Barray,*Carray,*btarray,*ctarray;
82: Mat_SeqDense *c=(Mat_SeqDense*)C->data;
83: Mat_MatTransMatMult *atb=c->atb;
84: Vec bt=atb->bt,ct=atb->ct;
87: /* create MAIJ matrix mA from A -- should be done in symbolic phase */
88: MatDestroy(&atb->mA);
89: MatCreateMAIJ(A,BN,&atb->mA);
91: /* transpose local arry of B, then copy it to vector bt */
92: MatDenseGetArray(B,&Barray);
93: VecGetArray(bt,&btarray);
95: k=0;
96: for (j=0; j<BN; j++) {
97: for (i=0; i<m; i++) btarray[i*BN + j] = Barray[k++];
98: }
99: VecRestoreArray(bt,&btarray);
100: MatDenseRestoreArray(B,&Barray);
101:
102: /* compute ct = mA^T * cb */
103: MatMultTranspose(atb->mA,bt,ct);
105: /* transpose local arry of ct to matrix C */
106: MatDenseGetArray(C,&Carray);
107: VecGetArray(ct,&ctarray);
108: k = 0;
109: for (j=0; j<BN; j++) {
110: for (i=0; i<n; i++) Carray[k++] = ctarray[i*BN + j];
111: }
112: VecRestoreArray(ct,&ctarray);
113: MatDenseRestoreArray(C,&Carray);
114: return(0);
115: }