1ceb9aaf7SBarry Smith 2af0996ceSBarry Smith #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/ 3ceb9aaf7SBarry Smith 4ceb9aaf7SBarry Smith typedef struct { 5*986c4d72SJose E. Roman Mat A; /* sparse matrix */ 6*986c4d72SJose E. Roman Mat U,V; /* dense tall-skinny matrices */ 7*986c4d72SJose E. Roman Vec c; /* sequential vector containing the diagonal of C */ 8*986c4d72SJose E. Roman Vec work1,work2; /* sequential vectors that hold partial products */ 9ab92ecdeSBarry Smith PetscMPIInt nwork; /* length of work vectors */ 10ab92ecdeSBarry Smith } Mat_LRC; 11ab92ecdeSBarry Smith 12ab92ecdeSBarry Smith 13cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatMultTranspose_SeqDense(Mat,Vec,Vec); 14cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatMultAdd_SeqDense(Mat,Vec,Vec,Vec); 1586eed97dSJose E. Roman PETSC_INTERN PetscErrorCode MatMult_SeqDense(Mat,Vec,Vec); 16ceb9aaf7SBarry Smith 17ceb9aaf7SBarry Smith #undef __FUNCT__ 18ab92ecdeSBarry Smith #define __FUNCT__ "MatMult_LRC" 19ab92ecdeSBarry Smith PetscErrorCode MatMult_LRC(Mat N,Vec x,Vec y) 20ceb9aaf7SBarry Smith { 21ab92ecdeSBarry Smith Mat_LRC *Na = (Mat_LRC*)N->data; 22ceb9aaf7SBarry Smith PetscErrorCode ierr; 23ab92ecdeSBarry Smith PetscScalar *w1,*w2; 24ceb9aaf7SBarry Smith 25ceb9aaf7SBarry Smith PetscFunctionBegin; 26ab92ecdeSBarry Smith /* multiply the local part of V with the local part of x */ 27ab92ecdeSBarry Smith /* note in this call x is treated as a sequential vector */ 28cc2e6a90SBarry Smith ierr = MatMultTranspose_SeqDense(Na->V,x,Na->work1);CHKERRQ(ierr); 29ab92ecdeSBarry Smith 3086eed97dSJose E. Roman /* form the sum of all the local multiplies: this is work2 = V'*x = 31ab92ecdeSBarry Smith sum_{all processors} work1 */ 32ab92ecdeSBarry Smith ierr = VecGetArray(Na->work1,&w1);CHKERRQ(ierr); 33ab92ecdeSBarry Smith ierr = VecGetArray(Na->work2,&w2);CHKERRQ(ierr); 34b2566f29SBarry Smith ierr = MPIU_Allreduce(w1,w2,Na->nwork,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)N));CHKERRQ(ierr); 35ab92ecdeSBarry Smith ierr = VecRestoreArray(Na->work1,&w1);CHKERRQ(ierr); 36ab92ecdeSBarry Smith ierr = VecRestoreArray(Na->work2,&w2);CHKERRQ(ierr); 37ab92ecdeSBarry Smith 38*986c4d72SJose E. Roman if (Na->c) { /* work2 = C*work2 */ 39*986c4d72SJose E. Roman ierr = VecPointwiseMult(Na->work2,Na->c,Na->work2);CHKERRQ(ierr); 40*986c4d72SJose E. Roman } 41*986c4d72SJose E. Roman 4286eed97dSJose E. Roman if (Na->A) { 4386eed97dSJose E. Roman /* form y = A*x */ 4486eed97dSJose E. Roman ierr = MatMult(Na->A,x,y);CHKERRQ(ierr); 4586eed97dSJose E. Roman /* multiply-add y = y + U*work2 */ 46ab92ecdeSBarry Smith /* note in this call y is treated as a sequential vector */ 47cc2e6a90SBarry Smith ierr = MatMultAdd_SeqDense(Na->U,Na->work2,y,y);CHKERRQ(ierr); 4886eed97dSJose E. Roman } else { 4986eed97dSJose E. Roman /* multiply y = U*work2 */ 5086eed97dSJose E. Roman ierr = MatMult_SeqDense(Na->U,Na->work2,y);CHKERRQ(ierr); 5186eed97dSJose E. Roman } 52ceb9aaf7SBarry Smith PetscFunctionReturn(0); 53ceb9aaf7SBarry Smith } 54ceb9aaf7SBarry Smith 55ceb9aaf7SBarry Smith #undef __FUNCT__ 56ab92ecdeSBarry Smith #define __FUNCT__ "MatDestroy_LRC" 57ab92ecdeSBarry Smith PetscErrorCode MatDestroy_LRC(Mat N) 58ceb9aaf7SBarry Smith { 59ab92ecdeSBarry Smith Mat_LRC *Na = (Mat_LRC*)N->data; 60ceb9aaf7SBarry Smith PetscErrorCode ierr; 61ceb9aaf7SBarry Smith 62ceb9aaf7SBarry Smith PetscFunctionBegin; 63bf0cc555SLisandro Dalcin ierr = MatDestroy(&Na->A);CHKERRQ(ierr); 64bf0cc555SLisandro Dalcin ierr = MatDestroy(&Na->U);CHKERRQ(ierr); 65bf0cc555SLisandro Dalcin ierr = MatDestroy(&Na->V);CHKERRQ(ierr); 66*986c4d72SJose E. Roman ierr = VecDestroy(&Na->c);CHKERRQ(ierr); 676bf464f9SBarry Smith ierr = VecDestroy(&Na->work1);CHKERRQ(ierr); 686bf464f9SBarry Smith ierr = VecDestroy(&Na->work2);CHKERRQ(ierr); 69bf0cc555SLisandro Dalcin ierr = PetscFree(N->data);CHKERRQ(ierr); 70ceb9aaf7SBarry Smith PetscFunctionReturn(0); 71ceb9aaf7SBarry Smith } 72ceb9aaf7SBarry Smith 73ceb9aaf7SBarry Smith #undef __FUNCT__ 74ab92ecdeSBarry Smith #define __FUNCT__ "MatCreateLRC" 75ceb9aaf7SBarry Smith /*@ 76*986c4d72SJose E. Roman MatCreateLRC - Creates a new matrix object that behaves like A + U*C*V' 77ceb9aaf7SBarry Smith 78ceb9aaf7SBarry Smith Collective on Mat 79ceb9aaf7SBarry Smith 809d9032efSJose E. Roman Input Parameters: 8186eed97dSJose E. Roman + A - the (sparse) matrix (can be NULL) 82*986c4d72SJose E. Roman . U, V - two dense rectangular (tall and skinny) matrices 83*986c4d72SJose E. Roman - c - a sequential vector containing the diagonal of C (can be NULL) 84ceb9aaf7SBarry Smith 85ceb9aaf7SBarry Smith Output Parameter: 86*986c4d72SJose E. Roman . N - the matrix that represents A + U*C*V' 87ceb9aaf7SBarry Smith 889d9032efSJose E. Roman Notes: 89*986c4d72SJose E. Roman The matrix A + U*C*V' is not formed! Rather the new matrix 90ceb9aaf7SBarry Smith object performs the matrix-vector product by first multiplying by 919d9032efSJose E. Roman A and then adding the other term. 929d9032efSJose E. Roman 93*986c4d72SJose E. Roman C is a diagonal matrix (represented as a vector) of order k, 94*986c4d72SJose E. Roman where k is the number of columns of both U and V. 9586eed97dSJose E. Roman 96*986c4d72SJose E. Roman If A is NULL then the new object behaves like a low-rank matrix U*C*V'. 97*986c4d72SJose E. Roman 98*986c4d72SJose E. Roman Use V=U (or V=NULL) for a symmetric low-rank correction, A + U*C*U'. 99*986c4d72SJose E. Roman 100*986c4d72SJose E. Roman If c is NULL then the low-rank correction is just U*V'. 10186eed97dSJose E. Roman 1029d9032efSJose E. Roman Level: intermediate 103ceb9aaf7SBarry Smith @*/ 104*986c4d72SJose E. Roman PetscErrorCode MatCreateLRC(Mat A,Mat U,Vec c,Mat V,Mat *N) 105ceb9aaf7SBarry Smith { 106ceb9aaf7SBarry Smith PetscErrorCode ierr; 107*986c4d72SJose E. Roman PetscBool match; 1089d9032efSJose E. Roman PetscInt m,n,k,m1,n1,k1; 109ab92ecdeSBarry Smith Mat_LRC *Na; 110ceb9aaf7SBarry Smith 111ceb9aaf7SBarry Smith PetscFunctionBegin; 11286eed97dSJose E. Roman if (A) PetscValidHeaderSpecific(A,MAT_CLASSID,1); 1139d9032efSJose E. Roman PetscValidHeaderSpecific(U,MAT_CLASSID,2); 114*986c4d72SJose E. Roman if (c) PetscValidHeaderSpecific(c,VEC_CLASSID,3); 115*986c4d72SJose E. Roman if (V) PetscValidHeaderSpecific(V,MAT_CLASSID,4); 11686eed97dSJose E. Roman else V=U; 11786eed97dSJose E. Roman if (A) PetscCheckSameComm(A,1,U,2); 118*986c4d72SJose E. Roman PetscCheckSameComm(U,2,V,4); 119*986c4d72SJose E. Roman 120*986c4d72SJose E. Roman ierr = PetscObjectTypeCompareAny((PetscObject)U,&match,MATSEQDENSE,MATMPIDENSE,"");CHKERRQ(ierr); 121*986c4d72SJose E. Roman if (!match) SETERRQ(PetscObjectComm((PetscObject)U),PETSC_ERR_SUP,"Matrix U must be of type dense"); 122*986c4d72SJose E. Roman if (V) { 123*986c4d72SJose E. Roman ierr = PetscObjectTypeCompareAny((PetscObject)V,&match,MATSEQDENSE,MATMPIDENSE,"");CHKERRQ(ierr); 124*986c4d72SJose E. Roman if (!match) SETERRQ(PetscObjectComm((PetscObject)V),PETSC_ERR_SUP,"Matrix V must be of type dense"); 125*986c4d72SJose E. Roman } 1269d9032efSJose E. Roman 1279d9032efSJose E. Roman ierr = MatGetSize(U,NULL,&k);CHKERRQ(ierr); 1289d9032efSJose E. Roman ierr = MatGetSize(V,NULL,&k1);CHKERRQ(ierr); 1299d9032efSJose E. Roman if (k!=k1) SETERRQ(PetscObjectComm((PetscObject)U),PETSC_ERR_ARG_INCOMP,"U and V have different number of columns"); 1309d9032efSJose E. Roman ierr = MatGetLocalSize(U,&m,NULL);CHKERRQ(ierr); 1319d9032efSJose E. Roman ierr = MatGetLocalSize(V,&n,NULL);CHKERRQ(ierr); 13286eed97dSJose E. Roman if (A) { 1339d9032efSJose E. Roman ierr = MatGetLocalSize(A,&m1,&n1);CHKERRQ(ierr); 1349d9032efSJose E. Roman if (m!=m1) SETERRQ(PetscObjectComm((PetscObject)U),PETSC_ERR_ARG_INCOMP,"Local dimensions of U and A do not match"); 1359d9032efSJose E. Roman if (n!=n1) SETERRQ(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_INCOMP,"Local dimensions of V and A do not match"); 13686eed97dSJose E. Roman } 137*986c4d72SJose E. Roman if (c) { 138*986c4d72SJose E. Roman ierr = VecGetSize(c,&k1);CHKERRQ(ierr); 139*986c4d72SJose E. Roman if (k!=k1) SETERRQ(PetscObjectComm((PetscObject)c),PETSC_ERR_ARG_INCOMP,"The length of c does not match the number of columns of U and V"); 140*986c4d72SJose E. Roman ierr = VecGetLocalSize(c,&k1);CHKERRQ(ierr); 141*986c4d72SJose E. Roman if (k!=k1) SETERRQ(PetscObjectComm((PetscObject)c),PETSC_ERR_ARG_INCOMP,"c must be a sequential vector"); 142*986c4d72SJose E. Roman } 1439d9032efSJose E. Roman 14486eed97dSJose E. Roman ierr = MatCreate(PetscObjectComm((PetscObject)U),N);CHKERRQ(ierr); 1459d9032efSJose E. Roman ierr = MatSetSizes(*N,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 146ab92ecdeSBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)*N,MATLRC);CHKERRQ(ierr); 147ceb9aaf7SBarry Smith 148b00a9115SJed Brown ierr = PetscNewLog(*N,&Na);CHKERRQ(ierr); 14938f2d2fdSLisandro Dalcin (*N)->data = (void*)Na; 150ceb9aaf7SBarry Smith Na->A = A; 151*986c4d72SJose E. Roman Na->c = c; 152ab92ecdeSBarry Smith 153ab92ecdeSBarry Smith ierr = MatDenseGetLocalMatrix(U,&Na->U);CHKERRQ(ierr); 154ab92ecdeSBarry Smith ierr = MatDenseGetLocalMatrix(V,&Na->V);CHKERRQ(ierr); 15586eed97dSJose E. Roman if (A) { ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); } 156ab92ecdeSBarry Smith ierr = PetscObjectReference((PetscObject)Na->U);CHKERRQ(ierr); 157ab92ecdeSBarry Smith ierr = PetscObjectReference((PetscObject)Na->V);CHKERRQ(ierr); 158*986c4d72SJose E. Roman if (c) { ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr); } 159ceb9aaf7SBarry Smith 160d0f46423SBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,U->cmap->N,&Na->work1);CHKERRQ(ierr); 161ab92ecdeSBarry Smith ierr = VecDuplicate(Na->work1,&Na->work2);CHKERRQ(ierr); 162*986c4d72SJose E. Roman ierr = PetscMPIIntCast(U->cmap->N,&Na->nwork);CHKERRQ(ierr); 163ab92ecdeSBarry Smith 164ab92ecdeSBarry Smith (*N)->ops->destroy = MatDestroy_LRC; 165ab92ecdeSBarry Smith (*N)->ops->mult = MatMult_LRC; 166ceb9aaf7SBarry Smith (*N)->assembled = PETSC_TRUE; 16726c5da1bSJose E. Roman (*N)->preallocated = PETSC_TRUE; 1689d9032efSJose E. Roman (*N)->cmap->N = V->rmap->N; 1699d9032efSJose E. Roman (*N)->rmap->N = U->rmap->N; 1709d9032efSJose E. Roman (*N)->cmap->n = V->rmap->n; 1719d9032efSJose E. Roman (*N)->rmap->n = U->rmap->n; 172ceb9aaf7SBarry Smith PetscFunctionReturn(0); 173ceb9aaf7SBarry Smith } 174ceb9aaf7SBarry Smith 175