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