1 2 /* 3 Defines matrix-matrix product routines 4 C = A^T * B 5 */ 6 7 #include <../src/mat/impls/aij/seq/aij.h> /*I "petscmat.h" I*/ 8 #include <../src/mat/impls/dense/seq/dense.h> 9 10 PetscErrorCode MatDestroy_SeqDense_MatTransMatMult(Mat A) 11 { 12 PetscErrorCode ierr; 13 Mat_SeqDense *a = (Mat_SeqDense*)A->data; 14 Mat_MatTransMatMult *atb = a->atb; 15 16 PetscFunctionBegin; 17 ierr = MatDestroy(&atb->mA);CHKERRQ(ierr); 18 ierr = VecDestroy(&atb->bt);CHKERRQ(ierr); 19 ierr = VecDestroy(&atb->ct);CHKERRQ(ierr); 20 ierr = (atb->destroy)(A);CHKERRQ(ierr); 21 ierr = PetscFree(atb);CHKERRQ(ierr); 22 PetscFunctionReturn(0); 23 } 24 25 PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 26 { 27 PetscErrorCode ierr; 28 29 PetscFunctionBegin; 30 if (scall == MAT_INITIAL_MATRIX) { 31 ierr = PetscLogEventBegin(MAT_TransposeMatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 32 ierr = MatTransposeMatMultSymbolic_SeqAIJ_SeqDense(A,B,fill,C);CHKERRQ(ierr); 33 ierr = PetscLogEventEnd(MAT_TransposeMatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 34 } 35 ierr = PetscLogEventBegin(MAT_TransposeMatMultNumeric,A,B,0,0);CHKERRQ(ierr); 36 ierr = MatTransposeMatMultNumeric_SeqAIJ_SeqDense(A,B,*C);CHKERRQ(ierr); 37 ierr = PetscLogEventEnd(MAT_TransposeMatMultNumeric,A,B,0,0);CHKERRQ(ierr); 38 PetscFunctionReturn(0); 39 } 40 41 PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C) 42 { 43 PetscErrorCode ierr; 44 PetscInt m=A->rmap->n,n=A->cmap->n,BN=B->cmap->N; 45 Mat_MatTransMatMult *atb; 46 Mat Cdense; 47 Vec bt,ct; 48 Mat_SeqDense *c; 49 50 PetscFunctionBegin; 51 ierr = PetscNew(&atb);CHKERRQ(ierr); 52 53 /* create output dense matrix C = A^T*B */ 54 ierr = MatCreate(PETSC_COMM_SELF,&Cdense);CHKERRQ(ierr); 55 ierr = MatSetSizes(Cdense,n,BN,n,BN);CHKERRQ(ierr); 56 ierr = MatSetType(Cdense,MATSEQDENSE);CHKERRQ(ierr); 57 ierr = MatSeqDenseSetPreallocation(Cdense,NULL);CHKERRQ(ierr); 58 59 /* create vectors bt and ct to hold locally transposed arrays of B and C */ 60 ierr = VecCreate(PETSC_COMM_SELF,&bt);CHKERRQ(ierr); 61 ierr = VecSetSizes(bt,m*BN,m*BN);CHKERRQ(ierr); 62 ierr = VecSetType(bt,VECSTANDARD);CHKERRQ(ierr); 63 ierr = VecCreate(PETSC_COMM_SELF,&ct);CHKERRQ(ierr); 64 ierr = VecSetSizes(ct,n*BN,n*BN);CHKERRQ(ierr); 65 ierr = VecSetType(ct,VECSTANDARD);CHKERRQ(ierr); 66 atb->bt = bt; 67 atb->ct = ct; 68 69 *C = Cdense; 70 c = (Mat_SeqDense*)Cdense->data; 71 c->atb = atb; 72 atb->destroy = Cdense->ops->destroy; 73 Cdense->ops->destroy = MatDestroy_SeqDense_MatTransMatMult; 74 PetscFunctionReturn(0); 75 } 76 77 PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C) 78 { 79 PetscErrorCode ierr; 80 PetscInt i,j,k,m=A->rmap->n,n=A->cmap->n,BN=B->cmap->N; 81 const PetscScalar *Barray,*ctarray; 82 PetscScalar *Carray,*btarray; 83 Mat_SeqDense *c=(Mat_SeqDense*)C->data; 84 Mat_MatTransMatMult *atb=c->atb; 85 Vec bt=atb->bt,ct=atb->ct; 86 87 PetscFunctionBegin; 88 /* create MAIJ matrix mA from A -- should be done in symbolic phase */ 89 ierr = MatDestroy(&atb->mA);CHKERRQ(ierr); 90 ierr = MatCreateMAIJ(A,BN,&atb->mA);CHKERRQ(ierr); 91 92 /* transpose local arry of B, then copy it to vector bt */ 93 ierr = MatDenseGetArrayRead(B,&Barray);CHKERRQ(ierr); 94 ierr = VecGetArray(bt,&btarray);CHKERRQ(ierr); 95 96 k=0; 97 for (j=0; j<BN; j++) { 98 for (i=0; i<m; i++) btarray[i*BN + j] = Barray[k++]; 99 } 100 ierr = VecRestoreArray(bt,&btarray);CHKERRQ(ierr); 101 ierr = MatDenseRestoreArrayRead(B,&Barray);CHKERRQ(ierr); 102 103 /* compute ct = mA^T * cb */ 104 ierr = MatMultTranspose(atb->mA,bt,ct);CHKERRQ(ierr); 105 106 /* transpose local arry of ct to matrix C */ 107 ierr = MatDenseGetArray(C,&Carray);CHKERRQ(ierr); 108 ierr = VecGetArrayRead(ct,&ctarray);CHKERRQ(ierr); 109 k = 0; 110 for (j=0; j<BN; j++) { 111 for (i=0; i<n; i++) Carray[k++] = ctarray[i*BN + j]; 112 } 113 ierr = VecRestoreArrayRead(ct,&ctarray);CHKERRQ(ierr); 114 ierr = MatDenseRestoreArray(C,&Carray);CHKERRQ(ierr); 115 PetscFunctionReturn(0); 116 } 117