xref: /petsc/src/mat/impls/aij/seq/mattransposematmult.c (revision 18992e5d44eadf48aacf2808ebf9acfbbacf9355)
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 
254222ddf1SHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat C)
26c608a817SHong Zhang {
27c608a817SHong Zhang   PetscErrorCode      ierr;
288949adfdSHong Zhang   PetscInt            m=A->rmap->n,n=A->cmap->n,BN=B->cmap->N;
297a94429cSHong Zhang   Mat_MatTransMatMult *atb;
307a94429cSHong Zhang   Vec                 bt,ct;
31afea5027SHong Zhang   Mat_SeqDense        *c;
327a94429cSHong Zhang 
337a94429cSHong Zhang   PetscFunctionBegin;
34b00a9115SJed Brown   ierr = PetscNew(&atb);CHKERRQ(ierr);
357a94429cSHong Zhang 
367a94429cSHong Zhang   /* create output dense matrix C = A^T*B */
374222ddf1SHong Zhang   ierr = MatSetSizes(C,n,BN,n,BN);CHKERRQ(ierr);
384222ddf1SHong Zhang   ierr = MatSetType(C,MATSEQDENSE);CHKERRQ(ierr);
39*18992e5dSStefano Zampini   ierr = MatSetUp(C);CHKERRQ(ierr);
407a94429cSHong Zhang 
417a94429cSHong Zhang   /* create vectors bt and ct to hold locally transposed arrays of B and C */
427a94429cSHong Zhang   ierr = VecCreate(PETSC_COMM_SELF,&bt);CHKERRQ(ierr);
43afea5027SHong Zhang   ierr = VecSetSizes(bt,m*BN,m*BN);CHKERRQ(ierr);
44c0dedaeaSBarry Smith   ierr = VecSetType(bt,VECSTANDARD);CHKERRQ(ierr);
457a94429cSHong Zhang   ierr = VecCreate(PETSC_COMM_SELF,&ct);CHKERRQ(ierr);
46afea5027SHong Zhang   ierr = VecSetSizes(ct,n*BN,n*BN);CHKERRQ(ierr);
47c0dedaeaSBarry Smith   ierr = VecSetType(ct,VECSTANDARD);CHKERRQ(ierr);
487a94429cSHong Zhang   atb->bt = bt;
497a94429cSHong Zhang   atb->ct = ct;
507a94429cSHong Zhang 
514222ddf1SHong Zhang   c                               = (Mat_SeqDense*)C->data;
52c608a817SHong Zhang   c->atb                          = atb;
534222ddf1SHong Zhang   atb->destroy                    = C->ops->destroy;
544222ddf1SHong Zhang   C->ops->destroy                 = MatDestroy_SeqDense_MatTransMatMult;
554222ddf1SHong Zhang   C->ops->transposematmultnumeric = MatTransposeMatMultNumeric_SeqAIJ_SeqDense;
56c608a817SHong Zhang   PetscFunctionReturn(0);
57c608a817SHong Zhang }
58c608a817SHong Zhang 
59c608a817SHong Zhang PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
60c608a817SHong Zhang {
61c608a817SHong Zhang   PetscErrorCode      ierr;
62c608a817SHong Zhang   PetscInt            i,j,k,m=A->rmap->n,n=A->cmap->n,BN=B->cmap->N;
631683a169SBarry Smith   const PetscScalar   *Barray,*ctarray;
641683a169SBarry Smith   PetscScalar         *Carray,*btarray;
65c608a817SHong Zhang   Mat_SeqDense        *c=(Mat_SeqDense*)C->data;
66c608a817SHong Zhang   Mat_MatTransMatMult *atb=c->atb;
67c608a817SHong Zhang   Vec                 bt=atb->bt,ct=atb->ct;
68c608a817SHong Zhang 
69c608a817SHong Zhang   PetscFunctionBegin;
70c608a817SHong Zhang   /* create MAIJ matrix mA from A -- should be done in symbolic phase */
71c608a817SHong Zhang   ierr = MatDestroy(&atb->mA);CHKERRQ(ierr);
72c608a817SHong Zhang   ierr = MatCreateMAIJ(A,BN,&atb->mA);CHKERRQ(ierr);
73c608a817SHong Zhang 
747a94429cSHong Zhang   /* transpose local arry of B, then copy it to vector bt */
751683a169SBarry Smith   ierr = MatDenseGetArrayRead(B,&Barray);CHKERRQ(ierr);
767a94429cSHong Zhang   ierr = VecGetArray(bt,&btarray);CHKERRQ(ierr);
777a94429cSHong Zhang 
787a94429cSHong Zhang   k=0;
79afea5027SHong Zhang   for (j=0; j<BN; j++) {
80afea5027SHong Zhang     for (i=0; i<m; i++) btarray[i*BN + j] = Barray[k++];
817a94429cSHong Zhang   }
827a94429cSHong Zhang   ierr = VecRestoreArray(bt,&btarray);CHKERRQ(ierr);
831683a169SBarry Smith   ierr = MatDenseRestoreArrayRead(B,&Barray);CHKERRQ(ierr);
847a94429cSHong Zhang 
857a94429cSHong Zhang   /* compute ct = mA^T * cb */
867a94429cSHong Zhang   ierr = MatMultTranspose(atb->mA,bt,ct);CHKERRQ(ierr);
877a94429cSHong Zhang 
887a94429cSHong Zhang   /* transpose local arry of ct to matrix C */
89c608a817SHong Zhang   ierr = MatDenseGetArray(C,&Carray);CHKERRQ(ierr);
901683a169SBarry Smith   ierr = VecGetArrayRead(ct,&ctarray);CHKERRQ(ierr);
917a94429cSHong Zhang   k = 0;
92afea5027SHong Zhang   for (j=0; j<BN; j++) {
93afea5027SHong Zhang     for (i=0; i<n; i++) Carray[k++] = ctarray[i*BN + j];
947a94429cSHong Zhang   }
951683a169SBarry Smith   ierr = VecRestoreArrayRead(ct,&ctarray);CHKERRQ(ierr);
96c608a817SHong Zhang   ierr = MatDenseRestoreArray(C,&Carray);CHKERRQ(ierr);
977a94429cSHong Zhang   PetscFunctionReturn(0);
987a94429cSHong Zhang }
99