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