TaoSetFromOptions#

Sets various Tao parameters from the options database

Synopsis#

#include "petsctao.h" 
PetscErrorCode TaoSetFromOptions(Tao tao)

Collective

Input Parameter#

  • tao - the Tao solver context

Options Database Keys#

  • -tao_type - The algorithm that Tao uses (lmvm, nls, etc.)

  • -tao_gatol - absolute error tolerance for ||gradient||

  • -tao_grtol - relative error tolerance for ||gradient||

  • -tao_gttol - reduction of ||gradient|| relative to initial gradient

  • -tao_max_it - sets maximum number of iterations

  • -tao_max_funcs - sets maximum number of function evaluations

  • -tao_fmin - stop if function value reaches fmin

  • -tao_steptol - stop if trust region radius less than

  • -tao_trust0 - initial trust region radius

  • -tao_view_solution - view the solution at the end of the optimization process

  • -tao_monitor - prints function value and residual norm at each iteration

  • -tao_monitor_short - same as -tao_monitor, but truncates very small values

  • -tao_monitor_constraint_norm - prints objective value, gradient, and constraint norm at each iteration

  • -tao_monitor_globalization - prints information about the globalization at each iteration

  • -tao_monitor_solution - prints solution vector at each iteration

  • -tao_monitor_ls_residual - prints least-squares residual vector at each iteration

  • -tao_monitor_step - prints step vector at each iteration

  • -tao_monitor_gradient - prints gradient vector at each iteration

  • -tao_monitor_solution_draw - graphically view solution vector at each iteration

  • -tao_monitor_step_draw - graphically view step vector at each iteration

  • -tao_monitor_gradient_draw - graphically view gradient at each iteration

  • -tao_monitor_cancel - cancels all monitors (except those set with command line)

  • -tao_fd_gradient - use gradient computed with finite differences

  • -tao_fd_hessian - use hessian computed with finite differences

  • -tao_mf_hessian - use matrix-free Hessian computed with finite differences

  • -tao_view - prints information about the Tao after solving

  • -tao_converged_reason - prints the reason Tao stopped iterating

Note#

To see all options, run your program with the -help option or consult the user’s manual. Should be called after TaoCreate() but before TaoSolve()

See Also#

TAO: Optimization Solvers, Tao, TaoCreate(), TaoSolve()

Level#

beginner

Location#

src/tao/interface/taosolver.c

Examples#

src/tao/unconstrained/tutorials/minsurf2.c
src/tao/complementarity/tutorials/minsurf1.c
src/tao/unconstrained/tutorials/rosenbrock3.c
src/tao/leastsquares/tutorials/tomography.c
src/tao/leastsquares/tutorials/chwirut1.c
src/tao/leastsquares/tutorials/chwirut2.c
src/tao/leastsquares/tutorials/cs1.c
src/tao/bound/tutorials/plate2.c
src/tao/complementarity/tutorials/blackscholes.c
src/tao/bound/tutorials/jbearing2.c

Implementations#

TaoSetFromOptions_BLMVM() in src/tao/bound/impls/blmvm/blmvm.c
TaoSetFromOptions_BNCG() in src/tao/bound/impls/bncg/bncg.c
TaoSetFromOptions_BNK() in src/tao/bound/impls/bnk/bnk.c
TaoSetFromOptions_BNTL() in src/tao/bound/impls/bnk/bntl.c
TaoSetFromOptions_BNTR() in src/tao/bound/impls/bnk/bntr.c
TaoSetFromOptions_BQNK() in src/tao/bound/impls/bqnk/bqnk.c
TaoSetFromOptions_BQNLS() in src/tao/bound/impls/bqnls/bqnls.c
TaoSetFromOptions_TRON() in src/tao/bound/impls/tron/tron.c
TaoSetFromOptions_SSLS() in src/tao/complementarity/impls/ssls/ssls.c
TaoSetFromOptions_ADMM() in src/tao/constrained/impls/admm/admm.c
TaoSetFromOptions_ALMM() in src/tao/constrained/impls/almm/almm.c
TaoSetFromOptions_IPM() in src/tao/constrained/impls/ipm/ipm.c
TaoSetFromOptions_PDIPM() in src/tao/constrained/impls/ipm/pdipm.c
TaoSetFromOptions_BRGN() in src/tao/leastsquares/impls/brgn/brgn.c
TaoSetFromOptions_POUNDERS() in src/tao/leastsquares/impls/pounders/pounders.c
TaoSetFromOptions_LCL() in src/tao/pde_constrained/impls/lcl/lcl.c
TaoSetFromOptions_BQPIP() in src/tao/quadratic/impls/bqpip/bqpip.c
TaoSetFromOptions_GPCG() in src/tao/quadratic/impls/gpcg/gpcg.c
TaoSetFromOptions_BMRM() in src/tao/unconstrained/impls/bmrm/bmrm.c
TaoSetFromOptions_CG() in src/tao/unconstrained/impls/cg/taocg.c
TaoSetFromOptions_LMVM() in src/tao/unconstrained/impls/lmvm/lmvm.c
TaoSetFromOptions_NM() in src/tao/unconstrained/impls/neldermead/neldermead.c
TaoSetFromOptions_NLS() in src/tao/unconstrained/impls/nls/nls.c
TaoSetFromOptions_NTL() in src/tao/unconstrained/impls/ntl/ntl.c
TaoSetFromOptions_NTR() in src/tao/unconstrained/impls/ntr/ntr.c
TaoSetFromOptions_OWLQN() in src/tao/unconstrained/impls/owlqn/owlqn.c


Index of all Tao routines
Table of Contents for all manual pages
Index of all manual pages