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