1ceb9aaf7SBarry Smith 2af0996ceSBarry Smith #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/ 3ceb9aaf7SBarry Smith 4ceb9aaf7SBarry Smith typedef struct { 5ab92ecdeSBarry Smith Mat A,U,V; 69d9032efSJose E. Roman Vec work1,work2; /* Sequential vectors that hold partial products */ 7ab92ecdeSBarry Smith PetscMPIInt nwork; /* length of work vectors */ 8ab92ecdeSBarry Smith } Mat_LRC; 9ab92ecdeSBarry Smith 10ab92ecdeSBarry Smith 11cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatMultTranspose_SeqDense(Mat,Vec,Vec); 12cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatMultAdd_SeqDense(Mat,Vec,Vec,Vec); 1386eed97dSJose E. Roman PETSC_INTERN PetscErrorCode MatMult_SeqDense(Mat,Vec,Vec); 14ceb9aaf7SBarry Smith 15ceb9aaf7SBarry Smith #undef __FUNCT__ 16ab92ecdeSBarry Smith #define __FUNCT__ "MatMult_LRC" 17ab92ecdeSBarry Smith PetscErrorCode MatMult_LRC(Mat N,Vec x,Vec y) 18ceb9aaf7SBarry Smith { 19ab92ecdeSBarry Smith Mat_LRC *Na = (Mat_LRC*)N->data; 20ceb9aaf7SBarry Smith PetscErrorCode ierr; 21ab92ecdeSBarry Smith PetscScalar *w1,*w2; 22ceb9aaf7SBarry Smith 23ceb9aaf7SBarry Smith PetscFunctionBegin; 24ab92ecdeSBarry Smith /* multiply the local part of V with the local part of x */ 25ab92ecdeSBarry Smith /* note in this call x is treated as a sequential vector */ 26cc2e6a90SBarry Smith ierr = MatMultTranspose_SeqDense(Na->V,x,Na->work1);CHKERRQ(ierr); 27ab92ecdeSBarry Smith 2886eed97dSJose E. Roman /* form the sum of all the local multiplies: this is work2 = V'*x = 29ab92ecdeSBarry Smith sum_{all processors} work1 */ 30ab92ecdeSBarry Smith ierr = VecGetArray(Na->work1,&w1);CHKERRQ(ierr); 31ab92ecdeSBarry Smith ierr = VecGetArray(Na->work2,&w2);CHKERRQ(ierr); 32b2566f29SBarry Smith ierr = MPIU_Allreduce(w1,w2,Na->nwork,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)N));CHKERRQ(ierr); 33ab92ecdeSBarry Smith ierr = VecRestoreArray(Na->work1,&w1);CHKERRQ(ierr); 34ab92ecdeSBarry Smith ierr = VecRestoreArray(Na->work2,&w2);CHKERRQ(ierr); 35ab92ecdeSBarry Smith 3686eed97dSJose E. Roman if (Na->A) { 3786eed97dSJose E. Roman /* form y = A*x */ 3886eed97dSJose E. Roman ierr = MatMult(Na->A,x,y);CHKERRQ(ierr); 3986eed97dSJose E. Roman /* multiply-add y = y + U*work2 */ 40ab92ecdeSBarry Smith /* note in this call y is treated as a sequential vector */ 41cc2e6a90SBarry Smith ierr = MatMultAdd_SeqDense(Na->U,Na->work2,y,y);CHKERRQ(ierr); 4286eed97dSJose E. Roman } else { 4386eed97dSJose E. Roman /* multiply y = U*work2 */ 4486eed97dSJose E. Roman ierr = MatMult_SeqDense(Na->U,Na->work2,y);CHKERRQ(ierr); 4586eed97dSJose E. Roman } 46ceb9aaf7SBarry Smith PetscFunctionReturn(0); 47ceb9aaf7SBarry Smith } 48ceb9aaf7SBarry Smith 49ceb9aaf7SBarry Smith #undef __FUNCT__ 50ab92ecdeSBarry Smith #define __FUNCT__ "MatDestroy_LRC" 51ab92ecdeSBarry Smith PetscErrorCode MatDestroy_LRC(Mat N) 52ceb9aaf7SBarry Smith { 53ab92ecdeSBarry Smith Mat_LRC *Na = (Mat_LRC*)N->data; 54ceb9aaf7SBarry Smith PetscErrorCode ierr; 55ceb9aaf7SBarry Smith 56ceb9aaf7SBarry Smith PetscFunctionBegin; 57bf0cc555SLisandro Dalcin ierr = MatDestroy(&Na->A);CHKERRQ(ierr); 58bf0cc555SLisandro Dalcin ierr = MatDestroy(&Na->U);CHKERRQ(ierr); 59bf0cc555SLisandro Dalcin ierr = MatDestroy(&Na->V);CHKERRQ(ierr); 606bf464f9SBarry Smith ierr = VecDestroy(&Na->work1);CHKERRQ(ierr); 616bf464f9SBarry Smith ierr = VecDestroy(&Na->work2);CHKERRQ(ierr); 62bf0cc555SLisandro Dalcin ierr = PetscFree(N->data);CHKERRQ(ierr); 63ceb9aaf7SBarry Smith PetscFunctionReturn(0); 64ceb9aaf7SBarry Smith } 65ceb9aaf7SBarry Smith 66ceb9aaf7SBarry Smith #undef __FUNCT__ 67ab92ecdeSBarry Smith #define __FUNCT__ "MatCreateLRC" 68ceb9aaf7SBarry Smith /*@ 69ab92ecdeSBarry Smith MatCreateLRC - Creates a new matrix object that behaves like A + U*V' 70ceb9aaf7SBarry Smith 71ceb9aaf7SBarry Smith Collective on Mat 72ceb9aaf7SBarry Smith 739d9032efSJose E. Roman Input Parameters: 7486eed97dSJose E. Roman + A - the (sparse) matrix (can be NULL) 759d9032efSJose E. Roman - U, V - two dense rectangular (tall and skinny) matrices 76ceb9aaf7SBarry Smith 77ceb9aaf7SBarry Smith Output Parameter: 78ab92ecdeSBarry Smith . N - the matrix that represents A + U*V' 79ceb9aaf7SBarry Smith 809d9032efSJose E. Roman Notes: 819d9032efSJose E. Roman The matrix A + U*V' is not formed! Rather the new matrix 82ceb9aaf7SBarry Smith object performs the matrix-vector product by first multiplying by 839d9032efSJose E. Roman A and then adding the other term. 849d9032efSJose E. Roman 8586eed97dSJose E. Roman If A is NULL then the new object behaves like a low-rank matrix U*V'. 8686eed97dSJose E. Roman 8786eed97dSJose E. Roman Use V=U (or V=NULL) for a symmetric low-rank correction, A + U*U'. 8886eed97dSJose E. Roman 899d9032efSJose E. Roman Level: intermediate 90ceb9aaf7SBarry Smith @*/ 917087cfbeSBarry Smith PetscErrorCode MatCreateLRC(Mat A,Mat U,Mat V,Mat *N) 92ceb9aaf7SBarry Smith { 93ceb9aaf7SBarry Smith PetscErrorCode ierr; 949d9032efSJose E. Roman PetscInt m,n,k,m1,n1,k1; 95ab92ecdeSBarry Smith Mat_LRC *Na; 96ceb9aaf7SBarry Smith 97ceb9aaf7SBarry Smith PetscFunctionBegin; 9886eed97dSJose E. Roman if (A) PetscValidHeaderSpecific(A,MAT_CLASSID,1); 999d9032efSJose E. Roman PetscValidHeaderSpecific(U,MAT_CLASSID,2); 10086eed97dSJose E. Roman if (V) PetscValidHeaderSpecific(V,MAT_CLASSID,3); 10186eed97dSJose E. Roman else V=U; 10286eed97dSJose E. Roman if (A) PetscCheckSameComm(A,1,U,2); 10386eed97dSJose E. Roman PetscCheckSameComm(U,2,V,3); 1049d9032efSJose E. Roman 1059d9032efSJose E. Roman ierr = MatGetSize(U,NULL,&k);CHKERRQ(ierr); 1069d9032efSJose E. Roman ierr = MatGetSize(V,NULL,&k1);CHKERRQ(ierr); 1079d9032efSJose E. Roman if (k!=k1) SETERRQ(PetscObjectComm((PetscObject)U),PETSC_ERR_ARG_INCOMP,"U and V have different number of columns"); 1089d9032efSJose E. Roman ierr = MatGetLocalSize(U,&m,NULL);CHKERRQ(ierr); 1099d9032efSJose E. Roman ierr = MatGetLocalSize(V,&n,NULL);CHKERRQ(ierr); 11086eed97dSJose E. Roman if (A) { 1119d9032efSJose E. Roman ierr = MatGetLocalSize(A,&m1,&n1);CHKERRQ(ierr); 1129d9032efSJose E. Roman if (m!=m1) SETERRQ(PetscObjectComm((PetscObject)U),PETSC_ERR_ARG_INCOMP,"Local dimensions of U and A do not match"); 1139d9032efSJose E. Roman if (n!=n1) SETERRQ(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_INCOMP,"Local dimensions of V and A do not match"); 11486eed97dSJose E. Roman } 1159d9032efSJose E. Roman 11686eed97dSJose E. Roman ierr = MatCreate(PetscObjectComm((PetscObject)U),N);CHKERRQ(ierr); 1179d9032efSJose E. Roman ierr = MatSetSizes(*N,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 118ab92ecdeSBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)*N,MATLRC);CHKERRQ(ierr); 119ceb9aaf7SBarry Smith 120b00a9115SJed Brown ierr = PetscNewLog(*N,&Na);CHKERRQ(ierr); 12138f2d2fdSLisandro Dalcin (*N)->data = (void*)Na; 122ceb9aaf7SBarry Smith Na->A = A; 123ab92ecdeSBarry Smith 124ab92ecdeSBarry Smith ierr = MatDenseGetLocalMatrix(U,&Na->U);CHKERRQ(ierr); 125ab92ecdeSBarry Smith ierr = MatDenseGetLocalMatrix(V,&Na->V);CHKERRQ(ierr); 12686eed97dSJose E. Roman if (A) { ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); } 127ab92ecdeSBarry Smith ierr = PetscObjectReference((PetscObject)Na->U);CHKERRQ(ierr); 128ab92ecdeSBarry Smith ierr = PetscObjectReference((PetscObject)Na->V);CHKERRQ(ierr); 129ceb9aaf7SBarry Smith 130d0f46423SBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,U->cmap->N,&Na->work1);CHKERRQ(ierr); 131ab92ecdeSBarry Smith ierr = VecDuplicate(Na->work1,&Na->work2);CHKERRQ(ierr); 132d0f46423SBarry Smith Na->nwork = U->cmap->N; 133ab92ecdeSBarry Smith 134ab92ecdeSBarry Smith (*N)->ops->destroy = MatDestroy_LRC; 135ab92ecdeSBarry Smith (*N)->ops->mult = MatMult_LRC; 136ceb9aaf7SBarry Smith (*N)->assembled = PETSC_TRUE; 137*26c5da1bSJose E. Roman (*N)->preallocated = PETSC_TRUE; 1389d9032efSJose E. Roman (*N)->cmap->N = V->rmap->N; 1399d9032efSJose E. Roman (*N)->rmap->N = U->rmap->N; 1409d9032efSJose E. Roman (*N)->cmap->n = V->rmap->n; 1419d9032efSJose E. Roman (*N)->rmap->n = U->rmap->n; 142ceb9aaf7SBarry Smith PetscFunctionReturn(0); 143ceb9aaf7SBarry Smith } 144ceb9aaf7SBarry Smith 145