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