17a94429cSHong Zhang 27a94429cSHong Zhang /* 32cff0574SHong Zhang Defines matrix-matrix product routines 47a94429cSHong Zhang C = A^T * B 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 10afea5027SHong Zhang #undef __FUNCT__ 11afea5027SHong Zhang #define __FUNCT__ "MatDestroy_SeqDense_MatTransMatMult" 12afea5027SHong Zhang PetscErrorCode MatDestroy_SeqDense_MatTransMatMult(Mat A) 13afea5027SHong Zhang { 14afea5027SHong Zhang PetscErrorCode ierr; 15afea5027SHong Zhang Mat_SeqDense *a = (Mat_SeqDense*)A->data; 16afea5027SHong Zhang Mat_MatTransMatMult *atb = a->atb; 17afea5027SHong Zhang 18afea5027SHong Zhang PetscFunctionBegin; 19afea5027SHong Zhang ierr = MatDestroy(&atb->mA);CHKERRQ(ierr); 20afea5027SHong Zhang ierr = VecDestroy(&atb->bt);CHKERRQ(ierr); 21afea5027SHong Zhang ierr = VecDestroy(&atb->ct);CHKERRQ(ierr); 22afea5027SHong Zhang ierr = (atb->destroy)(A);CHKERRQ(ierr); 23afea5027SHong Zhang ierr = PetscFree(atb);CHKERRQ(ierr); 24afea5027SHong Zhang PetscFunctionReturn(0); 25afea5027SHong Zhang } 267a94429cSHong Zhang 277a94429cSHong Zhang #undef __FUNCT__ 287a94429cSHong Zhang #define __FUNCT__ "MatTransposeMatMult_SeqAIJ_SeqDense" 297a94429cSHong Zhang PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 307a94429cSHong Zhang { 317a94429cSHong Zhang PetscErrorCode ierr; 32*c608a817SHong Zhang 33*c608a817SHong Zhang PetscFunctionBegin; 34*c608a817SHong Zhang if (scall == MAT_INITIAL_MATRIX) { 35*c608a817SHong Zhang ierr = PetscLogEventBegin(MAT_TransposeMatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 36*c608a817SHong Zhang ierr = MatTransposeMatMultSymbolic_SeqAIJ_SeqDense(A,B,fill,C);CHKERRQ(ierr); 37*c608a817SHong Zhang ierr = PetscLogEventEnd(MAT_TransposeMatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 38*c608a817SHong Zhang } 39*c608a817SHong Zhang ierr = PetscLogEventBegin(MAT_TransposeMatMultNumeric,A,B,0,0);CHKERRQ(ierr); 40*c608a817SHong Zhang ierr = MatTransposeMatMultNumeric_SeqAIJ_SeqDense(A,B,*C);CHKERRQ(ierr); 41*c608a817SHong Zhang ierr = PetscLogEventEnd(MAT_TransposeMatMultNumeric,A,B,0,0);CHKERRQ(ierr); 42*c608a817SHong Zhang PetscFunctionReturn(0); 43*c608a817SHong Zhang } 44*c608a817SHong Zhang 45*c608a817SHong Zhang #undef __FUNCT__ 46*c608a817SHong Zhang #define __FUNCT__ "MatTransposeMatMultSymbolic_SeqAIJ_SeqDense" 47*c608a817SHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C) 48*c608a817SHong Zhang { 49*c608a817SHong Zhang PetscErrorCode ierr; 50*c608a817SHong Zhang PetscInt m,n,BN=B->cmap->N; 517a94429cSHong Zhang Mat_MatTransMatMult *atb; 527a94429cSHong Zhang Mat Cdense; 537a94429cSHong Zhang Vec bt,ct; 54afea5027SHong Zhang Mat_SeqDense *c; 557a94429cSHong Zhang 567a94429cSHong Zhang PetscFunctionBegin; 577a94429cSHong Zhang ierr = PetscNew(Mat_MatTransMatMult,&atb);CHKERRQ(ierr); 587a94429cSHong Zhang 597a94429cSHong Zhang /* create output dense matrix C = A^T*B */ 607a94429cSHong Zhang ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr); 617a94429cSHong Zhang ierr = MatCreate(PetscObjectComm((PetscObject)A),&Cdense);CHKERRQ(ierr); 62*c608a817SHong Zhang ierr = MatSetSizes(Cdense,n,BN,n,BN);CHKERRQ(ierr); 637a94429cSHong Zhang ierr = MatSetType(Cdense,MATDENSE);CHKERRQ(ierr); 647a94429cSHong Zhang ierr = MatSeqDenseSetPreallocation(Cdense,NULL);CHKERRQ(ierr); 657a94429cSHong Zhang 667a94429cSHong Zhang /* create vectors bt and ct to hold locally transposed arrays of B and C */ 677a94429cSHong Zhang ierr = VecCreate(PETSC_COMM_SELF,&bt);CHKERRQ(ierr); 68afea5027SHong Zhang ierr = VecSetSizes(bt,m*BN,m*BN);CHKERRQ(ierr); 697a94429cSHong Zhang ierr = VecSetFromOptions(bt);CHKERRQ(ierr); 707a94429cSHong Zhang ierr = VecCreate(PETSC_COMM_SELF,&ct);CHKERRQ(ierr); 71afea5027SHong Zhang ierr = VecSetSizes(ct,n*BN,n*BN);CHKERRQ(ierr); 727a94429cSHong Zhang ierr = VecSetFromOptions(ct);CHKERRQ(ierr); 737a94429cSHong Zhang atb->bt = bt; 747a94429cSHong Zhang atb->ct = ct; 757a94429cSHong Zhang 76*c608a817SHong Zhang *C = Cdense; 77*c608a817SHong Zhang c = (Mat_SeqDense*)Cdense->data; 78*c608a817SHong Zhang c->atb = atb; 79*c608a817SHong Zhang atb->destroy = Cdense->ops->destroy; 80*c608a817SHong Zhang Cdense->ops->destroy = MatDestroy_SeqDense_MatTransMatMult; 81*c608a817SHong Zhang PetscFunctionReturn(0); 82*c608a817SHong Zhang } 83*c608a817SHong Zhang 84*c608a817SHong Zhang #undef __FUNCT__ 85*c608a817SHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_SeqAIJ_SeqDense" 86*c608a817SHong Zhang PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C) 87*c608a817SHong Zhang { 88*c608a817SHong Zhang PetscErrorCode ierr; 89*c608a817SHong Zhang PetscInt i,j,k,m=A->rmap->n,n=A->cmap->n,BN=B->cmap->N; 90*c608a817SHong Zhang PetscScalar *Barray,*Carray,*btarray,*ctarray; 91*c608a817SHong Zhang Mat_SeqDense *c=(Mat_SeqDense*)C->data; 92*c608a817SHong Zhang Mat_MatTransMatMult *atb=c->atb; 93*c608a817SHong Zhang Vec bt=atb->bt,ct=atb->ct; 94*c608a817SHong Zhang 95*c608a817SHong Zhang PetscFunctionBegin; 96*c608a817SHong Zhang //ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr); 97*c608a817SHong Zhang 98*c608a817SHong Zhang /* create MAIJ matrix mA from A -- should be done in symbolic phase */ 99*c608a817SHong Zhang ierr = MatDestroy(&atb->mA);CHKERRQ(ierr); 100*c608a817SHong Zhang ierr = MatCreateMAIJ(A,BN,&atb->mA);CHKERRQ(ierr); 101*c608a817SHong Zhang 1027a94429cSHong Zhang /* transpose local arry of B, then copy it to vector bt */ 1037a94429cSHong Zhang ierr = MatDenseGetArray(B,&Barray);CHKERRQ(ierr); 1047a94429cSHong Zhang ierr = VecGetArray(bt,&btarray);CHKERRQ(ierr); 1057a94429cSHong Zhang 1067a94429cSHong Zhang k=0; 107afea5027SHong Zhang for (j=0; j<BN; j++) { 108afea5027SHong Zhang for (i=0; i<m; i++) btarray[i*BN + j] = Barray[k++]; 1097a94429cSHong Zhang } 1107a94429cSHong Zhang ierr = VecRestoreArray(bt,&btarray);CHKERRQ(ierr); 1117a94429cSHong Zhang ierr = MatDenseRestoreArray(B,&Barray);CHKERRQ(ierr); 1127a94429cSHong Zhang 1137a94429cSHong Zhang /* compute ct = mA^T * cb */ 1147a94429cSHong Zhang ierr = MatMultTranspose(atb->mA,bt,ct);CHKERRQ(ierr); 1157a94429cSHong Zhang 1167a94429cSHong Zhang /* transpose local arry of ct to matrix C */ 117*c608a817SHong Zhang ierr = MatDenseGetArray(C,&Carray);CHKERRQ(ierr); 1187a94429cSHong Zhang ierr = VecGetArray(ct,&ctarray);CHKERRQ(ierr); 1197a94429cSHong Zhang k = 0; 120afea5027SHong Zhang for (j=0; j<BN; j++) { 121afea5027SHong Zhang for (i=0; i<n; i++) Carray[k++] = ctarray[i*BN + j]; 1227a94429cSHong Zhang } 1237a94429cSHong Zhang ierr = VecRestoreArray(ct,&ctarray);CHKERRQ(ierr); 124*c608a817SHong Zhang ierr = MatDenseRestoreArray(C,&Carray);CHKERRQ(ierr); 1257a94429cSHong Zhang PetscFunctionReturn(0); 1267a94429cSHong Zhang } 127