xref: /petsc/src/tao/interface/taosolver.c (revision e91ef88e73a631b88b8031ec3d9dd16cd2dcf0e1)
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 = (PetscViewer)ctx;
1541 
1542   PetscFunctionBegin;
1543   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
1544   its=tao->niter;
1545   fct=tao->fc;
1546   gnorm=tao->residual;
1547   ierr=PetscViewerASCIIPrintf(viewer,"iter = %3D,",its);CHKERRQ(ierr);
1548   ierr=PetscViewerASCIIPrintf(viewer," Function value: %g,",(double)fct);CHKERRQ(ierr);
1549   if (gnorm >= PETSC_INFINITY) {
1550     ierr=PetscViewerASCIIPrintf(viewer,"  Residual: Inf \n");CHKERRQ(ierr);
1551   } else {
1552     ierr=PetscViewerASCIIPrintf(viewer,"  Residual: %g \n",(double)gnorm);CHKERRQ(ierr);
1553   }
1554   PetscFunctionReturn(0);
1555 }
1556 
1557 #undef __FUNCT__
1558 #define __FUNCT__ "TaoDefaultSMonitor"
1559 /*@
1560    TaoDefaultSMonitor - Default routine for monitoring progress of the
1561    solver. Same as TaoDefaultMonitor() except
1562    it prints fewer digits of the residual as the residual gets smaller.
1563    This is because the later digits are meaningless and are often
1564    different on different machines; by using this routine different
1565    machines will usually generate the same output. It can be turned on
1566    by using the -tao_smonitor option
1567 
1568    Collective on Tao
1569 
1570    Input Parameters:
1571 +  tao - the Tao context
1572 -  ctx - PetscViewer context of type ASCII
1573 
1574    Options Database Keys:
1575 .  -tao_smonitor
1576 
1577    Level: advanced
1578 
1579 .seealso: TaoDefaultMonitor(), TaoSetMonitor()
1580 @*/
1581 PetscErrorCode TaoDefaultSMonitor(Tao tao, void *ctx)
1582 {
1583   PetscErrorCode ierr;
1584   PetscInt       its;
1585   PetscReal      fct,gnorm;
1586   PetscViewer    viewer = (PetscViewer)ctx;
1587 
1588   PetscFunctionBegin;
1589   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
1590   its=tao->niter;
1591   fct=tao->fc;
1592   gnorm=tao->residual;
1593   ierr=PetscViewerASCIIPrintf(viewer,"iter = %3D,",its);CHKERRQ(ierr);
1594   ierr=PetscViewerASCIIPrintf(viewer," Function value %g,",(double)fct);CHKERRQ(ierr);
1595   if (gnorm >= PETSC_INFINITY) {
1596     ierr=PetscViewerASCIIPrintf(viewer," Residual: Inf \n");CHKERRQ(ierr);
1597   } else if (gnorm > 1.e-6) {
1598     ierr=PetscViewerASCIIPrintf(viewer," Residual: %g \n",(double)gnorm);CHKERRQ(ierr);
1599   } else if (gnorm > 1.e-11) {
1600     ierr=PetscViewerASCIIPrintf(viewer," Residual: < 1.0e-6 \n");CHKERRQ(ierr);
1601   } else {
1602     ierr=PetscViewerASCIIPrintf(viewer," Residual: < 1.0e-11 \n");CHKERRQ(ierr);
1603   }
1604   PetscFunctionReturn(0);
1605 }
1606 
1607 #undef __FUNCT__
1608 #define __FUNCT__ "TaoDefaultCMonitor"
1609 /*@
1610    TaoDefaultCMonitor - same as TaoDefaultMonitor() except
1611    it prints the norm of the constraints function. It can be turned on
1612    from the command line using the -tao_cmonitor option
1613 
1614    Collective on Tao
1615 
1616    Input Parameters:
1617 +  tao - the Tao context
1618 -  ctx - PetscViewer context or NULL
1619 
1620    Options Database Keys:
1621 .  -tao_cmonitor
1622 
1623    Level: advanced
1624 
1625 .seealso: TaoDefaultMonitor(), TaoSetMonitor()
1626 @*/
1627 PetscErrorCode TaoDefaultCMonitor(Tao tao, void *ctx)
1628 {
1629   PetscErrorCode ierr;
1630   PetscInt       its;
1631   PetscReal      fct,gnorm;
1632   PetscViewer    viewer = (PetscViewer)ctx;
1633 
1634   PetscFunctionBegin;
1635   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
1636   its=tao->niter;
1637   fct=tao->fc;
1638   gnorm=tao->residual;
1639   ierr=PetscViewerASCIIPrintf(viewer,"iter = %D,",its);CHKERRQ(ierr);
1640   ierr=PetscViewerASCIIPrintf(viewer," Function value: %g,",(double)fct);CHKERRQ(ierr);
1641   ierr=PetscViewerASCIIPrintf(viewer,"  Residual: %g ",(double)gnorm);CHKERRQ(ierr);
1642   ierr = PetscViewerASCIIPrintf(viewer,"  Constraint: %g \n",(double)tao->cnorm);CHKERRQ(ierr);
1643   PetscFunctionReturn(0);
1644 }
1645 
1646 #undef __FUNCT__
1647 #define __FUNCT__ "TaoSolutionMonitor"
1648 /*@C
1649    TaoSolutionMonitor - Views the solution at each iteration
1650    It can be turned on from the command line using the
1651    -tao_view_solution option
1652 
1653    Collective on Tao
1654 
1655    Input Parameters:
1656 +  tao - the Tao context
1657 -  ctx - PetscViewer context or NULL
1658 
1659    Options Database Keys:
1660 .  -tao_view_solution
1661 
1662    Level: advanced
1663 
1664 .seealso: TaoDefaultSMonitor(), TaoSetMonitor()
1665 @*/
1666 PetscErrorCode TaoSolutionMonitor(Tao tao, void *ctx)
1667 {
1668   PetscErrorCode ierr;
1669   PetscViewer    viewer  = (PetscViewer)ctx;;
1670 
1671   PetscFunctionBegin;
1672   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
1673   ierr = VecView(tao->solution, viewer);CHKERRQ(ierr);
1674   PetscFunctionReturn(0);
1675 }
1676 
1677 #undef __FUNCT__
1678 #define __FUNCT__ "TaoGradientMonitor"
1679 /*@C
1680    TaoGradientMonitor - Views the gradient at each iteration
1681    It can be turned on from the command line using the
1682    -tao_view_gradient option
1683 
1684    Collective on Tao
1685 
1686    Input Parameters:
1687 +  tao - the Tao context
1688 -  ctx - PetscViewer context or NULL
1689 
1690    Options Database Keys:
1691 .  -tao_view_gradient
1692 
1693    Level: advanced
1694 
1695 .seealso: TaoDefaultSMonitor(), TaoSetMonitor()
1696 @*/
1697 PetscErrorCode TaoGradientMonitor(Tao tao, void *ctx)
1698 {
1699   PetscErrorCode ierr;
1700   PetscViewer    viewer = (PetscViewer)ctx;
1701 
1702   PetscFunctionBegin;
1703   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
1704   ierr = VecView(tao->gradient, viewer);CHKERRQ(ierr);
1705   PetscFunctionReturn(0);
1706 }
1707 
1708 #undef __FUNCT__
1709 #define __FUNCT__ "TaoStepDirectionMonitor"
1710 /*@C
1711    TaoStepDirectionMonitor - Views the gradient at each iteration
1712    It can be turned on from the command line using the
1713    -tao_view_gradient option
1714 
1715    Collective on Tao
1716 
1717    Input Parameters:
1718 +  tao - the Tao context
1719 -  ctx - PetscViewer context or NULL
1720 
1721    Options Database Keys:
1722 .  -tao_view_gradient
1723 
1724    Level: advanced
1725 
1726 .seealso: TaoDefaultSMonitor(), TaoSetMonitor()
1727 @*/
1728 PetscErrorCode TaoStepDirectionMonitor(Tao tao, void *ctx)
1729 {
1730   PetscErrorCode ierr;
1731   PetscViewer    viewer = (PetscViewer)ctx;
1732 
1733   PetscFunctionBegin;
1734   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
1735   ierr = VecView(tao->stepdirection, viewer);CHKERRQ(ierr);
1736   PetscFunctionReturn(0);
1737 }
1738 
1739 #undef __FUNCT__
1740 #define __FUNCT__ "TaoDrawSolutionMonitor"
1741 /*@C
1742    TaoDrawSolutionMonitor - Plots the solution at each iteration
1743    It can be turned on from the command line using the
1744    -tao_draw_solution option
1745 
1746    Collective on Tao
1747 
1748    Input Parameters:
1749 +  tao - the Tao context
1750 -  ctx - PetscViewer context
1751 
1752    Options Database Keys:
1753 .  -tao_draw_solution
1754 
1755    Level: advanced
1756 
1757 .seealso: TaoSolutionMonitor(), TaoSetMonitor(), TaoDrawGradientMonitor
1758 @*/
1759 PetscErrorCode TaoDrawSolutionMonitor(Tao tao, void *ctx)
1760 {
1761   PetscErrorCode ierr;
1762   PetscViewer    viewer = (PetscViewer) ctx;
1763 
1764   PetscFunctionBegin;
1765   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
1766   ierr = VecView(tao->solution, viewer);CHKERRQ(ierr);
1767   PetscFunctionReturn(0);
1768 }
1769 
1770 #undef __FUNCT__
1771 #define __FUNCT__ "TaoDrawGradientMonitor"
1772 /*@C
1773    TaoDrawGradientMonitor - Plots the gradient at each iteration
1774    It can be turned on from the command line using the
1775    -tao_draw_gradient option
1776 
1777    Collective on Tao
1778 
1779    Input Parameters:
1780 +  tao - the Tao context
1781 -  ctx - PetscViewer context
1782 
1783    Options Database Keys:
1784 .  -tao_draw_gradient
1785 
1786    Level: advanced
1787 
1788 .seealso: TaoGradientMonitor(), TaoSetMonitor(), TaoDrawSolutionMonitor
1789 @*/
1790 PetscErrorCode TaoDrawGradientMonitor(Tao tao, void *ctx)
1791 {
1792   PetscErrorCode ierr;
1793   PetscViewer    viewer = (PetscViewer)ctx;
1794 
1795   PetscFunctionBegin;
1796   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
1797   ierr = VecView(tao->gradient, viewer);CHKERRQ(ierr);
1798   PetscFunctionReturn(0);
1799 }
1800 
1801 #undef __FUNCT__
1802 #define __FUNCT__ "TaoDrawStepMonitor"
1803 /*@C
1804    TaoDrawStepMonitor - Plots the step direction at each iteration
1805    It can be turned on from the command line using the
1806    -tao_draw_step option
1807 
1808    Collective on Tao
1809 
1810    Input Parameters:
1811 +  tao - the Tao context
1812 -  ctx - PetscViewer context
1813 
1814    Options Database Keys:
1815 .  -tao_draw_step
1816 
1817    Level: advanced
1818 
1819 .seealso: TaoSetMonitor(), TaoDrawSolutionMonitor
1820 @*/
1821 PetscErrorCode TaoDrawStepMonitor(Tao tao, void *ctx)
1822 {
1823   PetscErrorCode ierr;
1824   PetscViewer    viewer = (PetscViewer)(ctx);
1825 
1826   PetscFunctionBegin;
1827   ierr = VecView(tao->stepdirection, viewer);CHKERRQ(ierr);
1828   PetscFunctionReturn(0);
1829 }
1830 
1831 #undef __FUNCT__
1832 #define __FUNCT__ "TaoSeparableObjectiveMonitor"
1833 /*@C
1834    TaoSeparableObjectiveMonitor - Views the separable objective function at each iteration
1835    It can be turned on from the command line using the
1836    -tao_view_separableobjective option
1837 
1838    Collective on Tao
1839 
1840    Input Parameters:
1841 +  tao - the Tao context
1842 -  ctx - PetscViewer context or NULL
1843 
1844    Options Database Keys:
1845 .  -tao_view_separableobjective
1846 
1847    Level: advanced
1848 
1849 .seealso: TaoDefaultSMonitor(), TaoSetMonitor()
1850 @*/
1851 PetscErrorCode TaoSeparableObjectiveMonitor(Tao tao, void *ctx)
1852 {
1853   PetscErrorCode ierr;
1854   PetscViewer    viewer  = (PetscViewer)ctx;
1855 
1856   PetscFunctionBegin;
1857   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
1858   ierr = VecView(tao->sep_objective,viewer);CHKERRQ(ierr);
1859   PetscFunctionReturn(0);
1860 }
1861 
1862 #undef __FUNCT__
1863 #define __FUNCT__ "TaoDefaultConvergenceTest"
1864 /*@
1865    TaoDefaultConvergenceTest - Determines whether the solver should continue iterating
1866    or terminate.
1867 
1868    Collective on Tao
1869 
1870    Input Parameters:
1871 +  tao - the Tao context
1872 -  dummy - unused dummy context
1873 
1874    Output Parameter:
1875 .  reason - for terminating
1876 
1877    Notes:
1878    This routine checks the residual in the optimality conditions, the
1879    relative residual in the optimity conditions, the number of function
1880    evaluations, and the function value to test convergence.  Some
1881    solvers may use different convergence routines.
1882 
1883    Level: developer
1884 
1885 .seealso: TaoSetTolerances(),TaoGetConvergedReason(),TaoSetConvergedReason()
1886 @*/
1887 
1888 PetscErrorCode TaoDefaultConvergenceTest(Tao tao,void *dummy)
1889 {
1890   PetscInt           niter=tao->niter, nfuncs=PetscMax(tao->nfuncs,tao->nfuncgrads);
1891   PetscInt           max_funcs=tao->max_funcs;
1892   PetscReal          gnorm=tao->residual, gnorm0=tao->gnorm0;
1893   PetscReal          f=tao->fc, steptol=tao->steptol,trradius=tao->step;
1894   PetscReal          gatol=tao->gatol,grtol=tao->grtol,gttol=tao->gttol;
1895   PetscReal          fatol=tao->fatol,frtol=tao->frtol,catol=tao->catol,crtol=tao->crtol;
1896   PetscReal          fmin=tao->fmin, cnorm=tao->cnorm, cnorm0=tao->cnorm0;
1897   PetscReal          gnorm2;
1898   TaoConvergedReason reason=tao->reason;
1899   PetscErrorCode     ierr;
1900 
1901   PetscFunctionBegin;
1902   PetscValidHeaderSpecific(tao, TAO_CLASSID,1);
1903   if (reason != TAO_CONTINUE_ITERATING) {
1904     PetscFunctionReturn(0);
1905   }
1906   if (gnorm >= PETSC_INFINITY) {
1907     gnorm2 = gnorm;
1908   } else {
1909     gnorm2=gnorm*gnorm;
1910   }
1911 
1912   if (PetscIsInfOrNanReal(f)) {
1913     ierr = PetscInfo(tao,"Failed to converged, function value is Inf or NaN\n");CHKERRQ(ierr);
1914     reason = TAO_DIVERGED_NAN;
1915   } else if (f <= fmin && cnorm <=catol) {
1916     ierr = PetscInfo2(tao,"Converged due to function value %g < minimum function value %g\n", (double)f,(double)fmin);CHKERRQ(ierr);
1917     reason = TAO_CONVERGED_MINF;
1918   } else if (gnorm2 <= fatol && cnorm <=catol) {
1919     ierr = PetscInfo2(tao,"Converged due to estimated f(X) - f(X*) = %g < %g\n",(double)gnorm2,(double)fatol);CHKERRQ(ierr);
1920     reason = TAO_CONVERGED_FATOL;
1921   } else if (f != 0 && gnorm2 / PetscAbsReal(f)<= frtol && cnorm/PetscMax(cnorm0,1.0) <= crtol) {
1922     ierr = PetscInfo2(tao,"Converged due to estimated |f(X)-f(X*)|/f(X) = %g < %g\n",(double)(gnorm2/PetscAbsReal(f)),(double)frtol);CHKERRQ(ierr);
1923     reason = TAO_CONVERGED_FRTOL;
1924   } else if (gnorm<= gatol && cnorm <=catol) {
1925     ierr = PetscInfo2(tao,"Converged due to residual norm ||g(X)||=%g < %g\n",(double)gnorm,(double)gatol);CHKERRQ(ierr);
1926     reason = TAO_CONVERGED_GATOL;
1927   } else if ( f!=0 && PetscAbsReal(gnorm/f) <= grtol && cnorm <= crtol) {
1928     ierr = PetscInfo2(tao,"Converged due to residual ||g(X)||/|f(X)| =%g < %g\n",(double)(gnorm/f),(double)grtol);CHKERRQ(ierr);
1929     reason = TAO_CONVERGED_GRTOL;
1930   } else if (gnorm0 != 0 && gnorm/gnorm0 <= gttol && cnorm <= crtol) {
1931     ierr = PetscInfo2(tao,"Converged due to relative residual norm ||g(X)||/||g(X0)|| = %g < %g\n",(double)(gnorm/gnorm0),(double)gttol);CHKERRQ(ierr);
1932     reason = TAO_CONVERGED_GTTOL;
1933   } else if (nfuncs > max_funcs){
1934     ierr = PetscInfo2(tao,"Exceeded maximum number of function evaluations: %D > %D\n", nfuncs,max_funcs);CHKERRQ(ierr);
1935     reason = TAO_DIVERGED_MAXFCN;
1936   } else if ( tao->lsflag != 0 ){
1937     ierr = PetscInfo(tao,"Tao Line Search failure.\n");CHKERRQ(ierr);
1938     reason = TAO_DIVERGED_LS_FAILURE;
1939   } else if (trradius < steptol && niter > 0){
1940     ierr = PetscInfo2(tao,"Trust region/step size too small: %g < %g\n", (double)trradius,(double)steptol);CHKERRQ(ierr);
1941     reason = TAO_CONVERGED_STEPTOL;
1942   } else if (niter > tao->max_it) {
1943     ierr = PetscInfo2(tao,"Exceeded maximum number of iterations: %D > %D\n",niter,tao->max_it);CHKERRQ(ierr);
1944     reason = TAO_DIVERGED_MAXITS;
1945   } else {
1946     reason = TAO_CONTINUE_ITERATING;
1947   }
1948   tao->reason = reason;
1949   PetscFunctionReturn(0);
1950 }
1951 
1952 #undef __FUNCT__
1953 #define __FUNCT__ "TaoSetOptionsPrefix"
1954 /*@C
1955    TaoSetOptionsPrefix - Sets the prefix used for searching for all
1956    TAO options in the database.
1957 
1958 
1959    Logically Collective on Tao
1960 
1961    Input Parameters:
1962 +  tao - the Tao context
1963 -  prefix - the prefix string to prepend to all TAO option requests
1964 
1965    Notes:
1966    A hyphen (-) must NOT be given at the beginning of the prefix name.
1967    The first character of all runtime options is AUTOMATICALLY the hyphen.
1968 
1969    For example, to distinguish between the runtime options for two
1970    different TAO solvers, one could call
1971 .vb
1972       TaoSetOptionsPrefix(tao1,"sys1_")
1973       TaoSetOptionsPrefix(tao2,"sys2_")
1974 .ve
1975 
1976    This would enable use of different options for each system, such as
1977 .vb
1978       -sys1_tao_method blmvm -sys1_tao_gtol 1.e-3
1979       -sys2_tao_method lmvm  -sys2_tao_gtol 1.e-4
1980 .ve
1981 
1982 
1983    Level: advanced
1984 
1985 .seealso: TaoAppendOptionsPrefix(), TaoGetOptionsPrefix()
1986 @*/
1987 
1988 PetscErrorCode TaoSetOptionsPrefix(Tao tao, const char p[])
1989 {
1990   PetscErrorCode ierr;
1991 
1992   PetscFunctionBegin;
1993   ierr = PetscObjectSetOptionsPrefix((PetscObject)tao,p);CHKERRQ(ierr);
1994   if (tao->linesearch) {
1995     ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,p);CHKERRQ(ierr);
1996   }
1997   if (tao->ksp) {
1998     ierr = KSPSetOptionsPrefix(tao->ksp,p);CHKERRQ(ierr);
1999   }
2000   PetscFunctionReturn(0);
2001 }
2002 
2003 #undef __FUNCT__
2004 #define __FUNCT__ "TaoAppendOptionsPrefix"
2005 /*@C
2006    TaoAppendOptionsPrefix - Appends to the prefix used for searching for all
2007    TAO options in the database.
2008 
2009 
2010    Logically Collective on Tao
2011 
2012    Input Parameters:
2013 +  tao - the Tao solver context
2014 -  prefix - the prefix string to prepend to all TAO option requests
2015 
2016    Notes:
2017    A hyphen (-) must NOT be given at the beginning of the prefix name.
2018    The first character of all runtime options is AUTOMATICALLY the hyphen.
2019 
2020 
2021    Level: advanced
2022 
2023 .seealso: TaoSetOptionsPrefix(), TaoGetOptionsPrefix()
2024 @*/
2025 PetscErrorCode TaoAppendOptionsPrefix(Tao tao, const char p[])
2026 {
2027   PetscErrorCode ierr;
2028 
2029   PetscFunctionBegin;
2030   ierr = PetscObjectAppendOptionsPrefix((PetscObject)tao,p);CHKERRQ(ierr);
2031   if (tao->linesearch) {
2032     ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,p);CHKERRQ(ierr);
2033   }
2034   if (tao->ksp) {
2035     ierr = KSPSetOptionsPrefix(tao->ksp,p);CHKERRQ(ierr);
2036   }
2037   PetscFunctionReturn(0);
2038 }
2039 
2040 #undef __FUNCT__
2041 #define __FUNCT__ "TaoGetOptionsPrefix"
2042 /*@C
2043   TaoGetOptionsPrefix - Gets the prefix used for searching for all
2044   TAO options in the database
2045 
2046   Not Collective
2047 
2048   Input Parameters:
2049 . tao - the Tao context
2050 
2051   Output Parameters:
2052 . prefix - pointer to the prefix string used is returned
2053 
2054   Notes: On the fortran side, the user should pass in a string 'prefix' of
2055   sufficient length to hold the prefix.
2056 
2057   Level: advanced
2058 
2059 .seealso: TaoSetOptionsPrefix(), TaoAppendOptionsPrefix()
2060 @*/
2061 PetscErrorCode TaoGetOptionsPrefix(Tao tao, const char *p[])
2062 {
2063    return PetscObjectGetOptionsPrefix((PetscObject)tao,p);
2064 }
2065 
2066 #undef __FUNCT__
2067 #define __FUNCT__ "TaoSetType"
2068 /*@C
2069    TaoSetType - Sets the method for the unconstrained minimization solver.
2070 
2071    Collective on Tao
2072 
2073    Input Parameters:
2074 +  solver - the Tao solver context
2075 -  type - a known method
2076 
2077    Options Database Key:
2078 .  -tao_type <type> - Sets the method; use -help for a list
2079    of available methods (for instance, "-tao_type lmvm" or "-tao_type tron")
2080 
2081    Available methods include:
2082 +    nls - Newton's method with line search for unconstrained minimization
2083 .    ntr - Newton's method with trust region for unconstrained minimization
2084 .    ntl - Newton's method with trust region, line search for unconstrained minimization
2085 .    lmvm - Limited memory variable metric method for unconstrained minimization
2086 .    cg - Nonlinear conjugate gradient method for unconstrained minimization
2087 .    nm - Nelder-Mead algorithm for derivate-free unconstrained minimization
2088 .    tron - Newton Trust Region method for bound constrained minimization
2089 .    gpcg - Newton Trust Region method for quadratic bound constrained minimization
2090 .    blmvm - Limited memory variable metric method for bound constrained minimization
2091 -    pounders - Model-based algorithm pounder extended for nonlinear least squares
2092 
2093   Level: intermediate
2094 
2095 .seealso: TaoCreate(), TaoGetType(), TaoType
2096 
2097 @*/
2098 PetscErrorCode TaoSetType(Tao tao, const TaoType type)
2099 {
2100   PetscErrorCode ierr;
2101   PetscErrorCode (*create_xxx)(Tao);
2102   PetscBool      issame;
2103 
2104   PetscFunctionBegin;
2105   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2106 
2107   ierr = PetscObjectTypeCompare((PetscObject)tao,type,&issame);CHKERRQ(ierr);
2108   if (issame) PetscFunctionReturn(0);
2109 
2110   ierr = PetscFunctionListFind(TaoList, type, (void(**)(void))&create_xxx);CHKERRQ(ierr);
2111   if (!create_xxx) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested Tao type %s",type);
2112 
2113   /* Destroy the existing solver information */
2114   if (tao->ops->destroy) {
2115     ierr = (*tao->ops->destroy)(tao);CHKERRQ(ierr);
2116   }
2117   ierr = KSPDestroy(&tao->ksp);CHKERRQ(ierr);
2118   ierr = TaoLineSearchDestroy(&tao->linesearch);CHKERRQ(ierr);
2119   ierr = VecDestroy(&tao->gradient);CHKERRQ(ierr);
2120   ierr = VecDestroy(&tao->stepdirection);CHKERRQ(ierr);
2121 
2122   tao->ops->setup = 0;
2123   tao->ops->solve = 0;
2124   tao->ops->view  = 0;
2125   tao->ops->setfromoptions = 0;
2126   tao->ops->destroy = 0;
2127 
2128   tao->setupcalled = PETSC_FALSE;
2129 
2130   ierr = (*create_xxx)(tao);CHKERRQ(ierr);
2131   ierr = PetscObjectChangeTypeName((PetscObject)tao,type);CHKERRQ(ierr);
2132   PetscFunctionReturn(0);
2133 }
2134 
2135 #undef __FUNCT__
2136 #define __FUNCT__ "TaoRegister"
2137 /*MC
2138    TaoRegister - Adds a method to the TAO package for unconstrained minimization.
2139 
2140    Synopsis:
2141    TaoRegister(char *name_solver,char *path,char *name_Create,int (*routine_Create)(Tao))
2142 
2143    Not collective
2144 
2145    Input Parameters:
2146 +  sname - name of a new user-defined solver
2147 -  func - routine to Create method context
2148 
2149    Notes:
2150    TaoRegister() may be called multiple times to add several user-defined solvers.
2151 
2152    Sample usage:
2153 .vb
2154    TaoRegister("my_solver",MySolverCreate);
2155 .ve
2156 
2157    Then, your solver can be chosen with the procedural interface via
2158 $     TaoSetType(tao,"my_solver")
2159    or at runtime via the option
2160 $     -tao_type my_solver
2161 
2162    Level: advanced
2163 
2164 .seealso: TaoRegisterAll(), TaoRegisterDestroy()
2165 M*/
2166 PetscErrorCode TaoRegister(const char sname[], PetscErrorCode (*func)(Tao))
2167 {
2168   PetscErrorCode ierr;
2169 
2170   PetscFunctionBegin;
2171   ierr = PetscFunctionListAdd(&TaoList,sname, (void (*)(void))func);CHKERRQ(ierr);
2172   PetscFunctionReturn(0);
2173 }
2174 
2175 #undef __FUNCT__
2176 #define __FUNCT__ "TaoRegisterDestroy"
2177 /*@C
2178    TaoRegisterDestroy - Frees the list of minimization solvers that were
2179    registered by TaoRegisterDynamic().
2180 
2181    Not Collective
2182 
2183    Level: advanced
2184 
2185 .seealso: TaoRegisterAll(), TaoRegister()
2186 @*/
2187 PetscErrorCode TaoRegisterDestroy(void)
2188 {
2189   PetscErrorCode ierr;
2190   PetscFunctionBegin;
2191   ierr = PetscFunctionListDestroy(&TaoList);CHKERRQ(ierr);
2192   TaoRegisterAllCalled = PETSC_FALSE;
2193   PetscFunctionReturn(0);
2194 }
2195 
2196 #undef __FUNCT__
2197 #define __FUNCT__ "TaoGetIterationNumber"
2198 /*@
2199    TaoGetIterationNumber - Gets the number of Tao iterations completed
2200    at this time.
2201 
2202    Not Collective
2203 
2204    Input Parameter:
2205 .  tao - Tao context
2206 
2207    Output Parameter:
2208 .  iter - iteration number
2209 
2210    Notes:
2211    For example, during the computation of iteration 2 this would return 1.
2212 
2213 
2214    Level: intermediate
2215 
2216 .keywords: Tao, nonlinear, get, iteration, number,
2217 
2218 .seealso:   TaoGetLinearSolveIterations()
2219 @*/
2220 PetscErrorCode  TaoGetIterationNumber(Tao tao,PetscInt *iter)
2221 {
2222   PetscFunctionBegin;
2223   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2224   PetscValidIntPointer(iter,2);
2225   *iter = tao->niter;
2226   PetscFunctionReturn(0);
2227 }
2228 
2229 #undef __FUNCT__
2230 #define __FUNCT__ "TaoSetIterationNumber"
2231 /*@
2232    TaoSetIterationNumber - Sets the current iteration number.
2233 
2234    Not Collective
2235 
2236    Input Parameter:
2237 .  tao - Tao context
2238 .  iter - iteration number
2239 
2240    Level: developer
2241 
2242 .keywords: Tao, nonlinear, set, iteration, number,
2243 
2244 .seealso:   TaoGetLinearSolveIterations()
2245 @*/
2246 PetscErrorCode  TaoSetIterationNumber(Tao tao,PetscInt iter)
2247 {
2248   PetscErrorCode ierr;
2249 
2250   PetscFunctionBegin;
2251   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2252   ierr       = PetscObjectSAWsTakeAccess((PetscObject)tao);CHKERRQ(ierr);
2253   tao->niter = iter;
2254   ierr       = PetscObjectSAWsGrantAccess((PetscObject)tao);CHKERRQ(ierr);
2255   PetscFunctionReturn(0);
2256 }
2257 
2258 #undef __FUNCT__
2259 #define __FUNCT__ "TaoGetTotalIterationNumber"
2260 /*@
2261    TaoGetTotalIterationNumber - Gets the total number of Tao iterations
2262    completed. This number keeps accumulating if multiple solves
2263    are called with the Tao object.
2264 
2265    Not Collective
2266 
2267    Input Parameter:
2268 .  tao - Tao context
2269 
2270    Output Parameter:
2271 .  iter - iteration number
2272 
2273    Notes:
2274    The total iteration count is updated after each solve, if there is a current
2275    TaoSolve() in progress then those iterations are not yet counted.
2276 
2277    Level: intermediate
2278 
2279 .keywords: Tao, nonlinear, get, iteration, number,
2280 
2281 .seealso:   TaoGetLinearSolveIterations()
2282 @*/
2283 PetscErrorCode  TaoGetTotalIterationNumber(Tao tao,PetscInt *iter)
2284 {
2285   PetscFunctionBegin;
2286   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2287   PetscValidIntPointer(iter,2);
2288   *iter = tao->ntotalits;
2289   PetscFunctionReturn(0);
2290 }
2291 
2292 #undef __FUNCT__
2293 #define __FUNCT__ "TaoSetTotalIterationNumber"
2294 /*@
2295    TaoSetTotalIterationNumber - Sets the current total iteration number.
2296 
2297    Not Collective
2298 
2299    Input Parameter:
2300 .  tao - Tao context
2301 .  iter - iteration number
2302 
2303    Level: developer
2304 
2305 .keywords: Tao, nonlinear, set, iteration, number,
2306 
2307 .seealso:   TaoGetLinearSolveIterations()
2308 @*/
2309 PetscErrorCode  TaoSetTotalIterationNumber(Tao tao,PetscInt iter)
2310 {
2311   PetscErrorCode ierr;
2312 
2313   PetscFunctionBegin;
2314   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2315   ierr       = PetscObjectSAWsTakeAccess((PetscObject)tao);CHKERRQ(ierr);
2316   tao->ntotalits = iter;
2317   ierr       = PetscObjectSAWsGrantAccess((PetscObject)tao);CHKERRQ(ierr);
2318   PetscFunctionReturn(0);
2319 }
2320 
2321 #undef __FUNCT__
2322 #define __FUNCT__ "TaoSetConvergedReason"
2323 /*@
2324   TaoSetConvergedReason - Sets the termination flag on a Tao object
2325 
2326   Logically Collective on Tao
2327 
2328   Input Parameters:
2329 + tao - the Tao context
2330 - reason - one of
2331 $     TAO_CONVERGED_ATOL (2),
2332 $     TAO_CONVERGED_RTOL (3),
2333 $     TAO_CONVERGED_STEPTOL (4),
2334 $     TAO_CONVERGED_MINF (5),
2335 $     TAO_CONVERGED_USER (6),
2336 $     TAO_DIVERGED_MAXITS (-2),
2337 $     TAO_DIVERGED_NAN (-4),
2338 $     TAO_DIVERGED_MAXFCN (-5),
2339 $     TAO_DIVERGED_LS_FAILURE (-6),
2340 $     TAO_DIVERGED_TR_REDUCTION (-7),
2341 $     TAO_DIVERGED_USER (-8),
2342 $     TAO_CONTINUE_ITERATING (0)
2343 
2344    Level: intermediate
2345 
2346 @*/
2347 PetscErrorCode TaoSetConvergedReason(Tao tao, TaoConvergedReason reason)
2348 {
2349   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2350   PetscFunctionBegin;
2351   tao->reason = reason;
2352   PetscFunctionReturn(0);
2353 }
2354 
2355 #undef __FUNCT__
2356 #define __FUNCT__ "TaoGetConvergedReason"
2357 /*@
2358    TaoGetConvergedReason - Gets the reason the Tao iteration was stopped.
2359 
2360    Not Collective
2361 
2362    Input Parameter:
2363 .  tao - the Tao solver context
2364 
2365    Output Parameter:
2366 .  reason - one of
2367 $  TAO_CONVERGED_FATOL (1)           f(X)-f(X*) <= fatol
2368 $  TAO_CONVERGED_FRTOL (2)           |f(X) - f(X*)|/|f(X)| < frtol
2369 $  TAO_CONVERGED_GATOL (3)           ||g(X)|| < gatol
2370 $  TAO_CONVERGED_GRTOL (4)           ||g(X)|| / f(X)  < grtol
2371 $  TAO_CONVERGED_GTTOL (5)           ||g(X)|| / ||g(X0)|| < gttol
2372 $  TAO_CONVERGED_STEPTOL (6)         step size small
2373 $  TAO_CONVERGED_MINF (7)            F < F_min
2374 $  TAO_CONVERGED_USER (8)            User defined
2375 $  TAO_DIVERGED_MAXITS (-2)          its > maxits
2376 $  TAO_DIVERGED_NAN (-4)             Numerical problems
2377 $  TAO_DIVERGED_MAXFCN (-5)          fevals > max_funcsals
2378 $  TAO_DIVERGED_LS_FAILURE (-6)      line search failure
2379 $  TAO_DIVERGED_TR_REDUCTION (-7)    trust region failure
2380 $  TAO_DIVERGED_USER(-8)             (user defined)
2381  $  TAO_CONTINUE_ITERATING (0)
2382 
2383    where
2384 +  X - current solution
2385 .  X0 - initial guess
2386 .  f(X) - current function value
2387 .  f(X*) - true solution (estimated)
2388 .  g(X) - current gradient
2389 .  its - current iterate number
2390 .  maxits - maximum number of iterates
2391 .  fevals - number of function evaluations
2392 -  max_funcsals - maximum number of function evaluations
2393 
2394    Level: intermediate
2395 
2396 .seealso: TaoSetConvergenceTest(), TaoSetTolerances()
2397 
2398 @*/
2399 PetscErrorCode TaoGetConvergedReason(Tao tao, TaoConvergedReason *reason)
2400 {
2401   PetscFunctionBegin;
2402   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2403   PetscValidPointer(reason,2);
2404   *reason = tao->reason;
2405   PetscFunctionReturn(0);
2406 }
2407 
2408 #undef __FUNCT__
2409 #define __FUNCT__ "TaoGetSolutionStatus"
2410 /*@
2411   TaoGetSolutionStatus - Get the current iterate, objective value,
2412   residual, infeasibility, and termination
2413 
2414   Not Collective
2415 
2416    Input Parameters:
2417 .  tao - the Tao context
2418 
2419    Output Parameters:
2420 +  iterate - the current iterate number (>=0)
2421 .  f - the current function value
2422 .  gnorm - the square of the gradient norm, duality gap, or other measure indicating distance from optimality.
2423 .  cnorm - the infeasibility of the current solution with regard to the constraints.
2424 .  xdiff - the step length or trust region radius of the most recent iterate.
2425 -  reason - The termination reason, which can equal TAO_CONTINUE_ITERATING
2426 
2427    Level: intermediate
2428 
2429    Note:
2430    TAO returns the values set by the solvers in the routine TaoMonitor().
2431 
2432    Note:
2433    If any of the output arguments are set to NULL, no corresponding value will be returned.
2434 
2435 .seealso: TaoMonitor(), TaoGetConvergedReason()
2436 @*/
2437 PetscErrorCode TaoGetSolutionStatus(Tao tao, PetscInt *its, PetscReal *f, PetscReal *gnorm, PetscReal *cnorm, PetscReal *xdiff, TaoConvergedReason *reason)
2438 {
2439   PetscFunctionBegin;
2440   if (its) *its=tao->niter;
2441   if (f) *f=tao->fc;
2442   if (gnorm) *gnorm=tao->residual;
2443   if (cnorm) *cnorm=tao->cnorm;
2444   if (reason) *reason=tao->reason;
2445   if (xdiff) *xdiff=tao->step;
2446   PetscFunctionReturn(0);
2447 }
2448 
2449 #undef __FUNCT__
2450 #define __FUNCT__ "TaoGetType"
2451 /*@C
2452    TaoGetType - Gets the current Tao algorithm.
2453 
2454    Not Collective
2455 
2456    Input Parameter:
2457 .  tao - the Tao solver context
2458 
2459    Output Parameter:
2460 .  type - Tao method
2461 
2462    Level: intermediate
2463 
2464 @*/
2465 PetscErrorCode TaoGetType(Tao tao, const TaoType *type)
2466 {
2467   PetscFunctionBegin;
2468   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2469   PetscValidPointer(type,2);
2470   *type=((PetscObject)tao)->type_name;
2471   PetscFunctionReturn(0);
2472 }
2473 
2474 #undef __FUNCT__
2475 #define __FUNCT__ "TaoMonitor"
2476 /*@C
2477   TaoMonitor - Monitor the solver and the current solution.  This
2478   routine will record the iteration number and residual statistics,
2479   call any monitors specified by the user, and calls the convergence-check routine.
2480 
2481    Input Parameters:
2482 +  tao - the Tao context
2483 .  its - the current iterate number (>=0)
2484 .  f - the current objective function value
2485 .  res - the gradient norm, square root of the duality gap, or other measure indicating distince from optimality.  This measure will be recorded and
2486           used for some termination tests.
2487 .  cnorm - the infeasibility of the current solution with regard to the constraints.
2488 -  steplength - multiple of the step direction added to the previous iterate.
2489 
2490    Output Parameters:
2491 .  reason - The termination reason, which can equal TAO_CONTINUE_ITERATING
2492 
2493    Options Database Key:
2494 .  -tao_monitor - Use the default monitor, which prints statistics to standard output
2495 
2496 .seealso TaoGetConvergedReason(), TaoDefaultMonitor(), TaoSetMonitor()
2497 
2498    Level: developer
2499 
2500 @*/
2501 PetscErrorCode TaoMonitor(Tao tao, PetscInt its, PetscReal f, PetscReal res, PetscReal cnorm, PetscReal steplength, TaoConvergedReason *reason)
2502 {
2503   PetscErrorCode ierr;
2504   PetscInt       i;
2505 
2506   PetscFunctionBegin;
2507   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2508   tao->fc = f;
2509   tao->residual = res;
2510   tao->cnorm = cnorm;
2511   tao->step = steplength;
2512   if (its == 0) {
2513     tao->cnorm0 = cnorm; tao->gnorm0 = res;
2514   }
2515   TaoLogConvergenceHistory(tao,f,res,cnorm,tao->ksp_its);
2516   if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(res)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
2517   if (tao->ops->convergencetest) {
2518     ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
2519   }
2520   for (i=0;i<tao->numbermonitors;i++) {
2521     ierr = (*tao->monitor[i])(tao,tao->monitorcontext[i]);CHKERRQ(ierr);
2522   }
2523   *reason = tao->reason;
2524   PetscFunctionReturn(0);
2525 }
2526 
2527 #undef __FUNCT__
2528 #define __FUNCT__ "TaoSetConvergenceHistory"
2529 /*@
2530    TaoSetConvergenceHistory - Sets the array used to hold the convergence history.
2531 
2532    Logically Collective on Tao
2533 
2534    Input Parameters:
2535 +  tao - the Tao solver context
2536 .  obj   - array to hold objective value history
2537 .  resid - array to hold residual history
2538 .  cnorm - array to hold constraint violation history
2539 .  lits - integer array holds the number of linear iterations for each Tao iteration
2540 .  na  - size of obj, resid, and cnorm
2541 -  reset - PetscTrue indicates each new minimization resets the history counter to zero,
2542            else it continues storing new values for new minimizations after the old ones
2543 
2544    Notes:
2545    If set, TAO will fill the given arrays with the indicated
2546    information at each iteration.  If 'obj','resid','cnorm','lits' are
2547    *all* NULL then space (using size na, or 1000 if na is PETSC_DECIDE or
2548    PETSC_DEFAULT) is allocated for the history.
2549    If not all are NULL, then only the non-NULL information categories
2550    will be stored, the others will be ignored.
2551 
2552    Any convergence information after iteration number 'na' will not be stored.
2553 
2554    This routine is useful, e.g., when running a code for purposes
2555    of accurate performance monitoring, when no I/O should be done
2556    during the section of code that is being timed.
2557 
2558    Level: intermediate
2559 
2560 .seealso: TaoGetConvergenceHistory()
2561 
2562 @*/
2563 PetscErrorCode TaoSetConvergenceHistory(Tao tao, PetscReal *obj, PetscReal *resid, PetscReal *cnorm, PetscInt *lits, PetscInt na,PetscBool reset)
2564 {
2565   PetscErrorCode ierr;
2566 
2567   PetscFunctionBegin;
2568   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2569   if (obj) PetscValidScalarPointer(obj,2);
2570   if (resid) PetscValidScalarPointer(resid,3);
2571   if (cnorm) PetscValidScalarPointer(cnorm,4);
2572   if (lits) PetscValidIntPointer(lits,5);
2573 
2574   if (na == PETSC_DECIDE || na == PETSC_DEFAULT) na = 1000;
2575   if (!obj && !resid && !cnorm && !lits) {
2576     ierr = PetscCalloc1(na,&obj);CHKERRQ(ierr);
2577     ierr = PetscCalloc1(na,&resid);CHKERRQ(ierr);
2578     ierr = PetscCalloc1(na,&cnorm);CHKERRQ(ierr);
2579     ierr = PetscCalloc1(na,&lits);CHKERRQ(ierr);
2580     tao->hist_malloc=PETSC_TRUE;
2581   }
2582 
2583   tao->hist_obj = obj;
2584   tao->hist_resid = resid;
2585   tao->hist_cnorm = cnorm;
2586   tao->hist_lits = lits;
2587   tao->hist_max   = na;
2588   tao->hist_reset = reset;
2589   tao->hist_len = 0;
2590   PetscFunctionReturn(0);
2591 }
2592 
2593 #undef __FUNCT__
2594 #define __FUNCT__ "TaoGetConvergenceHistory"
2595 /*@C
2596    TaoGetConvergenceHistory - Gets the arrays used to hold the convergence history.
2597 
2598    Collective on Tao
2599 
2600    Input Parameter:
2601 .  tao - the Tao context
2602 
2603    Output Parameters:
2604 +  obj   - array used to hold objective value history
2605 .  resid - array used to hold residual history
2606 .  cnorm - array used to hold constraint violation history
2607 .  lits  - integer array used to hold linear solver iteration count
2608 -  nhist  - size of obj, resid, cnorm, and lits (will be less than or equal to na given in TaoSetHistory)
2609 
2610    Notes:
2611     This routine must be preceded by calls to TaoSetConvergenceHistory()
2612     and TaoSolve(), otherwise it returns useless information.
2613 
2614     The calling sequence for this routine in Fortran is
2615 $   call TaoGetConvergenceHistory(Tao tao, PetscInt nhist, PetscErrorCode ierr)
2616 
2617    This routine is useful, e.g., when running a code for purposes
2618    of accurate performance monitoring, when no I/O should be done
2619    during the section of code that is being timed.
2620 
2621    Level: advanced
2622 
2623 .seealso: TaoSetConvergenceHistory()
2624 
2625 @*/
2626 PetscErrorCode TaoGetConvergenceHistory(Tao tao, PetscReal **obj, PetscReal **resid, PetscReal **cnorm, PetscInt **lits, PetscInt *nhist)
2627 {
2628   PetscFunctionBegin;
2629   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2630   if (obj)   *obj   = tao->hist_obj;
2631   if (cnorm) *cnorm = tao->hist_cnorm;
2632   if (resid) *resid = tao->hist_resid;
2633   if (nhist) *nhist   = tao->hist_len;
2634   PetscFunctionReturn(0);
2635 }
2636 
2637 #undef __FUNCT__
2638 #define __FUNCT__ "TaoSetApplicationContext"
2639 /*@
2640    TaoSetApplicationContext - Sets the optional user-defined context for
2641    a solver.
2642 
2643    Logically Collective on Tao
2644 
2645    Input Parameters:
2646 +  tao  - the Tao context
2647 -  usrP - optional user context
2648 
2649    Level: intermediate
2650 
2651 .seealso: TaoGetApplicationContext(), TaoSetApplicationContext()
2652 @*/
2653 PetscErrorCode  TaoSetApplicationContext(Tao tao,void *usrP)
2654 {
2655   PetscFunctionBegin;
2656   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2657   tao->user = usrP;
2658   PetscFunctionReturn(0);
2659 }
2660 
2661 #undef __FUNCT__
2662 #define __FUNCT__ "TaoGetApplicationContext"
2663 /*@
2664    TaoGetApplicationContext - Gets the user-defined context for a
2665    TAO solvers.
2666 
2667    Not Collective
2668 
2669    Input Parameter:
2670 .  tao  - Tao context
2671 
2672    Output Parameter:
2673 .  usrP - user context
2674 
2675    Level: intermediate
2676 
2677 .seealso: TaoSetApplicationContext()
2678 @*/
2679 PetscErrorCode  TaoGetApplicationContext(Tao tao,void *usrP)
2680 {
2681   PetscFunctionBegin;
2682   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2683   *(void**)usrP = tao->user;
2684   PetscFunctionReturn(0);
2685 }
2686 
2687 #undef __FUNCT__
2688 #define __FUNCT__ "TaoSetGradientNorm"
2689 /*@
2690    TaoSetGradientNorm - Sets the matrix used to define the inner product that measures the size of the gradient.
2691 
2692    Collective on tao
2693 
2694    Input Parameters:
2695 +  tao  - the Tao context
2696 -  M    - gradient norm
2697 
2698    Level: beginner
2699 
2700 .seealso: TaoGetGradientNorm(), TaoGradientNorm()
2701 @*/
2702 PetscErrorCode  TaoSetGradientNorm(Tao tao, Mat M)
2703 {
2704   PetscErrorCode ierr;
2705 
2706   PetscFunctionBegin;
2707   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2708 
2709   if (tao->gradient_norm) {
2710     ierr = PetscObjectDereference((PetscObject)tao->gradient_norm);CHKERRQ(ierr);
2711     ierr = VecDestroy(&tao->gradient_norm_tmp);CHKERRQ(ierr);
2712   }
2713 
2714   ierr = PetscObjectReference((PetscObject)M);CHKERRQ(ierr);
2715   tao->gradient_norm = M;
2716   ierr = MatCreateVecs(M, NULL, &tao->gradient_norm_tmp);CHKERRQ(ierr);
2717   PetscFunctionReturn(0);
2718 }
2719 
2720 #undef __FUNCT__
2721 #define __FUNCT__ "TaoGetGradientNorm"
2722 /*@
2723    TaoGetGradientNorm - Returns the matrix used to define the inner product for measuring the size of the gradient.
2724 
2725    Not Collective
2726 
2727    Input Parameter:
2728 .  tao  - Tao context
2729 
2730    Output Parameter:
2731 .  M - gradient norm
2732 
2733    Level: beginner
2734 
2735 .seealso: TaoSetGradientNorm(), TaoGradientNorm()
2736 @*/
2737 PetscErrorCode  TaoGetGradientNorm(Tao tao, Mat *M)
2738 {
2739   PetscFunctionBegin;
2740   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2741   *M = tao->gradient_norm;
2742   PetscFunctionReturn(0);
2743 }
2744 
2745 #undef __FUNCT__
2746 #define __FUNCT__ "TaoGradientNorm"
2747 /*c
2748    TaoGradientNorm - Compute the norm with respect to the inner product the user has set.
2749 
2750    Collective on tao
2751 
2752    Input Parameter:
2753 .  tao      - the Tao context
2754 .  gradient - the gradient to be computed
2755 .  norm     - the norm type
2756 
2757    Output Parameter:
2758 .  gnorm    - the gradient norm
2759 
2760    Level: developer
2761 
2762 .seealso: TaoSetGradientNorm(), TaoGetGradientNorm()
2763 @*/
2764 PetscErrorCode  TaoGradientNorm(Tao tao, Vec gradient, NormType type, PetscReal *gnorm)
2765 {
2766   PetscErrorCode ierr;
2767 
2768   PetscFunctionBegin;
2769   PetscValidHeaderSpecific(gradient,VEC_CLASSID,1);
2770 
2771   if (tao->gradient_norm) {
2772     PetscScalar gnorms;
2773 
2774     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.");
2775     ierr = MatMult(tao->gradient_norm, gradient, tao->gradient_norm_tmp);CHKERRQ(ierr);
2776     ierr = VecDot(gradient, tao->gradient_norm_tmp, &gnorms);CHKERRQ(ierr);
2777     *gnorm = PetscRealPart(PetscSqrtScalar(gnorms));
2778   } else {
2779     ierr = VecNorm(gradient, type, gnorm);CHKERRQ(ierr);
2780   }
2781   PetscFunctionReturn(0);
2782 }
2783 
2784 
2785