xref: /petsc/src/mat/impls/lrc/lrc.c (revision dfd57a172ac9fa6c7b5fe6de6ab5df85cefc2996)
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 PetscErrorCode MatMult_LRC(Mat N,Vec x,Vec y)
14 {
15   Mat_LRC           *Na = (Mat_LRC*)N->data;
16   Mat               Uloc,Vloc;
17   PetscErrorCode    ierr;
18   PetscScalar       *w1,*w2;
19   const PetscScalar *a;
20 
21   PetscFunctionBegin;
22   ierr = VecGetArrayRead(x,&a);CHKERRQ(ierr);
23   ierr = VecPlaceArray(Na->xl,a);CHKERRQ(ierr);
24   ierr = VecGetLocalVector(y,Na->yl);CHKERRQ(ierr);
25   ierr = MatDenseGetLocalMatrix(Na->U,&Uloc);CHKERRQ(ierr);
26   ierr = MatDenseGetLocalMatrix(Na->V,&Vloc);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(Vloc,Na->xl,Na->work1);CHKERRQ(ierr);
31 #else
32   ierr = MatMultTranspose(Vloc,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));CHKERRMPI(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(Uloc,Na->work2,Na->yl,Na->yl);CHKERRQ(ierr);
52   } else {
53     /* multiply y = U*work2 */
54     ierr = MatMult(Uloc,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 PetscErrorCode MatDestroy_LRC(Mat N)
64 {
65   Mat_LRC        *Na = (Mat_LRC*)N->data;
66   PetscErrorCode ierr;
67 
68   PetscFunctionBegin;
69   ierr = MatDestroy(&Na->A);CHKERRQ(ierr);
70   ierr = MatDestroy(&Na->U);CHKERRQ(ierr);
71   ierr = MatDestroy(&Na->V);CHKERRQ(ierr);
72   ierr = VecDestroy(&Na->c);CHKERRQ(ierr);
73   ierr = VecDestroy(&Na->work1);CHKERRQ(ierr);
74   ierr = VecDestroy(&Na->work2);CHKERRQ(ierr);
75   ierr = VecDestroy(&Na->xl);CHKERRQ(ierr);
76   ierr = VecDestroy(&Na->yl);CHKERRQ(ierr);
77   ierr = PetscFree(N->data);CHKERRQ(ierr);
78   ierr = PetscObjectComposeFunction((PetscObject)N,"MatLRCGetMats_C",0);CHKERRQ(ierr);
79   PetscFunctionReturn(0);
80 }
81 
82 PetscErrorCode MatLRCGetMats_LRC(Mat N,Mat *A,Mat *U,Vec *c,Mat *V)
83 {
84   Mat_LRC *Na = (Mat_LRC*)N->data;
85 
86   PetscFunctionBegin;
87   if (A) *A = Na->A;
88   if (U) *U = Na->U;
89   if (c) *c = Na->c;
90   if (V) *V = Na->V;
91   PetscFunctionReturn(0);
92 }
93 
94 /*@
95    MatLRCGetMats - Returns the constituents of an LRC matrix
96 
97    Collective on Mat
98 
99    Input Parameter:
100 .  N    - matrix of type LRC
101 
102    Output Parameters:
103 +  A    - the (sparse) matrix
104 .  U, V - two dense rectangular (tall and skinny) matrices
105 -  c    - a sequential vector containing the diagonal of C
106 
107    Note:
108    The returned matrices need not be destroyed by the caller.
109 
110    Level: intermediate
111 
112 .seealso: MatCreateLRC()
113 @*/
114 PetscErrorCode  MatLRCGetMats(Mat N,Mat *A,Mat *U,Vec *c,Mat *V)
115 {
116   PetscErrorCode ierr;
117 
118   PetscFunctionBegin;
119   ierr = PetscUseMethod(N,"MatLRCGetMats_C",(Mat,Mat*,Mat*,Vec*,Mat*),(N,A,U,c,V));CHKERRQ(ierr);
120   PetscFunctionReturn(0);
121 }
122 
123 /*@
124    MatCreateLRC - Creates a new matrix object that behaves like A + U*C*V'
125 
126    Collective on Mat
127 
128    Input Parameters:
129 +  A    - the (sparse) matrix (can be NULL)
130 .  U, V - two dense rectangular (tall and skinny) matrices
131 -  c    - a sequential vector containing the diagonal of C (can be NULL)
132 
133    Output Parameter:
134 .  N    - the matrix that represents A + U*C*V'
135 
136    Notes:
137    The matrix A + U*C*V' is not formed! Rather the new matrix
138    object performs the matrix-vector product by first multiplying by
139    A and then adding the other term.
140 
141    C is a diagonal matrix (represented as a vector) of order k,
142    where k is the number of columns of both U and V.
143 
144    If A is NULL then the new object behaves like a low-rank matrix U*C*V'.
145 
146    Use V=U (or V=NULL) for a symmetric low-rank correction, A + U*C*U'.
147 
148    If c is NULL then the low-rank correction is just U*V'.
149 
150    Level: intermediate
151 
152 .seealso: MatLRCGetMats()
153 @*/
154 PetscErrorCode MatCreateLRC(Mat A,Mat U,Vec c,Mat V,Mat *N)
155 {
156   PetscErrorCode ierr;
157   PetscBool      match;
158   PetscInt       m,n,k,m1,n1,k1;
159   Mat_LRC        *Na;
160 
161   PetscFunctionBegin;
162   if (A) PetscValidHeaderSpecific(A,MAT_CLASSID,1);
163   PetscValidHeaderSpecific(U,MAT_CLASSID,2);
164   if (c) PetscValidHeaderSpecific(c,VEC_CLASSID,3);
165   if (V) PetscValidHeaderSpecific(V,MAT_CLASSID,4);
166   else V=U;
167   if (A) PetscCheckSameComm(A,1,U,2);
168   PetscCheckSameComm(U,2,V,4);
169 
170   ierr = PetscObjectTypeCompareAny((PetscObject)U,&match,MATSEQDENSE,MATMPIDENSE,"");CHKERRQ(ierr);
171   if (!match) SETERRQ(PetscObjectComm((PetscObject)U),PETSC_ERR_SUP,"Matrix U must be of type dense");
172   if (V) {
173     ierr = PetscObjectTypeCompareAny((PetscObject)V,&match,MATSEQDENSE,MATMPIDENSE,"");CHKERRQ(ierr);
174     if (!match) SETERRQ(PetscObjectComm((PetscObject)V),PETSC_ERR_SUP,"Matrix V must be of type dense");
175   }
176 
177   ierr = MatGetSize(U,NULL,&k);CHKERRQ(ierr);
178   ierr = MatGetSize(V,NULL,&k1);CHKERRQ(ierr);
179   if (k!=k1) SETERRQ2(PetscObjectComm((PetscObject)U),PETSC_ERR_ARG_INCOMP,"U and V have different number of columns (%D vs %D)",k,k1);
180   ierr = MatGetLocalSize(U,&m,NULL);CHKERRQ(ierr);
181   ierr = MatGetLocalSize(V,&n,NULL);CHKERRQ(ierr);
182   if (A) {
183     ierr = MatGetLocalSize(A,&m1,&n1);CHKERRQ(ierr);
184     if (m!=m1) SETERRQ2(PetscObjectComm((PetscObject)U),PETSC_ERR_ARG_INCOMP,"Local dimensions of U %D and A %D do not match",m,m1);
185     if (n!=n1) SETERRQ2(PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_INCOMP,"Local dimensions of V %D and A %D do not match",n,n1);
186   }
187   if (c) {
188     ierr = VecGetSize(c,&k1);CHKERRQ(ierr);
189     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);
190     ierr = VecGetLocalSize(c,&k1);CHKERRQ(ierr);
191     if (k!=k1) SETERRQ(PetscObjectComm((PetscObject)c),PETSC_ERR_ARG_INCOMP,"c must be a sequential vector");
192   }
193 
194   ierr = MatCreate(PetscObjectComm((PetscObject)U),N);CHKERRQ(ierr);
195   ierr = MatSetSizes(*N,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
196   ierr = PetscObjectChangeTypeName((PetscObject)*N,MATLRC);CHKERRQ(ierr);
197 
198   ierr       = PetscNewLog(*N,&Na);CHKERRQ(ierr);
199   (*N)->data = (void*)Na;
200   Na->A      = A;
201   Na->U      = U;
202   Na->c      = c;
203   Na->V      = V;
204 
205   if (A) { ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); }
206   ierr = PetscObjectReference((PetscObject)Na->U);CHKERRQ(ierr);
207   ierr = PetscObjectReference((PetscObject)Na->V);CHKERRQ(ierr);
208   if (c) { ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr); }
209 
210   ierr = VecCreateSeq(PETSC_COMM_SELF,U->cmap->N,&Na->work1);CHKERRQ(ierr);
211   ierr = VecDuplicate(Na->work1,&Na->work2);CHKERRQ(ierr);
212   ierr = PetscMPIIntCast(U->cmap->N,&Na->nwork);CHKERRQ(ierr);
213 
214   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,V->rmap->n,NULL,&Na->xl);CHKERRQ(ierr);
215   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,U->rmap->n,NULL,&Na->yl);CHKERRQ(ierr);
216 
217   (*N)->ops->destroy = MatDestroy_LRC;
218   (*N)->ops->mult    = MatMult_LRC;
219   (*N)->assembled    = PETSC_TRUE;
220   (*N)->preallocated = PETSC_TRUE;
221   (*N)->cmap->N      = V->rmap->N;
222   (*N)->rmap->N      = U->rmap->N;
223   (*N)->cmap->n      = V->rmap->n;
224   (*N)->rmap->n      = U->rmap->n;
225 
226   ierr = PetscObjectComposeFunction((PetscObject)(*N),"MatLRCGetMats_C",MatLRCGetMats_LRC);CHKERRQ(ierr);
227   PetscFunctionReturn(0);
228 }
229