xref: /petsc/src/mat/impls/aij/seq/mattransposematmult.c (revision 1683a169f336407cb9d0d8662b3b43f2a09d25dc)
17a94429cSHong Zhang 
27a94429cSHong Zhang /*
32cff0574SHong Zhang   Defines matrix-matrix product routines
47a94429cSHong Zhang           C = A^T * B
57a94429cSHong Zhang */
67a94429cSHong Zhang 
77a94429cSHong Zhang #include <../src/mat/impls/aij/seq/aij.h> /*I "petscmat.h" I*/
8afea5027SHong Zhang #include <../src/mat/impls/dense/seq/dense.h>
9afea5027SHong Zhang 
10afea5027SHong Zhang PetscErrorCode MatDestroy_SeqDense_MatTransMatMult(Mat A)
11afea5027SHong Zhang {
12afea5027SHong Zhang   PetscErrorCode      ierr;
13afea5027SHong Zhang   Mat_SeqDense        *a = (Mat_SeqDense*)A->data;
14afea5027SHong Zhang   Mat_MatTransMatMult *atb = a->atb;
15afea5027SHong Zhang 
16afea5027SHong Zhang   PetscFunctionBegin;
17afea5027SHong Zhang   ierr = MatDestroy(&atb->mA);CHKERRQ(ierr);
18afea5027SHong Zhang   ierr = VecDestroy(&atb->bt);CHKERRQ(ierr);
19afea5027SHong Zhang   ierr = VecDestroy(&atb->ct);CHKERRQ(ierr);
20afea5027SHong Zhang   ierr = (atb->destroy)(A);CHKERRQ(ierr);
21afea5027SHong Zhang   ierr = PetscFree(atb);CHKERRQ(ierr);
22afea5027SHong Zhang   PetscFunctionReturn(0);
23afea5027SHong Zhang }
247a94429cSHong Zhang 
257a94429cSHong Zhang PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
267a94429cSHong Zhang {
277a94429cSHong Zhang   PetscErrorCode ierr;
28c608a817SHong Zhang 
29c608a817SHong Zhang   PetscFunctionBegin;
30c608a817SHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
31c608a817SHong Zhang     ierr = PetscLogEventBegin(MAT_TransposeMatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
32c608a817SHong Zhang     ierr = MatTransposeMatMultSymbolic_SeqAIJ_SeqDense(A,B,fill,C);CHKERRQ(ierr);
33c608a817SHong Zhang     ierr = PetscLogEventEnd(MAT_TransposeMatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
34c608a817SHong Zhang   }
35c608a817SHong Zhang   ierr = PetscLogEventBegin(MAT_TransposeMatMultNumeric,A,B,0,0);CHKERRQ(ierr);
36c608a817SHong Zhang   ierr = MatTransposeMatMultNumeric_SeqAIJ_SeqDense(A,B,*C);CHKERRQ(ierr);
37c608a817SHong Zhang   ierr = PetscLogEventEnd(MAT_TransposeMatMultNumeric,A,B,0,0);CHKERRQ(ierr);
38c608a817SHong Zhang   PetscFunctionReturn(0);
39c608a817SHong Zhang }
40c608a817SHong Zhang 
41c608a817SHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
42c608a817SHong Zhang {
43c608a817SHong Zhang   PetscErrorCode      ierr;
448949adfdSHong Zhang   PetscInt            m=A->rmap->n,n=A->cmap->n,BN=B->cmap->N;
457a94429cSHong Zhang   Mat_MatTransMatMult *atb;
467a94429cSHong Zhang   Mat                 Cdense;
477a94429cSHong Zhang   Vec                 bt,ct;
48afea5027SHong Zhang   Mat_SeqDense        *c;
497a94429cSHong Zhang 
507a94429cSHong Zhang   PetscFunctionBegin;
51b00a9115SJed Brown   ierr = PetscNew(&atb);CHKERRQ(ierr);
527a94429cSHong Zhang 
537a94429cSHong Zhang   /* create output dense matrix C = A^T*B */
548949adfdSHong Zhang   ierr = MatCreate(PETSC_COMM_SELF,&Cdense);CHKERRQ(ierr);
55c608a817SHong Zhang   ierr = MatSetSizes(Cdense,n,BN,n,BN);CHKERRQ(ierr);
568949adfdSHong Zhang   ierr = MatSetType(Cdense,MATSEQDENSE);CHKERRQ(ierr);
577a94429cSHong Zhang   ierr = MatSeqDenseSetPreallocation(Cdense,NULL);CHKERRQ(ierr);
587a94429cSHong Zhang 
597a94429cSHong Zhang   /* create vectors bt and ct to hold locally transposed arrays of B and C */
607a94429cSHong Zhang   ierr = VecCreate(PETSC_COMM_SELF,&bt);CHKERRQ(ierr);
61afea5027SHong Zhang   ierr = VecSetSizes(bt,m*BN,m*BN);CHKERRQ(ierr);
62c0dedaeaSBarry Smith   ierr = VecSetType(bt,VECSTANDARD);CHKERRQ(ierr);
637a94429cSHong Zhang   ierr = VecCreate(PETSC_COMM_SELF,&ct);CHKERRQ(ierr);
64afea5027SHong Zhang   ierr = VecSetSizes(ct,n*BN,n*BN);CHKERRQ(ierr);
65c0dedaeaSBarry Smith   ierr = VecSetType(ct,VECSTANDARD);CHKERRQ(ierr);
667a94429cSHong Zhang   atb->bt = bt;
677a94429cSHong Zhang   atb->ct = ct;
687a94429cSHong Zhang 
69c608a817SHong Zhang   *C                   = Cdense;
70c608a817SHong Zhang   c                    = (Mat_SeqDense*)Cdense->data;
71c608a817SHong Zhang   c->atb               = atb;
72c608a817SHong Zhang   atb->destroy         = Cdense->ops->destroy;
73c608a817SHong Zhang   Cdense->ops->destroy = MatDestroy_SeqDense_MatTransMatMult;
74c608a817SHong Zhang   PetscFunctionReturn(0);
75c608a817SHong Zhang }
76c608a817SHong Zhang 
77c608a817SHong Zhang PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
78c608a817SHong Zhang {
79c608a817SHong Zhang   PetscErrorCode      ierr;
80c608a817SHong Zhang   PetscInt            i,j,k,m=A->rmap->n,n=A->cmap->n,BN=B->cmap->N;
81*1683a169SBarry Smith   const PetscScalar   *Barray,*ctarray;
82*1683a169SBarry Smith   PetscScalar         *Carray,*btarray;
83c608a817SHong Zhang   Mat_SeqDense        *c=(Mat_SeqDense*)C->data;
84c608a817SHong Zhang   Mat_MatTransMatMult *atb=c->atb;
85c608a817SHong Zhang   Vec                 bt=atb->bt,ct=atb->ct;
86c608a817SHong Zhang 
87c608a817SHong Zhang   PetscFunctionBegin;
88c608a817SHong Zhang   /* create MAIJ matrix mA from A -- should be done in symbolic phase */
89c608a817SHong Zhang   ierr = MatDestroy(&atb->mA);CHKERRQ(ierr);
90c608a817SHong Zhang   ierr = MatCreateMAIJ(A,BN,&atb->mA);CHKERRQ(ierr);
91c608a817SHong Zhang 
927a94429cSHong Zhang   /* transpose local arry of B, then copy it to vector bt */
93*1683a169SBarry Smith   ierr = MatDenseGetArrayRead(B,&Barray);CHKERRQ(ierr);
947a94429cSHong Zhang   ierr = VecGetArray(bt,&btarray);CHKERRQ(ierr);
957a94429cSHong Zhang 
967a94429cSHong Zhang   k=0;
97afea5027SHong Zhang   for (j=0; j<BN; j++) {
98afea5027SHong Zhang     for (i=0; i<m; i++) btarray[i*BN + j] = Barray[k++];
997a94429cSHong Zhang   }
1007a94429cSHong Zhang   ierr = VecRestoreArray(bt,&btarray);CHKERRQ(ierr);
101*1683a169SBarry Smith   ierr = MatDenseRestoreArrayRead(B,&Barray);CHKERRQ(ierr);
1027a94429cSHong Zhang 
1037a94429cSHong Zhang   /* compute ct = mA^T * cb */
1047a94429cSHong Zhang   ierr = MatMultTranspose(atb->mA,bt,ct);CHKERRQ(ierr);
1057a94429cSHong Zhang 
1067a94429cSHong Zhang   /* transpose local arry of ct to matrix C */
107c608a817SHong Zhang   ierr = MatDenseGetArray(C,&Carray);CHKERRQ(ierr);
108*1683a169SBarry Smith   ierr = VecGetArrayRead(ct,&ctarray);CHKERRQ(ierr);
1097a94429cSHong Zhang   k = 0;
110afea5027SHong Zhang   for (j=0; j<BN; j++) {
111afea5027SHong Zhang     for (i=0; i<n; i++) Carray[k++] = ctarray[i*BN + j];
1127a94429cSHong Zhang   }
113*1683a169SBarry Smith   ierr = VecRestoreArrayRead(ct,&ctarray);CHKERRQ(ierr);
114c608a817SHong Zhang   ierr = MatDenseRestoreArray(C,&Carray);CHKERRQ(ierr);
1157a94429cSHong Zhang   PetscFunctionReturn(0);
1167a94429cSHong Zhang }
117