xref: /petsc/src/mat/impls/aij/seq/mattransposematmult.c (revision cdb0f33d09c128f365fdb48a6f07c56e211b6a43)
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 MatTransposeMatMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat C)
26 {
27   PetscErrorCode      ierr;
28   PetscInt            m=A->rmap->n,n=A->cmap->n,BN=B->cmap->N;
29   Mat_MatTransMatMult *atb;
30   Vec                 bt,ct;
31   Mat_SeqDense        *c;
32 
33   PetscFunctionBegin;
34   ierr = PetscNew(&atb);CHKERRQ(ierr);
35 
36   /* create output dense matrix C = A^T*B */
37   ierr = MatSetSizes(C,n,BN,n,BN);CHKERRQ(ierr);
38   ierr = MatSetType(C,MATSEQDENSE);CHKERRQ(ierr);
39   ierr = MatSetUp(C);CHKERRQ(ierr);
40 
41   /* create vectors bt and ct to hold locally transposed arrays of B and C */
42   ierr = VecCreate(PETSC_COMM_SELF,&bt);CHKERRQ(ierr);
43   ierr = VecSetSizes(bt,m*BN,m*BN);CHKERRQ(ierr);
44   ierr = VecSetType(bt,VECSTANDARD);CHKERRQ(ierr);
45   ierr = VecCreate(PETSC_COMM_SELF,&ct);CHKERRQ(ierr);
46   ierr = VecSetSizes(ct,n*BN,n*BN);CHKERRQ(ierr);
47   ierr = VecSetType(ct,VECSTANDARD);CHKERRQ(ierr);
48   atb->bt = bt;
49   atb->ct = ct;
50 
51   c                               = (Mat_SeqDense*)C->data;
52   c->atb                          = atb;
53   atb->destroy                    = C->ops->destroy;
54   C->ops->destroy                 = MatDestroy_SeqDense_MatTransMatMult;
55   C->ops->transposematmultnumeric = MatTransposeMatMultNumeric_SeqAIJ_SeqDense;
56   PetscFunctionReturn(0);
57 }
58 
59 PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
60 {
61   PetscErrorCode      ierr;
62   PetscInt            i,j,k,m=A->rmap->n,n=A->cmap->n,BN=B->cmap->N;
63   const PetscScalar   *Barray,*ctarray;
64   PetscScalar         *Carray,*btarray;
65   Mat_SeqDense        *c=(Mat_SeqDense*)C->data;
66   Mat_MatTransMatMult *atb=c->atb;
67   Vec                 bt=atb->bt,ct=atb->ct;
68 
69   PetscFunctionBegin;
70   /* create MAIJ matrix mA from A -- should be done in symbolic phase */
71   ierr = MatDestroy(&atb->mA);CHKERRQ(ierr);
72   ierr = MatCreateMAIJ(A,BN,&atb->mA);CHKERRQ(ierr);
73 
74   /* transpose local arry of B, then copy it to vector bt */
75   ierr = MatDenseGetArrayRead(B,&Barray);CHKERRQ(ierr);
76   ierr = VecGetArray(bt,&btarray);CHKERRQ(ierr);
77 
78   k=0;
79   for (j=0; j<BN; j++) {
80     for (i=0; i<m; i++) btarray[i*BN + j] = Barray[k++];
81   }
82   ierr = VecRestoreArray(bt,&btarray);CHKERRQ(ierr);
83   ierr = MatDenseRestoreArrayRead(B,&Barray);CHKERRQ(ierr);
84 
85   /* compute ct = mA^T * cb */
86   ierr = MatMultTranspose(atb->mA,bt,ct);CHKERRQ(ierr);
87 
88   /* transpose local arry of ct to matrix C */
89   ierr = MatDenseGetArray(C,&Carray);CHKERRQ(ierr);
90   ierr = VecGetArrayRead(ct,&ctarray);CHKERRQ(ierr);
91   k = 0;
92   for (j=0; j<BN; j++) {
93     for (i=0; i<n; i++) Carray[k++] = ctarray[i*BN + j];
94   }
95   ierr = VecRestoreArrayRead(ct,&ctarray);CHKERRQ(ierr);
96   ierr = MatDenseRestoreArray(C,&Carray);CHKERRQ(ierr);
97   PetscFunctionReturn(0);
98 }
99