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