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