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 119371c9d4SSatish Balay PetscErrorCode MatDestroy_SeqDense_MatTransMatMult(void *data) { 126718818eSStefano Zampini Mat_MatTransMatMult *atb = (Mat_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)); 19afea5027SHong Zhang PetscFunctionReturn(0); 20afea5027SHong Zhang } 217a94429cSHong Zhang 226718818eSStefano Zampini static PetscErrorCode MatTMatTMultNumeric_SeqAIJ_SeqDense(Mat, Mat, Mat); 236718818eSStefano Zampini 249371c9d4SSatish Balay PETSC_INTERN PetscErrorCode MatTMatTMultSymbolic_SeqAIJ_SeqDense(Mat A, Mat B, PetscReal fill, Mat C) { 257a94429cSHong Zhang Mat_MatTransMatMult *atb; 267a3c3d58SStefano Zampini PetscBool cisdense; 276718818eSStefano Zampini PetscInt dofm; 287a94429cSHong Zhang 297a94429cSHong Zhang PetscFunctionBegin; 306718818eSStefano Zampini MatCheckProduct(C, 4); 3128b400f6SJacob Faibussowitsch PetscCheck(!C->product->data, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Extra product struct not empty"); 32aed4548fSBarry 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]); 337a94429cSHong Zhang 346718818eSStefano Zampini /* create output dense matrix C */ 356718818eSStefano Zampini if (C->product->type == MATPRODUCT_AtB) { 369566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C, A->cmap->n, B->cmap->N, A->cmap->n, B->cmap->N)); 376718818eSStefano Zampini dofm = B->cmap->n; 386718818eSStefano Zampini } else { 399566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C, A->rmap->n, B->rmap->N, A->rmap->n, B->rmap->N)); 406718818eSStefano Zampini dofm = B->rmap->n; 416718818eSStefano Zampini } 429566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompareAny((PetscObject)C, &cisdense, MATSEQDENSE, MATSEQDENSECUDA, "")); 43*48a46eb9SPierre Jolivet if (!cisdense) PetscCall(MatSetType(C, ((PetscObject)B)->type_name)); 449566063dSJacob Faibussowitsch PetscCall(MatSetUp(C)); 457a94429cSHong Zhang 466718818eSStefano Zampini /* create additional data structure for the product */ 479566063dSJacob Faibussowitsch PetscCall(PetscNew(&atb)); 489566063dSJacob Faibussowitsch PetscCall(MatCreateMAIJ(A, dofm, &atb->mA)); 499566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(atb->mA, &atb->ct, &atb->bt)); 506718818eSStefano Zampini C->product->data = atb; 516718818eSStefano Zampini C->product->destroy = MatDestroy_SeqDense_MatTransMatMult; 527a94429cSHong Zhang 536718818eSStefano Zampini if (C->product->type == MATPRODUCT_AtB) { 546718818eSStefano Zampini C->ops->transposematmultnumeric = MatTMatTMultNumeric_SeqAIJ_SeqDense; 556718818eSStefano Zampini } else { 566718818eSStefano Zampini C->ops->mattransposemultnumeric = MatTMatTMultNumeric_SeqAIJ_SeqDense; 576718818eSStefano Zampini } 58c608a817SHong Zhang PetscFunctionReturn(0); 59c608a817SHong Zhang } 60c608a817SHong Zhang 619371c9d4SSatish Balay PetscErrorCode MatTMatTMultNumeric_SeqAIJ_SeqDense(Mat A, Mat B, Mat C) { 626718818eSStefano Zampini PetscInt i, j, m = A->rmap->n, n = A->cmap->n, blda, clda; 636718818eSStefano Zampini PetscInt mdof = C->cmap->N; 646718818eSStefano Zampini const PetscScalar *Barray; 656718818eSStefano Zampini PetscScalar *Carray; 666718818eSStefano Zampini Mat_MatTransMatMult *atb; 676718818eSStefano Zampini Vec bt, ct; 68c608a817SHong Zhang 69c608a817SHong Zhang PetscFunctionBegin; 706718818eSStefano Zampini MatCheckProduct(C, 3); 71aed4548fSBarry 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]); 726718818eSStefano Zampini atb = (Mat_MatTransMatMult *)C->product->data; 7328b400f6SJacob Faibussowitsch PetscCheck(atb, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Missing product struct"); 746718818eSStefano Zampini bt = atb->bt; 756718818eSStefano Zampini ct = atb->ct; 76c608a817SHong Zhang 779566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayRead(B, &Barray)); 789566063dSJacob Faibussowitsch PetscCall(MatDenseGetLDA(B, &blda)); 799566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayWrite(C, &Carray)); 809566063dSJacob Faibussowitsch PetscCall(MatDenseGetLDA(C, &clda)); 816718818eSStefano Zampini if (C->product->type == MATPRODUCT_AtB) { /* transpose local array of B, then copy it to vector bt */ 826718818eSStefano Zampini const PetscScalar *ctarray; 836718818eSStefano Zampini PetscScalar *btarray; 847a94429cSHong Zhang 859566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(bt, &btarray)); 866718818eSStefano Zampini for (j = 0; j < mdof; j++) { 876718818eSStefano Zampini for (i = 0; i < m; i++) btarray[i * mdof + j] = Barray[j * blda + i]; 887a94429cSHong Zhang } 899566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(bt, &btarray)); 907a94429cSHong Zhang 917a94429cSHong Zhang /* compute ct = mA^T * cb */ 929566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(atb->mA, bt, ct)); 937a94429cSHong Zhang 947a3c3d58SStefano Zampini /* transpose local array of ct to matrix C */ 959566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(ct, &ctarray)); 966718818eSStefano Zampini for (j = 0; j < mdof; j++) { 976718818eSStefano Zampini for (i = 0; i < n; i++) Carray[j * clda + i] = ctarray[i * mdof + j]; 987a94429cSHong Zhang } 999566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(ct, &ctarray)); 1006718818eSStefano Zampini } else { 1016718818eSStefano Zampini const PetscScalar *btarray; 1026718818eSStefano Zampini PetscScalar *ctarray; 1036718818eSStefano Zampini 1046718818eSStefano Zampini if (blda == B->rmap->n) { 1059566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(ct, Barray)); 1066718818eSStefano Zampini } else { 1076718818eSStefano Zampini PetscInt bn = B->cmap->n; 1086718818eSStefano Zampini PetscInt bm = B->rmap->n; 1096718818eSStefano Zampini 1109566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(ct, &ctarray)); 1116718818eSStefano Zampini for (j = 0; j < bn; j++) { 1126718818eSStefano Zampini for (i = 0; i < bm; i++) ctarray[j * bm + i] = Barray[j * blda + i]; 1136718818eSStefano Zampini } 1149566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(ct, &ctarray)); 1156718818eSStefano Zampini } 1166718818eSStefano Zampini 1179566063dSJacob Faibussowitsch PetscCall(MatMult(atb->mA, ct, bt)); 118*48a46eb9SPierre Jolivet if (blda == B->rmap->n) PetscCall(VecResetArray(ct)); 1199566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(bt, &btarray)); 1206718818eSStefano Zampini for (j = 0; j < mdof; j++) { 1216718818eSStefano Zampini for (i = 0; i < m; i++) Carray[j * clda + i] = btarray[i * mdof + j]; 1226718818eSStefano Zampini } 1239566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(bt, &btarray)); 1246718818eSStefano Zampini } 1259566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayRead(B, &Barray)); 1269566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArray(C, &Carray)); 1277a94429cSHong Zhang PetscFunctionReturn(0); 1287a94429cSHong Zhang } 129