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