xref: /petsc/src/mat/impls/aij/seq/mattransposematmult.c (revision 3ba1676111f5c958fe6c2729b46ca4d523958bb3)
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 
11d71ae5a4SJacob Faibussowitsch PetscErrorCode MatDestroy_SeqDense_MatTransMatMult(void *data)
12d71ae5a4SJacob Faibussowitsch {
136718818eSStefano Zampini   Mat_MatTransMatMult *atb = (Mat_MatTransMatMult *)data;
14afea5027SHong Zhang 
15afea5027SHong Zhang   PetscFunctionBegin;
169566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&atb->mA));
179566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&atb->bt));
189566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&atb->ct));
199566063dSJacob Faibussowitsch   PetscCall(PetscFree(atb));
20*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
21afea5027SHong Zhang }
227a94429cSHong Zhang 
236718818eSStefano Zampini static PetscErrorCode MatTMatTMultNumeric_SeqAIJ_SeqDense(Mat, Mat, Mat);
246718818eSStefano Zampini 
25d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatTMatTMultSymbolic_SeqAIJ_SeqDense(Mat A, Mat B, PetscReal fill, Mat C)
26d71ae5a4SJacob Faibussowitsch {
277a94429cSHong Zhang   Mat_MatTransMatMult *atb;
287a3c3d58SStefano Zampini   PetscBool            cisdense;
296718818eSStefano Zampini   PetscInt             dofm;
307a94429cSHong Zhang 
317a94429cSHong Zhang   PetscFunctionBegin;
326718818eSStefano Zampini   MatCheckProduct(C, 4);
3328b400f6SJacob Faibussowitsch   PetscCheck(!C->product->data, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Extra product struct not empty");
34aed4548fSBarry 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]);
357a94429cSHong Zhang 
366718818eSStefano Zampini   /* create output dense matrix C */
376718818eSStefano Zampini   if (C->product->type == MATPRODUCT_AtB) {
389566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(C, A->cmap->n, B->cmap->N, A->cmap->n, B->cmap->N));
396718818eSStefano Zampini     dofm = B->cmap->n;
406718818eSStefano Zampini   } else {
419566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(C, A->rmap->n, B->rmap->N, A->rmap->n, B->rmap->N));
426718818eSStefano Zampini     dofm = B->rmap->n;
436718818eSStefano Zampini   }
449566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompareAny((PetscObject)C, &cisdense, MATSEQDENSE, MATSEQDENSECUDA, ""));
4548a46eb9SPierre Jolivet   if (!cisdense) PetscCall(MatSetType(C, ((PetscObject)B)->type_name));
469566063dSJacob Faibussowitsch   PetscCall(MatSetUp(C));
477a94429cSHong Zhang 
486718818eSStefano Zampini   /* create additional data structure for the product */
499566063dSJacob Faibussowitsch   PetscCall(PetscNew(&atb));
509566063dSJacob Faibussowitsch   PetscCall(MatCreateMAIJ(A, dofm, &atb->mA));
519566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(atb->mA, &atb->ct, &atb->bt));
526718818eSStefano Zampini   C->product->data    = atb;
536718818eSStefano Zampini   C->product->destroy = MatDestroy_SeqDense_MatTransMatMult;
547a94429cSHong Zhang 
556718818eSStefano Zampini   if (C->product->type == MATPRODUCT_AtB) {
566718818eSStefano Zampini     C->ops->transposematmultnumeric = MatTMatTMultNumeric_SeqAIJ_SeqDense;
576718818eSStefano Zampini   } else {
586718818eSStefano Zampini     C->ops->mattransposemultnumeric = MatTMatTMultNumeric_SeqAIJ_SeqDense;
596718818eSStefano Zampini   }
60*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
61c608a817SHong Zhang }
62c608a817SHong Zhang 
63d71ae5a4SJacob Faibussowitsch PetscErrorCode MatTMatTMultNumeric_SeqAIJ_SeqDense(Mat A, Mat B, Mat C)
64d71ae5a4SJacob Faibussowitsch {
656718818eSStefano Zampini   PetscInt             i, j, m = A->rmap->n, n = A->cmap->n, blda, clda;
666718818eSStefano Zampini   PetscInt             mdof = C->cmap->N;
676718818eSStefano Zampini   const PetscScalar   *Barray;
686718818eSStefano Zampini   PetscScalar         *Carray;
696718818eSStefano Zampini   Mat_MatTransMatMult *atb;
706718818eSStefano Zampini   Vec                  bt, ct;
71c608a817SHong Zhang 
72c608a817SHong Zhang   PetscFunctionBegin;
736718818eSStefano Zampini   MatCheckProduct(C, 3);
74aed4548fSBarry 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]);
756718818eSStefano Zampini   atb = (Mat_MatTransMatMult *)C->product->data;
7628b400f6SJacob Faibussowitsch   PetscCheck(atb, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Missing product struct");
776718818eSStefano Zampini   bt = atb->bt;
786718818eSStefano Zampini   ct = atb->ct;
79c608a817SHong Zhang 
809566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayRead(B, &Barray));
819566063dSJacob Faibussowitsch   PetscCall(MatDenseGetLDA(B, &blda));
829566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayWrite(C, &Carray));
839566063dSJacob Faibussowitsch   PetscCall(MatDenseGetLDA(C, &clda));
846718818eSStefano Zampini   if (C->product->type == MATPRODUCT_AtB) { /* transpose local array of B, then copy it to vector bt */
856718818eSStefano Zampini     const PetscScalar *ctarray;
866718818eSStefano Zampini     PetscScalar       *btarray;
877a94429cSHong Zhang 
889566063dSJacob Faibussowitsch     PetscCall(VecGetArrayWrite(bt, &btarray));
896718818eSStefano Zampini     for (j = 0; j < mdof; j++) {
906718818eSStefano Zampini       for (i = 0; i < m; i++) btarray[i * mdof + j] = Barray[j * blda + i];
917a94429cSHong Zhang     }
929566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayWrite(bt, &btarray));
937a94429cSHong Zhang 
947a94429cSHong Zhang     /* compute ct = mA^T * cb */
959566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(atb->mA, bt, ct));
967a94429cSHong Zhang 
977a3c3d58SStefano Zampini     /* transpose local array of ct to matrix C */
989566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(ct, &ctarray));
996718818eSStefano Zampini     for (j = 0; j < mdof; j++) {
1006718818eSStefano Zampini       for (i = 0; i < n; i++) Carray[j * clda + i] = ctarray[i * mdof + j];
1017a94429cSHong Zhang     }
1029566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(ct, &ctarray));
1036718818eSStefano Zampini   } else {
1046718818eSStefano Zampini     const PetscScalar *btarray;
1056718818eSStefano Zampini     PetscScalar       *ctarray;
1066718818eSStefano Zampini 
1076718818eSStefano Zampini     if (blda == B->rmap->n) {
1089566063dSJacob Faibussowitsch       PetscCall(VecPlaceArray(ct, Barray));
1096718818eSStefano Zampini     } else {
1106718818eSStefano Zampini       PetscInt bn = B->cmap->n;
1116718818eSStefano Zampini       PetscInt bm = B->rmap->n;
1126718818eSStefano Zampini 
1139566063dSJacob Faibussowitsch       PetscCall(VecGetArrayWrite(ct, &ctarray));
1146718818eSStefano Zampini       for (j = 0; j < bn; j++) {
1156718818eSStefano Zampini         for (i = 0; i < bm; i++) ctarray[j * bm + i] = Barray[j * blda + i];
1166718818eSStefano Zampini       }
1179566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayWrite(ct, &ctarray));
1186718818eSStefano Zampini     }
1196718818eSStefano Zampini 
1209566063dSJacob Faibussowitsch     PetscCall(MatMult(atb->mA, ct, bt));
12148a46eb9SPierre Jolivet     if (blda == B->rmap->n) PetscCall(VecResetArray(ct));
1229566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(bt, &btarray));
1236718818eSStefano Zampini     for (j = 0; j < mdof; j++) {
1246718818eSStefano Zampini       for (i = 0; i < m; i++) Carray[j * clda + i] = btarray[i * mdof + j];
1256718818eSStefano Zampini     }
1269566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(bt, &btarray));
1276718818eSStefano Zampini   }
1289566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayRead(B, &Barray));
1299566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArray(C, &Carray));
130*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1317a94429cSHong Zhang }
132