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