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