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