Line data Source code
1 : /*
2 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3 : SLEPc - Scalable Library for Eigenvalue Problem Computations
4 : Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
5 :
6 : This file is part of SLEPc.
7 : SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9 : */
10 : /*
11 : SLEPc singular value solver: "lanczos"
12 :
13 : Method: Explicitly restarted Lanczos
14 :
15 : Algorithm:
16 :
17 : Golub-Kahan-Lanczos bidiagonalization with explicit restart.
18 :
19 : References:
20 :
21 : [1] G.H. Golub and W. Kahan, "Calculating the singular values
22 : and pseudo-inverse of a matrix", SIAM J. Numer. Anal. Ser.
23 : B 2:205-224, 1965.
24 :
25 : [2] V. Hernandez, J.E. Roman, and A. Tomas, "A robust and
26 : efficient parallel SVD solver based on restarted Lanczos
27 : bidiagonalization", Elec. Trans. Numer. Anal. 31:68-85,
28 : 2008.
29 : */
30 :
31 : #include <slepc/private/svdimpl.h> /*I "slepcsvd.h" I*/
32 :
33 : typedef struct {
34 : PetscBool oneside;
35 : } SVD_LANCZOS;
36 :
37 18 : static PetscErrorCode SVDSetUp_Lanczos(SVD svd)
38 : {
39 18 : SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;
40 18 : PetscInt N;
41 :
42 18 : PetscFunctionBegin;
43 18 : SVDCheckStandard(svd);
44 18 : SVDCheckDefinite(svd);
45 18 : PetscCall(MatGetSize(svd->A,NULL,&N));
46 18 : PetscCall(SVDSetDimensions_Default(svd));
47 18 : PetscCheck(svd->ncv<=svd->nsv+svd->mpd,PetscObjectComm((PetscObject)svd),PETSC_ERR_USER_INPUT,"The value of ncv must not be larger than nev+mpd");
48 18 : if (svd->max_it==PETSC_DEFAULT) svd->max_it = PetscMax(N/svd->ncv,100);
49 18 : svd->leftbasis = PetscNot(lanczos->oneside);
50 18 : PetscCall(SVDAllocateSolution(svd,1));
51 18 : PetscCall(DSSetType(svd->ds,DSSVD));
52 18 : PetscCall(DSSetCompact(svd->ds,PETSC_TRUE));
53 18 : PetscCall(DSSetExtraRow(svd->ds,PETSC_TRUE));
54 18 : PetscCall(DSAllocate(svd->ds,svd->ncv+1));
55 18 : PetscFunctionReturn(PETSC_SUCCESS);
56 : }
57 :
58 785 : PetscErrorCode SVDTwoSideLanczos(SVD svd,PetscReal *alpha,PetscReal *beta,BV V,BV U,PetscInt k,PetscInt *n,PetscBool *breakdown)
59 : {
60 785 : PetscInt i;
61 785 : Vec u,v;
62 785 : PetscBool lindep=PETSC_FALSE;
63 :
64 785 : PetscFunctionBegin;
65 785 : PetscCall(BVGetColumn(svd->V,k,&v));
66 785 : PetscCall(BVGetColumn(svd->U,k,&u));
67 785 : PetscCall(MatMult(svd->A,v,u));
68 785 : PetscCall(BVRestoreColumn(svd->V,k,&v));
69 785 : PetscCall(BVRestoreColumn(svd->U,k,&u));
70 785 : PetscCall(BVOrthonormalizeColumn(svd->U,k,PETSC_FALSE,alpha+k,&lindep));
71 785 : if (PetscUnlikely(lindep)) {
72 0 : *n = k;
73 0 : if (breakdown) *breakdown = lindep;
74 0 : PetscFunctionReturn(PETSC_SUCCESS);
75 : }
76 :
77 6356 : for (i=k+1;i<*n;i++) {
78 5571 : PetscCall(BVGetColumn(svd->V,i,&v));
79 5571 : PetscCall(BVGetColumn(svd->U,i-1,&u));
80 5571 : PetscCall(MatMult(svd->AT,u,v));
81 5571 : PetscCall(BVRestoreColumn(svd->V,i,&v));
82 5571 : PetscCall(BVRestoreColumn(svd->U,i-1,&u));
83 5571 : PetscCall(BVOrthonormalizeColumn(svd->V,i,PETSC_FALSE,beta+i-1,&lindep));
84 5571 : if (PetscUnlikely(lindep)) {
85 0 : *n = i;
86 0 : break;
87 : }
88 5571 : PetscCall(BVGetColumn(svd->V,i,&v));
89 5571 : PetscCall(BVGetColumn(svd->U,i,&u));
90 5571 : PetscCall(MatMult(svd->A,v,u));
91 5571 : PetscCall(BVRestoreColumn(svd->V,i,&v));
92 5571 : PetscCall(BVRestoreColumn(svd->U,i,&u));
93 5571 : PetscCall(BVOrthonormalizeColumn(svd->U,i,PETSC_FALSE,alpha+i,&lindep));
94 5571 : if (PetscUnlikely(lindep)) {
95 0 : *n = i;
96 0 : break;
97 : }
98 : }
99 :
100 785 : if (!lindep) {
101 785 : PetscCall(BVGetColumn(svd->V,*n,&v));
102 785 : PetscCall(BVGetColumn(svd->U,*n-1,&u));
103 785 : PetscCall(MatMult(svd->AT,u,v));
104 785 : PetscCall(BVRestoreColumn(svd->V,*n,&v));
105 785 : PetscCall(BVRestoreColumn(svd->U,*n-1,&u));
106 785 : PetscCall(BVOrthogonalizeColumn(svd->V,*n,NULL,beta+*n-1,&lindep));
107 : }
108 785 : if (breakdown) *breakdown = lindep;
109 785 : PetscFunctionReturn(PETSC_SUCCESS);
110 : }
111 :
112 54 : static PetscErrorCode SVDOneSideLanczos(SVD svd,PetscReal *alpha,PetscReal *beta,BV V,Vec u,Vec u_1,PetscInt k,PetscInt n,PetscScalar* work)
113 : {
114 54 : PetscInt i,bvl,bvk;
115 54 : PetscReal a,b;
116 54 : Vec z,temp;
117 :
118 54 : PetscFunctionBegin;
119 54 : PetscCall(BVGetActiveColumns(V,&bvl,&bvk));
120 54 : PetscCall(BVGetColumn(V,k,&z));
121 54 : PetscCall(MatMult(svd->A,z,u));
122 54 : PetscCall(BVRestoreColumn(V,k,&z));
123 :
124 468 : for (i=k+1;i<n;i++) {
125 414 : PetscCall(BVGetColumn(V,i,&z));
126 414 : PetscCall(MatMult(svd->AT,u,z));
127 414 : PetscCall(BVRestoreColumn(V,i,&z));
128 414 : PetscCall(VecNormBegin(u,NORM_2,&a));
129 414 : PetscCall(BVSetActiveColumns(V,0,i));
130 414 : PetscCall(BVDotColumnBegin(V,i,work));
131 414 : PetscCall(VecNormEnd(u,NORM_2,&a));
132 414 : PetscCall(BVDotColumnEnd(V,i,work));
133 414 : PetscCall(VecScale(u,1.0/a));
134 414 : PetscCall(BVMultColumn(V,-1.0/a,1.0/a,i,work));
135 :
136 : /* h = V^* z, z = z - V h */
137 414 : PetscCall(BVDotColumn(V,i,work));
138 414 : PetscCall(BVMultColumn(V,-1.0,1.0,i,work));
139 414 : PetscCall(BVNormColumn(V,i,NORM_2,&b));
140 414 : PetscCheck(PetscAbsReal(b)>10*PETSC_MACHINE_EPSILON,PetscObjectComm((PetscObject)svd),PETSC_ERR_PLIB,"Recurrence generated a zero vector; use a two-sided variant");
141 414 : PetscCall(BVScaleColumn(V,i,1.0/b));
142 :
143 414 : PetscCall(BVGetColumn(V,i,&z));
144 414 : PetscCall(MatMult(svd->A,z,u_1));
145 414 : PetscCall(BVRestoreColumn(V,i,&z));
146 414 : PetscCall(VecAXPY(u_1,-b,u));
147 414 : alpha[i-1] = a;
148 414 : beta[i-1] = b;
149 414 : temp = u;
150 414 : u = u_1;
151 414 : u_1 = temp;
152 : }
153 :
154 54 : PetscCall(BVGetColumn(V,n,&z));
155 54 : PetscCall(MatMult(svd->AT,u,z));
156 54 : PetscCall(BVRestoreColumn(V,n,&z));
157 54 : PetscCall(VecNormBegin(u,NORM_2,&a));
158 54 : PetscCall(BVDotColumnBegin(V,n,work));
159 54 : PetscCall(VecNormEnd(u,NORM_2,&a));
160 54 : PetscCall(BVDotColumnEnd(V,n,work));
161 54 : PetscCall(VecScale(u,1.0/a));
162 54 : PetscCall(BVMultColumn(V,-1.0/a,1.0/a,n,work));
163 :
164 : /* h = V^* z, z = z - V h */
165 54 : PetscCall(BVDotColumn(V,n,work));
166 54 : PetscCall(BVMultColumn(V,-1.0,1.0,n,work));
167 54 : PetscCall(BVNormColumn(V,i,NORM_2,&b));
168 :
169 54 : alpha[n-1] = a;
170 54 : beta[n-1] = b;
171 54 : PetscCall(BVSetActiveColumns(V,bvl,bvk));
172 54 : PetscFunctionReturn(PETSC_SUCCESS);
173 : }
174 :
175 : /*
176 : SVDKrylovConvergence - Implements the loop that checks for convergence
177 : in Krylov methods.
178 :
179 : Input Parameters:
180 : svd - the solver; some error estimates are updated in svd->errest
181 : getall - whether all residuals must be computed
182 : kini - initial value of k (the loop variable)
183 : nits - number of iterations of the loop
184 : normr - norm of triangular factor of qr([A;B]), used only in GSVD
185 :
186 : Output Parameter:
187 : kout - the first index where the convergence test failed
188 : */
189 1689 : PetscErrorCode SVDKrylovConvergence(SVD svd,PetscBool getall,PetscInt kini,PetscInt nits,PetscReal normr,PetscInt *kout)
190 : {
191 1689 : PetscInt k,marker,ld;
192 1689 : PetscReal *alpha,*beta,*betah,resnorm;
193 1689 : PetscBool extra;
194 :
195 1689 : PetscFunctionBegin;
196 1689 : if (PetscUnlikely(svd->conv == SVD_CONV_MAXIT && svd->its >= svd->max_it)) *kout = svd->nsv;
197 : else {
198 1689 : PetscCall(DSGetLeadingDimension(svd->ds,&ld));
199 1689 : PetscCall(DSGetExtraRow(svd->ds,&extra));
200 1689 : PetscCheck(extra,PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"Only implemented for DS with extra row");
201 1689 : marker = -1;
202 1689 : if (svd->trackall) getall = PETSC_TRUE;
203 1689 : PetscCall(DSGetArrayReal(svd->ds,DS_MAT_T,&alpha));
204 1689 : beta = alpha + ld;
205 1689 : betah = alpha + 2*ld;
206 2338 : for (k=kini;k<kini+nits;k++) {
207 2332 : if (svd->isgeneralized) resnorm = SlepcAbs(beta[k],betah[k])*normr;
208 1365 : else resnorm = PetscAbsReal(beta[k]);
209 2332 : PetscCall((*svd->converged)(svd,svd->sigma[k],resnorm,&svd->errest[k],svd->convergedctx));
210 2332 : if (marker==-1 && svd->errest[k] >= svd->tol) marker = k;
211 2332 : if (marker!=-1 && !getall) break;
212 : }
213 1689 : PetscCall(DSRestoreArrayReal(svd->ds,DS_MAT_T,&alpha));
214 1689 : if (marker!=-1) k = marker;
215 1689 : *kout = k;
216 : }
217 1689 : PetscFunctionReturn(PETSC_SUCCESS);
218 : }
219 :
220 18 : static PetscErrorCode SVDSolve_Lanczos(SVD svd)
221 : {
222 18 : SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;
223 18 : PetscReal *alpha,*beta;
224 18 : PetscScalar *swork,*w,*P;
225 18 : PetscInt i,k,j,nv,ld;
226 18 : Vec u=NULL,u_1=NULL;
227 18 : Mat U,V;
228 :
229 18 : PetscFunctionBegin;
230 : /* allocate working space */
231 18 : PetscCall(DSGetLeadingDimension(svd->ds,&ld));
232 18 : PetscCall(PetscMalloc2(ld,&w,svd->ncv,&swork));
233 :
234 18 : if (lanczos->oneside) {
235 2 : PetscCall(MatCreateVecs(svd->A,NULL,&u));
236 2 : PetscCall(MatCreateVecs(svd->A,NULL,&u_1));
237 : }
238 :
239 : /* normalize start vector */
240 18 : if (!svd->nini) {
241 14 : PetscCall(BVSetRandomColumn(svd->V,0));
242 14 : PetscCall(BVOrthonormalizeColumn(svd->V,0,PETSC_TRUE,NULL,NULL));
243 : }
244 :
245 589 : while (svd->reason == SVD_CONVERGED_ITERATING) {
246 571 : svd->its++;
247 :
248 : /* inner loop */
249 571 : nv = PetscMin(svd->nconv+svd->mpd,svd->ncv);
250 571 : PetscCall(BVSetActiveColumns(svd->V,svd->nconv,nv));
251 571 : PetscCall(DSGetArrayReal(svd->ds,DS_MAT_T,&alpha));
252 571 : beta = alpha + ld;
253 571 : if (lanczos->oneside) PetscCall(SVDOneSideLanczos(svd,alpha,beta,svd->V,u,u_1,svd->nconv,nv,swork));
254 : else {
255 517 : PetscCall(BVSetActiveColumns(svd->U,svd->nconv,nv));
256 517 : PetscCall(SVDTwoSideLanczos(svd,alpha,beta,svd->V,svd->U,svd->nconv,&nv,NULL));
257 : }
258 571 : PetscCall(DSRestoreArrayReal(svd->ds,DS_MAT_T,&alpha));
259 :
260 : /* compute SVD of bidiagonal matrix */
261 571 : PetscCall(DSSetDimensions(svd->ds,nv,svd->nconv,0));
262 571 : PetscCall(DSSVDSetDimensions(svd->ds,nv));
263 571 : PetscCall(DSSetState(svd->ds,DS_STATE_INTERMEDIATE));
264 571 : PetscCall(DSSolve(svd->ds,w,NULL));
265 571 : PetscCall(DSSort(svd->ds,w,NULL,NULL,NULL,NULL));
266 571 : PetscCall(DSUpdateExtraRow(svd->ds));
267 571 : PetscCall(DSSynchronize(svd->ds,w,NULL));
268 5920 : for (i=svd->nconv;i<nv;i++) svd->sigma[i] = PetscRealPart(w[i]);
269 :
270 : /* check convergence */
271 571 : PetscCall(SVDKrylovConvergence(svd,PETSC_FALSE,svd->nconv,nv-svd->nconv,1.0,&k));
272 571 : PetscCall((*svd->stopping)(svd,svd->its,svd->max_it,k,svd->nsv,&svd->reason,svd->stoppingctx));
273 :
274 : /* compute restart vector */
275 571 : if (svd->reason == SVD_CONVERGED_ITERATING) {
276 553 : if (k<nv) {
277 553 : PetscCall(DSGetArray(svd->ds,DS_MAT_V,&P));
278 5746 : for (j=svd->nconv;j<nv;j++) swork[j-svd->nconv] = PetscConj(P[j+k*ld]);
279 553 : PetscCall(DSRestoreArray(svd->ds,DS_MAT_V,&P));
280 553 : PetscCall(BVMultColumn(svd->V,1.0,0.0,nv,swork));
281 : } else {
282 : /* all approximations have converged, generate a new initial vector */
283 0 : PetscCall(BVSetRandomColumn(svd->V,nv));
284 0 : PetscCall(BVOrthonormalizeColumn(svd->V,nv,PETSC_FALSE,NULL,NULL));
285 : }
286 : }
287 :
288 : /* compute converged singular vectors */
289 571 : PetscCall(DSGetMat(svd->ds,DS_MAT_V,&V));
290 571 : PetscCall(BVMultInPlace(svd->V,V,svd->nconv,k));
291 571 : PetscCall(DSRestoreMat(svd->ds,DS_MAT_V,&V));
292 571 : if (!lanczos->oneside) {
293 517 : PetscCall(DSGetMat(svd->ds,DS_MAT_U,&U));
294 517 : PetscCall(BVMultInPlace(svd->U,U,svd->nconv,k));
295 517 : PetscCall(DSRestoreMat(svd->ds,DS_MAT_U,&U));
296 : }
297 :
298 : /* copy restart vector from the last column */
299 571 : if (svd->reason == SVD_CONVERGED_ITERATING) PetscCall(BVCopyColumn(svd->V,nv,k));
300 :
301 571 : svd->nconv = k;
302 589 : PetscCall(SVDMonitor(svd,svd->its,svd->nconv,svd->sigma,svd->errest,nv));
303 : }
304 :
305 : /* free working space */
306 18 : PetscCall(VecDestroy(&u));
307 18 : PetscCall(VecDestroy(&u_1));
308 18 : PetscCall(PetscFree2(w,swork));
309 18 : PetscFunctionReturn(PETSC_SUCCESS);
310 : }
311 :
312 12 : static PetscErrorCode SVDSetFromOptions_Lanczos(SVD svd,PetscOptionItems *PetscOptionsObject)
313 : {
314 12 : PetscBool set,val;
315 12 : SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;
316 :
317 12 : PetscFunctionBegin;
318 12 : PetscOptionsHeadBegin(PetscOptionsObject,"SVD Lanczos Options");
319 :
320 12 : PetscCall(PetscOptionsBool("-svd_lanczos_oneside","Use one-side reorthogonalization","SVDLanczosSetOneSide",lanczos->oneside,&val,&set));
321 12 : if (set) PetscCall(SVDLanczosSetOneSide(svd,val));
322 :
323 12 : PetscOptionsHeadEnd();
324 12 : PetscFunctionReturn(PETSC_SUCCESS);
325 : }
326 :
327 2 : static PetscErrorCode SVDLanczosSetOneSide_Lanczos(SVD svd,PetscBool oneside)
328 : {
329 2 : SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;
330 :
331 2 : PetscFunctionBegin;
332 2 : if (lanczos->oneside != oneside) {
333 2 : lanczos->oneside = oneside;
334 2 : svd->state = SVD_STATE_INITIAL;
335 : }
336 2 : PetscFunctionReturn(PETSC_SUCCESS);
337 : }
338 :
339 : /*@
340 : SVDLanczosSetOneSide - Indicate if the variant of the Lanczos method
341 : to be used is one-sided or two-sided.
342 :
343 : Logically Collective
344 :
345 : Input Parameters:
346 : + svd - singular value solver
347 : - oneside - boolean flag indicating if the method is one-sided or not
348 :
349 : Options Database Key:
350 : . -svd_lanczos_oneside <boolean> - Indicates the boolean flag
351 :
352 : Note:
353 : By default, a two-sided variant is selected, which is sometimes slightly
354 : more robust. However, the one-sided variant is faster because it avoids
355 : the orthogonalization associated to left singular vectors. It also saves
356 : the memory required for storing such vectors.
357 :
358 : Level: advanced
359 :
360 : .seealso: SVDTRLanczosSetOneSide()
361 : @*/
362 2 : PetscErrorCode SVDLanczosSetOneSide(SVD svd,PetscBool oneside)
363 : {
364 2 : PetscFunctionBegin;
365 2 : PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
366 8 : PetscValidLogicalCollectiveBool(svd,oneside,2);
367 2 : PetscTryMethod(svd,"SVDLanczosSetOneSide_C",(SVD,PetscBool),(svd,oneside));
368 2 : PetscFunctionReturn(PETSC_SUCCESS);
369 : }
370 :
371 2 : static PetscErrorCode SVDLanczosGetOneSide_Lanczos(SVD svd,PetscBool *oneside)
372 : {
373 2 : SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;
374 :
375 2 : PetscFunctionBegin;
376 2 : *oneside = lanczos->oneside;
377 2 : PetscFunctionReturn(PETSC_SUCCESS);
378 : }
379 :
380 : /*@
381 : SVDLanczosGetOneSide - Gets if the variant of the Lanczos method
382 : to be used is one-sided or two-sided.
383 :
384 : Not Collective
385 :
386 : Input Parameters:
387 : . svd - singular value solver
388 :
389 : Output Parameters:
390 : . oneside - boolean flag indicating if the method is one-sided or not
391 :
392 : Level: advanced
393 :
394 : .seealso: SVDLanczosSetOneSide()
395 : @*/
396 2 : PetscErrorCode SVDLanczosGetOneSide(SVD svd,PetscBool *oneside)
397 : {
398 2 : PetscFunctionBegin;
399 2 : PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
400 2 : PetscAssertPointer(oneside,2);
401 2 : PetscUseMethod(svd,"SVDLanczosGetOneSide_C",(SVD,PetscBool*),(svd,oneside));
402 2 : PetscFunctionReturn(PETSC_SUCCESS);
403 : }
404 :
405 14 : static PetscErrorCode SVDDestroy_Lanczos(SVD svd)
406 : {
407 14 : PetscFunctionBegin;
408 14 : PetscCall(PetscFree(svd->data));
409 14 : PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDLanczosSetOneSide_C",NULL));
410 14 : PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDLanczosGetOneSide_C",NULL));
411 14 : PetscFunctionReturn(PETSC_SUCCESS);
412 : }
413 :
414 0 : static PetscErrorCode SVDView_Lanczos(SVD svd,PetscViewer viewer)
415 : {
416 0 : SVD_LANCZOS *lanczos = (SVD_LANCZOS*)svd->data;
417 0 : PetscBool isascii;
418 :
419 0 : PetscFunctionBegin;
420 0 : PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
421 0 : if (isascii) PetscCall(PetscViewerASCIIPrintf(viewer," %s-sided reorthogonalization\n",lanczos->oneside? "one": "two"));
422 0 : PetscFunctionReturn(PETSC_SUCCESS);
423 : }
424 :
425 14 : SLEPC_EXTERN PetscErrorCode SVDCreate_Lanczos(SVD svd)
426 : {
427 14 : SVD_LANCZOS *ctx;
428 :
429 14 : PetscFunctionBegin;
430 14 : PetscCall(PetscNew(&ctx));
431 14 : svd->data = (void*)ctx;
432 :
433 14 : svd->ops->setup = SVDSetUp_Lanczos;
434 14 : svd->ops->solve = SVDSolve_Lanczos;
435 14 : svd->ops->destroy = SVDDestroy_Lanczos;
436 14 : svd->ops->setfromoptions = SVDSetFromOptions_Lanczos;
437 14 : svd->ops->view = SVDView_Lanczos;
438 14 : svd->ops->computevectors = SVDComputeVectors_Left;
439 14 : PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDLanczosSetOneSide_C",SVDLanczosSetOneSide_Lanczos));
440 14 : PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDLanczosGetOneSide_C",SVDLanczosGetOneSide_Lanczos));
441 14 : PetscFunctionReturn(PETSC_SUCCESS);
442 : }
|