Actual source code: svdopts.c

slepc-3.21.0 2024-03-30
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_DEFAULT for either argument to assign a reasonably good value.

 99:    Level: intermediate

101: .seealso: SVDGetTolerances()
102: @*/
103: PetscErrorCode SVDSetTolerances(SVD svd,PetscReal tol,PetscInt maxits)
104: {
105:   PetscFunctionBegin;
109:   if (tol == (PetscReal)PETSC_DEFAULT) {
110:     svd->tol   = PETSC_DEFAULT;
111:     svd->state = SVD_STATE_INITIAL;
112:   } else {
113:     PetscCheck(tol>0.0,PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of tol. Must be > 0");
114:     svd->tol = tol;
115:   }
116:   if (maxits == PETSC_DEFAULT || maxits == PETSC_DECIDE) {
117:     svd->max_it = PETSC_DEFAULT;
118:     svd->state  = SVD_STATE_INITIAL;
119:   } else {
120:     PetscCheck(maxits>0,PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of maxits. Must be > 0");
121:     svd->max_it = maxits;
122:   }
123:   PetscFunctionReturn(PETSC_SUCCESS);
124: }

126: /*@C
127:    SVDGetTolerances - Gets the tolerance and maximum
128:    iteration count used by the default SVD convergence tests.

130:    Not Collective

132:    Input Parameter:
133: .  svd - the singular value solver context

135:    Output Parameters:
136: +  tol - the convergence tolerance
137: -  maxits - maximum number of iterations

139:    Notes:
140:    The user can specify NULL for any parameter that is not needed.

142:    Level: intermediate

144: .seealso: SVDSetTolerances()
145: @*/
146: PetscErrorCode SVDGetTolerances(SVD svd,PetscReal *tol,PetscInt *maxits)
147: {
148:   PetscFunctionBegin;
150:   if (tol)    *tol    = svd->tol;
151:   if (maxits) *maxits = svd->max_it;
152:   PetscFunctionReturn(PETSC_SUCCESS);
153: }

155: /*@
156:    SVDSetDimensions - Sets the number of singular values to compute
157:    and the dimension of the subspace.

159:    Logically Collective

161:    Input Parameters:
162: +  svd - the singular value solver context
163: .  nsv - number of singular values to compute
164: .  ncv - the maximum dimension of the subspace to be used by the solver
165: -  mpd - the maximum dimension allowed for the projected problem

167:    Options Database Keys:
168: +  -svd_nsv <nsv> - Sets the number of singular values
169: .  -svd_ncv <ncv> - Sets the dimension of the subspace
170: -  -svd_mpd <mpd> - Sets the maximum projected dimension

172:    Notes:
173:    Use PETSC_DEFAULT for ncv and mpd to assign a reasonably good value, which is
174:    dependent on the solution method and the number of singular values required.

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

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

185:    Level: intermediate

187: .seealso: SVDGetDimensions()
188: @*/
189: PetscErrorCode SVDSetDimensions(SVD svd,PetscInt nsv,PetscInt ncv,PetscInt mpd)
190: {
191:   PetscFunctionBegin;
196:   PetscCheck(nsv>0,PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of nsv. Must be > 0");
197:   svd->nsv = nsv;
198:   if (ncv == PETSC_DEFAULT || ncv == PETSC_DECIDE) {
199:     svd->ncv = PETSC_DEFAULT;
200:   } else {
201:     PetscCheck(ncv>0,PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of ncv. Must be > 0");
202:     svd->ncv = ncv;
203:   }
204:   if (mpd == PETSC_DECIDE || mpd == PETSC_DEFAULT) {
205:     svd->mpd = PETSC_DEFAULT;
206:   } else {
207:     PetscCheck(mpd>0,PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of mpd. Must be > 0");
208:     svd->mpd = mpd;
209:   }
210:   svd->state = SVD_STATE_INITIAL;
211:   PetscFunctionReturn(PETSC_SUCCESS);
212: }

214: /*@C
215:    SVDGetDimensions - Gets the number of singular values to compute
216:    and the dimension of the subspace.

218:    Not Collective

220:    Input Parameter:
221: .  svd - the singular value context

223:    Output Parameters:
224: +  nsv - number of singular values to compute
225: .  ncv - the maximum dimension of the subspace to be used by the solver
226: -  mpd - the maximum dimension allowed for the projected problem

228:    Notes:
229:    The user can specify NULL for any parameter that is not needed.

231:    Level: intermediate

233: .seealso: SVDSetDimensions()
234: @*/
235: PetscErrorCode SVDGetDimensions(SVD svd,PetscInt *nsv,PetscInt *ncv,PetscInt *mpd)
236: {
237:   PetscFunctionBegin;
239:   if (nsv) *nsv = svd->nsv;
240:   if (ncv) *ncv = svd->ncv;
241:   if (mpd) *mpd = svd->mpd;
242:   PetscFunctionReturn(PETSC_SUCCESS);
243: }

245: /*@
246:     SVDSetWhichSingularTriplets - Specifies which singular triplets are
247:     to be sought.

249:     Logically Collective

251:     Input Parameter:
252: .   svd - singular value solver context obtained from SVDCreate()

254:     Output Parameter:
255: .   which - which singular triplets are to be sought

257:     Options Database Keys:
258: +   -svd_largest  - Sets largest singular values
259: -   -svd_smallest - Sets smallest singular values

261:     Notes:
262:     The parameter 'which' can have one of these values

264: +     SVD_LARGEST  - largest singular values
265: -     SVD_SMALLEST - smallest singular values

267:     Level: intermediate

269: .seealso: SVDGetWhichSingularTriplets(), SVDWhich
270: @*/
271: PetscErrorCode SVDSetWhichSingularTriplets(SVD svd,SVDWhich which)
272: {
273:   PetscFunctionBegin;
276:   switch (which) {
277:     case SVD_LARGEST:
278:     case SVD_SMALLEST:
279:       if (svd->which != which) {
280:         svd->state = SVD_STATE_INITIAL;
281:         svd->which = which;
282:       }
283:       break;
284:   default:
285:     SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'which' parameter");
286:   }
287:   PetscFunctionReturn(PETSC_SUCCESS);
288: }

290: /*@
291:     SVDGetWhichSingularTriplets - Returns which singular triplets are
292:     to be sought.

294:     Not Collective

296:     Input Parameter:
297: .   svd - singular value solver context obtained from SVDCreate()

299:     Output Parameter:
300: .   which - which singular triplets are to be sought

302:     Notes:
303:     See SVDSetWhichSingularTriplets() for possible values of which

305:     Level: intermediate

307: .seealso: SVDSetWhichSingularTriplets(), SVDWhich
308: @*/
309: PetscErrorCode SVDGetWhichSingularTriplets(SVD svd,SVDWhich *which)
310: {
311:   PetscFunctionBegin;
313:   PetscAssertPointer(which,2);
314:   *which = svd->which;
315:   PetscFunctionReturn(PETSC_SUCCESS);
316: }

318: /*@C
319:    SVDSetConvergenceTestFunction - Sets a function to compute the error estimate
320:    used in the convergence test.

322:    Logically Collective

324:    Input Parameters:
325: +  svd     - singular value solver context obtained from SVDCreate()
326: .  conv    - a pointer to the convergence test function
327: .  ctx     - context for private data for the convergence routine (may be null)
328: -  destroy - a routine for destroying the context (may be null)

330:    Calling sequence of conv:
331: $  PetscErrorCode conv(SVD svd,PetscReal sigma,PetscReal res,PetscReal *errest,void *ctx)
332: +   svd    - singular value solver context obtained from SVDCreate()
333: .   sigma  - computed singular value
334: .   res    - residual norm associated to the singular triplet
335: .   errest - (output) computed error estimate
336: -   ctx    - optional context, as set by SVDSetConvergenceTestFunction()

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,PetscErrorCode (*conv)(SVD svd,PetscReal sigma,PetscReal res,PetscReal *errest,void *ctx),void* ctx,PetscErrorCode (*destroy)(void*))
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    - pointer to the stopping test function
451: .  ctx     - context for private data for the stopping routine (may be null)
452: -  destroy - a routine for destroying the context (may be null)

454:    Calling sequence of stop:
455: $  PetscErrorCode stop(SVD svd,PetscInt its,PetscInt max_it,PetscInt nconv,PetscInt nsv,SVDConvergedReason *reason,void *ctx)
456: +   svd    - singular value solver context obtained from SVDCreate()
457: .   its    - current number of iterations
458: .   max_it - maximum number of iterations
459: .   nconv  - number of currently converged singular triplets
460: .   nsv    - number of requested singular triplets
461: .   reason - (output) result of the stopping test
462: -   ctx    - optional context, as set by SVDSetStoppingTestFunction()

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

470:    Level: advanced

472: .seealso: SVDSetStoppingTest(), SVDStoppingBasic()
473: @*/
474: PetscErrorCode SVDSetStoppingTestFunction(SVD svd,PetscErrorCode (*stop)(SVD svd,PetscInt its,PetscInt max_it,PetscInt nconv,PetscInt nsv,SVDConvergedReason *reason,void *ctx),void* ctx,PetscErrorCode (*destroy)(void*))
475: {
476:   PetscFunctionBegin;
478:   if (svd->stoppingdestroy) PetscCall((*svd->stoppingdestroy)(svd->stoppingctx));
479:   svd->stoppinguser    = stop;
480:   svd->stoppingdestroy = destroy;
481:   svd->stoppingctx     = ctx;
482:   if (stop == SVDStoppingBasic) svd->stop = SVD_STOP_BASIC;
483:   else {
484:     svd->stop     = SVD_STOP_USER;
485:     svd->stopping = svd->stoppinguser;
486:   }
487:   PetscFunctionReturn(PETSC_SUCCESS);
488: }

490: /*@
491:    SVDSetStoppingTest - Specifies how to decide the termination of the outer
492:    loop of the singular value solver.

494:    Logically Collective

496:    Input Parameters:
497: +  svd  - singular value solver context obtained from SVDCreate()
498: -  stop - the type of stopping test

500:    Options Database Keys:
501: +  -svd_stop_basic - Sets the default stopping test
502: -  -svd_stop_user  - Selects the user-defined stopping test

504:    Note:
505:    The parameter 'stop' can have one of these values
506: +     SVD_STOP_BASIC - default stopping test
507: -     SVD_STOP_USER  - function set by SVDSetStoppingTestFunction()

509:    Level: advanced

511: .seealso: SVDGetStoppingTest(), SVDSetStoppingTestFunction(), SVDSetConvergenceTest(), SVDStop
512: @*/
513: PetscErrorCode SVDSetStoppingTest(SVD svd,SVDStop stop)
514: {
515:   PetscFunctionBegin;
518:   switch (stop) {
519:     case SVD_STOP_BASIC: svd->stopping = SVDStoppingBasic; break;
520:     case SVD_STOP_USER:
521:       PetscCheck(svd->stoppinguser,PetscObjectComm((PetscObject)svd),PETSC_ERR_ORDER,"Must call SVDSetStoppingTestFunction() first");
522:       svd->stopping = svd->stoppinguser;
523:       break;
524:     default:
525:       SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'stop' value");
526:   }
527:   svd->stop = stop;
528:   PetscFunctionReturn(PETSC_SUCCESS);
529: }

531: /*@
532:    SVDGetStoppingTest - Gets the method used to decide the termination of the outer
533:    loop of the singular value solver.

535:    Not Collective

537:    Input Parameters:
538: .  svd   - singular value solver context obtained from SVDCreate()

540:    Output Parameters:
541: .  stop  - the type of stopping test

543:    Level: advanced

545: .seealso: SVDSetStoppingTest(), SVDStop
546: @*/
547: PetscErrorCode SVDGetStoppingTest(SVD svd,SVDStop *stop)
548: {
549:   PetscFunctionBegin;
551:   PetscAssertPointer(stop,2);
552:   *stop = svd->stop;
553:   PetscFunctionReturn(PETSC_SUCCESS);
554: }

556: /*@C
557:    SVDMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type
558:    indicated by the user.

560:    Collective

562:    Input Parameters:
563: +  svd      - the singular value solver context
564: .  opt      - the command line option for this monitor
565: .  name     - the monitor type one is seeking
566: .  ctx      - an optional user context for the monitor, or NULL
567: -  trackall - whether this monitor tracks all singular values or not

569:    Level: developer

571: .seealso: SVDMonitorSet(), SVDSetTrackAll()
572: @*/
573: PetscErrorCode SVDMonitorSetFromOptions(SVD svd,const char opt[],const char name[],void *ctx,PetscBool trackall)
574: {
575:   PetscErrorCode       (*mfunc)(SVD,PetscInt,PetscInt,PetscReal*,PetscReal*,PetscInt,void*);
576:   PetscErrorCode       (*cfunc)(PetscViewer,PetscViewerFormat,void*,PetscViewerAndFormat**);
577:   PetscErrorCode       (*dfunc)(PetscViewerAndFormat**);
578:   PetscViewerAndFormat *vf;
579:   PetscViewer          viewer;
580:   PetscViewerFormat    format;
581:   PetscViewerType      vtype;
582:   char                 key[PETSC_MAX_PATH_LEN];
583:   PetscBool            flg;

585:   PetscFunctionBegin;
586:   PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)svd),((PetscObject)svd)->options,((PetscObject)svd)->prefix,opt,&viewer,&format,&flg));
587:   if (!flg) PetscFunctionReturn(PETSC_SUCCESS);

589:   PetscCall(PetscViewerGetType(viewer,&vtype));
590:   PetscCall(SlepcMonitorMakeKey_Internal(name,vtype,format,key));
591:   PetscCall(PetscFunctionListFind(SVDMonitorList,key,&mfunc));
592:   PetscCheck(mfunc,PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"Specified viewer and format not supported");
593:   PetscCall(PetscFunctionListFind(SVDMonitorCreateList,key,&cfunc));
594:   PetscCall(PetscFunctionListFind(SVDMonitorDestroyList,key,&dfunc));
595:   if (!cfunc) cfunc = PetscViewerAndFormatCreate_Internal;
596:   if (!dfunc) dfunc = PetscViewerAndFormatDestroy;

598:   PetscCall((*cfunc)(viewer,format,ctx,&vf));
599:   PetscCall(PetscOptionsRestoreViewer(&viewer));
600:   PetscCall(SVDMonitorSet(svd,mfunc,vf,(PetscErrorCode(*)(void **))dfunc));
601:   if (trackall) PetscCall(SVDSetTrackAll(svd,PETSC_TRUE));
602:   PetscFunctionReturn(PETSC_SUCCESS);
603: }

605: /*@
606:    SVDSetFromOptions - Sets SVD options from the options database.
607:    This routine must be called before SVDSetUp() if the user is to be
608:    allowed to set the solver type.

610:    Collective

612:    Input Parameters:
613: .  svd - the singular value solver context

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

618:    Level: beginner

620: .seealso: SVDSetOptionsPrefix()
621: @*/
622: PetscErrorCode SVDSetFromOptions(SVD svd)
623: {
624:   char           type[256];
625:   PetscBool      set,flg,val,flg1,flg2,flg3;
626:   PetscInt       i,j,k;
627:   PetscReal      r;

629:   PetscFunctionBegin;
631:   PetscCall(SVDRegisterAll());
632:   PetscObjectOptionsBegin((PetscObject)svd);
633:     PetscCall(PetscOptionsFList("-svd_type","SVD solver method","SVDSetType",SVDList,(char*)(((PetscObject)svd)->type_name?((PetscObject)svd)->type_name:SVDCROSS),type,sizeof(type),&flg));
634:     if (flg) PetscCall(SVDSetType(svd,type));
635:     else if (!((PetscObject)svd)->type_name) PetscCall(SVDSetType(svd,SVDCROSS));

637:     PetscCall(PetscOptionsBoolGroupBegin("-svd_standard","Singular value decomposition (SVD)","SVDSetProblemType",&flg));
638:     if (flg) PetscCall(SVDSetProblemType(svd,SVD_STANDARD));
639:     PetscCall(PetscOptionsBoolGroup("-svd_generalized","Generalized singular value decomposition (GSVD)","SVDSetProblemType",&flg));
640:     if (flg) PetscCall(SVDSetProblemType(svd,SVD_GENERALIZED));
641:     PetscCall(PetscOptionsBoolGroupEnd("-svd_hyperbolic","Hyperbolic singular value decomposition (HSVD)","SVDSetProblemType",&flg));
642:     if (flg) PetscCall(SVDSetProblemType(svd,SVD_HYPERBOLIC));

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

647:     i = svd->max_it;
648:     PetscCall(PetscOptionsInt("-svd_max_it","Maximum number of iterations","SVDSetTolerances",svd->max_it,&i,&flg1));
649:     r = svd->tol;
650:     PetscCall(PetscOptionsReal("-svd_tol","Tolerance","SVDSetTolerances",SlepcDefaultTol(svd->tol),&r,&flg2));
651:     if (flg1 || flg2) PetscCall(SVDSetTolerances(svd,r,i));

653:     PetscCall(PetscOptionsBoolGroupBegin("-svd_conv_abs","Absolute error convergence test","SVDSetConvergenceTest",&flg));
654:     if (flg) PetscCall(SVDSetConvergenceTest(svd,SVD_CONV_ABS));
655:     PetscCall(PetscOptionsBoolGroup("-svd_conv_rel","Relative error convergence test","SVDSetConvergenceTest",&flg));
656:     if (flg) PetscCall(SVDSetConvergenceTest(svd,SVD_CONV_REL));
657:     PetscCall(PetscOptionsBoolGroup("-svd_conv_norm","Convergence test relative to the matrix norms","SVDSetConvergenceTest",&flg));
658:     if (flg) PetscCall(SVDSetConvergenceTest(svd,SVD_CONV_NORM));
659:     PetscCall(PetscOptionsBoolGroup("-svd_conv_maxit","Maximum iterations convergence test","SVDSetConvergenceTest",&flg));
660:     if (flg) PetscCall(SVDSetConvergenceTest(svd,SVD_CONV_MAXIT));
661:     PetscCall(PetscOptionsBoolGroupEnd("-svd_conv_user","User-defined convergence test","SVDSetConvergenceTest",&flg));
662:     if (flg) PetscCall(SVDSetConvergenceTest(svd,SVD_CONV_USER));

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

669:     i = svd->nsv;
670:     PetscCall(PetscOptionsInt("-svd_nsv","Number of singular values to compute","SVDSetDimensions",svd->nsv,&i,&flg1));
671:     j = svd->ncv;
672:     PetscCall(PetscOptionsInt("-svd_ncv","Number of basis vectors","SVDSetDimensions",svd->ncv,&j,&flg2));
673:     k = svd->mpd;
674:     PetscCall(PetscOptionsInt("-svd_mpd","Maximum dimension of projected problem","SVDSetDimensions",svd->mpd,&k,&flg3));
675:     if (flg1 || flg2 || flg3) PetscCall(SVDSetDimensions(svd,i,j,k));

677:     PetscCall(PetscOptionsBoolGroupBegin("-svd_largest","Compute largest singular values","SVDSetWhichSingularTriplets",&flg));
678:     if (flg) PetscCall(SVDSetWhichSingularTriplets(svd,SVD_LARGEST));
679:     PetscCall(PetscOptionsBoolGroupEnd("-svd_smallest","Compute smallest singular values","SVDSetWhichSingularTriplets",&flg));
680:     if (flg) PetscCall(SVDSetWhichSingularTriplets(svd,SVD_SMALLEST));

682:     /* -----------------------------------------------------------------------*/
683:     /*
684:       Cancels all monitors hardwired into code before call to SVDSetFromOptions()
685:     */
686:     PetscCall(PetscOptionsBool("-svd_monitor_cancel","Remove any hardwired monitor routines","SVDMonitorCancel",PETSC_FALSE,&flg,&set));
687:     if (set && flg) PetscCall(SVDMonitorCancel(svd));
688:     PetscCall(SVDMonitorSetFromOptions(svd,"-svd_monitor","first_approximation",NULL,PETSC_FALSE));
689:     PetscCall(SVDMonitorSetFromOptions(svd,"-svd_monitor_all","all_approximations",NULL,PETSC_TRUE));
690:     PetscCall(SVDMonitorSetFromOptions(svd,"-svd_monitor_conv","convergence_history",NULL,PETSC_FALSE));
691:     PetscCall(SVDMonitorSetFromOptions(svd,"-svd_monitor_conditioning","conditioning",NULL,PETSC_FALSE));

693:     /* -----------------------------------------------------------------------*/
694:     PetscCall(PetscOptionsName("-svd_view","Print detailed information on solver used","SVDView",&set));
695:     PetscCall(PetscOptionsName("-svd_view_vectors","View computed singular vectors","SVDVectorsView",&set));
696:     PetscCall(PetscOptionsName("-svd_view_values","View computed singular values","SVDValuesView",&set));
697:     PetscCall(PetscOptionsName("-svd_converged_reason","Print reason for convergence, and number of iterations","SVDConvergedReasonView",&set));
698:     PetscCall(PetscOptionsName("-svd_error_absolute","Print absolute errors of each singular triplet","SVDErrorView",&set));
699:     PetscCall(PetscOptionsName("-svd_error_relative","Print relative errors of each singular triplet","SVDErrorView",&set));
700:     PetscCall(PetscOptionsName("-svd_error_norm","Print errors relative to the matrix norms of each singular triplet","SVDErrorView",&set));

702:     PetscTryTypeMethod(svd,setfromoptions,PetscOptionsObject);
703:     PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)svd,PetscOptionsObject));
704:   PetscOptionsEnd();

706:   if (!svd->V) PetscCall(SVDGetBV(svd,&svd->V,NULL));
707:   PetscCall(BVSetFromOptions(svd->V));
708:   if (!svd->U) PetscCall(SVDGetBV(svd,NULL,&svd->U));
709:   PetscCall(BVSetFromOptions(svd->U));
710:   if (!svd->ds) PetscCall(SVDGetDS(svd,&svd->ds));
711:   PetscCall(SVDSetDSType(svd));
712:   PetscCall(DSSetFromOptions(svd->ds));
713:   PetscFunctionReturn(PETSC_SUCCESS);
714: }

716: /*@
717:    SVDSetProblemType - Specifies the type of the singular value problem.

719:    Logically Collective

721:    Input Parameters:
722: +  svd  - the singular value solver context
723: -  type - a known type of singular value problem

725:    Options Database Keys:
726: +  -svd_standard    - standard singular value decomposition (SVD)
727: .  -svd_generalized - generalized singular value problem (GSVD)
728: -  -svd_hyperbolic  - hyperbolic singular value problem (HSVD)

730:    Notes:
731:    The GSVD requires that two matrices have been passed via SVDSetOperators().
732:    The HSVD requires that a signature matrix has been passed via SVDSetSignature().

734:    Level: intermediate

736: .seealso: SVDSetOperators(), SVDSetSignature(), SVDSetType(), SVDGetProblemType(), SVDProblemType
737: @*/
738: PetscErrorCode SVDSetProblemType(SVD svd,SVDProblemType type)
739: {
740:   PetscFunctionBegin;
743:   if (type == svd->problem_type) PetscFunctionReturn(PETSC_SUCCESS);
744:   switch (type) {
745:     case SVD_STANDARD:
746:       svd->isgeneralized = PETSC_FALSE;
747:       svd->ishyperbolic  = PETSC_FALSE;
748:       break;
749:     case SVD_GENERALIZED:
750:       svd->isgeneralized = PETSC_TRUE;
751:       svd->ishyperbolic  = PETSC_FALSE;
752:       break;
753:     case SVD_HYPERBOLIC:
754:       svd->isgeneralized = PETSC_FALSE;
755:       svd->ishyperbolic  = PETSC_TRUE;
756:       break;
757:     default:
758:       SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_WRONG,"Unknown singular value problem type");
759:   }
760:   svd->problem_type = type;
761:   svd->state = SVD_STATE_INITIAL;
762:   PetscFunctionReturn(PETSC_SUCCESS);
763: }

765: /*@
766:    SVDGetProblemType - Gets the problem type from the SVD object.

768:    Not Collective

770:    Input Parameter:
771: .  svd - the singular value solver context

773:    Output Parameter:
774: .  type - the problem type

776:    Level: intermediate

778: .seealso: SVDSetProblemType(), SVDProblemType
779: @*/
780: PetscErrorCode SVDGetProblemType(SVD svd,SVDProblemType *type)
781: {
782:   PetscFunctionBegin;
784:   PetscAssertPointer(type,2);
785:   *type = svd->problem_type;
786:   PetscFunctionReturn(PETSC_SUCCESS);
787: }

789: /*@
790:    SVDIsGeneralized - Ask if the SVD object corresponds to a generalized
791:    singular value problem.

793:    Not Collective

795:    Input Parameter:
796: .  svd - the singular value solver context

798:    Output Parameter:
799: .  is - the answer

801:    Level: intermediate

803: .seealso: SVDIsHyperbolic()
804: @*/
805: PetscErrorCode SVDIsGeneralized(SVD svd,PetscBool* is)
806: {
807:   PetscFunctionBegin;
809:   PetscAssertPointer(is,2);
810:   *is = svd->isgeneralized;
811:   PetscFunctionReturn(PETSC_SUCCESS);
812: }

814: /*@
815:    SVDIsHyperbolic - Ask if the SVD object corresponds to a hyperbolic
816:    singular value problem.

818:    Not Collective

820:    Input Parameter:
821: .  svd - the singular value solver context

823:    Output Parameter:
824: .  is - the answer

826:    Level: intermediate

828: .seealso: SVDIsGeneralized()
829: @*/
830: PetscErrorCode SVDIsHyperbolic(SVD svd,PetscBool* is)
831: {
832:   PetscFunctionBegin;
834:   PetscAssertPointer(is,2);
835:   *is = svd->ishyperbolic;
836:   PetscFunctionReturn(PETSC_SUCCESS);
837: }

839: /*@
840:    SVDSetTrackAll - Specifies if the solver must compute the residual norm of all
841:    approximate singular value or not.

843:    Logically Collective

845:    Input Parameters:
846: +  svd      - the singular value solver context
847: -  trackall - whether to compute all residuals or not

849:    Notes:
850:    If the user sets trackall=PETSC_TRUE then the solver computes (or estimates)
851:    the residual norm for each singular value approximation. Computing the residual is
852:    usually an expensive operation and solvers commonly compute only the residual
853:    associated to the first unconverged singular value.

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

857:    Level: developer

859: .seealso: SVDGetTrackAll()
860: @*/
861: PetscErrorCode SVDSetTrackAll(SVD svd,PetscBool trackall)
862: {
863:   PetscFunctionBegin;
866:   svd->trackall = trackall;
867:   PetscFunctionReturn(PETSC_SUCCESS);
868: }

870: /*@
871:    SVDGetTrackAll - Returns the flag indicating whether all residual norms must
872:    be computed or not.

874:    Not Collective

876:    Input Parameter:
877: .  svd - the singular value solver context

879:    Output Parameter:
880: .  trackall - the returned flag

882:    Level: developer

884: .seealso: SVDSetTrackAll()
885: @*/
886: PetscErrorCode SVDGetTrackAll(SVD svd,PetscBool *trackall)
887: {
888:   PetscFunctionBegin;
890:   PetscAssertPointer(trackall,2);
891:   *trackall = svd->trackall;
892:   PetscFunctionReturn(PETSC_SUCCESS);
893: }

895: /*@C
896:    SVDSetOptionsPrefix - Sets the prefix used for searching for all
897:    SVD options in the database.

899:    Logically Collective

901:    Input Parameters:
902: +  svd - the singular value solver context
903: -  prefix - the prefix string to prepend to all SVD option requests

905:    Notes:
906:    A hyphen (-) must NOT be given at the beginning of the prefix name.
907:    The first character of all runtime options is AUTOMATICALLY the
908:    hyphen.

910:    For example, to distinguish between the runtime options for two
911:    different SVD contexts, one could call
912: .vb
913:       SVDSetOptionsPrefix(svd1,"svd1_")
914:       SVDSetOptionsPrefix(svd2,"svd2_")
915: .ve

917:    Level: advanced

919: .seealso: SVDAppendOptionsPrefix(), SVDGetOptionsPrefix()
920: @*/
921: PetscErrorCode SVDSetOptionsPrefix(SVD svd,const char *prefix)
922: {
923:   PetscFunctionBegin;
925:   if (!svd->V) PetscCall(SVDGetBV(svd,&svd->V,&svd->U));
926:   PetscCall(BVSetOptionsPrefix(svd->V,prefix));
927:   PetscCall(BVSetOptionsPrefix(svd->U,prefix));
928:   if (!svd->ds) PetscCall(SVDGetDS(svd,&svd->ds));
929:   PetscCall(DSSetOptionsPrefix(svd->ds,prefix));
930:   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)svd,prefix));
931:   PetscFunctionReturn(PETSC_SUCCESS);
932: }

934: /*@C
935:    SVDAppendOptionsPrefix - Appends to the prefix used for searching for all
936:    SVD options in the database.

938:    Logically Collective

940:    Input Parameters:
941: +  svd - the singular value solver context
942: -  prefix - the prefix string to prepend to all SVD option requests

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

948:    Level: advanced

950: .seealso: SVDSetOptionsPrefix(), SVDGetOptionsPrefix()
951: @*/
952: PetscErrorCode SVDAppendOptionsPrefix(SVD svd,const char *prefix)
953: {
954:   PetscFunctionBegin;
956:   if (!svd->V) PetscCall(SVDGetBV(svd,&svd->V,&svd->U));
957:   PetscCall(BVAppendOptionsPrefix(svd->V,prefix));
958:   PetscCall(BVAppendOptionsPrefix(svd->U,prefix));
959:   if (!svd->ds) PetscCall(SVDGetDS(svd,&svd->ds));
960:   PetscCall(DSAppendOptionsPrefix(svd->ds,prefix));
961:   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)svd,prefix));
962:   PetscFunctionReturn(PETSC_SUCCESS);
963: }

965: /*@C
966:    SVDGetOptionsPrefix - Gets the prefix used for searching for all
967:    SVD options in the database.

969:    Not Collective

971:    Input Parameters:
972: .  svd - the singular value solver context

974:    Output Parameters:
975: .  prefix - pointer to the prefix string used is returned

977:    Note:
978:    On the Fortran side, the user should pass in a string 'prefix' of
979:    sufficient length to hold the prefix.

981:    Level: advanced

983: .seealso: SVDSetOptionsPrefix(), SVDAppendOptionsPrefix()
984: @*/
985: PetscErrorCode SVDGetOptionsPrefix(SVD svd,const char *prefix[])
986: {
987:   PetscFunctionBegin;
989:   PetscAssertPointer(prefix,2);
990:   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)svd,prefix));
991:   PetscFunctionReturn(PETSC_SUCCESS);
992: }