17a94429cSHong Zhang 27a94429cSHong Zhang /* 3*6718818eSStefano Zampini Defines matrix-matrix product routines for 4*6718818eSStefano Zampini C = A^T * B and C = A * B^t 5*6718818eSStefano 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 11*6718818eSStefano Zampini PetscErrorCode MatDestroy_SeqDense_MatTransMatMult(void *data) 12afea5027SHong Zhang { 13afea5027SHong Zhang PetscErrorCode ierr; 14*6718818eSStefano 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 24*6718818eSStefano Zampini static PetscErrorCode MatTMatTMultNumeric_SeqAIJ_SeqDense(Mat,Mat,Mat); 25*6718818eSStefano Zampini 26*6718818eSStefano 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; 31*6718818eSStefano Zampini PetscInt dofm; 327a94429cSHong Zhang 337a94429cSHong Zhang PetscFunctionBegin; 34*6718818eSStefano Zampini MatCheckProduct(C,4); 35*6718818eSStefano Zampini if (C->product->data) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Extra product struct not empty"); 36*6718818eSStefano Zampini if (C->product->type != MATPRODUCT_ABt && C->product->type != MATPRODUCT_AtB) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Not for product type %s",MatProductTypes[C->product->type]); 377a94429cSHong Zhang 38*6718818eSStefano Zampini /* create output dense matrix C */ 39*6718818eSStefano Zampini if (C->product->type == MATPRODUCT_AtB) { 40*6718818eSStefano Zampini ierr = MatSetSizes(C,A->cmap->n,B->cmap->N,A->cmap->n,B->cmap->N);CHKERRQ(ierr); 41*6718818eSStefano Zampini dofm = B->cmap->n; 42*6718818eSStefano Zampini } else { 43*6718818eSStefano Zampini ierr = MatSetSizes(C,A->rmap->n,B->rmap->N,A->rmap->n,B->rmap->N);CHKERRQ(ierr); 44*6718818eSStefano Zampini dofm = B->rmap->n; 45*6718818eSStefano 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 52*6718818eSStefano Zampini /* create additional data structure for the product */ 53*6718818eSStefano Zampini ierr = PetscNew(&atb);CHKERRQ(ierr); 54*6718818eSStefano Zampini ierr = MatCreateMAIJ(A,dofm,&atb->mA);CHKERRQ(ierr); 55*6718818eSStefano Zampini ierr = MatCreateVecs(atb->mA,&atb->ct,&atb->bt);CHKERRQ(ierr); 56*6718818eSStefano Zampini C->product->data = atb; 57*6718818eSStefano Zampini C->product->destroy = MatDestroy_SeqDense_MatTransMatMult; 587a94429cSHong Zhang 59*6718818eSStefano Zampini if (C->product->type == MATPRODUCT_AtB) { 60*6718818eSStefano Zampini C->ops->transposematmultnumeric = MatTMatTMultNumeric_SeqAIJ_SeqDense; 61*6718818eSStefano Zampini } else { 62*6718818eSStefano Zampini C->ops->mattransposemultnumeric = MatTMatTMultNumeric_SeqAIJ_SeqDense; 63*6718818eSStefano Zampini } 64c608a817SHong Zhang PetscFunctionReturn(0); 65c608a817SHong Zhang } 66c608a817SHong Zhang 67*6718818eSStefano Zampini PetscErrorCode MatTMatTMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C) 68c608a817SHong Zhang { 69c608a817SHong Zhang PetscErrorCode ierr; 70*6718818eSStefano Zampini PetscInt i,j,m=A->rmap->n,n=A->cmap->n,blda,clda; 71*6718818eSStefano Zampini PetscInt mdof = C->cmap->N; 72*6718818eSStefano Zampini const PetscScalar *Barray; 73*6718818eSStefano Zampini PetscScalar *Carray; 74*6718818eSStefano Zampini Mat_MatTransMatMult *atb; 75*6718818eSStefano Zampini Vec bt,ct; 76c608a817SHong Zhang 77c608a817SHong Zhang PetscFunctionBegin; 78*6718818eSStefano Zampini MatCheckProduct(C,3); 79*6718818eSStefano Zampini if (C->product->type != MATPRODUCT_ABt && C->product->type != MATPRODUCT_AtB) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Not for product type %s",MatProductTypes[C->product->type]); 80*6718818eSStefano Zampini atb = (Mat_MatTransMatMult *)C->product->data; 81*6718818eSStefano Zampini if (!atb) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing product struct"); 82*6718818eSStefano Zampini bt = atb->bt; 83*6718818eSStefano Zampini ct = atb->ct; 84c608a817SHong Zhang 851683a169SBarry Smith ierr = MatDenseGetArrayRead(B,&Barray);CHKERRQ(ierr); 86*6718818eSStefano Zampini ierr = MatDenseGetLDA(B,&blda);CHKERRQ(ierr); 87*6718818eSStefano Zampini ierr = MatDenseGetArrayWrite(C,&Carray);CHKERRQ(ierr); 88*6718818eSStefano Zampini ierr = MatDenseGetLDA(C,&clda);CHKERRQ(ierr); 89*6718818eSStefano Zampini if (C->product->type == MATPRODUCT_AtB) { /* transpose local array of B, then copy it to vector bt */ 90*6718818eSStefano Zampini const PetscScalar *ctarray; 91*6718818eSStefano Zampini PetscScalar *btarray; 927a94429cSHong Zhang 93*6718818eSStefano Zampini ierr = VecGetArrayWrite(bt,&btarray);CHKERRQ(ierr); 94*6718818eSStefano Zampini for (j=0; j<mdof; j++) { 95*6718818eSStefano Zampini for (i=0; i<m; i++) btarray[i*mdof + j] = Barray[j*blda + i]; 967a94429cSHong Zhang } 97*6718818eSStefano 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); 104*6718818eSStefano Zampini for (j=0; j<mdof; j++) { 105*6718818eSStefano 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); 108*6718818eSStefano Zampini } else { 109*6718818eSStefano Zampini const PetscScalar *btarray; 110*6718818eSStefano Zampini PetscScalar *ctarray; 111*6718818eSStefano Zampini 112*6718818eSStefano Zampini if (blda == B->rmap->n) { 113*6718818eSStefano Zampini ierr = VecPlaceArray(ct,Barray);CHKERRQ(ierr); 114*6718818eSStefano Zampini } else { 115*6718818eSStefano Zampini PetscInt bn = B->cmap->n; 116*6718818eSStefano Zampini PetscInt bm = B->rmap->n; 117*6718818eSStefano Zampini 118*6718818eSStefano Zampini ierr = VecGetArrayWrite(ct,&ctarray);CHKERRQ(ierr); 119*6718818eSStefano Zampini for (j=0; j<bn; j++) { 120*6718818eSStefano Zampini for (i=0; i<bm; i++) ctarray[j*bm + i] = Barray[j*blda + i]; 121*6718818eSStefano Zampini } 122*6718818eSStefano Zampini ierr = VecRestoreArrayWrite(ct,&ctarray);CHKERRQ(ierr); 123*6718818eSStefano Zampini } 124*6718818eSStefano Zampini 125*6718818eSStefano Zampini ierr = MatMult(atb->mA,ct,bt);CHKERRQ(ierr); 126*6718818eSStefano Zampini if (blda == B->rmap->n) { 127*6718818eSStefano Zampini ierr = VecResetArray(ct);CHKERRQ(ierr); 128*6718818eSStefano Zampini } 129*6718818eSStefano Zampini ierr = VecGetArrayRead(bt,&btarray);CHKERRQ(ierr); 130*6718818eSStefano Zampini for (j=0; j<mdof; j++) { 131*6718818eSStefano Zampini for (i=0; i<m; i++) Carray[j*clda + i] = btarray[i*mdof + j]; 132*6718818eSStefano Zampini } 133*6718818eSStefano Zampini ierr = VecRestoreArrayRead(bt,&btarray);CHKERRQ(ierr); 134*6718818eSStefano Zampini } 135*6718818eSStefano Zampini ierr = MatDenseRestoreArrayRead(B,&Barray);CHKERRQ(ierr); 136c608a817SHong Zhang ierr = MatDenseRestoreArray(C,&Carray);CHKERRQ(ierr); 1377a94429cSHong Zhang PetscFunctionReturn(0); 1387a94429cSHong Zhang } 139