xref: /petsc/src/mat/impls/lrc/lrc.c (revision 2bf68e3e0f2a61f71e7c65bee250bfa1c8ce0cdb)
1 
2 #include <petsc/private/matimpl.h>          /*I "petscmat.h" I*/
3 
4 typedef struct {
5   Mat         A;           /* sparse matrix */
6   Mat         U,V;         /* dense tall-skinny matrices */
7   Vec         c;           /* sequential vector containing the diagonal of C */
8   Vec         work1,work2; /* sequential vectors that hold partial products */
9   PetscMPIInt nwork;       /* length of work vectors */
10   Vec         xl,yl;       /* auxiliary sequential vectors for matmult operation */
11 } Mat_LRC;
12 
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   const PetscScalar *a;
22 
23   PetscFunctionBegin;
24   ierr = VecGetArrayRead(x,&a);CHKERRQ(ierr);
25   ierr = VecPlaceArray(Na->xl,a);CHKERRQ(ierr);
26   ierr = VecGetLocalVector(y,Na->yl);CHKERRQ(ierr);
27 
28   /* multiply the local part of V with the local part of x */
29 #if defined(PETSC_USE_COMPLEX)
30   ierr = MatMultHermitianTranspose(Na->V,Na->xl,Na->work1);CHKERRQ(ierr);
31 #else
32   ierr = MatMultTranspose(Na->V,Na->xl,Na->work1);CHKERRQ(ierr);
33 #endif
34 
35   /* form the sum of all the local multiplies: this is work2 = V'*x =
36      sum_{all processors} work1 */
37   ierr = VecGetArray(Na->work1,&w1);CHKERRQ(ierr);
38   ierr = VecGetArray(Na->work2,&w2);CHKERRQ(ierr);
39   ierr = MPIU_Allreduce(w1,w2,Na->nwork,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)N));CHKERRQ(ierr);
40   ierr = VecRestoreArray(Na->work1,&w1);CHKERRQ(ierr);
41   ierr = VecRestoreArray(Na->work2,&w2);CHKERRQ(ierr);
42 
43   if (Na->c) {  /* work2 = C*work2 */
44     ierr = VecPointwiseMult(Na->work2,Na->c,Na->work2);CHKERRQ(ierr);
45   }
46 
47   if (Na->A) {
48     /* form y = A*x */
49     ierr = MatMult(Na->A,x,y);CHKERRQ(ierr);
50     /* multiply-add y = y + U*work2 */
51     ierr = MatMultAdd(Na->U,Na->work2,Na->yl,Na->yl);CHKERRQ(ierr);
52   } else {
53     /* multiply y = U*work2 */
54     ierr = MatMult(Na->U,Na->work2,Na->yl);CHKERRQ(ierr);
55   }
56 
57   ierr = VecRestoreArrayRead(x,&a);CHKERRQ(ierr);
58   ierr = VecResetArray(Na->xl);CHKERRQ(ierr);
59   ierr = VecRestoreLocalVector(y,Na->yl);CHKERRQ(ierr);
60   PetscFunctionReturn(0);
61 }
62 
63 #undef __FUNCT__
64 #define __FUNCT__ "MatDestroy_LRC"
65 PetscErrorCode MatDestroy_LRC(Mat N)
66 {
67   Mat_LRC        *Na = (Mat_LRC*)N->data;
68   PetscErrorCode ierr;
69 
70   PetscFunctionBegin;
71   ierr = MatDestroy(&Na->A);CHKERRQ(ierr);
72   ierr = MatDestroy(&Na->U);CHKERRQ(ierr);
73   ierr = MatDestroy(&Na->V);CHKERRQ(ierr);
74   ierr = VecDestroy(&Na->c);CHKERRQ(ierr);
75   ierr = VecDestroy(&Na->work1);CHKERRQ(ierr);
76   ierr = VecDestroy(&Na->work2);CHKERRQ(ierr);
77   ierr = VecDestroy(&Na->xl);CHKERRQ(ierr);
78   ierr = VecDestroy(&Na->yl);CHKERRQ(ierr);
79   ierr = PetscFree(N->data);CHKERRQ(ierr);
80   ierr = PetscObjectComposeFunction((PetscObject)N,"MatLRCGetMats_C",0);CHKERRQ(ierr);
81   PetscFunctionReturn(0);
82 }
83 
84 #undef __FUNCT__
85 #define __FUNCT__ "MatLRCGetMats_LRC"
86 PetscErrorCode MatLRCGetMats_LRC(Mat N,Mat *A,Mat *U,Vec *c,Mat *V)
87 {
88   Mat_LRC *Na = (Mat_LRC*)N->data;
89 
90   PetscFunctionBegin;
91   if (A) *A = Na->A;
92   if (U) *U = Na->U;
93   if (c) *c = Na->c;
94   if (V) *V = Na->V;
95   PetscFunctionReturn(0);
96 }
97 
98 #undef __FUNCT__
99 #define __FUNCT__ "MatLRCGetMats"
100 /*@
101    MatLRCGetMats - Returns the constituents of an LRC matrix
102 
103    Collective on Mat
104 
105    Input Parameter:
106 .  N    - matrix of type LRC
107 
108    Output Parameters:
109 +  A    - the (sparse) matrix
110 .  U, V - two dense rectangular (tall and skinny) matrices
111 -  c    - a sequential vector containing the diagonal of C
112 
113    Note:
114    The returned matrices need not be destroyed by the caller.
115 
116    Level: intermediate
117 
118 .seealso: MatCreateLRC()
119 @*/
120 PetscErrorCode  MatLRCGetMats(Mat N,Mat *A,Mat *U,Vec *c,Mat *V)
121 {
122   PetscErrorCode ierr;
123 
124   PetscFunctionBegin;
125   ierr = PetscUseMethod(N,"MatLRCGetMats_C",(Mat,Mat*,Mat*,Vec*,Mat*),(N,A,U,c,V));CHKERRQ(ierr);
126   PetscFunctionReturn(0);
127 }
128 
129 #undef __FUNCT__
130 #define __FUNCT__ "MatCreateLRC"
131 /*@
132    MatCreateLRC - Creates a new matrix object that behaves like A + U*C*V'
133 
134    Collective on Mat
135 
136    Input Parameters:
137 +  A    - the (sparse) matrix (can be NULL)
138 .  U, V - two dense rectangular (tall and skinny) matrices
139 -  c    - a sequential vector containing the diagonal of C (can be NULL)
140 
141    Output Parameter:
142 .  N    - the matrix that represents A + U*C*V'
143 
144    Notes:
145    The matrix A + U*C*V' is not formed! Rather the new matrix
146    object performs the matrix-vector product by first multiplying by
147    A and then adding the other term.
148 
149    C is a diagonal matrix (represented as a vector) of order k,
150    where k is the number of columns of both U and V.
151 
152    If A is NULL then the new object behaves like a low-rank matrix U*C*V'.
153 
154    Use V=U (or V=NULL) for a symmetric low-rank correction, A + U*C*U'.
155 
156    If c is NULL then the low-rank correction is just U*V'.
157 
158    Level: intermediate
159 
160 .seealso: MatLRCGetMats()
161 @*/
162 PetscErrorCode MatCreateLRC(Mat A,Mat U,Vec c,Mat V,Mat *N)
163 {
164   PetscErrorCode ierr;
165   PetscBool      match;
166   PetscInt       m,n,k,m1,n1,k1;
167   Mat_LRC        *Na;
168 
169   PetscFunctionBegin;
170   if (A) PetscValidHeaderSpecific(A,MAT_CLASSID,1);
171   PetscValidHeaderSpecific(U,MAT_CLASSID,2);
172   if (c) PetscValidHeaderSpecific(c,VEC_CLASSID,3);
173   if (V) PetscValidHeaderSpecific(V,MAT_CLASSID,4);
174   else V=U;
175   if (A) PetscCheckSameComm(A,1,U,2);
176   PetscCheckSameComm(U,2,V,4);
177 
178   ierr = PetscObjectTypeCompareAny((PetscObject)U,&match,MATSEQDENSE,MATMPIDENSE,"");CHKERRQ(ierr);
179   if (!match) SETERRQ(PetscObjectComm((PetscObject)U),PETSC_ERR_SUP,"Matrix U must be of type dense");
180   if (V) {
181     ierr = PetscObjectTypeCompareAny((PetscObject)V,&match,MATSEQDENSE,MATMPIDENSE,"");CHKERRQ(ierr);
182     if (!match) SETERRQ(PetscObjectComm((PetscObject)V),PETSC_ERR_SUP,"Matrix V must be of type dense");
183   }
184 
185   ierr = MatGetSize(U,NULL,&k);CHKERRQ(ierr);
186   ierr = MatGetSize(V,NULL,&k1);CHKERRQ(ierr);
187   if (k!=k1) SETERRQ2(PetscObjectComm((PetscObject)U),PETSC_ERR_ARG_INCOMP,"U and V have different number of columns (%D vs %D)",k,k1);
188   ierr = MatGetLocalSize(U,&m,NULL);CHKERRQ(ierr);
189   ierr = MatGetLocalSize(V,&n,NULL);CHKERRQ(ierr);
190   if (A) {
191     ierr = MatGetLocalSize(A,&m1,&n1);CHKERRQ(ierr);
192     if (m!=m1) SETERRQ2(PetscObjectComm((PetscObject)U),PETSC_ERR_ARG_INCOMP,"Local dimensions of U %D and A %D do not match",m,m1);
193     if (n!=n1) SETERRQ2(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_INCOMP,"Local dimensions of V %D and A %D do not match",n,n1);
194   }
195   if (c) {
196     ierr = VecGetSize(c,&k1);CHKERRQ(ierr);
197     if (k!=k1) SETERRQ2(PetscObjectComm((PetscObject)c),PETSC_ERR_ARG_INCOMP,"The length of c %D does not match the number of columns of U and V (%D)",k1,k);
198     ierr = VecGetLocalSize(c,&k1);CHKERRQ(ierr);
199     if (k!=k1) SETERRQ(PetscObjectComm((PetscObject)c),PETSC_ERR_ARG_INCOMP,"c must be a sequential vector");
200   }
201 
202   ierr = MatCreate(PetscObjectComm((PetscObject)U),N);CHKERRQ(ierr);
203   ierr = MatSetSizes(*N,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
204   ierr = PetscObjectChangeTypeName((PetscObject)*N,MATLRC);CHKERRQ(ierr);
205 
206   ierr       = PetscNewLog(*N,&Na);CHKERRQ(ierr);
207   (*N)->data = (void*)Na;
208   Na->A      = A;
209   Na->c      = c;
210 
211   ierr = MatDenseGetLocalMatrix(U,&Na->U);CHKERRQ(ierr);
212   ierr = MatDenseGetLocalMatrix(V,&Na->V);CHKERRQ(ierr);
213   if (A) { ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); }
214   ierr = PetscObjectReference((PetscObject)Na->U);CHKERRQ(ierr);
215   ierr = PetscObjectReference((PetscObject)Na->V);CHKERRQ(ierr);
216   if (c) { ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr); }
217 
218   ierr = VecCreateSeq(PETSC_COMM_SELF,U->cmap->N,&Na->work1);CHKERRQ(ierr);
219   ierr = VecDuplicate(Na->work1,&Na->work2);CHKERRQ(ierr);
220   ierr = PetscMPIIntCast(U->cmap->N,&Na->nwork);CHKERRQ(ierr);
221 
222   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,V->rmap->n,NULL,&Na->xl);CHKERRQ(ierr);
223   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,U->rmap->n,NULL,&Na->yl);CHKERRQ(ierr);
224 
225   (*N)->ops->destroy = MatDestroy_LRC;
226   (*N)->ops->mult    = MatMult_LRC;
227   (*N)->assembled    = PETSC_TRUE;
228   (*N)->preallocated = PETSC_TRUE;
229   (*N)->cmap->N      = V->rmap->N;
230   (*N)->rmap->N      = U->rmap->N;
231   (*N)->cmap->n      = V->rmap->n;
232   (*N)->rmap->n      = U->rmap->n;
233 
234   ierr = PetscObjectComposeFunction((PetscObject)(*N),"MatLRCGetMats_C",MatLRCGetMats_LRC);CHKERRQ(ierr);
235   PetscFunctionReturn(0);
236 }
237 
238