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