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