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