Actual source code: mattransposematmult.c
petsc-3.6.4 2016-04-12
2: /*
3: Defines matrix-matrix product routines
4: C = A^T * B
5: */
7: #include <../src/mat/impls/aij/seq/aij.h> /*I "petscmat.h" I*/
8: #include <../src/mat/impls/dense/seq/dense.h>
12: PetscErrorCode MatDestroy_SeqDense_MatTransMatMult(Mat A)
13: {
14: PetscErrorCode ierr;
15: Mat_SeqDense *a = (Mat_SeqDense*)A->data;
16: Mat_MatTransMatMult *atb = a->atb;
19: MatDestroy(&atb->mA);
20: VecDestroy(&atb->bt);
21: VecDestroy(&atb->ct);
22: (atb->destroy)(A);
23: PetscFree(atb);
24: return(0);
25: }
29: PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
30: {
32:
34: if (scall == MAT_INITIAL_MATRIX) {
35: PetscLogEventBegin(MAT_TransposeMatMultSymbolic,A,B,0,0);
36: MatTransposeMatMultSymbolic_SeqAIJ_SeqDense(A,B,fill,C);
37: PetscLogEventEnd(MAT_TransposeMatMultSymbolic,A,B,0,0);
38: }
39: PetscLogEventBegin(MAT_TransposeMatMultNumeric,A,B,0,0);
40: MatTransposeMatMultNumeric_SeqAIJ_SeqDense(A,B,*C);
41: PetscLogEventEnd(MAT_TransposeMatMultNumeric,A,B,0,0);
42: return(0);
43: }
47: PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
48: {
49: PetscErrorCode ierr;
50: PetscInt m=A->rmap->n,n=A->cmap->n,BN=B->cmap->N;
51: Mat_MatTransMatMult *atb;
52: Mat Cdense;
53: Vec bt,ct;
54: Mat_SeqDense *c;
57: PetscNew(&atb);
59: /* create output dense matrix C = A^T*B */
60: MatCreate(PETSC_COMM_SELF,&Cdense);
61: MatSetSizes(Cdense,n,BN,n,BN);
62: MatSetType(Cdense,MATSEQDENSE);
63: MatSeqDenseSetPreallocation(Cdense,NULL);
65: /* create vectors bt and ct to hold locally transposed arrays of B and C */
66: VecCreate(PETSC_COMM_SELF,&bt);
67: VecSetSizes(bt,m*BN,m*BN);
68: VecSetType(bt,VECSTANDARD);
69: VecCreate(PETSC_COMM_SELF,&ct);
70: VecSetSizes(ct,n*BN,n*BN);
71: VecSetType(ct,VECSTANDARD);
72: atb->bt = bt;
73: atb->ct = ct;
75: *C = Cdense;
76: c = (Mat_SeqDense*)Cdense->data;
77: c->atb = atb;
78: atb->destroy = Cdense->ops->destroy;
79: Cdense->ops->destroy = MatDestroy_SeqDense_MatTransMatMult;
80: return(0);
81: }
85: PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
86: {
87: PetscErrorCode ierr;
88: PetscInt i,j,k,m=A->rmap->n,n=A->cmap->n,BN=B->cmap->N;
89: PetscScalar *Barray,*Carray,*btarray,*ctarray;
90: Mat_SeqDense *c=(Mat_SeqDense*)C->data;
91: Mat_MatTransMatMult *atb=c->atb;
92: Vec bt=atb->bt,ct=atb->ct;
95: /* create MAIJ matrix mA from A -- should be done in symbolic phase */
96: MatDestroy(&atb->mA);
97: MatCreateMAIJ(A,BN,&atb->mA);
99: /* transpose local arry of B, then copy it to vector bt */
100: MatDenseGetArray(B,&Barray);
101: VecGetArray(bt,&btarray);
103: k=0;
104: for (j=0; j<BN; j++) {
105: for (i=0; i<m; i++) btarray[i*BN + j] = Barray[k++];
106: }
107: VecRestoreArray(bt,&btarray);
108: MatDenseRestoreArray(B,&Barray);
109:
110: /* compute ct = mA^T * cb */
111: MatMultTranspose(atb->mA,bt,ct);
113: /* transpose local arry of ct to matrix C */
114: MatDenseGetArray(C,&Carray);
115: VecGetArray(ct,&ctarray);
116: k = 0;
117: for (j=0; j<BN; j++) {
118: for (i=0; i<n; i++) Carray[k++] = ctarray[i*BN + j];
119: }
120: VecRestoreArray(ct,&ctarray);
121: MatDenseRestoreArray(C,&Carray);
122: return(0);
123: }