Actual source code: svdopts.c

slepc-3.15.1 2021-05-28
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2021, 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 on svd

 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: {
 49:   if (svd->impltrans!=impl) {
 50:     svd->impltrans = impl;
 51:     svd->state     = SVD_STATE_INITIAL;
 52:   }
 53:   return(0);
 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: {
 77:   *impl = svd->impltrans;
 78:   return(0);
 79: }

 81: /*@
 82:    SVDSetTolerances - Sets the tolerance and maximum
 83:    iteration count used by the default SVD convergence testers.

 85:    Logically Collective on svd

 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: {
109:   if (tol == PETSC_DEFAULT) {
110:     svd->tol   = PETSC_DEFAULT;
111:     svd->state = SVD_STATE_INITIAL;
112:   } else {
113:     if (tol <= 0.0) SETERRQ(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:     if (maxits <= 0) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of maxits. Must be > 0");
121:     svd->max_it = maxits;
122:   }
123:   return(0);
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: {
150:   if (tol)    *tol    = svd->tol;
151:   if (maxits) *maxits = svd->max_it;
152:   return(0);
153: }

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

159:    Logically Collective on svd

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: {
196:   if (nsv<1) SETERRQ(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:     if (ncv<1) SETERRQ(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:     if (mpd<1) SETERRQ(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:   return(0);
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: {
239:   if (nsv) *nsv = svd->nsv;
240:   if (ncv) *ncv = svd->ncv;
241:   if (mpd) *mpd = svd->mpd;
242:   return(0);
243: }

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

249:     Logically Collective on svd

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:     Possible values:
258:     The parameter 'which' can have one of these values

260: +     SVD_LARGEST  - largest singular values
261: -     SVD_SMALLEST - smallest singular values

263:     Options Database Keys:
264: +   -svd_largest  - Sets largest singular values
265: -   -svd_smallest - Sets smallest singular values

267:     Level: intermediate

269: .seealso: SVDGetWhichSingularTriplets(), SVDWhich
270: @*/
271: PetscErrorCode SVDSetWhichSingularTriplets(SVD svd,SVDWhich which)
272: {
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:   return(0);
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: {
314:   *which = svd->which;
315:   return(0);
316: }

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

322:    Logically Collective on svd

324:    Input Parameters:
325: +  svd     - singular value solver context obtained from SVDCreate()
326: .  func    - 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 func:
331: $   func(SVD svd,PetscReal sigma,PetscReal res,PetscReal *errest,void *ctx)

333: +   svd    - singular value solver context obtained from SVDCreate()
334: .   sigma  - computed singular value
335: .   res    - residual norm associated to the singular triplet
336: .   errest - (output) computed error estimate
337: -   ctx    - optional context, as set by SVDSetConvergenceTestFunction()

339:    Note:
340:    If the error estimate returned by the convergence test function is less than
341:    the tolerance, then the singular value is accepted as converged.

343:    Level: advanced

345: .seealso: SVDSetConvergenceTest(), SVDSetTolerances()
346: @*/
347: PetscErrorCode SVDSetConvergenceTestFunction(SVD svd,PetscErrorCode (*func)(SVD,PetscReal,PetscReal,PetscReal*,void*),void* ctx,PetscErrorCode (*destroy)(void*))
348: {

353:   if (svd->convergeddestroy) {
354:     (*svd->convergeddestroy)(svd->convergedctx);
355:   }
356:   svd->convergeduser    = func;
357:   svd->convergeddestroy = destroy;
358:   svd->convergedctx     = ctx;
359:   if (func == SVDConvergedAbsolute) svd->conv = SVD_CONV_ABS;
360:   else if (func == SVDConvergedRelative) svd->conv = SVD_CONV_REL;
361:   else if (func == SVDConvergedMaxIt) svd->conv = SVD_CONV_MAXIT;
362:   else {
363:     svd->conv      = SVD_CONV_USER;
364:     svd->converged = svd->convergeduser;
365:   }
366:   return(0);
367: }

369: /*@
370:    SVDSetConvergenceTest - Specifies how to compute the error estimate
371:    used in the convergence test.

373:    Logically Collective on svd

375:    Input Parameters:
376: +  svd  - singular value solver context obtained from SVDCreate()
377: -  conv - the type of convergence test

379:    Options Database Keys:
380: +  -svd_conv_abs   - Sets the absolute convergence test
381: .  -svd_conv_rel   - Sets the convergence test relative to the singular value
382: .  -svd_conv_maxit - Forces the maximum number of iterations as set by -svd_max_it
383: -  -svd_conv_user  - Selects the user-defined convergence test

385:    Note:
386:    The parameter 'conv' can have one of these values
387: +     SVD_CONV_ABS   - absolute error ||r||
388: .     SVD_CONV_REL   - error relative to the singular value l, ||r||/sigma
389: .     SVD_CONV_MAXIT - no convergence until maximum number of iterations has been reached
390: -     SVD_CONV_USER  - function set by SVDSetConvergenceTestFunction()

392:    Level: intermediate

394: .seealso: SVDGetConvergenceTest(), SVDSetConvergenceTestFunction(), SVDSetStoppingTest(), SVDConv
395: @*/
396: PetscErrorCode SVDSetConvergenceTest(SVD svd,SVDConv conv)
397: {
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_MAXIT: svd->converged = SVDConvergedMaxIt; break;
405:     case SVD_CONV_USER:
406:       if (!svd->convergeduser) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ORDER,"Must call SVDSetConvergenceTestFunction() first");
407:       svd->converged = svd->convergeduser;
408:       break;
409:     default:
410:       SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'conv' value");
411:   }
412:   svd->conv = conv;
413:   return(0);
414: }

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

420:    Not Collective

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

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

428:    Level: intermediate

430: .seealso: SVDSetConvergenceTest(), SVDConv
431: @*/
432: PetscErrorCode SVDGetConvergenceTest(SVD svd,SVDConv *conv)
433: {
437:   *conv = svd->conv;
438:   return(0);
439: }

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

445:    Logically Collective on svd

447:    Input Parameters:
448: +  svd     - singular value solver context obtained from SVDCreate()
449: .  func    - pointer to the stopping test function
450: .  ctx     - context for private data for the stopping routine (may be null)
451: -  destroy - a routine for destroying the context (may be null)

453:    Calling Sequence of func:
454: $   func(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 (*func)(SVD,PetscInt,PetscInt,PetscInt,PetscInt,SVDConvergedReason*,void*),void* ctx,PetscErrorCode (*destroy)(void*))
475: {

480:   if (svd->stoppingdestroy) {
481:     (*svd->stoppingdestroy)(svd->stoppingctx);
482:   }
483:   svd->stoppinguser    = func;
484:   svd->stoppingdestroy = destroy;
485:   svd->stoppingctx     = ctx;
486:   if (func == SVDStoppingBasic) svd->stop = SVD_STOP_BASIC;
487:   else {
488:     svd->stop     = SVD_STOP_USER;
489:     svd->stopping = svd->stoppinguser;
490:   }
491:   return(0);
492: }

494: /*@
495:    SVDSetStoppingTest - Specifies how to decide the termination of the outer
496:    loop of the singular value solver.

498:    Logically Collective on svd

500:    Input Parameters:
501: +  svd  - singular value solver context obtained from SVDCreate()
502: -  stop - the type of stopping test

504:    Options Database Keys:
505: +  -svd_stop_basic - Sets the default stopping test
506: -  -svd_stop_user  - Selects the user-defined stopping test

508:    Note:
509:    The parameter 'stop' can have one of these values
510: +     SVD_STOP_BASIC - default stopping test
511: -     SVD_STOP_USER  - function set by SVDSetStoppingTestFunction()

513:    Level: advanced

515: .seealso: SVDGetStoppingTest(), SVDSetStoppingTestFunction(), SVDSetConvergenceTest(), SVDStop
516: @*/
517: PetscErrorCode SVDSetStoppingTest(SVD svd,SVDStop stop)
518: {
522:   switch (stop) {
523:     case SVD_STOP_BASIC: svd->stopping = SVDStoppingBasic; break;
524:     case SVD_STOP_USER:
525:       if (!svd->stoppinguser) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ORDER,"Must call SVDSetStoppingTestFunction() first");
526:       svd->stopping = svd->stoppinguser;
527:       break;
528:     default:
529:       SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'stop' value");
530:   }
531:   svd->stop = stop;
532:   return(0);
533: }

535: /*@
536:    SVDGetStoppingTest - Gets the method used to decide the termination of the outer
537:    loop of the singular value solver.

539:    Not Collective

541:    Input Parameters:
542: .  svd   - singular value solver context obtained from SVDCreate()

544:    Output Parameters:
545: .  stop  - the type of stopping test

547:    Level: advanced

549: .seealso: SVDSetStoppingTest(), SVDStop
550: @*/
551: PetscErrorCode SVDGetStoppingTest(SVD svd,SVDStop *stop)
552: {
556:   *stop = svd->stop;
557:   return(0);
558: }

560: /*@C
561:    SVDMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type
562:    indicated by the user.

564:    Collective on svd

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

573:    Level: developer

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

591:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)svd),((PetscObject)svd)->options,((PetscObject)svd)->prefix,opt,&viewer,&format,&flg);
592:   if (!flg) return(0);

594:   PetscViewerGetType(viewer,&vtype);
595:   SlepcMonitorMakeKey_Internal(name,vtype,format,key);
596:   PetscFunctionListFind(SVDMonitorList,key,&mfunc);
597:   PetscFunctionListFind(SVDMonitorCreateList,key,&cfunc);
598:   PetscFunctionListFind(SVDMonitorDestroyList,key,&dfunc);
599:   if (!cfunc) cfunc = PetscViewerAndFormatCreate_Internal;
600:   if (!dfunc) dfunc = PetscViewerAndFormatDestroy;

602:   (*cfunc)(viewer,format,ctx,&vf);
603:   PetscObjectDereference((PetscObject)viewer);
604:   SVDMonitorSet(svd,mfunc,vf,(PetscErrorCode(*)(void **))dfunc);
605:   if (trackall) {
606:     SVDSetTrackAll(svd,PETSC_TRUE);
607:   }
608:   return(0);
609: }

611: /*@
612:    SVDSetFromOptions - Sets SVD options from the options database.
613:    This routine must be called before SVDSetUp() if the user is to be
614:    allowed to set the solver type.

616:    Collective on svd

618:    Input Parameters:
619: .  svd - the singular value solver context

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

624:    Level: beginner

626: .seealso:
627: @*/
628: PetscErrorCode SVDSetFromOptions(SVD svd)
629: {
631:   char           type[256];
632:   PetscBool      set,flg,val,flg1,flg2,flg3;
633:   PetscInt       i,j,k;
634:   PetscReal      r;

638:   SVDRegisterAll();
639:   PetscObjectOptionsBegin((PetscObject)svd);
640:     PetscOptionsFList("-svd_type","SVD solver method","SVDSetType",SVDList,(char*)(((PetscObject)svd)->type_name?((PetscObject)svd)->type_name:SVDCROSS),type,sizeof(type),&flg);
641:     if (flg) {
642:       SVDSetType(svd,type);
643:     } else if (!((PetscObject)svd)->type_name) {
644:       SVDSetType(svd,SVDCROSS);
645:     }

647:     PetscOptionsBoolGroupBegin("-svd_standard","Singular value decomposition (SVD)","SVDSetProblemType",&flg);
648:     if (flg) { SVDSetProblemType(svd,SVD_STANDARD); }
649:     PetscOptionsBoolGroupEnd("-svd_generalized","Generalized singular value decomposition (GSVD)","SVDSetProblemType",&flg);
650:     if (flg) { SVDSetProblemType(svd,SVD_GENERALIZED); }

652:     PetscOptionsBool("-svd_implicittranspose","Handle matrix transpose implicitly","SVDSetImplicitTranspose",svd->impltrans,&val,&flg);
653:     if (flg) { SVDSetImplicitTranspose(svd,val); }

655:     i = svd->max_it;
656:     PetscOptionsInt("-svd_max_it","Maximum number of iterations","SVDSetTolerances",svd->max_it,&i,&flg1);
657:     r = svd->tol;
658:     PetscOptionsReal("-svd_tol","Tolerance","SVDSetTolerances",svd->tol==PETSC_DEFAULT?SLEPC_DEFAULT_TOL:svd->tol,&r,&flg2);
659:     if (flg1 || flg2) { SVDSetTolerances(svd,r,i); }

661:     PetscOptionsBoolGroupBegin("-svd_conv_abs","Absolute error convergence test","SVDSetConvergenceTest",&flg);
662:     if (flg) { SVDSetConvergenceTest(svd,SVD_CONV_ABS); }
663:     PetscOptionsBoolGroup("-svd_conv_rel","Relative error convergence test","SVDSetConvergenceTest",&flg);
664:     if (flg) { SVDSetConvergenceTest(svd,SVD_CONV_REL); }
665:     PetscOptionsBoolGroup("-svd_conv_maxit","Maximum iterations convergence test","SVDSetConvergenceTest",&flg);
666:     if (flg) { SVDSetConvergenceTest(svd,SVD_CONV_MAXIT); }
667:     PetscOptionsBoolGroupEnd("-svd_conv_user","User-defined convergence test","SVDSetConvergenceTest",&flg);
668:     if (flg) { SVDSetConvergenceTest(svd,SVD_CONV_USER); }

670:     PetscOptionsBoolGroupBegin("-svd_stop_basic","Stop iteration if all singular values converged or max_it reached","SVDSetStoppingTest",&flg);
671:     if (flg) { SVDSetStoppingTest(svd,SVD_STOP_BASIC); }
672:     PetscOptionsBoolGroupEnd("-svd_stop_user","User-defined stopping test","SVDSetStoppingTest",&flg);
673:     if (flg) { SVDSetStoppingTest(svd,SVD_STOP_USER); }

675:     i = svd->nsv;
676:     PetscOptionsInt("-svd_nsv","Number of singular values to compute","SVDSetDimensions",svd->nsv,&i,&flg1);
677:     j = svd->ncv;
678:     PetscOptionsInt("-svd_ncv","Number of basis vectors","SVDSetDimensions",svd->ncv,&j,&flg2);
679:     k = svd->mpd;
680:     PetscOptionsInt("-svd_mpd","Maximum dimension of projected problem","SVDSetDimensions",svd->mpd,&k,&flg3);
681:     if (flg1 || flg2 || flg3) { SVDSetDimensions(svd,i,j,k); }

683:     PetscOptionsBoolGroupBegin("-svd_largest","Compute largest singular values","SVDSetWhichSingularTriplets",&flg);
684:     if (flg) { SVDSetWhichSingularTriplets(svd,SVD_LARGEST); }
685:     PetscOptionsBoolGroupEnd("-svd_smallest","Compute smallest singular values","SVDSetWhichSingularTriplets",&flg);
686:     if (flg) { SVDSetWhichSingularTriplets(svd,SVD_SMALLEST); }

688:     /* -----------------------------------------------------------------------*/
689:     /*
690:       Cancels all monitors hardwired into code before call to SVDSetFromOptions()
691:     */
692:     PetscOptionsBool("-svd_monitor_cancel","Remove any hardwired monitor routines","SVDMonitorCancel",PETSC_FALSE,&flg,&set);
693:     if (set && flg) { SVDMonitorCancel(svd); }
694:     SVDMonitorSetFromOptions(svd,"-svd_monitor","first_approximation",NULL,PETSC_FALSE);
695:     SVDMonitorSetFromOptions(svd,"-svd_monitor_all","all_approximations",NULL,PETSC_TRUE);
696:     SVDMonitorSetFromOptions(svd,"-svd_monitor_conv","convergence_history",NULL,PETSC_FALSE);

698:     /* -----------------------------------------------------------------------*/
699:     PetscOptionsName("-svd_view","Print detailed information on solver used","SVDView",NULL);
700:     PetscOptionsName("-svd_view_vectors","View computed singular vectors","SVDVectorsView",NULL);
701:     PetscOptionsName("-svd_view_values","View computed singular values","SVDValuesView",NULL);
702:     PetscOptionsName("-svd_converged_reason","Print reason for convergence, and number of iterations","SVDConvergedReasonView",NULL);
703:     PetscOptionsName("-svd_error_absolute","Print absolute errors of each singular triplet","SVDErrorView",NULL);
704:     PetscOptionsName("-svd_error_relative","Print relative errors of each singular triplet","SVDErrorView",NULL);

706:     if (svd->ops->setfromoptions) {
707:       (*svd->ops->setfromoptions)(PetscOptionsObject,svd);
708:     }
709:     PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)svd);
710:   PetscOptionsEnd();

712:   if (!svd->V) { SVDGetBV(svd,&svd->V,NULL); }
713:   BVSetFromOptions(svd->V);
714:   if (!svd->U) { SVDGetBV(svd,NULL,&svd->U); }
715:   BVSetFromOptions(svd->U);
716:   if (!svd->ds) { SVDGetDS(svd,&svd->ds); }
717:   DSSetFromOptions(svd->ds);
718:   return(0);
719: }

721: /*@
722:    SVDSetProblemType - Specifies the type of the singular value problem.

724:    Logically Collective on svd

726:    Input Parameters:
727: +  svd  - the singular value solver context
728: -  type - a known type of singular value problem

730:    Options Database Keys:
731: +  -svd_standard    - standard singular value decomposition (SVD)
732: -  -svd_generalized - generalized singular value problem (GSVD)

734:    Notes:
735:    The GSVD requires that two matrices have been passed via SVDSetOperators().

737:    Level: intermediate

739: .seealso: SVDSetOperators(), SVDSetType(), SVDGetProblemType(), SVDProblemType
740: @*/
741: PetscErrorCode SVDSetProblemType(SVD svd,SVDProblemType type)
742: {
746:   if (type == svd->problem_type) return(0);
747:   switch (type) {
748:     case SVD_STANDARD:
749:       svd->isgeneralized = PETSC_FALSE;
750:       break;
751:     case SVD_GENERALIZED:
752:       svd->isgeneralized = PETSC_TRUE;
753:       break;
754:     default:
755:       SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_WRONG,"Unknown singular value problem type");
756:   }
757:   svd->problem_type = type;
758:   svd->state = SVD_STATE_INITIAL;
759:   return(0);
760: }

762: /*@
763:    SVDGetProblemType - Gets the problem type from the SVD object.

765:    Not Collective

767:    Input Parameter:
768: .  svd - the singular value solver context

770:    Output Parameter:
771: .  type - the problem type

773:    Level: intermediate

775: .seealso: SVDSetProblemType(), SVDProblemType
776: @*/
777: PetscErrorCode SVDGetProblemType(SVD svd,SVDProblemType *type)
778: {
782:   *type = svd->problem_type;
783:   return(0);
784: }

786: /*@
787:    SVDIsGeneralized - Ask if the SVD object corresponds to a generalized
788:    singular value problem.

790:    Not collective

792:    Input Parameter:
793: .  svd - the singular value solver context

795:    Output Parameter:
796: .  is - the answer

798:    Level: intermediate

800: .seealso: SVDIsHermitian(), SVDIsPositive()
801: @*/
802: PetscErrorCode SVDIsGeneralized(SVD svd,PetscBool* is)
803: {
807:   *is = svd->isgeneralized;
808:   return(0);
809: }

811: /*@
812:    SVDSetTrackAll - Specifies if the solver must compute the residual norm of all
813:    approximate singular value or not.

815:    Logically Collective on svd

817:    Input Parameters:
818: +  svd      - the singular value solver context
819: -  trackall - whether to compute all residuals or not

821:    Notes:
822:    If the user sets trackall=PETSC_TRUE then the solver computes (or estimates)
823:    the residual norm for each singular value approximation. Computing the residual is
824:    usually an expensive operation and solvers commonly compute only the residual
825:    associated to the first unconverged singular value.

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

829:    Level: developer

831: .seealso: SVDGetTrackAll()
832: @*/
833: PetscErrorCode SVDSetTrackAll(SVD svd,PetscBool trackall)
834: {
838:   svd->trackall = trackall;
839:   return(0);
840: }

842: /*@
843:    SVDGetTrackAll - Returns the flag indicating whether all residual norms must
844:    be computed or not.

846:    Not Collective

848:    Input Parameter:
849: .  svd - the singular value solver context

851:    Output Parameter:
852: .  trackall - the returned flag

854:    Level: developer

856: .seealso: SVDSetTrackAll()
857: @*/
858: PetscErrorCode SVDGetTrackAll(SVD svd,PetscBool *trackall)
859: {
863:   *trackall = svd->trackall;
864:   return(0);
865: }


868: /*@C
869:    SVDSetOptionsPrefix - Sets the prefix used for searching for all
870:    SVD options in the database.

872:    Logically Collective on svd

874:    Input Parameters:
875: +  svd - the singular value solver context
876: -  prefix - the prefix string to prepend to all SVD option requests

878:    Notes:
879:    A hyphen (-) must NOT be given at the beginning of the prefix name.
880:    The first character of all runtime options is AUTOMATICALLY the
881:    hyphen.

883:    For example, to distinguish between the runtime options for two
884:    different SVD contexts, one could call
885: .vb
886:       SVDSetOptionsPrefix(svd1,"svd1_")
887:       SVDSetOptionsPrefix(svd2,"svd2_")
888: .ve

890:    Level: advanced

892: .seealso: SVDAppendOptionsPrefix(), SVDGetOptionsPrefix()
893: @*/
894: PetscErrorCode SVDSetOptionsPrefix(SVD svd,const char *prefix)
895: {

900:   if (!svd->V) { SVDGetBV(svd,&svd->V,&svd->U); }
901:   BVSetOptionsPrefix(svd->V,prefix);
902:   BVSetOptionsPrefix(svd->U,prefix);
903:   if (!svd->ds) { SVDGetDS(svd,&svd->ds); }
904:   DSSetOptionsPrefix(svd->ds,prefix);
905:   PetscObjectSetOptionsPrefix((PetscObject)svd,prefix);
906:   return(0);
907: }

909: /*@C
910:    SVDAppendOptionsPrefix - Appends to the prefix used for searching for all
911:    SVD options in the database.

913:    Logically Collective on svd

915:    Input Parameters:
916: +  svd - the singular value solver context
917: -  prefix - the prefix string to prepend to all SVD option requests

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

923:    Level: advanced

925: .seealso: SVDSetOptionsPrefix(), SVDGetOptionsPrefix()
926: @*/
927: PetscErrorCode SVDAppendOptionsPrefix(SVD svd,const char *prefix)
928: {

933:   if (!svd->V) { SVDGetBV(svd,&svd->V,&svd->U); }
934:   BVAppendOptionsPrefix(svd->V,prefix);
935:   BVAppendOptionsPrefix(svd->U,prefix);
936:   if (!svd->ds) { SVDGetDS(svd,&svd->ds); }
937:   DSAppendOptionsPrefix(svd->ds,prefix);
938:   PetscObjectAppendOptionsPrefix((PetscObject)svd,prefix);
939:   return(0);
940: }

942: /*@C
943:    SVDGetOptionsPrefix - Gets the prefix used for searching for all
944:    SVD options in the database.

946:    Not Collective

948:    Input Parameters:
949: .  svd - the singular value solver context

951:    Output Parameters:
952: .  prefix - pointer to the prefix string used is returned

954:    Note:
955:    On the Fortran side, the user should pass in a string 'prefix' of
956:    sufficient length to hold the prefix.

958:    Level: advanced

960: .seealso: SVDSetOptionsPrefix(), SVDAppendOptionsPrefix()
961: @*/
962: PetscErrorCode SVDGetOptionsPrefix(SVD svd,const char *prefix[])
963: {

969:   PetscObjectGetOptionsPrefix((PetscObject)svd,prefix);
970:   return(0);
971: }