xref: /petsc/src/mat/impls/lrc/lrc.c (revision 267ca84cfdc7de74c397a288718dfdf4253e4090)
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 
14ceb9aaf7SBarry Smith #undef __FUNCT__
15ab92ecdeSBarry Smith #define __FUNCT__ "MatMult_LRC"
16ab92ecdeSBarry Smith PetscErrorCode MatMult_LRC(Mat N,Vec x,Vec y)
17ceb9aaf7SBarry Smith {
18ab92ecdeSBarry Smith   Mat_LRC           *Na = (Mat_LRC*)N->data;
19ceb9aaf7SBarry Smith   PetscErrorCode    ierr;
20ab92ecdeSBarry Smith   PetscScalar       *w1,*w2;
210575051bSJose E. Roman   const PetscScalar *a;
22ceb9aaf7SBarry Smith 
23ceb9aaf7SBarry Smith   PetscFunctionBegin;
240575051bSJose E. Roman   ierr = VecGetArrayRead(x,&a);CHKERRQ(ierr);
250575051bSJose E. Roman   ierr = VecPlaceArray(Na->xl,a);CHKERRQ(ierr);
260575051bSJose E. Roman   ierr = VecGetLocalVector(y,Na->yl);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)
300575051bSJose E. Roman   ierr = MatMultHermitianTranspose(Na->V,Na->xl,Na->work1);CHKERRQ(ierr);
310575051bSJose E. Roman #else
320575051bSJose E. Roman   ierr = MatMultTranspose(Na->V,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);
39b2566f29SBarry Smith   ierr = MPIU_Allreduce(w1,w2,Na->nwork,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)N));CHKERRQ(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 */
510575051bSJose E. Roman     ierr = MatMultAdd(Na->U,Na->work2,Na->yl,Na->yl);CHKERRQ(ierr);
5286eed97dSJose E. Roman   } else {
5386eed97dSJose E. Roman     /* multiply y = U*work2 */
540575051bSJose E. Roman     ierr = MatMult(Na->U,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 
63ceb9aaf7SBarry Smith #undef __FUNCT__
64ab92ecdeSBarry Smith #define __FUNCT__ "MatDestroy_LRC"
65ab92ecdeSBarry Smith PetscErrorCode MatDestroy_LRC(Mat N)
66ceb9aaf7SBarry Smith {
67ab92ecdeSBarry Smith   Mat_LRC        *Na = (Mat_LRC*)N->data;
68ceb9aaf7SBarry Smith   PetscErrorCode ierr;
69ceb9aaf7SBarry Smith 
70ceb9aaf7SBarry Smith   PetscFunctionBegin;
71bf0cc555SLisandro Dalcin   ierr = MatDestroy(&Na->A);CHKERRQ(ierr);
72bf0cc555SLisandro Dalcin   ierr = MatDestroy(&Na->U);CHKERRQ(ierr);
73bf0cc555SLisandro Dalcin   ierr = MatDestroy(&Na->V);CHKERRQ(ierr);
74986c4d72SJose E. Roman   ierr = VecDestroy(&Na->c);CHKERRQ(ierr);
756bf464f9SBarry Smith   ierr = VecDestroy(&Na->work1);CHKERRQ(ierr);
766bf464f9SBarry Smith   ierr = VecDestroy(&Na->work2);CHKERRQ(ierr);
770575051bSJose E. Roman   ierr = VecDestroy(&Na->xl);CHKERRQ(ierr);
780575051bSJose E. Roman   ierr = VecDestroy(&Na->yl);CHKERRQ(ierr);
79bf0cc555SLisandro Dalcin   ierr = PetscFree(N->data);CHKERRQ(ierr);
80*267ca84cSJose E. Roman   ierr = PetscObjectComposeFunction((PetscObject)N,"MatLRCGetMats_C",0);CHKERRQ(ierr);
81*267ca84cSJose E. Roman   PetscFunctionReturn(0);
82*267ca84cSJose E. Roman }
83*267ca84cSJose E. Roman 
84*267ca84cSJose E. Roman #undef __FUNCT__
85*267ca84cSJose E. Roman #define __FUNCT__ "MatLRCGetMats_LRC"
86*267ca84cSJose E. Roman PetscErrorCode MatLRCGetMats_LRC(Mat N,Mat *A,Mat *U,Vec *c,Mat *V)
87*267ca84cSJose E. Roman {
88*267ca84cSJose E. Roman   Mat_LRC *Na = (Mat_LRC*)N->data;
89*267ca84cSJose E. Roman 
90*267ca84cSJose E. Roman   PetscFunctionBegin;
91*267ca84cSJose E. Roman   if (A) *A = Na->A;
92*267ca84cSJose E. Roman   if (U) *U = Na->U;
93*267ca84cSJose E. Roman   if (c) *c = Na->c;
94*267ca84cSJose E. Roman   if (V) *V = Na->V;
95*267ca84cSJose E. Roman   PetscFunctionReturn(0);
96*267ca84cSJose E. Roman }
97*267ca84cSJose E. Roman 
98*267ca84cSJose E. Roman #undef __FUNCT__
99*267ca84cSJose E. Roman #define __FUNCT__ "MatLRCGetMats"
100*267ca84cSJose E. Roman /*@
101*267ca84cSJose E. Roman    MatLRCGetMats - Returns the constituents of an LRC matrix
102*267ca84cSJose E. Roman 
103*267ca84cSJose E. Roman    Collective on Mat
104*267ca84cSJose E. Roman 
105*267ca84cSJose E. Roman    Input Parameter:
106*267ca84cSJose E. Roman .  N    - matrix of type LRC
107*267ca84cSJose E. Roman 
108*267ca84cSJose E. Roman    Output Parameters:
109*267ca84cSJose E. Roman +  A    - the (sparse) matrix
110*267ca84cSJose E. Roman .  U, V - two dense rectangular (tall and skinny) matrices
111*267ca84cSJose E. Roman -  c    - a sequential vector containing the diagonal of C
112*267ca84cSJose E. Roman 
113*267ca84cSJose E. Roman    Note:
114*267ca84cSJose E. Roman    The returned matrices need not be destroyed by the caller.
115*267ca84cSJose E. Roman 
116*267ca84cSJose E. Roman    Level: intermediate
117*267ca84cSJose E. Roman 
118*267ca84cSJose E. Roman .seealso: MatCreateLRC()
119*267ca84cSJose E. Roman @*/
120*267ca84cSJose E. Roman PetscErrorCode  MatLRCGetMats(Mat N,Mat *A,Mat *U,Vec *c,Mat *V)
121*267ca84cSJose E. Roman {
122*267ca84cSJose E. Roman   PetscErrorCode ierr;
123*267ca84cSJose E. Roman 
124*267ca84cSJose E. Roman   PetscFunctionBegin;
125*267ca84cSJose E. Roman   ierr = PetscUseMethod(N,"MatLRCGetMats_C",(Mat,Mat*,Mat*,Vec*,Mat*),(N,A,U,c,V));CHKERRQ(ierr);
126ceb9aaf7SBarry Smith   PetscFunctionReturn(0);
127ceb9aaf7SBarry Smith }
128ceb9aaf7SBarry Smith 
129ceb9aaf7SBarry Smith #undef __FUNCT__
130ab92ecdeSBarry Smith #define __FUNCT__ "MatCreateLRC"
131ceb9aaf7SBarry Smith /*@
132986c4d72SJose E. Roman    MatCreateLRC - Creates a new matrix object that behaves like A + U*C*V'
133ceb9aaf7SBarry Smith 
134ceb9aaf7SBarry Smith    Collective on Mat
135ceb9aaf7SBarry Smith 
1369d9032efSJose E. Roman    Input Parameters:
13786eed97dSJose E. Roman +  A    - the (sparse) matrix (can be NULL)
138986c4d72SJose E. Roman .  U, V - two dense rectangular (tall and skinny) matrices
139986c4d72SJose E. Roman -  c    - a sequential vector containing the diagonal of C (can be NULL)
140ceb9aaf7SBarry Smith 
141ceb9aaf7SBarry Smith    Output Parameter:
142986c4d72SJose E. Roman .  N    - the matrix that represents A + U*C*V'
143ceb9aaf7SBarry Smith 
1449d9032efSJose E. Roman    Notes:
145986c4d72SJose E. Roman    The matrix A + U*C*V' is not formed! Rather the new matrix
146ceb9aaf7SBarry Smith    object performs the matrix-vector product by first multiplying by
1479d9032efSJose E. Roman    A and then adding the other term.
1489d9032efSJose E. Roman 
149986c4d72SJose E. Roman    C is a diagonal matrix (represented as a vector) of order k,
150986c4d72SJose E. Roman    where k is the number of columns of both U and V.
15186eed97dSJose E. Roman 
152986c4d72SJose E. Roman    If A is NULL then the new object behaves like a low-rank matrix U*C*V'.
153986c4d72SJose E. Roman 
154986c4d72SJose E. Roman    Use V=U (or V=NULL) for a symmetric low-rank correction, A + U*C*U'.
155986c4d72SJose E. Roman 
156986c4d72SJose E. Roman    If c is NULL then the low-rank correction is just U*V'.
15786eed97dSJose E. Roman 
1589d9032efSJose E. Roman    Level: intermediate
159*267ca84cSJose E. Roman 
160*267ca84cSJose E. Roman .seealso: MatLRCGetMats()
161ceb9aaf7SBarry Smith @*/
162986c4d72SJose E. Roman PetscErrorCode MatCreateLRC(Mat A,Mat U,Vec c,Mat V,Mat *N)
163ceb9aaf7SBarry Smith {
164ceb9aaf7SBarry Smith   PetscErrorCode ierr;
165986c4d72SJose E. Roman   PetscBool      match;
1669d9032efSJose E. Roman   PetscInt       m,n,k,m1,n1,k1;
167ab92ecdeSBarry Smith   Mat_LRC        *Na;
168ceb9aaf7SBarry Smith 
169ceb9aaf7SBarry Smith   PetscFunctionBegin;
17086eed97dSJose E. Roman   if (A) PetscValidHeaderSpecific(A,MAT_CLASSID,1);
1719d9032efSJose E. Roman   PetscValidHeaderSpecific(U,MAT_CLASSID,2);
172986c4d72SJose E. Roman   if (c) PetscValidHeaderSpecific(c,VEC_CLASSID,3);
173986c4d72SJose E. Roman   if (V) PetscValidHeaderSpecific(V,MAT_CLASSID,4);
17486eed97dSJose E. Roman   else V=U;
17586eed97dSJose E. Roman   if (A) PetscCheckSameComm(A,1,U,2);
176986c4d72SJose E. Roman   PetscCheckSameComm(U,2,V,4);
177986c4d72SJose E. Roman 
178986c4d72SJose E. Roman   ierr = PetscObjectTypeCompareAny((PetscObject)U,&match,MATSEQDENSE,MATMPIDENSE,"");CHKERRQ(ierr);
179986c4d72SJose E. Roman   if (!match) SETERRQ(PetscObjectComm((PetscObject)U),PETSC_ERR_SUP,"Matrix U must be of type dense");
180986c4d72SJose E. Roman   if (V) {
181986c4d72SJose E. Roman     ierr = PetscObjectTypeCompareAny((PetscObject)V,&match,MATSEQDENSE,MATMPIDENSE,"");CHKERRQ(ierr);
182986c4d72SJose E. Roman     if (!match) SETERRQ(PetscObjectComm((PetscObject)V),PETSC_ERR_SUP,"Matrix V must be of type dense");
183986c4d72SJose E. Roman   }
1849d9032efSJose E. Roman 
1859d9032efSJose E. Roman   ierr = MatGetSize(U,NULL,&k);CHKERRQ(ierr);
1869d9032efSJose E. Roman   ierr = MatGetSize(V,NULL,&k1);CHKERRQ(ierr);
1879d9032efSJose E. Roman   if (k!=k1) SETERRQ(PetscObjectComm((PetscObject)U),PETSC_ERR_ARG_INCOMP,"U and V have different number of columns");
1889d9032efSJose E. Roman   ierr = MatGetLocalSize(U,&m,NULL);CHKERRQ(ierr);
1899d9032efSJose E. Roman   ierr = MatGetLocalSize(V,&n,NULL);CHKERRQ(ierr);
19086eed97dSJose E. Roman   if (A) {
1919d9032efSJose E. Roman     ierr = MatGetLocalSize(A,&m1,&n1);CHKERRQ(ierr);
1929d9032efSJose E. Roman     if (m!=m1) SETERRQ(PetscObjectComm((PetscObject)U),PETSC_ERR_ARG_INCOMP,"Local dimensions of U and A do not match");
1939d9032efSJose E. Roman     if (n!=n1) SETERRQ(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_INCOMP,"Local dimensions of V and A do not match");
19486eed97dSJose E. Roman   }
195986c4d72SJose E. Roman   if (c) {
196986c4d72SJose E. Roman     ierr = VecGetSize(c,&k1);CHKERRQ(ierr);
197986c4d72SJose 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");
198986c4d72SJose E. Roman     ierr = VecGetLocalSize(c,&k1);CHKERRQ(ierr);
199986c4d72SJose E. Roman     if (k!=k1) SETERRQ(PetscObjectComm((PetscObject)c),PETSC_ERR_ARG_INCOMP,"c must be a sequential vector");
200986c4d72SJose E. Roman   }
2019d9032efSJose E. Roman 
20286eed97dSJose E. Roman   ierr = MatCreate(PetscObjectComm((PetscObject)U),N);CHKERRQ(ierr);
2039d9032efSJose E. Roman   ierr = MatSetSizes(*N,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
204ab92ecdeSBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)*N,MATLRC);CHKERRQ(ierr);
205ceb9aaf7SBarry Smith 
206b00a9115SJed Brown   ierr       = PetscNewLog(*N,&Na);CHKERRQ(ierr);
20738f2d2fdSLisandro Dalcin   (*N)->data = (void*)Na;
208ceb9aaf7SBarry Smith   Na->A      = A;
209986c4d72SJose E. Roman   Na->c      = c;
210ab92ecdeSBarry Smith 
211ab92ecdeSBarry Smith   ierr = MatDenseGetLocalMatrix(U,&Na->U);CHKERRQ(ierr);
212ab92ecdeSBarry Smith   ierr = MatDenseGetLocalMatrix(V,&Na->V);CHKERRQ(ierr);
21386eed97dSJose E. Roman   if (A) { ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); }
214ab92ecdeSBarry Smith   ierr = PetscObjectReference((PetscObject)Na->U);CHKERRQ(ierr);
215ab92ecdeSBarry Smith   ierr = PetscObjectReference((PetscObject)Na->V);CHKERRQ(ierr);
216986c4d72SJose E. Roman   if (c) { ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr); }
217ceb9aaf7SBarry Smith 
218d0f46423SBarry Smith   ierr = VecCreateSeq(PETSC_COMM_SELF,U->cmap->N,&Na->work1);CHKERRQ(ierr);
219ab92ecdeSBarry Smith   ierr = VecDuplicate(Na->work1,&Na->work2);CHKERRQ(ierr);
220986c4d72SJose E. Roman   ierr = PetscMPIIntCast(U->cmap->N,&Na->nwork);CHKERRQ(ierr);
221ab92ecdeSBarry Smith 
2220575051bSJose E. Roman   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,V->rmap->n,NULL,&Na->xl);CHKERRQ(ierr);
2230575051bSJose E. Roman   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,U->rmap->n,NULL,&Na->yl);CHKERRQ(ierr);
2240575051bSJose E. Roman 
225ab92ecdeSBarry Smith   (*N)->ops->destroy = MatDestroy_LRC;
226ab92ecdeSBarry Smith   (*N)->ops->mult    = MatMult_LRC;
227ceb9aaf7SBarry Smith   (*N)->assembled    = PETSC_TRUE;
22826c5da1bSJose E. Roman   (*N)->preallocated = PETSC_TRUE;
2299d9032efSJose E. Roman   (*N)->cmap->N      = V->rmap->N;
2309d9032efSJose E. Roman   (*N)->rmap->N      = U->rmap->N;
2319d9032efSJose E. Roman   (*N)->cmap->n      = V->rmap->n;
2329d9032efSJose E. Roman   (*N)->rmap->n      = U->rmap->n;
233*267ca84cSJose E. Roman 
234*267ca84cSJose E. Roman   ierr = PetscObjectComposeFunction((PetscObject)(*N),"MatLRCGetMats_C",MatLRCGetMats_LRC);CHKERRQ(ierr);
235ceb9aaf7SBarry Smith   PetscFunctionReturn(0);
236ceb9aaf7SBarry Smith }
237ceb9aaf7SBarry Smith 
238