xref: /petsc/src/mat/impls/aij/seq/mattransposematmult.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
17a94429cSHong Zhang /*
26718818eSStefano Zampini   Defines matrix-matrix product routines for
36718818eSStefano Zampini           C = A^T * B and C = A * B^t
46718818eSStefano Zampini   with A SeqAIJ and B SeqDense
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 
MatProductCtxDestroy_SeqDense_MatTransMatMult(PetscCtxRt data)10*2a8381b2SBarry Smith static PetscErrorCode MatProductCtxDestroy_SeqDense_MatTransMatMult(PetscCtxRt data)
11d71ae5a4SJacob Faibussowitsch {
12cc1eb50dSBarry Smith   MatProductCtx_MatTransMatMult *atb = *(MatProductCtx_MatTransMatMult **)data;
13afea5027SHong Zhang 
14afea5027SHong Zhang   PetscFunctionBegin;
159566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&atb->mA));
169566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&atb->bt));
179566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&atb->ct));
189566063dSJacob Faibussowitsch   PetscCall(PetscFree(atb));
193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
20afea5027SHong Zhang }
217a94429cSHong Zhang 
226718818eSStefano Zampini static PetscErrorCode MatTMatTMultNumeric_SeqAIJ_SeqDense(Mat, Mat, Mat);
236718818eSStefano Zampini 
MatTMatTMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat C)24d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatTMatTMultSymbolic_SeqAIJ_SeqDense(Mat A, Mat B, PetscReal fill, Mat C)
25d71ae5a4SJacob Faibussowitsch {
26cc1eb50dSBarry Smith   MatProductCtx_MatTransMatMult *atb;
277a3c3d58SStefano Zampini   PetscBool                      cisdense;
286718818eSStefano Zampini   PetscInt                       dofm;
297a94429cSHong Zhang 
307a94429cSHong Zhang   PetscFunctionBegin;
316718818eSStefano Zampini   MatCheckProduct(C, 4);
3228b400f6SJacob Faibussowitsch   PetscCheck(!C->product->data, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Extra product struct not empty");
33aed4548fSBarry Smith   PetscCheck(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]);
347a94429cSHong Zhang 
356718818eSStefano Zampini   /* create output dense matrix C */
366718818eSStefano Zampini   if (C->product->type == MATPRODUCT_AtB) {
379566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(C, A->cmap->n, B->cmap->N, A->cmap->n, B->cmap->N));
386718818eSStefano Zampini     dofm = B->cmap->n;
396718818eSStefano Zampini   } else {
409566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(C, A->rmap->n, B->rmap->N, A->rmap->n, B->rmap->N));
416718818eSStefano Zampini     dofm = B->rmap->n;
426718818eSStefano Zampini   }
439566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompareAny((PetscObject)C, &cisdense, MATSEQDENSE, MATSEQDENSECUDA, ""));
4448a46eb9SPierre Jolivet   if (!cisdense) PetscCall(MatSetType(C, ((PetscObject)B)->type_name));
459566063dSJacob Faibussowitsch   PetscCall(MatSetUp(C));
467a94429cSHong Zhang 
476718818eSStefano Zampini   /* create additional data structure for the product */
489566063dSJacob Faibussowitsch   PetscCall(PetscNew(&atb));
499566063dSJacob Faibussowitsch   PetscCall(MatCreateMAIJ(A, dofm, &atb->mA));
509566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(atb->mA, &atb->ct, &atb->bt));
516718818eSStefano Zampini   C->product->data    = atb;
52cc1eb50dSBarry Smith   C->product->destroy = MatProductCtxDestroy_SeqDense_MatTransMatMult;
537a94429cSHong Zhang 
546718818eSStefano Zampini   if (C->product->type == MATPRODUCT_AtB) {
556718818eSStefano Zampini     C->ops->transposematmultnumeric = MatTMatTMultNumeric_SeqAIJ_SeqDense;
566718818eSStefano Zampini   } else {
576718818eSStefano Zampini     C->ops->mattransposemultnumeric = MatTMatTMultNumeric_SeqAIJ_SeqDense;
586718818eSStefano Zampini   }
593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
60c608a817SHong Zhang }
61c608a817SHong Zhang 
MatTMatTMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)6266976f2fSJacob Faibussowitsch static PetscErrorCode MatTMatTMultNumeric_SeqAIJ_SeqDense(Mat A, Mat B, Mat C)
63d71ae5a4SJacob Faibussowitsch {
646718818eSStefano Zampini   PetscInt                       i, j, m = A->rmap->n, n = A->cmap->n, blda, clda;
656718818eSStefano Zampini   PetscInt                       mdof = C->cmap->N;
666718818eSStefano Zampini   const PetscScalar             *Barray;
676718818eSStefano Zampini   PetscScalar                   *Carray;
68cc1eb50dSBarry Smith   MatProductCtx_MatTransMatMult *atb;
696718818eSStefano Zampini   Vec                            bt, ct;
70c608a817SHong Zhang 
71c608a817SHong Zhang   PetscFunctionBegin;
726718818eSStefano Zampini   MatCheckProduct(C, 3);
73aed4548fSBarry Smith   PetscCheck(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]);
74cc1eb50dSBarry Smith   atb = (MatProductCtx_MatTransMatMult *)C->product->data;
7528b400f6SJacob Faibussowitsch   PetscCheck(atb, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Missing product struct");
766718818eSStefano Zampini   bt = atb->bt;
776718818eSStefano Zampini   ct = atb->ct;
78c608a817SHong Zhang 
799566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayRead(B, &Barray));
809566063dSJacob Faibussowitsch   PetscCall(MatDenseGetLDA(B, &blda));
819566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayWrite(C, &Carray));
829566063dSJacob Faibussowitsch   PetscCall(MatDenseGetLDA(C, &clda));
836718818eSStefano Zampini   if (C->product->type == MATPRODUCT_AtB) { /* transpose local array of B, then copy it to vector bt */
846718818eSStefano Zampini     const PetscScalar *ctarray;
856718818eSStefano Zampini     PetscScalar       *btarray;
867a94429cSHong Zhang 
879566063dSJacob Faibussowitsch     PetscCall(VecGetArrayWrite(bt, &btarray));
886718818eSStefano Zampini     for (j = 0; j < mdof; j++) {
896718818eSStefano Zampini       for (i = 0; i < m; i++) btarray[i * mdof + j] = Barray[j * blda + i];
907a94429cSHong Zhang     }
919566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayWrite(bt, &btarray));
927a94429cSHong Zhang 
937a94429cSHong Zhang     /* compute ct = mA^T * cb */
949566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(atb->mA, bt, ct));
957a94429cSHong Zhang 
967a3c3d58SStefano Zampini     /* transpose local array of ct to matrix C */
979566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(ct, &ctarray));
986718818eSStefano Zampini     for (j = 0; j < mdof; j++) {
996718818eSStefano Zampini       for (i = 0; i < n; i++) Carray[j * clda + i] = ctarray[i * mdof + j];
1007a94429cSHong Zhang     }
1019566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(ct, &ctarray));
1026718818eSStefano Zampini   } else {
1036718818eSStefano Zampini     const PetscScalar *btarray;
1046718818eSStefano Zampini     PetscScalar       *ctarray;
1056718818eSStefano Zampini 
1066718818eSStefano Zampini     if (blda == B->rmap->n) {
1079566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(ct, Barray));
1086718818eSStefano Zampini     } else {
1096718818eSStefano Zampini       PetscInt bn = B->cmap->n;
1106718818eSStefano Zampini       PetscInt bm = B->rmap->n;
1116718818eSStefano Zampini 
1129566063dSJacob Faibussowitsch       PetscCall(VecGetArrayWrite(ct, &ctarray));
1136718818eSStefano Zampini       for (j = 0; j < bn; j++) {
1146718818eSStefano Zampini         for (i = 0; i < bm; i++) ctarray[j * bm + i] = Barray[j * blda + i];
1156718818eSStefano Zampini       }
1169566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayWrite(ct, &ctarray));
1176718818eSStefano Zampini     }
1186718818eSStefano Zampini 
1199566063dSJacob Faibussowitsch     PetscCall(MatMult(atb->mA, ct, bt));
12048a46eb9SPierre Jolivet     if (blda == B->rmap->n) PetscCall(VecResetArray(ct));
1219566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(bt, &btarray));
1226718818eSStefano Zampini     for (j = 0; j < mdof; j++) {
1236718818eSStefano Zampini       for (i = 0; i < m; i++) Carray[j * clda + i] = btarray[i * mdof + j];
1246718818eSStefano Zampini     }
1259566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(bt, &btarray));
1266718818eSStefano Zampini   }
1279566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayRead(B, &Barray));
1289566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArray(C, &Carray));
1293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1307a94429cSHong Zhang }
131