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