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