Actual source code: svdopts.c

slepc-main 2024-12-17
Report Typos and Errors
  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 - Activate 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).

 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.

 40:    Level: advanced

 42: .seealso: 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: 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
 83:    iteration count used by the default SVD convergence testers.

 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_UMLIMITED to indicate there is no upper bound on this value.

102:    Level: intermediate

104: .seealso: 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: 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.

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).

186:    The test against the threshold is done for converged singular values, which
187:    implies that the final number of converged singular values will be at least
188:    one more than the actual number of values below/above the threshold.

190:    Since the number of computed singular values is not known a priori, the solver
191:    will need to reallocate the basis of vectors internally, to have enough room
192:    to accommodate all the singular vectors. Hence, this option must be used with
193:    caution to avoid out-of-memory problems. The recommendation is to set the value
194:    of ncv to be larger than the estimated number of singular values, to minimize
195:    the number of reallocations.

197:    This functionality is most useful when computing largest singular values. A
198:    typical use case is to compute a low rank approximation of a matrix. Suppose
199:    we know that singular values decay abruptly around a certain index k, which
200:    is unknown. Then using a small relative threshold such as 0.2 will guarantee that
201:    the computed singular vectors capture the numerical rank k. However, if the matrix
202:    does not have low rank, i.e., singular values decay progressively, then a
203:    value of 0.2 will imply a very high cost, both computationally and in memory.

205:    If a number of wanted singular values has been set with SVDSetDimensions()
206:    it is also taken into account and the solver will stop when one of the two
207:    conditions (threshold or number of converged values) is met.

209:    Use SVDSetStoppingTest() to return to the usual computation of a fixed number
210:    of singular values.

212:    Level: advanced

214: .seealso: SVDGetThreshold(), SVDSetStoppingTest(), SVDSetDimensions(), SVDSetWhichSingularTriplets()
215: @*/
216: PetscErrorCode SVDSetThreshold(SVD svd,PetscReal thres,PetscBool rel)
217: {
218:   PetscFunctionBegin;
222:   PetscCheck(thres>0.0,PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of the threshold. Must be > 0");
223:   if (svd->thres != thres || svd->threlative != rel) {
224:     svd->thres = thres;
225:     svd->threlative = rel;
226:     svd->state = SVD_STATE_INITIAL;
227:     PetscCall(SVDSetStoppingTest(svd,SVD_STOP_THRESHOLD));
228:   }
229:   PetscFunctionReturn(PETSC_SUCCESS);
230: }

232: /*@
233:    SVDGetThreshold - Gets the threshold used by the threshold stopping test.

235:    Not Collective

237:    Input Parameter:
238: .  svd - the singular value solver context

240:    Output Parameters:
241: +  thres - the threshold
242: -  rel   - whether the threshold is relative or not

244:    Level: advanced

246: .seealso: SVDSetThreshold()
247: @*/
248: PetscErrorCode SVDGetThreshold(SVD svd,PetscReal *thres,PetscBool *rel)
249: {
250:   PetscFunctionBegin;
252:   if (thres) *thres = svd->thres;
253:   if (rel)   *rel   = svd->threlative;
254:   PetscFunctionReturn(PETSC_SUCCESS);
255: }

257: /*@
258:    SVDSetDimensions - Sets the number of singular values to compute
259:    and the dimension of the subspace.

261:    Logically Collective

263:    Input Parameters:
264: +  svd - the singular value solver context
265: .  nsv - number of singular values to compute
266: .  ncv - the maximum dimension of the subspace to be used by the solver
267: -  mpd - the maximum dimension allowed for the projected problem

269:    Options Database Keys:
270: +  -svd_nsv <nsv> - Sets the number of singular values
271: .  -svd_ncv <ncv> - Sets the dimension of the subspace
272: -  -svd_mpd <mpd> - Sets the maximum projected dimension

274:    Notes:
275:    Use PETSC_DETERMINE for ncv and mpd to assign a reasonably good value, which is
276:    dependent on the solution method and the number of singular values required. For
277:    any of the arguments, use PETSC_CURRENT to preserve the current value.

279:    The parameters ncv and mpd are intimately related, so that the user is advised
280:    to set one of them at most. Normal usage is that
281:    (a) in cases where nsv is small, the user sets ncv (a reasonable default is 2*nsv); and
282:    (b) in cases where nsv is large, the user sets mpd.

284:    The value of ncv should always be between nsv and (nsv+mpd), typically
285:    ncv=nsv+mpd. If nsv is not too large, mpd=nsv is a reasonable choice, otherwise
286:    a smaller value should be used.

288:    Level: intermediate

290: .seealso: SVDGetDimensions()
291: @*/
292: PetscErrorCode SVDSetDimensions(SVD svd,PetscInt nsv,PetscInt ncv,PetscInt mpd)
293: {
294:   PetscFunctionBegin;
299:   if (nsv != PETSC_CURRENT) {
300:     PetscCheck(nsv>0,PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of nsv. Must be > 0");
301:     svd->nsv = nsv;
302:   }
303:   if (ncv == PETSC_DETERMINE) {
304:     svd->ncv = PETSC_DETERMINE;
305:   } else if (ncv != PETSC_CURRENT) {
306:     PetscCheck(ncv>0,PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of ncv. Must be > 0");
307:     svd->ncv = ncv;
308:   }
309:   if (mpd == PETSC_DETERMINE) {
310:     svd->mpd = PETSC_DETERMINE;
311:   } else if (mpd != PETSC_CURRENT) {
312:     PetscCheck(mpd>0,PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of mpd. Must be > 0");
313:     svd->mpd = mpd;
314:   }
315:   svd->state = SVD_STATE_INITIAL;
316:   PetscFunctionReturn(PETSC_SUCCESS);
317: }

319: /*@
320:    SVDGetDimensions - Gets the number of singular values to compute
321:    and the dimension of the subspace.

323:    Not Collective

325:    Input Parameter:
326: .  svd - the singular value context

328:    Output Parameters:
329: +  nsv - number of singular values to compute
330: .  ncv - the maximum dimension of the subspace to be used by the solver
331: -  mpd - the maximum dimension allowed for the projected problem

333:    Notes:
334:    The user can specify NULL for any parameter that is not needed.

336:    Level: intermediate

338: .seealso: SVDSetDimensions()
339: @*/
340: PetscErrorCode SVDGetDimensions(SVD svd,PetscInt *nsv,PetscInt *ncv,PetscInt *mpd)
341: {
342:   PetscFunctionBegin;
344:   if (nsv) *nsv = svd->nsv? svd->nsv: 1;
345:   if (ncv) *ncv = svd->ncv;
346:   if (mpd) *mpd = svd->mpd;
347:   PetscFunctionReturn(PETSC_SUCCESS);
348: }

350: /*@
351:     SVDSetWhichSingularTriplets - Specifies which singular triplets are
352:     to be sought.

354:     Logically Collective

356:     Input Parameter:
357: .   svd - singular value solver context obtained from SVDCreate()

359:     Output Parameter:
360: .   which - which singular triplets are to be sought

362:     Options Database Keys:
363: +   -svd_largest  - Sets largest singular values
364: -   -svd_smallest - Sets smallest singular values

366:     Notes:
367:     The parameter 'which' can have one of these values

369: +     SVD_LARGEST  - largest singular values
370: -     SVD_SMALLEST - smallest singular values

372:     Level: intermediate

374: .seealso: SVDGetWhichSingularTriplets(), SVDWhich
375: @*/
376: PetscErrorCode SVDSetWhichSingularTriplets(SVD svd,SVDWhich which)
377: {
378:   PetscFunctionBegin;
381:   switch (which) {
382:     case SVD_LARGEST:
383:     case SVD_SMALLEST:
384:       if (svd->which != which) {
385:         svd->state = SVD_STATE_INITIAL;
386:         svd->which = which;
387:       }
388:       break;
389:   default:
390:     SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'which' parameter");
391:   }
392:   PetscFunctionReturn(PETSC_SUCCESS);
393: }

395: /*@
396:     SVDGetWhichSingularTriplets - Returns which singular triplets are
397:     to be sought.

399:     Not Collective

401:     Input Parameter:
402: .   svd - singular value solver context obtained from SVDCreate()

404:     Output Parameter:
405: .   which - which singular triplets are to be sought

407:     Notes:
408:     See SVDSetWhichSingularTriplets() for possible values of which

410:     Level: intermediate

412: .seealso: SVDSetWhichSingularTriplets(), SVDWhich
413: @*/
414: PetscErrorCode SVDGetWhichSingularTriplets(SVD svd,SVDWhich *which)
415: {
416:   PetscFunctionBegin;
418:   PetscAssertPointer(which,2);
419:   *which = svd->which;
420:   PetscFunctionReturn(PETSC_SUCCESS);
421: }

423: /*@C
424:    SVDSetConvergenceTestFunction - Sets a function to compute the error estimate
425:    used in the convergence test.

427:    Logically Collective

429:    Input Parameters:
430: +  svd     - singular value solver context obtained from SVDCreate()
431: .  conv    - the convergence test function, see SVDConvergenceTestFn for the calling sequence
432: .  ctx     - context for private data for the convergence routine (may be NULL)
433: -  destroy - a routine for destroying the context (may be NULL), see PetscCtxDestroyFn for the calling sequence

435:    Note:
436:    If the error estimate returned by the convergence test function is less than
437:    the tolerance, then the singular value is accepted as converged.

439:    Level: advanced

441: .seealso: SVDSetConvergenceTest(), SVDSetTolerances()
442: @*/
443: PetscErrorCode SVDSetConvergenceTestFunction(SVD svd,SVDConvergenceTestFn *conv,void* ctx,PetscCtxDestroyFn *destroy)
444: {
445:   PetscFunctionBegin;
447:   if (svd->convergeddestroy) PetscCall((*svd->convergeddestroy)(&svd->convergedctx));
448:   svd->convergeduser    = conv;
449:   svd->convergeddestroy = destroy;
450:   svd->convergedctx     = ctx;
451:   if (conv == SVDConvergedAbsolute) svd->conv = SVD_CONV_ABS;
452:   else if (conv == SVDConvergedRelative) svd->conv = SVD_CONV_REL;
453:   else if (conv == SVDConvergedNorm) svd->conv = SVD_CONV_NORM;
454:   else if (conv == SVDConvergedMaxIt) svd->conv = SVD_CONV_MAXIT;
455:   else {
456:     svd->conv      = SVD_CONV_USER;
457:     svd->converged = svd->convergeduser;
458:   }
459:   PetscFunctionReturn(PETSC_SUCCESS);
460: }

462: /*@
463:    SVDSetConvergenceTest - Specifies how to compute the error estimate
464:    used in the convergence test.

466:    Logically Collective

468:    Input Parameters:
469: +  svd  - singular value solver context obtained from SVDCreate()
470: -  conv - the type of convergence test

472:    Options Database Keys:
473: +  -svd_conv_abs   - Sets the absolute convergence test
474: .  -svd_conv_rel   - Sets the convergence test relative to the singular value
475: .  -svd_conv_norm  - Sets the convergence test relative to the matrix norm
476: .  -svd_conv_maxit - Forces the maximum number of iterations as set by -svd_max_it
477: -  -svd_conv_user  - Selects the user-defined convergence test

479:    Notes:
480:    The parameter 'conv' can have one of these values
481: +     SVD_CONV_ABS   - absolute error ||r||
482: .     SVD_CONV_REL   - error relative to the singular value sigma, ||r||/sigma
483: .     SVD_CONV_NORM  - error relative to the matrix norms, ||r||/||Z||, with Z=A or Z=[A;B]
484: .     SVD_CONV_MAXIT - no convergence until maximum number of iterations has been reached
485: -     SVD_CONV_USER  - function set by SVDSetConvergenceTestFunction()

487:    The default in standard SVD is SVD_CONV_REL, while in GSVD the default is SVD_CONV_NORM.

489:    Level: intermediate

491: .seealso: SVDGetConvergenceTest(), SVDSetConvergenceTestFunction(), SVDSetStoppingTest(), SVDConv
492: @*/
493: PetscErrorCode SVDSetConvergenceTest(SVD svd,SVDConv conv)
494: {
495:   PetscFunctionBegin;
498:   switch (conv) {
499:     case SVD_CONV_ABS:   svd->converged = SVDConvergedAbsolute; break;
500:     case SVD_CONV_REL:   svd->converged = SVDConvergedRelative; break;
501:     case SVD_CONV_NORM:  svd->converged = SVDConvergedNorm; break;
502:     case SVD_CONV_MAXIT: svd->converged = SVDConvergedMaxIt; break;
503:     case SVD_CONV_USER:
504:       PetscCheck(svd->convergeduser,PetscObjectComm((PetscObject)svd),PETSC_ERR_ORDER,"Must call SVDSetConvergenceTestFunction() first");
505:       svd->converged = svd->convergeduser;
506:       break;
507:     default:
508:       SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'conv' value");
509:   }
510:   svd->conv = conv;
511:   PetscFunctionReturn(PETSC_SUCCESS);
512: }

514: /*@
515:    SVDGetConvergenceTest - Gets the method used to compute the error estimate
516:    used in the convergence test.

518:    Not Collective

520:    Input Parameters:
521: .  svd   - singular value solver context obtained from SVDCreate()

523:    Output Parameters:
524: .  conv  - the type of convergence test

526:    Level: intermediate

528: .seealso: SVDSetConvergenceTest(), SVDConv
529: @*/
530: PetscErrorCode SVDGetConvergenceTest(SVD svd,SVDConv *conv)
531: {
532:   PetscFunctionBegin;
534:   PetscAssertPointer(conv,2);
535:   *conv = svd->conv;
536:   PetscFunctionReturn(PETSC_SUCCESS);
537: }

539: /*@C
540:    SVDSetStoppingTestFunction - Sets a function to decide when to stop the outer
541:    iteration of the singular value solver.

543:    Logically Collective

545:    Input Parameters:
546: +  svd     - singular value solver context obtained from SVDCreate()
547: .  stop    - the stopping test function, see SVDStoppingTestFn for the calling sequence
548: .  ctx     - context for private data for the stopping routine (may be NULL)
549: -  destroy - a routine for destroying the context (may be NULL), see PetscCtxDestroyFn for the calling sequence

551:    Note:
552:    Normal usage is to first call the default routine SVDStoppingBasic() and then
553:    set reason to SVD_CONVERGED_USER if some user-defined conditions have been
554:    met. To let the singular value solver continue iterating, the result must be
555:    left as SVD_CONVERGED_ITERATING.

557:    Level: advanced

559: .seealso: SVDSetStoppingTest(), SVDStoppingBasic()
560: @*/
561: PetscErrorCode SVDSetStoppingTestFunction(SVD svd,SVDStoppingTestFn *stop,void* ctx,PetscCtxDestroyFn *destroy)
562: {
563:   PetscFunctionBegin;
565:   if (svd->stoppingdestroy) PetscCall((*svd->stoppingdestroy)(&svd->stoppingctx));
566:   svd->stoppinguser    = stop;
567:   svd->stoppingdestroy = destroy;
568:   svd->stoppingctx     = ctx;
569:   if (stop == SVDStoppingBasic) PetscCall(SVDSetStoppingTest(svd,SVD_STOP_BASIC));
570:   else if (stop == SVDStoppingThreshold) PetscCall(SVDSetStoppingTest(svd,SVD_STOP_THRESHOLD));
571:   else {
572:     svd->stop     = SVD_STOP_USER;
573:     svd->stopping = svd->stoppinguser;
574:   }
575:   PetscFunctionReturn(PETSC_SUCCESS);
576: }

578: /*@
579:    SVDSetStoppingTest - Specifies how to decide the termination of the outer
580:    loop of the singular value solver.

582:    Logically Collective

584:    Input Parameters:
585: +  svd  - singular value solver context obtained from SVDCreate()
586: -  stop - the type of stopping test

588:    Options Database Keys:
589: +  -svd_stop_basic     - Sets the default stopping test
590: .  -svd_stop_threshold - Sets the threshold stopping test
591: -  -svd_stop_user      - Selects the user-defined stopping test

593:    Note:
594:    The parameter 'stop' can have one of these values
595: +     SVD_STOP_BASIC     - default stopping test
596: .     SVD_STOP_THRESHOLD - threshold stopping test
597: -     SVD_STOP_USER      - function set by SVDSetStoppingTestFunction()

599:    Level: advanced

601: .seealso: SVDGetStoppingTest(), SVDSetStoppingTestFunction(), SVDSetConvergenceTest(), SVDStop
602: @*/
603: PetscErrorCode SVDSetStoppingTest(SVD svd,SVDStop stop)
604: {
605:   PetscFunctionBegin;
608:   switch (stop) {
609:     case SVD_STOP_BASIC: svd->stopping = SVDStoppingBasic; break;
610:     case SVD_STOP_THRESHOLD: svd->stopping = SVDStoppingThreshold; break;
611:     case SVD_STOP_USER:
612:       PetscCheck(svd->stoppinguser,PetscObjectComm((PetscObject)svd),PETSC_ERR_ORDER,"Must call SVDSetStoppingTestFunction() first");
613:       svd->stopping = svd->stoppinguser;
614:       break;
615:     default:
616:       SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'stop' value");
617:   }
618:   svd->stop = stop;
619:   PetscFunctionReturn(PETSC_SUCCESS);
620: }

622: /*@
623:    SVDGetStoppingTest - Gets the method used to decide the termination of the outer
624:    loop of the singular value solver.

626:    Not Collective

628:    Input Parameters:
629: .  svd   - singular value solver context obtained from SVDCreate()

631:    Output Parameters:
632: .  stop  - the type of stopping test

634:    Level: advanced

636: .seealso: SVDSetStoppingTest(), SVDStop
637: @*/
638: PetscErrorCode SVDGetStoppingTest(SVD svd,SVDStop *stop)
639: {
640:   PetscFunctionBegin;
642:   PetscAssertPointer(stop,2);
643:   *stop = svd->stop;
644:   PetscFunctionReturn(PETSC_SUCCESS);
645: }

647: /*@C
648:    SVDMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type
649:    indicated by the user.

651:    Collective

653:    Input Parameters:
654: +  svd      - the singular value solver context
655: .  opt      - the command line option for this monitor
656: .  name     - the monitor type one is seeking
657: .  ctx      - an optional user context for the monitor, or NULL
658: -  trackall - whether this monitor tracks all singular values or not

660:    Level: developer

662: .seealso: SVDMonitorSet(), SVDSetTrackAll()
663: @*/
664: PetscErrorCode SVDMonitorSetFromOptions(SVD svd,const char opt[],const char name[],void *ctx,PetscBool trackall)
665: {
666:   PetscErrorCode       (*mfunc)(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,void*);
667:   PetscErrorCode       (*cfunc)(PetscViewer,PetscViewerFormat,void*,PetscViewerAndFormat**);
668:   PetscErrorCode       (*dfunc)(PetscViewerAndFormat**);
669:   PetscViewerAndFormat *vf;
670:   PetscViewer          viewer;
671:   PetscViewerFormat    format;
672:   PetscViewerType      vtype;
673:   char                 key[PETSC_MAX_PATH_LEN];
674:   PetscBool            flg;

676:   PetscFunctionBegin;
677:   PetscCall(PetscOptionsCreateViewer(PetscObjectComm((PetscObject)svd),((PetscObject)svd)->options,((PetscObject)svd)->prefix,opt,&viewer,&format,&flg));
678:   if (!flg) PetscFunctionReturn(PETSC_SUCCESS);

680:   PetscCall(PetscViewerGetType(viewer,&vtype));
681:   PetscCall(SlepcMonitorMakeKey_Internal(name,vtype,format,key));
682:   PetscCall(PetscFunctionListFind(SVDMonitorList,key,&mfunc));
683:   PetscCheck(mfunc,PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"Specified viewer and format not supported");
684:   PetscCall(PetscFunctionListFind(SVDMonitorCreateList,key,&cfunc));
685:   PetscCall(PetscFunctionListFind(SVDMonitorDestroyList,key,&dfunc));
686:   if (!cfunc) cfunc = PetscViewerAndFormatCreate_Internal;
687:   if (!dfunc) dfunc = PetscViewerAndFormatDestroy;

689:   PetscCall((*cfunc)(viewer,format,ctx,&vf));
690:   PetscCall(PetscViewerDestroy(&viewer));
691:   PetscCall(SVDMonitorSet(svd,mfunc,vf,(PetscCtxDestroyFn*)dfunc));
692:   if (trackall) PetscCall(SVDSetTrackAll(svd,PETSC_TRUE));
693:   PetscFunctionReturn(PETSC_SUCCESS);
694: }

696: /*@
697:    SVDSetFromOptions - Sets SVD options from the options database.
698:    This routine must be called before SVDSetUp() if the user is to be
699:    allowed to set the solver type.

701:    Collective

703:    Input Parameters:
704: .  svd - the singular value solver context

706:    Notes:
707:    To see all options, run your program with the -help option.

709:    Level: beginner

711: .seealso: SVDSetOptionsPrefix()
712: @*/
713: PetscErrorCode SVDSetFromOptions(SVD svd)
714: {
715:   char           type[256];
716:   PetscBool      set,flg,val,flg1,flg2,flg3;
717:   PetscInt       i,j,k;
718:   PetscReal      r;

720:   PetscFunctionBegin;
722:   PetscCall(SVDRegisterAll());
723:   PetscObjectOptionsBegin((PetscObject)svd);
724:     PetscCall(PetscOptionsFList("-svd_type","SVD solver method","SVDSetType",SVDList,(char*)(((PetscObject)svd)->type_name?((PetscObject)svd)->type_name:SVDCROSS),type,sizeof(type),&flg));
725:     if (flg) PetscCall(SVDSetType(svd,type));
726:     else if (!((PetscObject)svd)->type_name) PetscCall(SVDSetType(svd,SVDCROSS));

728:     PetscCall(PetscOptionsBoolGroupBegin("-svd_standard","Singular value decomposition (SVD)","SVDSetProblemType",&flg));
729:     if (flg) PetscCall(SVDSetProblemType(svd,SVD_STANDARD));
730:     PetscCall(PetscOptionsBoolGroup("-svd_generalized","Generalized singular value decomposition (GSVD)","SVDSetProblemType",&flg));
731:     if (flg) PetscCall(SVDSetProblemType(svd,SVD_GENERALIZED));
732:     PetscCall(PetscOptionsBoolGroupEnd("-svd_hyperbolic","Hyperbolic singular value decomposition (HSVD)","SVDSetProblemType",&flg));
733:     if (flg) PetscCall(SVDSetProblemType(svd,SVD_HYPERBOLIC));

735:     PetscCall(PetscOptionsBool("-svd_implicittranspose","Handle matrix transpose implicitly","SVDSetImplicitTranspose",svd->impltrans,&val,&flg));
736:     if (flg) PetscCall(SVDSetImplicitTranspose(svd,val));

738:     i = svd->max_it;
739:     PetscCall(PetscOptionsInt("-svd_max_it","Maximum number of iterations","SVDSetTolerances",svd->max_it,&i,&flg1));
740:     r = svd->tol;
741:     PetscCall(PetscOptionsReal("-svd_tol","Tolerance","SVDSetTolerances",SlepcDefaultTol(svd->tol),&r,&flg2));
742:     if (flg1 || flg2) PetscCall(SVDSetTolerances(svd,r,i));

744:     r = svd->thres;
745:     PetscCall(PetscOptionsReal("-svd_threshold_absolute","Absolute threshold","SVDSetThreshold",r,&r,&flg));
746:     if (flg) PetscCall(SVDSetThreshold(svd,r,PETSC_FALSE));
747:     PetscCall(PetscOptionsReal("-svd_threshold_relative","Relative threshold","SVDSetThreshold",r,&r,&flg));
748:     if (flg) PetscCall(SVDSetThreshold(svd,r,PETSC_TRUE));

750:     PetscCall(PetscOptionsBoolGroupBegin("-svd_conv_abs","Absolute error convergence test","SVDSetConvergenceTest",&flg));
751:     if (flg) PetscCall(SVDSetConvergenceTest(svd,SVD_CONV_ABS));
752:     PetscCall(PetscOptionsBoolGroup("-svd_conv_rel","Relative error convergence test","SVDSetConvergenceTest",&flg));
753:     if (flg) PetscCall(SVDSetConvergenceTest(svd,SVD_CONV_REL));
754:     PetscCall(PetscOptionsBoolGroup("-svd_conv_norm","Convergence test relative to the matrix norms","SVDSetConvergenceTest",&flg));
755:     if (flg) PetscCall(SVDSetConvergenceTest(svd,SVD_CONV_NORM));
756:     PetscCall(PetscOptionsBoolGroup("-svd_conv_maxit","Maximum iterations convergence test","SVDSetConvergenceTest",&flg));
757:     if (flg) PetscCall(SVDSetConvergenceTest(svd,SVD_CONV_MAXIT));
758:     PetscCall(PetscOptionsBoolGroupEnd("-svd_conv_user","User-defined convergence test","SVDSetConvergenceTest",&flg));
759:     if (flg) PetscCall(SVDSetConvergenceTest(svd,SVD_CONV_USER));

761:     PetscCall(PetscOptionsBoolGroupBegin("-svd_stop_basic","Stop iteration if all singular values converged or max_it reached","SVDSetStoppingTest",&flg));
762:     if (flg) PetscCall(SVDSetStoppingTest(svd,SVD_STOP_BASIC));
763:     PetscCall(PetscOptionsBoolGroup("-svd_stop_threshold","Stop iteration if a converged singular value is below/above the threshold","SVDSetStoppingTest",&flg));
764:     if (flg) PetscCall(SVDSetStoppingTest(svd,SVD_STOP_THRESHOLD));
765:     PetscCall(PetscOptionsBoolGroupEnd("-svd_stop_user","User-defined stopping test","SVDSetStoppingTest",&flg));
766:     if (flg) PetscCall(SVDSetStoppingTest(svd,SVD_STOP_USER));

768:     i = svd->nsv;
769:     PetscCall(PetscOptionsInt("-svd_nsv","Number of singular values to compute","SVDSetDimensions",svd->nsv,&i,&flg1));
770:     if (!flg1) i = PETSC_CURRENT;
771:     j = svd->ncv;
772:     PetscCall(PetscOptionsInt("-svd_ncv","Number of basis vectors","SVDSetDimensions",svd->ncv,&j,&flg2));
773:     k = svd->mpd;
774:     PetscCall(PetscOptionsInt("-svd_mpd","Maximum dimension of projected problem","SVDSetDimensions",svd->mpd,&k,&flg3));
775:     if (flg1 || flg2 || flg3) PetscCall(SVDSetDimensions(svd,i,j,k));

777:     PetscCall(PetscOptionsBoolGroupBegin("-svd_largest","Compute largest singular values","SVDSetWhichSingularTriplets",&flg));
778:     if (flg) PetscCall(SVDSetWhichSingularTriplets(svd,SVD_LARGEST));
779:     PetscCall(PetscOptionsBoolGroupEnd("-svd_smallest","Compute smallest singular values","SVDSetWhichSingularTriplets",&flg));
780:     if (flg) PetscCall(SVDSetWhichSingularTriplets(svd,SVD_SMALLEST));

782:     /* -----------------------------------------------------------------------*/
783:     /*
784:       Cancels all monitors hardwired into code before call to SVDSetFromOptions()
785:     */
786:     PetscCall(PetscOptionsBool("-svd_monitor_cancel","Remove any hardwired monitor routines","SVDMonitorCancel",PETSC_FALSE,&flg,&set));
787:     if (set && flg) PetscCall(SVDMonitorCancel(svd));
788:     PetscCall(SVDMonitorSetFromOptions(svd,"-svd_monitor","first_approximation",NULL,PETSC_FALSE));
789:     PetscCall(SVDMonitorSetFromOptions(svd,"-svd_monitor_all","all_approximations",NULL,PETSC_TRUE));
790:     PetscCall(SVDMonitorSetFromOptions(svd,"-svd_monitor_conv","convergence_history",NULL,PETSC_FALSE));
791:     PetscCall(SVDMonitorSetFromOptions(svd,"-svd_monitor_conditioning","conditioning",NULL,PETSC_FALSE));

793:     /* -----------------------------------------------------------------------*/
794:     PetscCall(PetscOptionsName("-svd_view","Print detailed information on solver used","SVDView",&set));
795:     PetscCall(PetscOptionsName("-svd_view_vectors","View computed singular vectors","SVDVectorsView",&set));
796:     PetscCall(PetscOptionsName("-svd_view_values","View computed singular values","SVDValuesView",&set));
797:     PetscCall(PetscOptionsName("-svd_converged_reason","Print reason for convergence, and number of iterations","SVDConvergedReasonView",&set));
798:     PetscCall(PetscOptionsName("-svd_error_absolute","Print absolute errors of each singular triplet","SVDErrorView",&set));
799:     PetscCall(PetscOptionsName("-svd_error_relative","Print relative errors of each singular triplet","SVDErrorView",&set));
800:     PetscCall(PetscOptionsName("-svd_error_norm","Print errors relative to the matrix norms of each singular triplet","SVDErrorView",&set));

802:     PetscTryTypeMethod(svd,setfromoptions,PetscOptionsObject);
803:     PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)svd,PetscOptionsObject));
804:   PetscOptionsEnd();

806:   if (!svd->V) PetscCall(SVDGetBV(svd,&svd->V,NULL));
807:   PetscCall(BVSetFromOptions(svd->V));
808:   if (!svd->U) PetscCall(SVDGetBV(svd,NULL,&svd->U));
809:   PetscCall(BVSetFromOptions(svd->U));
810:   if (!svd->ds) PetscCall(SVDGetDS(svd,&svd->ds));
811:   PetscCall(SVDSetDSType(svd));
812:   PetscCall(DSSetFromOptions(svd->ds));
813:   PetscFunctionReturn(PETSC_SUCCESS);
814: }

816: /*@
817:    SVDSetProblemType - Specifies the type of the singular value problem.

819:    Logically Collective

821:    Input Parameters:
822: +  svd  - the singular value solver context
823: -  type - a known type of singular value problem

825:    Options Database Keys:
826: +  -svd_standard    - standard singular value decomposition (SVD)
827: .  -svd_generalized - generalized singular value problem (GSVD)
828: -  -svd_hyperbolic  - hyperbolic singular value problem (HSVD)

830:    Notes:
831:    The GSVD requires that two matrices have been passed via SVDSetOperators().
832:    The HSVD requires that a signature matrix has been passed via SVDSetSignature().

834:    Level: intermediate

836: .seealso: SVDSetOperators(), SVDSetSignature(), SVDSetType(), SVDGetProblemType(), SVDProblemType
837: @*/
838: PetscErrorCode SVDSetProblemType(SVD svd,SVDProblemType type)
839: {
840:   PetscFunctionBegin;
843:   if (type == svd->problem_type) PetscFunctionReturn(PETSC_SUCCESS);
844:   switch (type) {
845:     case SVD_STANDARD:
846:       svd->isgeneralized = PETSC_FALSE;
847:       svd->ishyperbolic  = PETSC_FALSE;
848:       break;
849:     case SVD_GENERALIZED:
850:       svd->isgeneralized = PETSC_TRUE;
851:       svd->ishyperbolic  = PETSC_FALSE;
852:       break;
853:     case SVD_HYPERBOLIC:
854:       svd->isgeneralized = PETSC_FALSE;
855:       svd->ishyperbolic  = PETSC_TRUE;
856:       break;
857:     default:
858:       SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_WRONG,"Unknown singular value problem type");
859:   }
860:   svd->problem_type = type;
861:   svd->state = SVD_STATE_INITIAL;
862:   PetscFunctionReturn(PETSC_SUCCESS);
863: }

865: /*@
866:    SVDGetProblemType - Gets the problem type from the SVD object.

868:    Not Collective

870:    Input Parameter:
871: .  svd - the singular value solver context

873:    Output Parameter:
874: .  type - the problem type

876:    Level: intermediate

878: .seealso: SVDSetProblemType(), SVDProblemType
879: @*/
880: PetscErrorCode SVDGetProblemType(SVD svd,SVDProblemType *type)
881: {
882:   PetscFunctionBegin;
884:   PetscAssertPointer(type,2);
885:   *type = svd->problem_type;
886:   PetscFunctionReturn(PETSC_SUCCESS);
887: }

889: /*@
890:    SVDIsGeneralized - Ask if the SVD object corresponds to a generalized
891:    singular value problem.

893:    Not Collective

895:    Input Parameter:
896: .  svd - the singular value solver context

898:    Output Parameter:
899: .  is - the answer

901:    Level: intermediate

903: .seealso: SVDIsHyperbolic()
904: @*/
905: PetscErrorCode SVDIsGeneralized(SVD svd,PetscBool* is)
906: {
907:   PetscFunctionBegin;
909:   PetscAssertPointer(is,2);
910:   *is = svd->isgeneralized;
911:   PetscFunctionReturn(PETSC_SUCCESS);
912: }

914: /*@
915:    SVDIsHyperbolic - Ask if the SVD object corresponds to a hyperbolic
916:    singular value problem.

918:    Not Collective

920:    Input Parameter:
921: .  svd - the singular value solver context

923:    Output Parameter:
924: .  is - the answer

926:    Level: intermediate

928: .seealso: SVDIsGeneralized()
929: @*/
930: PetscErrorCode SVDIsHyperbolic(SVD svd,PetscBool* is)
931: {
932:   PetscFunctionBegin;
934:   PetscAssertPointer(is,2);
935:   *is = svd->ishyperbolic;
936:   PetscFunctionReturn(PETSC_SUCCESS);
937: }

939: /*@
940:    SVDSetTrackAll - Specifies if the solver must compute the residual norm of all
941:    approximate singular value or not.

943:    Logically Collective

945:    Input Parameters:
946: +  svd      - the singular value solver context
947: -  trackall - whether to compute all residuals or not

949:    Notes:
950:    If the user sets trackall=PETSC_TRUE then the solver computes (or estimates)
951:    the residual norm for each singular value approximation. Computing the residual is
952:    usually an expensive operation and solvers commonly compute only the residual
953:    associated to the first unconverged singular value.

955:    The option '-svd_monitor_all' automatically activates this option.

957:    Level: developer

959: .seealso: SVDGetTrackAll()
960: @*/
961: PetscErrorCode SVDSetTrackAll(SVD svd,PetscBool trackall)
962: {
963:   PetscFunctionBegin;
966:   svd->trackall = trackall;
967:   PetscFunctionReturn(PETSC_SUCCESS);
968: }

970: /*@
971:    SVDGetTrackAll - Returns the flag indicating whether all residual norms must
972:    be computed or not.

974:    Not Collective

976:    Input Parameter:
977: .  svd - the singular value solver context

979:    Output Parameter:
980: .  trackall - the returned flag

982:    Level: developer

984: .seealso: SVDSetTrackAll()
985: @*/
986: PetscErrorCode SVDGetTrackAll(SVD svd,PetscBool *trackall)
987: {
988:   PetscFunctionBegin;
990:   PetscAssertPointer(trackall,2);
991:   *trackall = svd->trackall;
992:   PetscFunctionReturn(PETSC_SUCCESS);
993: }

995: /*@
996:    SVDSetOptionsPrefix - Sets the prefix used for searching for all
997:    SVD options in the database.

999:    Logically Collective

1001:    Input Parameters:
1002: +  svd - the singular value solver context
1003: -  prefix - the prefix string to prepend to all SVD option requests

1005:    Notes:
1006:    A hyphen (-) must NOT be given at the beginning of the prefix name.
1007:    The first character of all runtime options is AUTOMATICALLY the
1008:    hyphen.

1010:    For example, to distinguish between the runtime options for two
1011:    different SVD contexts, one could call
1012: .vb
1013:       SVDSetOptionsPrefix(svd1,"svd1_")
1014:       SVDSetOptionsPrefix(svd2,"svd2_")
1015: .ve

1017:    Level: advanced

1019: .seealso: SVDAppendOptionsPrefix(), SVDGetOptionsPrefix()
1020: @*/
1021: PetscErrorCode SVDSetOptionsPrefix(SVD svd,const char *prefix)
1022: {
1023:   PetscFunctionBegin;
1025:   if (!svd->V) PetscCall(SVDGetBV(svd,&svd->V,&svd->U));
1026:   PetscCall(BVSetOptionsPrefix(svd->V,prefix));
1027:   PetscCall(BVSetOptionsPrefix(svd->U,prefix));
1028:   if (!svd->ds) PetscCall(SVDGetDS(svd,&svd->ds));
1029:   PetscCall(DSSetOptionsPrefix(svd->ds,prefix));
1030:   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)svd,prefix));
1031:   PetscFunctionReturn(PETSC_SUCCESS);
1032: }

1034: /*@
1035:    SVDAppendOptionsPrefix - Appends to the prefix used for searching for all
1036:    SVD options in the database.

1038:    Logically Collective

1040:    Input Parameters:
1041: +  svd - the singular value solver context
1042: -  prefix - the prefix string to prepend to all SVD option requests

1044:    Notes:
1045:    A hyphen (-) must NOT be given at the beginning of the prefix name.
1046:    The first character of all runtime options is AUTOMATICALLY the hyphen.

1048:    Level: advanced

1050: .seealso: SVDSetOptionsPrefix(), SVDGetOptionsPrefix()
1051: @*/
1052: PetscErrorCode SVDAppendOptionsPrefix(SVD svd,const char *prefix)
1053: {
1054:   PetscFunctionBegin;
1056:   if (!svd->V) PetscCall(SVDGetBV(svd,&svd->V,&svd->U));
1057:   PetscCall(BVAppendOptionsPrefix(svd->V,prefix));
1058:   PetscCall(BVAppendOptionsPrefix(svd->U,prefix));
1059:   if (!svd->ds) PetscCall(SVDGetDS(svd,&svd->ds));
1060:   PetscCall(DSAppendOptionsPrefix(svd->ds,prefix));
1061:   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)svd,prefix));
1062:   PetscFunctionReturn(PETSC_SUCCESS);
1063: }

1065: /*@
1066:    SVDGetOptionsPrefix - Gets the prefix used for searching for all
1067:    SVD options in the database.

1069:    Not Collective

1071:    Input Parameters:
1072: .  svd - the singular value solver context

1074:    Output Parameters:
1075: .  prefix - pointer to the prefix string used is returned

1077:    Note:
1078:    On the Fortran side, the user should pass in a string 'prefix' of
1079:    sufficient length to hold the prefix.

1081:    Level: advanced

1083: .seealso: SVDSetOptionsPrefix(), SVDAppendOptionsPrefix()
1084: @*/
1085: PetscErrorCode SVDGetOptionsPrefix(SVD svd,const char *prefix[])
1086: {
1087:   PetscFunctionBegin;
1089:   PetscAssertPointer(prefix,2);
1090:   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)svd,prefix));
1091:   PetscFunctionReturn(PETSC_SUCCESS);
1092: }