xref: /petsc/src/mat/impls/aij/seq/mattransposematmult.c (revision 2c71b3e237ead271e4f3aa1505f92bf476e3413d)
17a94429cSHong Zhang 
27a94429cSHong Zhang /*
36718818eSStefano Zampini   Defines matrix-matrix product routines for
46718818eSStefano Zampini           C = A^T * B and C = A * B^t
56718818eSStefano Zampini   with A SeqAIJ and B SeqDense
67a94429cSHong Zhang */
77a94429cSHong Zhang 
87a94429cSHong Zhang #include <../src/mat/impls/aij/seq/aij.h> /*I "petscmat.h" I*/
9afea5027SHong Zhang #include <../src/mat/impls/dense/seq/dense.h>
10afea5027SHong Zhang 
116718818eSStefano Zampini PetscErrorCode MatDestroy_SeqDense_MatTransMatMult(void *data)
12afea5027SHong Zhang {
13afea5027SHong Zhang   PetscErrorCode      ierr;
146718818eSStefano Zampini   Mat_MatTransMatMult *atb = (Mat_MatTransMatMult *)data;
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 = PetscFree(atb);CHKERRQ(ierr);
21afea5027SHong Zhang   PetscFunctionReturn(0);
22afea5027SHong Zhang }
237a94429cSHong Zhang 
246718818eSStefano Zampini static PetscErrorCode MatTMatTMultNumeric_SeqAIJ_SeqDense(Mat,Mat,Mat);
256718818eSStefano Zampini 
266718818eSStefano Zampini PETSC_INTERN PetscErrorCode MatTMatTMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat C)
27c608a817SHong Zhang {
28c608a817SHong Zhang   PetscErrorCode      ierr;
297a94429cSHong Zhang   Mat_MatTransMatMult *atb;
307a3c3d58SStefano Zampini   PetscBool           cisdense;
316718818eSStefano Zampini   PetscInt            dofm;
327a94429cSHong Zhang 
337a94429cSHong Zhang   PetscFunctionBegin;
346718818eSStefano Zampini   MatCheckProduct(C,4);
35*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(C->product->data,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Extra product struct not empty");
36*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(C->product->type != MATPRODUCT_ABt && C->product->type != MATPRODUCT_AtB,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Not for product type %s",MatProductTypes[C->product->type]);
377a94429cSHong Zhang 
386718818eSStefano Zampini   /* create output dense matrix C */
396718818eSStefano Zampini   if (C->product->type == MATPRODUCT_AtB) {
406718818eSStefano Zampini     ierr = MatSetSizes(C,A->cmap->n,B->cmap->N,A->cmap->n,B->cmap->N);CHKERRQ(ierr);
416718818eSStefano Zampini     dofm = B->cmap->n;
426718818eSStefano Zampini   } else {
436718818eSStefano Zampini     ierr = MatSetSizes(C,A->rmap->n,B->rmap->N,A->rmap->n,B->rmap->N);CHKERRQ(ierr);
446718818eSStefano Zampini     dofm = B->rmap->n;
456718818eSStefano Zampini   }
467a3c3d58SStefano Zampini   ierr = PetscObjectTypeCompareAny((PetscObject)C,&cisdense,MATSEQDENSE,MATSEQDENSECUDA,"");CHKERRQ(ierr);
477a3c3d58SStefano Zampini   if (!cisdense) {
487a3c3d58SStefano Zampini     ierr = MatSetType(C,((PetscObject)B)->type_name);CHKERRQ(ierr);
497a3c3d58SStefano Zampini   }
5018992e5dSStefano Zampini   ierr = MatSetUp(C);CHKERRQ(ierr);
517a94429cSHong Zhang 
526718818eSStefano Zampini   /* create additional data structure for the product */
536718818eSStefano Zampini   ierr = PetscNew(&atb);CHKERRQ(ierr);
546718818eSStefano Zampini   ierr = MatCreateMAIJ(A,dofm,&atb->mA);CHKERRQ(ierr);
556718818eSStefano Zampini   ierr = MatCreateVecs(atb->mA,&atb->ct,&atb->bt);CHKERRQ(ierr);
566718818eSStefano Zampini   C->product->data    = atb;
576718818eSStefano Zampini   C->product->destroy = MatDestroy_SeqDense_MatTransMatMult;
587a94429cSHong Zhang 
596718818eSStefano Zampini   if (C->product->type == MATPRODUCT_AtB) {
606718818eSStefano Zampini     C->ops->transposematmultnumeric = MatTMatTMultNumeric_SeqAIJ_SeqDense;
616718818eSStefano Zampini   } else {
626718818eSStefano Zampini     C->ops->mattransposemultnumeric = MatTMatTMultNumeric_SeqAIJ_SeqDense;
636718818eSStefano Zampini   }
64c608a817SHong Zhang   PetscFunctionReturn(0);
65c608a817SHong Zhang }
66c608a817SHong Zhang 
676718818eSStefano Zampini PetscErrorCode MatTMatTMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
68c608a817SHong Zhang {
69c608a817SHong Zhang   PetscErrorCode      ierr;
706718818eSStefano Zampini   PetscInt            i,j,m=A->rmap->n,n=A->cmap->n,blda,clda;
716718818eSStefano Zampini   PetscInt            mdof = C->cmap->N;
726718818eSStefano Zampini   const PetscScalar   *Barray;
736718818eSStefano Zampini   PetscScalar         *Carray;
746718818eSStefano Zampini   Mat_MatTransMatMult *atb;
756718818eSStefano Zampini   Vec                 bt,ct;
76c608a817SHong Zhang 
77c608a817SHong Zhang   PetscFunctionBegin;
786718818eSStefano Zampini   MatCheckProduct(C,3);
79*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(C->product->type != MATPRODUCT_ABt && C->product->type != MATPRODUCT_AtB,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Not for product type %s",MatProductTypes[C->product->type]);
806718818eSStefano Zampini   atb = (Mat_MatTransMatMult *)C->product->data;
81*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!atb,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing product struct");
826718818eSStefano Zampini   bt = atb->bt;
836718818eSStefano Zampini   ct = atb->ct;
84c608a817SHong Zhang 
851683a169SBarry Smith   ierr = MatDenseGetArrayRead(B,&Barray);CHKERRQ(ierr);
866718818eSStefano Zampini   ierr = MatDenseGetLDA(B,&blda);CHKERRQ(ierr);
876718818eSStefano Zampini   ierr = MatDenseGetArrayWrite(C,&Carray);CHKERRQ(ierr);
886718818eSStefano Zampini   ierr = MatDenseGetLDA(C,&clda);CHKERRQ(ierr);
896718818eSStefano Zampini   if (C->product->type == MATPRODUCT_AtB) { /* transpose local array of B, then copy it to vector bt */
906718818eSStefano Zampini     const PetscScalar *ctarray;
916718818eSStefano Zampini     PetscScalar       *btarray;
927a94429cSHong Zhang 
936718818eSStefano Zampini     ierr = VecGetArrayWrite(bt,&btarray);CHKERRQ(ierr);
946718818eSStefano Zampini     for (j=0; j<mdof; j++) {
956718818eSStefano Zampini       for (i=0; i<m; i++) btarray[i*mdof + j] = Barray[j*blda + i];
967a94429cSHong Zhang     }
976718818eSStefano Zampini     ierr = VecRestoreArrayWrite(bt,&btarray);CHKERRQ(ierr);
987a94429cSHong Zhang 
997a94429cSHong Zhang     /* compute ct = mA^T * cb */
1007a94429cSHong Zhang     ierr = MatMultTranspose(atb->mA,bt,ct);CHKERRQ(ierr);
1017a94429cSHong Zhang 
1027a3c3d58SStefano Zampini     /* transpose local array of ct to matrix C */
1031683a169SBarry Smith     ierr = VecGetArrayRead(ct,&ctarray);CHKERRQ(ierr);
1046718818eSStefano Zampini     for (j=0; j<mdof; j++) {
1056718818eSStefano Zampini       for (i=0; i<n; i++) Carray[j*clda + i] = ctarray[i*mdof + j];
1067a94429cSHong Zhang     }
1071683a169SBarry Smith     ierr = VecRestoreArrayRead(ct,&ctarray);CHKERRQ(ierr);
1086718818eSStefano Zampini   } else {
1096718818eSStefano Zampini     const PetscScalar *btarray;
1106718818eSStefano Zampini     PetscScalar       *ctarray;
1116718818eSStefano Zampini 
1126718818eSStefano Zampini     if (blda == B->rmap->n) {
1136718818eSStefano Zampini       ierr = VecPlaceArray(ct,Barray);CHKERRQ(ierr);
1146718818eSStefano Zampini     } else {
1156718818eSStefano Zampini       PetscInt bn = B->cmap->n;
1166718818eSStefano Zampini       PetscInt bm = B->rmap->n;
1176718818eSStefano Zampini 
1186718818eSStefano Zampini       ierr = VecGetArrayWrite(ct,&ctarray);CHKERRQ(ierr);
1196718818eSStefano Zampini       for (j=0; j<bn; j++) {
1206718818eSStefano Zampini         for (i=0; i<bm; i++) ctarray[j*bm + i] = Barray[j*blda + i];
1216718818eSStefano Zampini       }
1226718818eSStefano Zampini       ierr = VecRestoreArrayWrite(ct,&ctarray);CHKERRQ(ierr);
1236718818eSStefano Zampini     }
1246718818eSStefano Zampini 
1256718818eSStefano Zampini     ierr = MatMult(atb->mA,ct,bt);CHKERRQ(ierr);
1266718818eSStefano Zampini     if (blda == B->rmap->n) {
1276718818eSStefano Zampini       ierr = VecResetArray(ct);CHKERRQ(ierr);
1286718818eSStefano Zampini     }
1296718818eSStefano Zampini     ierr = VecGetArrayRead(bt,&btarray);CHKERRQ(ierr);
1306718818eSStefano Zampini     for (j=0; j<mdof; j++) {
1316718818eSStefano Zampini       for (i=0; i<m; i++) Carray[j*clda + i] = btarray[i*mdof + j];
1326718818eSStefano Zampini     }
1336718818eSStefano Zampini     ierr = VecRestoreArrayRead(bt,&btarray);CHKERRQ(ierr);
1346718818eSStefano Zampini   }
1356718818eSStefano Zampini   ierr = MatDenseRestoreArrayRead(B,&Barray);CHKERRQ(ierr);
136c608a817SHong Zhang   ierr = MatDenseRestoreArray(C,&Carray);CHKERRQ(ierr);
1377a94429cSHong Zhang   PetscFunctionReturn(0);
1387a94429cSHong Zhang }
139