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 : SVD routines for setting solver options
12 : */
13 :
14 : #include <slepc/private/svdimpl.h> /*I "slepcsvd.h" I*/
15 : #include <petscdraw.h>
16 :
17 : /*@
18 : SVDSetImplicitTranspose - Indicates how to handle the transpose of the matrix
19 : associated with the singular value problem.
20 :
21 : Logically Collective
22 :
23 : Input Parameters:
24 : + svd - the singular value solver context
25 : - impl - how to handle the transpose (implicitly or not)
26 :
27 : Options Database Key:
28 : . -svd_implicittranspose - Activate the implicit transpose mode.
29 :
30 : Notes:
31 : By default, the transpose of the matrix is explicitly built (if the matrix
32 : has defined the MatTranspose operation).
33 :
34 : If this flag is set to true, the solver does not build the transpose, but
35 : handles it implicitly via MatMultTranspose() (or MatMultHermitianTranspose()
36 : in the complex case) operations. This is likely to be more inefficient
37 : than the default behaviour, both in sequential and in parallel, but
38 : requires less storage.
39 :
40 : Level: advanced
41 :
42 : .seealso: SVDGetImplicitTranspose(), SVDSolve(), SVDSetOperators()
43 : @*/
44 11 : PetscErrorCode SVDSetImplicitTranspose(SVD svd,PetscBool impl)
45 : {
46 11 : PetscFunctionBegin;
47 11 : PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
48 33 : PetscValidLogicalCollectiveBool(svd,impl,2);
49 11 : if (svd->impltrans!=impl) {
50 11 : svd->impltrans = impl;
51 11 : svd->state = SVD_STATE_INITIAL;
52 : }
53 11 : PetscFunctionReturn(PETSC_SUCCESS);
54 : }
55 :
56 : /*@
57 : SVDGetImplicitTranspose - Gets the mode used to handle the transpose
58 : of the matrix associated with the singular value problem.
59 :
60 : Not Collective
61 :
62 : Input Parameter:
63 : . svd - the singular value solver context
64 :
65 : Output Parameter:
66 : . impl - how to handle the transpose (implicitly or not)
67 :
68 : Level: advanced
69 :
70 : .seealso: SVDSetImplicitTranspose(), SVDSolve(), SVDSetOperators()
71 : @*/
72 11 : PetscErrorCode SVDGetImplicitTranspose(SVD svd,PetscBool *impl)
73 : {
74 11 : PetscFunctionBegin;
75 11 : PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
76 11 : PetscAssertPointer(impl,2);
77 11 : *impl = svd->impltrans;
78 11 : PetscFunctionReturn(PETSC_SUCCESS);
79 : }
80 :
81 : /*@
82 : SVDSetTolerances - Sets the tolerance and maximum
83 : iteration count used by the default SVD convergence testers.
84 :
85 : Logically Collective
86 :
87 : Input Parameters:
88 : + svd - the singular value solver context
89 : . tol - the convergence tolerance
90 : - maxits - maximum number of iterations to use
91 :
92 : Options Database Keys:
93 : + -svd_tol <tol> - Sets the convergence tolerance
94 : - -svd_max_it <maxits> - Sets the maximum number of iterations allowed
95 :
96 : Note:
97 : Use PETSC_CURRENT to retain the current value of any of the parameters.
98 : Use PETSC_DETERMINE for either argument to assign a default value computed
99 : internally (may be different in each solver).
100 : For maxits use PETSC_UMLIMITED to indicate there is no upper bound on this value.
101 :
102 : Level: intermediate
103 :
104 : .seealso: SVDGetTolerances()
105 : @*/
106 56 : PetscErrorCode SVDSetTolerances(SVD svd,PetscReal tol,PetscInt maxits)
107 : {
108 56 : PetscFunctionBegin;
109 56 : PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
110 168 : PetscValidLogicalCollectiveReal(svd,tol,2);
111 168 : PetscValidLogicalCollectiveInt(svd,maxits,3);
112 56 : if (tol == (PetscReal)PETSC_DETERMINE) {
113 2 : svd->tol = PETSC_DETERMINE;
114 2 : svd->state = SVD_STATE_INITIAL;
115 54 : } else if (tol != (PetscReal)PETSC_CURRENT) {
116 39 : PetscCheck(tol>0.0,PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of tol. Must be > 0");
117 39 : svd->tol = tol;
118 : }
119 56 : if (maxits == PETSC_DETERMINE) {
120 10 : svd->max_it = PETSC_DETERMINE;
121 10 : svd->state = SVD_STATE_INITIAL;
122 46 : } else if (maxits == PETSC_UNLIMITED) {
123 0 : svd->max_it = PETSC_INT_MAX;
124 46 : } else if (maxits != PETSC_CURRENT) {
125 38 : PetscCheck(maxits>0,PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of maxits. Must be > 0");
126 38 : svd->max_it = maxits;
127 : }
128 56 : PetscFunctionReturn(PETSC_SUCCESS);
129 : }
130 :
131 : /*@
132 : SVDGetTolerances - Gets the tolerance and maximum
133 : iteration count used by the default SVD convergence tests.
134 :
135 : Not Collective
136 :
137 : Input Parameter:
138 : . svd - the singular value solver context
139 :
140 : Output Parameters:
141 : + tol - the convergence tolerance
142 : - maxits - maximum number of iterations
143 :
144 : Notes:
145 : The user can specify NULL for any parameter that is not needed.
146 :
147 : Level: intermediate
148 :
149 : .seealso: SVDSetTolerances()
150 : @*/
151 69 : PetscErrorCode SVDGetTolerances(SVD svd,PetscReal *tol,PetscInt *maxits)
152 : {
153 69 : PetscFunctionBegin;
154 69 : PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
155 69 : if (tol) *tol = svd->tol;
156 69 : if (maxits) *maxits = svd->max_it;
157 69 : PetscFunctionReturn(PETSC_SUCCESS);
158 : }
159 :
160 : /*@
161 : SVDSetDimensions - Sets the number of singular values to compute
162 : and the dimension of the subspace.
163 :
164 : Logically Collective
165 :
166 : Input Parameters:
167 : + svd - the singular value solver context
168 : . nsv - number of singular values to compute
169 : . ncv - the maximum dimension of the subspace to be used by the solver
170 : - mpd - the maximum dimension allowed for the projected problem
171 :
172 : Options Database Keys:
173 : + -svd_nsv <nsv> - Sets the number of singular values
174 : . -svd_ncv <ncv> - Sets the dimension of the subspace
175 : - -svd_mpd <mpd> - Sets the maximum projected dimension
176 :
177 : Notes:
178 : Use PETSC_DETERMINE for ncv and mpd to assign a reasonably good value, which is
179 : dependent on the solution method and the number of singular values required. For
180 : any of the arguments, use PETSC_CURRENT to preserve the current value.
181 :
182 : The parameters ncv and mpd are intimately related, so that the user is advised
183 : to set one of them at most. Normal usage is that
184 : (a) in cases where nsv is small, the user sets ncv (a reasonable default is 2*nsv); and
185 : (b) in cases where nsv is large, the user sets mpd.
186 :
187 : The value of ncv should always be between nsv and (nsv+mpd), typically
188 : ncv=nsv+mpd. If nsv is not too large, mpd=nsv is a reasonable choice, otherwise
189 : a smaller value should be used.
190 :
191 : Level: intermediate
192 :
193 : .seealso: SVDGetDimensions()
194 : @*/
195 183 : PetscErrorCode SVDSetDimensions(SVD svd,PetscInt nsv,PetscInt ncv,PetscInt mpd)
196 : {
197 183 : PetscFunctionBegin;
198 183 : PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
199 549 : PetscValidLogicalCollectiveInt(svd,nsv,2);
200 549 : PetscValidLogicalCollectiveInt(svd,ncv,3);
201 549 : PetscValidLogicalCollectiveInt(svd,mpd,4);
202 183 : if (nsv != PETSC_CURRENT) {
203 183 : PetscCheck(nsv>0,PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of nsv. Must be > 0");
204 183 : svd->nsv = nsv;
205 : }
206 183 : if (ncv == PETSC_DETERMINE) {
207 126 : svd->ncv = PETSC_DETERMINE;
208 57 : } else if (ncv != PETSC_CURRENT) {
209 57 : PetscCheck(ncv>0,PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of ncv. Must be > 0");
210 57 : svd->ncv = ncv;
211 : }
212 183 : if (mpd == PETSC_DETERMINE) {
213 181 : svd->mpd = PETSC_DETERMINE;
214 2 : } else if (mpd != PETSC_CURRENT) {
215 2 : PetscCheck(mpd>0,PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of mpd. Must be > 0");
216 2 : svd->mpd = mpd;
217 : }
218 183 : svd->state = SVD_STATE_INITIAL;
219 183 : PetscFunctionReturn(PETSC_SUCCESS);
220 : }
221 :
222 : /*@
223 : SVDGetDimensions - Gets the number of singular values to compute
224 : and the dimension of the subspace.
225 :
226 : Not Collective
227 :
228 : Input Parameter:
229 : . svd - the singular value context
230 :
231 : Output Parameters:
232 : + nsv - number of singular values to compute
233 : . ncv - the maximum dimension of the subspace to be used by the solver
234 : - mpd - the maximum dimension allowed for the projected problem
235 :
236 : Notes:
237 : The user can specify NULL for any parameter that is not needed.
238 :
239 : Level: intermediate
240 :
241 : .seealso: SVDSetDimensions()
242 : @*/
243 75 : PetscErrorCode SVDGetDimensions(SVD svd,PetscInt *nsv,PetscInt *ncv,PetscInt *mpd)
244 : {
245 75 : PetscFunctionBegin;
246 75 : PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
247 75 : if (nsv) *nsv = svd->nsv;
248 75 : if (ncv) *ncv = svd->ncv;
249 75 : if (mpd) *mpd = svd->mpd;
250 75 : PetscFunctionReturn(PETSC_SUCCESS);
251 : }
252 :
253 : /*@
254 : SVDSetWhichSingularTriplets - Specifies which singular triplets are
255 : to be sought.
256 :
257 : Logically Collective
258 :
259 : Input Parameter:
260 : . svd - singular value solver context obtained from SVDCreate()
261 :
262 : Output Parameter:
263 : . which - which singular triplets are to be sought
264 :
265 : Options Database Keys:
266 : + -svd_largest - Sets largest singular values
267 : - -svd_smallest - Sets smallest singular values
268 :
269 : Notes:
270 : The parameter 'which' can have one of these values
271 :
272 : + SVD_LARGEST - largest singular values
273 : - SVD_SMALLEST - smallest singular values
274 :
275 : Level: intermediate
276 :
277 : .seealso: SVDGetWhichSingularTriplets(), SVDWhich
278 : @*/
279 41 : PetscErrorCode SVDSetWhichSingularTriplets(SVD svd,SVDWhich which)
280 : {
281 41 : PetscFunctionBegin;
282 41 : PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
283 123 : PetscValidLogicalCollectiveEnum(svd,which,2);
284 41 : switch (which) {
285 41 : case SVD_LARGEST:
286 : case SVD_SMALLEST:
287 41 : if (svd->which != which) {
288 25 : svd->state = SVD_STATE_INITIAL;
289 25 : svd->which = which;
290 : }
291 41 : break;
292 0 : default:
293 0 : SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'which' parameter");
294 : }
295 41 : PetscFunctionReturn(PETSC_SUCCESS);
296 : }
297 :
298 : /*@
299 : SVDGetWhichSingularTriplets - Returns which singular triplets are
300 : to be sought.
301 :
302 : Not Collective
303 :
304 : Input Parameter:
305 : . svd - singular value solver context obtained from SVDCreate()
306 :
307 : Output Parameter:
308 : . which - which singular triplets are to be sought
309 :
310 : Notes:
311 : See SVDSetWhichSingularTriplets() for possible values of which
312 :
313 : Level: intermediate
314 :
315 : .seealso: SVDSetWhichSingularTriplets(), SVDWhich
316 : @*/
317 11 : PetscErrorCode SVDGetWhichSingularTriplets(SVD svd,SVDWhich *which)
318 : {
319 11 : PetscFunctionBegin;
320 11 : PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
321 11 : PetscAssertPointer(which,2);
322 11 : *which = svd->which;
323 11 : PetscFunctionReturn(PETSC_SUCCESS);
324 : }
325 :
326 : /*@C
327 : SVDSetConvergenceTestFunction - Sets a function to compute the error estimate
328 : used in the convergence test.
329 :
330 : Logically Collective
331 :
332 : Input Parameters:
333 : + svd - singular value solver context obtained from SVDCreate()
334 : . conv - the convergence test function, see SVDConvergenceTestFn for the calling sequence
335 : . ctx - context for private data for the convergence routine (may be NULL)
336 : - destroy - a routine for destroying the context (may be NULL), see PetscCtxDestroyFn for the calling sequence
337 :
338 : Note:
339 : If the error estimate returned by the convergence test function is less than
340 : the tolerance, then the singular value is accepted as converged.
341 :
342 : Level: advanced
343 :
344 : .seealso: SVDSetConvergenceTest(), SVDSetTolerances()
345 : @*/
346 1 : PetscErrorCode SVDSetConvergenceTestFunction(SVD svd,SVDConvergenceTestFn *conv,void* ctx,PetscCtxDestroyFn *destroy)
347 : {
348 1 : PetscFunctionBegin;
349 1 : PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
350 1 : if (svd->convergeddestroy) PetscCall((*svd->convergeddestroy)(&svd->convergedctx));
351 1 : svd->convergeduser = conv;
352 1 : svd->convergeddestroy = destroy;
353 1 : svd->convergedctx = ctx;
354 1 : if (conv == SVDConvergedAbsolute) svd->conv = SVD_CONV_ABS;
355 1 : else if (conv == SVDConvergedRelative) svd->conv = SVD_CONV_REL;
356 1 : else if (conv == SVDConvergedNorm) svd->conv = SVD_CONV_NORM;
357 1 : else if (conv == SVDConvergedMaxIt) svd->conv = SVD_CONV_MAXIT;
358 : else {
359 1 : svd->conv = SVD_CONV_USER;
360 1 : svd->converged = svd->convergeduser;
361 : }
362 1 : PetscFunctionReturn(PETSC_SUCCESS);
363 : }
364 :
365 : /*@
366 : SVDSetConvergenceTest - Specifies how to compute the error estimate
367 : used in the convergence test.
368 :
369 : Logically Collective
370 :
371 : Input Parameters:
372 : + svd - singular value solver context obtained from SVDCreate()
373 : - conv - the type of convergence test
374 :
375 : Options Database Keys:
376 : + -svd_conv_abs - Sets the absolute convergence test
377 : . -svd_conv_rel - Sets the convergence test relative to the singular value
378 : . -svd_conv_norm - Sets the convergence test relative to the matrix norm
379 : . -svd_conv_maxit - Forces the maximum number of iterations as set by -svd_max_it
380 : - -svd_conv_user - Selects the user-defined convergence test
381 :
382 : Notes:
383 : The parameter 'conv' can have one of these values
384 : + SVD_CONV_ABS - absolute error ||r||
385 : . SVD_CONV_REL - error relative to the singular value sigma, ||r||/sigma
386 : . SVD_CONV_NORM - error relative to the matrix norms, ||r||/||Z||, with Z=A or Z=[A;B]
387 : . SVD_CONV_MAXIT - no convergence until maximum number of iterations has been reached
388 : - SVD_CONV_USER - function set by SVDSetConvergenceTestFunction()
389 :
390 : The default in standard SVD is SVD_CONV_REL, while in GSVD the default is SVD_CONV_NORM.
391 :
392 : Level: intermediate
393 :
394 : .seealso: SVDGetConvergenceTest(), SVDSetConvergenceTestFunction(), SVDSetStoppingTest(), SVDConv
395 : @*/
396 198 : PetscErrorCode SVDSetConvergenceTest(SVD svd,SVDConv conv)
397 : {
398 198 : PetscFunctionBegin;
399 198 : PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
400 594 : PetscValidLogicalCollectiveEnum(svd,conv,2);
401 198 : switch (conv) {
402 11 : case SVD_CONV_ABS: svd->converged = SVDConvergedAbsolute; break;
403 110 : case SVD_CONV_REL: svd->converged = SVDConvergedRelative; break;
404 75 : case SVD_CONV_NORM: svd->converged = SVDConvergedNorm; break;
405 2 : case SVD_CONV_MAXIT: svd->converged = SVDConvergedMaxIt; break;
406 0 : case SVD_CONV_USER:
407 0 : PetscCheck(svd->convergeduser,PetscObjectComm((PetscObject)svd),PETSC_ERR_ORDER,"Must call SVDSetConvergenceTestFunction() first");
408 0 : svd->converged = svd->convergeduser;
409 0 : break;
410 0 : default:
411 0 : SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'conv' value");
412 : }
413 198 : svd->conv = conv;
414 198 : PetscFunctionReturn(PETSC_SUCCESS);
415 : }
416 :
417 : /*@
418 : SVDGetConvergenceTest - Gets the method used to compute the error estimate
419 : used in the convergence test.
420 :
421 : Not Collective
422 :
423 : Input Parameters:
424 : . svd - singular value solver context obtained from SVDCreate()
425 :
426 : Output Parameters:
427 : . conv - the type of convergence test
428 :
429 : Level: intermediate
430 :
431 : .seealso: SVDSetConvergenceTest(), SVDConv
432 : @*/
433 11 : PetscErrorCode SVDGetConvergenceTest(SVD svd,SVDConv *conv)
434 : {
435 11 : PetscFunctionBegin;
436 11 : PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
437 11 : PetscAssertPointer(conv,2);
438 11 : *conv = svd->conv;
439 11 : PetscFunctionReturn(PETSC_SUCCESS);
440 : }
441 :
442 : /*@C
443 : SVDSetStoppingTestFunction - Sets a function to decide when to stop the outer
444 : iteration of the singular value solver.
445 :
446 : Logically Collective
447 :
448 : Input Parameters:
449 : + svd - singular value solver context obtained from SVDCreate()
450 : . stop - the stopping test function, see SVDStoppingTestFn for the calling sequence
451 : . ctx - context for private data for the stopping routine (may be NULL)
452 : - destroy - a routine for destroying the context (may be NULL), see PetscCtxDestroyFn for the calling sequence
453 :
454 : Note:
455 : Normal usage is to first call the default routine SVDStoppingBasic() and then
456 : set reason to SVD_CONVERGED_USER if some user-defined conditions have been
457 : met. To let the singular value solver continue iterating, the result must be
458 : left as SVD_CONVERGED_ITERATING.
459 :
460 : Level: advanced
461 :
462 : .seealso: SVDSetStoppingTest(), SVDStoppingBasic()
463 : @*/
464 1 : PetscErrorCode SVDSetStoppingTestFunction(SVD svd,SVDStoppingTestFn *stop,void* ctx,PetscCtxDestroyFn *destroy)
465 : {
466 1 : PetscFunctionBegin;
467 1 : PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
468 1 : if (svd->stoppingdestroy) PetscCall((*svd->stoppingdestroy)(&svd->stoppingctx));
469 1 : svd->stoppinguser = stop;
470 1 : svd->stoppingdestroy = destroy;
471 1 : svd->stoppingctx = ctx;
472 1 : if (stop == SVDStoppingBasic) svd->stop = SVD_STOP_BASIC;
473 : else {
474 1 : svd->stop = SVD_STOP_USER;
475 1 : svd->stopping = svd->stoppinguser;
476 : }
477 1 : PetscFunctionReturn(PETSC_SUCCESS);
478 : }
479 :
480 : /*@
481 : SVDSetStoppingTest - Specifies how to decide the termination of the outer
482 : loop of the singular value solver.
483 :
484 : Logically Collective
485 :
486 : Input Parameters:
487 : + svd - singular value solver context obtained from SVDCreate()
488 : - stop - the type of stopping test
489 :
490 : Options Database Keys:
491 : + -svd_stop_basic - Sets the default stopping test
492 : - -svd_stop_user - Selects the user-defined stopping test
493 :
494 : Note:
495 : The parameter 'stop' can have one of these values
496 : + SVD_STOP_BASIC - default stopping test
497 : - SVD_STOP_USER - function set by SVDSetStoppingTestFunction()
498 :
499 : Level: advanced
500 :
501 : .seealso: SVDGetStoppingTest(), SVDSetStoppingTestFunction(), SVDSetConvergenceTest(), SVDStop
502 : @*/
503 11 : PetscErrorCode SVDSetStoppingTest(SVD svd,SVDStop stop)
504 : {
505 11 : PetscFunctionBegin;
506 11 : PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
507 33 : PetscValidLogicalCollectiveEnum(svd,stop,2);
508 11 : switch (stop) {
509 11 : case SVD_STOP_BASIC: svd->stopping = SVDStoppingBasic; break;
510 0 : case SVD_STOP_USER:
511 0 : PetscCheck(svd->stoppinguser,PetscObjectComm((PetscObject)svd),PETSC_ERR_ORDER,"Must call SVDSetStoppingTestFunction() first");
512 0 : svd->stopping = svd->stoppinguser;
513 0 : break;
514 0 : default:
515 0 : SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'stop' value");
516 : }
517 11 : svd->stop = stop;
518 11 : PetscFunctionReturn(PETSC_SUCCESS);
519 : }
520 :
521 : /*@
522 : SVDGetStoppingTest - Gets the method used to decide the termination of the outer
523 : loop of the singular value solver.
524 :
525 : Not Collective
526 :
527 : Input Parameters:
528 : . svd - singular value solver context obtained from SVDCreate()
529 :
530 : Output Parameters:
531 : . stop - the type of stopping test
532 :
533 : Level: advanced
534 :
535 : .seealso: SVDSetStoppingTest(), SVDStop
536 : @*/
537 11 : PetscErrorCode SVDGetStoppingTest(SVD svd,SVDStop *stop)
538 : {
539 11 : PetscFunctionBegin;
540 11 : PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
541 11 : PetscAssertPointer(stop,2);
542 11 : *stop = svd->stop;
543 11 : PetscFunctionReturn(PETSC_SUCCESS);
544 : }
545 :
546 : /*@C
547 : SVDMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type
548 : indicated by the user.
549 :
550 : Collective
551 :
552 : Input Parameters:
553 : + svd - the singular value solver context
554 : . opt - the command line option for this monitor
555 : . name - the monitor type one is seeking
556 : . ctx - an optional user context for the monitor, or NULL
557 : - trackall - whether this monitor tracks all singular values or not
558 :
559 : Level: developer
560 :
561 : .seealso: SVDMonitorSet(), SVDSetTrackAll()
562 : @*/
563 760 : PetscErrorCode SVDMonitorSetFromOptions(SVD svd,const char opt[],const char name[],void *ctx,PetscBool trackall)
564 : {
565 760 : PetscErrorCode (*mfunc)(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,void*);
566 760 : PetscErrorCode (*cfunc)(PetscViewer,PetscViewerFormat,void*,PetscViewerAndFormat**);
567 760 : PetscErrorCode (*dfunc)(PetscViewerAndFormat**);
568 760 : PetscViewerAndFormat *vf;
569 760 : PetscViewer viewer;
570 760 : PetscViewerFormat format;
571 760 : PetscViewerType vtype;
572 760 : char key[PETSC_MAX_PATH_LEN];
573 760 : PetscBool flg;
574 :
575 760 : PetscFunctionBegin;
576 760 : PetscCall(PetscOptionsCreateViewer(PetscObjectComm((PetscObject)svd),((PetscObject)svd)->options,((PetscObject)svd)->prefix,opt,&viewer,&format,&flg));
577 760 : if (!flg) PetscFunctionReturn(PETSC_SUCCESS);
578 :
579 6 : PetscCall(PetscViewerGetType(viewer,&vtype));
580 6 : PetscCall(SlepcMonitorMakeKey_Internal(name,vtype,format,key));
581 6 : PetscCall(PetscFunctionListFind(SVDMonitorList,key,&mfunc));
582 6 : PetscCheck(mfunc,PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"Specified viewer and format not supported");
583 6 : PetscCall(PetscFunctionListFind(SVDMonitorCreateList,key,&cfunc));
584 6 : PetscCall(PetscFunctionListFind(SVDMonitorDestroyList,key,&dfunc));
585 6 : if (!cfunc) cfunc = PetscViewerAndFormatCreate_Internal;
586 6 : if (!dfunc) dfunc = PetscViewerAndFormatDestroy;
587 :
588 6 : PetscCall((*cfunc)(viewer,format,ctx,&vf));
589 6 : PetscCall(PetscViewerDestroy(&viewer));
590 6 : PetscCall(SVDMonitorSet(svd,mfunc,vf,(PetscCtxDestroyFn*)dfunc));
591 6 : if (trackall) PetscCall(SVDSetTrackAll(svd,PETSC_TRUE));
592 6 : PetscFunctionReturn(PETSC_SUCCESS);
593 : }
594 :
595 : /*@
596 : SVDSetFromOptions - Sets SVD options from the options database.
597 : This routine must be called before SVDSetUp() if the user is to be
598 : allowed to set the solver type.
599 :
600 : Collective
601 :
602 : Input Parameters:
603 : . svd - the singular value solver context
604 :
605 : Notes:
606 : To see all options, run your program with the -help option.
607 :
608 : Level: beginner
609 :
610 : .seealso: SVDSetOptionsPrefix()
611 : @*/
612 190 : PetscErrorCode SVDSetFromOptions(SVD svd)
613 : {
614 190 : char type[256];
615 190 : PetscBool set,flg,val,flg1,flg2,flg3;
616 190 : PetscInt i,j,k;
617 190 : PetscReal r;
618 :
619 190 : PetscFunctionBegin;
620 190 : PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
621 190 : PetscCall(SVDRegisterAll());
622 570 : PetscObjectOptionsBegin((PetscObject)svd);
623 214 : PetscCall(PetscOptionsFList("-svd_type","SVD solver method","SVDSetType",SVDList,(char*)(((PetscObject)svd)->type_name?((PetscObject)svd)->type_name:SVDCROSS),type,sizeof(type),&flg));
624 190 : if (flg) PetscCall(SVDSetType(svd,type));
625 27 : else if (!((PetscObject)svd)->type_name) PetscCall(SVDSetType(svd,SVDCROSS));
626 :
627 190 : PetscCall(PetscOptionsBoolGroupBegin("-svd_standard","Singular value decomposition (SVD)","SVDSetProblemType",&flg));
628 190 : if (flg) PetscCall(SVDSetProblemType(svd,SVD_STANDARD));
629 190 : PetscCall(PetscOptionsBoolGroup("-svd_generalized","Generalized singular value decomposition (GSVD)","SVDSetProblemType",&flg));
630 190 : if (flg) PetscCall(SVDSetProblemType(svd,SVD_GENERALIZED));
631 190 : PetscCall(PetscOptionsBoolGroupEnd("-svd_hyperbolic","Hyperbolic singular value decomposition (HSVD)","SVDSetProblemType",&flg));
632 190 : if (flg) PetscCall(SVDSetProblemType(svd,SVD_HYPERBOLIC));
633 :
634 190 : PetscCall(PetscOptionsBool("-svd_implicittranspose","Handle matrix transpose implicitly","SVDSetImplicitTranspose",svd->impltrans,&val,&flg));
635 190 : if (flg) PetscCall(SVDSetImplicitTranspose(svd,val));
636 :
637 190 : i = svd->max_it;
638 190 : PetscCall(PetscOptionsInt("-svd_max_it","Maximum number of iterations","SVDSetTolerances",svd->max_it,&i,&flg1));
639 190 : r = svd->tol;
640 359 : PetscCall(PetscOptionsReal("-svd_tol","Tolerance","SVDSetTolerances",SlepcDefaultTol(svd->tol),&r,&flg2));
641 190 : if (flg1 || flg2) PetscCall(SVDSetTolerances(svd,r,i));
642 :
643 190 : PetscCall(PetscOptionsBoolGroupBegin("-svd_conv_abs","Absolute error convergence test","SVDSetConvergenceTest",&flg));
644 190 : if (flg) PetscCall(SVDSetConvergenceTest(svd,SVD_CONV_ABS));
645 190 : PetscCall(PetscOptionsBoolGroup("-svd_conv_rel","Relative error convergence test","SVDSetConvergenceTest",&flg));
646 190 : if (flg) PetscCall(SVDSetConvergenceTest(svd,SVD_CONV_REL));
647 190 : PetscCall(PetscOptionsBoolGroup("-svd_conv_norm","Convergence test relative to the matrix norms","SVDSetConvergenceTest",&flg));
648 190 : if (flg) PetscCall(SVDSetConvergenceTest(svd,SVD_CONV_NORM));
649 190 : PetscCall(PetscOptionsBoolGroup("-svd_conv_maxit","Maximum iterations convergence test","SVDSetConvergenceTest",&flg));
650 190 : if (flg) PetscCall(SVDSetConvergenceTest(svd,SVD_CONV_MAXIT));
651 190 : PetscCall(PetscOptionsBoolGroupEnd("-svd_conv_user","User-defined convergence test","SVDSetConvergenceTest",&flg));
652 190 : if (flg) PetscCall(SVDSetConvergenceTest(svd,SVD_CONV_USER));
653 :
654 190 : PetscCall(PetscOptionsBoolGroupBegin("-svd_stop_basic","Stop iteration if all singular values converged or max_it reached","SVDSetStoppingTest",&flg));
655 190 : if (flg) PetscCall(SVDSetStoppingTest(svd,SVD_STOP_BASIC));
656 190 : PetscCall(PetscOptionsBoolGroupEnd("-svd_stop_user","User-defined stopping test","SVDSetStoppingTest",&flg));
657 190 : if (flg) PetscCall(SVDSetStoppingTest(svd,SVD_STOP_USER));
658 :
659 190 : i = svd->nsv;
660 190 : PetscCall(PetscOptionsInt("-svd_nsv","Number of singular values to compute","SVDSetDimensions",svd->nsv,&i,&flg1));
661 190 : j = svd->ncv;
662 190 : PetscCall(PetscOptionsInt("-svd_ncv","Number of basis vectors","SVDSetDimensions",svd->ncv,&j,&flg2));
663 190 : k = svd->mpd;
664 190 : PetscCall(PetscOptionsInt("-svd_mpd","Maximum dimension of projected problem","SVDSetDimensions",svd->mpd,&k,&flg3));
665 190 : if (flg1 || flg2 || flg3) PetscCall(SVDSetDimensions(svd,i,j,k));
666 :
667 190 : PetscCall(PetscOptionsBoolGroupBegin("-svd_largest","Compute largest singular values","SVDSetWhichSingularTriplets",&flg));
668 190 : if (flg) PetscCall(SVDSetWhichSingularTriplets(svd,SVD_LARGEST));
669 190 : PetscCall(PetscOptionsBoolGroupEnd("-svd_smallest","Compute smallest singular values","SVDSetWhichSingularTriplets",&flg));
670 190 : if (flg) PetscCall(SVDSetWhichSingularTriplets(svd,SVD_SMALLEST));
671 :
672 : /* -----------------------------------------------------------------------*/
673 : /*
674 : Cancels all monitors hardwired into code before call to SVDSetFromOptions()
675 : */
676 190 : PetscCall(PetscOptionsBool("-svd_monitor_cancel","Remove any hardwired monitor routines","SVDMonitorCancel",PETSC_FALSE,&flg,&set));
677 190 : if (set && flg) PetscCall(SVDMonitorCancel(svd));
678 190 : PetscCall(SVDMonitorSetFromOptions(svd,"-svd_monitor","first_approximation",NULL,PETSC_FALSE));
679 190 : PetscCall(SVDMonitorSetFromOptions(svd,"-svd_monitor_all","all_approximations",NULL,PETSC_TRUE));
680 190 : PetscCall(SVDMonitorSetFromOptions(svd,"-svd_monitor_conv","convergence_history",NULL,PETSC_FALSE));
681 190 : PetscCall(SVDMonitorSetFromOptions(svd,"-svd_monitor_conditioning","conditioning",NULL,PETSC_FALSE));
682 :
683 : /* -----------------------------------------------------------------------*/
684 190 : PetscCall(PetscOptionsName("-svd_view","Print detailed information on solver used","SVDView",&set));
685 190 : PetscCall(PetscOptionsName("-svd_view_vectors","View computed singular vectors","SVDVectorsView",&set));
686 190 : PetscCall(PetscOptionsName("-svd_view_values","View computed singular values","SVDValuesView",&set));
687 190 : PetscCall(PetscOptionsName("-svd_converged_reason","Print reason for convergence, and number of iterations","SVDConvergedReasonView",&set));
688 190 : PetscCall(PetscOptionsName("-svd_error_absolute","Print absolute errors of each singular triplet","SVDErrorView",&set));
689 190 : PetscCall(PetscOptionsName("-svd_error_relative","Print relative errors of each singular triplet","SVDErrorView",&set));
690 190 : PetscCall(PetscOptionsName("-svd_error_norm","Print errors relative to the matrix norms of each singular triplet","SVDErrorView",&set));
691 :
692 190 : PetscTryTypeMethod(svd,setfromoptions,PetscOptionsObject);
693 190 : PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)svd,PetscOptionsObject));
694 190 : PetscOptionsEnd();
695 :
696 190 : if (!svd->V) PetscCall(SVDGetBV(svd,&svd->V,NULL));
697 190 : PetscCall(BVSetFromOptions(svd->V));
698 190 : if (!svd->U) PetscCall(SVDGetBV(svd,NULL,&svd->U));
699 190 : PetscCall(BVSetFromOptions(svd->U));
700 190 : if (!svd->ds) PetscCall(SVDGetDS(svd,&svd->ds));
701 190 : PetscCall(SVDSetDSType(svd));
702 190 : PetscCall(DSSetFromOptions(svd->ds));
703 190 : PetscFunctionReturn(PETSC_SUCCESS);
704 : }
705 :
706 : /*@
707 : SVDSetProblemType - Specifies the type of the singular value problem.
708 :
709 : Logically Collective
710 :
711 : Input Parameters:
712 : + svd - the singular value solver context
713 : - type - a known type of singular value problem
714 :
715 : Options Database Keys:
716 : + -svd_standard - standard singular value decomposition (SVD)
717 : . -svd_generalized - generalized singular value problem (GSVD)
718 : - -svd_hyperbolic - hyperbolic singular value problem (HSVD)
719 :
720 : Notes:
721 : The GSVD requires that two matrices have been passed via SVDSetOperators().
722 : The HSVD requires that a signature matrix has been passed via SVDSetSignature().
723 :
724 : Level: intermediate
725 :
726 : .seealso: SVDSetOperators(), SVDSetSignature(), SVDSetType(), SVDGetProblemType(), SVDProblemType
727 : @*/
728 198 : PetscErrorCode SVDSetProblemType(SVD svd,SVDProblemType type)
729 : {
730 198 : PetscFunctionBegin;
731 198 : PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
732 594 : PetscValidLogicalCollectiveEnum(svd,type,2);
733 198 : if (type == svd->problem_type) PetscFunctionReturn(PETSC_SUCCESS);
734 198 : switch (type) {
735 98 : case SVD_STANDARD:
736 98 : svd->isgeneralized = PETSC_FALSE;
737 98 : svd->ishyperbolic = PETSC_FALSE;
738 98 : break;
739 74 : case SVD_GENERALIZED:
740 74 : svd->isgeneralized = PETSC_TRUE;
741 74 : svd->ishyperbolic = PETSC_FALSE;
742 74 : break;
743 26 : case SVD_HYPERBOLIC:
744 26 : svd->isgeneralized = PETSC_FALSE;
745 26 : svd->ishyperbolic = PETSC_TRUE;
746 26 : break;
747 0 : default:
748 0 : SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_WRONG,"Unknown singular value problem type");
749 : }
750 198 : svd->problem_type = type;
751 198 : svd->state = SVD_STATE_INITIAL;
752 198 : PetscFunctionReturn(PETSC_SUCCESS);
753 : }
754 :
755 : /*@
756 : SVDGetProblemType - Gets the problem type from the SVD object.
757 :
758 : Not Collective
759 :
760 : Input Parameter:
761 : . svd - the singular value solver context
762 :
763 : Output Parameter:
764 : . type - the problem type
765 :
766 : Level: intermediate
767 :
768 : .seealso: SVDSetProblemType(), SVDProblemType
769 : @*/
770 11 : PetscErrorCode SVDGetProblemType(SVD svd,SVDProblemType *type)
771 : {
772 11 : PetscFunctionBegin;
773 11 : PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
774 11 : PetscAssertPointer(type,2);
775 11 : *type = svd->problem_type;
776 11 : PetscFunctionReturn(PETSC_SUCCESS);
777 : }
778 :
779 : /*@
780 : SVDIsGeneralized - Ask if the SVD object corresponds to a generalized
781 : singular value problem.
782 :
783 : Not Collective
784 :
785 : Input Parameter:
786 : . svd - the singular value solver context
787 :
788 : Output Parameter:
789 : . is - the answer
790 :
791 : Level: intermediate
792 :
793 : .seealso: SVDIsHyperbolic()
794 : @*/
795 11 : PetscErrorCode SVDIsGeneralized(SVD svd,PetscBool* is)
796 : {
797 11 : PetscFunctionBegin;
798 11 : PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
799 11 : PetscAssertPointer(is,2);
800 11 : *is = svd->isgeneralized;
801 11 : PetscFunctionReturn(PETSC_SUCCESS);
802 : }
803 :
804 : /*@
805 : SVDIsHyperbolic - Ask if the SVD object corresponds to a hyperbolic
806 : singular value problem.
807 :
808 : Not Collective
809 :
810 : Input Parameter:
811 : . svd - the singular value solver context
812 :
813 : Output Parameter:
814 : . is - the answer
815 :
816 : Level: intermediate
817 :
818 : .seealso: SVDIsGeneralized()
819 : @*/
820 13 : PetscErrorCode SVDIsHyperbolic(SVD svd,PetscBool* is)
821 : {
822 13 : PetscFunctionBegin;
823 13 : PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
824 13 : PetscAssertPointer(is,2);
825 13 : *is = svd->ishyperbolic;
826 13 : PetscFunctionReturn(PETSC_SUCCESS);
827 : }
828 :
829 : /*@
830 : SVDSetTrackAll - Specifies if the solver must compute the residual norm of all
831 : approximate singular value or not.
832 :
833 : Logically Collective
834 :
835 : Input Parameters:
836 : + svd - the singular value solver context
837 : - trackall - whether to compute all residuals or not
838 :
839 : Notes:
840 : If the user sets trackall=PETSC_TRUE then the solver computes (or estimates)
841 : the residual norm for each singular value approximation. Computing the residual is
842 : usually an expensive operation and solvers commonly compute only the residual
843 : associated to the first unconverged singular value.
844 :
845 : The option '-svd_monitor_all' automatically activates this option.
846 :
847 : Level: developer
848 :
849 : .seealso: SVDGetTrackAll()
850 : @*/
851 2 : PetscErrorCode SVDSetTrackAll(SVD svd,PetscBool trackall)
852 : {
853 2 : PetscFunctionBegin;
854 2 : PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
855 6 : PetscValidLogicalCollectiveBool(svd,trackall,2);
856 2 : svd->trackall = trackall;
857 2 : PetscFunctionReturn(PETSC_SUCCESS);
858 : }
859 :
860 : /*@
861 : SVDGetTrackAll - Returns the flag indicating whether all residual norms must
862 : be computed or not.
863 :
864 : Not Collective
865 :
866 : Input Parameter:
867 : . svd - the singular value solver context
868 :
869 : Output Parameter:
870 : . trackall - the returned flag
871 :
872 : Level: developer
873 :
874 : .seealso: SVDSetTrackAll()
875 : @*/
876 90 : PetscErrorCode SVDGetTrackAll(SVD svd,PetscBool *trackall)
877 : {
878 90 : PetscFunctionBegin;
879 90 : PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
880 90 : PetscAssertPointer(trackall,2);
881 90 : *trackall = svd->trackall;
882 90 : PetscFunctionReturn(PETSC_SUCCESS);
883 : }
884 :
885 : /*@
886 : SVDSetOptionsPrefix - Sets the prefix used for searching for all
887 : SVD options in the database.
888 :
889 : Logically Collective
890 :
891 : Input Parameters:
892 : + svd - the singular value solver context
893 : - prefix - the prefix string to prepend to all SVD option requests
894 :
895 : Notes:
896 : A hyphen (-) must NOT be given at the beginning of the prefix name.
897 : The first character of all runtime options is AUTOMATICALLY the
898 : hyphen.
899 :
900 : For example, to distinguish between the runtime options for two
901 : different SVD contexts, one could call
902 : .vb
903 : SVDSetOptionsPrefix(svd1,"svd1_")
904 : SVDSetOptionsPrefix(svd2,"svd2_")
905 : .ve
906 :
907 : Level: advanced
908 :
909 : .seealso: SVDAppendOptionsPrefix(), SVDGetOptionsPrefix()
910 : @*/
911 3 : PetscErrorCode SVDSetOptionsPrefix(SVD svd,const char *prefix)
912 : {
913 3 : PetscFunctionBegin;
914 3 : PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
915 3 : if (!svd->V) PetscCall(SVDGetBV(svd,&svd->V,&svd->U));
916 3 : PetscCall(BVSetOptionsPrefix(svd->V,prefix));
917 3 : PetscCall(BVSetOptionsPrefix(svd->U,prefix));
918 3 : if (!svd->ds) PetscCall(SVDGetDS(svd,&svd->ds));
919 3 : PetscCall(DSSetOptionsPrefix(svd->ds,prefix));
920 3 : PetscCall(PetscObjectSetOptionsPrefix((PetscObject)svd,prefix));
921 3 : PetscFunctionReturn(PETSC_SUCCESS);
922 : }
923 :
924 : /*@
925 : SVDAppendOptionsPrefix - Appends to the prefix used for searching for all
926 : SVD options in the database.
927 :
928 : Logically Collective
929 :
930 : Input Parameters:
931 : + svd - the singular value solver context
932 : - prefix - the prefix string to prepend to all SVD option requests
933 :
934 : Notes:
935 : A hyphen (-) must NOT be given at the beginning of the prefix name.
936 : The first character of all runtime options is AUTOMATICALLY the hyphen.
937 :
938 : Level: advanced
939 :
940 : .seealso: SVDSetOptionsPrefix(), SVDGetOptionsPrefix()
941 : @*/
942 3 : PetscErrorCode SVDAppendOptionsPrefix(SVD svd,const char *prefix)
943 : {
944 3 : PetscFunctionBegin;
945 3 : PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
946 3 : if (!svd->V) PetscCall(SVDGetBV(svd,&svd->V,&svd->U));
947 3 : PetscCall(BVAppendOptionsPrefix(svd->V,prefix));
948 3 : PetscCall(BVAppendOptionsPrefix(svd->U,prefix));
949 3 : if (!svd->ds) PetscCall(SVDGetDS(svd,&svd->ds));
950 3 : PetscCall(DSAppendOptionsPrefix(svd->ds,prefix));
951 3 : PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)svd,prefix));
952 3 : PetscFunctionReturn(PETSC_SUCCESS);
953 : }
954 :
955 : /*@
956 : SVDGetOptionsPrefix - Gets the prefix used for searching for all
957 : SVD options in the database.
958 :
959 : Not Collective
960 :
961 : Input Parameters:
962 : . svd - the singular value solver context
963 :
964 : Output Parameters:
965 : . prefix - pointer to the prefix string used is returned
966 :
967 : Note:
968 : On the Fortran side, the user should pass in a string 'prefix' of
969 : sufficient length to hold the prefix.
970 :
971 : Level: advanced
972 :
973 : .seealso: SVDSetOptionsPrefix(), SVDAppendOptionsPrefix()
974 : @*/
975 3 : PetscErrorCode SVDGetOptionsPrefix(SVD svd,const char *prefix[])
976 : {
977 3 : PetscFunctionBegin;
978 3 : PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
979 3 : PetscAssertPointer(prefix,2);
980 3 : PetscCall(PetscObjectGetOptionsPrefix((PetscObject)svd,prefix));
981 3 : PetscFunctionReturn(PETSC_SUCCESS);
982 : }
|