xref: /petsc/src/mat/impls/aij/seq/mattransposematmult.c (revision 28b400f66ebc7ae0049166a2294dfcd3df27e64b)
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 {
136718818eSStefano Zampini   Mat_MatTransMatMult *atb = (Mat_MatTransMatMult *)data;
14afea5027SHong Zhang 
15afea5027SHong Zhang   PetscFunctionBegin;
165f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&atb->mA));
175f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&atb->bt));
185f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&atb->ct));
195f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(atb));
20afea5027SHong Zhang   PetscFunctionReturn(0);
21afea5027SHong Zhang }
227a94429cSHong Zhang 
236718818eSStefano Zampini static PetscErrorCode MatTMatTMultNumeric_SeqAIJ_SeqDense(Mat,Mat,Mat);
246718818eSStefano Zampini 
256718818eSStefano Zampini PETSC_INTERN PetscErrorCode MatTMatTMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat C)
26c608a817SHong Zhang {
277a94429cSHong Zhang   Mat_MatTransMatMult *atb;
287a3c3d58SStefano Zampini   PetscBool           cisdense;
296718818eSStefano Zampini   PetscInt            dofm;
307a94429cSHong Zhang 
317a94429cSHong Zhang   PetscFunctionBegin;
326718818eSStefano Zampini   MatCheckProduct(C,4);
33*28b400f6SJacob Faibussowitsch   PetscCheck(!C->product->data,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Extra product struct not empty");
342c71b3e2SJacob 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]);
357a94429cSHong Zhang 
366718818eSStefano Zampini   /* create output dense matrix C */
376718818eSStefano Zampini   if (C->product->type == MATPRODUCT_AtB) {
385f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetSizes(C,A->cmap->n,B->cmap->N,A->cmap->n,B->cmap->N));
396718818eSStefano Zampini     dofm = B->cmap->n;
406718818eSStefano Zampini   } else {
415f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetSizes(C,A->rmap->n,B->rmap->N,A->rmap->n,B->rmap->N));
426718818eSStefano Zampini     dofm = B->rmap->n;
436718818eSStefano Zampini   }
445f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectTypeCompareAny((PetscObject)C,&cisdense,MATSEQDENSE,MATSEQDENSECUDA,""));
457a3c3d58SStefano Zampini   if (!cisdense) {
465f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetType(C,((PetscObject)B)->type_name));
477a3c3d58SStefano Zampini   }
485f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(C));
497a94429cSHong Zhang 
506718818eSStefano Zampini   /* create additional data structure for the product */
515f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscNew(&atb));
525f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateMAIJ(A,dofm,&atb->mA));
535f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(atb->mA,&atb->ct,&atb->bt));
546718818eSStefano Zampini   C->product->data    = atb;
556718818eSStefano Zampini   C->product->destroy = MatDestroy_SeqDense_MatTransMatMult;
567a94429cSHong Zhang 
576718818eSStefano Zampini   if (C->product->type == MATPRODUCT_AtB) {
586718818eSStefano Zampini     C->ops->transposematmultnumeric = MatTMatTMultNumeric_SeqAIJ_SeqDense;
596718818eSStefano Zampini   } else {
606718818eSStefano Zampini     C->ops->mattransposemultnumeric = MatTMatTMultNumeric_SeqAIJ_SeqDense;
616718818eSStefano Zampini   }
62c608a817SHong Zhang   PetscFunctionReturn(0);
63c608a817SHong Zhang }
64c608a817SHong Zhang 
656718818eSStefano Zampini PetscErrorCode MatTMatTMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
66c608a817SHong Zhang {
676718818eSStefano Zampini   PetscInt            i,j,m=A->rmap->n,n=A->cmap->n,blda,clda;
686718818eSStefano Zampini   PetscInt            mdof = C->cmap->N;
696718818eSStefano Zampini   const PetscScalar   *Barray;
706718818eSStefano Zampini   PetscScalar         *Carray;
716718818eSStefano Zampini   Mat_MatTransMatMult *atb;
726718818eSStefano Zampini   Vec                 bt,ct;
73c608a817SHong Zhang 
74c608a817SHong Zhang   PetscFunctionBegin;
756718818eSStefano Zampini   MatCheckProduct(C,3);
762c71b3e2SJacob 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]);
776718818eSStefano Zampini   atb = (Mat_MatTransMatMult *)C->product->data;
78*28b400f6SJacob Faibussowitsch   PetscCheck(atb,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing product struct");
796718818eSStefano Zampini   bt = atb->bt;
806718818eSStefano Zampini   ct = atb->ct;
81c608a817SHong Zhang 
825f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDenseGetArrayRead(B,&Barray));
835f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDenseGetLDA(B,&blda));
845f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDenseGetArrayWrite(C,&Carray));
855f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDenseGetLDA(C,&clda));
866718818eSStefano Zampini   if (C->product->type == MATPRODUCT_AtB) { /* transpose local array of B, then copy it to vector bt */
876718818eSStefano Zampini     const PetscScalar *ctarray;
886718818eSStefano Zampini     PetscScalar       *btarray;
897a94429cSHong Zhang 
905f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetArrayWrite(bt,&btarray));
916718818eSStefano Zampini     for (j=0; j<mdof; j++) {
926718818eSStefano Zampini       for (i=0; i<m; i++) btarray[i*mdof + j] = Barray[j*blda + i];
937a94429cSHong Zhang     }
945f80ce2aSJacob Faibussowitsch     CHKERRQ(VecRestoreArrayWrite(bt,&btarray));
957a94429cSHong Zhang 
967a94429cSHong Zhang     /* compute ct = mA^T * cb */
975f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMultTranspose(atb->mA,bt,ct));
987a94429cSHong Zhang 
997a3c3d58SStefano Zampini     /* transpose local array of ct to matrix C */
1005f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetArrayRead(ct,&ctarray));
1016718818eSStefano Zampini     for (j=0; j<mdof; j++) {
1026718818eSStefano Zampini       for (i=0; i<n; i++) Carray[j*clda + i] = ctarray[i*mdof + j];
1037a94429cSHong Zhang     }
1045f80ce2aSJacob Faibussowitsch     CHKERRQ(VecRestoreArrayRead(ct,&ctarray));
1056718818eSStefano Zampini   } else {
1066718818eSStefano Zampini     const PetscScalar *btarray;
1076718818eSStefano Zampini     PetscScalar       *ctarray;
1086718818eSStefano Zampini 
1096718818eSStefano Zampini     if (blda == B->rmap->n) {
1105f80ce2aSJacob Faibussowitsch       CHKERRQ(VecPlaceArray(ct,Barray));
1116718818eSStefano Zampini     } else {
1126718818eSStefano Zampini       PetscInt bn = B->cmap->n;
1136718818eSStefano Zampini       PetscInt bm = B->rmap->n;
1146718818eSStefano Zampini 
1155f80ce2aSJacob Faibussowitsch       CHKERRQ(VecGetArrayWrite(ct,&ctarray));
1166718818eSStefano Zampini       for (j=0; j<bn; j++) {
1176718818eSStefano Zampini         for (i=0; i<bm; i++) ctarray[j*bm + i] = Barray[j*blda + i];
1186718818eSStefano Zampini       }
1195f80ce2aSJacob Faibussowitsch       CHKERRQ(VecRestoreArrayWrite(ct,&ctarray));
1206718818eSStefano Zampini     }
1216718818eSStefano Zampini 
1225f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMult(atb->mA,ct,bt));
1236718818eSStefano Zampini     if (blda == B->rmap->n) {
1245f80ce2aSJacob Faibussowitsch       CHKERRQ(VecResetArray(ct));
1256718818eSStefano Zampini     }
1265f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetArrayRead(bt,&btarray));
1276718818eSStefano Zampini     for (j=0; j<mdof; j++) {
1286718818eSStefano Zampini       for (i=0; i<m; i++) Carray[j*clda + i] = btarray[i*mdof + j];
1296718818eSStefano Zampini     }
1305f80ce2aSJacob Faibussowitsch     CHKERRQ(VecRestoreArrayRead(bt,&btarray));
1316718818eSStefano Zampini   }
1325f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDenseRestoreArrayRead(B,&Barray));
1335f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDenseRestoreArray(C,&Carray));
1347a94429cSHong Zhang   PetscFunctionReturn(0);
1357a94429cSHong Zhang }
136