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