Actual source code: svdopts.c
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
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: */
14: #include <slepc/private/svdimpl.h>
15: #include <petscdraw.h>
17: /*@
18: SVDSetImplicitTranspose - Indicates how to handle the transpose of the matrix
19: associated with the singular value problem.
21: Logically Collective
23: Input Parameters:
24: + svd - the singular value solver context
25: - impl - how to handle the transpose (implicitly or not)
27: Options Database Key:
28: . -svd_implicittranspose - enable the implicit transpose mode
30: Notes:
31: By default, the transpose of the matrix is explicitly built (if the matrix
32: has defined the `MatTranspose()` operation, `MATOP_TRANSPOSE`).
34: If this flag is set to `PETSC_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 behavior, both sequentially and in parallel, but
38: requires less storage.
40: Level: advanced
42: .seealso: [](ch:svd), `SVDGetImplicitTranspose()`, `SVDSolve()`, `SVDSetOperators()`
43: @*/
44: PetscErrorCode SVDSetImplicitTranspose(SVD svd,PetscBool impl)
45: {
46: PetscFunctionBegin;
49: if (svd->impltrans!=impl) {
50: svd->impltrans = impl;
51: svd->state = SVD_STATE_INITIAL;
52: }
53: PetscFunctionReturn(PETSC_SUCCESS);
54: }
56: /*@
57: SVDGetImplicitTranspose - Gets the mode used to handle the transpose
58: of the matrix associated with the singular value problem.
60: Not Collective
62: Input Parameter:
63: . svd - the singular value solver context
65: Output Parameter:
66: . impl - how to handle the transpose (implicitly or not)
68: Level: advanced
70: .seealso: [](ch:svd), `SVDSetImplicitTranspose()`, `SVDSolve()`, `SVDSetOperators()`
71: @*/
72: PetscErrorCode SVDGetImplicitTranspose(SVD svd,PetscBool *impl)
73: {
74: PetscFunctionBegin;
76: PetscAssertPointer(impl,2);
77: *impl = svd->impltrans;
78: PetscFunctionReturn(PETSC_SUCCESS);
79: }
81: /*@
82: SVDSetTolerances - Sets the tolerance and maximum iteration count used
83: by the `SVD` convergence tests.
85: Logically Collective
87: Input Parameters:
88: + svd - the singular value solver context
89: . tol - the convergence tolerance
90: - maxits - maximum number of iterations to use
92: Options Database Keys:
93: + -svd_tol \<tol\> - sets the convergence tolerance
94: - -svd_max_it \<maxits\> - sets the maximum number of iterations allowed
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_UNLIMITED` to indicate there is no upper bound on this value.
102: Level: intermediate
104: .seealso: [](ch:svd), `SVDGetTolerances()`
105: @*/
106: PetscErrorCode SVDSetTolerances(SVD svd,PetscReal tol,PetscInt maxits)
107: {
108: PetscFunctionBegin;
112: if (tol == (PetscReal)PETSC_DETERMINE) {
113: svd->tol = PETSC_DETERMINE;
114: svd->state = SVD_STATE_INITIAL;
115: } else if (tol != (PetscReal)PETSC_CURRENT) {
116: PetscCheck(tol>0.0,PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of tol. Must be > 0");
117: svd->tol = tol;
118: }
119: if (maxits == PETSC_DETERMINE) {
120: svd->max_it = PETSC_DETERMINE;
121: svd->state = SVD_STATE_INITIAL;
122: } else if (maxits == PETSC_UNLIMITED) {
123: svd->max_it = PETSC_INT_MAX;
124: } else if (maxits != PETSC_CURRENT) {
125: PetscCheck(maxits>0,PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of maxits. Must be > 0");
126: svd->max_it = maxits;
127: }
128: PetscFunctionReturn(PETSC_SUCCESS);
129: }
131: /*@
132: SVDGetTolerances - Gets the tolerance and maximum
133: iteration count used by the default SVD convergence tests.
135: Not Collective
137: Input Parameter:
138: . svd - the singular value solver context
140: Output Parameters:
141: + tol - the convergence tolerance
142: - maxits - maximum number of iterations
144: Notes:
145: The user can specify `NULL` for any parameter that is not needed.
147: Level: intermediate
149: .seealso: [](ch:svd), `SVDSetTolerances()`
150: @*/
151: PetscErrorCode SVDGetTolerances(SVD svd,PetscReal *tol,PetscInt *maxits)
152: {
153: PetscFunctionBegin;
155: if (tol) *tol = svd->tol;
156: if (maxits) *maxits = svd->max_it;
157: PetscFunctionReturn(PETSC_SUCCESS);
158: }
160: /*@
161: SVDSetThreshold - Sets the threshold used in the threshold stopping test.
163: Logically Collective
165: Input Parameters:
166: + svd - the singular value solver context
167: . thres - the threshold value
168: - rel - whether the threshold is relative or not
170: Options Database Keys:
171: + -svd_threshold_absolute \<thres\> - sets an absolute threshold
172: - -svd_threshold_relative \<thres\> - sets a relative threshold
174: Notes:
175: This function internally calls `SVDSetStoppingTest()` to set a special stopping
176: test based on the threshold, where singular values are computed in sequence
177: until one of the computed singular values is below the threshold `thres`.
179: If the solver is configured to compute smallest singular values, then the
180: threshold must be interpreted in the opposite direction, i.e., the computation
181: will stop when one of the computed singular values is above the threshold.
183: In the case of largest singular values, the threshold can be made relative
184: with respect to the largest singular value (i.e., the matrix norm). Otherwise,
185: the argument `rel` should be `PETSC_FALSE`.
187: The test against the threshold is done for converged singular values, which
188: implies that the final number of converged singular values will be at least
189: one more than the actual number of values below/above the threshold.
191: Since the number of computed singular values is not known a priori, the solver
192: will need to reallocate the basis of vectors internally, to have enough room
193: to accommodate all the singular vectors. Hence, this option must be used with
194: caution to avoid out-of-memory problems. The recommendation is to set the value
195: of `ncv` to be larger than the estimated number of singular values, to minimize
196: the number of reallocations.
198: This functionality is most useful when computing largest singular values. A
199: typical use case is to compute a low rank approximation of a matrix. Suppose
200: we know that singular values decay abruptly around a certain index $k$, which
201: is unknown. Then using a small relative threshold such as 0.2 will guarantee that
202: the computed singular vectors capture the numerical rank $k$. However, if the matrix
203: does not have low rank, i.e., singular values decay progressively, then a
204: value of 0.2 will imply a very high cost, both computationally and in memory.
206: If a number of wanted singular values has been set with `SVDSetDimensions()`
207: it is also taken into account and the solver will stop when one of the two
208: conditions (threshold or number of converged values) is met.
210: Use `SVDSetStoppingTest()` to return to the usual computation of a fixed number
211: of singular values.
213: Level: advanced
215: .seealso: [](ch:svd), `SVDGetThreshold()`, `SVDSetStoppingTest()`, `SVDSetDimensions()`, `SVDSetWhichSingularTriplets()`
216: @*/
217: PetscErrorCode SVDSetThreshold(SVD svd,PetscReal thres,PetscBool rel)
218: {
219: PetscFunctionBegin;
223: PetscCheck(thres>0.0,PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of the threshold. Must be > 0");
224: if (svd->thres != thres || svd->threlative != rel) {
225: svd->thres = thres;
226: svd->threlative = rel;
227: svd->state = SVD_STATE_INITIAL;
228: PetscCall(SVDSetStoppingTest(svd,SVD_STOP_THRESHOLD));
229: }
230: PetscFunctionReturn(PETSC_SUCCESS);
231: }
233: /*@
234: SVDGetThreshold - Gets the threshold used by the threshold stopping test.
236: Not Collective
238: Input Parameter:
239: . svd - the singular value solver context
241: Output Parameters:
242: + thres - the threshold
243: - rel - whether the threshold is relative or not
245: Level: advanced
247: .seealso: [](ch:svd), `SVDSetThreshold()`
248: @*/
249: PetscErrorCode SVDGetThreshold(SVD svd,PetscReal *thres,PetscBool *rel)
250: {
251: PetscFunctionBegin;
253: if (thres) *thres = svd->thres;
254: if (rel) *rel = svd->threlative;
255: PetscFunctionReturn(PETSC_SUCCESS);
256: }
258: /*@
259: SVDSetDimensions - Sets the number of singular values to compute
260: and the dimension of the subspace.
262: Logically Collective
264: Input Parameters:
265: + svd - the singular value solver context
266: . nsv - number of singular values to compute
267: . ncv - the maximum dimension of the subspace to be used by the solver
268: - mpd - the maximum dimension allowed for the projected problem
270: Options Database Keys:
271: + -svd_nsv \<nsv\> - sets the number of singular values
272: . -svd_ncv \<ncv\> - sets the dimension of the subspace
273: - -svd_mpd \<mpd\> - sets the maximum projected dimension
275: Notes:
276: Use `PETSC_DETERMINE` for `ncv` and `mpd` to assign a reasonably good value, which is
277: dependent on the solution method and the number of singular values required. For
278: any of the arguments, use `PETSC_CURRENT` to preserve the current value.
280: The parameters `ncv` and `mpd` are intimately related, so that the user is advised
281: to set one of them at most. Normal usage is\:
283: 1. In cases where `nsv` is small, the user sets `ncv` (a reasonable default is `2*nsv`).
284: 1. In cases where `nsv` is large, the user sets `mpd`.
286: The value of `ncv` should always be between `nsv` and `(nsv+mpd)`, typically
287: `ncv=nsv+mpd`. If `nsv` is not too large, `mpd=nsv` is a reasonable choice, otherwise
288: a smaller value should be used.
290: Level: intermediate
292: .seealso: [](ch:svd), `SVDGetDimensions()`
293: @*/
294: PetscErrorCode SVDSetDimensions(SVD svd,PetscInt nsv,PetscInt ncv,PetscInt mpd)
295: {
296: PetscFunctionBegin;
301: if (nsv != PETSC_CURRENT) {
302: PetscCheck(nsv>0,PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of nsv. Must be > 0");
303: svd->nsv = nsv;
304: }
305: if (ncv == PETSC_DETERMINE) {
306: svd->ncv = PETSC_DETERMINE;
307: } else if (ncv != PETSC_CURRENT) {
308: PetscCheck(ncv>0,PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of ncv. Must be > 0");
309: svd->ncv = ncv;
310: }
311: if (mpd == PETSC_DETERMINE) {
312: svd->mpd = PETSC_DETERMINE;
313: } else if (mpd != PETSC_CURRENT) {
314: PetscCheck(mpd>0,PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of mpd. Must be > 0");
315: svd->mpd = mpd;
316: }
317: svd->state = SVD_STATE_INITIAL;
318: PetscFunctionReturn(PETSC_SUCCESS);
319: }
321: /*@
322: SVDGetDimensions - Gets the number of singular values to compute
323: and the dimension of the subspace.
325: Not Collective
327: Input Parameter:
328: . svd - the singular value solver context
330: Output Parameters:
331: + nsv - number of singular values to compute
332: . ncv - the maximum dimension of the subspace to be used by the solver
333: - mpd - the maximum dimension allowed for the projected problem
335: Notes:
336: The user can specify `NULL` for any parameter that is not needed.
338: Level: intermediate
340: .seealso: [](ch:svd), `SVDSetDimensions()`
341: @*/
342: PetscErrorCode SVDGetDimensions(SVD svd,PetscInt *nsv,PetscInt *ncv,PetscInt *mpd)
343: {
344: PetscFunctionBegin;
346: if (nsv) *nsv = svd->nsv? svd->nsv: 1;
347: if (ncv) *ncv = svd->ncv;
348: if (mpd) *mpd = svd->mpd;
349: PetscFunctionReturn(PETSC_SUCCESS);
350: }
352: /*@
353: SVDSetWhichSingularTriplets - Specifies which singular triplets are
354: to be sought.
356: Logically Collective
358: Input Parameter:
359: . svd - the singular value solver context
361: Output Parameter:
362: . which - which singular triplets are to be sought, see `SVDWhich` for possible values
364: Options Database Keys:
365: + -svd_largest - sets largest singular values
366: - -svd_smallest - sets smallest singular values
368: Level: intermediate
370: .seealso: [](ch:svd), `SVDGetWhichSingularTriplets()`, `SVDWhich`
371: @*/
372: PetscErrorCode SVDSetWhichSingularTriplets(SVD svd,SVDWhich which)
373: {
374: PetscFunctionBegin;
377: switch (which) {
378: case SVD_LARGEST:
379: case SVD_SMALLEST:
380: if (svd->which != which) {
381: svd->state = SVD_STATE_INITIAL;
382: svd->which = which;
383: }
384: break;
385: default:
386: SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'which' parameter");
387: }
388: PetscFunctionReturn(PETSC_SUCCESS);
389: }
391: /*@
392: SVDGetWhichSingularTriplets - Returns which singular triplets are
393: to be sought.
395: Not Collective
397: Input Parameter:
398: . svd - the singular value solver context
400: Output Parameter:
401: . which - which singular triplets are to be sought
403: Level: intermediate
405: .seealso: [](ch:svd), `SVDSetWhichSingularTriplets()`, `SVDWhich`
406: @*/
407: PetscErrorCode SVDGetWhichSingularTriplets(SVD svd,SVDWhich *which)
408: {
409: PetscFunctionBegin;
411: PetscAssertPointer(which,2);
412: *which = svd->which;
413: PetscFunctionReturn(PETSC_SUCCESS);
414: }
416: /*@C
417: SVDSetConvergenceTestFunction - Sets a function to compute the error estimate
418: used in the convergence test.
420: Logically Collective
422: Input Parameters:
423: + svd - the singular value solver context
424: . conv - the convergence test function, see `SVDConvergenceTestFn` for the calling sequence
425: . ctx - context for private data for the convergence routine (may be `NULL`)
426: - destroy - a routine for destroying the context (may be `NULL`), see `PetscCtxDestroyFn`
427: for the calling sequence
429: Notes:
430: When this is called with a user-defined function, then the convergence
431: criterion is set to `SVD_CONV_USER`, see `SVDSetConvergenceTest()`.
433: If the error estimate returned by the convergence test function is less than
434: the tolerance, then the singular value is accepted as converged.
436: Level: advanced
438: .seealso: [](ch:svd), `SVDSetConvergenceTest()`, `SVDSetTolerances()`
439: @*/
440: PetscErrorCode SVDSetConvergenceTestFunction(SVD svd,SVDConvergenceTestFn *conv,void *ctx,PetscCtxDestroyFn *destroy)
441: {
442: PetscFunctionBegin;
444: if (svd->convergeddestroy) PetscCall((*svd->convergeddestroy)(&svd->convergedctx));
445: svd->convergeduser = conv;
446: svd->convergeddestroy = destroy;
447: svd->convergedctx = ctx;
448: if (conv == SVDConvergedAbsolute) svd->conv = SVD_CONV_ABS;
449: else if (conv == SVDConvergedRelative) svd->conv = SVD_CONV_REL;
450: else if (conv == SVDConvergedNorm) svd->conv = SVD_CONV_NORM;
451: else if (conv == SVDConvergedMaxIt) svd->conv = SVD_CONV_MAXIT;
452: else {
453: svd->conv = SVD_CONV_USER;
454: svd->converged = svd->convergeduser;
455: }
456: PetscFunctionReturn(PETSC_SUCCESS);
457: }
459: /*@
460: SVDSetConvergenceTest - Specifies how to compute the error estimate
461: used in the convergence test.
463: Logically Collective
465: Input Parameters:
466: + svd - the singular value solver context
467: - conv - the type of convergence test, see `SVDConv` for possible values
469: Options Database Keys:
470: + -svd_conv_abs - sets the absolute convergence test
471: . -svd_conv_rel - sets the convergence test relative to the singular value
472: . -svd_conv_norm - sets the convergence test relative to the matrix norms
473: . -svd_conv_maxit - forces the maximum number of iterations as set by `-svd_max_it`
474: - -svd_conv_user - selects the user-defined convergence test
476: Note:
477: The default in standard SVD is `SVD_CONV_REL`, while in GSVD the default is `SVD_CONV_NORM`.
479: Level: intermediate
481: .seealso: [](ch:svd), `SVDGetConvergenceTest()`, `SVDSetConvergenceTestFunction()`, `SVDSetStoppingTest()`, `SVDConv`
482: @*/
483: PetscErrorCode SVDSetConvergenceTest(SVD svd,SVDConv conv)
484: {
485: PetscFunctionBegin;
488: switch (conv) {
489: case SVD_CONV_ABS: svd->converged = SVDConvergedAbsolute; break;
490: case SVD_CONV_REL: svd->converged = SVDConvergedRelative; break;
491: case SVD_CONV_NORM: svd->converged = SVDConvergedNorm; break;
492: case SVD_CONV_MAXIT: svd->converged = SVDConvergedMaxIt; break;
493: case SVD_CONV_USER:
494: PetscCheck(svd->convergeduser,PetscObjectComm((PetscObject)svd),PETSC_ERR_ORDER,"Must call SVDSetConvergenceTestFunction() first");
495: svd->converged = svd->convergeduser;
496: break;
497: default:
498: SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'conv' value");
499: }
500: svd->conv = conv;
501: PetscFunctionReturn(PETSC_SUCCESS);
502: }
504: /*@
505: SVDGetConvergenceTest - Gets the method used to compute the error estimate
506: used in the convergence test.
508: Not Collective
510: Input Parameter:
511: . svd - the singular value solver context
513: Output Parameter:
514: . conv - the type of convergence test
516: Level: intermediate
518: .seealso: [](ch:svd), `SVDSetConvergenceTest()`, `SVDConv`
519: @*/
520: PetscErrorCode SVDGetConvergenceTest(SVD svd,SVDConv *conv)
521: {
522: PetscFunctionBegin;
524: PetscAssertPointer(conv,2);
525: *conv = svd->conv;
526: PetscFunctionReturn(PETSC_SUCCESS);
527: }
529: /*@C
530: SVDSetStoppingTestFunction - Sets a function to decide when to stop the outer
531: iteration of the singular value solver.
533: Logically Collective
535: Input Parameters:
536: + svd - the singular value solver context
537: . stop - the stopping test function, see `SVDStoppingTestFn` for the calling sequence
538: . ctx - context for private data for the stopping routine (may be `NULL`)
539: - destroy - a routine for destroying the context (may be `NULL`), see `PetscCtxDestroyFn`
540: for the calling sequence
542: Note:
543: When implementing a function for this, normal usage is to first call the
544: default routine `SVDStoppingBasic()` and then set `reason` to `SVD_CONVERGED_USER`
545: if some user-defined conditions have been met. To let the singular value solver
546: continue iterating, the result must be left as `SVD_CONVERGED_ITERATING`.
548: Level: advanced
550: .seealso: [](ch:svd), `SVDSetStoppingTest()`, `SVDStoppingBasic()`
551: @*/
552: PetscErrorCode SVDSetStoppingTestFunction(SVD svd,SVDStoppingTestFn *stop,void *ctx,PetscCtxDestroyFn *destroy)
553: {
554: PetscFunctionBegin;
556: if (svd->stoppingdestroy) PetscCall((*svd->stoppingdestroy)(&svd->stoppingctx));
557: svd->stoppinguser = stop;
558: svd->stoppingdestroy = destroy;
559: svd->stoppingctx = ctx;
560: if (stop == SVDStoppingBasic) PetscCall(SVDSetStoppingTest(svd,SVD_STOP_BASIC));
561: else if (stop == SVDStoppingThreshold) PetscCall(SVDSetStoppingTest(svd,SVD_STOP_THRESHOLD));
562: else {
563: svd->stop = SVD_STOP_USER;
564: svd->stopping = svd->stoppinguser;
565: }
566: PetscFunctionReturn(PETSC_SUCCESS);
567: }
569: /*@
570: SVDSetStoppingTest - Specifies how to decide the termination of the outer
571: loop of the singular value solver.
573: Logically Collective
575: Input Parameters:
576: + svd - the singular value solver context
577: - stop - the type of stopping test, see `SVDStop`
579: Options Database Keys:
580: + -svd_stop_basic - sets the default stopping test
581: . -svd_stop_threshold - sets the threshold stopping test
582: - -svd_stop_user - selects the user-defined stopping test
584: Level: advanced
586: .seealso: [](ch:svd), `SVDGetStoppingTest()`, `SVDSetStoppingTestFunction()`, `SVDSetConvergenceTest()`, `SVDStop`
587: @*/
588: PetscErrorCode SVDSetStoppingTest(SVD svd,SVDStop stop)
589: {
590: PetscFunctionBegin;
593: switch (stop) {
594: case SVD_STOP_BASIC: svd->stopping = SVDStoppingBasic; break;
595: case SVD_STOP_THRESHOLD: svd->stopping = SVDStoppingThreshold; break;
596: case SVD_STOP_USER:
597: PetscCheck(svd->stoppinguser,PetscObjectComm((PetscObject)svd),PETSC_ERR_ORDER,"Must call SVDSetStoppingTestFunction() first");
598: svd->stopping = svd->stoppinguser;
599: break;
600: default:
601: SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'stop' value");
602: }
603: svd->stop = stop;
604: PetscFunctionReturn(PETSC_SUCCESS);
605: }
607: /*@
608: SVDGetStoppingTest - Gets the method used to decide the termination of the outer
609: loop of the singular value solver.
611: Not Collective
613: Input Parameter:
614: . svd - the singular value solver context
616: Output Parameter:
617: . stop - the type of stopping test
619: Level: advanced
621: .seealso: [](ch:svd), `SVDSetStoppingTest()`, `SVDStop`
622: @*/
623: PetscErrorCode SVDGetStoppingTest(SVD svd,SVDStop *stop)
624: {
625: PetscFunctionBegin;
627: PetscAssertPointer(stop,2);
628: *stop = svd->stop;
629: PetscFunctionReturn(PETSC_SUCCESS);
630: }
632: /*@C
633: SVDMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type
634: indicated by the user.
636: Collective
638: Input Parameters:
639: + svd - the singular value solver context
640: . opt - the command line option for this monitor
641: . name - the monitor type one is seeking
642: . ctx - an optional user context for the monitor, or `NULL`
643: - trackall - whether this monitor tracks all singular values or not
645: Level: developer
647: .seealso: [](ch:svd), `SVDMonitorSet()`, `SVDSetTrackAll()`
648: @*/
649: PetscErrorCode SVDMonitorSetFromOptions(SVD svd,const char opt[],const char name[],void *ctx,PetscBool trackall)
650: {
651: PetscErrorCode (*mfunc)(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,void*);
652: PetscErrorCode (*cfunc)(PetscViewer,PetscViewerFormat,void*,PetscViewerAndFormat**);
653: PetscErrorCode (*dfunc)(PetscViewerAndFormat**);
654: PetscViewerAndFormat *vf;
655: PetscViewer viewer;
656: PetscViewerFormat format;
657: PetscViewerType vtype;
658: char key[PETSC_MAX_PATH_LEN];
659: PetscBool flg;
661: PetscFunctionBegin;
662: PetscCall(PetscOptionsCreateViewer(PetscObjectComm((PetscObject)svd),((PetscObject)svd)->options,((PetscObject)svd)->prefix,opt,&viewer,&format,&flg));
663: if (!flg) PetscFunctionReturn(PETSC_SUCCESS);
665: PetscCall(PetscViewerGetType(viewer,&vtype));
666: PetscCall(SlepcMonitorMakeKey_Internal(name,vtype,format,key));
667: PetscCall(PetscFunctionListFind(SVDMonitorList,key,&mfunc));
668: PetscCheck(mfunc,PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"Specified viewer and format not supported");
669: PetscCall(PetscFunctionListFind(SVDMonitorCreateList,key,&cfunc));
670: PetscCall(PetscFunctionListFind(SVDMonitorDestroyList,key,&dfunc));
671: if (!cfunc) cfunc = PetscViewerAndFormatCreate_Internal;
672: if (!dfunc) dfunc = PetscViewerAndFormatDestroy;
674: PetscCall((*cfunc)(viewer,format,ctx,&vf));
675: PetscCall(PetscViewerDestroy(&viewer));
676: PetscCall(SVDMonitorSet(svd,mfunc,vf,(PetscCtxDestroyFn*)dfunc));
677: if (trackall) PetscCall(SVDSetTrackAll(svd,PETSC_TRUE));
678: PetscFunctionReturn(PETSC_SUCCESS);
679: }
681: /*@
682: SVDSetFromOptions - Sets `SVD` options from the options database.
683: This routine must be called before `SVDSetUp()` if the user is to be
684: allowed to configure the solver.
686: Collective
688: Input Parameter:
689: . svd - the singular value solver context
691: Note:
692: To see all options, run your program with the `-help` option.
694: Level: beginner
696: .seealso: [](ch:svd), `SVDSetOptionsPrefix()`
697: @*/
698: PetscErrorCode SVDSetFromOptions(SVD svd)
699: {
700: char type[256];
701: PetscBool set,flg,val,flg1,flg2,flg3;
702: PetscInt i,j,k;
703: PetscReal r;
705: PetscFunctionBegin;
707: PetscCall(SVDRegisterAll());
708: PetscObjectOptionsBegin((PetscObject)svd);
709: PetscCall(PetscOptionsFList("-svd_type","SVD solver method","SVDSetType",SVDList,(char*)(((PetscObject)svd)->type_name?((PetscObject)svd)->type_name:SVDCROSS),type,sizeof(type),&flg));
710: if (flg) PetscCall(SVDSetType(svd,type));
711: else if (!((PetscObject)svd)->type_name) PetscCall(SVDSetType(svd,SVDCROSS));
713: PetscCall(PetscOptionsBoolGroupBegin("-svd_standard","Singular value decomposition (SVD)","SVDSetProblemType",&flg));
714: if (flg) PetscCall(SVDSetProblemType(svd,SVD_STANDARD));
715: PetscCall(PetscOptionsBoolGroup("-svd_generalized","Generalized singular value decomposition (GSVD)","SVDSetProblemType",&flg));
716: if (flg) PetscCall(SVDSetProblemType(svd,SVD_GENERALIZED));
717: PetscCall(PetscOptionsBoolGroupEnd("-svd_hyperbolic","Hyperbolic singular value decomposition (HSVD)","SVDSetProblemType",&flg));
718: if (flg) PetscCall(SVDSetProblemType(svd,SVD_HYPERBOLIC));
720: PetscCall(PetscOptionsBool("-svd_implicittranspose","Handle matrix transpose implicitly","SVDSetImplicitTranspose",svd->impltrans,&val,&flg));
721: if (flg) PetscCall(SVDSetImplicitTranspose(svd,val));
723: i = svd->max_it;
724: PetscCall(PetscOptionsInt("-svd_max_it","Maximum number of iterations","SVDSetTolerances",svd->max_it,&i,&flg1));
725: r = svd->tol;
726: PetscCall(PetscOptionsReal("-svd_tol","Tolerance","SVDSetTolerances",SlepcDefaultTol(svd->tol),&r,&flg2));
727: if (flg1 || flg2) PetscCall(SVDSetTolerances(svd,r,i));
729: r = svd->thres;
730: PetscCall(PetscOptionsReal("-svd_threshold_absolute","Absolute threshold","SVDSetThreshold",r,&r,&flg));
731: if (flg) PetscCall(SVDSetThreshold(svd,r,PETSC_FALSE));
732: PetscCall(PetscOptionsReal("-svd_threshold_relative","Relative threshold","SVDSetThreshold",r,&r,&flg));
733: if (flg) PetscCall(SVDSetThreshold(svd,r,PETSC_TRUE));
735: PetscCall(PetscOptionsBoolGroupBegin("-svd_conv_abs","Absolute error convergence test","SVDSetConvergenceTest",&flg));
736: if (flg) PetscCall(SVDSetConvergenceTest(svd,SVD_CONV_ABS));
737: PetscCall(PetscOptionsBoolGroup("-svd_conv_rel","Relative error convergence test","SVDSetConvergenceTest",&flg));
738: if (flg) PetscCall(SVDSetConvergenceTest(svd,SVD_CONV_REL));
739: PetscCall(PetscOptionsBoolGroup("-svd_conv_norm","Convergence test relative to the matrix norms","SVDSetConvergenceTest",&flg));
740: if (flg) PetscCall(SVDSetConvergenceTest(svd,SVD_CONV_NORM));
741: PetscCall(PetscOptionsBoolGroup("-svd_conv_maxit","Maximum iterations convergence test","SVDSetConvergenceTest",&flg));
742: if (flg) PetscCall(SVDSetConvergenceTest(svd,SVD_CONV_MAXIT));
743: PetscCall(PetscOptionsBoolGroupEnd("-svd_conv_user","User-defined convergence test","SVDSetConvergenceTest",&flg));
744: if (flg) PetscCall(SVDSetConvergenceTest(svd,SVD_CONV_USER));
746: PetscCall(PetscOptionsBoolGroupBegin("-svd_stop_basic","Stop iteration if all singular values converged or max_it reached","SVDSetStoppingTest",&flg));
747: if (flg) PetscCall(SVDSetStoppingTest(svd,SVD_STOP_BASIC));
748: PetscCall(PetscOptionsBoolGroup("-svd_stop_threshold","Stop iteration if a converged singular value is below/above the threshold","SVDSetStoppingTest",&flg));
749: if (flg) PetscCall(SVDSetStoppingTest(svd,SVD_STOP_THRESHOLD));
750: PetscCall(PetscOptionsBoolGroupEnd("-svd_stop_user","User-defined stopping test","SVDSetStoppingTest",&flg));
751: if (flg) PetscCall(SVDSetStoppingTest(svd,SVD_STOP_USER));
753: i = svd->nsv;
754: PetscCall(PetscOptionsInt("-svd_nsv","Number of singular values to compute","SVDSetDimensions",svd->nsv,&i,&flg1));
755: if (!flg1) i = PETSC_CURRENT;
756: j = svd->ncv;
757: PetscCall(PetscOptionsInt("-svd_ncv","Number of basis vectors","SVDSetDimensions",svd->ncv,&j,&flg2));
758: k = svd->mpd;
759: PetscCall(PetscOptionsInt("-svd_mpd","Maximum dimension of projected problem","SVDSetDimensions",svd->mpd,&k,&flg3));
760: if (flg1 || flg2 || flg3) PetscCall(SVDSetDimensions(svd,i,j,k));
762: PetscCall(PetscOptionsBoolGroupBegin("-svd_largest","Compute largest singular values","SVDSetWhichSingularTriplets",&flg));
763: if (flg) PetscCall(SVDSetWhichSingularTriplets(svd,SVD_LARGEST));
764: PetscCall(PetscOptionsBoolGroupEnd("-svd_smallest","Compute smallest singular values","SVDSetWhichSingularTriplets",&flg));
765: if (flg) PetscCall(SVDSetWhichSingularTriplets(svd,SVD_SMALLEST));
767: /* -----------------------------------------------------------------------*/
768: /*
769: Cancels all monitors hardwired into code before call to SVDSetFromOptions()
770: */
771: PetscCall(PetscOptionsBool("-svd_monitor_cancel","Remove any hardwired monitor routines","SVDMonitorCancel",PETSC_FALSE,&flg,&set));
772: if (set && flg) PetscCall(SVDMonitorCancel(svd));
773: PetscCall(SVDMonitorSetFromOptions(svd,"-svd_monitor","first_approximation",NULL,PETSC_FALSE));
774: PetscCall(SVDMonitorSetFromOptions(svd,"-svd_monitor_all","all_approximations",NULL,PETSC_TRUE));
775: PetscCall(SVDMonitorSetFromOptions(svd,"-svd_monitor_conv","convergence_history",NULL,PETSC_FALSE));
776: PetscCall(SVDMonitorSetFromOptions(svd,"-svd_monitor_conditioning","conditioning",NULL,PETSC_FALSE));
778: /* -----------------------------------------------------------------------*/
779: PetscCall(PetscOptionsName("-svd_view","Print detailed information on solver used","SVDView",&set));
780: PetscCall(PetscOptionsName("-svd_view_vectors","View computed singular vectors","SVDVectorsView",&set));
781: PetscCall(PetscOptionsName("-svd_view_values","View computed singular values","SVDValuesView",&set));
782: PetscCall(PetscOptionsName("-svd_converged_reason","Print reason for convergence, and number of iterations","SVDConvergedReasonView",&set));
783: PetscCall(PetscOptionsName("-svd_error_absolute","Print absolute errors of each singular triplet","SVDErrorView",&set));
784: PetscCall(PetscOptionsName("-svd_error_relative","Print relative errors of each singular triplet","SVDErrorView",&set));
785: PetscCall(PetscOptionsName("-svd_error_norm","Print errors relative to the matrix norms of each singular triplet","SVDErrorView",&set));
787: PetscTryTypeMethod(svd,setfromoptions,PetscOptionsObject);
788: PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)svd,PetscOptionsObject));
789: PetscOptionsEnd();
791: if (!svd->V) PetscCall(SVDGetBV(svd,&svd->V,NULL));
792: PetscCall(BVSetFromOptions(svd->V));
793: if (!svd->U) PetscCall(SVDGetBV(svd,NULL,&svd->U));
794: PetscCall(BVSetFromOptions(svd->U));
795: if (!svd->ds) PetscCall(SVDGetDS(svd,&svd->ds));
796: PetscCall(SVDSetDSType(svd));
797: PetscCall(DSSetFromOptions(svd->ds));
798: PetscFunctionReturn(PETSC_SUCCESS);
799: }
801: /*@
802: SVDSetProblemType - Specifies the type of the singular value problem.
804: Logically Collective
806: Input Parameters:
807: + svd - the singular value solver context
808: - type - a known type of singular value problem
810: Options Database Keys:
811: + -svd_standard - standard singular value decomposition (SVD)
812: . -svd_generalized - generalized singular value problem (GSVD)
813: - -svd_hyperbolic - hyperbolic singular value problem (HSVD)
815: Note:
816: The GSVD requires that two matrices have been passed via `SVDSetOperators()`.
817: The HSVD requires that a signature matrix has been passed via `SVDSetSignature()`.
819: Level: intermediate
821: .seealso: [](ch:svd), `SVDSetOperators()`, `SVDSetSignature()`, `SVDSetType()`, `SVDGetProblemType()`, `SVDProblemType`
822: @*/
823: PetscErrorCode SVDSetProblemType(SVD svd,SVDProblemType type)
824: {
825: PetscFunctionBegin;
828: if (type == svd->problem_type) PetscFunctionReturn(PETSC_SUCCESS);
829: switch (type) {
830: case SVD_STANDARD:
831: svd->isgeneralized = PETSC_FALSE;
832: svd->ishyperbolic = PETSC_FALSE;
833: break;
834: case SVD_GENERALIZED:
835: svd->isgeneralized = PETSC_TRUE;
836: svd->ishyperbolic = PETSC_FALSE;
837: break;
838: case SVD_HYPERBOLIC:
839: svd->isgeneralized = PETSC_FALSE;
840: svd->ishyperbolic = PETSC_TRUE;
841: break;
842: default:
843: SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_WRONG,"Unknown singular value problem type");
844: }
845: svd->problem_type = type;
846: svd->state = SVD_STATE_INITIAL;
847: PetscFunctionReturn(PETSC_SUCCESS);
848: }
850: /*@
851: SVDGetProblemType - Gets the problem type from the SVD object.
853: Not Collective
855: Input Parameter:
856: . svd - the singular value solver context
858: Output Parameter:
859: . type - the problem type
861: Level: intermediate
863: .seealso: [](ch:svd), `SVDSetProblemType()`, `SVDProblemType`
864: @*/
865: PetscErrorCode SVDGetProblemType(SVD svd,SVDProblemType *type)
866: {
867: PetscFunctionBegin;
869: PetscAssertPointer(type,2);
870: *type = svd->problem_type;
871: PetscFunctionReturn(PETSC_SUCCESS);
872: }
874: /*@
875: SVDIsGeneralized - Ask if the SVD object corresponds to a generalized
876: singular value problem.
878: Not Collective
880: Input Parameter:
881: . svd - the singular value solver context
883: Output Parameter:
884: . is - the answer
886: Level: intermediate
888: .seealso: [](ch:svd), `SVDIsHyperbolic()`
889: @*/
890: PetscErrorCode SVDIsGeneralized(SVD svd,PetscBool* is)
891: {
892: PetscFunctionBegin;
894: PetscAssertPointer(is,2);
895: *is = svd->isgeneralized;
896: PetscFunctionReturn(PETSC_SUCCESS);
897: }
899: /*@
900: SVDIsHyperbolic - Ask if the SVD object corresponds to a hyperbolic
901: singular value problem.
903: Not Collective
905: Input Parameter:
906: . svd - the singular value solver context
908: Output Parameter:
909: . is - the answer
911: Level: intermediate
913: .seealso: [](ch:svd), `SVDIsGeneralized()`
914: @*/
915: PetscErrorCode SVDIsHyperbolic(SVD svd,PetscBool* is)
916: {
917: PetscFunctionBegin;
919: PetscAssertPointer(is,2);
920: *is = svd->ishyperbolic;
921: PetscFunctionReturn(PETSC_SUCCESS);
922: }
924: /*@
925: SVDSetTrackAll - Specifies if the solver must compute the residual norm of all
926: approximate singular value or not.
928: Logically Collective
930: Input Parameters:
931: + svd - the singular value solver context
932: - trackall - whether to compute all residuals or not
934: Notes:
935: If the user sets `trackall`=`PETSC_TRUE` then the solver computes (or estimates)
936: the residual norm for each singular value approximation. Computing the residual is
937: usually an expensive operation and solvers commonly compute only the residual
938: associated to the first unconverged singular value.
940: The option `-svd_monitor_all` automatically activates this option.
942: Level: developer
944: .seealso: [](ch:svd), `SVDGetTrackAll()`
945: @*/
946: PetscErrorCode SVDSetTrackAll(SVD svd,PetscBool trackall)
947: {
948: PetscFunctionBegin;
951: svd->trackall = trackall;
952: PetscFunctionReturn(PETSC_SUCCESS);
953: }
955: /*@
956: SVDGetTrackAll - Returns the flag indicating whether all residual norms must
957: be computed or not.
959: Not Collective
961: Input Parameter:
962: . svd - the singular value solver context
964: Output Parameter:
965: . trackall - the returned flag
967: Level: developer
969: .seealso: [](ch:svd), `SVDSetTrackAll()`
970: @*/
971: PetscErrorCode SVDGetTrackAll(SVD svd,PetscBool *trackall)
972: {
973: PetscFunctionBegin;
975: PetscAssertPointer(trackall,2);
976: *trackall = svd->trackall;
977: PetscFunctionReturn(PETSC_SUCCESS);
978: }
980: /*@
981: SVDSetOptionsPrefix - Sets the prefix used for searching for all
982: `SVD` options in the database.
984: Logically Collective
986: Input Parameters:
987: + svd - the singular value solver context
988: - prefix - the prefix string to prepend to all `SVD` option requests
990: Notes:
991: A hyphen (-) must NOT be given at the beginning of the prefix name.
992: The first character of all runtime options is AUTOMATICALLY the
993: hyphen.
995: For example, to distinguish between the runtime options for two
996: different `SVD` contexts, one could call
997: .vb
998: SVDSetOptionsPrefix(svd1,"svd1_")
999: SVDSetOptionsPrefix(svd2,"svd2_")
1000: .ve
1002: Level: advanced
1004: .seealso: [](ch:svd), `SVDAppendOptionsPrefix()`, `SVDGetOptionsPrefix()`
1005: @*/
1006: PetscErrorCode SVDSetOptionsPrefix(SVD svd,const char prefix[])
1007: {
1008: PetscFunctionBegin;
1010: if (!svd->V) PetscCall(SVDGetBV(svd,&svd->V,&svd->U));
1011: PetscCall(BVSetOptionsPrefix(svd->V,prefix));
1012: PetscCall(BVSetOptionsPrefix(svd->U,prefix));
1013: if (!svd->ds) PetscCall(SVDGetDS(svd,&svd->ds));
1014: PetscCall(DSSetOptionsPrefix(svd->ds,prefix));
1015: PetscCall(PetscObjectSetOptionsPrefix((PetscObject)svd,prefix));
1016: PetscFunctionReturn(PETSC_SUCCESS);
1017: }
1019: /*@
1020: SVDAppendOptionsPrefix - Appends to the prefix used for searching for all
1021: `SVD` options in the database.
1023: Logically Collective
1025: Input Parameters:
1026: + svd - the singular value solver context
1027: - prefix - the prefix string to prepend to all `SVD` option requests
1029: Notes:
1030: A hyphen (-) must NOT be given at the beginning of the prefix name.
1031: The first character of all runtime options is AUTOMATICALLY the hyphen.
1033: Level: advanced
1035: .seealso: [](ch:svd), `SVDSetOptionsPrefix()`, `SVDGetOptionsPrefix()`
1036: @*/
1037: PetscErrorCode SVDAppendOptionsPrefix(SVD svd,const char prefix[])
1038: {
1039: PetscFunctionBegin;
1041: if (!svd->V) PetscCall(SVDGetBV(svd,&svd->V,&svd->U));
1042: PetscCall(BVAppendOptionsPrefix(svd->V,prefix));
1043: PetscCall(BVAppendOptionsPrefix(svd->U,prefix));
1044: if (!svd->ds) PetscCall(SVDGetDS(svd,&svd->ds));
1045: PetscCall(DSAppendOptionsPrefix(svd->ds,prefix));
1046: PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)svd,prefix));
1047: PetscFunctionReturn(PETSC_SUCCESS);
1048: }
1050: /*@
1051: SVDGetOptionsPrefix - Gets the prefix used for searching for all
1052: `SVD` options in the database.
1054: Not Collective
1056: Input Parameter:
1057: . svd - the singular value solver context
1059: Output Parameter:
1060: . prefix - pointer to the prefix string used is returned
1062: Level: advanced
1064: .seealso: [](ch:svd), `SVDSetOptionsPrefix()`, `SVDAppendOptionsPrefix()`
1065: @*/
1066: PetscErrorCode SVDGetOptionsPrefix(SVD svd,const char *prefix[])
1067: {
1068: PetscFunctionBegin;
1070: PetscAssertPointer(prefix,2);
1071: PetscCall(PetscObjectGetOptionsPrefix((PetscObject)svd,prefix));
1072: PetscFunctionReturn(PETSC_SUCCESS);
1073: }