1 2 #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/ 3 4 typedef struct { 5 Mat A,U,V; 6 Vec work1,work2; /* Sequential vectors that hold partial products */ 7 PetscMPIInt nwork; /* length of work vectors */ 8 } Mat_LRC; 9 10 11 PETSC_INTERN PetscErrorCode MatMultTranspose_SeqDense(Mat,Vec,Vec); 12 PETSC_INTERN PetscErrorCode MatMultAdd_SeqDense(Mat,Vec,Vec,Vec); 13 14 #undef __FUNCT__ 15 #define __FUNCT__ "MatMult_LRC" 16 PetscErrorCode MatMult_LRC(Mat N,Vec x,Vec y) 17 { 18 Mat_LRC *Na = (Mat_LRC*)N->data; 19 PetscErrorCode ierr; 20 PetscScalar *w1,*w2; 21 22 PetscFunctionBegin; 23 ierr = MatMult(Na->A,x,y);CHKERRQ(ierr); 24 25 /* multiply the local part of V with the local part of x */ 26 /* note in this call x is treated as a sequential vector */ 27 ierr = MatMultTranspose_SeqDense(Na->V,x,Na->work1);CHKERRQ(ierr); 28 29 /* Form the sum of all the local multiplies : this is work2 = V'*x = 30 sum_{all processors} work1 */ 31 32 ierr = VecGetArray(Na->work1,&w1);CHKERRQ(ierr); 33 ierr = VecGetArray(Na->work2,&w2);CHKERRQ(ierr); 34 ierr = MPIU_Allreduce(w1,w2,Na->nwork,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)N));CHKERRQ(ierr); 35 ierr = VecRestoreArray(Na->work1,&w1);CHKERRQ(ierr); 36 ierr = VecRestoreArray(Na->work2,&w2);CHKERRQ(ierr); 37 38 /* multiply-sub y = y + U*work2 */ 39 /* note in this call y is treated as a sequential vector */ 40 ierr = MatMultAdd_SeqDense(Na->U,Na->work2,y,y);CHKERRQ(ierr); 41 PetscFunctionReturn(0); 42 } 43 44 #undef __FUNCT__ 45 #define __FUNCT__ "MatDestroy_LRC" 46 PetscErrorCode MatDestroy_LRC(Mat N) 47 { 48 Mat_LRC *Na = (Mat_LRC*)N->data; 49 PetscErrorCode ierr; 50 51 PetscFunctionBegin; 52 ierr = MatDestroy(&Na->A);CHKERRQ(ierr); 53 ierr = MatDestroy(&Na->U);CHKERRQ(ierr); 54 ierr = MatDestroy(&Na->V);CHKERRQ(ierr); 55 ierr = VecDestroy(&Na->work1);CHKERRQ(ierr); 56 ierr = VecDestroy(&Na->work2);CHKERRQ(ierr); 57 ierr = PetscFree(N->data);CHKERRQ(ierr); 58 PetscFunctionReturn(0); 59 } 60 61 62 #undef __FUNCT__ 63 #define __FUNCT__ "MatCreateLRC" 64 /*@ 65 MatCreateLRC - Creates a new matrix object that behaves like A + U*V' 66 67 Collective on Mat 68 69 Input Parameters: 70 + A - the (sparse) matrix 71 - U, V - two dense rectangular (tall and skinny) matrices 72 73 Output Parameter: 74 . N - the matrix that represents A + U*V' 75 76 Notes: 77 The matrix A + U*V' is not formed! Rather the new matrix 78 object performs the matrix-vector product by first multiplying by 79 A and then adding the other term. 80 81 Level: intermediate 82 @*/ 83 PetscErrorCode MatCreateLRC(Mat A,Mat U,Mat V,Mat *N) 84 { 85 PetscErrorCode ierr; 86 PetscInt m,n,k,m1,n1,k1; 87 Mat_LRC *Na; 88 89 PetscFunctionBegin; 90 PetscValidHeaderSpecific(A,MAT_CLASSID,1); 91 PetscValidHeaderSpecific(U,MAT_CLASSID,2); 92 PetscValidHeaderSpecific(V,MAT_CLASSID,3); 93 94 ierr = MatGetSize(U,NULL,&k);CHKERRQ(ierr); 95 ierr = MatGetSize(V,NULL,&k1);CHKERRQ(ierr); 96 if (k!=k1) SETERRQ(PetscObjectComm((PetscObject)U),PETSC_ERR_ARG_INCOMP,"U and V have different number of columns"); 97 ierr = MatGetLocalSize(U,&m,NULL);CHKERRQ(ierr); 98 ierr = MatGetLocalSize(V,&n,NULL);CHKERRQ(ierr); 99 ierr = MatGetLocalSize(A,&m1,&n1);CHKERRQ(ierr); 100 if (m!=m1) SETERRQ(PetscObjectComm((PetscObject)U),PETSC_ERR_ARG_INCOMP,"Local dimensions of U and A do not match"); 101 if (n!=n1) SETERRQ(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_INCOMP,"Local dimensions of V and A do not match"); 102 103 ierr = MatCreate(PetscObjectComm((PetscObject)A),N);CHKERRQ(ierr); 104 ierr = MatSetSizes(*N,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 105 ierr = PetscObjectChangeTypeName((PetscObject)*N,MATLRC);CHKERRQ(ierr); 106 107 ierr = PetscNewLog(*N,&Na);CHKERRQ(ierr); 108 (*N)->data = (void*)Na; 109 Na->A = A; 110 111 ierr = MatDenseGetLocalMatrix(U,&Na->U);CHKERRQ(ierr); 112 ierr = MatDenseGetLocalMatrix(V,&Na->V);CHKERRQ(ierr); 113 ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 114 ierr = PetscObjectReference((PetscObject)Na->U);CHKERRQ(ierr); 115 ierr = PetscObjectReference((PetscObject)Na->V);CHKERRQ(ierr); 116 117 ierr = VecCreateSeq(PETSC_COMM_SELF,U->cmap->N,&Na->work1);CHKERRQ(ierr); 118 ierr = VecDuplicate(Na->work1,&Na->work2);CHKERRQ(ierr); 119 Na->nwork = U->cmap->N; 120 121 (*N)->ops->destroy = MatDestroy_LRC; 122 (*N)->ops->mult = MatMult_LRC; 123 (*N)->assembled = PETSC_TRUE; 124 (*N)->cmap->N = V->rmap->N; 125 (*N)->rmap->N = U->rmap->N; 126 (*N)->cmap->n = V->rmap->n; 127 (*N)->rmap->n = U->rmap->n; 128 PetscFunctionReturn(0); 129 } 130 131