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