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