1ceb9aaf7SBarry Smith 2af0996ceSBarry Smith #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/ 3ceb9aaf7SBarry Smith 4ceb9aaf7SBarry Smith typedef struct { 5986c4d72SJose E. Roman Mat A; /* sparse matrix */ 6986c4d72SJose E. Roman Mat U,V; /* dense tall-skinny matrices */ 7986c4d72SJose E. Roman Vec c; /* sequential vector containing the diagonal of C */ 8986c4d72SJose E. Roman Vec work1,work2; /* sequential vectors that hold partial products */ 9ab92ecdeSBarry Smith PetscMPIInt nwork; /* length of work vectors */ 100575051bSJose E. Roman Vec xl,yl; /* auxiliary sequential vectors for matmult operation */ 11ab92ecdeSBarry Smith } Mat_LRC; 12ab92ecdeSBarry Smith 13ab92ecdeSBarry Smith PetscErrorCode MatMult_LRC(Mat N,Vec x,Vec y) 14ceb9aaf7SBarry Smith { 15ab92ecdeSBarry Smith Mat_LRC *Na = (Mat_LRC*)N->data; 16a639b4dcSJose E. Roman Mat Uloc,Vloc; 17ceb9aaf7SBarry Smith PetscErrorCode ierr; 18ab92ecdeSBarry Smith PetscScalar *w1,*w2; 190575051bSJose E. Roman const PetscScalar *a; 20ceb9aaf7SBarry Smith 21ceb9aaf7SBarry Smith PetscFunctionBegin; 220575051bSJose E. Roman ierr = VecGetArrayRead(x,&a);CHKERRQ(ierr); 230575051bSJose E. Roman ierr = VecPlaceArray(Na->xl,a);CHKERRQ(ierr); 240575051bSJose E. Roman ierr = VecGetLocalVector(y,Na->yl);CHKERRQ(ierr); 25a639b4dcSJose E. Roman ierr = MatDenseGetLocalMatrix(Na->U,&Uloc);CHKERRQ(ierr); 26a639b4dcSJose E. Roman ierr = MatDenseGetLocalMatrix(Na->V,&Vloc);CHKERRQ(ierr); 270575051bSJose E. Roman 28ab92ecdeSBarry Smith /* multiply the local part of V with the local part of x */ 290575051bSJose E. Roman #if defined(PETSC_USE_COMPLEX) 30a639b4dcSJose E. Roman ierr = MatMultHermitianTranspose(Vloc,Na->xl,Na->work1);CHKERRQ(ierr); 310575051bSJose E. Roman #else 32a639b4dcSJose E. Roman ierr = MatMultTranspose(Vloc,Na->xl,Na->work1);CHKERRQ(ierr); 330575051bSJose E. Roman #endif 34ab92ecdeSBarry Smith 3586eed97dSJose E. Roman /* form the sum of all the local multiplies: this is work2 = V'*x = 36ab92ecdeSBarry Smith sum_{all processors} work1 */ 37ab92ecdeSBarry Smith ierr = VecGetArray(Na->work1,&w1);CHKERRQ(ierr); 38ab92ecdeSBarry Smith ierr = VecGetArray(Na->work2,&w2);CHKERRQ(ierr); 39820f2d46SBarry Smith ierr = MPIU_Allreduce(w1,w2,Na->nwork,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)N));CHKERRMPI(ierr); 40ab92ecdeSBarry Smith ierr = VecRestoreArray(Na->work1,&w1);CHKERRQ(ierr); 41ab92ecdeSBarry Smith ierr = VecRestoreArray(Na->work2,&w2);CHKERRQ(ierr); 42ab92ecdeSBarry Smith 43986c4d72SJose E. Roman if (Na->c) { /* work2 = C*work2 */ 44986c4d72SJose E. Roman ierr = VecPointwiseMult(Na->work2,Na->c,Na->work2);CHKERRQ(ierr); 45986c4d72SJose E. Roman } 46986c4d72SJose E. Roman 4786eed97dSJose E. Roman if (Na->A) { 4886eed97dSJose E. Roman /* form y = A*x */ 4986eed97dSJose E. Roman ierr = MatMult(Na->A,x,y);CHKERRQ(ierr); 5086eed97dSJose E. Roman /* multiply-add y = y + U*work2 */ 51a639b4dcSJose E. Roman ierr = MatMultAdd(Uloc,Na->work2,Na->yl,Na->yl);CHKERRQ(ierr); 5286eed97dSJose E. Roman } else { 5386eed97dSJose E. Roman /* multiply y = U*work2 */ 54a639b4dcSJose E. Roman ierr = MatMult(Uloc,Na->work2,Na->yl);CHKERRQ(ierr); 5586eed97dSJose E. Roman } 560575051bSJose E. Roman 570575051bSJose E. Roman ierr = VecRestoreArrayRead(x,&a);CHKERRQ(ierr); 580575051bSJose E. Roman ierr = VecResetArray(Na->xl);CHKERRQ(ierr); 590575051bSJose E. Roman ierr = VecRestoreLocalVector(y,Na->yl);CHKERRQ(ierr); 60ceb9aaf7SBarry Smith PetscFunctionReturn(0); 61ceb9aaf7SBarry Smith } 62ceb9aaf7SBarry Smith 63ab92ecdeSBarry Smith PetscErrorCode MatDestroy_LRC(Mat N) 64ceb9aaf7SBarry Smith { 65ab92ecdeSBarry Smith Mat_LRC *Na = (Mat_LRC*)N->data; 66ceb9aaf7SBarry Smith PetscErrorCode ierr; 67ceb9aaf7SBarry Smith 68ceb9aaf7SBarry Smith PetscFunctionBegin; 69bf0cc555SLisandro Dalcin ierr = MatDestroy(&Na->A);CHKERRQ(ierr); 70bf0cc555SLisandro Dalcin ierr = MatDestroy(&Na->U);CHKERRQ(ierr); 71bf0cc555SLisandro Dalcin ierr = MatDestroy(&Na->V);CHKERRQ(ierr); 72986c4d72SJose E. Roman ierr = VecDestroy(&Na->c);CHKERRQ(ierr); 736bf464f9SBarry Smith ierr = VecDestroy(&Na->work1);CHKERRQ(ierr); 746bf464f9SBarry Smith ierr = VecDestroy(&Na->work2);CHKERRQ(ierr); 750575051bSJose E. Roman ierr = VecDestroy(&Na->xl);CHKERRQ(ierr); 760575051bSJose E. Roman ierr = VecDestroy(&Na->yl);CHKERRQ(ierr); 77bf0cc555SLisandro Dalcin ierr = PetscFree(N->data);CHKERRQ(ierr); 78267ca84cSJose E. Roman ierr = PetscObjectComposeFunction((PetscObject)N,"MatLRCGetMats_C",0);CHKERRQ(ierr); 79267ca84cSJose E. Roman PetscFunctionReturn(0); 80267ca84cSJose E. Roman } 81267ca84cSJose E. Roman 82267ca84cSJose E. Roman PetscErrorCode MatLRCGetMats_LRC(Mat N,Mat *A,Mat *U,Vec *c,Mat *V) 83267ca84cSJose E. Roman { 84267ca84cSJose E. Roman Mat_LRC *Na = (Mat_LRC*)N->data; 85267ca84cSJose E. Roman 86267ca84cSJose E. Roman PetscFunctionBegin; 87267ca84cSJose E. Roman if (A) *A = Na->A; 88267ca84cSJose E. Roman if (U) *U = Na->U; 89267ca84cSJose E. Roman if (c) *c = Na->c; 90267ca84cSJose E. Roman if (V) *V = Na->V; 91267ca84cSJose E. Roman PetscFunctionReturn(0); 92267ca84cSJose E. Roman } 93267ca84cSJose E. Roman 94267ca84cSJose E. Roman /*@ 95267ca84cSJose E. Roman MatLRCGetMats - Returns the constituents of an LRC matrix 96267ca84cSJose E. Roman 97267ca84cSJose E. Roman Collective on Mat 98267ca84cSJose E. Roman 99267ca84cSJose E. Roman Input Parameter: 100267ca84cSJose E. Roman . N - matrix of type LRC 101267ca84cSJose E. Roman 102267ca84cSJose E. Roman Output Parameters: 103267ca84cSJose E. Roman + A - the (sparse) matrix 1046b867d5aSJose E. Roman . U - first dense rectangular (tall and skinny) matrix 1056b867d5aSJose E. Roman . c - a sequential vector containing the diagonal of C 1066b867d5aSJose E. Roman - V - second dense rectangular (tall and skinny) matrix 107267ca84cSJose E. Roman 108267ca84cSJose E. Roman Note: 109267ca84cSJose E. Roman The returned matrices need not be destroyed by the caller. 110267ca84cSJose E. Roman 111267ca84cSJose E. Roman Level: intermediate 112267ca84cSJose E. Roman 113267ca84cSJose E. Roman .seealso: MatCreateLRC() 114267ca84cSJose E. Roman @*/ 115267ca84cSJose E. Roman PetscErrorCode MatLRCGetMats(Mat N,Mat *A,Mat *U,Vec *c,Mat *V) 116267ca84cSJose E. Roman { 117267ca84cSJose E. Roman PetscErrorCode ierr; 118267ca84cSJose E. Roman 119267ca84cSJose E. Roman PetscFunctionBegin; 120267ca84cSJose E. Roman ierr = PetscUseMethod(N,"MatLRCGetMats_C",(Mat,Mat*,Mat*,Vec*,Mat*),(N,A,U,c,V));CHKERRQ(ierr); 121ceb9aaf7SBarry Smith PetscFunctionReturn(0); 122ceb9aaf7SBarry Smith } 123ceb9aaf7SBarry Smith 124ceb9aaf7SBarry Smith /*@ 125986c4d72SJose E. Roman MatCreateLRC - Creates a new matrix object that behaves like A + U*C*V' 126ceb9aaf7SBarry Smith 127ceb9aaf7SBarry Smith Collective on Mat 128ceb9aaf7SBarry Smith 1299d9032efSJose E. Roman Input Parameters: 13086eed97dSJose E. Roman + A - the (sparse) matrix (can be NULL) 131986c4d72SJose E. Roman . U, V - two dense rectangular (tall and skinny) matrices 132986c4d72SJose E. Roman - c - a sequential vector containing the diagonal of C (can be NULL) 133ceb9aaf7SBarry Smith 134ceb9aaf7SBarry Smith Output Parameter: 135986c4d72SJose E. Roman . N - the matrix that represents A + U*C*V' 136ceb9aaf7SBarry Smith 1379d9032efSJose E. Roman Notes: 138986c4d72SJose E. Roman The matrix A + U*C*V' is not formed! Rather the new matrix 139ceb9aaf7SBarry Smith object performs the matrix-vector product by first multiplying by 1409d9032efSJose E. Roman A and then adding the other term. 1419d9032efSJose E. Roman 142986c4d72SJose E. Roman C is a diagonal matrix (represented as a vector) of order k, 143986c4d72SJose E. Roman where k is the number of columns of both U and V. 14486eed97dSJose E. Roman 145986c4d72SJose E. Roman If A is NULL then the new object behaves like a low-rank matrix U*C*V'. 146986c4d72SJose E. Roman 147986c4d72SJose E. Roman Use V=U (or V=NULL) for a symmetric low-rank correction, A + U*C*U'. 148986c4d72SJose E. Roman 149986c4d72SJose E. Roman If c is NULL then the low-rank correction is just U*V'. 15086eed97dSJose E. Roman 1519d9032efSJose E. Roman Level: intermediate 152267ca84cSJose E. Roman 153267ca84cSJose E. Roman .seealso: MatLRCGetMats() 154ceb9aaf7SBarry Smith @*/ 155986c4d72SJose E. Roman PetscErrorCode MatCreateLRC(Mat A,Mat U,Vec c,Mat V,Mat *N) 156ceb9aaf7SBarry Smith { 157ceb9aaf7SBarry Smith PetscErrorCode ierr; 158986c4d72SJose E. Roman PetscBool match; 1599d9032efSJose E. Roman PetscInt m,n,k,m1,n1,k1; 160ab92ecdeSBarry Smith Mat_LRC *Na; 161ceb9aaf7SBarry Smith 162ceb9aaf7SBarry Smith PetscFunctionBegin; 16386eed97dSJose E. Roman if (A) PetscValidHeaderSpecific(A,MAT_CLASSID,1); 1649d9032efSJose E. Roman PetscValidHeaderSpecific(U,MAT_CLASSID,2); 165986c4d72SJose E. Roman if (c) PetscValidHeaderSpecific(c,VEC_CLASSID,3); 166986c4d72SJose E. Roman if (V) PetscValidHeaderSpecific(V,MAT_CLASSID,4); 16786eed97dSJose E. Roman else V=U; 16886eed97dSJose E. Roman if (A) PetscCheckSameComm(A,1,U,2); 169986c4d72SJose E. Roman PetscCheckSameComm(U,2,V,4); 170986c4d72SJose E. Roman 171986c4d72SJose E. Roman ierr = PetscObjectTypeCompareAny((PetscObject)U,&match,MATSEQDENSE,MATMPIDENSE,"");CHKERRQ(ierr); 172*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(!match,PetscObjectComm((PetscObject)U),PETSC_ERR_SUP,"Matrix U must be of type dense"); 173986c4d72SJose E. Roman if (V) { 174986c4d72SJose E. Roman ierr = PetscObjectTypeCompareAny((PetscObject)V,&match,MATSEQDENSE,MATMPIDENSE,"");CHKERRQ(ierr); 175*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(!match,PetscObjectComm((PetscObject)V),PETSC_ERR_SUP,"Matrix V must be of type dense"); 176986c4d72SJose E. Roman } 1779d9032efSJose E. Roman 1789d9032efSJose E. Roman ierr = MatGetSize(U,NULL,&k);CHKERRQ(ierr); 1799d9032efSJose E. Roman ierr = MatGetSize(V,NULL,&k1);CHKERRQ(ierr); 180*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(k!=k1,PetscObjectComm((PetscObject)U),PETSC_ERR_ARG_INCOMP,"U and V have different number of columns (%" PetscInt_FMT " vs %" PetscInt_FMT ")",k,k1); 1819d9032efSJose E. Roman ierr = MatGetLocalSize(U,&m,NULL);CHKERRQ(ierr); 1829d9032efSJose E. Roman ierr = MatGetLocalSize(V,&n,NULL);CHKERRQ(ierr); 18386eed97dSJose E. Roman if (A) { 1849d9032efSJose E. Roman ierr = MatGetLocalSize(A,&m1,&n1);CHKERRQ(ierr); 185*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(m!=m1,PetscObjectComm((PetscObject)U),PETSC_ERR_ARG_INCOMP,"Local dimensions of U %" PetscInt_FMT " and A %" PetscInt_FMT " do not match",m,m1); 186*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(n!=n1,PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_INCOMP,"Local dimensions of V %" PetscInt_FMT " and A %" PetscInt_FMT " do not match",n,n1); 18786eed97dSJose E. Roman } 188986c4d72SJose E. Roman if (c) { 189986c4d72SJose E. Roman ierr = VecGetSize(c,&k1);CHKERRQ(ierr); 190*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(k!=k1,PetscObjectComm((PetscObject)c),PETSC_ERR_ARG_INCOMP,"The length of c %" PetscInt_FMT " does not match the number of columns of U and V (%" PetscInt_FMT ")",k1,k); 191986c4d72SJose E. Roman ierr = VecGetLocalSize(c,&k1);CHKERRQ(ierr); 192*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(k!=k1,PetscObjectComm((PetscObject)c),PETSC_ERR_ARG_INCOMP,"c must be a sequential vector"); 193986c4d72SJose E. Roman } 1949d9032efSJose E. Roman 19586eed97dSJose E. Roman ierr = MatCreate(PetscObjectComm((PetscObject)U),N);CHKERRQ(ierr); 1969d9032efSJose E. Roman ierr = MatSetSizes(*N,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 197ab92ecdeSBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)*N,MATLRC);CHKERRQ(ierr); 198ceb9aaf7SBarry Smith 199b00a9115SJed Brown ierr = PetscNewLog(*N,&Na);CHKERRQ(ierr); 20038f2d2fdSLisandro Dalcin (*N)->data = (void*)Na; 201ceb9aaf7SBarry Smith Na->A = A; 202a639b4dcSJose E. Roman Na->U = U; 203986c4d72SJose E. Roman Na->c = c; 204a639b4dcSJose E. Roman Na->V = V; 205ab92ecdeSBarry Smith 20686eed97dSJose E. Roman if (A) { ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); } 207ab92ecdeSBarry Smith ierr = PetscObjectReference((PetscObject)Na->U);CHKERRQ(ierr); 208ab92ecdeSBarry Smith ierr = PetscObjectReference((PetscObject)Na->V);CHKERRQ(ierr); 209986c4d72SJose E. Roman if (c) { ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr); } 210ceb9aaf7SBarry Smith 211d0f46423SBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,U->cmap->N,&Na->work1);CHKERRQ(ierr); 212ab92ecdeSBarry Smith ierr = VecDuplicate(Na->work1,&Na->work2);CHKERRQ(ierr); 213986c4d72SJose E. Roman ierr = PetscMPIIntCast(U->cmap->N,&Na->nwork);CHKERRQ(ierr); 214ab92ecdeSBarry Smith 2150575051bSJose E. Roman ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,V->rmap->n,NULL,&Na->xl);CHKERRQ(ierr); 2160575051bSJose E. Roman ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,U->rmap->n,NULL,&Na->yl);CHKERRQ(ierr); 2170575051bSJose E. Roman 218ab92ecdeSBarry Smith (*N)->ops->destroy = MatDestroy_LRC; 219ab92ecdeSBarry Smith (*N)->ops->mult = MatMult_LRC; 220ceb9aaf7SBarry Smith (*N)->assembled = PETSC_TRUE; 22126c5da1bSJose E. Roman (*N)->preallocated = PETSC_TRUE; 2229d9032efSJose E. Roman (*N)->cmap->N = V->rmap->N; 2239d9032efSJose E. Roman (*N)->rmap->N = U->rmap->N; 2249d9032efSJose E. Roman (*N)->cmap->n = V->rmap->n; 2259d9032efSJose E. Roman (*N)->rmap->n = U->rmap->n; 226267ca84cSJose E. Roman 227267ca84cSJose E. Roman ierr = PetscObjectComposeFunction((PetscObject)(*N),"MatLRCGetMats_C",MatLRCGetMats_LRC);CHKERRQ(ierr); 228ceb9aaf7SBarry Smith PetscFunctionReturn(0); 229ceb9aaf7SBarry Smith } 230