Actual source code: svdopts.c

slepc-main 2024-11-09
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:    SVDSetDimensions - Sets the number of singular values to compute
162:    and the dimension of the subspace.

164:    Logically Collective

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

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

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.

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.

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.

191:    Level: intermediate

193: .seealso: SVDGetDimensions()
194: @*/
195: PetscErrorCode SVDSetDimensions(SVD svd,PetscInt nsv,PetscInt ncv,PetscInt mpd)
196: {
197:   PetscFunctionBegin;
202:   if (nsv != PETSC_CURRENT) {
203:     PetscCheck(nsv>0,PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of nsv. Must be > 0");
204:     svd->nsv = nsv;
205:   }
206:   if (ncv == PETSC_DETERMINE) {
207:     svd->ncv = PETSC_DETERMINE;
208:   } else if (ncv != PETSC_CURRENT) {
209:     PetscCheck(ncv>0,PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of ncv. Must be > 0");
210:     svd->ncv = ncv;
211:   }
212:   if (mpd == PETSC_DETERMINE) {
213:     svd->mpd = PETSC_DETERMINE;
214:   } else if (mpd != PETSC_CURRENT) {
215:     PetscCheck(mpd>0,PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of mpd. Must be > 0");
216:     svd->mpd = mpd;
217:   }
218:   svd->state = SVD_STATE_INITIAL;
219:   PetscFunctionReturn(PETSC_SUCCESS);
220: }

222: /*@
223:    SVDGetDimensions - Gets the number of singular values to compute
224:    and the dimension of the subspace.

226:    Not Collective

228:    Input Parameter:
229: .  svd - the singular value context

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

236:    Notes:
237:    The user can specify NULL for any parameter that is not needed.

239:    Level: intermediate

241: .seealso: SVDSetDimensions()
242: @*/
243: PetscErrorCode SVDGetDimensions(SVD svd,PetscInt *nsv,PetscInt *ncv,PetscInt *mpd)
244: {
245:   PetscFunctionBegin;
247:   if (nsv) *nsv = svd->nsv;
248:   if (ncv) *ncv = svd->ncv;
249:   if (mpd) *mpd = svd->mpd;
250:   PetscFunctionReturn(PETSC_SUCCESS);
251: }

253: /*@
254:     SVDSetWhichSingularTriplets - Specifies which singular triplets are
255:     to be sought.

257:     Logically Collective

259:     Input Parameter:
260: .   svd - singular value solver context obtained from SVDCreate()

262:     Output Parameter:
263: .   which - which singular triplets are to be sought

265:     Options Database Keys:
266: +   -svd_largest  - Sets largest singular values
267: -   -svd_smallest - Sets smallest singular values

269:     Notes:
270:     The parameter 'which' can have one of these values

272: +     SVD_LARGEST  - largest singular values
273: -     SVD_SMALLEST - smallest singular values

275:     Level: intermediate

277: .seealso: SVDGetWhichSingularTriplets(), SVDWhich
278: @*/
279: PetscErrorCode SVDSetWhichSingularTriplets(SVD svd,SVDWhich which)
280: {
281:   PetscFunctionBegin;
284:   switch (which) {
285:     case SVD_LARGEST:
286:     case SVD_SMALLEST:
287:       if (svd->which != which) {
288:         svd->state = SVD_STATE_INITIAL;
289:         svd->which = which;
290:       }
291:       break;
292:   default:
293:     SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'which' parameter");
294:   }
295:   PetscFunctionReturn(PETSC_SUCCESS);
296: }

298: /*@
299:     SVDGetWhichSingularTriplets - Returns which singular triplets are
300:     to be sought.

302:     Not Collective

304:     Input Parameter:
305: .   svd - singular value solver context obtained from SVDCreate()

307:     Output Parameter:
308: .   which - which singular triplets are to be sought

310:     Notes:
311:     See SVDSetWhichSingularTriplets() for possible values of which

313:     Level: intermediate

315: .seealso: SVDSetWhichSingularTriplets(), SVDWhich
316: @*/
317: PetscErrorCode SVDGetWhichSingularTriplets(SVD svd,SVDWhich *which)
318: {
319:   PetscFunctionBegin;
321:   PetscAssertPointer(which,2);
322:   *which = svd->which;
323:   PetscFunctionReturn(PETSC_SUCCESS);
324: }

326: /*@C
327:    SVDSetConvergenceTestFunction - Sets a function to compute the error estimate
328:    used in the convergence test.

330:    Logically Collective

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

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.

342:    Level: advanced

344: .seealso: SVDSetConvergenceTest(), SVDSetTolerances()
345: @*/
346: PetscErrorCode SVDSetConvergenceTestFunction(SVD svd,SVDConvergenceTestFn *conv,void* ctx,PetscCtxDestroyFn *destroy)
347: {
348:   PetscFunctionBegin;
350:   if (svd->convergeddestroy) PetscCall((*svd->convergeddestroy)(&svd->convergedctx));
351:   svd->convergeduser    = conv;
352:   svd->convergeddestroy = destroy;
353:   svd->convergedctx     = ctx;
354:   if (conv == SVDConvergedAbsolute) svd->conv = SVD_CONV_ABS;
355:   else if (conv == SVDConvergedRelative) svd->conv = SVD_CONV_REL;
356:   else if (conv == SVDConvergedNorm) svd->conv = SVD_CONV_NORM;
357:   else if (conv == SVDConvergedMaxIt) svd->conv = SVD_CONV_MAXIT;
358:   else {
359:     svd->conv      = SVD_CONV_USER;
360:     svd->converged = svd->convergeduser;
361:   }
362:   PetscFunctionReturn(PETSC_SUCCESS);
363: }

365: /*@
366:    SVDSetConvergenceTest - Specifies how to compute the error estimate
367:    used in the convergence test.

369:    Logically Collective

371:    Input Parameters:
372: +  svd  - singular value solver context obtained from SVDCreate()
373: -  conv - the type of convergence test

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

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

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

392:    Level: intermediate

394: .seealso: SVDGetConvergenceTest(), SVDSetConvergenceTestFunction(), SVDSetStoppingTest(), SVDConv
395: @*/
396: PetscErrorCode SVDSetConvergenceTest(SVD svd,SVDConv conv)
397: {
398:   PetscFunctionBegin;
401:   switch (conv) {
402:     case SVD_CONV_ABS:   svd->converged = SVDConvergedAbsolute; break;
403:     case SVD_CONV_REL:   svd->converged = SVDConvergedRelative; break;
404:     case SVD_CONV_NORM:  svd->converged = SVDConvergedNorm; break;
405:     case SVD_CONV_MAXIT: svd->converged = SVDConvergedMaxIt; break;
406:     case SVD_CONV_USER:
407:       PetscCheck(svd->convergeduser,PetscObjectComm((PetscObject)svd),PETSC_ERR_ORDER,"Must call SVDSetConvergenceTestFunction() first");
408:       svd->converged = svd->convergeduser;
409:       break;
410:     default:
411:       SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'conv' value");
412:   }
413:   svd->conv = conv;
414:   PetscFunctionReturn(PETSC_SUCCESS);
415: }

417: /*@
418:    SVDGetConvergenceTest - Gets the method used to compute the error estimate
419:    used in the convergence test.

421:    Not Collective

423:    Input Parameters:
424: .  svd   - singular value solver context obtained from SVDCreate()

426:    Output Parameters:
427: .  conv  - the type of convergence test

429:    Level: intermediate

431: .seealso: SVDSetConvergenceTest(), SVDConv
432: @*/
433: PetscErrorCode SVDGetConvergenceTest(SVD svd,SVDConv *conv)
434: {
435:   PetscFunctionBegin;
437:   PetscAssertPointer(conv,2);
438:   *conv = svd->conv;
439:   PetscFunctionReturn(PETSC_SUCCESS);
440: }

442: /*@C
443:    SVDSetStoppingTestFunction - Sets a function to decide when to stop the outer
444:    iteration of the singular value solver.

446:    Logically Collective

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

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.

460:    Level: advanced

462: .seealso: SVDSetStoppingTest(), SVDStoppingBasic()
463: @*/
464: PetscErrorCode SVDSetStoppingTestFunction(SVD svd,SVDStoppingTestFn *stop,void* ctx,PetscCtxDestroyFn *destroy)
465: {
466:   PetscFunctionBegin;
468:   if (svd->stoppingdestroy) PetscCall((*svd->stoppingdestroy)(&svd->stoppingctx));
469:   svd->stoppinguser    = stop;
470:   svd->stoppingdestroy = destroy;
471:   svd->stoppingctx     = ctx;
472:   if (stop == SVDStoppingBasic) svd->stop = SVD_STOP_BASIC;
473:   else {
474:     svd->stop     = SVD_STOP_USER;
475:     svd->stopping = svd->stoppinguser;
476:   }
477:   PetscFunctionReturn(PETSC_SUCCESS);
478: }

480: /*@
481:    SVDSetStoppingTest - Specifies how to decide the termination of the outer
482:    loop of the singular value solver.

484:    Logically Collective

486:    Input Parameters:
487: +  svd  - singular value solver context obtained from SVDCreate()
488: -  stop - the type of stopping test

490:    Options Database Keys:
491: +  -svd_stop_basic - Sets the default stopping test
492: -  -svd_stop_user  - Selects the user-defined stopping test

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

499:    Level: advanced

501: .seealso: SVDGetStoppingTest(), SVDSetStoppingTestFunction(), SVDSetConvergenceTest(), SVDStop
502: @*/
503: PetscErrorCode SVDSetStoppingTest(SVD svd,SVDStop stop)
504: {
505:   PetscFunctionBegin;
508:   switch (stop) {
509:     case SVD_STOP_BASIC: svd->stopping = SVDStoppingBasic; break;
510:     case SVD_STOP_USER:
511:       PetscCheck(svd->stoppinguser,PetscObjectComm((PetscObject)svd),PETSC_ERR_ORDER,"Must call SVDSetStoppingTestFunction() first");
512:       svd->stopping = svd->stoppinguser;
513:       break;
514:     default:
515:       SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'stop' value");
516:   }
517:   svd->stop = stop;
518:   PetscFunctionReturn(PETSC_SUCCESS);
519: }

521: /*@
522:    SVDGetStoppingTest - Gets the method used to decide the termination of the outer
523:    loop of the singular value solver.

525:    Not Collective

527:    Input Parameters:
528: .  svd   - singular value solver context obtained from SVDCreate()

530:    Output Parameters:
531: .  stop  - the type of stopping test

533:    Level: advanced

535: .seealso: SVDSetStoppingTest(), SVDStop
536: @*/
537: PetscErrorCode SVDGetStoppingTest(SVD svd,SVDStop *stop)
538: {
539:   PetscFunctionBegin;
541:   PetscAssertPointer(stop,2);
542:   *stop = svd->stop;
543:   PetscFunctionReturn(PETSC_SUCCESS);
544: }

546: /*@C
547:    SVDMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type
548:    indicated by the user.

550:    Collective

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

559:    Level: developer

561: .seealso: SVDMonitorSet(), SVDSetTrackAll()
562: @*/
563: PetscErrorCode SVDMonitorSetFromOptions(SVD svd,const char opt[],const char name[],void *ctx,PetscBool trackall)
564: {
565:   PetscErrorCode       (*mfunc)(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,void*);
566:   PetscErrorCode       (*cfunc)(PetscViewer,PetscViewerFormat,void*,PetscViewerAndFormat**);
567:   PetscErrorCode       (*dfunc)(PetscViewerAndFormat**);
568:   PetscViewerAndFormat *vf;
569:   PetscViewer          viewer;
570:   PetscViewerFormat    format;
571:   PetscViewerType      vtype;
572:   char                 key[PETSC_MAX_PATH_LEN];
573:   PetscBool            flg;

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

579:   PetscCall(PetscViewerGetType(viewer,&vtype));
580:   PetscCall(SlepcMonitorMakeKey_Internal(name,vtype,format,key));
581:   PetscCall(PetscFunctionListFind(SVDMonitorList,key,&mfunc));
582:   PetscCheck(mfunc,PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"Specified viewer and format not supported");
583:   PetscCall(PetscFunctionListFind(SVDMonitorCreateList,key,&cfunc));
584:   PetscCall(PetscFunctionListFind(SVDMonitorDestroyList,key,&dfunc));
585:   if (!cfunc) cfunc = PetscViewerAndFormatCreate_Internal;
586:   if (!dfunc) dfunc = PetscViewerAndFormatDestroy;

588:   PetscCall((*cfunc)(viewer,format,ctx,&vf));
589:   PetscCall(PetscViewerDestroy(&viewer));
590:   PetscCall(SVDMonitorSet(svd,mfunc,vf,(PetscCtxDestroyFn*)dfunc));
591:   if (trackall) PetscCall(SVDSetTrackAll(svd,PETSC_TRUE));
592:   PetscFunctionReturn(PETSC_SUCCESS);
593: }

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.

600:    Collective

602:    Input Parameters:
603: .  svd - the singular value solver context

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

608:    Level: beginner

610: .seealso: SVDSetOptionsPrefix()
611: @*/
612: PetscErrorCode SVDSetFromOptions(SVD svd)
613: {
614:   char           type[256];
615:   PetscBool      set,flg,val,flg1,flg2,flg3;
616:   PetscInt       i,j,k;
617:   PetscReal      r;

619:   PetscFunctionBegin;
621:   PetscCall(SVDRegisterAll());
622:   PetscObjectOptionsBegin((PetscObject)svd);
623:     PetscCall(PetscOptionsFList("-svd_type","SVD solver method","SVDSetType",SVDList,(char*)(((PetscObject)svd)->type_name?((PetscObject)svd)->type_name:SVDCROSS),type,sizeof(type),&flg));
624:     if (flg) PetscCall(SVDSetType(svd,type));
625:     else if (!((PetscObject)svd)->type_name) PetscCall(SVDSetType(svd,SVDCROSS));

627:     PetscCall(PetscOptionsBoolGroupBegin("-svd_standard","Singular value decomposition (SVD)","SVDSetProblemType",&flg));
628:     if (flg) PetscCall(SVDSetProblemType(svd,SVD_STANDARD));
629:     PetscCall(PetscOptionsBoolGroup("-svd_generalized","Generalized singular value decomposition (GSVD)","SVDSetProblemType",&flg));
630:     if (flg) PetscCall(SVDSetProblemType(svd,SVD_GENERALIZED));
631:     PetscCall(PetscOptionsBoolGroupEnd("-svd_hyperbolic","Hyperbolic singular value decomposition (HSVD)","SVDSetProblemType",&flg));
632:     if (flg) PetscCall(SVDSetProblemType(svd,SVD_HYPERBOLIC));

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

637:     i = svd->max_it;
638:     PetscCall(PetscOptionsInt("-svd_max_it","Maximum number of iterations","SVDSetTolerances",svd->max_it,&i,&flg1));
639:     r = svd->tol;
640:     PetscCall(PetscOptionsReal("-svd_tol","Tolerance","SVDSetTolerances",SlepcDefaultTol(svd->tol),&r,&flg2));
641:     if (flg1 || flg2) PetscCall(SVDSetTolerances(svd,r,i));

643:     PetscCall(PetscOptionsBoolGroupBegin("-svd_conv_abs","Absolute error convergence test","SVDSetConvergenceTest",&flg));
644:     if (flg) PetscCall(SVDSetConvergenceTest(svd,SVD_CONV_ABS));
645:     PetscCall(PetscOptionsBoolGroup("-svd_conv_rel","Relative error convergence test","SVDSetConvergenceTest",&flg));
646:     if (flg) PetscCall(SVDSetConvergenceTest(svd,SVD_CONV_REL));
647:     PetscCall(PetscOptionsBoolGroup("-svd_conv_norm","Convergence test relative to the matrix norms","SVDSetConvergenceTest",&flg));
648:     if (flg) PetscCall(SVDSetConvergenceTest(svd,SVD_CONV_NORM));
649:     PetscCall(PetscOptionsBoolGroup("-svd_conv_maxit","Maximum iterations convergence test","SVDSetConvergenceTest",&flg));
650:     if (flg) PetscCall(SVDSetConvergenceTest(svd,SVD_CONV_MAXIT));
651:     PetscCall(PetscOptionsBoolGroupEnd("-svd_conv_user","User-defined convergence test","SVDSetConvergenceTest",&flg));
652:     if (flg) PetscCall(SVDSetConvergenceTest(svd,SVD_CONV_USER));

654:     PetscCall(PetscOptionsBoolGroupBegin("-svd_stop_basic","Stop iteration if all singular values converged or max_it reached","SVDSetStoppingTest",&flg));
655:     if (flg) PetscCall(SVDSetStoppingTest(svd,SVD_STOP_BASIC));
656:     PetscCall(PetscOptionsBoolGroupEnd("-svd_stop_user","User-defined stopping test","SVDSetStoppingTest",&flg));
657:     if (flg) PetscCall(SVDSetStoppingTest(svd,SVD_STOP_USER));

659:     i = svd->nsv;
660:     PetscCall(PetscOptionsInt("-svd_nsv","Number of singular values to compute","SVDSetDimensions",svd->nsv,&i,&flg1));
661:     j = svd->ncv;
662:     PetscCall(PetscOptionsInt("-svd_ncv","Number of basis vectors","SVDSetDimensions",svd->ncv,&j,&flg2));
663:     k = svd->mpd;
664:     PetscCall(PetscOptionsInt("-svd_mpd","Maximum dimension of projected problem","SVDSetDimensions",svd->mpd,&k,&flg3));
665:     if (flg1 || flg2 || flg3) PetscCall(SVDSetDimensions(svd,i,j,k));

667:     PetscCall(PetscOptionsBoolGroupBegin("-svd_largest","Compute largest singular values","SVDSetWhichSingularTriplets",&flg));
668:     if (flg) PetscCall(SVDSetWhichSingularTriplets(svd,SVD_LARGEST));
669:     PetscCall(PetscOptionsBoolGroupEnd("-svd_smallest","Compute smallest singular values","SVDSetWhichSingularTriplets",&flg));
670:     if (flg) PetscCall(SVDSetWhichSingularTriplets(svd,SVD_SMALLEST));

672:     /* -----------------------------------------------------------------------*/
673:     /*
674:       Cancels all monitors hardwired into code before call to SVDSetFromOptions()
675:     */
676:     PetscCall(PetscOptionsBool("-svd_monitor_cancel","Remove any hardwired monitor routines","SVDMonitorCancel",PETSC_FALSE,&flg,&set));
677:     if (set && flg) PetscCall(SVDMonitorCancel(svd));
678:     PetscCall(SVDMonitorSetFromOptions(svd,"-svd_monitor","first_approximation",NULL,PETSC_FALSE));
679:     PetscCall(SVDMonitorSetFromOptions(svd,"-svd_monitor_all","all_approximations",NULL,PETSC_TRUE));
680:     PetscCall(SVDMonitorSetFromOptions(svd,"-svd_monitor_conv","convergence_history",NULL,PETSC_FALSE));
681:     PetscCall(SVDMonitorSetFromOptions(svd,"-svd_monitor_conditioning","conditioning",NULL,PETSC_FALSE));

683:     /* -----------------------------------------------------------------------*/
684:     PetscCall(PetscOptionsName("-svd_view","Print detailed information on solver used","SVDView",&set));
685:     PetscCall(PetscOptionsName("-svd_view_vectors","View computed singular vectors","SVDVectorsView",&set));
686:     PetscCall(PetscOptionsName("-svd_view_values","View computed singular values","SVDValuesView",&set));
687:     PetscCall(PetscOptionsName("-svd_converged_reason","Print reason for convergence, and number of iterations","SVDConvergedReasonView",&set));
688:     PetscCall(PetscOptionsName("-svd_error_absolute","Print absolute errors of each singular triplet","SVDErrorView",&set));
689:     PetscCall(PetscOptionsName("-svd_error_relative","Print relative errors of each singular triplet","SVDErrorView",&set));
690:     PetscCall(PetscOptionsName("-svd_error_norm","Print errors relative to the matrix norms of each singular triplet","SVDErrorView",&set));

692:     PetscTryTypeMethod(svd,setfromoptions,PetscOptionsObject);
693:     PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)svd,PetscOptionsObject));
694:   PetscOptionsEnd();

696:   if (!svd->V) PetscCall(SVDGetBV(svd,&svd->V,NULL));
697:   PetscCall(BVSetFromOptions(svd->V));
698:   if (!svd->U) PetscCall(SVDGetBV(svd,NULL,&svd->U));
699:   PetscCall(BVSetFromOptions(svd->U));
700:   if (!svd->ds) PetscCall(SVDGetDS(svd,&svd->ds));
701:   PetscCall(SVDSetDSType(svd));
702:   PetscCall(DSSetFromOptions(svd->ds));
703:   PetscFunctionReturn(PETSC_SUCCESS);
704: }

706: /*@
707:    SVDSetProblemType - Specifies the type of the singular value problem.

709:    Logically Collective

711:    Input Parameters:
712: +  svd  - the singular value solver context
713: -  type - a known type of singular value problem

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)

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

724:    Level: intermediate

726: .seealso: SVDSetOperators(), SVDSetSignature(), SVDSetType(), SVDGetProblemType(), SVDProblemType
727: @*/
728: PetscErrorCode SVDSetProblemType(SVD svd,SVDProblemType type)
729: {
730:   PetscFunctionBegin;
733:   if (type == svd->problem_type) PetscFunctionReturn(PETSC_SUCCESS);
734:   switch (type) {
735:     case SVD_STANDARD:
736:       svd->isgeneralized = PETSC_FALSE;
737:       svd->ishyperbolic  = PETSC_FALSE;
738:       break;
739:     case SVD_GENERALIZED:
740:       svd->isgeneralized = PETSC_TRUE;
741:       svd->ishyperbolic  = PETSC_FALSE;
742:       break;
743:     case SVD_HYPERBOLIC:
744:       svd->isgeneralized = PETSC_FALSE;
745:       svd->ishyperbolic  = PETSC_TRUE;
746:       break;
747:     default:
748:       SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_WRONG,"Unknown singular value problem type");
749:   }
750:   svd->problem_type = type;
751:   svd->state = SVD_STATE_INITIAL;
752:   PetscFunctionReturn(PETSC_SUCCESS);
753: }

755: /*@
756:    SVDGetProblemType - Gets the problem type from the SVD object.

758:    Not Collective

760:    Input Parameter:
761: .  svd - the singular value solver context

763:    Output Parameter:
764: .  type - the problem type

766:    Level: intermediate

768: .seealso: SVDSetProblemType(), SVDProblemType
769: @*/
770: PetscErrorCode SVDGetProblemType(SVD svd,SVDProblemType *type)
771: {
772:   PetscFunctionBegin;
774:   PetscAssertPointer(type,2);
775:   *type = svd->problem_type;
776:   PetscFunctionReturn(PETSC_SUCCESS);
777: }

779: /*@
780:    SVDIsGeneralized - Ask if the SVD object corresponds to a generalized
781:    singular value problem.

783:    Not Collective

785:    Input Parameter:
786: .  svd - the singular value solver context

788:    Output Parameter:
789: .  is - the answer

791:    Level: intermediate

793: .seealso: SVDIsHyperbolic()
794: @*/
795: PetscErrorCode SVDIsGeneralized(SVD svd,PetscBool* is)
796: {
797:   PetscFunctionBegin;
799:   PetscAssertPointer(is,2);
800:   *is = svd->isgeneralized;
801:   PetscFunctionReturn(PETSC_SUCCESS);
802: }

804: /*@
805:    SVDIsHyperbolic - Ask if the SVD object corresponds to a hyperbolic
806:    singular value problem.

808:    Not Collective

810:    Input Parameter:
811: .  svd - the singular value solver context

813:    Output Parameter:
814: .  is - the answer

816:    Level: intermediate

818: .seealso: SVDIsGeneralized()
819: @*/
820: PetscErrorCode SVDIsHyperbolic(SVD svd,PetscBool* is)
821: {
822:   PetscFunctionBegin;
824:   PetscAssertPointer(is,2);
825:   *is = svd->ishyperbolic;
826:   PetscFunctionReturn(PETSC_SUCCESS);
827: }

829: /*@
830:    SVDSetTrackAll - Specifies if the solver must compute the residual norm of all
831:    approximate singular value or not.

833:    Logically Collective

835:    Input Parameters:
836: +  svd      - the singular value solver context
837: -  trackall - whether to compute all residuals or not

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.

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

847:    Level: developer

849: .seealso: SVDGetTrackAll()
850: @*/
851: PetscErrorCode SVDSetTrackAll(SVD svd,PetscBool trackall)
852: {
853:   PetscFunctionBegin;
856:   svd->trackall = trackall;
857:   PetscFunctionReturn(PETSC_SUCCESS);
858: }

860: /*@
861:    SVDGetTrackAll - Returns the flag indicating whether all residual norms must
862:    be computed or not.

864:    Not Collective

866:    Input Parameter:
867: .  svd - the singular value solver context

869:    Output Parameter:
870: .  trackall - the returned flag

872:    Level: developer

874: .seealso: SVDSetTrackAll()
875: @*/
876: PetscErrorCode SVDGetTrackAll(SVD svd,PetscBool *trackall)
877: {
878:   PetscFunctionBegin;
880:   PetscAssertPointer(trackall,2);
881:   *trackall = svd->trackall;
882:   PetscFunctionReturn(PETSC_SUCCESS);
883: }

885: /*@
886:    SVDSetOptionsPrefix - Sets the prefix used for searching for all
887:    SVD options in the database.

889:    Logically Collective

891:    Input Parameters:
892: +  svd - the singular value solver context
893: -  prefix - the prefix string to prepend to all SVD option requests

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.

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

907:    Level: advanced

909: .seealso: SVDAppendOptionsPrefix(), SVDGetOptionsPrefix()
910: @*/
911: PetscErrorCode SVDSetOptionsPrefix(SVD svd,const char *prefix)
912: {
913:   PetscFunctionBegin;
915:   if (!svd->V) PetscCall(SVDGetBV(svd,&svd->V,&svd->U));
916:   PetscCall(BVSetOptionsPrefix(svd->V,prefix));
917:   PetscCall(BVSetOptionsPrefix(svd->U,prefix));
918:   if (!svd->ds) PetscCall(SVDGetDS(svd,&svd->ds));
919:   PetscCall(DSSetOptionsPrefix(svd->ds,prefix));
920:   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)svd,prefix));
921:   PetscFunctionReturn(PETSC_SUCCESS);
922: }

924: /*@
925:    SVDAppendOptionsPrefix - Appends to the prefix used for searching for all
926:    SVD options in the database.

928:    Logically Collective

930:    Input Parameters:
931: +  svd - the singular value solver context
932: -  prefix - the prefix string to prepend to all SVD option requests

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.

938:    Level: advanced

940: .seealso: SVDSetOptionsPrefix(), SVDGetOptionsPrefix()
941: @*/
942: PetscErrorCode SVDAppendOptionsPrefix(SVD svd,const char *prefix)
943: {
944:   PetscFunctionBegin;
946:   if (!svd->V) PetscCall(SVDGetBV(svd,&svd->V,&svd->U));
947:   PetscCall(BVAppendOptionsPrefix(svd->V,prefix));
948:   PetscCall(BVAppendOptionsPrefix(svd->U,prefix));
949:   if (!svd->ds) PetscCall(SVDGetDS(svd,&svd->ds));
950:   PetscCall(DSAppendOptionsPrefix(svd->ds,prefix));
951:   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)svd,prefix));
952:   PetscFunctionReturn(PETSC_SUCCESS);
953: }

955: /*@
956:    SVDGetOptionsPrefix - Gets the prefix used for searching for all
957:    SVD options in the database.

959:    Not Collective

961:    Input Parameters:
962: .  svd - the singular value solver context

964:    Output Parameters:
965: .  prefix - pointer to the prefix string used is returned

967:    Note:
968:    On the Fortran side, the user should pass in a string 'prefix' of
969:    sufficient length to hold the prefix.

971:    Level: advanced

973: .seealso: SVDSetOptionsPrefix(), SVDAppendOptionsPrefix()
974: @*/
975: PetscErrorCode SVDGetOptionsPrefix(SVD svd,const char *prefix[])
976: {
977:   PetscFunctionBegin;
979:   PetscAssertPointer(prefix,2);
980:   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)svd,prefix));
981:   PetscFunctionReturn(PETSC_SUCCESS);
982: }