xref: /petsc/src/mat/impls/aij/seq/mattransposematmult.c (revision afea5027d9b82d962ae48d03fe31bc0ccbaef054)
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*/
8*afea5027SHong Zhang #include <../src/mat/impls/dense/seq/dense.h>
9*afea5027SHong Zhang 
10*afea5027SHong Zhang #undef __FUNCT__
11*afea5027SHong Zhang #define __FUNCT__ "MatDestroy_SeqDense_MatTransMatMult"
12*afea5027SHong Zhang PetscErrorCode MatDestroy_SeqDense_MatTransMatMult(Mat A)
13*afea5027SHong Zhang {
14*afea5027SHong Zhang   PetscErrorCode      ierr;
15*afea5027SHong Zhang   Mat_SeqDense        *a = (Mat_SeqDense*)A->data;
16*afea5027SHong Zhang   Mat_MatTransMatMult *atb = a->atb;
17*afea5027SHong Zhang 
18*afea5027SHong Zhang   PetscFunctionBegin;
19*afea5027SHong Zhang   ierr = MatDestroy(&atb->mA);CHKERRQ(ierr);
20*afea5027SHong Zhang   ierr = VecDestroy(&atb->bt);CHKERRQ(ierr);
21*afea5027SHong Zhang   ierr = VecDestroy(&atb->ct);CHKERRQ(ierr);
22*afea5027SHong Zhang   ierr = (atb->destroy)(A);CHKERRQ(ierr);
23*afea5027SHong Zhang   ierr = PetscFree(atb);CHKERRQ(ierr);
24*afea5027SHong Zhang   PetscFunctionReturn(0);
25*afea5027SHong Zhang }
267a94429cSHong Zhang 
277a94429cSHong Zhang #undef __FUNCT__
287a94429cSHong Zhang #define __FUNCT__ "MatTransposeMatMult_SeqAIJ_SeqDense"
297a94429cSHong Zhang PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
307a94429cSHong Zhang {
317a94429cSHong Zhang   PetscErrorCode      ierr;
32*afea5027SHong Zhang   PetscInt            i,j,k,m,n,BN=B->cmap->N;
337a94429cSHong Zhang   Mat_MatTransMatMult *atb;
347a94429cSHong Zhang   Mat                 Cdense;
357a94429cSHong Zhang   PetscScalar         *Barray,*Carray,*btarray,*ctarray;
367a94429cSHong Zhang   Vec                 bt,ct;
37*afea5027SHong Zhang   Mat_SeqDense        *c;
387a94429cSHong Zhang 
397a94429cSHong Zhang   PetscFunctionBegin;
407a94429cSHong Zhang   ierr = PetscNew(Mat_MatTransMatMult,&atb);CHKERRQ(ierr);
417a94429cSHong Zhang 
427a94429cSHong Zhang   /* create MAIJ matrix mA from A */
43*afea5027SHong Zhang   ierr = MatCreateMAIJ(A,BN,&atb->mA);CHKERRQ(ierr);
447a94429cSHong Zhang 
457a94429cSHong Zhang   /* create output dense matrix C = A^T*B */
467a94429cSHong Zhang   ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr);
477a94429cSHong Zhang   ierr = MatCreate(PetscObjectComm((PetscObject)A),&Cdense);CHKERRQ(ierr);
48*afea5027SHong Zhang   ierr = MatSetSizes(Cdense,n,BN,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
497a94429cSHong Zhang   ierr = MatSetType(Cdense,MATDENSE);CHKERRQ(ierr);
507a94429cSHong Zhang   ierr = MatSeqDenseSetPreallocation(Cdense,NULL);CHKERRQ(ierr);
517a94429cSHong Zhang 
527a94429cSHong Zhang   /* create vectors bt and ct to hold locally transposed arrays of B and C */
537a94429cSHong Zhang   ierr = VecCreate(PETSC_COMM_SELF,&bt);CHKERRQ(ierr);
54*afea5027SHong Zhang   ierr = VecSetSizes(bt,m*BN,m*BN);CHKERRQ(ierr);
557a94429cSHong Zhang   ierr = VecSetFromOptions(bt);CHKERRQ(ierr);
567a94429cSHong Zhang   ierr = VecCreate(PETSC_COMM_SELF,&ct);CHKERRQ(ierr);
57*afea5027SHong Zhang   ierr = VecSetSizes(ct,n*BN,n*BN);CHKERRQ(ierr);
587a94429cSHong Zhang   ierr = VecSetFromOptions(ct);CHKERRQ(ierr);
597a94429cSHong Zhang   atb->bt = bt;
607a94429cSHong Zhang   atb->ct = ct;
617a94429cSHong Zhang 
627a94429cSHong Zhang   /* transpose local arry of B, then copy it to vector bt */
637a94429cSHong Zhang   ierr = MatDenseGetArray(B,&Barray);CHKERRQ(ierr);
647a94429cSHong Zhang   ierr = VecGetArray(bt,&btarray);CHKERRQ(ierr);
657a94429cSHong Zhang 
667a94429cSHong Zhang   k=0;
67*afea5027SHong Zhang   for (j=0; j<BN; j++) {
68*afea5027SHong Zhang     for (i=0; i<m; i++) btarray[i*BN + j] = Barray[k++];
697a94429cSHong Zhang   }
707a94429cSHong Zhang   ierr = VecRestoreArray(bt,&btarray);CHKERRQ(ierr);
717a94429cSHong Zhang   ierr = MatDenseRestoreArray(B,&Barray);CHKERRQ(ierr);
727a94429cSHong Zhang 
737a94429cSHong Zhang   /* compute ct = mA^T * cb */
747a94429cSHong Zhang   ierr = MatMultTranspose(atb->mA,bt,ct);CHKERRQ(ierr);
757a94429cSHong Zhang 
767a94429cSHong Zhang   /* transpose local arry of ct to matrix C */
777a94429cSHong Zhang   ierr = MatDenseGetArray(Cdense,&Carray);CHKERRQ(ierr);
787a94429cSHong Zhang   ierr = VecGetArray(ct,&ctarray);CHKERRQ(ierr);
797a94429cSHong Zhang   k = 0;
80*afea5027SHong Zhang   for (j=0; j<BN; j++) {
81*afea5027SHong Zhang     for (i=0; i<n; i++) Carray[k++] = ctarray[i*BN + j];
827a94429cSHong Zhang   }
837a94429cSHong Zhang   ierr = VecRestoreArray(ct,&ctarray);CHKERRQ(ierr);
847a94429cSHong Zhang   ierr = MatDenseRestoreArray(Cdense,&Carray);CHKERRQ(ierr);
857a94429cSHong Zhang 
867a94429cSHong Zhang   *C                   = Cdense;
87*afea5027SHong Zhang   c                    = (Mat_SeqDense*)Cdense->data;
88*afea5027SHong Zhang   c->atb               = atb;
89*afea5027SHong Zhang   atb->destroy         = Cdense->ops->destroy;
90*afea5027SHong Zhang   Cdense->ops->destroy = MatDestroy_SeqDense_MatTransMatMult;
917a94429cSHong Zhang   PetscFunctionReturn(0);
927a94429cSHong Zhang }
937a94429cSHong Zhang 
947a94429cSHong Zhang #undef __FUNCT__
957a94429cSHong Zhang #define __FUNCT__ "MatTransposeMatMultSymbolic_SeqAIJ_SeqDense"
967a94429cSHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
977a94429cSHong Zhang {
987a94429cSHong Zhang   //PetscErrorCode ierr;
997a94429cSHong Zhang 
1007a94429cSHong Zhang   PetscFunctionBegin;
1017a94429cSHong Zhang   PetscFunctionReturn(0);
1027a94429cSHong Zhang }
1037a94429cSHong Zhang 
1047a94429cSHong Zhang #undef __FUNCT__
1057a94429cSHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_SeqAIJ_SeqDense"
1067a94429cSHong Zhang PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
1077a94429cSHong Zhang {
1087a94429cSHong Zhang     //PetscErrorCode ierr;
1097a94429cSHong Zhang 
1107a94429cSHong Zhang   PetscFunctionBegin;
1117a94429cSHong Zhang   PetscFunctionReturn(0);
1127a94429cSHong Zhang }
113