xref: /petsc/src/mat/impls/aij/seq/mattransposematmult.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1 
2 /*
3   Defines matrix-matrix product routines
4           C = A^T * B
5 */
6 
7 #include <../src/mat/impls/aij/seq/aij.h> /*I "petscmat.h" I*/
8 #include <../src/mat/impls/dense/seq/dense.h>
9 
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;
15 
16   PetscFunctionBegin;
17   ierr = MatDestroy(&atb->mA);CHKERRQ(ierr);
18   ierr = VecDestroy(&atb->bt);CHKERRQ(ierr);
19   ierr = VecDestroy(&atb->ct);CHKERRQ(ierr);
20   ierr = (atb->destroy)(A);CHKERRQ(ierr);
21   ierr = PetscFree(atb);CHKERRQ(ierr);
22   PetscFunctionReturn(0);
23 }
24 
25 PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
26 {
27   PetscErrorCode ierr;
28 
29   PetscFunctionBegin;
30   if (scall == MAT_INITIAL_MATRIX) {
31     ierr = PetscLogEventBegin(MAT_TransposeMatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
32     ierr = MatTransposeMatMultSymbolic_SeqAIJ_SeqDense(A,B,fill,C);CHKERRQ(ierr);
33     ierr = PetscLogEventEnd(MAT_TransposeMatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
34   }
35   ierr = PetscLogEventBegin(MAT_TransposeMatMultNumeric,A,B,0,0);CHKERRQ(ierr);
36   ierr = MatTransposeMatMultNumeric_SeqAIJ_SeqDense(A,B,*C);CHKERRQ(ierr);
37   ierr = PetscLogEventEnd(MAT_TransposeMatMultNumeric,A,B,0,0);CHKERRQ(ierr);
38   PetscFunctionReturn(0);
39 }
40 
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;
49 
50   PetscFunctionBegin;
51   ierr = PetscNew(&atb);CHKERRQ(ierr);
52 
53   /* create output dense matrix C = A^T*B */
54   ierr = MatCreate(PETSC_COMM_SELF,&Cdense);CHKERRQ(ierr);
55   ierr = MatSetSizes(Cdense,n,BN,n,BN);CHKERRQ(ierr);
56   ierr = MatSetType(Cdense,MATSEQDENSE);CHKERRQ(ierr);
57   ierr = MatSeqDenseSetPreallocation(Cdense,NULL);CHKERRQ(ierr);
58 
59   /* create vectors bt and ct to hold locally transposed arrays of B and C */
60   ierr = VecCreate(PETSC_COMM_SELF,&bt);CHKERRQ(ierr);
61   ierr = VecSetSizes(bt,m*BN,m*BN);CHKERRQ(ierr);
62   ierr = VecSetType(bt,VECSTANDARD);CHKERRQ(ierr);
63   ierr = VecCreate(PETSC_COMM_SELF,&ct);CHKERRQ(ierr);
64   ierr = VecSetSizes(ct,n*BN,n*BN);CHKERRQ(ierr);
65   ierr = VecSetType(ct,VECSTANDARD);CHKERRQ(ierr);
66   atb->bt = bt;
67   atb->ct = ct;
68 
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   Cdense->ops->transposematmultnumeric = MatTransposeMatMultNumeric_SeqAIJ_SeqDense;
75   PetscFunctionReturn(0);
76 }
77 
78 PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
79 {
80   PetscErrorCode      ierr;
81   PetscInt            i,j,k,m=A->rmap->n,n=A->cmap->n,BN=B->cmap->N;
82   const PetscScalar   *Barray,*ctarray;
83   PetscScalar         *Carray,*btarray;
84   Mat_SeqDense        *c=(Mat_SeqDense*)C->data;
85   Mat_MatTransMatMult *atb=c->atb;
86   Vec                 bt=atb->bt,ct=atb->ct;
87 
88   PetscFunctionBegin;
89   /* create MAIJ matrix mA from A -- should be done in symbolic phase */
90   ierr = MatDestroy(&atb->mA);CHKERRQ(ierr);
91   ierr = MatCreateMAIJ(A,BN,&atb->mA);CHKERRQ(ierr);
92 
93   /* transpose local arry of B, then copy it to vector bt */
94   ierr = MatDenseGetArrayRead(B,&Barray);CHKERRQ(ierr);
95   ierr = VecGetArray(bt,&btarray);CHKERRQ(ierr);
96 
97   k=0;
98   for (j=0; j<BN; j++) {
99     for (i=0; i<m; i++) btarray[i*BN + j] = Barray[k++];
100   }
101   ierr = VecRestoreArray(bt,&btarray);CHKERRQ(ierr);
102   ierr = MatDenseRestoreArrayRead(B,&Barray);CHKERRQ(ierr);
103 
104   /* compute ct = mA^T * cb */
105   ierr = MatMultTranspose(atb->mA,bt,ct);CHKERRQ(ierr);
106 
107   /* transpose local arry of ct to matrix C */
108   ierr = MatDenseGetArray(C,&Carray);CHKERRQ(ierr);
109   ierr = VecGetArrayRead(ct,&ctarray);CHKERRQ(ierr);
110   k = 0;
111   for (j=0; j<BN; j++) {
112     for (i=0; i<n; i++) Carray[k++] = ctarray[i*BN + j];
113   }
114   ierr = VecRestoreArrayRead(ct,&ctarray);CHKERRQ(ierr);
115   ierr = MatDenseRestoreArray(C,&Carray);CHKERRQ(ierr);
116   PetscFunctionReturn(0);
117 }
118