xref: /petsc/src/tao/interface/taosolver.c (revision c3b5f7ba6bc5ce25a01a67bb37ba5d34b02bbbd7)
1 #define TAO_DLL
2 
3 #include <petsc/private/taoimpl.h> /*I "petsctao.h" I*/
4 
5 PetscBool TaoRegisterAllCalled = PETSC_FALSE;
6 PetscFunctionList TaoList = NULL;
7 
8 PetscClassId TAO_CLASSID;
9 
10 PetscLogEvent TAO_Solve;
11 PetscLogEvent TAO_ObjectiveEval;
12 PetscLogEvent TAO_GradientEval;
13 PetscLogEvent TAO_ObjGradEval;
14 PetscLogEvent TAO_HessianEval;
15 PetscLogEvent TAO_JacobianEval;
16 PetscLogEvent TAO_ConstraintsEval;
17 
18 const char *TaoSubSetTypes[] = {"subvec","mask","matrixfree","TaoSubSetType","TAO_SUBSET_",NULL};
19 
20 struct _n_TaoMonitorDrawCtx {
21   PetscViewer viewer;
22   PetscInt    howoften;  /* when > 0 uses iteration % howoften, when negative only final solution plotted */
23 };
24 
25 /*@
26   TaoCreate - Creates a TAO solver
27 
28   Collective
29 
30   Input Parameter:
31 . comm - MPI communicator
32 
33   Output Parameter:
34 . newtao - the new Tao context
35 
36   Available methods include:
37 +    nls - Newton's method with line search for unconstrained minimization
38 .    ntr - Newton's method with trust region for unconstrained minimization
39 .    ntl - Newton's method with trust region, line search for unconstrained minimization
40 .    lmvm - Limited memory variable metric method for unconstrained minimization
41 .    cg - Nonlinear conjugate gradient method for unconstrained minimization
42 .    nm - Nelder-Mead algorithm for derivate-free unconstrained minimization
43 .    tron - Newton Trust Region method for bound constrained minimization
44 .    gpcg - Newton Trust Region method for quadratic bound constrained minimization
45 .    blmvm - Limited memory variable metric method for bound constrained minimization
46 .    lcl - Linearly constrained Lagrangian method for pde-constrained minimization
47 -    pounders - Model-based algorithm for nonlinear least squares
48 
49    Options Database Keys:
50 .   -tao_type - select which method TAO should use
51 
52    Level: beginner
53 
54 .seealso: TaoSolve(), TaoDestroy()
55 @*/
56 PetscErrorCode TaoCreate(MPI_Comm comm, Tao *newtao)
57 {
58   Tao            tao;
59 
60   PetscFunctionBegin;
61   PetscValidPointer(newtao,2);
62   PetscCall(TaoInitializePackage());
63   PetscCall(TaoLineSearchInitializePackage());
64   PetscCall(PetscHeaderCreate(tao,TAO_CLASSID,"Tao","Optimization solver","Tao",comm,TaoDestroy,TaoView));
65 
66   /* Set non-NULL defaults */
67   tao->ops->convergencetest = TaoDefaultConvergenceTest;
68 
69   tao->max_it    = 10000;
70   tao->max_funcs = -1;
71 #if defined(PETSC_USE_REAL_SINGLE)
72   tao->gatol     = 1e-5;
73   tao->grtol     = 1e-5;
74   tao->crtol     = 1e-5;
75   tao->catol     = 1e-5;
76 #else
77   tao->gatol     = 1e-8;
78   tao->grtol     = 1e-8;
79   tao->crtol     = 1e-8;
80   tao->catol     = 1e-8;
81 #endif
82   tao->gttol     = 0.0;
83   tao->steptol   = 0.0;
84   tao->trust0    = PETSC_INFINITY;
85   tao->fmin      = PETSC_NINFINITY;
86 
87   tao->hist_reset  = PETSC_TRUE;
88 
89   PetscCall(TaoResetStatistics(tao));
90   *newtao = tao;
91   PetscFunctionReturn(0);
92 }
93 
94 /*@
95   TaoSolve - Solves an optimization problem min F(x) s.t. l <= x <= u
96 
97   Collective on Tao
98 
99   Input Parameters:
100 . tao - the Tao context
101 
102   Notes:
103   The user must set up the Tao with calls to TaoSetSolution(),
104   TaoSetObjective(),
105   TaoSetGradient(), and (if using 2nd order method) TaoSetHessian().
106 
107   You should call TaoGetConvergedReason() or run with -tao_converged_reason to determine if the optimization algorithm actually succeeded or
108   why it failed.
109 
110   Level: beginner
111 
112 .seealso: TaoCreate(), TaoSetObjective(), TaoSetGradient(), TaoSetHessian(), TaoGetConvergedReason()
113  @*/
114 PetscErrorCode TaoSolve(Tao tao)
115 {
116   PetscErrorCode   ierr;
117   static PetscBool set = PETSC_FALSE;
118 
119   PetscFunctionBegin;
120   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
121   ierr = PetscCitationsRegister("@TechReport{tao-user-ref,\n"
122                                 "title   = {Toolkit for Advanced Optimization (TAO) Users Manual},\n"
123                                 "author  = {Todd Munson and Jason Sarich and Stefan Wild and Steve Benson and Lois Curfman McInnes},\n"
124                                 "Institution = {Argonne National Laboratory},\n"
125                                 "Year   = 2014,\n"
126                                 "Number = {ANL/MCS-TM-322 - Revision 3.5},\n"
127                                 "url    = {https://www.mcs.anl.gov/research/projects/tao/}\n}\n",&set);PetscCall(ierr);
128   tao->header_printed = PETSC_FALSE;
129   PetscCall(TaoSetUp(tao));
130   PetscCall(TaoResetStatistics(tao));
131   if (tao->linesearch) {
132     PetscCall(TaoLineSearchReset(tao->linesearch));
133   }
134 
135   PetscCall(PetscLogEventBegin(TAO_Solve,tao,0,0,0));
136   if (tao->ops->solve) PetscCall((*tao->ops->solve)(tao));
137   PetscCall(PetscLogEventEnd(TAO_Solve,tao,0,0,0));
138 
139   PetscCall(VecViewFromOptions(tao->solution,(PetscObject)tao,"-tao_view_solution"));
140 
141   tao->ntotalits += tao->niter;
142   PetscCall(TaoViewFromOptions(tao,NULL,"-tao_view"));
143 
144   if (tao->printreason) {
145     if (tao->reason > 0) {
146       PetscCall(PetscPrintf(((PetscObject)tao)->comm,"TAO solve converged due to %s iterations %D\n",TaoConvergedReasons[tao->reason],tao->niter));
147     } else {
148       PetscCall(PetscPrintf(((PetscObject)tao)->comm,"TAO solve did not converge due to %s iteration %D\n",TaoConvergedReasons[tao->reason],tao->niter));
149     }
150   }
151   PetscFunctionReturn(0);
152 }
153 
154 /*@
155   TaoSetUp - Sets up the internal data structures for the later use
156   of a Tao solver
157 
158   Collective on tao
159 
160   Input Parameters:
161 . tao - the TAO context
162 
163   Notes:
164   The user will not need to explicitly call TaoSetUp(), as it will
165   automatically be called in TaoSolve().  However, if the user
166   desires to call it explicitly, it should come after TaoCreate()
167   and any TaoSetSomething() routines, but before TaoSolve().
168 
169   Level: advanced
170 
171 .seealso: TaoCreate(), TaoSolve()
172 @*/
173 PetscErrorCode TaoSetUp(Tao tao)
174 {
175   PetscFunctionBegin;
176   PetscValidHeaderSpecific(tao, TAO_CLASSID,1);
177   if (tao->setupcalled) PetscFunctionReturn(0);
178   PetscCheck(tao->solution,PetscObjectComm((PetscObject)tao),PETSC_ERR_ARG_WRONGSTATE,"Must call TaoSetSolution");
179   if (tao->ops->setup) {
180     PetscCall((*tao->ops->setup)(tao));
181   }
182   tao->setupcalled = PETSC_TRUE;
183   PetscFunctionReturn(0);
184 }
185 
186 /*@C
187   TaoDestroy - Destroys the TAO context that was created with
188   TaoCreate()
189 
190   Collective on Tao
191 
192   Input Parameter:
193 . tao - the Tao context
194 
195   Level: beginner
196 
197 .seealso: TaoCreate(), TaoSolve()
198 @*/
199 PetscErrorCode TaoDestroy(Tao *tao)
200 {
201   PetscFunctionBegin;
202   if (!*tao) PetscFunctionReturn(0);
203   PetscValidHeaderSpecific(*tao,TAO_CLASSID,1);
204   if (--((PetscObject)*tao)->refct > 0) {*tao = NULL;PetscFunctionReturn(0);}
205 
206   if ((*tao)->ops->destroy) {
207     PetscCall((*((*tao))->ops->destroy)(*tao));
208   }
209   PetscCall(KSPDestroy(&(*tao)->ksp));
210   PetscCall(TaoLineSearchDestroy(&(*tao)->linesearch));
211 
212   if ((*tao)->ops->convergencedestroy) {
213     PetscCall((*(*tao)->ops->convergencedestroy)((*tao)->cnvP));
214     if ((*tao)->jacobian_state_inv) {
215       PetscCall(MatDestroy(&(*tao)->jacobian_state_inv));
216     }
217   }
218   PetscCall(VecDestroy(&(*tao)->solution));
219   PetscCall(VecDestroy(&(*tao)->gradient));
220   PetscCall(VecDestroy(&(*tao)->ls_res));
221 
222   if ((*tao)->gradient_norm) {
223     PetscCall(PetscObjectDereference((PetscObject)(*tao)->gradient_norm));
224     PetscCall(VecDestroy(&(*tao)->gradient_norm_tmp));
225   }
226 
227   PetscCall(VecDestroy(&(*tao)->XL));
228   PetscCall(VecDestroy(&(*tao)->XU));
229   PetscCall(VecDestroy(&(*tao)->IL));
230   PetscCall(VecDestroy(&(*tao)->IU));
231   PetscCall(VecDestroy(&(*tao)->DE));
232   PetscCall(VecDestroy(&(*tao)->DI));
233   PetscCall(VecDestroy(&(*tao)->constraints));
234   PetscCall(VecDestroy(&(*tao)->constraints_equality));
235   PetscCall(VecDestroy(&(*tao)->constraints_inequality));
236   PetscCall(VecDestroy(&(*tao)->stepdirection));
237   PetscCall(MatDestroy(&(*tao)->hessian_pre));
238   PetscCall(MatDestroy(&(*tao)->hessian));
239   PetscCall(MatDestroy(&(*tao)->ls_jac));
240   PetscCall(MatDestroy(&(*tao)->ls_jac_pre));
241   PetscCall(MatDestroy(&(*tao)->jacobian_pre));
242   PetscCall(MatDestroy(&(*tao)->jacobian));
243   PetscCall(MatDestroy(&(*tao)->jacobian_state_pre));
244   PetscCall(MatDestroy(&(*tao)->jacobian_state));
245   PetscCall(MatDestroy(&(*tao)->jacobian_state_inv));
246   PetscCall(MatDestroy(&(*tao)->jacobian_design));
247   PetscCall(MatDestroy(&(*tao)->jacobian_equality));
248   PetscCall(MatDestroy(&(*tao)->jacobian_equality_pre));
249   PetscCall(MatDestroy(&(*tao)->jacobian_inequality));
250   PetscCall(MatDestroy(&(*tao)->jacobian_inequality_pre));
251   PetscCall(ISDestroy(&(*tao)->state_is));
252   PetscCall(ISDestroy(&(*tao)->design_is));
253   PetscCall(VecDestroy(&(*tao)->res_weights_v));
254   PetscCall(TaoCancelMonitors(*tao));
255   if ((*tao)->hist_malloc) {
256     PetscCall(PetscFree4((*tao)->hist_obj,(*tao)->hist_resid,(*tao)->hist_cnorm,(*tao)->hist_lits));
257   }
258   if ((*tao)->res_weights_n) {
259     PetscCall(PetscFree((*tao)->res_weights_rows));
260     PetscCall(PetscFree((*tao)->res_weights_cols));
261     PetscCall(PetscFree((*tao)->res_weights_w));
262   }
263   PetscCall(PetscHeaderDestroy(tao));
264   PetscFunctionReturn(0);
265 }
266 
267 /*@
268   TaoSetFromOptions - Sets various Tao parameters from user
269   options.
270 
271   Collective on Tao
272 
273   Input Parameter:
274 . tao - the Tao solver context
275 
276   options Database Keys:
277 + -tao_type <type> - The algorithm that TAO uses (lmvm, nls, etc.)
278 . -tao_gatol <gatol> - absolute error tolerance for ||gradient||
279 . -tao_grtol <grtol> - relative error tolerance for ||gradient||
280 . -tao_gttol <gttol> - reduction of ||gradient|| relative to initial gradient
281 . -tao_max_it <max> - sets maximum number of iterations
282 . -tao_max_funcs <max> - sets maximum number of function evaluations
283 . -tao_fmin <fmin> - stop if function value reaches fmin
284 . -tao_steptol <tol> - stop if trust region radius less than <tol>
285 . -tao_trust0 <t> - initial trust region radius
286 . -tao_monitor - prints function value and residual at each iteration
287 . -tao_smonitor - same as tao_monitor, but truncates very small values
288 . -tao_cmonitor - prints function value, residual, and constraint norm at each iteration
289 . -tao_view_solution - prints solution vector at each iteration
290 . -tao_view_ls_residual - prints least-squares residual vector at each iteration
291 . -tao_view_stepdirection - prints step direction vector at each iteration
292 . -tao_view_gradient - prints gradient vector at each iteration
293 . -tao_draw_solution - graphically view solution vector at each iteration
294 . -tao_draw_step - graphically view step vector at each iteration
295 . -tao_draw_gradient - graphically view gradient at each iteration
296 . -tao_fd_gradient - use gradient computed with finite differences
297 . -tao_fd_hessian - use hessian computed with finite differences
298 . -tao_mf_hessian - use matrix-free hessian computed with finite differences
299 . -tao_cancelmonitors - cancels all monitors (except those set with command line)
300 . -tao_view - prints information about the Tao after solving
301 - -tao_converged_reason - prints the reason TAO stopped iterating
302 
303   Notes:
304   To see all options, run your program with the -help option or consult the
305   user's manual. Should be called after TaoCreate() but before TaoSolve()
306 
307   Level: beginner
308 @*/
309 PetscErrorCode TaoSetFromOptions(Tao tao)
310 {
311   PetscErrorCode ierr;
312   TaoType        default_type = TAOLMVM;
313   char           type[256], monfilename[PETSC_MAX_PATH_LEN];
314   PetscViewer    monviewer;
315   PetscBool      flg;
316   MPI_Comm       comm;
317 
318   PetscFunctionBegin;
319   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
320   PetscCall(PetscObjectGetComm((PetscObject)tao,&comm));
321 
322   /* So no warnings are given about unused options */
323   PetscCall(PetscOptionsHasName(((PetscObject)tao)->options,((PetscObject)tao)->prefix,"-tao_ls_type",&flg));
324 
325   ierr = PetscObjectOptionsBegin((PetscObject)tao);PetscCall(ierr);
326   {
327     if (((PetscObject)tao)->type_name) default_type = ((PetscObject)tao)->type_name;
328     /* Check for type from options */
329     PetscCall(PetscOptionsFList("-tao_type","Tao Solver type","TaoSetType",TaoList,default_type,type,256,&flg));
330     if (flg) {
331       PetscCall(TaoSetType(tao,type));
332     } else if (!((PetscObject)tao)->type_name) {
333       PetscCall(TaoSetType(tao,default_type));
334     }
335 
336     PetscCall(PetscOptionsReal("-tao_catol","Stop if constraints violations within","TaoSetConstraintTolerances",tao->catol,&tao->catol,&flg));
337     if (flg) tao->catol_changed = PETSC_TRUE;
338     PetscCall(PetscOptionsReal("-tao_crtol","Stop if relative constraint violations within","TaoSetConstraintTolerances",tao->crtol,&tao->crtol,&flg));
339     if (flg) tao->crtol_changed = PETSC_TRUE;
340     PetscCall(PetscOptionsReal("-tao_gatol","Stop if norm of gradient less than","TaoSetTolerances",tao->gatol,&tao->gatol,&flg));
341     if (flg) tao->gatol_changed = PETSC_TRUE;
342     PetscCall(PetscOptionsReal("-tao_grtol","Stop if norm of gradient divided by the function value is less than","TaoSetTolerances",tao->grtol,&tao->grtol,&flg));
343     if (flg) tao->grtol_changed = PETSC_TRUE;
344     PetscCall(PetscOptionsReal("-tao_gttol","Stop if the norm of the gradient is less than the norm of the initial gradient times tol","TaoSetTolerances",tao->gttol,&tao->gttol,&flg));
345     if (flg) tao->gttol_changed = PETSC_TRUE;
346     PetscCall(PetscOptionsInt("-tao_max_it","Stop if iteration number exceeds","TaoSetMaximumIterations",tao->max_it,&tao->max_it,&flg));
347     if (flg) tao->max_it_changed = PETSC_TRUE;
348     PetscCall(PetscOptionsInt("-tao_max_funcs","Stop if number of function evaluations exceeds","TaoSetMaximumFunctionEvaluations",tao->max_funcs,&tao->max_funcs,&flg));
349     if (flg) tao->max_funcs_changed = PETSC_TRUE;
350     PetscCall(PetscOptionsReal("-tao_fmin","Stop if function less than","TaoSetFunctionLowerBound",tao->fmin,&tao->fmin,&flg));
351     if (flg) tao->fmin_changed = PETSC_TRUE;
352     PetscCall(PetscOptionsReal("-tao_steptol","Stop if step size or trust region radius less than","",tao->steptol,&tao->steptol,&flg));
353     if (flg) tao->steptol_changed = PETSC_TRUE;
354     PetscCall(PetscOptionsReal("-tao_trust0","Initial trust region radius","TaoSetTrustRegionRadius",tao->trust0,&tao->trust0,&flg));
355     if (flg) tao->trust0_changed = PETSC_TRUE;
356     PetscCall(PetscOptionsString("-tao_view_solution","view solution vector after each evaluation","TaoSetMonitor","stdout",monfilename,sizeof(monfilename),&flg));
357     if (flg) {
358       PetscCall(PetscViewerASCIIOpen(comm,monfilename,&monviewer));
359       PetscCall(TaoSetMonitor(tao,TaoSolutionMonitor,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy));
360     }
361 
362     PetscCall(PetscOptionsBool("-tao_converged_reason","Print reason for TAO converged","TaoSolve",tao->printreason,&tao->printreason,NULL));
363     PetscCall(PetscOptionsString("-tao_view_gradient","view gradient vector after each evaluation","TaoSetMonitor","stdout",monfilename,sizeof(monfilename),&flg));
364     if (flg) {
365       PetscCall(PetscViewerASCIIOpen(comm,monfilename,&monviewer));
366       PetscCall(TaoSetMonitor(tao,TaoGradientMonitor,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy));
367     }
368 
369     PetscCall(PetscOptionsString("-tao_view_stepdirection","view step direction vector after each iteration","TaoSetMonitor","stdout",monfilename,sizeof(monfilename),&flg));
370     if (flg) {
371       PetscCall(PetscViewerASCIIOpen(comm,monfilename,&monviewer));
372       PetscCall(TaoSetMonitor(tao,TaoStepDirectionMonitor,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy));
373     }
374 
375     PetscCall(PetscOptionsString("-tao_view_residual","view least-squares residual vector after each evaluation","TaoSetMonitor","stdout",monfilename,sizeof(monfilename),&flg));
376     if (flg) {
377       PetscCall(PetscViewerASCIIOpen(comm,monfilename,&monviewer));
378       PetscCall(TaoSetMonitor(tao,TaoResidualMonitor,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy));
379     }
380 
381     PetscCall(PetscOptionsString("-tao_monitor","Use the default convergence monitor","TaoSetMonitor","stdout",monfilename,sizeof(monfilename),&flg));
382     if (flg) {
383       PetscCall(PetscViewerASCIIOpen(comm,monfilename,&monviewer));
384       PetscCall(TaoSetMonitor(tao,TaoMonitorDefault,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy));
385     }
386 
387     PetscCall(PetscOptionsString("-tao_gmonitor","Use the convergence monitor with extra globalization info","TaoSetMonitor","stdout",monfilename,sizeof(monfilename),&flg));
388     if (flg) {
389       PetscCall(PetscViewerASCIIOpen(comm,monfilename,&monviewer));
390       PetscCall(TaoSetMonitor(tao,TaoDefaultGMonitor,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy));
391     }
392 
393     PetscCall(PetscOptionsString("-tao_smonitor","Use the short convergence monitor","TaoSetMonitor","stdout",monfilename,sizeof(monfilename),&flg));
394     if (flg) {
395       PetscCall(PetscViewerASCIIOpen(comm,monfilename,&monviewer));
396       PetscCall(TaoSetMonitor(tao,TaoDefaultSMonitor,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy));
397     }
398 
399     PetscCall(PetscOptionsString("-tao_cmonitor","Use the default convergence monitor with constraint norm","TaoSetMonitor","stdout",monfilename,sizeof(monfilename),&flg));
400     if (flg) {
401       PetscCall(PetscViewerASCIIOpen(comm,monfilename,&monviewer));
402       PetscCall(TaoSetMonitor(tao,TaoDefaultCMonitor,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy));
403     }
404 
405     flg = PETSC_FALSE;
406     PetscCall(PetscOptionsBool("-tao_cancelmonitors","cancel all monitors and call any registered destroy routines","TaoCancelMonitors",flg,&flg,NULL));
407     if (flg) PetscCall(TaoCancelMonitors(tao));
408 
409     flg = PETSC_FALSE;
410     PetscCall(PetscOptionsBool("-tao_draw_solution","Plot solution vector at each iteration","TaoSetMonitor",flg,&flg,NULL));
411     if (flg) {
412       TaoMonitorDrawCtx drawctx;
413       PetscInt          howoften = 1;
414       PetscCall(TaoMonitorDrawCtxCreate(PetscObjectComm((PetscObject)tao),NULL,NULL,PETSC_DECIDE,PETSC_DECIDE,300,300,howoften,&drawctx));
415       PetscCall(TaoSetMonitor(tao,TaoDrawSolutionMonitor,drawctx,(PetscErrorCode (*)(void**))TaoMonitorDrawCtxDestroy));
416     }
417 
418     flg = PETSC_FALSE;
419     PetscCall(PetscOptionsBool("-tao_draw_step","plots step direction at each iteration","TaoSetMonitor",flg,&flg,NULL));
420     if (flg) {
421       PetscCall(TaoSetMonitor(tao,TaoDrawStepMonitor,NULL,NULL));
422     }
423 
424     flg = PETSC_FALSE;
425     PetscCall(PetscOptionsBool("-tao_draw_gradient","plots gradient at each iteration","TaoSetMonitor",flg,&flg,NULL));
426     if (flg) {
427       TaoMonitorDrawCtx drawctx;
428       PetscInt          howoften = 1;
429       PetscCall(TaoMonitorDrawCtxCreate(PetscObjectComm((PetscObject)tao),NULL,NULL,PETSC_DECIDE,PETSC_DECIDE,300,300,howoften,&drawctx));
430       PetscCall(TaoSetMonitor(tao,TaoDrawGradientMonitor,drawctx,(PetscErrorCode (*)(void**))TaoMonitorDrawCtxDestroy));
431     }
432     flg = PETSC_FALSE;
433     PetscCall(PetscOptionsBool("-tao_fd_gradient","compute gradient using finite differences","TaoDefaultComputeGradient",flg,&flg,NULL));
434     if (flg) {
435       PetscCall(TaoSetGradient(tao,NULL,TaoDefaultComputeGradient,NULL));
436     }
437     flg = PETSC_FALSE;
438     PetscCall(PetscOptionsBool("-tao_fd_hessian","compute hessian using finite differences","TaoDefaultComputeHessian",flg,&flg,NULL));
439     if (flg) {
440       Mat H;
441 
442       PetscCall(MatCreate(PetscObjectComm((PetscObject)tao),&H));
443       PetscCall(MatSetType(H,MATAIJ));
444       PetscCall(TaoSetHessian(tao,H,H,TaoDefaultComputeHessian,NULL));
445       PetscCall(MatDestroy(&H));
446     }
447     flg = PETSC_FALSE;
448     PetscCall(PetscOptionsBool("-tao_mf_hessian","compute matrix-free hessian using finite differences","TaoDefaultComputeHessianMFFD",flg,&flg,NULL));
449     if (flg) {
450       Mat H;
451 
452       PetscCall(MatCreate(PetscObjectComm((PetscObject)tao),&H));
453       PetscCall(TaoSetHessian(tao,H,H,TaoDefaultComputeHessianMFFD,NULL));
454       PetscCall(MatDestroy(&H));
455     }
456     flg = PETSC_FALSE;
457     PetscCall(PetscOptionsBool("-tao_recycle_history","enable recycling/re-using information from the previous TaoSolve() call for some algorithms","TaoSetRecycleHistory",flg,&flg,NULL));
458     if (flg) {
459       PetscCall(TaoSetRecycleHistory(tao,PETSC_TRUE));
460     }
461     PetscCall(PetscOptionsEnum("-tao_subset_type","subset type","",TaoSubSetTypes,(PetscEnum)tao->subset_type,(PetscEnum*)&tao->subset_type,NULL));
462 
463     if (tao->linesearch) {
464       PetscCall(TaoLineSearchSetFromOptions(tao->linesearch));
465     }
466 
467     if (tao->ops->setfromoptions) {
468       PetscCall((*tao->ops->setfromoptions)(PetscOptionsObject,tao));
469     }
470   }
471   ierr = PetscOptionsEnd();PetscCall(ierr);
472   PetscFunctionReturn(0);
473 }
474 
475 /*@C
476    TaoViewFromOptions - View from Options
477 
478    Collective on Tao
479 
480    Input Parameters:
481 +  A - the  Tao context
482 .  obj - Optional object
483 -  name - command line option
484 
485    Level: intermediate
486 .seealso:  Tao, TaoView, PetscObjectViewFromOptions(), TaoCreate()
487 @*/
488 PetscErrorCode  TaoViewFromOptions(Tao A,PetscObject obj,const char name[])
489 {
490   PetscFunctionBegin;
491   PetscValidHeaderSpecific(A,TAO_CLASSID,1);
492   PetscCall(PetscObjectViewFromOptions((PetscObject)A,obj,name));
493   PetscFunctionReturn(0);
494 }
495 
496 /*@C
497   TaoView - Prints information about the Tao
498 
499   Collective on Tao
500 
501   InputParameters:
502 + tao - the Tao context
503 - viewer - visualization context
504 
505   Options Database Key:
506 . -tao_view - Calls TaoView() at the end of TaoSolve()
507 
508   Notes:
509   The available visualization contexts include
510 +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
511 -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
512          output where only the first processor opens
513          the file.  All other processors send their
514          data to the first processor to print.
515 
516   Level: beginner
517 
518 .seealso: PetscViewerASCIIOpen()
519 @*/
520 PetscErrorCode TaoView(Tao tao, PetscViewer viewer)
521 {
522   PetscBool           isascii,isstring;
523   TaoType             type;
524 
525   PetscFunctionBegin;
526   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
527   if (!viewer) {
528     PetscCall(PetscViewerASCIIGetStdout(((PetscObject)tao)->comm,&viewer));
529   }
530   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
531   PetscCheckSameComm(tao,1,viewer,2);
532 
533   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
534   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring));
535   if (isascii) {
536     PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)tao,viewer));
537 
538     if (tao->ops->view) {
539       PetscCall(PetscViewerASCIIPushTab(viewer));
540       PetscCall((*tao->ops->view)(tao,viewer));
541       PetscCall(PetscViewerASCIIPopTab(viewer));
542     }
543     if (tao->linesearch) {
544       PetscCall(PetscViewerASCIIPushTab(viewer));
545       PetscCall(TaoLineSearchView(tao->linesearch,viewer));
546       PetscCall(PetscViewerASCIIPopTab(viewer));
547     }
548     if (tao->ksp) {
549       PetscCall(PetscViewerASCIIPushTab(viewer));
550       PetscCall(KSPView(tao->ksp,viewer));
551       PetscCall(PetscViewerASCIIPrintf(viewer,"total KSP iterations: %D\n",tao->ksp_tot_its));
552       PetscCall(PetscViewerASCIIPopTab(viewer));
553     }
554 
555     PetscCall(PetscViewerASCIIPushTab(viewer));
556 
557     if (tao->XL || tao->XU) {
558       PetscCall(PetscViewerASCIIPrintf(viewer,"Active Set subset type: %s\n",TaoSubSetTypes[tao->subset_type]));
559     }
560 
561     PetscCall(PetscViewerASCIIPrintf(viewer,"convergence tolerances: gatol=%g,",(double)tao->gatol));
562     PetscCall(PetscViewerASCIIPrintf(viewer," steptol=%g,",(double)tao->steptol));
563     PetscCall(PetscViewerASCIIPrintf(viewer," gttol=%g\n",(double)tao->gttol));
564     PetscCall(PetscViewerASCIIPrintf(viewer,"Residual in Function/Gradient:=%g\n",(double)tao->residual));
565 
566     if (tao->constrained) {
567       PetscCall(PetscViewerASCIIPrintf(viewer,"convergence tolerances:"));
568       PetscCall(PetscViewerASCIIPrintf(viewer," catol=%g,",(double)tao->catol));
569       PetscCall(PetscViewerASCIIPrintf(viewer," crtol=%g\n",(double)tao->crtol));
570       PetscCall(PetscViewerASCIIPrintf(viewer,"Residual in Constraints:=%g\n",(double)tao->cnorm));
571     }
572 
573     if (tao->trust < tao->steptol) {
574       PetscCall(PetscViewerASCIIPrintf(viewer,"convergence tolerances: steptol=%g\n",(double)tao->steptol));
575       PetscCall(PetscViewerASCIIPrintf(viewer,"Final trust region radius:=%g\n",(double)tao->trust));
576     }
577 
578     if (tao->fmin>-1.e25) {
579       PetscCall(PetscViewerASCIIPrintf(viewer,"convergence tolerances: function minimum=%g\n",(double)tao->fmin));
580     }
581     PetscCall(PetscViewerASCIIPrintf(viewer,"Objective value=%g\n",(double)tao->fc));
582 
583     PetscCall(PetscViewerASCIIPrintf(viewer,"total number of iterations=%D,          ",tao->niter));
584     PetscCall(PetscViewerASCIIPrintf(viewer,"              (max: %D)\n",tao->max_it));
585 
586     if (tao->nfuncs>0) {
587       PetscCall(PetscViewerASCIIPrintf(viewer,"total number of function evaluations=%D,",tao->nfuncs));
588       PetscCall(PetscViewerASCIIPrintf(viewer,"                max: %D\n",tao->max_funcs));
589     }
590     if (tao->ngrads>0) {
591       PetscCall(PetscViewerASCIIPrintf(viewer,"total number of gradient evaluations=%D,",tao->ngrads));
592       PetscCall(PetscViewerASCIIPrintf(viewer,"                max: %D\n",tao->max_funcs));
593     }
594     if (tao->nfuncgrads>0) {
595       PetscCall(PetscViewerASCIIPrintf(viewer,"total number of function/gradient evaluations=%D,",tao->nfuncgrads));
596       PetscCall(PetscViewerASCIIPrintf(viewer,"    (max: %D)\n",tao->max_funcs));
597     }
598     if (tao->nhess>0) {
599       PetscCall(PetscViewerASCIIPrintf(viewer,"total number of Hessian evaluations=%D\n",tao->nhess));
600     }
601     if (tao->nconstraints>0) {
602       PetscCall(PetscViewerASCIIPrintf(viewer,"total number of constraint function evaluations=%D\n",tao->nconstraints));
603     }
604     if (tao->njac>0) {
605       PetscCall(PetscViewerASCIIPrintf(viewer,"total number of Jacobian evaluations=%D\n",tao->njac));
606     }
607 
608     if (tao->reason>0) {
609       PetscCall(PetscViewerASCIIPrintf(viewer,    "Solution converged: "));
610       switch (tao->reason) {
611       case TAO_CONVERGED_GATOL:
612         PetscCall(PetscViewerASCIIPrintf(viewer," ||g(X)|| <= gatol\n"));
613         break;
614       case TAO_CONVERGED_GRTOL:
615         PetscCall(PetscViewerASCIIPrintf(viewer," ||g(X)||/|f(X)| <= grtol\n"));
616         break;
617       case TAO_CONVERGED_GTTOL:
618         PetscCall(PetscViewerASCIIPrintf(viewer," ||g(X)||/||g(X0)|| <= gttol\n"));
619         break;
620       case TAO_CONVERGED_STEPTOL:
621         PetscCall(PetscViewerASCIIPrintf(viewer," Steptol -- step size small\n"));
622         break;
623       case TAO_CONVERGED_MINF:
624         PetscCall(PetscViewerASCIIPrintf(viewer," Minf --  f < fmin\n"));
625         break;
626       case TAO_CONVERGED_USER:
627         PetscCall(PetscViewerASCIIPrintf(viewer," User Terminated\n"));
628         break;
629       default:
630         PetscCall(PetscViewerASCIIPrintf(viewer,"\n"));
631         break;
632       }
633     } else {
634       PetscCall(PetscViewerASCIIPrintf(viewer,"Solver terminated: %d",tao->reason));
635       switch (tao->reason) {
636       case TAO_DIVERGED_MAXITS:
637         PetscCall(PetscViewerASCIIPrintf(viewer," Maximum Iterations\n"));
638         break;
639       case TAO_DIVERGED_NAN:
640         PetscCall(PetscViewerASCIIPrintf(viewer," NAN or Inf encountered\n"));
641         break;
642       case TAO_DIVERGED_MAXFCN:
643         PetscCall(PetscViewerASCIIPrintf(viewer," Maximum Function Evaluations\n"));
644         break;
645       case TAO_DIVERGED_LS_FAILURE:
646         PetscCall(PetscViewerASCIIPrintf(viewer," Line Search Failure\n"));
647         break;
648       case TAO_DIVERGED_TR_REDUCTION:
649         PetscCall(PetscViewerASCIIPrintf(viewer," Trust Region too small\n"));
650         break;
651       case TAO_DIVERGED_USER:
652         PetscCall(PetscViewerASCIIPrintf(viewer," User Terminated\n"));
653         break;
654       default:
655         PetscCall(PetscViewerASCIIPrintf(viewer,"\n"));
656         break;
657       }
658     }
659     PetscCall(PetscViewerASCIIPopTab(viewer));
660   } else if (isstring) {
661     PetscCall(TaoGetType(tao,&type));
662     PetscCall(PetscViewerStringSPrintf(viewer," %-3.3s",type));
663   }
664   PetscFunctionReturn(0);
665 }
666 
667 /*@
668   TaoSetRecycleHistory - Sets the boolean flag to enable/disable re-using
669   iterate information from the previous TaoSolve(). This feature is disabled by
670   default.
671 
672   For conjugate gradient methods (BNCG), this re-uses the latest search direction
673   from the previous TaoSolve() call when computing the first search direction in a
674   new solution. By default, CG methods set the first search direction to the
675   negative gradient.
676 
677   For quasi-Newton family of methods (BQNLS, BQNKLS, BQNKTR, BQNKTL), this re-uses
678   the accumulated quasi-Newton Hessian approximation from the previous TaoSolve()
679   call. By default, QN family of methods reset the initial Hessian approximation to
680   the identity matrix.
681 
682   For any other algorithm, this setting has no effect.
683 
684   Logically collective on Tao
685 
686   Input Parameters:
687 + tao - the Tao context
688 - recycle - boolean flag
689 
690   Options Database Keys:
691 . -tao_recycle_history <true,false> - reuse the history
692 
693   Level: intermediate
694 
695 .seealso: TaoSetRecycleHistory(), TAOBNCG, TAOBQNLS, TAOBQNKLS, TAOBQNKTR, TAOBQNKTL
696 
697 @*/
698 PetscErrorCode TaoSetRecycleHistory(Tao tao, PetscBool recycle)
699 {
700   PetscFunctionBegin;
701   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
702   PetscValidLogicalCollectiveBool(tao,recycle,2);
703   tao->recycle = recycle;
704   PetscFunctionReturn(0);
705 }
706 
707 /*@
708   TaoGetRecycleHistory - Retrieve the boolean flag for re-using iterate information
709   from the previous TaoSolve(). This feature is disabled by default.
710 
711   Logically collective on Tao
712 
713   Input Parameters:
714 . tao - the Tao context
715 
716   Output Parameters:
717 . recycle - boolean flag
718 
719   Level: intermediate
720 
721 .seealso: TaoGetRecycleHistory(), TAOBNCG, TAOBQNLS, TAOBQNKLS, TAOBQNKTR, TAOBQNKTL
722 
723 @*/
724 PetscErrorCode TaoGetRecycleHistory(Tao tao, PetscBool *recycle)
725 {
726   PetscFunctionBegin;
727   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
728   PetscValidBoolPointer(recycle,2);
729   *recycle = tao->recycle;
730   PetscFunctionReturn(0);
731 }
732 
733 /*@
734   TaoSetTolerances - Sets parameters used in TAO convergence tests
735 
736   Logically collective on Tao
737 
738   Input Parameters:
739 + tao - the Tao context
740 . gatol - stop if norm of gradient is less than this
741 . grtol - stop if relative norm of gradient is less than this
742 - gttol - stop if norm of gradient is reduced by this factor
743 
744   Options Database Keys:
745 + -tao_gatol <gatol> - Sets gatol
746 . -tao_grtol <grtol> - Sets grtol
747 - -tao_gttol <gttol> - Sets gttol
748 
749   Stopping Criteria:
750 $ ||g(X)||                            <= gatol
751 $ ||g(X)|| / |f(X)|                   <= grtol
752 $ ||g(X)|| / ||g(X0)||                <= gttol
753 
754   Notes:
755   Use PETSC_DEFAULT to leave one or more tolerances unchanged.
756 
757   Level: beginner
758 
759 .seealso: TaoGetTolerances()
760 
761 @*/
762 PetscErrorCode TaoSetTolerances(Tao tao, PetscReal gatol, PetscReal grtol, PetscReal gttol)
763 {
764   PetscFunctionBegin;
765   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
766   PetscValidLogicalCollectiveReal(tao,gatol,2);
767   PetscValidLogicalCollectiveReal(tao,grtol,3);
768   PetscValidLogicalCollectiveReal(tao,gttol,4);
769 
770   if (gatol != PETSC_DEFAULT) {
771     if (gatol<0) {
772       PetscCall(PetscInfo(tao,"Tried to set negative gatol -- ignored.\n"));
773     } else {
774       tao->gatol = PetscMax(0,gatol);
775       tao->gatol_changed = PETSC_TRUE;
776     }
777   }
778 
779   if (grtol != PETSC_DEFAULT) {
780     if (grtol<0) {
781       PetscCall(PetscInfo(tao,"Tried to set negative grtol -- ignored.\n"));
782     } else {
783       tao->grtol = PetscMax(0,grtol);
784       tao->grtol_changed = PETSC_TRUE;
785     }
786   }
787 
788   if (gttol != PETSC_DEFAULT) {
789     if (gttol<0) {
790       PetscCall(PetscInfo(tao,"Tried to set negative gttol -- ignored.\n"));
791     } else {
792       tao->gttol = PetscMax(0,gttol);
793       tao->gttol_changed = PETSC_TRUE;
794     }
795   }
796   PetscFunctionReturn(0);
797 }
798 
799 /*@
800   TaoSetConstraintTolerances - Sets constraint tolerance parameters used in TAO  convergence tests
801 
802   Logically collective on Tao
803 
804   Input Parameters:
805 + tao - the Tao context
806 . catol - absolute constraint tolerance, constraint norm must be less than catol for used for gatol convergence criteria
807 - crtol - relative constraint tolerance, constraint norm must be less than crtol for used for gatol, gttol convergence criteria
808 
809   Options Database Keys:
810 + -tao_catol <catol> - Sets catol
811 - -tao_crtol <crtol> - Sets crtol
812 
813   Notes:
814   Use PETSC_DEFAULT to leave any tolerance unchanged.
815 
816   Level: intermediate
817 
818 .seealso: TaoGetTolerances(), TaoGetConstraintTolerances(), TaoSetTolerances()
819 
820 @*/
821 PetscErrorCode TaoSetConstraintTolerances(Tao tao, PetscReal catol, PetscReal crtol)
822 {
823   PetscFunctionBegin;
824   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
825   PetscValidLogicalCollectiveReal(tao,catol,2);
826   PetscValidLogicalCollectiveReal(tao,crtol,3);
827 
828   if (catol != PETSC_DEFAULT) {
829     if (catol<0) {
830       PetscCall(PetscInfo(tao,"Tried to set negative catol -- ignored.\n"));
831     } else {
832       tao->catol = PetscMax(0,catol);
833       tao->catol_changed = PETSC_TRUE;
834     }
835   }
836 
837   if (crtol != PETSC_DEFAULT) {
838     if (crtol<0) {
839       PetscCall(PetscInfo(tao,"Tried to set negative crtol -- ignored.\n"));
840     } else {
841       tao->crtol = PetscMax(0,crtol);
842       tao->crtol_changed = PETSC_TRUE;
843     }
844   }
845   PetscFunctionReturn(0);
846 }
847 
848 /*@
849   TaoGetConstraintTolerances - Gets constraint tolerance parameters used in TAO  convergence tests
850 
851   Not ollective
852 
853   Input Parameter:
854 . tao - the Tao context
855 
856   Output Parameters:
857 + catol - absolute constraint tolerance, constraint norm must be less than catol for used for gatol convergence criteria
858 - crtol - relative constraint tolerance, constraint norm must be less than crtol for used for gatol, gttol convergence criteria
859 
860   Level: intermediate
861 
862 .seealso: TaoGetTolerances(), TaoSetTolerances(), TaoSetConstraintTolerances()
863 
864 @*/
865 PetscErrorCode TaoGetConstraintTolerances(Tao tao, PetscReal *catol, PetscReal *crtol)
866 {
867   PetscFunctionBegin;
868   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
869   if (catol) *catol = tao->catol;
870   if (crtol) *crtol = tao->crtol;
871   PetscFunctionReturn(0);
872 }
873 
874 /*@
875    TaoSetFunctionLowerBound - Sets a bound on the solution objective value.
876    When an approximate solution with an objective value below this number
877    has been found, the solver will terminate.
878 
879    Logically Collective on Tao
880 
881    Input Parameters:
882 +  tao - the Tao solver context
883 -  fmin - the tolerance
884 
885    Options Database Keys:
886 .    -tao_fmin <fmin> - sets the minimum function value
887 
888    Level: intermediate
889 
890 .seealso: TaoSetTolerances()
891 @*/
892 PetscErrorCode TaoSetFunctionLowerBound(Tao tao,PetscReal fmin)
893 {
894   PetscFunctionBegin;
895   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
896   PetscValidLogicalCollectiveReal(tao,fmin,2);
897   tao->fmin = fmin;
898   tao->fmin_changed = PETSC_TRUE;
899   PetscFunctionReturn(0);
900 }
901 
902 /*@
903    TaoGetFunctionLowerBound - Gets the bound on the solution objective value.
904    When an approximate solution with an objective value below this number
905    has been found, the solver will terminate.
906 
907    Not collective on Tao
908 
909    Input Parameters:
910 .  tao - the Tao solver context
911 
912    OutputParameters:
913 .  fmin - the minimum function value
914 
915    Level: intermediate
916 
917 .seealso: TaoSetFunctionLowerBound()
918 @*/
919 PetscErrorCode TaoGetFunctionLowerBound(Tao tao,PetscReal *fmin)
920 {
921   PetscFunctionBegin;
922   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
923   PetscValidRealPointer(fmin,2);
924   *fmin = tao->fmin;
925   PetscFunctionReturn(0);
926 }
927 
928 /*@
929    TaoSetMaximumFunctionEvaluations - Sets a maximum number of
930    function evaluations.
931 
932    Logically Collective on Tao
933 
934    Input Parameters:
935 +  tao - the Tao solver context
936 -  nfcn - the maximum number of function evaluations (>=0)
937 
938    Options Database Keys:
939 .    -tao_max_funcs <nfcn> - sets the maximum number of function evaluations
940 
941    Level: intermediate
942 
943 .seealso: TaoSetTolerances(), TaoSetMaximumIterations()
944 @*/
945 
946 PetscErrorCode TaoSetMaximumFunctionEvaluations(Tao tao,PetscInt nfcn)
947 {
948   PetscFunctionBegin;
949   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
950   PetscValidLogicalCollectiveInt(tao,nfcn,2);
951   if (nfcn >= 0) { tao->max_funcs = PetscMax(0,nfcn); }
952   else { tao->max_funcs = -1; }
953   tao->max_funcs_changed = PETSC_TRUE;
954   PetscFunctionReturn(0);
955 }
956 
957 /*@
958    TaoGetMaximumFunctionEvaluations - Sets a maximum number of
959    function evaluations.
960 
961    Not Collective
962 
963    Input Parameters:
964 .  tao - the Tao solver context
965 
966    Output Parameters:
967 .  nfcn - the maximum number of function evaluations
968 
969    Level: intermediate
970 
971 .seealso: TaoSetMaximumFunctionEvaluations(), TaoGetMaximumIterations()
972 @*/
973 
974 PetscErrorCode TaoGetMaximumFunctionEvaluations(Tao tao,PetscInt *nfcn)
975 {
976   PetscFunctionBegin;
977   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
978   PetscValidIntPointer(nfcn,2);
979   *nfcn = tao->max_funcs;
980   PetscFunctionReturn(0);
981 }
982 
983 /*@
984    TaoGetCurrentFunctionEvaluations - Get current number of
985    function evaluations.
986 
987    Not Collective
988 
989    Input Parameters:
990 .  tao - the Tao solver context
991 
992    Output Parameters:
993 .  nfuncs - the current number of function evaluations (maximum between gradient and function evaluations)
994 
995    Level: intermediate
996 
997 .seealso: TaoSetMaximumFunctionEvaluations(), TaoGetMaximumFunctionEvaluations(), TaoGetMaximumIterations()
998 @*/
999 
1000 PetscErrorCode TaoGetCurrentFunctionEvaluations(Tao tao,PetscInt *nfuncs)
1001 {
1002   PetscFunctionBegin;
1003   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1004   PetscValidIntPointer(nfuncs,2);
1005   *nfuncs = PetscMax(tao->nfuncs,tao->nfuncgrads);
1006   PetscFunctionReturn(0);
1007 }
1008 
1009 /*@
1010    TaoSetMaximumIterations - Sets a maximum number of iterates.
1011 
1012    Logically Collective on Tao
1013 
1014    Input Parameters:
1015 +  tao - the Tao solver context
1016 -  maxits - the maximum number of iterates (>=0)
1017 
1018    Options Database Keys:
1019 .    -tao_max_it <its> - sets the maximum number of iterations
1020 
1021    Level: intermediate
1022 
1023 .seealso: TaoSetTolerances(), TaoSetMaximumFunctionEvaluations()
1024 @*/
1025 PetscErrorCode TaoSetMaximumIterations(Tao tao,PetscInt maxits)
1026 {
1027   PetscFunctionBegin;
1028   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1029   PetscValidLogicalCollectiveInt(tao,maxits,2);
1030   tao->max_it = PetscMax(0,maxits);
1031   tao->max_it_changed = PETSC_TRUE;
1032   PetscFunctionReturn(0);
1033 }
1034 
1035 /*@
1036    TaoGetMaximumIterations - Sets a maximum number of iterates.
1037 
1038    Not Collective
1039 
1040    Input Parameters:
1041 .  tao - the Tao solver context
1042 
1043    Output Parameters:
1044 .  maxits - the maximum number of iterates
1045 
1046    Level: intermediate
1047 
1048 .seealso: TaoSetMaximumIterations(), TaoGetMaximumFunctionEvaluations()
1049 @*/
1050 PetscErrorCode TaoGetMaximumIterations(Tao tao,PetscInt *maxits)
1051 {
1052   PetscFunctionBegin;
1053   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1054   PetscValidIntPointer(maxits,2);
1055   *maxits = tao->max_it;
1056   PetscFunctionReturn(0);
1057 }
1058 
1059 /*@
1060    TaoSetInitialTrustRegionRadius - Sets the initial trust region radius.
1061 
1062    Logically collective on Tao
1063 
1064    Input Parameters:
1065 +  tao - a TAO optimization solver
1066 -  radius - the trust region radius
1067 
1068    Level: intermediate
1069 
1070    Options Database Key:
1071 .  -tao_trust0 <t0> - sets initial trust region radius
1072 
1073 .seealso: TaoGetTrustRegionRadius(), TaoSetTrustRegionTolerance()
1074 @*/
1075 PetscErrorCode TaoSetInitialTrustRegionRadius(Tao tao, PetscReal radius)
1076 {
1077   PetscFunctionBegin;
1078   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1079   PetscValidLogicalCollectiveReal(tao,radius,2);
1080   tao->trust0 = PetscMax(0.0,radius);
1081   tao->trust0_changed = PETSC_TRUE;
1082   PetscFunctionReturn(0);
1083 }
1084 
1085 /*@
1086    TaoGetInitialTrustRegionRadius - Sets the initial trust region radius.
1087 
1088    Not Collective
1089 
1090    Input Parameter:
1091 .  tao - a TAO optimization solver
1092 
1093    Output Parameter:
1094 .  radius - the trust region radius
1095 
1096    Level: intermediate
1097 
1098 .seealso: TaoSetInitialTrustRegionRadius(), TaoGetCurrentTrustRegionRadius()
1099 @*/
1100 PetscErrorCode TaoGetInitialTrustRegionRadius(Tao tao, PetscReal *radius)
1101 {
1102   PetscFunctionBegin;
1103   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1104   PetscValidRealPointer(radius,2);
1105   *radius = tao->trust0;
1106   PetscFunctionReturn(0);
1107 }
1108 
1109 /*@
1110    TaoGetCurrentTrustRegionRadius - Gets the current trust region radius.
1111 
1112    Not Collective
1113 
1114    Input Parameter:
1115 .  tao - a TAO optimization solver
1116 
1117    Output Parameter:
1118 .  radius - the trust region radius
1119 
1120    Level: intermediate
1121 
1122 .seealso: TaoSetInitialTrustRegionRadius(), TaoGetInitialTrustRegionRadius()
1123 @*/
1124 PetscErrorCode TaoGetCurrentTrustRegionRadius(Tao tao, PetscReal *radius)
1125 {
1126   PetscFunctionBegin;
1127   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1128   PetscValidRealPointer(radius,2);
1129   *radius = tao->trust;
1130   PetscFunctionReturn(0);
1131 }
1132 
1133 /*@
1134   TaoGetTolerances - gets the current values of tolerances
1135 
1136   Not Collective
1137 
1138   Input Parameter:
1139 . tao - the Tao context
1140 
1141   Output Parameters:
1142 + gatol - stop if norm of gradient is less than this
1143 . grtol - stop if relative norm of gradient is less than this
1144 - gttol - stop if norm of gradient is reduced by a this factor
1145 
1146   Note: NULL can be used as an argument if not all tolerances values are needed
1147 
1148 .seealso TaoSetTolerances()
1149 
1150   Level: intermediate
1151 @*/
1152 PetscErrorCode TaoGetTolerances(Tao tao, PetscReal *gatol, PetscReal *grtol, PetscReal *gttol)
1153 {
1154   PetscFunctionBegin;
1155   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1156   if (gatol) *gatol = tao->gatol;
1157   if (grtol) *grtol = tao->grtol;
1158   if (gttol) *gttol = tao->gttol;
1159   PetscFunctionReturn(0);
1160 }
1161 
1162 /*@
1163   TaoGetKSP - Gets the linear solver used by the optimization solver.
1164   Application writers should use TaoGetKSP if they need direct access
1165   to the PETSc KSP object.
1166 
1167   Not Collective
1168 
1169    Input Parameters:
1170 .  tao - the TAO solver
1171 
1172    Output Parameters:
1173 .  ksp - the KSP linear solver used in the optimization solver
1174 
1175    Level: intermediate
1176 
1177 @*/
1178 PetscErrorCode TaoGetKSP(Tao tao, KSP *ksp)
1179 {
1180   PetscFunctionBegin;
1181   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1182   PetscValidPointer(ksp,2);
1183   *ksp = tao->ksp;
1184   PetscFunctionReturn(0);
1185 }
1186 
1187 /*@
1188    TaoGetLinearSolveIterations - Gets the total number of linear iterations
1189    used by the TAO solver
1190 
1191    Not Collective
1192 
1193    Input Parameter:
1194 .  tao - TAO context
1195 
1196    Output Parameter:
1197 .  lits - number of linear iterations
1198 
1199    Notes:
1200    This counter is reset to zero for each successive call to TaoSolve()
1201 
1202    Level: intermediate
1203 
1204 .seealso:  TaoGetKSP()
1205 @*/
1206 PetscErrorCode TaoGetLinearSolveIterations(Tao tao, PetscInt *lits)
1207 {
1208   PetscFunctionBegin;
1209   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1210   PetscValidIntPointer(lits,2);
1211   *lits = tao->ksp_tot_its;
1212   PetscFunctionReturn(0);
1213 }
1214 
1215 /*@
1216   TaoGetLineSearch - Gets the line search used by the optimization solver.
1217   Application writers should use TaoGetLineSearch if they need direct access
1218   to the TaoLineSearch object.
1219 
1220   Not Collective
1221 
1222    Input Parameters:
1223 .  tao - the TAO solver
1224 
1225    Output Parameters:
1226 .  ls - the line search used in the optimization solver
1227 
1228    Level: intermediate
1229 
1230 @*/
1231 PetscErrorCode TaoGetLineSearch(Tao tao, TaoLineSearch *ls)
1232 {
1233   PetscFunctionBegin;
1234   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1235   PetscValidPointer(ls,2);
1236   *ls = tao->linesearch;
1237   PetscFunctionReturn(0);
1238 }
1239 
1240 /*@
1241   TaoAddLineSearchCounts - Adds the number of function evaluations spent
1242   in the line search to the running total.
1243 
1244    Input Parameters:
1245 +  tao - the TAO solver
1246 -  ls - the line search used in the optimization solver
1247 
1248    Level: developer
1249 
1250 .seealso: TaoLineSearchApply()
1251 @*/
1252 PetscErrorCode TaoAddLineSearchCounts(Tao tao)
1253 {
1254   PetscBool      flg;
1255   PetscInt       nfeval,ngeval,nfgeval;
1256 
1257   PetscFunctionBegin;
1258   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1259   if (tao->linesearch) {
1260     PetscCall(TaoLineSearchIsUsingTaoRoutines(tao->linesearch,&flg));
1261     if (!flg) {
1262       PetscCall(TaoLineSearchGetNumberFunctionEvaluations(tao->linesearch,&nfeval,&ngeval,&nfgeval));
1263       tao->nfuncs += nfeval;
1264       tao->ngrads += ngeval;
1265       tao->nfuncgrads += nfgeval;
1266     }
1267   }
1268   PetscFunctionReturn(0);
1269 }
1270 
1271 /*@
1272   TaoGetSolution - Returns the vector with the current TAO solution
1273 
1274   Not Collective
1275 
1276   Input Parameter:
1277 . tao - the Tao context
1278 
1279   Output Parameter:
1280 . X - the current solution
1281 
1282   Level: intermediate
1283 
1284   Note:  The returned vector will be the same object that was passed into TaoSetSolution()
1285 @*/
1286 PetscErrorCode TaoGetSolution(Tao tao, Vec *X)
1287 {
1288   PetscFunctionBegin;
1289   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1290   PetscValidPointer(X,2);
1291   *X = tao->solution;
1292   PetscFunctionReturn(0);
1293 }
1294 
1295 /*@
1296    TaoResetStatistics - Initialize the statistics used by TAO for all of the solvers.
1297    These statistics include the iteration number, residual norms, and convergence status.
1298    This routine gets called before solving each optimization problem.
1299 
1300    Collective on Tao
1301 
1302    Input Parameters:
1303 .  solver - the Tao context
1304 
1305    Level: developer
1306 
1307 .seealso: TaoCreate(), TaoSolve()
1308 @*/
1309 PetscErrorCode TaoResetStatistics(Tao tao)
1310 {
1311   PetscFunctionBegin;
1312   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1313   tao->niter        = 0;
1314   tao->nfuncs       = 0;
1315   tao->nfuncgrads   = 0;
1316   tao->ngrads       = 0;
1317   tao->nhess        = 0;
1318   tao->njac         = 0;
1319   tao->nconstraints = 0;
1320   tao->ksp_its      = 0;
1321   tao->ksp_tot_its  = 0;
1322   tao->reason       = TAO_CONTINUE_ITERATING;
1323   tao->residual     = 0.0;
1324   tao->cnorm        = 0.0;
1325   tao->step         = 0.0;
1326   tao->lsflag       = PETSC_FALSE;
1327   if (tao->hist_reset) tao->hist_len = 0;
1328   PetscFunctionReturn(0);
1329 }
1330 
1331 /*@C
1332   TaoSetUpdate - Sets the general-purpose update function called
1333   at the beginning of every iteration of the nonlinear solve. Specifically
1334   it is called at the top of every iteration, after the new solution and the gradient
1335   is determined, but before the Hessian is computed (if applicable).
1336 
1337   Logically Collective on Tao
1338 
1339   Input Parameters:
1340 + tao - The tao solver context
1341 - func - The function
1342 
1343   Calling sequence of func:
1344 $ func (Tao tao, PetscInt step);
1345 
1346 . step - The current step of the iteration
1347 
1348   Level: advanced
1349 
1350 .seealso TaoSolve()
1351 @*/
1352 PetscErrorCode TaoSetUpdate(Tao tao, PetscErrorCode (*func)(Tao, PetscInt, void*), void *ctx)
1353 {
1354   PetscFunctionBegin;
1355   PetscValidHeaderSpecific(tao, TAO_CLASSID,1);
1356   tao->ops->update = func;
1357   tao->user_update = ctx;
1358   PetscFunctionReturn(0);
1359 }
1360 
1361 /*@C
1362   TaoSetConvergenceTest - Sets the function that is to be used to test
1363   for convergence o fthe iterative minimization solution.  The new convergence
1364   testing routine will replace TAO's default convergence test.
1365 
1366   Logically Collective on Tao
1367 
1368   Input Parameters:
1369 + tao - the Tao object
1370 . conv - the routine to test for convergence
1371 - ctx - [optional] context for private data for the convergence routine
1372         (may be NULL)
1373 
1374   Calling sequence of conv:
1375 $   PetscErrorCode conv(Tao tao, void *ctx)
1376 
1377 + tao - the Tao object
1378 - ctx - [optional] convergence context
1379 
1380   Note: The new convergence testing routine should call TaoSetConvergedReason().
1381 
1382   Level: advanced
1383 
1384 .seealso: TaoSetConvergedReason(), TaoGetSolutionStatus(), TaoGetTolerances(), TaoSetMonitor
1385 
1386 @*/
1387 PetscErrorCode TaoSetConvergenceTest(Tao tao, PetscErrorCode (*conv)(Tao, void*), void *ctx)
1388 {
1389   PetscFunctionBegin;
1390   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1391   tao->ops->convergencetest = conv;
1392   tao->cnvP = ctx;
1393   PetscFunctionReturn(0);
1394 }
1395 
1396 /*@C
1397    TaoSetMonitor - Sets an ADDITIONAL function that is to be used at every
1398    iteration of the solver to display the iteration's
1399    progress.
1400 
1401    Logically Collective on Tao
1402 
1403    Input Parameters:
1404 +  tao - the Tao solver context
1405 .  mymonitor - monitoring routine
1406 -  mctx - [optional] user-defined context for private data for the
1407           monitor routine (may be NULL)
1408 
1409    Calling sequence of mymonitor:
1410 .vb
1411      PetscErrorCode mymonitor(Tao tao,void *mctx)
1412 .ve
1413 
1414 +    tao - the Tao solver context
1415 -    mctx - [optional] monitoring context
1416 
1417    Options Database Keys:
1418 +    -tao_monitor        - sets TaoMonitorDefault()
1419 .    -tao_smonitor       - sets short monitor
1420 .    -tao_cmonitor       - same as smonitor plus constraint norm
1421 .    -tao_view_solution   - view solution at each iteration
1422 .    -tao_view_gradient   - view gradient at each iteration
1423 .    -tao_view_ls_residual - view least-squares residual vector at each iteration
1424 -    -tao_cancelmonitors - cancels all monitors that have been hardwired into a code by calls to TaoSetMonitor(), but does not cancel those set via the options database.
1425 
1426    Notes:
1427    Several different monitoring routines may be set by calling
1428    TaoSetMonitor() multiple times; all will be called in the
1429    order in which they were set.
1430 
1431    Fortran Notes:
1432     Only one monitor function may be set
1433 
1434    Level: intermediate
1435 
1436 .seealso: TaoMonitorDefault(), TaoCancelMonitors(),  TaoSetDestroyRoutine()
1437 @*/
1438 PetscErrorCode TaoSetMonitor(Tao tao, PetscErrorCode (*func)(Tao, void*), void *ctx,PetscErrorCode (*dest)(void**))
1439 {
1440   PetscInt       i;
1441   PetscBool      identical;
1442 
1443   PetscFunctionBegin;
1444   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1445   PetscCheck(tao->numbermonitors < MAXTAOMONITORS,PetscObjectComm((PetscObject)tao),PETSC_ERR_SUP,"Cannot attach another monitor -- max=%d",MAXTAOMONITORS);
1446 
1447   for (i=0; i<tao->numbermonitors;i++) {
1448     PetscCall(PetscMonitorCompare((PetscErrorCode (*)(void))func,ctx,dest,(PetscErrorCode (*)(void))tao->monitor[i],tao->monitorcontext[i],tao->monitordestroy[i],&identical));
1449     if (identical) PetscFunctionReturn(0);
1450   }
1451   tao->monitor[tao->numbermonitors] = func;
1452   tao->monitorcontext[tao->numbermonitors] = (void*)ctx;
1453   tao->monitordestroy[tao->numbermonitors] = dest;
1454   ++tao->numbermonitors;
1455   PetscFunctionReturn(0);
1456 }
1457 
1458 /*@
1459    TaoCancelMonitors - Clears all the monitor functions for a Tao object.
1460 
1461    Logically Collective on Tao
1462 
1463    Input Parameters:
1464 .  tao - the Tao solver context
1465 
1466    Options Database:
1467 .  -tao_cancelmonitors - cancels all monitors that have been hardwired
1468     into a code by calls to TaoSetMonitor(), but does not cancel those
1469     set via the options database
1470 
1471    Notes:
1472    There is no way to clear one specific monitor from a Tao object.
1473 
1474    Level: advanced
1475 
1476 .seealso: TaoMonitorDefault(), TaoSetMonitor()
1477 @*/
1478 PetscErrorCode TaoCancelMonitors(Tao tao)
1479 {
1480   PetscInt       i;
1481 
1482   PetscFunctionBegin;
1483   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1484   for (i=0;i<tao->numbermonitors;i++) {
1485     if (tao->monitordestroy[i]) {
1486       PetscCall((*tao->monitordestroy[i])(&tao->monitorcontext[i]));
1487     }
1488   }
1489   tao->numbermonitors = 0;
1490   PetscFunctionReturn(0);
1491 }
1492 
1493 /*@
1494    TaoMonitorDefault - Default routine for monitoring progress of the
1495    Tao solvers (default).  This monitor prints the function value and gradient
1496    norm at each iteration.  It can be turned on from the command line using the
1497    -tao_monitor option
1498 
1499    Collective on Tao
1500 
1501    Input Parameters:
1502 +  tao - the Tao context
1503 -  ctx - PetscViewer context or NULL
1504 
1505    Options Database Keys:
1506 .  -tao_monitor - turn on default monitoring
1507 
1508    Level: advanced
1509 
1510 .seealso: TaoDefaultSMonitor(), TaoSetMonitor()
1511 @*/
1512 PetscErrorCode TaoMonitorDefault(Tao tao, void *ctx)
1513 {
1514   PetscInt       its, tabs;
1515   PetscReal      fct,gnorm;
1516   PetscViewer    viewer = (PetscViewer)ctx;
1517 
1518   PetscFunctionBegin;
1519   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1520   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
1521   its = tao->niter;
1522   fct = tao->fc;
1523   gnorm = tao->residual;
1524   PetscCall(PetscViewerASCIIGetTab(viewer, &tabs));
1525   PetscCall(PetscViewerASCIISetTab(viewer, ((PetscObject)tao)->tablevel));
1526   if (its == 0 && ((PetscObject)tao)->prefix && !tao->header_printed) {
1527      PetscCall(PetscViewerASCIIPrintf(viewer,"  Iteration information for %s solve.\n",((PetscObject)tao)->prefix));
1528      tao->header_printed = PETSC_TRUE;
1529    }
1530   PetscCall(PetscViewerASCIIPrintf(viewer,"%3D TAO,",its));
1531   PetscCall(PetscViewerASCIIPrintf(viewer,"  Function value: %g,",(double)fct));
1532   if (gnorm >= PETSC_INFINITY) {
1533     PetscCall(PetscViewerASCIIPrintf(viewer,"  Residual: Inf \n"));
1534   } else {
1535     PetscCall(PetscViewerASCIIPrintf(viewer,"  Residual: %g \n",(double)gnorm));
1536   }
1537   PetscCall(PetscViewerASCIISetTab(viewer, tabs));
1538   PetscFunctionReturn(0);
1539 }
1540 
1541 /*@
1542    TaoDefaultGMonitor - Default routine for monitoring progress of the
1543    Tao solvers (default) with extra detail on the globalization method.
1544    This monitor prints the function value and gradient norm at each
1545    iteration, as well as the step size and trust radius. Note that the
1546    step size and trust radius may be the same for some algorithms.
1547    It can be turned on from the command line using the
1548    -tao_gmonitor option
1549 
1550    Collective on Tao
1551 
1552    Input Parameters:
1553 +  tao - the Tao context
1554 -  ctx - PetscViewer context or NULL
1555 
1556    Options Database Keys:
1557 .  -tao_gmonitor - turn on monitoring with globalization information
1558 
1559    Level: advanced
1560 
1561 .seealso: TaoDefaultSMonitor(), TaoSetMonitor()
1562 @*/
1563 PetscErrorCode TaoDefaultGMonitor(Tao tao, void *ctx)
1564 {
1565   PetscInt       its, tabs;
1566   PetscReal      fct,gnorm,stp,tr;
1567   PetscViewer    viewer = (PetscViewer)ctx;
1568 
1569   PetscFunctionBegin;
1570   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1571   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
1572   its = tao->niter;
1573   fct = tao->fc;
1574   gnorm = tao->residual;
1575   stp = tao->step;
1576   tr = tao->trust;
1577   PetscCall(PetscViewerASCIIGetTab(viewer, &tabs));
1578   PetscCall(PetscViewerASCIISetTab(viewer, ((PetscObject)tao)->tablevel));
1579   if (its == 0 && ((PetscObject)tao)->prefix && !tao->header_printed) {
1580      PetscCall(PetscViewerASCIIPrintf(viewer,"  Iteration information for %s solve.\n",((PetscObject)tao)->prefix));
1581      tao->header_printed = PETSC_TRUE;
1582    }
1583   PetscCall(PetscViewerASCIIPrintf(viewer,"%3D TAO,",its));
1584   PetscCall(PetscViewerASCIIPrintf(viewer,"  Function value: %g,",(double)fct));
1585   if (gnorm >= PETSC_INFINITY) {
1586     PetscCall(PetscViewerASCIIPrintf(viewer,"  Residual: Inf,"));
1587   } else {
1588     PetscCall(PetscViewerASCIIPrintf(viewer,"  Residual: %g,",(double)gnorm));
1589   }
1590   PetscCall(PetscViewerASCIIPrintf(viewer,"  Step: %g,  Trust: %g\n",(double)stp,(double)tr));
1591   PetscCall(PetscViewerASCIISetTab(viewer, tabs));
1592   PetscFunctionReturn(0);
1593 }
1594 
1595 /*@
1596    TaoDefaultSMonitor - Default routine for monitoring progress of the
1597    solver. Same as TaoMonitorDefault() except
1598    it prints fewer digits of the residual as the residual gets smaller.
1599    This is because the later digits are meaningless and are often
1600    different on different machines; by using this routine different
1601    machines will usually generate the same output. It can be turned on
1602    by using the -tao_smonitor option
1603 
1604    Collective on Tao
1605 
1606    Input Parameters:
1607 +  tao - the Tao context
1608 -  ctx - PetscViewer context of type ASCII
1609 
1610    Options Database Keys:
1611 .  -tao_smonitor - turn on default short monitoring
1612 
1613    Level: advanced
1614 
1615 .seealso: TaoMonitorDefault(), TaoSetMonitor()
1616 @*/
1617 PetscErrorCode TaoDefaultSMonitor(Tao tao, void *ctx)
1618 {
1619   PetscInt       its, tabs;
1620   PetscReal      fct,gnorm;
1621   PetscViewer    viewer = (PetscViewer)ctx;
1622 
1623   PetscFunctionBegin;
1624   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1625   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
1626   its = tao->niter;
1627   fct = tao->fc;
1628   gnorm = tao->residual;
1629   PetscCall(PetscViewerASCIIGetTab(viewer, &tabs));
1630   PetscCall(PetscViewerASCIISetTab(viewer, ((PetscObject)tao)->tablevel));
1631   PetscCall(PetscViewerASCIIPrintf(viewer,"iter = %3D,",its));
1632   PetscCall(PetscViewerASCIIPrintf(viewer," Function value %g,",(double)fct));
1633   if (gnorm >= PETSC_INFINITY) {
1634     PetscCall(PetscViewerASCIIPrintf(viewer," Residual: Inf \n"));
1635   } else if (gnorm > 1.e-6) {
1636     PetscCall(PetscViewerASCIIPrintf(viewer," Residual: %g \n",(double)gnorm));
1637   } else if (gnorm > 1.e-11) {
1638     PetscCall(PetscViewerASCIIPrintf(viewer," Residual: < 1.0e-6 \n"));
1639   } else {
1640     PetscCall(PetscViewerASCIIPrintf(viewer," Residual: < 1.0e-11 \n"));
1641   }
1642   PetscCall(PetscViewerASCIISetTab(viewer, tabs));
1643   PetscFunctionReturn(0);
1644 }
1645 
1646 /*@
1647    TaoDefaultCMonitor - same as TaoMonitorDefault() except
1648    it prints the norm of the constraints function. It can be turned on
1649    from the command line using the -tao_cmonitor option
1650 
1651    Collective on Tao
1652 
1653    Input Parameters:
1654 +  tao - the Tao context
1655 -  ctx - PetscViewer context or NULL
1656 
1657    Options Database Keys:
1658 .  -tao_cmonitor - monitor the constraints
1659 
1660    Level: advanced
1661 
1662 .seealso: TaoMonitorDefault(), TaoSetMonitor()
1663 @*/
1664 PetscErrorCode TaoDefaultCMonitor(Tao tao, void *ctx)
1665 {
1666   PetscInt       its, tabs;
1667   PetscReal      fct,gnorm;
1668   PetscViewer    viewer = (PetscViewer)ctx;
1669 
1670   PetscFunctionBegin;
1671   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1672   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
1673   its = tao->niter;
1674   fct = tao->fc;
1675   gnorm = tao->residual;
1676   PetscCall(PetscViewerASCIIGetTab(viewer, &tabs));
1677   PetscCall(PetscViewerASCIISetTab(viewer, ((PetscObject)tao)->tablevel));
1678   PetscCall(PetscViewerASCIIPrintf(viewer,"iter = %D,",its));
1679   PetscCall(PetscViewerASCIIPrintf(viewer," Function value: %g,",(double)fct));
1680   PetscCall(PetscViewerASCIIPrintf(viewer,"  Residual: %g ",(double)gnorm));
1681   PetscCall(PetscViewerASCIIPrintf(viewer,"  Constraint: %g \n",(double)tao->cnorm));
1682   PetscCall(PetscViewerASCIISetTab(viewer, tabs));
1683   PetscFunctionReturn(0);
1684 }
1685 
1686 /*@C
1687    TaoSolutionMonitor - Views the solution at each iteration
1688    It can be turned on from the command line using the
1689    -tao_view_solution option
1690 
1691    Collective on Tao
1692 
1693    Input Parameters:
1694 +  tao - the Tao context
1695 -  ctx - PetscViewer context or NULL
1696 
1697    Options Database Keys:
1698 .  -tao_view_solution - view the solution
1699 
1700    Level: advanced
1701 
1702 .seealso: TaoDefaultSMonitor(), TaoSetMonitor()
1703 @*/
1704 PetscErrorCode TaoSolutionMonitor(Tao tao, void *ctx)
1705 {
1706   PetscViewer    viewer  = (PetscViewer)ctx;
1707 
1708   PetscFunctionBegin;
1709   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1710   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
1711   PetscCall(VecView(tao->solution,viewer));
1712   PetscFunctionReturn(0);
1713 }
1714 
1715 /*@C
1716    TaoGradientMonitor - Views the gradient at each iteration
1717    It can be turned on from the command line using the
1718    -tao_view_gradient option
1719 
1720    Collective on Tao
1721 
1722    Input Parameters:
1723 +  tao - the Tao context
1724 -  ctx - PetscViewer context or NULL
1725 
1726    Options Database Keys:
1727 .  -tao_view_gradient - view the gradient at each iteration
1728 
1729    Level: advanced
1730 
1731 .seealso: TaoDefaultSMonitor(), TaoSetMonitor()
1732 @*/
1733 PetscErrorCode TaoGradientMonitor(Tao tao, void *ctx)
1734 {
1735   PetscViewer    viewer = (PetscViewer)ctx;
1736 
1737   PetscFunctionBegin;
1738   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1739   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
1740   PetscCall(VecView(tao->gradient,viewer));
1741   PetscFunctionReturn(0);
1742 }
1743 
1744 /*@C
1745    TaoStepDirectionMonitor - Views the step-direction at each iteration
1746 
1747    Collective on Tao
1748 
1749    Input Parameters:
1750 +  tao - the Tao context
1751 -  ctx - PetscViewer context or NULL
1752 
1753    Options Database Keys:
1754 .  -tao_view_gradient - view the gradient at each iteration
1755 
1756    Level: advanced
1757 
1758 .seealso: TaoDefaultSMonitor(), TaoSetMonitor()
1759 @*/
1760 PetscErrorCode TaoStepDirectionMonitor(Tao tao, void *ctx)
1761 {
1762   PetscViewer    viewer = (PetscViewer)ctx;
1763 
1764   PetscFunctionBegin;
1765   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1766   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
1767   PetscCall(VecView(tao->stepdirection,viewer));
1768   PetscFunctionReturn(0);
1769 }
1770 
1771 /*@C
1772    TaoDrawSolutionMonitor - Plots the solution at each iteration
1773    It can be turned on from the command line using the
1774    -tao_draw_solution option
1775 
1776    Collective on Tao
1777 
1778    Input Parameters:
1779 +  tao - the Tao context
1780 -  ctx - TaoMonitorDraw context
1781 
1782    Options Database Keys:
1783 .  -tao_draw_solution - draw the solution at each iteration
1784 
1785    Level: advanced
1786 
1787 .seealso: TaoSolutionMonitor(), TaoSetMonitor(), TaoDrawGradientMonitor
1788 @*/
1789 PetscErrorCode TaoDrawSolutionMonitor(Tao tao, void *ctx)
1790 {
1791   TaoMonitorDrawCtx ictx = (TaoMonitorDrawCtx)ctx;
1792 
1793   PetscFunctionBegin;
1794   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1795   if (!(((ictx->howoften > 0) && (!(tao->niter % ictx->howoften))) || ((ictx->howoften == -1) && tao->reason))) PetscFunctionReturn(0);
1796   PetscCall(VecView(tao->solution,ictx->viewer));
1797   PetscFunctionReturn(0);
1798 }
1799 
1800 /*@C
1801    TaoDrawGradientMonitor - Plots the gradient at each iteration
1802    It can be turned on from the command line using the
1803    -tao_draw_gradient option
1804 
1805    Collective on Tao
1806 
1807    Input Parameters:
1808 +  tao - the Tao context
1809 -  ctx - PetscViewer context
1810 
1811    Options Database Keys:
1812 .  -tao_draw_gradient - draw the gradient at each iteration
1813 
1814    Level: advanced
1815 
1816 .seealso: TaoGradientMonitor(), TaoSetMonitor(), TaoDrawSolutionMonitor
1817 @*/
1818 PetscErrorCode TaoDrawGradientMonitor(Tao tao, void *ctx)
1819 {
1820   TaoMonitorDrawCtx ictx = (TaoMonitorDrawCtx)ctx;
1821 
1822   PetscFunctionBegin;
1823   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1824   if (!(((ictx->howoften > 0) && (!(tao->niter % ictx->howoften))) || ((ictx->howoften == -1) && tao->reason))) PetscFunctionReturn(0);
1825   PetscCall(VecView(tao->gradient,ictx->viewer));
1826   PetscFunctionReturn(0);
1827 }
1828 
1829 /*@C
1830    TaoDrawStepMonitor - Plots the step direction at each iteration
1831 
1832    Collective on Tao
1833 
1834    Input Parameters:
1835 +  tao - the Tao context
1836 -  ctx - PetscViewer context
1837 
1838    Options Database Keys:
1839 .  -tao_draw_step - draw the step direction at each iteration
1840 
1841    Level: advanced
1842 
1843 .seealso: TaoSetMonitor(), TaoDrawSolutionMonitor
1844 @*/
1845 PetscErrorCode TaoDrawStepMonitor(Tao tao, void *ctx)
1846 {
1847   PetscViewer    viewer = (PetscViewer)ctx;
1848 
1849   PetscFunctionBegin;
1850   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1851   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
1852   PetscCall(VecView(tao->stepdirection,viewer));
1853   PetscFunctionReturn(0);
1854 }
1855 
1856 /*@C
1857    TaoResidualMonitor - Views the least-squares residual at each iteration
1858 
1859    Collective on Tao
1860 
1861    Input Parameters:
1862 +  tao - the Tao context
1863 -  ctx - PetscViewer context or NULL
1864 
1865    Options Database Keys:
1866 .  -tao_view_ls_residual - view the least-squares residual at each iteration
1867 
1868    Level: advanced
1869 
1870 .seealso: TaoDefaultSMonitor(), TaoSetMonitor()
1871 @*/
1872 PetscErrorCode TaoResidualMonitor(Tao tao, void *ctx)
1873 {
1874   PetscViewer    viewer  = (PetscViewer)ctx;
1875 
1876   PetscFunctionBegin;
1877   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1878   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
1879   PetscCall(VecView(tao->ls_res,viewer));
1880   PetscFunctionReturn(0);
1881 }
1882 
1883 /*@
1884    TaoDefaultConvergenceTest - Determines whether the solver should continue iterating
1885    or terminate.
1886 
1887    Collective on Tao
1888 
1889    Input Parameters:
1890 +  tao - the Tao context
1891 -  dummy - unused dummy context
1892 
1893    Output Parameter:
1894 .  reason - for terminating
1895 
1896    Notes:
1897    This routine checks the residual in the optimality conditions, the
1898    relative residual in the optimity conditions, the number of function
1899    evaluations, and the function value to test convergence.  Some
1900    solvers may use different convergence routines.
1901 
1902    Level: developer
1903 
1904 .seealso: TaoSetTolerances(),TaoGetConvergedReason(),TaoSetConvergedReason()
1905 @*/
1906 
1907 PetscErrorCode TaoDefaultConvergenceTest(Tao tao,void *dummy)
1908 {
1909   PetscInt           niter=tao->niter, nfuncs=PetscMax(tao->nfuncs,tao->nfuncgrads);
1910   PetscInt           max_funcs=tao->max_funcs;
1911   PetscReal          gnorm=tao->residual, gnorm0=tao->gnorm0;
1912   PetscReal          f=tao->fc, steptol=tao->steptol,trradius=tao->step;
1913   PetscReal          gatol=tao->gatol,grtol=tao->grtol,gttol=tao->gttol;
1914   PetscReal          catol=tao->catol,crtol=tao->crtol;
1915   PetscReal          fmin=tao->fmin, cnorm=tao->cnorm;
1916   TaoConvergedReason reason=tao->reason;
1917 
1918   PetscFunctionBegin;
1919   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1920   if (reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
1921 
1922   if (PetscIsInfOrNanReal(f)) {
1923     PetscCall(PetscInfo(tao,"Failed to converged, function value is Inf or NaN\n"));
1924     reason = TAO_DIVERGED_NAN;
1925   } else if (f <= fmin && cnorm <=catol) {
1926     PetscCall(PetscInfo(tao,"Converged due to function value %g < minimum function value %g\n", (double)f,(double)fmin));
1927     reason = TAO_CONVERGED_MINF;
1928   } else if (gnorm<= gatol && cnorm <=catol) {
1929     PetscCall(PetscInfo(tao,"Converged due to residual norm ||g(X)||=%g < %g\n",(double)gnorm,(double)gatol));
1930     reason = TAO_CONVERGED_GATOL;
1931   } else if (f!=0 && PetscAbsReal(gnorm/f) <= grtol && cnorm <= crtol) {
1932     PetscCall(PetscInfo(tao,"Converged due to residual ||g(X)||/|f(X)| =%g < %g\n",(double)(gnorm/f),(double)grtol));
1933     reason = TAO_CONVERGED_GRTOL;
1934   } else if (gnorm0 != 0 && ((gttol == 0 && gnorm == 0) || gnorm/gnorm0 < gttol) && cnorm <= crtol) {
1935     PetscCall(PetscInfo(tao,"Converged due to relative residual norm ||g(X)||/||g(X0)|| = %g < %g\n",(double)(gnorm/gnorm0),(double)gttol));
1936     reason = TAO_CONVERGED_GTTOL;
1937   } else if (max_funcs >=0 && nfuncs > max_funcs) {
1938     PetscCall(PetscInfo(tao,"Exceeded maximum number of function evaluations: %" PetscInt_FMT " > %" PetscInt_FMT "\n", nfuncs,max_funcs));
1939     reason = TAO_DIVERGED_MAXFCN;
1940   } else if (tao->lsflag != 0) {
1941     PetscCall(PetscInfo(tao,"Tao Line Search failure.\n"));
1942     reason = TAO_DIVERGED_LS_FAILURE;
1943   } else if (trradius < steptol && niter > 0) {
1944     PetscCall(PetscInfo(tao,"Trust region/step size too small: %g < %g\n", (double)trradius,(double)steptol));
1945     reason = TAO_CONVERGED_STEPTOL;
1946   } else if (niter >= tao->max_it) {
1947     PetscCall(PetscInfo(tao,"Exceeded maximum number of iterations: %D > %D\n",niter,tao->max_it));
1948     reason = TAO_DIVERGED_MAXITS;
1949   } else {
1950     reason = TAO_CONTINUE_ITERATING;
1951   }
1952   tao->reason = reason;
1953   PetscFunctionReturn(0);
1954 }
1955 
1956 /*@C
1957    TaoSetOptionsPrefix - Sets the prefix used for searching for all
1958    TAO options in the database.
1959 
1960    Logically Collective on Tao
1961 
1962    Input Parameters:
1963 +  tao - the Tao context
1964 -  prefix - the prefix string to prepend to all TAO option requests
1965 
1966    Notes:
1967    A hyphen (-) must NOT be given at the beginning of the prefix name.
1968    The first character of all runtime options is AUTOMATICALLY the hyphen.
1969 
1970    For example, to distinguish between the runtime options for two
1971    different TAO solvers, one could call
1972 .vb
1973       TaoSetOptionsPrefix(tao1,"sys1_")
1974       TaoSetOptionsPrefix(tao2,"sys2_")
1975 .ve
1976 
1977    This would enable use of different options for each system, such as
1978 .vb
1979       -sys1_tao_method blmvm -sys1_tao_grtol 1.e-3
1980       -sys2_tao_method lmvm  -sys2_tao_grtol 1.e-4
1981 .ve
1982 
1983    Level: advanced
1984 
1985 .seealso: TaoAppendOptionsPrefix(), TaoGetOptionsPrefix()
1986 @*/
1987 
1988 PetscErrorCode TaoSetOptionsPrefix(Tao tao, const char p[])
1989 {
1990   PetscFunctionBegin;
1991   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1992   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)tao,p));
1993   if (tao->linesearch) {
1994     PetscCall(TaoLineSearchSetOptionsPrefix(tao->linesearch,p));
1995   }
1996   if (tao->ksp) {
1997     PetscCall(KSPSetOptionsPrefix(tao->ksp,p));
1998   }
1999   PetscFunctionReturn(0);
2000 }
2001 
2002 /*@C
2003    TaoAppendOptionsPrefix - Appends to the prefix used for searching for all
2004    TAO options in the database.
2005 
2006    Logically Collective on Tao
2007 
2008    Input Parameters:
2009 +  tao - the Tao solver context
2010 -  prefix - the prefix string to prepend to all TAO option requests
2011 
2012    Notes:
2013    A hyphen (-) must NOT be given at the beginning of the prefix name.
2014    The first character of all runtime options is AUTOMATICALLY the hyphen.
2015 
2016    Level: advanced
2017 
2018 .seealso: TaoSetOptionsPrefix(), TaoGetOptionsPrefix()
2019 @*/
2020 PetscErrorCode TaoAppendOptionsPrefix(Tao tao, const char p[])
2021 {
2022   PetscFunctionBegin;
2023   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2024   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)tao,p));
2025   if (tao->linesearch) {
2026     PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)tao->linesearch,p));
2027   }
2028   if (tao->ksp) {
2029     PetscCall(KSPAppendOptionsPrefix(tao->ksp,p));
2030   }
2031   PetscFunctionReturn(0);
2032 }
2033 
2034 /*@C
2035   TaoGetOptionsPrefix - Gets the prefix used for searching for all
2036   TAO options in the database
2037 
2038   Not Collective
2039 
2040   Input Parameters:
2041 . tao - the Tao context
2042 
2043   Output Parameters:
2044 . prefix - pointer to the prefix string used is returned
2045 
2046   Notes:
2047     On the fortran side, the user should pass in a string 'prefix' of
2048   sufficient length to hold the prefix.
2049 
2050   Level: advanced
2051 
2052 .seealso: TaoSetOptionsPrefix(), TaoAppendOptionsPrefix()
2053 @*/
2054 PetscErrorCode TaoGetOptionsPrefix(Tao tao, const char *p[])
2055 {
2056   PetscFunctionBegin;
2057   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2058   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)tao,p));
2059   PetscFunctionReturn(0);
2060 }
2061 
2062 /*@C
2063    TaoSetType - Sets the method for the unconstrained minimization solver.
2064 
2065    Collective on Tao
2066 
2067    Input Parameters:
2068 +  solver - the Tao solver context
2069 -  type - a known method
2070 
2071    Options Database Key:
2072 .  -tao_type <type> - Sets the method; use -help for a list
2073    of available methods (for instance, "-tao_type lmvm" or "-tao_type tron")
2074 
2075    Available methods include:
2076 +    nls - Newton's method with line search for unconstrained minimization
2077 .    ntr - Newton's method with trust region for unconstrained minimization
2078 .    ntl - Newton's method with trust region, line search for unconstrained minimization
2079 .    lmvm - Limited memory variable metric method for unconstrained minimization
2080 .    cg - Nonlinear conjugate gradient method for unconstrained minimization
2081 .    nm - Nelder-Mead algorithm for derivate-free unconstrained minimization
2082 .    tron - Newton Trust Region method for bound constrained minimization
2083 .    gpcg - Newton Trust Region method for quadratic bound constrained minimization
2084 .    blmvm - Limited memory variable metric method for bound constrained minimization
2085 -    pounders - Model-based algorithm pounder extended for nonlinear least squares
2086 
2087   Level: intermediate
2088 
2089 .seealso: TaoCreate(), TaoGetType(), TaoType
2090 
2091 @*/
2092 PetscErrorCode TaoSetType(Tao tao, TaoType type)
2093 {
2094   PetscErrorCode (*create_xxx)(Tao);
2095   PetscBool      issame;
2096 
2097   PetscFunctionBegin;
2098   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2099 
2100   PetscCall(PetscObjectTypeCompare((PetscObject)tao,type,&issame));
2101   if (issame) PetscFunctionReturn(0);
2102 
2103   PetscCall(PetscFunctionListFind(TaoList, type, (void(**)(void))&create_xxx));
2104   PetscCheck(create_xxx,PetscObjectComm((PetscObject)tao),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested Tao type %s",type);
2105 
2106   /* Destroy the existing solver information */
2107   if (tao->ops->destroy) {
2108     PetscCall((*tao->ops->destroy)(tao));
2109   }
2110   PetscCall(KSPDestroy(&tao->ksp));
2111   PetscCall(TaoLineSearchDestroy(&tao->linesearch));
2112   PetscCall(VecDestroy(&tao->gradient));
2113   PetscCall(VecDestroy(&tao->stepdirection));
2114 
2115   tao->ops->setup = NULL;
2116   tao->ops->solve = NULL;
2117   tao->ops->view  = NULL;
2118   tao->ops->setfromoptions = NULL;
2119   tao->ops->destroy = NULL;
2120 
2121   tao->setupcalled = PETSC_FALSE;
2122 
2123   PetscCall((*create_xxx)(tao));
2124   PetscCall(PetscObjectChangeTypeName((PetscObject)tao,type));
2125   PetscFunctionReturn(0);
2126 }
2127 
2128 /*MC
2129    TaoRegister - Adds a method to the TAO package for unconstrained minimization.
2130 
2131    Synopsis:
2132    TaoRegister(char *name_solver,char *path,char *name_Create,PetscErrorCode (*routine_Create)(Tao))
2133 
2134    Not collective
2135 
2136    Input Parameters:
2137 +  sname - name of a new user-defined solver
2138 -  func - routine to Create method context
2139 
2140    Notes:
2141    TaoRegister() may be called multiple times to add several user-defined solvers.
2142 
2143    Sample usage:
2144 .vb
2145    TaoRegister("my_solver",MySolverCreate);
2146 .ve
2147 
2148    Then, your solver can be chosen with the procedural interface via
2149 $     TaoSetType(tao,"my_solver")
2150    or at runtime via the option
2151 $     -tao_type my_solver
2152 
2153    Level: advanced
2154 
2155 .seealso: TaoRegisterAll(), TaoRegisterDestroy()
2156 M*/
2157 PetscErrorCode TaoRegister(const char sname[], PetscErrorCode (*func)(Tao))
2158 {
2159   PetscFunctionBegin;
2160   PetscCall(TaoInitializePackage());
2161   PetscCall(PetscFunctionListAdd(&TaoList,sname,(void (*)(void))func));
2162   PetscFunctionReturn(0);
2163 }
2164 
2165 /*@C
2166    TaoRegisterDestroy - Frees the list of minimization solvers that were
2167    registered by TaoRegisterDynamic().
2168 
2169    Not Collective
2170 
2171    Level: advanced
2172 
2173 .seealso: TaoRegisterAll(), TaoRegister()
2174 @*/
2175 PetscErrorCode TaoRegisterDestroy(void)
2176 {
2177   PetscFunctionBegin;
2178   PetscCall(PetscFunctionListDestroy(&TaoList));
2179   TaoRegisterAllCalled = PETSC_FALSE;
2180   PetscFunctionReturn(0);
2181 }
2182 
2183 /*@
2184    TaoGetIterationNumber - Gets the number of Tao iterations completed
2185    at this time.
2186 
2187    Not Collective
2188 
2189    Input Parameter:
2190 .  tao - Tao context
2191 
2192    Output Parameter:
2193 .  iter - iteration number
2194 
2195    Notes:
2196    For example, during the computation of iteration 2 this would return 1.
2197 
2198    Level: intermediate
2199 
2200 .seealso:   TaoGetLinearSolveIterations(), TaoGetResidualNorm(), TaoGetObjective()
2201 @*/
2202 PetscErrorCode  TaoGetIterationNumber(Tao tao,PetscInt *iter)
2203 {
2204   PetscFunctionBegin;
2205   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2206   PetscValidIntPointer(iter,2);
2207   *iter = tao->niter;
2208   PetscFunctionReturn(0);
2209 }
2210 
2211 /*@
2212    TaoGetResidualNorm - Gets the current value of the norm of the residual
2213    at this time.
2214 
2215    Not Collective
2216 
2217    Input Parameter:
2218 .  tao - Tao context
2219 
2220    Output Parameter:
2221 .  value - the current value
2222 
2223    Level: intermediate
2224 
2225    Developer Note: This is the 2-norm of the residual, we cannot use TaoGetGradientNorm() because that has
2226                    a different meaning. For some reason Tao sometimes calls the gradient the residual.
2227 
2228 .seealso:   TaoGetLinearSolveIterations(), TaoGetIterationNumber(), TaoGetObjective()
2229 @*/
2230 PetscErrorCode TaoGetResidualNorm(Tao tao,PetscReal *value)
2231 {
2232   PetscFunctionBegin;
2233   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2234   PetscValidRealPointer(value,2);
2235   *value = tao->residual;
2236   PetscFunctionReturn(0);
2237 }
2238 
2239 /*@
2240    TaoSetIterationNumber - Sets the current iteration number.
2241 
2242    Logically Collective on Tao
2243 
2244    Input Parameters:
2245 +  tao - Tao context
2246 -  iter - iteration number
2247 
2248    Level: developer
2249 
2250 .seealso:   TaoGetLinearSolveIterations()
2251 @*/
2252 PetscErrorCode  TaoSetIterationNumber(Tao tao,PetscInt iter)
2253 {
2254   PetscFunctionBegin;
2255   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2256   PetscValidLogicalCollectiveInt(tao,iter,2);
2257   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)tao));
2258   tao->niter = iter;
2259   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)tao));
2260   PetscFunctionReturn(0);
2261 }
2262 
2263 /*@
2264    TaoGetTotalIterationNumber - Gets the total number of Tao iterations
2265    completed. This number keeps accumulating if multiple solves
2266    are called with the Tao object.
2267 
2268    Not Collective
2269 
2270    Input Parameter:
2271 .  tao - Tao context
2272 
2273    Output Parameter:
2274 .  iter - iteration number
2275 
2276    Notes:
2277    The total iteration count is updated after each solve, if there is a current
2278    TaoSolve() in progress then those iterations are not yet counted.
2279 
2280    Level: intermediate
2281 
2282 .seealso:   TaoGetLinearSolveIterations()
2283 @*/
2284 PetscErrorCode  TaoGetTotalIterationNumber(Tao tao,PetscInt *iter)
2285 {
2286   PetscFunctionBegin;
2287   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2288   PetscValidIntPointer(iter,2);
2289   *iter = tao->ntotalits;
2290   PetscFunctionReturn(0);
2291 }
2292 
2293 /*@
2294    TaoSetTotalIterationNumber - Sets the current total iteration number.
2295 
2296    Logically Collective on Tao
2297 
2298    Input Parameters:
2299 +  tao - Tao context
2300 -  iter - iteration number
2301 
2302    Level: developer
2303 
2304 .seealso:   TaoGetLinearSolveIterations()
2305 @*/
2306 PetscErrorCode  TaoSetTotalIterationNumber(Tao tao,PetscInt iter)
2307 {
2308   PetscFunctionBegin;
2309   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2310   PetscValidLogicalCollectiveInt(tao,iter,2);
2311   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)tao));
2312   tao->ntotalits = iter;
2313   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)tao));
2314   PetscFunctionReturn(0);
2315 }
2316 
2317 /*@
2318   TaoSetConvergedReason - Sets the termination flag on a Tao object
2319 
2320   Logically Collective on Tao
2321 
2322   Input Parameters:
2323 + tao - the Tao context
2324 - reason - one of
2325 $     TAO_CONVERGED_ATOL (2),
2326 $     TAO_CONVERGED_RTOL (3),
2327 $     TAO_CONVERGED_STEPTOL (4),
2328 $     TAO_CONVERGED_MINF (5),
2329 $     TAO_CONVERGED_USER (6),
2330 $     TAO_DIVERGED_MAXITS (-2),
2331 $     TAO_DIVERGED_NAN (-4),
2332 $     TAO_DIVERGED_MAXFCN (-5),
2333 $     TAO_DIVERGED_LS_FAILURE (-6),
2334 $     TAO_DIVERGED_TR_REDUCTION (-7),
2335 $     TAO_DIVERGED_USER (-8),
2336 $     TAO_CONTINUE_ITERATING (0)
2337 
2338    Level: intermediate
2339 
2340 @*/
2341 PetscErrorCode TaoSetConvergedReason(Tao tao, TaoConvergedReason reason)
2342 {
2343   PetscFunctionBegin;
2344   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2345   PetscValidLogicalCollectiveEnum(tao,reason,2);
2346   tao->reason = reason;
2347   PetscFunctionReturn(0);
2348 }
2349 
2350 /*@
2351    TaoGetConvergedReason - Gets the reason the Tao iteration was stopped.
2352 
2353    Not Collective
2354 
2355    Input Parameter:
2356 .  tao - the Tao solver context
2357 
2358    Output Parameter:
2359 .  reason - one of
2360 $  TAO_CONVERGED_GATOL (3)           ||g(X)|| < gatol
2361 $  TAO_CONVERGED_GRTOL (4)           ||g(X)|| / f(X)  < grtol
2362 $  TAO_CONVERGED_GTTOL (5)           ||g(X)|| / ||g(X0)|| < gttol
2363 $  TAO_CONVERGED_STEPTOL (6)         step size small
2364 $  TAO_CONVERGED_MINF (7)            F < F_min
2365 $  TAO_CONVERGED_USER (8)            User defined
2366 $  TAO_DIVERGED_MAXITS (-2)          its > maxits
2367 $  TAO_DIVERGED_NAN (-4)             Numerical problems
2368 $  TAO_DIVERGED_MAXFCN (-5)          fevals > max_funcsals
2369 $  TAO_DIVERGED_LS_FAILURE (-6)      line search failure
2370 $  TAO_DIVERGED_TR_REDUCTION (-7)    trust region failure
2371 $  TAO_DIVERGED_USER(-8)             (user defined)
2372 $  TAO_CONTINUE_ITERATING (0)
2373 
2374    where
2375 +  X - current solution
2376 .  X0 - initial guess
2377 .  f(X) - current function value
2378 .  f(X*) - true solution (estimated)
2379 .  g(X) - current gradient
2380 .  its - current iterate number
2381 .  maxits - maximum number of iterates
2382 .  fevals - number of function evaluations
2383 -  max_funcsals - maximum number of function evaluations
2384 
2385    Level: intermediate
2386 
2387 .seealso: TaoSetConvergenceTest(), TaoSetTolerances()
2388 
2389 @*/
2390 PetscErrorCode TaoGetConvergedReason(Tao tao, TaoConvergedReason *reason)
2391 {
2392   PetscFunctionBegin;
2393   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2394   PetscValidPointer(reason,2);
2395   *reason = tao->reason;
2396   PetscFunctionReturn(0);
2397 }
2398 
2399 /*@
2400    TaoGetSolutionStatus - Get the current iterate, objective value,
2401    residual, infeasibility, and termination
2402 
2403    Not Collective
2404 
2405    Input Parameter:
2406 .  tao - the Tao context
2407 
2408    Output Parameters:
2409 +  iterate - the current iterate number (>=0)
2410 .  f - the current function value
2411 .  gnorm - the square of the gradient norm, duality gap, or other measure indicating distance from optimality.
2412 .  cnorm - the infeasibility of the current solution with regard to the constraints.
2413 .  xdiff - the step length or trust region radius of the most recent iterate.
2414 -  reason - The termination reason, which can equal TAO_CONTINUE_ITERATING
2415 
2416    Level: intermediate
2417 
2418    Note:
2419    TAO returns the values set by the solvers in the routine TaoMonitor().
2420 
2421    Note:
2422    If any of the output arguments are set to NULL, no corresponding value will be returned.
2423 
2424 .seealso: TaoMonitor(), TaoGetConvergedReason()
2425 @*/
2426 PetscErrorCode TaoGetSolutionStatus(Tao tao, PetscInt *its, PetscReal *f, PetscReal *gnorm, PetscReal *cnorm, PetscReal *xdiff, TaoConvergedReason *reason)
2427 {
2428   PetscFunctionBegin;
2429   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2430   if (its) *its = tao->niter;
2431   if (f) *f = tao->fc;
2432   if (gnorm) *gnorm = tao->residual;
2433   if (cnorm) *cnorm = tao->cnorm;
2434   if (reason) *reason = tao->reason;
2435   if (xdiff) *xdiff = tao->step;
2436   PetscFunctionReturn(0);
2437 }
2438 
2439 /*@C
2440    TaoGetType - Gets the current Tao algorithm.
2441 
2442    Not Collective
2443 
2444    Input Parameter:
2445 .  tao - the Tao solver context
2446 
2447    Output Parameter:
2448 .  type - Tao method
2449 
2450    Level: intermediate
2451 
2452 @*/
2453 PetscErrorCode TaoGetType(Tao tao,TaoType *type)
2454 {
2455   PetscFunctionBegin;
2456   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2457   PetscValidPointer(type,2);
2458   *type = ((PetscObject)tao)->type_name;
2459   PetscFunctionReturn(0);
2460 }
2461 
2462 /*@C
2463   TaoMonitor - Monitor the solver and the current solution.  This
2464   routine will record the iteration number and residual statistics,
2465   call any monitors specified by the user, and calls the convergence-check routine.
2466 
2467    Input Parameters:
2468 +  tao - the Tao context
2469 .  its - the current iterate number (>=0)
2470 .  f - the current objective function value
2471 .  res - the gradient norm, square root of the duality gap, or other measure indicating distince from optimality.  This measure will be recorded and
2472           used for some termination tests.
2473 .  cnorm - the infeasibility of the current solution with regard to the constraints.
2474 -  steplength - multiple of the step direction added to the previous iterate.
2475 
2476    Output Parameters:
2477 .  reason - The termination reason, which can equal TAO_CONTINUE_ITERATING
2478 
2479    Options Database Key:
2480 .  -tao_monitor - Use the default monitor, which prints statistics to standard output
2481 
2482 .seealso TaoGetConvergedReason(), TaoMonitorDefault(), TaoSetMonitor()
2483 
2484    Level: developer
2485 
2486 @*/
2487 PetscErrorCode TaoMonitor(Tao tao, PetscInt its, PetscReal f, PetscReal res, PetscReal cnorm, PetscReal steplength)
2488 {
2489   PetscInt       i;
2490 
2491   PetscFunctionBegin;
2492   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2493   tao->fc = f;
2494   tao->residual = res;
2495   tao->cnorm = cnorm;
2496   tao->step = steplength;
2497   if (!its) {
2498     tao->cnorm0 = cnorm;
2499     tao->gnorm0 = res;
2500   }
2501   PetscCheck(!PetscIsInfOrNanReal(f) && !PetscIsInfOrNanReal(res),PetscObjectComm((PetscObject)tao),PETSC_ERR_USER, "User provided compute function generated Inf or NaN");
2502   for (i=0;i<tao->numbermonitors;i++) {
2503     PetscCall((*tao->monitor[i])(tao,tao->monitorcontext[i]));
2504   }
2505   PetscFunctionReturn(0);
2506 }
2507 
2508 /*@
2509    TaoSetConvergenceHistory - Sets the array used to hold the convergence history.
2510 
2511    Logically Collective on Tao
2512 
2513    Input Parameters:
2514 +  tao - the Tao solver context
2515 .  obj   - array to hold objective value history
2516 .  resid - array to hold residual history
2517 .  cnorm - array to hold constraint violation history
2518 .  lits - integer array holds the number of linear iterations for each Tao iteration
2519 .  na  - size of obj, resid, and cnorm
2520 -  reset - PetscTrue indicates each new minimization resets the history counter to zero,
2521            else it continues storing new values for new minimizations after the old ones
2522 
2523    Notes:
2524    If set, TAO will fill the given arrays with the indicated
2525    information at each iteration.  If 'obj','resid','cnorm','lits' are
2526    *all* NULL then space (using size na, or 1000 if na is PETSC_DECIDE or
2527    PETSC_DEFAULT) is allocated for the history.
2528    If not all are NULL, then only the non-NULL information categories
2529    will be stored, the others will be ignored.
2530 
2531    Any convergence information after iteration number 'na' will not be stored.
2532 
2533    This routine is useful, e.g., when running a code for purposes
2534    of accurate performance monitoring, when no I/O should be done
2535    during the section of code that is being timed.
2536 
2537    Level: intermediate
2538 
2539 .seealso: TaoGetConvergenceHistory()
2540 
2541 @*/
2542 PetscErrorCode TaoSetConvergenceHistory(Tao tao, PetscReal obj[], PetscReal resid[], PetscReal cnorm[], PetscInt lits[], PetscInt na,PetscBool reset)
2543 {
2544   PetscFunctionBegin;
2545   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2546   if (obj) PetscValidRealPointer(obj,2);
2547   if (resid) PetscValidRealPointer(resid,3);
2548   if (cnorm) PetscValidRealPointer(cnorm,4);
2549   if (lits) PetscValidIntPointer(lits,5);
2550 
2551   if (na == PETSC_DECIDE || na == PETSC_DEFAULT) na = 1000;
2552   if (!obj && !resid && !cnorm && !lits) {
2553     PetscCall(PetscCalloc4(na,&obj,na,&resid,na,&cnorm,na,&lits));
2554     tao->hist_malloc = PETSC_TRUE;
2555   }
2556 
2557   tao->hist_obj = obj;
2558   tao->hist_resid = resid;
2559   tao->hist_cnorm = cnorm;
2560   tao->hist_lits = lits;
2561   tao->hist_max   = na;
2562   tao->hist_reset = reset;
2563   tao->hist_len = 0;
2564   PetscFunctionReturn(0);
2565 }
2566 
2567 /*@C
2568    TaoGetConvergenceHistory - Gets the arrays used to hold the convergence history.
2569 
2570    Collective on Tao
2571 
2572    Input Parameter:
2573 .  tao - the Tao context
2574 
2575    Output Parameters:
2576 +  obj   - array used to hold objective value history
2577 .  resid - array used to hold residual history
2578 .  cnorm - array used to hold constraint violation history
2579 .  lits  - integer array used to hold linear solver iteration count
2580 -  nhist  - size of obj, resid, cnorm, and lits
2581 
2582    Notes:
2583     This routine must be preceded by calls to TaoSetConvergenceHistory()
2584     and TaoSolve(), otherwise it returns useless information.
2585 
2586     The calling sequence for this routine in Fortran is
2587 $   call TaoGetConvergenceHistory(Tao tao, PetscInt nhist, PetscErrorCode ierr)
2588 
2589    This routine is useful, e.g., when running a code for purposes
2590    of accurate performance monitoring, when no I/O should be done
2591    during the section of code that is being timed.
2592 
2593    Level: advanced
2594 
2595 .seealso: TaoSetConvergenceHistory()
2596 
2597 @*/
2598 PetscErrorCode TaoGetConvergenceHistory(Tao tao, PetscReal **obj, PetscReal **resid, PetscReal **cnorm, PetscInt **lits, PetscInt *nhist)
2599 {
2600   PetscFunctionBegin;
2601   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2602   if (obj)   *obj   = tao->hist_obj;
2603   if (cnorm) *cnorm = tao->hist_cnorm;
2604   if (resid) *resid = tao->hist_resid;
2605   if (nhist) *nhist = tao->hist_len;
2606   PetscFunctionReturn(0);
2607 }
2608 
2609 /*@
2610    TaoSetApplicationContext - Sets the optional user-defined context for
2611    a solver.
2612 
2613    Logically Collective on Tao
2614 
2615    Input Parameters:
2616 +  tao  - the Tao context
2617 -  usrP - optional user context
2618 
2619    Level: intermediate
2620 
2621 .seealso: TaoGetApplicationContext(), TaoSetApplicationContext()
2622 @*/
2623 PetscErrorCode  TaoSetApplicationContext(Tao tao,void *usrP)
2624 {
2625   PetscFunctionBegin;
2626   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2627   tao->user = usrP;
2628   PetscFunctionReturn(0);
2629 }
2630 
2631 /*@
2632    TaoGetApplicationContext - Gets the user-defined context for a
2633    TAO solvers.
2634 
2635    Not Collective
2636 
2637    Input Parameter:
2638 .  tao  - Tao context
2639 
2640    Output Parameter:
2641 .  usrP - user context
2642 
2643    Level: intermediate
2644 
2645 .seealso: TaoSetApplicationContext()
2646 @*/
2647 PetscErrorCode  TaoGetApplicationContext(Tao tao,void *usrP)
2648 {
2649   PetscFunctionBegin;
2650   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2651   PetscValidPointer(usrP,2);
2652   *(void**)usrP = tao->user;
2653   PetscFunctionReturn(0);
2654 }
2655 
2656 /*@
2657    TaoSetGradientNorm - Sets the matrix used to define the inner product that measures the size of the gradient.
2658 
2659    Collective on tao
2660 
2661    Input Parameters:
2662 +  tao  - the Tao context
2663 -  M    - gradient norm
2664 
2665    Level: beginner
2666 
2667 .seealso: TaoGetGradientNorm(), TaoGradientNorm()
2668 @*/
2669 PetscErrorCode  TaoSetGradientNorm(Tao tao, Mat M)
2670 {
2671   PetscFunctionBegin;
2672   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2673   PetscValidHeaderSpecific(M,MAT_CLASSID,2);
2674   PetscCall(PetscObjectReference((PetscObject)M));
2675   PetscCall(MatDestroy(&tao->gradient_norm));
2676   PetscCall(VecDestroy(&tao->gradient_norm_tmp));
2677   tao->gradient_norm = M;
2678   PetscCall(MatCreateVecs(M, NULL, &tao->gradient_norm_tmp));
2679   PetscFunctionReturn(0);
2680 }
2681 
2682 /*@
2683    TaoGetGradientNorm - Returns the matrix used to define the inner product for measuring the size of the gradient.
2684 
2685    Not Collective
2686 
2687    Input Parameter:
2688 .  tao  - Tao context
2689 
2690    Output Parameter:
2691 .  M - gradient norm
2692 
2693    Level: beginner
2694 
2695 .seealso: TaoSetGradientNorm(), TaoGradientNorm()
2696 @*/
2697 PetscErrorCode  TaoGetGradientNorm(Tao tao, Mat *M)
2698 {
2699   PetscFunctionBegin;
2700   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2701   PetscValidPointer(M,2);
2702   *M = tao->gradient_norm;
2703   PetscFunctionReturn(0);
2704 }
2705 
2706 /*@C
2707    TaoGradientNorm - Compute the norm with respect to the inner product the user has set.
2708 
2709    Collective on tao
2710 
2711    Input Parameters:
2712 +  tao      - the Tao context
2713 .  gradient - the gradient to be computed
2714 -  norm     - the norm type
2715 
2716    Output Parameter:
2717 .  gnorm    - the gradient norm
2718 
2719    Level: developer
2720 
2721 .seealso: TaoSetGradientNorm(), TaoGetGradientNorm()
2722 @*/
2723 PetscErrorCode  TaoGradientNorm(Tao tao, Vec gradient, NormType type, PetscReal *gnorm)
2724 {
2725   PetscFunctionBegin;
2726   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2727   PetscValidHeaderSpecific(gradient,VEC_CLASSID,2);
2728   PetscValidLogicalCollectiveEnum(tao,type,3);
2729   PetscValidRealPointer(gnorm,4);
2730   if (tao->gradient_norm) {
2731     PetscScalar gnorms;
2732 
2733     PetscCheck(type == NORM_2,PetscObjectComm((PetscObject)gradient), PETSC_ERR_ARG_WRONG, "Norm type must be NORM_2 if an inner product for the gradient norm is set.");
2734     PetscCall(MatMult(tao->gradient_norm, gradient, tao->gradient_norm_tmp));
2735     PetscCall(VecDot(gradient, tao->gradient_norm_tmp, &gnorms));
2736     *gnorm = PetscRealPart(PetscSqrtScalar(gnorms));
2737   } else {
2738     PetscCall(VecNorm(gradient, type, gnorm));
2739   }
2740   PetscFunctionReturn(0);
2741 }
2742 
2743 /*@C
2744    TaoMonitorDrawCtxCreate - Creates the monitor context for TaoMonitorDrawCtx
2745 
2746    Collective on Tao
2747 
2748    Output Patameter:
2749 .    ctx - the monitor context
2750 
2751    Options Database:
2752 .   -tao_draw_solution_initial - show initial guess as well as current solution
2753 
2754    Level: intermediate
2755 
2756 .seealso: TaoMonitorSet(), TaoMonitorDefault(), VecView(), TaoMonitorDrawCtx()
2757 @*/
2758 PetscErrorCode  TaoMonitorDrawCtxCreate(MPI_Comm comm,const char host[],const char label[],int x,int y,int m,int n,PetscInt howoften,TaoMonitorDrawCtx *ctx)
2759 {
2760   PetscFunctionBegin;
2761   PetscCall(PetscNew(ctx));
2762   PetscCall(PetscViewerDrawOpen(comm,host,label,x,y,m,n,&(*ctx)->viewer));
2763   PetscCall(PetscViewerSetFromOptions((*ctx)->viewer));
2764   (*ctx)->howoften = howoften;
2765   PetscFunctionReturn(0);
2766 }
2767 
2768 /*@C
2769    TaoMonitorDrawCtxDestroy - Destroys the monitor context for TaoMonitorDrawSolution()
2770 
2771    Collective on Tao
2772 
2773    Input Parameters:
2774 .    ctx - the monitor context
2775 
2776    Level: intermediate
2777 
2778 .seealso: TaoMonitorSet(), TaoMonitorDefault(), VecView(), TaoMonitorDrawSolution()
2779 @*/
2780 PetscErrorCode  TaoMonitorDrawCtxDestroy(TaoMonitorDrawCtx *ictx)
2781 {
2782   PetscFunctionBegin;
2783   PetscCall(PetscViewerDestroy(&(*ictx)->viewer));
2784   PetscCall(PetscFree(*ictx));
2785   PetscFunctionReturn(0);
2786 }
2787