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