xref: /petsc/src/tao/interface/taosolver.c (revision 125ddc8ac7d0bab4f90817717775c95e255c43d8)
1441846f8SBarry Smith #define TAO_DLL
2a7e14dcfSSatish Balay 
3af0996ceSBarry Smith #include <petsc/private/taoimpl.h> /*I "petsctao.h" I*/
4a7e14dcfSSatish Balay 
5441846f8SBarry Smith PetscBool TaoRegisterAllCalled = PETSC_FALSE;
6441846f8SBarry Smith PetscFunctionList TaoList = NULL;
7a7e14dcfSSatish Balay 
8441846f8SBarry Smith PetscClassId TAO_CLASSID;
9441846f8SBarry Smith PetscLogEvent Tao_Solve, Tao_ObjectiveEval, Tao_GradientEval, Tao_ObjGradientEval, Tao_HessianEval, Tao_ConstraintsEval, Tao_JacobianEval;
10a7e14dcfSSatish Balay 
115b0d99aeSBarry Smith const char *TaoSubSetTypes[] = {  "subvec","mask","matrixfree","TaoSubSetType","TAO_SUBSET_",0};
12a7e14dcfSSatish Balay 
13a7e14dcfSSatish Balay #undef __FUNCT__
14a7e14dcfSSatish Balay #define __FUNCT__ "TaoCreate"
15a7e14dcfSSatish Balay /*@
16a7e14dcfSSatish Balay   TaoCreate - Creates a TAO solver
17a7e14dcfSSatish Balay 
18a7e14dcfSSatish Balay   Collective on MPI_Comm
19a7e14dcfSSatish Balay 
20a7e14dcfSSatish Balay   Input Parameter:
21a7e14dcfSSatish Balay . comm - MPI communicator
22a7e14dcfSSatish Balay 
23a7e14dcfSSatish Balay   Output Parameter:
24441846f8SBarry Smith . newtao - the new Tao context
25a7e14dcfSSatish Balay 
26a7e14dcfSSatish Balay   Available methods include:
2758417fe7SBarry Smith +    nls - Newton's method with line search for unconstrained minimization
2858417fe7SBarry Smith .    ntr - Newton's method with trust region for unconstrained minimization
2958417fe7SBarry Smith .    ntl - Newton's method with trust region, line search for unconstrained minimization
3058417fe7SBarry Smith .    lmvm - Limited memory variable metric method for unconstrained minimization
3158417fe7SBarry Smith .    cg - Nonlinear conjugate gradient method for unconstrained minimization
3258417fe7SBarry Smith .    nm - Nelder-Mead algorithm for derivate-free unconstrained minimization
3358417fe7SBarry Smith .    tron - Newton Trust Region method for bound constrained minimization
3458417fe7SBarry Smith .    gpcg - Newton Trust Region method for quadratic bound constrained minimization
3558417fe7SBarry Smith .    blmvm - Limited memory variable metric method for bound constrained minimization
3658417fe7SBarry Smith .    lcl - Linearly constrained Lagrangian method for pde-constrained minimization
3758417fe7SBarry Smith -    pounders - Model-based algorithm for nonlinear least squares
38a7e14dcfSSatish Balay 
39a7e14dcfSSatish Balay    Options Database Keys:
4058417fe7SBarry Smith .   -tao_type - select which method TAO should use
41a7e14dcfSSatish Balay 
42a7e14dcfSSatish Balay    Level: beginner
43a7e14dcfSSatish Balay 
44a7e14dcfSSatish Balay .seealso: TaoSolve(), TaoDestroy()
45a7e14dcfSSatish Balay @*/
46441846f8SBarry Smith PetscErrorCode TaoCreate(MPI_Comm comm, Tao *newtao)
47a7e14dcfSSatish Balay {
48a7e14dcfSSatish Balay   PetscErrorCode ierr;
49441846f8SBarry Smith   Tao            tao;
50a7e14dcfSSatish Balay 
51a7e14dcfSSatish Balay   PetscFunctionBegin;
52a7e14dcfSSatish Balay   PetscValidPointer(newtao,2);
5347a47007SBarry Smith   *newtao = NULL;
54a7e14dcfSSatish Balay 
55a7e14dcfSSatish Balay   ierr = TaoInitializePackage();CHKERRQ(ierr);
56a7e14dcfSSatish Balay   ierr = TaoLineSearchInitializePackage();CHKERRQ(ierr);
57a7e14dcfSSatish Balay 
5873107ff1SLisandro Dalcin   ierr = PetscHeaderCreate(tao,TAO_CLASSID,"Tao","Optimization solver","Tao",comm,TaoDestroy,TaoView);CHKERRQ(ierr);
59a7e14dcfSSatish Balay   tao->ops->computeobjective=0;
60a7e14dcfSSatish Balay   tao->ops->computeobjectiveandgradient=0;
61a7e14dcfSSatish Balay   tao->ops->computegradient=0;
62a7e14dcfSSatish Balay   tao->ops->computehessian=0;
63a7e14dcfSSatish Balay   tao->ops->computeseparableobjective=0;
64a7e14dcfSSatish Balay   tao->ops->computeconstraints=0;
65a7e14dcfSSatish Balay   tao->ops->computejacobian=0;
66a7e14dcfSSatish Balay   tao->ops->computejacobianequality=0;
67a7e14dcfSSatish Balay   tao->ops->computejacobianinequality=0;
68a7e14dcfSSatish Balay   tao->ops->computeequalityconstraints=0;
69a7e14dcfSSatish Balay   tao->ops->computeinequalityconstraints=0;
70a7e14dcfSSatish Balay   tao->ops->convergencetest=TaoDefaultConvergenceTest;
71a7e14dcfSSatish Balay   tao->ops->convergencedestroy=0;
72a7e14dcfSSatish Balay   tao->ops->computedual=0;
73a7e14dcfSSatish Balay   tao->ops->setup=0;
74a7e14dcfSSatish Balay   tao->ops->solve=0;
75a7e14dcfSSatish Balay   tao->ops->view=0;
76a7e14dcfSSatish Balay   tao->ops->setfromoptions=0;
77a7e14dcfSSatish Balay   tao->ops->destroy=0;
78a7e14dcfSSatish Balay 
7947a47007SBarry Smith   tao->solution=NULL;
8047a47007SBarry Smith   tao->gradient=NULL;
8147a47007SBarry Smith   tao->sep_objective = NULL;
8247a47007SBarry Smith   tao->constraints=NULL;
8347a47007SBarry Smith   tao->constraints_equality=NULL;
8447a47007SBarry Smith   tao->constraints_inequality=NULL;
85*125ddc8aSJason Sarich   tao->sep_weights_v=NULL;
86*125ddc8aSJason Sarich   tao->sep_weights_w=NULL;
8747a47007SBarry Smith   tao->stepdirection=NULL;
888931d482SJason Sarich   tao->niter=0;
898931d482SJason Sarich   tao->ntotalits=0;
9047a47007SBarry Smith   tao->XL = NULL;
9147a47007SBarry Smith   tao->XU = NULL;
9247a47007SBarry Smith   tao->IL = NULL;
9347a47007SBarry Smith   tao->IU = NULL;
9447a47007SBarry Smith   tao->DI = NULL;
9547a47007SBarry Smith   tao->DE = NULL;
96a9603a14SPatrick Farrell   tao->gradient_norm = NULL;
97a9603a14SPatrick Farrell   tao->gradient_norm_tmp = NULL;
9847a47007SBarry Smith   tao->hessian = NULL;
9947a47007SBarry Smith   tao->hessian_pre = NULL;
10047a47007SBarry Smith   tao->jacobian = NULL;
10147a47007SBarry Smith   tao->jacobian_pre = NULL;
10247a47007SBarry Smith   tao->jacobian_state = NULL;
10347a47007SBarry Smith   tao->jacobian_state_pre = NULL;
10447a47007SBarry Smith   tao->jacobian_state_inv = NULL;
10547a47007SBarry Smith   tao->jacobian_design = NULL;
10647a47007SBarry Smith   tao->jacobian_design_pre = NULL;
10747a47007SBarry Smith   tao->jacobian_equality = NULL;
10847a47007SBarry Smith   tao->jacobian_equality_pre = NULL;
10947a47007SBarry Smith   tao->jacobian_inequality = NULL;
11047a47007SBarry Smith   tao->jacobian_inequality_pre = NULL;
11147a47007SBarry Smith   tao->state_is = NULL;
11247a47007SBarry Smith   tao->design_is = NULL;
113a7e14dcfSSatish Balay 
114a7e14dcfSSatish Balay   tao->max_it     = 10000;
115a7e14dcfSSatish Balay   tao->max_funcs   = 10000;
1166f4723b1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE)
1176f4723b1SBarry Smith   tao->gatol       = 1e-5;
1186f4723b1SBarry Smith   tao->grtol       = 1e-5;
1196f4723b1SBarry Smith #else
120a7e14dcfSSatish Balay   tao->gatol       = 1e-8;
121a7e14dcfSSatish Balay   tao->grtol       = 1e-8;
1226f4723b1SBarry Smith #endif
123a7e14dcfSSatish Balay   tao->crtol       = 0.0;
1246552cf8aSJason Sarich   tao->catol       = 0.0;
125a7e14dcfSSatish Balay   tao->steptol     = 0.0;
1266552cf8aSJason Sarich   tao->gttol       = 0.0;
127e270355aSBarry Smith   tao->trust0      = PETSC_INFINITY;
1286f4723b1SBarry Smith   tao->fmin        = PETSC_NINFINITY;
129ae93cb3cSJason Sarich   tao->hist_malloc = PETSC_FALSE;
130a7e14dcfSSatish Balay   tao->hist_reset = PETSC_TRUE;
131a7e14dcfSSatish Balay   tao->hist_max = 0;
132a7e14dcfSSatish Balay   tao->hist_len = 0;
13347a47007SBarry Smith   tao->hist_obj = NULL;
13447a47007SBarry Smith   tao->hist_resid = NULL;
13547a47007SBarry Smith   tao->hist_cnorm = NULL;
136ae93cb3cSJason Sarich   tao->hist_lits = NULL;
137a7e14dcfSSatish Balay 
138a7e14dcfSSatish Balay   tao->numbermonitors=0;
139a7e14dcfSSatish Balay   tao->viewsolution=PETSC_FALSE;
140a7e14dcfSSatish Balay   tao->viewhessian=PETSC_FALSE;
141a7e14dcfSSatish Balay   tao->viewgradient=PETSC_FALSE;
142a7e14dcfSSatish Balay   tao->viewjacobian=PETSC_FALSE;
143a7e14dcfSSatish Balay   tao->viewconstraints = PETSC_FALSE;
144a7e14dcfSSatish Balay 
1456552cf8aSJason Sarich   /* These flags prevents algorithms from overriding user options */
1466552cf8aSJason Sarich   tao->max_it_changed   =PETSC_FALSE;
1476552cf8aSJason Sarich   tao->max_funcs_changed=PETSC_FALSE;
1486552cf8aSJason Sarich   tao->gatol_changed    =PETSC_FALSE;
1496552cf8aSJason Sarich   tao->grtol_changed    =PETSC_FALSE;
1506552cf8aSJason Sarich   tao->gttol_changed    =PETSC_FALSE;
1516552cf8aSJason Sarich   tao->steptol_changed  =PETSC_FALSE;
1526552cf8aSJason Sarich   tao->trust0_changed   =PETSC_FALSE;
1536552cf8aSJason Sarich   tao->fmin_changed     =PETSC_FALSE;
1546552cf8aSJason Sarich   tao->catol_changed    =PETSC_FALSE;
1556552cf8aSJason Sarich   tao->crtol_changed    =PETSC_FALSE;
156a7e14dcfSSatish Balay   ierr = TaoResetStatistics(tao);CHKERRQ(ierr);
157a7e14dcfSSatish Balay   *newtao = tao;
158a7e14dcfSSatish Balay   PetscFunctionReturn(0);
159a7e14dcfSSatish Balay }
160a7e14dcfSSatish Balay 
161a7e14dcfSSatish Balay #undef __FUNCT__
162a7e14dcfSSatish Balay #define __FUNCT__ "TaoSolve"
163a7e14dcfSSatish Balay /*@
164a7e14dcfSSatish Balay   TaoSolve - Solves an optimization problem min F(x) s.t. l <= x <= u
165a7e14dcfSSatish Balay 
166441846f8SBarry Smith   Collective on Tao
167a7e14dcfSSatish Balay 
168a7e14dcfSSatish Balay   Input Parameters:
169441846f8SBarry Smith . tao - the Tao context
170a7e14dcfSSatish Balay 
171a7e14dcfSSatish Balay   Notes:
172441846f8SBarry Smith   The user must set up the Tao with calls to TaoSetInitialVector(),
173a7e14dcfSSatish Balay   TaoSetObjectiveRoutine(),
174a7e14dcfSSatish Balay   TaoSetGradientRoutine(), and (if using 2nd order method) TaoSetHessianRoutine().
175a7e14dcfSSatish Balay 
176a7e14dcfSSatish Balay   Level: beginner
177a7e14dcfSSatish Balay 
178a7e14dcfSSatish Balay .seealso: TaoCreate(), TaoSetObjectiveRoutine(), TaoSetGradientRoutine(), TaoSetHessianRoutine()
179a7e14dcfSSatish Balay  @*/
180441846f8SBarry Smith PetscErrorCode TaoSolve(Tao tao)
181a7e14dcfSSatish Balay {
182a7e14dcfSSatish Balay   PetscErrorCode   ierr;
183e2379f4fSBarry Smith   static PetscBool set = PETSC_FALSE;
18447a47007SBarry Smith 
185a7e14dcfSSatish Balay   PetscFunctionBegin;
186441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
187e2379f4fSBarry Smith   ierr = PetscCitationsRegister("@TechReport{tao-user-ref,\n"
188e2379f4fSBarry Smith                                 "title   = {Toolkit for Advanced Optimization (TAO) Users Manual},\n"
189e2379f4fSBarry Smith                                 "author  = {Todd Munson and Jason Sarich and Stefan Wild and Steve Benson and Lois Curfman McInnes},\n"
190e2379f4fSBarry Smith                                 "Institution = {Argonne National Laboratory},\n"
191e2379f4fSBarry Smith                                 "Year   = 2014,\n"
192e2379f4fSBarry Smith                                 "Number = {ANL/MCS-TM-322 - Revision 3.5},\n"
193302440fdSBarry Smith                                 "url    = {http://www.mcs.anl.gov/tao}\n}\n",&set);CHKERRQ(ierr);
194e2379f4fSBarry Smith 
195a7e14dcfSSatish Balay   ierr = TaoSetUp(tao);CHKERRQ(ierr);
196a7e14dcfSSatish Balay   ierr = TaoResetStatistics(tao);CHKERRQ(ierr);
197a7e14dcfSSatish Balay   if (tao->linesearch) {
198a7e14dcfSSatish Balay     ierr = TaoLineSearchReset(tao->linesearch);CHKERRQ(ierr);
199a7e14dcfSSatish Balay   }
200a7e14dcfSSatish Balay 
201441846f8SBarry Smith   ierr = PetscLogEventBegin(Tao_Solve,tao,0,0,0);CHKERRQ(ierr);
202a7e14dcfSSatish Balay   if (tao->ops->solve){ ierr = (*tao->ops->solve)(tao);CHKERRQ(ierr); }
203441846f8SBarry Smith   ierr = PetscLogEventEnd(Tao_Solve,tao,0,0,0);CHKERRQ(ierr);
204a7e14dcfSSatish Balay 
2058931d482SJason Sarich   tao->ntotalits += tao->niter;
206fbe0838dSJason Sarich   ierr = TaoViewFromOptions(tao,NULL,"-tao_view");CHKERRQ(ierr);
207a7e14dcfSSatish Balay 
208a7e14dcfSSatish Balay   if (tao->printreason) {
209a7e14dcfSSatish Balay     if (tao->reason > 0) {
2108931d482SJason Sarich       ierr = PetscPrintf(((PetscObject)tao)->comm,"TAO solve converged due to %s iterations %D\n",TaoConvergedReasons[tao->reason],tao->niter);CHKERRQ(ierr);
211a7e14dcfSSatish Balay     } else {
2128931d482SJason Sarich       ierr = PetscPrintf(((PetscObject)tao)->comm,"TAO solve did not converge due to %s iteration %D\n",TaoConvergedReasons[tao->reason],tao->niter);CHKERRQ(ierr);
213a7e14dcfSSatish Balay     }
214a7e14dcfSSatish Balay   }
215a7e14dcfSSatish Balay   PetscFunctionReturn(0);
216a7e14dcfSSatish Balay }
217a7e14dcfSSatish Balay 
218a7e14dcfSSatish Balay #undef __FUNCT__
219a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetUp"
220a7e14dcfSSatish Balay /*@
221a7e14dcfSSatish Balay   TaoSetUp - Sets up the internal data structures for the later use
222a7e14dcfSSatish Balay   of a Tao solver
223a7e14dcfSSatish Balay 
224a7e14dcfSSatish Balay   Collective on tao
225a7e14dcfSSatish Balay 
226a7e14dcfSSatish Balay   Input Parameters:
227a7e14dcfSSatish Balay . tao - the TAO context
228a7e14dcfSSatish Balay 
229a7e14dcfSSatish Balay   Notes:
230a7e14dcfSSatish Balay   The user will not need to explicitly call TaoSetUp(), as it will
231a7e14dcfSSatish Balay   automatically be called in TaoSolve().  However, if the user
232a7e14dcfSSatish Balay   desires to call it explicitly, it should come after TaoCreate()
233a7e14dcfSSatish Balay   and any TaoSetSomething() routines, but before TaoSolve().
234a7e14dcfSSatish Balay 
235a7e14dcfSSatish Balay   Level: advanced
236a7e14dcfSSatish Balay 
237a7e14dcfSSatish Balay .seealso: TaoCreate(), TaoSolve()
238a7e14dcfSSatish Balay @*/
239441846f8SBarry Smith PetscErrorCode TaoSetUp(Tao tao)
240a7e14dcfSSatish Balay {
241a7e14dcfSSatish Balay   PetscErrorCode ierr;
24247a47007SBarry Smith 
243a7e14dcfSSatish Balay   PetscFunctionBegin;
244441846f8SBarry Smith   PetscValidHeaderSpecific(tao, TAO_CLASSID,1);
245a7e14dcfSSatish Balay   if (tao->setupcalled) PetscFunctionReturn(0);
246a7e14dcfSSatish Balay 
24747a47007SBarry Smith   if (!tao->solution) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TaoSetInitialVector");
248a7e14dcfSSatish Balay   if (tao->ops->setup) {
249a7e14dcfSSatish Balay     ierr = (*tao->ops->setup)(tao);CHKERRQ(ierr);
250a7e14dcfSSatish Balay   }
251a7e14dcfSSatish Balay   tao->setupcalled = PETSC_TRUE;
252a7e14dcfSSatish Balay   PetscFunctionReturn(0);
253a7e14dcfSSatish Balay }
254a7e14dcfSSatish Balay 
255a7e14dcfSSatish Balay #undef __FUNCT__
256a7e14dcfSSatish Balay #define __FUNCT__ "TaoDestroy"
257a7e14dcfSSatish Balay /*@
258a7e14dcfSSatish Balay   TaoDestroy - Destroys the TAO context that was created with
259a7e14dcfSSatish Balay   TaoCreate()
260a7e14dcfSSatish Balay 
261441846f8SBarry Smith   Collective on Tao
262a7e14dcfSSatish Balay 
263a7e14dcfSSatish Balay   Input Parameter:
264441846f8SBarry Smith . tao - the Tao context
265a7e14dcfSSatish Balay 
266a7e14dcfSSatish Balay   Level: beginner
267a7e14dcfSSatish Balay 
268a7e14dcfSSatish Balay .seealso: TaoCreate(), TaoSolve()
269a7e14dcfSSatish Balay @*/
270441846f8SBarry Smith PetscErrorCode TaoDestroy(Tao *tao)
271a7e14dcfSSatish Balay {
272a7e14dcfSSatish Balay   PetscErrorCode ierr;
27347a47007SBarry Smith 
274a7e14dcfSSatish Balay   PetscFunctionBegin;
275a7e14dcfSSatish Balay   if (!*tao) PetscFunctionReturn(0);
276441846f8SBarry Smith   PetscValidHeaderSpecific(*tao,TAO_CLASSID,1);
277a7e14dcfSSatish Balay   if (--((PetscObject)*tao)->refct > 0) {*tao=0;PetscFunctionReturn(0);}
278a7e14dcfSSatish Balay 
279a7e14dcfSSatish Balay   if ((*tao)->ops->destroy) {
280a7e14dcfSSatish Balay     ierr = (*((*tao))->ops->destroy)(*tao);CHKERRQ(ierr);
281a7e14dcfSSatish Balay   }
282a7e14dcfSSatish Balay   ierr = KSPDestroy(&(*tao)->ksp);CHKERRQ(ierr);
283a7e14dcfSSatish Balay   ierr = TaoLineSearchDestroy(&(*tao)->linesearch);CHKERRQ(ierr);
284a7e14dcfSSatish Balay 
285a7e14dcfSSatish Balay   if ((*tao)->ops->convergencedestroy) {
286a7e14dcfSSatish Balay     ierr = (*(*tao)->ops->convergencedestroy)((*tao)->cnvP);CHKERRQ(ierr);
287a7e14dcfSSatish Balay     if ((*tao)->jacobian_state_inv) {
288a7e14dcfSSatish Balay       ierr = MatDestroy(&(*tao)->jacobian_state_inv);CHKERRQ(ierr);
289a7e14dcfSSatish Balay     }
290a7e14dcfSSatish Balay   }
291a7e14dcfSSatish Balay   ierr = VecDestroy(&(*tao)->solution);CHKERRQ(ierr);
292a7e14dcfSSatish Balay   ierr = VecDestroy(&(*tao)->gradient);CHKERRQ(ierr);
293a7e14dcfSSatish Balay 
294a9603a14SPatrick Farrell   if ((*tao)->gradient_norm) {
295a9603a14SPatrick Farrell     ierr = PetscObjectDereference((PetscObject)(*tao)->gradient_norm);CHKERRQ(ierr);
296a9603a14SPatrick Farrell     ierr = VecDestroy(&(*tao)->gradient_norm_tmp);CHKERRQ(ierr);
297a9603a14SPatrick Farrell   }
298a9603a14SPatrick Farrell 
299a7e14dcfSSatish Balay   ierr = VecDestroy(&(*tao)->XL);CHKERRQ(ierr);
300a7e14dcfSSatish Balay   ierr = VecDestroy(&(*tao)->XU);CHKERRQ(ierr);
301a7e14dcfSSatish Balay   ierr = VecDestroy(&(*tao)->IL);CHKERRQ(ierr);
302a7e14dcfSSatish Balay   ierr = VecDestroy(&(*tao)->IU);CHKERRQ(ierr);
303a7e14dcfSSatish Balay   ierr = VecDestroy(&(*tao)->DE);CHKERRQ(ierr);
304a7e14dcfSSatish Balay   ierr = VecDestroy(&(*tao)->DI);CHKERRQ(ierr);
305a7e14dcfSSatish Balay   ierr = VecDestroy(&(*tao)->constraints_equality);CHKERRQ(ierr);
306a7e14dcfSSatish Balay   ierr = VecDestroy(&(*tao)->constraints_inequality);CHKERRQ(ierr);
307a7e14dcfSSatish Balay   ierr = VecDestroy(&(*tao)->stepdirection);CHKERRQ(ierr);
308a7e14dcfSSatish Balay   ierr = MatDestroy(&(*tao)->hessian_pre);CHKERRQ(ierr);
309a7e14dcfSSatish Balay   ierr = MatDestroy(&(*tao)->hessian);CHKERRQ(ierr);
310a7e14dcfSSatish Balay   ierr = MatDestroy(&(*tao)->jacobian_pre);CHKERRQ(ierr);
311a7e14dcfSSatish Balay   ierr = MatDestroy(&(*tao)->jacobian);CHKERRQ(ierr);
312a7e14dcfSSatish Balay   ierr = MatDestroy(&(*tao)->jacobian_state_pre);CHKERRQ(ierr);
313a7e14dcfSSatish Balay   ierr = MatDestroy(&(*tao)->jacobian_state);CHKERRQ(ierr);
314a7e14dcfSSatish Balay   ierr = MatDestroy(&(*tao)->jacobian_state_inv);CHKERRQ(ierr);
315067f330cSJason Sarich   ierr = MatDestroy(&(*tao)->jacobian_design);CHKERRQ(ierr);
316a7e14dcfSSatish Balay   ierr = MatDestroy(&(*tao)->jacobian_equality);CHKERRQ(ierr);
317a7e14dcfSSatish Balay   ierr = MatDestroy(&(*tao)->jacobian_equality_pre);CHKERRQ(ierr);
318a7e14dcfSSatish Balay   ierr = MatDestroy(&(*tao)->jacobian_inequality);CHKERRQ(ierr);
319a7e14dcfSSatish Balay   ierr = MatDestroy(&(*tao)->jacobian_inequality_pre);CHKERRQ(ierr);
320a7e14dcfSSatish Balay   ierr = ISDestroy(&(*tao)->state_is);CHKERRQ(ierr);
321a7e14dcfSSatish Balay   ierr = ISDestroy(&(*tao)->design_is);CHKERRQ(ierr);
322*125ddc8aSJason Sarich   ierr = VecDestroy(&(*tao)->sep_weights_v);CHKERRQ(ierr);
323a7e14dcfSSatish Balay   ierr = TaoCancelMonitors(*tao);CHKERRQ(ierr);
324ae93cb3cSJason Sarich   if ((*tao)->hist_malloc) {
325ae93cb3cSJason Sarich     ierr = PetscFree((*tao)->hist_obj);CHKERRQ(ierr);
326ae93cb3cSJason Sarich     ierr = PetscFree((*tao)->hist_resid);CHKERRQ(ierr);
327ae93cb3cSJason Sarich     ierr = PetscFree((*tao)->hist_cnorm);CHKERRQ(ierr);
328ae93cb3cSJason Sarich     ierr = PetscFree((*tao)->hist_lits);CHKERRQ(ierr);
329ae93cb3cSJason Sarich   }
330*125ddc8aSJason Sarich   if ((*tao)->sep_weights_n) {
331*125ddc8aSJason Sarich     ierr = PetscFree((*tao)->sep_weights_rows);CHKERRQ(ierr);
332*125ddc8aSJason Sarich     ierr = PetscFree((*tao)->sep_weights_cols);CHKERRQ(ierr);
333*125ddc8aSJason Sarich     ierr = PetscFree((*tao)->sep_weights_w);CHKERRQ(ierr);
334*125ddc8aSJason Sarich   }
335a7e14dcfSSatish Balay   ierr = PetscHeaderDestroy(tao);CHKERRQ(ierr);
336a7e14dcfSSatish Balay   PetscFunctionReturn(0);
337a7e14dcfSSatish Balay }
338a7e14dcfSSatish Balay 
339a7e14dcfSSatish Balay #undef __FUNCT__
340a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetFromOptions"
341a7e14dcfSSatish Balay /*@
342441846f8SBarry Smith   TaoSetFromOptions - Sets various Tao parameters from user
343a7e14dcfSSatish Balay   options.
344a7e14dcfSSatish Balay 
345441846f8SBarry Smith   Collective on Tao
346a7e14dcfSSatish Balay 
347a7e14dcfSSatish Balay   Input Paremeter:
348441846f8SBarry Smith . tao - the Tao solver context
349a7e14dcfSSatish Balay 
350a7e14dcfSSatish Balay   options Database Keys:
35158417fe7SBarry Smith + -tao_type <type> - The algorithm that TAO uses (lmvm, nls, etc.)
352a7e14dcfSSatish Balay . -tao_gatol <gatol> - absolute error tolerance for ||gradient||
353a7e14dcfSSatish Balay . -tao_grtol <grtol> - relative error tolerance for ||gradient||
354a7e14dcfSSatish Balay . -tao_gttol <gttol> - reduction of ||gradient|| relative to initial gradient
355a7e14dcfSSatish Balay . -tao_max_it <max> - sets maximum number of iterations
356a7e14dcfSSatish Balay . -tao_max_funcs <max> - sets maximum number of function evaluations
357a7e14dcfSSatish Balay . -tao_fmin <fmin> - stop if function value reaches fmin
358a7e14dcfSSatish Balay . -tao_steptol <tol> - stop if trust region radius less than <tol>
359a7e14dcfSSatish Balay . -tao_trust0 <t> - initial trust region radius
360a7e14dcfSSatish Balay . -tao_monitor - prints function value and residual at each iteration
361a7e14dcfSSatish Balay . -tao_smonitor - same as tao_monitor, but truncates very small values
362a7e14dcfSSatish Balay . -tao_cmonitor - prints function value, residual, and constraint norm at each iteration
363a7e14dcfSSatish Balay . -tao_view_solution - prints solution vector at each iteration
364a7e14dcfSSatish Balay . -tao_view_separableobjective - prints separable objective vector at each iteration
365a7e14dcfSSatish Balay . -tao_view_step - prints step direction vector at each iteration
366a7e14dcfSSatish Balay . -tao_view_gradient - prints gradient vector at each iteration
367a7e14dcfSSatish Balay . -tao_draw_solution - graphically view solution vector at each iteration
368a7e14dcfSSatish Balay . -tao_draw_step - graphically view step vector at each iteration
369a7e14dcfSSatish Balay . -tao_draw_gradient - graphically view gradient at each iteration
370a7e14dcfSSatish Balay . -tao_fd_gradient - use gradient computed with finite differences
371a7e14dcfSSatish Balay . -tao_cancelmonitors - cancels all monitors (except those set with command line)
372441846f8SBarry Smith . -tao_view - prints information about the Tao after solving
373a7e14dcfSSatish Balay - -tao_converged_reason - prints the reason TAO stopped iterating
374a7e14dcfSSatish Balay 
375a7e14dcfSSatish Balay   Notes:
376a7e14dcfSSatish Balay   To see all options, run your program with the -help option or consult the
377a7e14dcfSSatish Balay   user's manual. Should be called after TaoCreate() but before TaoSolve()
378a7e14dcfSSatish Balay 
379a7e14dcfSSatish Balay   Level: beginner
380a7e14dcfSSatish Balay @*/
381441846f8SBarry Smith PetscErrorCode TaoSetFromOptions(Tao tao)
382a7e14dcfSSatish Balay {
383a7e14dcfSSatish Balay   PetscErrorCode ierr;
38458417fe7SBarry Smith   const TaoType  default_type = TAOLMVM;
385a7e14dcfSSatish Balay   char           type[256], monfilename[PETSC_MAX_PATH_LEN];
386a7e14dcfSSatish Balay   PetscViewer    monviewer;
387a7e14dcfSSatish Balay   PetscBool      flg;
388a7e14dcfSSatish Balay   MPI_Comm       comm;
389a7e14dcfSSatish Balay 
390a7e14dcfSSatish Balay   PetscFunctionBegin;
391441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
392a7e14dcfSSatish Balay   ierr = PetscObjectGetComm((PetscObject)tao,&comm);CHKERRQ(ierr);
393*125ddc8aSJason Sarich 
394a7e14dcfSSatish Balay   /* So no warnings are given about unused options */
395c5929fdfSBarry Smith   ierr = PetscOptionsHasName(((PetscObject)tao)->options,((PetscObject)tao)->prefix,"-tao_ls_type",&flg);CHKERRQ(ierr);
396a7e14dcfSSatish Balay 
397a7e14dcfSSatish Balay   ierr = PetscObjectOptionsBegin((PetscObject)tao);CHKERRQ(ierr);
398a7e14dcfSSatish Balay   {
399441846f8SBarry Smith     ierr = TaoRegisterAll();CHKERRQ(ierr);
400a7e14dcfSSatish Balay     if (((PetscObject)tao)->type_name) {
401a7e14dcfSSatish Balay       default_type = ((PetscObject)tao)->type_name;
402a7e14dcfSSatish Balay     }
403a7e14dcfSSatish Balay     /* Check for type from options */
404441846f8SBarry Smith     ierr = PetscOptionsFList("-tao_type","Tao Solver type","TaoSetType",TaoList,default_type,type,256,&flg);CHKERRQ(ierr);
405a7e14dcfSSatish Balay     if (flg) {
406a7e14dcfSSatish Balay       ierr = TaoSetType(tao,type);CHKERRQ(ierr);
407a7e14dcfSSatish Balay     } else if (!((PetscObject)tao)->type_name) {
408302440fdSBarry Smith       ierr = TaoSetType(tao,default_type);CHKERRQ(ierr);
409a7e14dcfSSatish Balay     }
410a7e14dcfSSatish Balay 
4116552cf8aSJason Sarich     ierr = PetscOptionsReal("-tao_catol","Stop if constraints violations within","TaoSetConstraintTolerances",tao->catol,&tao->catol,&flg);CHKERRQ(ierr);
4126552cf8aSJason Sarich     if (flg) tao->catol_changed=PETSC_TRUE;
4136552cf8aSJason Sarich     ierr = PetscOptionsReal("-tao_crtol","Stop if relative contraint violations within","TaoSetConstraintTolerances",tao->crtol,&tao->crtol,&flg);CHKERRQ(ierr);
4146552cf8aSJason Sarich     if (flg) tao->crtol_changed=PETSC_TRUE;
4156552cf8aSJason Sarich     ierr = PetscOptionsReal("-tao_gatol","Stop if norm of gradient less than","TaoSetTolerances",tao->gatol,&tao->gatol,&flg);CHKERRQ(ierr);
4166552cf8aSJason Sarich     if (flg) tao->gatol_changed=PETSC_TRUE;
4176552cf8aSJason Sarich     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);
4186552cf8aSJason Sarich     if (flg) tao->grtol_changed=PETSC_TRUE;
4196552cf8aSJason Sarich     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);
4206552cf8aSJason Sarich     if (flg) tao->gttol_changed=PETSC_TRUE;
4216552cf8aSJason Sarich     ierr = PetscOptionsInt("-tao_max_it","Stop if iteration number exceeds","TaoSetMaximumIterations",tao->max_it,&tao->max_it,&flg);CHKERRQ(ierr);
4226552cf8aSJason Sarich     if (flg) tao->max_it_changed=PETSC_TRUE;
4236552cf8aSJason Sarich     ierr = PetscOptionsInt("-tao_max_funcs","Stop if number of function evaluations exceeds","TaoSetMaximumFunctionEvaluations",tao->max_funcs,&tao->max_funcs,&flg);CHKERRQ(ierr);
4246552cf8aSJason Sarich     if (flg) tao->max_funcs_changed=PETSC_TRUE;
4256552cf8aSJason Sarich     ierr = PetscOptionsReal("-tao_fmin","Stop if function less than","TaoSetFunctionLowerBound",tao->fmin,&tao->fmin,&flg);CHKERRQ(ierr);
4266552cf8aSJason Sarich     if (flg) tao->fmin_changed=PETSC_TRUE;
4276552cf8aSJason Sarich     ierr = PetscOptionsReal("-tao_steptol","Stop if step size or trust region radius less than","",tao->steptol,&tao->steptol,&flg);CHKERRQ(ierr);
4286552cf8aSJason Sarich     if (flg) tao->steptol_changed=PETSC_TRUE;
4296552cf8aSJason Sarich     ierr = PetscOptionsReal("-tao_trust0","Initial trust region radius","TaoSetTrustRegionRadius",tao->trust0,&tao->trust0,&flg);CHKERRQ(ierr);
4306552cf8aSJason Sarich     if (flg) tao->trust0_changed=PETSC_TRUE;
431a7e14dcfSSatish Balay     ierr = PetscOptionsString("-tao_view_solution","view solution vector after each evaluation","TaoSetMonitor","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
432a7e14dcfSSatish Balay     if (flg) {
433a7e14dcfSSatish Balay       ierr = PetscViewerASCIIOpen(comm,monfilename,&monviewer);CHKERRQ(ierr);
434a7e14dcfSSatish Balay       ierr = TaoSetMonitor(tao,TaoSolutionMonitor,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
435a7e14dcfSSatish Balay     }
436a7e14dcfSSatish Balay 
4378afaa268SBarry Smith     ierr = PetscOptionsBool("-tao_converged_reason","Print reason for TAO converged","TaoSolve",tao->printreason,&tao->printreason,NULL);CHKERRQ(ierr);
438a7e14dcfSSatish Balay     ierr = PetscOptionsString("-tao_view_gradient","view gradient vector after each evaluation","TaoSetMonitor","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
439a7e14dcfSSatish Balay     if (flg) {
440a7e14dcfSSatish Balay       ierr = PetscViewerASCIIOpen(comm,monfilename,&monviewer);CHKERRQ(ierr);
441a7e14dcfSSatish Balay       ierr = TaoSetMonitor(tao,TaoGradientMonitor,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
442a7e14dcfSSatish Balay     }
443a7e14dcfSSatish Balay 
444a7e14dcfSSatish Balay     ierr = PetscOptionsString("-tao_view_stepdirection","view step direction vector after each iteration","TaoSetMonitor","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
445a7e14dcfSSatish Balay     if (flg) {
446a7e14dcfSSatish Balay       ierr = PetscViewerASCIIOpen(comm,monfilename,&monviewer);CHKERRQ(ierr);
447a7e14dcfSSatish Balay       ierr = TaoSetMonitor(tao,TaoStepDirectionMonitor,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
448a7e14dcfSSatish Balay     }
449a7e14dcfSSatish Balay 
450a7e14dcfSSatish Balay     ierr = PetscOptionsString("-tao_view_separableobjective","view separable objective vector after each evaluation","TaoSetMonitor","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
451a7e14dcfSSatish Balay     if (flg) {
452a7e14dcfSSatish Balay       ierr = PetscViewerASCIIOpen(comm,monfilename,&monviewer);CHKERRQ(ierr);
453a7e14dcfSSatish Balay       ierr = TaoSetMonitor(tao,TaoSeparableObjectiveMonitor,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
454a7e14dcfSSatish Balay     }
455a7e14dcfSSatish Balay 
456a7e14dcfSSatish Balay     ierr = PetscOptionsString("-tao_monitor","Use the default convergence monitor","TaoSetMonitor","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
457a7e14dcfSSatish Balay     if (flg) {
458a7e14dcfSSatish Balay       ierr = PetscViewerASCIIOpen(comm,monfilename,&monviewer);CHKERRQ(ierr);
459a7e14dcfSSatish Balay       ierr = TaoSetMonitor(tao,TaoDefaultMonitor,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
460a7e14dcfSSatish Balay     }
461a7e14dcfSSatish Balay 
462a7e14dcfSSatish Balay     ierr = PetscOptionsString("-tao_smonitor","Use the short convergence monitor","TaoSetMonitor","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
463a7e14dcfSSatish Balay     if (flg) {
464a7e14dcfSSatish Balay       ierr = PetscViewerASCIIOpen(comm,monfilename,&monviewer);CHKERRQ(ierr);
465a7e14dcfSSatish Balay       ierr = TaoSetMonitor(tao,TaoDefaultSMonitor,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
466a7e14dcfSSatish Balay     }
467a7e14dcfSSatish Balay 
468a7e14dcfSSatish Balay     ierr = PetscOptionsString("-tao_cmonitor","Use the default convergence monitor with constraint norm","TaoSetMonitor","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
469a7e14dcfSSatish Balay     if (flg) {
470a7e14dcfSSatish Balay       ierr = PetscViewerASCIIOpen(comm,monfilename,&monviewer);CHKERRQ(ierr);
471a7e14dcfSSatish Balay       ierr = TaoSetMonitor(tao,TaoDefaultCMonitor,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
472a7e14dcfSSatish Balay     }
473a7e14dcfSSatish Balay 
474a7e14dcfSSatish Balay 
4758afaa268SBarry Smith     flg = PETSC_FALSE;
4768afaa268SBarry Smith     ierr = PetscOptionsBool("-tao_cancelmonitors","cancel all monitors and call any registered destroy routines","TaoCancelMonitors",flg,&flg,NULL);CHKERRQ(ierr);
477a7e14dcfSSatish Balay     if (flg) {ierr = TaoCancelMonitors(tao);CHKERRQ(ierr);}
478a7e14dcfSSatish Balay 
4798afaa268SBarry Smith     flg = PETSC_FALSE;
4808afaa268SBarry Smith     ierr = PetscOptionsBool("-tao_draw_solution","Plot solution vector at each iteration","TaoSetMonitor",flg,&flg,NULL);CHKERRQ(ierr);
481a7e14dcfSSatish Balay     if (flg) {
48247a47007SBarry Smith       ierr = TaoSetMonitor(tao,TaoDrawSolutionMonitor,NULL,NULL);CHKERRQ(ierr);
483a7e14dcfSSatish Balay     }
484a7e14dcfSSatish Balay 
4858afaa268SBarry Smith     flg = PETSC_FALSE;
4868afaa268SBarry Smith     ierr = PetscOptionsBool("-tao_draw_step","plots step direction at each iteration","TaoSetMonitor",flg,&flg,NULL);CHKERRQ(ierr);
487a7e14dcfSSatish Balay     if (flg) {
48847a47007SBarry Smith       ierr = TaoSetMonitor(tao,TaoDrawStepMonitor,NULL,NULL);CHKERRQ(ierr);
489a7e14dcfSSatish Balay     }
490a7e14dcfSSatish Balay 
4918afaa268SBarry Smith     flg = PETSC_FALSE;
4928afaa268SBarry Smith     ierr = PetscOptionsBool("-tao_draw_gradient","plots gradient at each iteration","TaoSetMonitor",flg,&flg,NULL);CHKERRQ(ierr);
493a7e14dcfSSatish Balay     if (flg) {
49447a47007SBarry Smith       ierr = TaoSetMonitor(tao,TaoDrawGradientMonitor,NULL,NULL);CHKERRQ(ierr);
495a7e14dcfSSatish Balay     }
4968afaa268SBarry Smith     flg = PETSC_FALSE;
4978afaa268SBarry Smith     ierr = PetscOptionsBool("-tao_fd_gradient","compute gradient using finite differences","TaoDefaultComputeGradient",flg,&flg,NULL);CHKERRQ(ierr);
498a7e14dcfSSatish Balay     if (flg) {
49947a47007SBarry Smith       ierr = TaoSetGradientRoutine(tao,TaoDefaultComputeGradient,NULL);CHKERRQ(ierr);
500a7e14dcfSSatish Balay     }
501e1cc520bSBarry Smith     ierr = PetscOptionsEnum("-tao_subset_type","subset type", "", TaoSubSetTypes,(PetscEnum)tao->subset_type, (PetscEnum*)&tao->subset_type, 0);CHKERRQ(ierr);
502a7e14dcfSSatish Balay 
503a7e14dcfSSatish Balay     if (tao->ops->setfromoptions) {
5041a1499c8SBarry Smith       ierr = (*tao->ops->setfromoptions)(PetscOptionsObject,tao);CHKERRQ(ierr);
505a7e14dcfSSatish Balay     }
506a7e14dcfSSatish Balay   }
507a7e14dcfSSatish Balay   ierr = PetscOptionsEnd();CHKERRQ(ierr);
508a7e14dcfSSatish Balay   PetscFunctionReturn(0);
509a7e14dcfSSatish Balay }
510a7e14dcfSSatish Balay 
511a7e14dcfSSatish Balay #undef __FUNCT__
512a7e14dcfSSatish Balay #define __FUNCT__ "TaoView"
513a7e14dcfSSatish Balay /*@C
514441846f8SBarry Smith   TaoView - Prints information about the Tao
515a7e14dcfSSatish Balay 
516441846f8SBarry Smith   Collective on Tao
517a7e14dcfSSatish Balay 
518a7e14dcfSSatish Balay   InputParameters:
519441846f8SBarry Smith + tao - the Tao context
520a7e14dcfSSatish Balay - viewer - visualization context
521a7e14dcfSSatish Balay 
522a7e14dcfSSatish Balay   Options Database Key:
523a7e14dcfSSatish Balay . -tao_view - Calls TaoView() at the end of TaoSolve()
524a7e14dcfSSatish Balay 
525a7e14dcfSSatish Balay   Notes:
526a7e14dcfSSatish Balay   The available visualization contexts include
527a7e14dcfSSatish Balay +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
528a7e14dcfSSatish Balay -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
529a7e14dcfSSatish Balay          output where only the first processor opens
530a7e14dcfSSatish Balay          the file.  All other processors send their
531a7e14dcfSSatish Balay          data to the first processor to print.
532a7e14dcfSSatish Balay 
533a7e14dcfSSatish Balay   Level: beginner
534a7e14dcfSSatish Balay 
535a7e14dcfSSatish Balay .seealso: PetscViewerASCIIOpen()
536a7e14dcfSSatish Balay @*/
537441846f8SBarry Smith PetscErrorCode TaoView(Tao tao, PetscViewer viewer)
538a7e14dcfSSatish Balay {
539a7e14dcfSSatish Balay   PetscErrorCode      ierr;
540a7e14dcfSSatish Balay   PetscBool           isascii,isstring;
541441846f8SBarry Smith   const TaoType type;
54247a47007SBarry Smith 
543a7e14dcfSSatish Balay   PetscFunctionBegin;
544441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
545a7e14dcfSSatish Balay   if (!viewer) {
546a7e14dcfSSatish Balay     ierr = PetscViewerASCIIGetStdout(((PetscObject)tao)->comm,&viewer);CHKERRQ(ierr);
547a7e14dcfSSatish Balay   }
548a7e14dcfSSatish Balay   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
549a7e14dcfSSatish Balay   PetscCheckSameComm(tao,1,viewer,2);
550a7e14dcfSSatish Balay 
551a7e14dcfSSatish Balay   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
552a7e14dcfSSatish Balay   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
553a7e14dcfSSatish Balay   if (isascii) {
554a7e14dcfSSatish Balay     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)tao,viewer);CHKERRQ(ierr);
555a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
556a7e14dcfSSatish Balay 
557a7e14dcfSSatish Balay     if (tao->ops->view) {
558a7e14dcfSSatish Balay       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
559a7e14dcfSSatish Balay       ierr = (*tao->ops->view)(tao,viewer);CHKERRQ(ierr);
560a7e14dcfSSatish Balay       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
561a7e14dcfSSatish Balay     }
562a7e14dcfSSatish Balay     if (tao->linesearch) {
563a7e14dcfSSatish Balay       ierr = PetscObjectPrintClassNamePrefixType((PetscObject)(tao->linesearch),viewer);CHKERRQ(ierr);
564a7e14dcfSSatish Balay     }
565a7e14dcfSSatish Balay     if (tao->ksp) {
566a7e14dcfSSatish Balay       ierr = PetscObjectPrintClassNamePrefixType((PetscObject)(tao->ksp),viewer);CHKERRQ(ierr);
567ae93cb3cSJason Sarich       ierr = PetscViewerASCIIPrintf(viewer,"total KSP iterations: %D\n",tao->ksp_tot_its);CHKERRQ(ierr);
568a7e14dcfSSatish Balay     }
569a7e14dcfSSatish Balay     if (tao->XL || tao->XU) {
570302440fdSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"Active Set subset type: %s\n",TaoSubSetTypes[tao->subset_type]);CHKERRQ(ierr);
571a7e14dcfSSatish Balay     }
572a7e14dcfSSatish Balay 
57347a47007SBarry Smith     ierr=PetscViewerASCIIPrintf(viewer,"convergence tolerances: gatol=%g,",(double)tao->gatol);CHKERRQ(ierr);
57447a47007SBarry Smith     ierr=PetscViewerASCIIPrintf(viewer," steptol=%g,",(double)tao->steptol);CHKERRQ(ierr);
57547a47007SBarry Smith     ierr=PetscViewerASCIIPrintf(viewer," gttol=%g\n",(double)tao->gttol);CHKERRQ(ierr);
576a7e14dcfSSatish Balay 
57747a47007SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"Residual in Function/Gradient:=%g\n",(double)tao->residual);CHKERRQ(ierr);
578a7e14dcfSSatish Balay 
579a7e14dcfSSatish Balay     if (tao->cnorm>0 || tao->catol>0 || tao->crtol>0){
580a7e14dcfSSatish Balay       ierr=PetscViewerASCIIPrintf(viewer,"convergence tolerances:");CHKERRQ(ierr);
58147a47007SBarry Smith       ierr=PetscViewerASCIIPrintf(viewer," catol=%g,",(double)tao->catol);CHKERRQ(ierr);
58247a47007SBarry Smith       ierr=PetscViewerASCIIPrintf(viewer," crtol=%g\n",(double)tao->crtol);CHKERRQ(ierr);
58347a47007SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"Residual in Constraints:=%g\n",(double)tao->cnorm);CHKERRQ(ierr);
584a7e14dcfSSatish Balay     }
585a7e14dcfSSatish Balay 
586a7e14dcfSSatish Balay     if (tao->trust < tao->steptol){
58747a47007SBarry Smith       ierr=PetscViewerASCIIPrintf(viewer,"convergence tolerances: steptol=%g\n",(double)tao->steptol);CHKERRQ(ierr);
58847a47007SBarry Smith       ierr=PetscViewerASCIIPrintf(viewer,"Final trust region radius:=%g\n",(double)tao->trust);CHKERRQ(ierr);
589a7e14dcfSSatish Balay     }
590a7e14dcfSSatish Balay 
591a7e14dcfSSatish Balay     if (tao->fmin>-1.e25){
59247a47007SBarry Smith       ierr=PetscViewerASCIIPrintf(viewer,"convergence tolerances: function minimum=%g\n",(double)tao->fmin);CHKERRQ(ierr);
593a7e14dcfSSatish Balay     }
59447a47007SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"Objective value=%g\n",(double)tao->fc);CHKERRQ(ierr);
595a7e14dcfSSatish Balay 
59647a47007SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"total number of iterations=%D,          ",tao->niter);CHKERRQ(ierr);
597a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"              (max: %D)\n",tao->max_it);CHKERRQ(ierr);
598a7e14dcfSSatish Balay 
599a7e14dcfSSatish Balay     if (tao->nfuncs>0){
60047a47007SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"total number of function evaluations=%D,",tao->nfuncs);CHKERRQ(ierr);
60147a47007SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"                max: %D\n",tao->max_funcs);CHKERRQ(ierr);
602a7e14dcfSSatish Balay     }
603a7e14dcfSSatish Balay     if (tao->ngrads>0){
60447a47007SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"total number of gradient evaluations=%D,",tao->ngrads);CHKERRQ(ierr);
60547a47007SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"                max: %D\n",tao->max_funcs);CHKERRQ(ierr);
606a7e14dcfSSatish Balay     }
607a7e14dcfSSatish Balay     if (tao->nfuncgrads>0){
60847a47007SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"total number of function/gradient evaluations=%D,",tao->nfuncgrads);CHKERRQ(ierr);
60947a47007SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    (max: %D)\n",tao->max_funcs);CHKERRQ(ierr);
610a7e14dcfSSatish Balay     }
611a7e14dcfSSatish Balay     if (tao->nhess>0){
61247a47007SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"total number of Hessian evaluations=%D\n",tao->nhess);CHKERRQ(ierr);
613a7e14dcfSSatish Balay     }
614a7e14dcfSSatish Balay     /*  if (tao->linear_its>0){
61547a47007SBarry Smith      ierr = PetscViewerASCIIPrintf(viewer,"  total Krylov method iterations=%D\n",tao->linear_its);CHKERRQ(ierr);
616a7e14dcfSSatish Balay      }*/
617a7e14dcfSSatish Balay     if (tao->nconstraints>0){
61847a47007SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"total number of constraint function evaluations=%D\n",tao->nconstraints);CHKERRQ(ierr);
619a7e14dcfSSatish Balay     }
620a7e14dcfSSatish Balay     if (tao->njac>0){
62147a47007SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"total number of Jacobian evaluations=%D\n",tao->njac);CHKERRQ(ierr);
622a7e14dcfSSatish Balay     }
623a7e14dcfSSatish Balay 
624a7e14dcfSSatish Balay     if (tao->reason>0){
625a7e14dcfSSatish Balay       ierr = PetscViewerASCIIPrintf(viewer,    "Solution converged: ");CHKERRQ(ierr);
626a7e14dcfSSatish Balay       switch (tao->reason) {
627a7e14dcfSSatish Balay       case TAO_CONVERGED_GATOL:
628a7e14dcfSSatish Balay         ierr = PetscViewerASCIIPrintf(viewer," ||g(X)|| <= gatol\n");CHKERRQ(ierr);
629a7e14dcfSSatish Balay         break;
630a7e14dcfSSatish Balay       case TAO_CONVERGED_GRTOL:
631a7e14dcfSSatish Balay         ierr = PetscViewerASCIIPrintf(viewer," ||g(X)||/|f(X)| <= grtol\n");CHKERRQ(ierr);
632a7e14dcfSSatish Balay         break;
633a7e14dcfSSatish Balay       case TAO_CONVERGED_GTTOL:
634a7e14dcfSSatish Balay         ierr = PetscViewerASCIIPrintf(viewer," ||g(X)||/||g(X0)|| <= gttol\n");CHKERRQ(ierr);
635a7e14dcfSSatish Balay         break;
636a7e14dcfSSatish Balay       case TAO_CONVERGED_STEPTOL:
637a7e14dcfSSatish Balay         ierr = PetscViewerASCIIPrintf(viewer," Steptol -- step size small\n");CHKERRQ(ierr);
638a7e14dcfSSatish Balay         break;
639a7e14dcfSSatish Balay       case TAO_CONVERGED_MINF:
640a7e14dcfSSatish Balay         ierr = PetscViewerASCIIPrintf(viewer," Minf --  f < fmin\n");CHKERRQ(ierr);
641a7e14dcfSSatish Balay         break;
642a7e14dcfSSatish Balay       case TAO_CONVERGED_USER:
643a7e14dcfSSatish Balay         ierr = PetscViewerASCIIPrintf(viewer," User Terminated\n");CHKERRQ(ierr);
644a7e14dcfSSatish Balay         break;
645a7e14dcfSSatish Balay       default:
646a7e14dcfSSatish Balay         ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
647a7e14dcfSSatish Balay         break;
648a7e14dcfSSatish Balay       }
649a7e14dcfSSatish Balay 
650a7e14dcfSSatish Balay     } else {
65147cfc4faSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"Solver terminated: %d",tao->reason);CHKERRQ(ierr);
652a7e14dcfSSatish Balay       switch (tao->reason) {
653a7e14dcfSSatish Balay       case TAO_DIVERGED_MAXITS:
65447a47007SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer," Maximum Iterations\n");CHKERRQ(ierr);
655a7e14dcfSSatish Balay         break;
656a7e14dcfSSatish Balay       case TAO_DIVERGED_NAN:
65747a47007SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer," NAN or Inf encountered\n");CHKERRQ(ierr);
658a7e14dcfSSatish Balay         break;
659a7e14dcfSSatish Balay       case TAO_DIVERGED_MAXFCN:
66047a47007SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer," Maximum Function Evaluations\n");CHKERRQ(ierr);
661a7e14dcfSSatish Balay         break;
662a7e14dcfSSatish Balay       case TAO_DIVERGED_LS_FAILURE:
66347a47007SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer," Line Search Failure\n");CHKERRQ(ierr);
664a7e14dcfSSatish Balay         break;
665a7e14dcfSSatish Balay       case TAO_DIVERGED_TR_REDUCTION:
66647a47007SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer," Trust Region too small\n");CHKERRQ(ierr);
667a7e14dcfSSatish Balay         break;
668a7e14dcfSSatish Balay       case TAO_DIVERGED_USER:
66947a47007SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer," User Terminated\n");CHKERRQ(ierr);
670a7e14dcfSSatish Balay         break;
671a7e14dcfSSatish Balay       default:
672a7e14dcfSSatish Balay         ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
673a7e14dcfSSatish Balay         break;
674a7e14dcfSSatish Balay       }
675a7e14dcfSSatish Balay     }
676a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
677a7e14dcfSSatish Balay   } else if (isstring) {
678a7e14dcfSSatish Balay     ierr = TaoGetType(tao,&type);CHKERRQ(ierr);
679a7e14dcfSSatish Balay     ierr = PetscViewerStringSPrintf(viewer," %-3.3s",type);CHKERRQ(ierr);
680a7e14dcfSSatish Balay   }
681a7e14dcfSSatish Balay   PetscFunctionReturn(0);
682a7e14dcfSSatish Balay }
683a7e14dcfSSatish Balay 
684a7e14dcfSSatish Balay #undef __FUNCT__
685a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetTolerances"
686a7e14dcfSSatish Balay /*@
687a7e14dcfSSatish Balay   TaoSetTolerances - Sets parameters used in TAO convergence tests
688a7e14dcfSSatish Balay 
689441846f8SBarry Smith   Logically collective on Tao
690a7e14dcfSSatish Balay 
691a7e14dcfSSatish Balay   Input Parameters:
692441846f8SBarry Smith + tao - the Tao context
693a7e14dcfSSatish Balay . gatol - stop if norm of gradient is less than this
694a7e14dcfSSatish Balay . grtol - stop if relative norm of gradient is less than this
695a7e14dcfSSatish Balay - gttol - stop if norm of gradient is reduced by this factor
696a7e14dcfSSatish Balay 
697a7e14dcfSSatish Balay   Options Database Keys:
698e52336cbSBarry Smith + -tao_gatol <gatol> - Sets gatol
699a7e14dcfSSatish Balay . -tao_grtol <grtol> - Sets grtol
700a7e14dcfSSatish Balay - -tao_gttol <gttol> - Sets gttol
701a7e14dcfSSatish Balay 
702a7e14dcfSSatish Balay   Stopping Criteria:
703a7e14dcfSSatish Balay $ ||g(X)||                            <= gatol
704a7e14dcfSSatish Balay $ ||g(X)|| / |f(X)|                   <= grtol
705a7e14dcfSSatish Balay $ ||g(X)|| / ||g(X0)||                <= gttol
706a7e14dcfSSatish Balay 
707a7e14dcfSSatish Balay   Notes:
708a7e14dcfSSatish Balay   Use PETSC_DEFAULT to leave one or more tolerances unchanged.
709a7e14dcfSSatish Balay 
710a7e14dcfSSatish Balay   Level: beginner
711a7e14dcfSSatish Balay 
712a7e14dcfSSatish Balay .seealso: TaoGetTolerances()
713a7e14dcfSSatish Balay 
714a7e14dcfSSatish Balay @*/
715e52336cbSBarry Smith PetscErrorCode TaoSetTolerances(Tao tao, PetscReal gatol, PetscReal grtol, PetscReal gttol)
716a7e14dcfSSatish Balay {
717a7e14dcfSSatish Balay   PetscErrorCode ierr;
71847a47007SBarry Smith 
719a7e14dcfSSatish Balay   PetscFunctionBegin;
720441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
721a7e14dcfSSatish Balay 
722a7e14dcfSSatish Balay   if (gatol != PETSC_DEFAULT) {
723a7e14dcfSSatish Balay     if (gatol<0) {
724955c1f14SBarry Smith       ierr = PetscInfo(tao,"Tried to set negative gatol -- ignored.\n");CHKERRQ(ierr);
725a7e14dcfSSatish Balay     } else {
726a7e14dcfSSatish Balay       tao->gatol = PetscMax(0,gatol);
7276552cf8aSJason Sarich       tao->gatol_changed=PETSC_TRUE;
728a7e14dcfSSatish Balay     }
729a7e14dcfSSatish Balay   }
730a7e14dcfSSatish Balay 
731a7e14dcfSSatish Balay   if (grtol != PETSC_DEFAULT) {
732a7e14dcfSSatish Balay     if (grtol<0) {
733955c1f14SBarry Smith       ierr = PetscInfo(tao,"Tried to set negative grtol -- ignored.\n");CHKERRQ(ierr);
734a7e14dcfSSatish Balay     } else {
735a7e14dcfSSatish Balay       tao->grtol = PetscMax(0,grtol);
7366552cf8aSJason Sarich       tao->grtol_changed=PETSC_TRUE;
737a7e14dcfSSatish Balay     }
738a7e14dcfSSatish Balay   }
739a7e14dcfSSatish Balay 
740a7e14dcfSSatish Balay   if (gttol != PETSC_DEFAULT) {
741a7e14dcfSSatish Balay     if (gttol<0) {
742955c1f14SBarry Smith       ierr = PetscInfo(tao,"Tried to set negative gttol -- ignored.\n");CHKERRQ(ierr);
743a7e14dcfSSatish Balay     } else {
744a7e14dcfSSatish Balay       tao->gttol = PetscMax(0,gttol);
7456552cf8aSJason Sarich       tao->gttol_changed=PETSC_TRUE;
746a7e14dcfSSatish Balay     }
747a7e14dcfSSatish Balay   }
748a7e14dcfSSatish Balay   PetscFunctionReturn(0);
749a7e14dcfSSatish Balay }
750a7e14dcfSSatish Balay 
751a7e14dcfSSatish Balay #undef __FUNCT__
752a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetConstraintTolerances"
753a7e14dcfSSatish Balay /*@
75458417fe7SBarry Smith   TaoSetConstraintTolerances - Sets constraint tolerance parameters used in TAO  convergence tests
755a7e14dcfSSatish Balay 
756441846f8SBarry Smith   Logically collective on Tao
757a7e14dcfSSatish Balay 
758a7e14dcfSSatish Balay   Input Parameters:
759441846f8SBarry Smith + tao - the Tao context
760e52336cbSBarry Smith . catol - absolute constraint tolerance, constraint norm must be less than catol for used for gatol convergence criteria
761e52336cbSBarry Smith - crtol - relative contraint tolerance, constraint norm must be less than crtol for used for gatol, gttol convergence criteria
762a7e14dcfSSatish Balay 
763a7e14dcfSSatish Balay   Options Database Keys:
764a7e14dcfSSatish Balay + -tao_catol <catol> - Sets catol
765a7e14dcfSSatish Balay - -tao_crtol <crtol> - Sets crtol
766a7e14dcfSSatish Balay 
7676552cf8aSJason Sarich   Notes:
7686552cf8aSJason Sarich   Use PETSC_DEFAULT to leave any tolerance unchanged.
7696552cf8aSJason Sarich 
770a7e14dcfSSatish Balay   Level: intermediate
771a7e14dcfSSatish Balay 
77258417fe7SBarry Smith .seealso: TaoGetTolerances(), TaoGetConstraintTolerances(), TaoSetTolerances()
773a7e14dcfSSatish Balay 
774a7e14dcfSSatish Balay @*/
775441846f8SBarry Smith PetscErrorCode TaoSetConstraintTolerances(Tao tao, PetscReal catol, PetscReal crtol)
776a7e14dcfSSatish Balay {
777a7e14dcfSSatish Balay   PetscErrorCode ierr;
77847a47007SBarry Smith 
779a7e14dcfSSatish Balay   PetscFunctionBegin;
780441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
781a7e14dcfSSatish Balay 
782a7e14dcfSSatish Balay   if (catol != PETSC_DEFAULT) {
783a7e14dcfSSatish Balay     if (catol<0) {
784955c1f14SBarry Smith       ierr = PetscInfo(tao,"Tried to set negative catol -- ignored.\n");CHKERRQ(ierr);
785a7e14dcfSSatish Balay     } else {
786a7e14dcfSSatish Balay       tao->catol = PetscMax(0,catol);
7876552cf8aSJason Sarich       tao->catol_changed=PETSC_TRUE;
788a7e14dcfSSatish Balay     }
789a7e14dcfSSatish Balay   }
790a7e14dcfSSatish Balay 
791a7e14dcfSSatish Balay   if (crtol != PETSC_DEFAULT) {
792a7e14dcfSSatish Balay     if (crtol<0) {
793955c1f14SBarry Smith       ierr = PetscInfo(tao,"Tried to set negative crtol -- ignored.\n");CHKERRQ(ierr);
794a7e14dcfSSatish Balay     } else {
795a7e14dcfSSatish Balay       tao->crtol = PetscMax(0,crtol);
7966552cf8aSJason Sarich       tao->crtol_changed=PETSC_TRUE;
797a7e14dcfSSatish Balay     }
798a7e14dcfSSatish Balay   }
799a7e14dcfSSatish Balay   PetscFunctionReturn(0);
800a7e14dcfSSatish Balay }
801a7e14dcfSSatish Balay 
802a7e14dcfSSatish Balay #undef __FUNCT__
80358417fe7SBarry Smith #define __FUNCT__ "TaoGetConstraintTolerances"
80458417fe7SBarry Smith /*@
80558417fe7SBarry Smith   TaoGetConstraintTolerances - Gets constraint tolerance parameters used in TAO  convergence tests
80658417fe7SBarry Smith 
80758417fe7SBarry Smith   Not ollective
80858417fe7SBarry Smith 
80958417fe7SBarry Smith   Input Parameter:
81058417fe7SBarry Smith . tao - the Tao context
81158417fe7SBarry Smith 
81258417fe7SBarry Smith   Output Parameter:
813e52336cbSBarry Smith + catol - absolute constraint tolerance, constraint norm must be less than catol for used for gatol convergence criteria
814e52336cbSBarry Smith - crtol - relative contraint tolerance, constraint norm must be less than crtol for used for gatol, gttol convergence criteria
81558417fe7SBarry Smith 
81658417fe7SBarry Smith   Level: intermediate
81758417fe7SBarry Smith 
81858417fe7SBarry Smith .seealso: TaoGetTolerances(), TaoSetTolerances(), TaoSetConstraintTolerances()
81958417fe7SBarry Smith 
82058417fe7SBarry Smith @*/
82158417fe7SBarry Smith PetscErrorCode TaoGetConstraintTolerances(Tao tao, PetscReal *catol, PetscReal *crtol)
82258417fe7SBarry Smith {
82358417fe7SBarry Smith   PetscFunctionBegin;
82458417fe7SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
82558417fe7SBarry Smith   if (catol) *catol = tao->catol;
82658417fe7SBarry Smith   if (crtol) *crtol = tao->crtol;
82758417fe7SBarry Smith   PetscFunctionReturn(0);
82858417fe7SBarry Smith }
82958417fe7SBarry Smith 
83058417fe7SBarry Smith #undef __FUNCT__
831a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetFunctionLowerBound"
832a7e14dcfSSatish Balay /*@
833a7e14dcfSSatish Balay    TaoSetFunctionLowerBound - Sets a bound on the solution objective value.
834a7e14dcfSSatish Balay    When an approximate solution with an objective value below this number
835a7e14dcfSSatish Balay    has been found, the solver will terminate.
836a7e14dcfSSatish Balay 
837441846f8SBarry Smith    Logically Collective on Tao
838a7e14dcfSSatish Balay 
839a7e14dcfSSatish Balay    Input Parameters:
840441846f8SBarry Smith +  tao - the Tao solver context
841a7e14dcfSSatish Balay -  fmin - the tolerance
842a7e14dcfSSatish Balay 
843a7e14dcfSSatish Balay    Options Database Keys:
844a7e14dcfSSatish Balay .    -tao_fmin <fmin> - sets the minimum function value
845a7e14dcfSSatish Balay 
846a7e14dcfSSatish Balay    Level: intermediate
847a7e14dcfSSatish Balay 
848a7e14dcfSSatish Balay .seealso: TaoSetTolerances()
849a7e14dcfSSatish Balay @*/
850441846f8SBarry Smith PetscErrorCode TaoSetFunctionLowerBound(Tao tao,PetscReal fmin)
851a7e14dcfSSatish Balay {
852a7e14dcfSSatish Balay   PetscFunctionBegin;
853441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
854a7e14dcfSSatish Balay   tao->fmin = fmin;
8556552cf8aSJason Sarich   tao->fmin_changed=PETSC_TRUE;
856a7e14dcfSSatish Balay   PetscFunctionReturn(0);
857a7e14dcfSSatish Balay }
858a7e14dcfSSatish Balay 
859a7e14dcfSSatish Balay #undef __FUNCT__
860a7e14dcfSSatish Balay #define __FUNCT__ "TaoGetFunctionLowerBound"
861a7e14dcfSSatish Balay /*@
8628e96d397SJason Sarich    TaoGetFunctionLowerBound - Gets the bound on the solution objective value.
863a7e14dcfSSatish Balay    When an approximate solution with an objective value below this number
864a7e14dcfSSatish Balay    has been found, the solver will terminate.
865a7e14dcfSSatish Balay 
866441846f8SBarry Smith    Not collective on Tao
867a7e14dcfSSatish Balay 
868a7e14dcfSSatish Balay    Input Parameters:
869441846f8SBarry Smith .  tao - the Tao solver context
870a7e14dcfSSatish Balay 
871a7e14dcfSSatish Balay    OutputParameters:
872a7e14dcfSSatish Balay .  fmin - the minimum function value
873a7e14dcfSSatish Balay 
874a7e14dcfSSatish Balay    Level: intermediate
875a7e14dcfSSatish Balay 
876a7e14dcfSSatish Balay .seealso: TaoSetFunctionLowerBound()
877a7e14dcfSSatish Balay @*/
878441846f8SBarry Smith PetscErrorCode TaoGetFunctionLowerBound(Tao tao,PetscReal *fmin)
879a7e14dcfSSatish Balay {
880a7e14dcfSSatish Balay   PetscFunctionBegin;
881441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
882a7e14dcfSSatish Balay   *fmin = tao->fmin;
883a7e14dcfSSatish Balay   PetscFunctionReturn(0);
884a7e14dcfSSatish Balay }
885a7e14dcfSSatish Balay 
886a7e14dcfSSatish Balay #undef __FUNCT__
887a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetMaximumFunctionEvaluations"
888a7e14dcfSSatish Balay /*@
889a7e14dcfSSatish Balay    TaoSetMaximumFunctionEvaluations - Sets a maximum number of
890a7e14dcfSSatish Balay    function evaluations.
891a7e14dcfSSatish Balay 
892441846f8SBarry Smith    Logically Collective on Tao
893a7e14dcfSSatish Balay 
894a7e14dcfSSatish Balay    Input Parameters:
895441846f8SBarry Smith +  tao - the Tao solver context
896a7e14dcfSSatish Balay -  nfcn - the maximum number of function evaluations (>=0)
897a7e14dcfSSatish Balay 
898a7e14dcfSSatish Balay    Options Database Keys:
899a7e14dcfSSatish Balay .    -tao_max_funcs <nfcn> - sets the maximum number of function evaluations
900a7e14dcfSSatish Balay 
901a7e14dcfSSatish Balay    Level: intermediate
902a7e14dcfSSatish Balay 
903a7e14dcfSSatish Balay .seealso: TaoSetTolerances(), TaoSetMaximumIterations()
904a7e14dcfSSatish Balay @*/
905a7e14dcfSSatish Balay 
906441846f8SBarry Smith PetscErrorCode TaoSetMaximumFunctionEvaluations(Tao tao,PetscInt nfcn)
907a7e14dcfSSatish Balay {
908a7e14dcfSSatish Balay   PetscFunctionBegin;
909441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
91047a47007SBarry Smith   tao->max_funcs = PetscMax(0,nfcn);
9116552cf8aSJason Sarich   tao->max_funcs_changed=PETSC_TRUE;
912a7e14dcfSSatish Balay   PetscFunctionReturn(0);
913a7e14dcfSSatish Balay }
914a7e14dcfSSatish Balay 
915a7e14dcfSSatish Balay #undef __FUNCT__
916a7e14dcfSSatish Balay #define __FUNCT__ "TaoGetMaximumFunctionEvaluations"
917a7e14dcfSSatish Balay /*@
918a7e14dcfSSatish Balay    TaoGetMaximumFunctionEvaluations - Sets a maximum number of
919a7e14dcfSSatish Balay    function evaluations.
920a7e14dcfSSatish Balay 
921a7e14dcfSSatish Balay    Not Collective
922a7e14dcfSSatish Balay 
923a7e14dcfSSatish Balay    Input Parameters:
924441846f8SBarry Smith .  tao - the Tao solver context
925a7e14dcfSSatish Balay 
926a7e14dcfSSatish Balay    Output Parameters:
927a7e14dcfSSatish Balay .  nfcn - the maximum number of function evaluations
928a7e14dcfSSatish Balay 
929a7e14dcfSSatish Balay    Level: intermediate
930a7e14dcfSSatish Balay 
931a7e14dcfSSatish Balay .seealso: TaoSetMaximumFunctionEvaluations(), TaoGetMaximumIterations()
932a7e14dcfSSatish Balay @*/
933a7e14dcfSSatish Balay 
934441846f8SBarry Smith PetscErrorCode TaoGetMaximumFunctionEvaluations(Tao tao,PetscInt *nfcn)
935a7e14dcfSSatish Balay {
936a7e14dcfSSatish Balay   PetscFunctionBegin;
937441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
938a7e14dcfSSatish Balay   *nfcn = tao->max_funcs;
939a7e14dcfSSatish Balay   PetscFunctionReturn(0);
940a7e14dcfSSatish Balay }
941a7e14dcfSSatish Balay 
942a7e14dcfSSatish Balay #undef __FUNCT__
943770232b9SCe Qin #define __FUNCT__ "TaoGetCurrentFunctionEvaluations"
944770232b9SCe Qin /*@
945770232b9SCe Qin    TaoGetCurrentFunctionEvaluations - Get current number of
946770232b9SCe Qin    function evaluations.
947770232b9SCe Qin 
948770232b9SCe Qin    Not Collective
949770232b9SCe Qin 
950770232b9SCe Qin    Input Parameters:
951770232b9SCe Qin .  tao - the Tao solver context
952770232b9SCe Qin 
953770232b9SCe Qin    Output Parameters:
954770232b9SCe Qin .  nfuncs - the current number of function evaluations
955770232b9SCe Qin 
956770232b9SCe Qin    Level: intermediate
957770232b9SCe Qin 
958770232b9SCe Qin .seealso: TaoSetMaximumFunctionEvaluations(), TaoGetMaximumFunctionEvaluations(), TaoGetMaximumIterations()
959770232b9SCe Qin @*/
960770232b9SCe Qin 
961770232b9SCe Qin PetscErrorCode TaoGetCurrentFunctionEvaluations(Tao tao,PetscInt *nfuncs)
962770232b9SCe Qin {
963770232b9SCe Qin   PetscFunctionBegin;
964770232b9SCe Qin   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
965770232b9SCe Qin   *nfuncs=PetscMax(tao->nfuncs,tao->nfuncgrads);
966770232b9SCe Qin   PetscFunctionReturn(0);
967770232b9SCe Qin }
968770232b9SCe Qin 
969770232b9SCe Qin #undef __FUNCT__
970a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetMaximumIterations"
971a7e14dcfSSatish Balay /*@
972a7e14dcfSSatish Balay    TaoSetMaximumIterations - Sets a maximum number of iterates.
973a7e14dcfSSatish Balay 
974441846f8SBarry Smith    Logically Collective on Tao
975a7e14dcfSSatish Balay 
976a7e14dcfSSatish Balay    Input Parameters:
977441846f8SBarry Smith +  tao - the Tao solver context
978a7e14dcfSSatish Balay -  maxits - the maximum number of iterates (>=0)
979a7e14dcfSSatish Balay 
980a7e14dcfSSatish Balay    Options Database Keys:
981a7e14dcfSSatish Balay .    -tao_max_it <its> - sets the maximum number of iterations
982a7e14dcfSSatish Balay 
983a7e14dcfSSatish Balay    Level: intermediate
984a7e14dcfSSatish Balay 
985a7e14dcfSSatish Balay .seealso: TaoSetTolerances(), TaoSetMaximumFunctionEvaluations()
986a7e14dcfSSatish Balay @*/
987441846f8SBarry Smith PetscErrorCode TaoSetMaximumIterations(Tao tao,PetscInt maxits)
988a7e14dcfSSatish Balay {
989a7e14dcfSSatish Balay   PetscFunctionBegin;
990441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
99147a47007SBarry Smith   tao->max_it = PetscMax(0,maxits);
9926552cf8aSJason Sarich   tao->max_it_changed=PETSC_TRUE;
993a7e14dcfSSatish Balay   PetscFunctionReturn(0);
994a7e14dcfSSatish Balay }
995a7e14dcfSSatish Balay 
996a7e14dcfSSatish Balay #undef __FUNCT__
997a7e14dcfSSatish Balay #define __FUNCT__ "TaoGetMaximumIterations"
998a7e14dcfSSatish Balay /*@
999a7e14dcfSSatish Balay    TaoGetMaximumIterations - Sets a maximum number of iterates.
1000a7e14dcfSSatish Balay 
1001a7e14dcfSSatish Balay    Not Collective
1002a7e14dcfSSatish Balay 
1003a7e14dcfSSatish Balay    Input Parameters:
1004441846f8SBarry Smith .  tao - the Tao solver context
1005a7e14dcfSSatish Balay 
1006a7e14dcfSSatish Balay    Output Parameters:
1007a7e14dcfSSatish Balay .  maxits - the maximum number of iterates
1008a7e14dcfSSatish Balay 
1009a7e14dcfSSatish Balay    Level: intermediate
1010a7e14dcfSSatish Balay 
1011a7e14dcfSSatish Balay .seealso: TaoSetMaximumIterations(), TaoGetMaximumFunctionEvaluations()
1012a7e14dcfSSatish Balay @*/
1013441846f8SBarry Smith PetscErrorCode TaoGetMaximumIterations(Tao tao,PetscInt *maxits)
1014a7e14dcfSSatish Balay {
1015a7e14dcfSSatish Balay   PetscFunctionBegin;
1016441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1017a7e14dcfSSatish Balay   *maxits = tao->max_it;
1018a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1019a7e14dcfSSatish Balay }
1020a7e14dcfSSatish Balay 
1021a7e14dcfSSatish Balay #undef __FUNCT__
1022a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetInitialTrustRegionRadius"
1023a7e14dcfSSatish Balay /*@
1024a7e14dcfSSatish Balay    TaoSetInitialTrustRegionRadius - Sets the initial trust region radius.
1025a7e14dcfSSatish Balay 
1026441846f8SBarry Smith    Logically collective on Tao
1027a7e14dcfSSatish Balay 
1028a7e14dcfSSatish Balay    Input Parameter:
1029a7e14dcfSSatish Balay +  tao - a TAO optimization solver
1030a7e14dcfSSatish Balay -  radius - the trust region radius
1031a7e14dcfSSatish Balay 
1032a7e14dcfSSatish Balay    Level: intermediate
1033a7e14dcfSSatish Balay 
1034a7e14dcfSSatish Balay    Options Database Key:
1035a7e14dcfSSatish Balay .  -tao_trust0 <t0> - sets initial trust region radius
1036a7e14dcfSSatish Balay 
1037a7e14dcfSSatish Balay .seealso: TaoGetTrustRegionRadius(), TaoSetTrustRegionTolerance()
1038a7e14dcfSSatish Balay @*/
1039441846f8SBarry Smith PetscErrorCode TaoSetInitialTrustRegionRadius(Tao tao, PetscReal radius)
1040a7e14dcfSSatish Balay {
1041a7e14dcfSSatish Balay   PetscFunctionBegin;
1042441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1043a7e14dcfSSatish Balay   tao->trust0 = PetscMax(0.0,radius);
10446552cf8aSJason Sarich   tao->trust0_changed=PETSC_TRUE;
1045a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1046a7e14dcfSSatish Balay }
1047a7e14dcfSSatish Balay 
1048a7e14dcfSSatish Balay #undef __FUNCT__
1049a7e14dcfSSatish Balay #define __FUNCT__ "TaoGetInitialTrustRegionRadius"
1050a7e14dcfSSatish Balay /*@
1051a7e14dcfSSatish Balay    TaoGetInitialTrustRegionRadius - Sets the initial trust region radius.
1052a7e14dcfSSatish Balay 
1053a7e14dcfSSatish Balay    Not Collective
1054a7e14dcfSSatish Balay 
1055a7e14dcfSSatish Balay    Input Parameter:
1056a7e14dcfSSatish Balay .  tao - a TAO optimization solver
1057a7e14dcfSSatish Balay 
1058a7e14dcfSSatish Balay    Output Parameter:
1059a7e14dcfSSatish Balay .  radius - the trust region radius
1060a7e14dcfSSatish Balay 
1061a7e14dcfSSatish Balay    Level: intermediate
1062a7e14dcfSSatish Balay 
1063a7e14dcfSSatish Balay .seealso: TaoSetInitialTrustRegionRadius(), TaoGetCurrentTrustRegionRadius()
1064a7e14dcfSSatish Balay @*/
1065441846f8SBarry Smith PetscErrorCode TaoGetInitialTrustRegionRadius(Tao tao, PetscReal *radius)
1066a7e14dcfSSatish Balay {
1067a7e14dcfSSatish Balay   PetscFunctionBegin;
1068441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1069a7e14dcfSSatish Balay   *radius = tao->trust0;
1070a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1071a7e14dcfSSatish Balay }
1072a7e14dcfSSatish Balay 
1073a7e14dcfSSatish Balay #undef __FUNCT__
1074a7e14dcfSSatish Balay #define __FUNCT__ "TaoGetCurrentTrustRegionRadius"
1075a7e14dcfSSatish Balay /*@
1076a7e14dcfSSatish Balay    TaoGetCurrentTrustRegionRadius - Gets the current trust region radius.
1077a7e14dcfSSatish Balay 
1078a7e14dcfSSatish Balay    Not Collective
1079a7e14dcfSSatish Balay 
1080a7e14dcfSSatish Balay    Input Parameter:
1081a7e14dcfSSatish Balay .  tao - a TAO optimization solver
1082a7e14dcfSSatish Balay 
1083a7e14dcfSSatish Balay    Output Parameter:
1084a7e14dcfSSatish Balay .  radius - the trust region radius
1085a7e14dcfSSatish Balay 
1086a7e14dcfSSatish Balay    Level: intermediate
1087a7e14dcfSSatish Balay 
1088a7e14dcfSSatish Balay .seealso: TaoSetInitialTrustRegionRadius(), TaoGetInitialTrustRegionRadius()
1089a7e14dcfSSatish Balay @*/
1090441846f8SBarry Smith PetscErrorCode TaoGetCurrentTrustRegionRadius(Tao tao, PetscReal *radius)
1091a7e14dcfSSatish Balay {
1092a7e14dcfSSatish Balay   PetscFunctionBegin;
1093441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1094a7e14dcfSSatish Balay   *radius = tao->trust;
1095a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1096a7e14dcfSSatish Balay }
1097a7e14dcfSSatish Balay 
1098a7e14dcfSSatish Balay #undef __FUNCT__
1099a7e14dcfSSatish Balay #define __FUNCT__ "TaoGetTolerances"
1100a7e14dcfSSatish Balay /*@
1101a7e14dcfSSatish Balay   TaoGetTolerances - gets the current values of tolerances
1102a7e14dcfSSatish Balay 
1103a7e14dcfSSatish Balay   Not Collective
1104a7e14dcfSSatish Balay 
1105a7e14dcfSSatish Balay   Input Parameters:
1106441846f8SBarry Smith . tao - the Tao context
1107a7e14dcfSSatish Balay 
1108a7e14dcfSSatish Balay   Output Parameters:
1109e52336cbSBarry Smith + gatol - stop if norm of gradient is less than this
1110a7e14dcfSSatish Balay . grtol - stop if relative norm of gradient is less than this
1111a7e14dcfSSatish Balay - gttol - stop if norm of gradient is reduced by a this factor
1112a7e14dcfSSatish Balay 
111347a47007SBarry Smith   Note: NULL can be used as an argument if not all tolerances values are needed
1114a7e14dcfSSatish Balay 
1115a7e14dcfSSatish Balay .seealso TaoSetTolerances()
1116a7e14dcfSSatish Balay 
1117a7e14dcfSSatish Balay   Level: intermediate
1118a7e14dcfSSatish Balay @*/
1119e52336cbSBarry Smith PetscErrorCode TaoGetTolerances(Tao tao, PetscReal *gatol, PetscReal *grtol, PetscReal *gttol)
1120a7e14dcfSSatish Balay {
1121a7e14dcfSSatish Balay   PetscFunctionBegin;
1122441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1123a7e14dcfSSatish Balay   if (gatol) *gatol=tao->gatol;
1124a7e14dcfSSatish Balay   if (grtol) *grtol=tao->grtol;
1125a7e14dcfSSatish Balay   if (gttol) *gttol=tao->gttol;
1126a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1127a7e14dcfSSatish Balay }
1128a7e14dcfSSatish Balay 
1129a7e14dcfSSatish Balay #undef __FUNCT__
1130a7e14dcfSSatish Balay #define __FUNCT__ "TaoGetKSP"
1131a7e14dcfSSatish Balay /*@
1132a7e14dcfSSatish Balay   TaoGetKSP - Gets the linear solver used by the optimization solver.
1133a7e14dcfSSatish Balay   Application writers should use TaoGetKSP if they need direct access
1134a7e14dcfSSatish Balay   to the PETSc KSP object.
1135a7e14dcfSSatish Balay 
1136a7e14dcfSSatish Balay   Not Collective
1137a7e14dcfSSatish Balay 
1138a7e14dcfSSatish Balay    Input Parameters:
1139a7e14dcfSSatish Balay .  tao - the TAO solver
1140a7e14dcfSSatish Balay 
1141a7e14dcfSSatish Balay    Output Parameters:
1142a7e14dcfSSatish Balay .  ksp - the KSP linear solver used in the optimization solver
1143a7e14dcfSSatish Balay 
1144a7e14dcfSSatish Balay    Level: intermediate
1145a7e14dcfSSatish Balay 
1146a7e14dcfSSatish Balay @*/
1147441846f8SBarry Smith PetscErrorCode TaoGetKSP(Tao tao, KSP *ksp)
114847a47007SBarry Smith {
1149a7e14dcfSSatish Balay   PetscFunctionBegin;
1150a7e14dcfSSatish Balay   *ksp = tao->ksp;
1151a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1152a7e14dcfSSatish Balay }
1153a7e14dcfSSatish Balay 
1154a7e14dcfSSatish Balay #undef __FUNCT__
1155025e9500SJason Sarich #define __FUNCT__ "TaoGetLinearSolveIterations"
1156025e9500SJason Sarich /*@
1157025e9500SJason Sarich    TaoGetLinearSolveIterations - Gets the total number of linear iterations
1158025e9500SJason Sarich    used by the TAO solver
1159025e9500SJason Sarich 
1160025e9500SJason Sarich    Not Collective
1161025e9500SJason Sarich 
1162025e9500SJason Sarich    Input Parameter:
1163025e9500SJason Sarich .  tao - TAO context
1164025e9500SJason Sarich 
1165025e9500SJason Sarich    Output Parameter:
1166025e9500SJason Sarich .  lits - number of linear iterations
1167025e9500SJason Sarich 
1168025e9500SJason Sarich    Notes:
1169025e9500SJason Sarich    This counter is reset to zero for each successive call to TaoSolve()
1170025e9500SJason Sarich 
1171025e9500SJason Sarich    Level: intermediate
1172025e9500SJason Sarich 
1173025e9500SJason Sarich .keywords: TAO
1174025e9500SJason Sarich 
1175025e9500SJason Sarich .seealso:  TaoGetKSP()
1176025e9500SJason Sarich @*/
1177025e9500SJason Sarich PetscErrorCode  TaoGetLinearSolveIterations(Tao tao,PetscInt *lits)
1178025e9500SJason Sarich {
1179025e9500SJason Sarich   PetscFunctionBegin;
1180025e9500SJason Sarich   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1181025e9500SJason Sarich   PetscValidIntPointer(lits,2);
11822d9aa51bSJason Sarich   *lits = tao->ksp_tot_its;
1183025e9500SJason Sarich   PetscFunctionReturn(0);
1184025e9500SJason Sarich }
1185025e9500SJason Sarich 
1186025e9500SJason Sarich #undef __FUNCT__
1187a7e14dcfSSatish Balay #define __FUNCT__ "TaoGetLineSearch"
1188a7e14dcfSSatish Balay /*@
1189a7e14dcfSSatish Balay   TaoGetLineSearch - Gets the line search used by the optimization solver.
1190a7e14dcfSSatish Balay   Application writers should use TaoGetLineSearch if they need direct access
1191a7e14dcfSSatish Balay   to the TaoLineSearch object.
1192a7e14dcfSSatish Balay 
1193a7e14dcfSSatish Balay   Not Collective
1194a7e14dcfSSatish Balay 
1195a7e14dcfSSatish Balay    Input Parameters:
1196a7e14dcfSSatish Balay .  tao - the TAO solver
1197a7e14dcfSSatish Balay 
1198a7e14dcfSSatish Balay    Output Parameters:
1199a7e14dcfSSatish Balay .  ls - the line search used in the optimization solver
1200a7e14dcfSSatish Balay 
1201a7e14dcfSSatish Balay    Level: intermediate
1202a7e14dcfSSatish Balay 
1203a7e14dcfSSatish Balay @*/
1204441846f8SBarry Smith PetscErrorCode TaoGetLineSearch(Tao tao, TaoLineSearch *ls)
120547a47007SBarry Smith {
1206a7e14dcfSSatish Balay   PetscFunctionBegin;
1207a7e14dcfSSatish Balay   *ls = tao->linesearch;
1208a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1209a7e14dcfSSatish Balay }
1210a7e14dcfSSatish Balay 
1211a7e14dcfSSatish Balay #undef __FUNCT__
1212a7e14dcfSSatish Balay #define __FUNCT__ "TaoAddLineSearchCounts"
1213a7e14dcfSSatish Balay /*@
1214a7e14dcfSSatish Balay   TaoAddLineSearchCounts - Adds the number of function evaluations spent
1215a7e14dcfSSatish Balay   in the line search to the running total.
1216a7e14dcfSSatish Balay 
1217a7e14dcfSSatish Balay    Input Parameters:
1218a7e14dcfSSatish Balay +  tao - the TAO solver
1219a7e14dcfSSatish Balay -  ls - the line search used in the optimization solver
1220a7e14dcfSSatish Balay 
1221a7e14dcfSSatish Balay    Level: developer
1222a7e14dcfSSatish Balay 
1223a7e14dcfSSatish Balay .seealso: TaoLineSearchApply()
1224a7e14dcfSSatish Balay @*/
1225441846f8SBarry Smith PetscErrorCode TaoAddLineSearchCounts(Tao tao)
122647a47007SBarry Smith {
1227a7e14dcfSSatish Balay   PetscErrorCode ierr;
1228a7e14dcfSSatish Balay   PetscBool      flg;
1229a7e14dcfSSatish Balay   PetscInt       nfeval,ngeval,nfgeval;
123047a47007SBarry Smith 
1231a7e14dcfSSatish Balay   PetscFunctionBegin;
1232441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1233a7e14dcfSSatish Balay   if (tao->linesearch) {
1234302440fdSBarry Smith     ierr = TaoLineSearchIsUsingTaoRoutines(tao->linesearch,&flg);CHKERRQ(ierr);
123594ae4db5SBarry Smith     if (!flg) {
123647a47007SBarry Smith       ierr = TaoLineSearchGetNumberFunctionEvaluations(tao->linesearch,&nfeval,&ngeval,&nfgeval);CHKERRQ(ierr);
1237a7e14dcfSSatish Balay       tao->nfuncs+=nfeval;
1238a7e14dcfSSatish Balay       tao->ngrads+=ngeval;
1239a7e14dcfSSatish Balay       tao->nfuncgrads+=nfgeval;
1240a7e14dcfSSatish Balay     }
1241a7e14dcfSSatish Balay   }
1242a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1243a7e14dcfSSatish Balay }
1244a7e14dcfSSatish Balay 
1245a7e14dcfSSatish Balay #undef __FUNCT__
1246a7e14dcfSSatish Balay #define __FUNCT__ "TaoGetSolutionVector"
1247a7e14dcfSSatish Balay /*@
1248a7e14dcfSSatish Balay   TaoGetSolutionVector - Returns the vector with the current TAO solution
1249a7e14dcfSSatish Balay 
1250a7e14dcfSSatish Balay   Not Collective
1251a7e14dcfSSatish Balay 
1252a7e14dcfSSatish Balay   Input Parameter:
1253441846f8SBarry Smith . tao - the Tao context
1254a7e14dcfSSatish Balay 
1255a7e14dcfSSatish Balay   Output Parameter:
1256a7e14dcfSSatish Balay . X - the current solution
1257a7e14dcfSSatish Balay 
1258a7e14dcfSSatish Balay   Level: intermediate
1259a7e14dcfSSatish Balay 
1260a7e14dcfSSatish Balay   Note:  The returned vector will be the same object that was passed into TaoSetInitialVector()
1261a7e14dcfSSatish Balay @*/
1262441846f8SBarry Smith PetscErrorCode TaoGetSolutionVector(Tao tao, Vec *X)
1263a7e14dcfSSatish Balay {
1264a7e14dcfSSatish Balay   PetscFunctionBegin;
1265441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1266a7e14dcfSSatish Balay   *X = tao->solution;
1267a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1268a7e14dcfSSatish Balay }
1269a7e14dcfSSatish Balay 
1270a7e14dcfSSatish Balay #undef __FUNCT__
1271a7e14dcfSSatish Balay #define __FUNCT__ "TaoGetGradientVector"
1272a7e14dcfSSatish Balay /*@
1273a7e14dcfSSatish Balay   TaoGetGradientVector - Returns the vector with the current TAO gradient
1274a7e14dcfSSatish Balay 
1275a7e14dcfSSatish Balay   Not Collective
1276a7e14dcfSSatish Balay 
1277a7e14dcfSSatish Balay   Input Parameter:
1278441846f8SBarry Smith . tao - the Tao context
1279a7e14dcfSSatish Balay 
1280a7e14dcfSSatish Balay   Output Parameter:
1281a7e14dcfSSatish Balay . G - the current solution
1282a7e14dcfSSatish Balay 
1283a7e14dcfSSatish Balay   Level: intermediate
1284a7e14dcfSSatish Balay @*/
1285441846f8SBarry Smith PetscErrorCode TaoGetGradientVector(Tao tao, Vec *G)
1286a7e14dcfSSatish Balay {
1287a7e14dcfSSatish Balay   PetscFunctionBegin;
1288441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1289a7e14dcfSSatish Balay   *G = tao->gradient;
1290a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1291a7e14dcfSSatish Balay }
1292a7e14dcfSSatish Balay 
1293a7e14dcfSSatish Balay #undef __FUNCT__
1294a7e14dcfSSatish Balay #define __FUNCT__ "TaoResetStatistics"
1295a7e14dcfSSatish Balay /*@
1296a7e14dcfSSatish Balay    TaoResetStatistics - Initialize the statistics used by TAO for all of the solvers.
1297a7e14dcfSSatish Balay    These statistics include the iteration number, residual norms, and convergence status.
1298a7e14dcfSSatish Balay    This routine gets called before solving each optimization problem.
1299a7e14dcfSSatish Balay 
1300441846f8SBarry Smith    Collective on Tao
1301a7e14dcfSSatish Balay 
1302a7e14dcfSSatish Balay    Input Parameters:
1303441846f8SBarry Smith .  solver - the Tao context
1304a7e14dcfSSatish Balay 
1305a7e14dcfSSatish Balay    Level: developer
1306a7e14dcfSSatish Balay 
1307a7e14dcfSSatish Balay .seealso: TaoCreate(), TaoSolve()
1308a7e14dcfSSatish Balay @*/
1309441846f8SBarry Smith PetscErrorCode TaoResetStatistics(Tao tao)
1310a7e14dcfSSatish Balay {
1311a7e14dcfSSatish Balay   PetscFunctionBegin;
1312441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1313a7e14dcfSSatish Balay   tao->niter        = 0;
1314a7e14dcfSSatish Balay   tao->nfuncs       = 0;
1315a7e14dcfSSatish Balay   tao->nfuncgrads   = 0;
1316a7e14dcfSSatish Balay   tao->ngrads       = 0;
1317a7e14dcfSSatish Balay   tao->nhess        = 0;
1318a7e14dcfSSatish Balay   tao->njac         = 0;
1319a7e14dcfSSatish Balay   tao->nconstraints = 0;
1320a7e14dcfSSatish Balay   tao->ksp_its      = 0;
1321ae93cb3cSJason Sarich   tao->ksp_tot_its      = 0;
1322a7e14dcfSSatish Balay   tao->reason       = TAO_CONTINUE_ITERATING;
1323a7e14dcfSSatish Balay   tao->residual     = 0.0;
1324a7e14dcfSSatish Balay   tao->cnorm        = 0.0;
1325a7e14dcfSSatish Balay   tao->step         = 0.0;
1326a7e14dcfSSatish Balay   tao->lsflag       = PETSC_FALSE;
1327a7e14dcfSSatish Balay   if (tao->hist_reset) tao->hist_len=0;
1328a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1329a7e14dcfSSatish Balay }
1330a7e14dcfSSatish Balay 
1331a7e14dcfSSatish Balay #undef __FUNCT__
1332a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetConvergenceTest"
1333a7e14dcfSSatish Balay /*@C
1334a7e14dcfSSatish Balay   TaoSetConvergenceTest - Sets the function that is to be used to test
1335a7e14dcfSSatish Balay   for convergence o fthe iterative minimization solution.  The new convergence
1336a7e14dcfSSatish Balay   testing routine will replace TAO's default convergence test.
1337a7e14dcfSSatish Balay 
1338441846f8SBarry Smith   Logically Collective on Tao
1339a7e14dcfSSatish Balay 
1340a7e14dcfSSatish Balay   Input Parameters:
1341441846f8SBarry Smith + tao - the Tao object
1342a7e14dcfSSatish Balay . conv - the routine to test for convergence
1343a7e14dcfSSatish Balay - ctx - [optional] context for private data for the convergence routine
134447a47007SBarry Smith         (may be NULL)
1345a7e14dcfSSatish Balay 
1346a7e14dcfSSatish Balay   Calling sequence of conv:
1347441846f8SBarry Smith $   PetscErrorCode conv(Tao tao, void *ctx)
1348a7e14dcfSSatish Balay 
1349441846f8SBarry Smith + tao - the Tao object
1350a7e14dcfSSatish Balay - ctx - [optional] convergence context
1351a7e14dcfSSatish Balay 
1352e4cb33bbSBarry Smith   Note: The new convergence testing routine should call TaoSetConvergedReason().
1353a7e14dcfSSatish Balay 
1354a7e14dcfSSatish Balay   Level: advanced
1355a7e14dcfSSatish Balay 
1356e4cb33bbSBarry Smith .seealso: TaoSetConvergedReason(), TaoGetSolutionStatus(), TaoGetTolerances(), TaoSetMonitor
1357a7e14dcfSSatish Balay 
1358a7e14dcfSSatish Balay @*/
1359441846f8SBarry Smith PetscErrorCode TaoSetConvergenceTest(Tao tao, PetscErrorCode (*conv)(Tao,void*), void *ctx)
1360a7e14dcfSSatish Balay {
1361a7e14dcfSSatish Balay   PetscFunctionBegin;
1362441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1363a7e14dcfSSatish Balay   (tao)->ops->convergencetest = conv;
1364a7e14dcfSSatish Balay   (tao)->cnvP = ctx;
1365a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1366a7e14dcfSSatish Balay }
1367a7e14dcfSSatish Balay 
1368a7e14dcfSSatish Balay #undef __FUNCT__
1369a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetMonitor"
1370a7e14dcfSSatish Balay /*@C
1371a7e14dcfSSatish Balay    TaoSetMonitor - Sets an ADDITIONAL function that is to be used at every
1372a7e14dcfSSatish Balay    iteration of the solver to display the iteration's
1373a7e14dcfSSatish Balay    progress.
1374a7e14dcfSSatish Balay 
1375441846f8SBarry Smith    Logically Collective on Tao
1376a7e14dcfSSatish Balay 
1377a7e14dcfSSatish Balay    Input Parameters:
1378441846f8SBarry Smith +  tao - the Tao solver context
1379a7e14dcfSSatish Balay .  mymonitor - monitoring routine
1380a7e14dcfSSatish Balay -  mctx - [optional] user-defined context for private data for the
138147a47007SBarry Smith           monitor routine (may be NULL)
1382a7e14dcfSSatish Balay 
1383a7e14dcfSSatish Balay    Calling sequence of mymonitor:
1384441846f8SBarry Smith $     int mymonitor(Tao tao,void *mctx)
1385a7e14dcfSSatish Balay 
1386441846f8SBarry Smith +    tao - the Tao solver context
1387a7e14dcfSSatish Balay -    mctx - [optional] monitoring context
1388a7e14dcfSSatish Balay 
1389a7e14dcfSSatish Balay 
1390a7e14dcfSSatish Balay    Options Database Keys:
1391a7e14dcfSSatish Balay +    -tao_monitor        - sets TaoDefaultMonitor()
1392a7e14dcfSSatish Balay .    -tao_smonitor       - sets short monitor
1393a7e14dcfSSatish Balay .    -tao_cmonitor       - same as smonitor plus constraint norm
1394a7e14dcfSSatish Balay .    -tao_view_solution   - view solution at each iteration
1395a7e14dcfSSatish Balay .    -tao_view_gradient   - view gradient at each iteration
1396a7e14dcfSSatish Balay .    -tao_view_separableobjective - view separable objective function at each iteration
1397a7e14dcfSSatish Balay -    -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.
1398a7e14dcfSSatish Balay 
1399a7e14dcfSSatish Balay 
1400a7e14dcfSSatish Balay    Notes:
1401a7e14dcfSSatish Balay    Several different monitoring routines may be set by calling
1402a7e14dcfSSatish Balay    TaoSetMonitor() multiple times; all will be called in the
1403a7e14dcfSSatish Balay    order in which they were set.
1404a7e14dcfSSatish Balay 
1405a7e14dcfSSatish Balay    Fortran Notes: Only one monitor function may be set
1406a7e14dcfSSatish Balay 
1407a7e14dcfSSatish Balay    Level: intermediate
1408a7e14dcfSSatish Balay 
1409a7e14dcfSSatish Balay .seealso: TaoDefaultMonitor(), TaoCancelMonitors(),  TaoSetDestroyRoutine()
1410a7e14dcfSSatish Balay @*/
1411441846f8SBarry Smith PetscErrorCode TaoSetMonitor(Tao tao, PetscErrorCode (*func)(Tao, void*), void *ctx,PetscErrorCode (*dest)(void**))
1412a7e14dcfSSatish Balay {
141308d19d1fSJason Sarich   PetscErrorCode ierr;
141408d19d1fSJason Sarich   PetscInt       i;
141508d19d1fSJason Sarich 
1416a7e14dcfSSatish Balay   PetscFunctionBegin;
1417441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
141847a47007SBarry Smith   if (tao->numbermonitors >= MAXTAOMONITORS) SETERRQ1(PETSC_COMM_SELF,1,"Cannot attach another monitor -- max=",MAXTAOMONITORS);
141908d19d1fSJason Sarich 
142008d19d1fSJason Sarich   for (i=0; i<tao->numbermonitors;i++) {
142108d19d1fSJason Sarich     if (func == tao->monitor[i] && dest == tao->monitordestroy[i] && ctx == tao->monitorcontext[i]) {
142208d19d1fSJason Sarich       if (dest) {
142308d19d1fSJason Sarich         ierr = (*dest)(&ctx);CHKERRQ(ierr);
142408d19d1fSJason Sarich       }
142508d19d1fSJason Sarich       PetscFunctionReturn(0);
142608d19d1fSJason Sarich     }
142708d19d1fSJason Sarich   }
1428a7e14dcfSSatish Balay   tao->monitor[tao->numbermonitors] = func;
1429a7e14dcfSSatish Balay   tao->monitorcontext[tao->numbermonitors] = ctx;
1430a7e14dcfSSatish Balay   tao->monitordestroy[tao->numbermonitors] = dest;
1431a7e14dcfSSatish Balay   ++tao->numbermonitors;
1432a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1433a7e14dcfSSatish Balay }
1434a7e14dcfSSatish Balay 
1435a7e14dcfSSatish Balay #undef __FUNCT__
1436a7e14dcfSSatish Balay #define __FUNCT__ "TaoCancelMonitors"
1437a7e14dcfSSatish Balay /*@
1438441846f8SBarry Smith    TaoCancelMonitors - Clears all the monitor functions for a Tao object.
1439a7e14dcfSSatish Balay 
1440441846f8SBarry Smith    Logically Collective on Tao
1441a7e14dcfSSatish Balay 
1442a7e14dcfSSatish Balay    Input Parameters:
1443441846f8SBarry Smith .  tao - the Tao solver context
1444a7e14dcfSSatish Balay 
1445a7e14dcfSSatish Balay    Options Database:
1446a7e14dcfSSatish Balay .  -tao_cancelmonitors - cancels all monitors that have been hardwired
1447a7e14dcfSSatish Balay     into a code by calls to TaoSetMonitor(), but does not cancel those
1448a7e14dcfSSatish Balay     set via the options database
1449a7e14dcfSSatish Balay 
1450a7e14dcfSSatish Balay    Notes:
1451441846f8SBarry Smith    There is no way to clear one specific monitor from a Tao object.
1452a7e14dcfSSatish Balay 
1453a7e14dcfSSatish Balay    Level: advanced
1454a7e14dcfSSatish Balay 
1455a7e14dcfSSatish Balay .seealso: TaoDefaultMonitor(), TaoSetMonitor()
1456a7e14dcfSSatish Balay @*/
1457441846f8SBarry Smith PetscErrorCode TaoCancelMonitors(Tao tao)
1458a7e14dcfSSatish Balay {
1459a7e14dcfSSatish Balay   PetscInt       i;
1460a7e14dcfSSatish Balay   PetscErrorCode ierr;
146147a47007SBarry Smith 
1462a7e14dcfSSatish Balay   PetscFunctionBegin;
1463441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1464a7e14dcfSSatish Balay   for (i=0;i<tao->numbermonitors;i++) {
1465a7e14dcfSSatish Balay     if (tao->monitordestroy[i]) {
1466a7e14dcfSSatish Balay       ierr = (*tao->monitordestroy[i])(&tao->monitorcontext[i]);CHKERRQ(ierr);
1467a7e14dcfSSatish Balay     }
1468a7e14dcfSSatish Balay   }
1469a7e14dcfSSatish Balay   tao->numbermonitors=0;
1470a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1471a7e14dcfSSatish Balay }
1472a7e14dcfSSatish Balay 
1473a7e14dcfSSatish Balay #undef __FUNCT__
1474a7e14dcfSSatish Balay #define __FUNCT__ "TaoDefaultMonitor"
14757fab98ebSJason Sarich /*@
1476a7e14dcfSSatish Balay    TaoDefaultMonitor - Default routine for monitoring progress of the
1477441846f8SBarry Smith    Tao solvers (default).  This monitor prints the function value and gradient
1478a7e14dcfSSatish Balay    norm at each iteration.  It can be turned on from the command line using the
1479a7e14dcfSSatish Balay    -tao_monitor option
1480a7e14dcfSSatish Balay 
1481441846f8SBarry Smith    Collective on Tao
1482a7e14dcfSSatish Balay 
1483a7e14dcfSSatish Balay    Input Parameters:
1484441846f8SBarry Smith +  tao - the Tao context
148547a47007SBarry Smith -  ctx - PetscViewer context or NULL
1486a7e14dcfSSatish Balay 
1487a7e14dcfSSatish Balay    Options Database Keys:
1488a7e14dcfSSatish Balay .  -tao_monitor
1489a7e14dcfSSatish Balay 
1490a7e14dcfSSatish Balay    Level: advanced
1491a7e14dcfSSatish Balay 
1492a7e14dcfSSatish Balay .seealso: TaoDefaultSMonitor(), TaoSetMonitor()
1493a7e14dcfSSatish Balay @*/
1494441846f8SBarry Smith PetscErrorCode TaoDefaultMonitor(Tao tao, void *ctx)
1495a7e14dcfSSatish Balay {
1496a7e14dcfSSatish Balay   PetscErrorCode ierr;
1497a7e14dcfSSatish Balay   PetscInt       its;
1498a7e14dcfSSatish Balay   PetscReal      fct,gnorm;
14998163d661SBarry Smith   PetscViewer    viewer = (PetscViewer)ctx;
150047a47007SBarry Smith 
1501a7e14dcfSSatish Balay   PetscFunctionBegin;
15028163d661SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
1503a7e14dcfSSatish Balay   its=tao->niter;
1504a7e14dcfSSatish Balay   fct=tao->fc;
1505a7e14dcfSSatish Balay   gnorm=tao->residual;
1506a7e14dcfSSatish Balay   ierr=PetscViewerASCIIPrintf(viewer,"iter = %3D,",its);CHKERRQ(ierr);
150747a47007SBarry Smith   ierr=PetscViewerASCIIPrintf(viewer," Function value: %g,",(double)fct);CHKERRQ(ierr);
150808d19d1fSJason Sarich   if (gnorm >= PETSC_INFINITY) {
150908d19d1fSJason Sarich     ierr=PetscViewerASCIIPrintf(viewer,"  Residual: Inf \n");CHKERRQ(ierr);
151008d19d1fSJason Sarich   } else {
151147a47007SBarry Smith     ierr=PetscViewerASCIIPrintf(viewer,"  Residual: %g \n",(double)gnorm);CHKERRQ(ierr);
151208d19d1fSJason Sarich   }
1513a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1514a7e14dcfSSatish Balay }
1515a7e14dcfSSatish Balay 
1516a7e14dcfSSatish Balay #undef __FUNCT__
1517a7e14dcfSSatish Balay #define __FUNCT__ "TaoDefaultSMonitor"
15187fab98ebSJason Sarich /*@
1519a7e14dcfSSatish Balay    TaoDefaultSMonitor - Default routine for monitoring progress of the
1520a7e14dcfSSatish Balay    solver. Same as TaoDefaultMonitor() except
1521a7e14dcfSSatish Balay    it prints fewer digits of the residual as the residual gets smaller.
1522a7e14dcfSSatish Balay    This is because the later digits are meaningless and are often
1523a7e14dcfSSatish Balay    different on different machines; by using this routine different
1524a7e14dcfSSatish Balay    machines will usually generate the same output. It can be turned on
1525a7e14dcfSSatish Balay    by using the -tao_smonitor option
1526a7e14dcfSSatish Balay 
1527441846f8SBarry Smith    Collective on Tao
1528a7e14dcfSSatish Balay 
1529a7e14dcfSSatish Balay    Input Parameters:
1530441846f8SBarry Smith +  tao - the Tao context
15314d4332d5SBarry Smith -  ctx - PetscViewer context of type ASCII
1532a7e14dcfSSatish Balay 
1533a7e14dcfSSatish Balay    Options Database Keys:
1534a7e14dcfSSatish Balay .  -tao_smonitor
1535a7e14dcfSSatish Balay 
1536a7e14dcfSSatish Balay    Level: advanced
1537a7e14dcfSSatish Balay 
1538a7e14dcfSSatish Balay .seealso: TaoDefaultMonitor(), TaoSetMonitor()
1539a7e14dcfSSatish Balay @*/
1540441846f8SBarry Smith PetscErrorCode TaoDefaultSMonitor(Tao tao, void *ctx)
1541a7e14dcfSSatish Balay {
1542a7e14dcfSSatish Balay   PetscErrorCode ierr;
1543a7e14dcfSSatish Balay   PetscInt       its;
1544a7e14dcfSSatish Balay   PetscReal      fct,gnorm;
15454d4332d5SBarry Smith   PetscViewer    viewer = (PetscViewer)ctx;
154647a47007SBarry Smith 
1547a7e14dcfSSatish Balay   PetscFunctionBegin;
15484d4332d5SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
1549a7e14dcfSSatish Balay   its=tao->niter;
1550a7e14dcfSSatish Balay   fct=tao->fc;
1551a7e14dcfSSatish Balay   gnorm=tao->residual;
1552a7e14dcfSSatish Balay   ierr=PetscViewerASCIIPrintf(viewer,"iter = %3D,",its);CHKERRQ(ierr);
155347a47007SBarry Smith   ierr=PetscViewerASCIIPrintf(viewer," Function value %g,",(double)fct);CHKERRQ(ierr);
1554d393f493SBarry Smith   if (gnorm >= PETSC_INFINITY) {
155508d19d1fSJason Sarich     ierr=PetscViewerASCIIPrintf(viewer," Residual: Inf \n");CHKERRQ(ierr);
155608d19d1fSJason Sarich   } else if (gnorm > 1.e-6) {
155747a47007SBarry Smith     ierr=PetscViewerASCIIPrintf(viewer," Residual: %g \n",(double)gnorm);CHKERRQ(ierr);
1558a7e14dcfSSatish Balay   } else if (gnorm > 1.e-11) {
1559a7e14dcfSSatish Balay     ierr=PetscViewerASCIIPrintf(viewer," Residual: < 1.0e-6 \n");CHKERRQ(ierr);
1560a7e14dcfSSatish Balay   } else {
1561a7e14dcfSSatish Balay     ierr=PetscViewerASCIIPrintf(viewer," Residual: < 1.0e-11 \n");CHKERRQ(ierr);
1562a7e14dcfSSatish Balay   }
1563a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1564a7e14dcfSSatish Balay }
1565a7e14dcfSSatish Balay 
1566a7e14dcfSSatish Balay #undef __FUNCT__
1567a7e14dcfSSatish Balay #define __FUNCT__ "TaoDefaultCMonitor"
15687fab98ebSJason Sarich /*@
1569a7e14dcfSSatish Balay    TaoDefaultCMonitor - same as TaoDefaultMonitor() except
1570a7e14dcfSSatish Balay    it prints the norm of the constraints function. It can be turned on
1571a7e14dcfSSatish Balay    from the command line using the -tao_cmonitor option
1572a7e14dcfSSatish Balay 
1573441846f8SBarry Smith    Collective on Tao
1574a7e14dcfSSatish Balay 
1575a7e14dcfSSatish Balay    Input Parameters:
1576441846f8SBarry Smith +  tao - the Tao context
157747a47007SBarry Smith -  ctx - PetscViewer context or NULL
1578a7e14dcfSSatish Balay 
1579a7e14dcfSSatish Balay    Options Database Keys:
1580a7e14dcfSSatish Balay .  -tao_cmonitor
1581a7e14dcfSSatish Balay 
1582a7e14dcfSSatish Balay    Level: advanced
1583a7e14dcfSSatish Balay 
1584a7e14dcfSSatish Balay .seealso: TaoDefaultMonitor(), TaoSetMonitor()
1585a7e14dcfSSatish Balay @*/
1586441846f8SBarry Smith PetscErrorCode TaoDefaultCMonitor(Tao tao, void *ctx)
1587a7e14dcfSSatish Balay {
1588a7e14dcfSSatish Balay   PetscErrorCode ierr;
1589a7e14dcfSSatish Balay   PetscInt       its;
1590a7e14dcfSSatish Balay   PetscReal      fct,gnorm;
1591d393f493SBarry Smith   PetscViewer    viewer = (PetscViewer)ctx;
1592a7e14dcfSSatish Balay 
1593a7e14dcfSSatish Balay   PetscFunctionBegin;
1594d393f493SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
1595a7e14dcfSSatish Balay   its=tao->niter;
1596a7e14dcfSSatish Balay   fct=tao->fc;
1597a7e14dcfSSatish Balay   gnorm=tao->residual;
1598a7e14dcfSSatish Balay   ierr=PetscViewerASCIIPrintf(viewer,"iter = %D,",its);CHKERRQ(ierr);
159947a47007SBarry Smith   ierr=PetscViewerASCIIPrintf(viewer," Function value: %g,",(double)fct);CHKERRQ(ierr);
160047a47007SBarry Smith   ierr=PetscViewerASCIIPrintf(viewer,"  Residual: %g ",(double)gnorm);CHKERRQ(ierr);
160147a47007SBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"  Constraint: %g \n",(double)tao->cnorm);CHKERRQ(ierr);
1602a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1603a7e14dcfSSatish Balay }
1604a7e14dcfSSatish Balay 
1605a7e14dcfSSatish Balay #undef __FUNCT__
1606a7e14dcfSSatish Balay #define __FUNCT__ "TaoSolutionMonitor"
1607a7e14dcfSSatish Balay /*@C
1608a7e14dcfSSatish Balay    TaoSolutionMonitor - Views the solution at each iteration
1609a7e14dcfSSatish Balay    It can be turned on from the command line using the
1610a7e14dcfSSatish Balay    -tao_view_solution option
1611a7e14dcfSSatish Balay 
1612441846f8SBarry Smith    Collective on Tao
1613a7e14dcfSSatish Balay 
1614a7e14dcfSSatish Balay    Input Parameters:
1615441846f8SBarry Smith +  tao - the Tao context
161647a47007SBarry Smith -  ctx - PetscViewer context or NULL
1617a7e14dcfSSatish Balay 
1618a7e14dcfSSatish Balay    Options Database Keys:
1619a7e14dcfSSatish Balay .  -tao_view_solution
1620a7e14dcfSSatish Balay 
1621a7e14dcfSSatish Balay    Level: advanced
1622a7e14dcfSSatish Balay 
1623a7e14dcfSSatish Balay .seealso: TaoDefaultSMonitor(), TaoSetMonitor()
1624a7e14dcfSSatish Balay @*/
1625441846f8SBarry Smith PetscErrorCode TaoSolutionMonitor(Tao tao, void *ctx)
1626a7e14dcfSSatish Balay {
1627a7e14dcfSSatish Balay   PetscErrorCode ierr;
1628d393f493SBarry Smith   PetscViewer    viewer  = (PetscViewer)ctx;;
162947a47007SBarry Smith 
1630a7e14dcfSSatish Balay   PetscFunctionBegin;
1631d393f493SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
1632a7e14dcfSSatish Balay   ierr = VecView(tao->solution, viewer);CHKERRQ(ierr);
1633a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1634a7e14dcfSSatish Balay }
1635a7e14dcfSSatish Balay 
1636a7e14dcfSSatish Balay #undef __FUNCT__
1637a7e14dcfSSatish Balay #define __FUNCT__ "TaoGradientMonitor"
1638a7e14dcfSSatish Balay /*@C
1639a7e14dcfSSatish Balay    TaoGradientMonitor - Views the gradient at each iteration
1640a7e14dcfSSatish Balay    It can be turned on from the command line using the
1641a7e14dcfSSatish Balay    -tao_view_gradient option
1642a7e14dcfSSatish Balay 
1643441846f8SBarry Smith    Collective on Tao
1644a7e14dcfSSatish Balay 
1645a7e14dcfSSatish Balay    Input Parameters:
1646441846f8SBarry Smith +  tao - the Tao context
164747a47007SBarry Smith -  ctx - PetscViewer context or NULL
1648a7e14dcfSSatish Balay 
1649a7e14dcfSSatish Balay    Options Database Keys:
1650a7e14dcfSSatish Balay .  -tao_view_gradient
1651a7e14dcfSSatish Balay 
1652a7e14dcfSSatish Balay    Level: advanced
1653a7e14dcfSSatish Balay 
1654a7e14dcfSSatish Balay .seealso: TaoDefaultSMonitor(), TaoSetMonitor()
1655a7e14dcfSSatish Balay @*/
1656441846f8SBarry Smith PetscErrorCode TaoGradientMonitor(Tao tao, void *ctx)
1657a7e14dcfSSatish Balay {
1658a7e14dcfSSatish Balay   PetscErrorCode ierr;
1659d393f493SBarry Smith   PetscViewer    viewer = (PetscViewer)ctx;
166047a47007SBarry Smith 
1661a7e14dcfSSatish Balay   PetscFunctionBegin;
1662d393f493SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
1663a7e14dcfSSatish Balay   ierr = VecView(tao->gradient, viewer);CHKERRQ(ierr);
1664a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1665a7e14dcfSSatish Balay }
1666a7e14dcfSSatish Balay 
1667a7e14dcfSSatish Balay #undef __FUNCT__
1668a7e14dcfSSatish Balay #define __FUNCT__ "TaoStepDirectionMonitor"
1669a7e14dcfSSatish Balay /*@C
1670a7e14dcfSSatish Balay    TaoStepDirectionMonitor - Views the gradient at each iteration
1671a7e14dcfSSatish Balay    It can be turned on from the command line using the
1672a7e14dcfSSatish Balay    -tao_view_gradient option
1673a7e14dcfSSatish Balay 
1674441846f8SBarry Smith    Collective on Tao
1675a7e14dcfSSatish Balay 
1676a7e14dcfSSatish Balay    Input Parameters:
1677441846f8SBarry Smith +  tao - the Tao context
167847a47007SBarry Smith -  ctx - PetscViewer context or NULL
1679a7e14dcfSSatish Balay 
1680a7e14dcfSSatish Balay    Options Database Keys:
1681a7e14dcfSSatish Balay .  -tao_view_gradient
1682a7e14dcfSSatish Balay 
1683a7e14dcfSSatish Balay    Level: advanced
1684a7e14dcfSSatish Balay 
1685a7e14dcfSSatish Balay .seealso: TaoDefaultSMonitor(), TaoSetMonitor()
1686a7e14dcfSSatish Balay @*/
1687441846f8SBarry Smith PetscErrorCode TaoStepDirectionMonitor(Tao tao, void *ctx)
1688a7e14dcfSSatish Balay {
1689a7e14dcfSSatish Balay   PetscErrorCode ierr;
1690d393f493SBarry Smith   PetscViewer    viewer = (PetscViewer)ctx;
1691d393f493SBarry Smith 
1692a7e14dcfSSatish Balay   PetscFunctionBegin;
1693d393f493SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
1694a7e14dcfSSatish Balay   ierr = VecView(tao->stepdirection, viewer);CHKERRQ(ierr);
1695a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1696a7e14dcfSSatish Balay }
1697a7e14dcfSSatish Balay 
1698a7e14dcfSSatish Balay #undef __FUNCT__
1699a7e14dcfSSatish Balay #define __FUNCT__ "TaoDrawSolutionMonitor"
1700a7e14dcfSSatish Balay /*@C
1701a7e14dcfSSatish Balay    TaoDrawSolutionMonitor - Plots the solution at each iteration
1702a7e14dcfSSatish Balay    It can be turned on from the command line using the
1703a7e14dcfSSatish Balay    -tao_draw_solution option
1704a7e14dcfSSatish Balay 
1705441846f8SBarry Smith    Collective on Tao
1706a7e14dcfSSatish Balay 
1707a7e14dcfSSatish Balay    Input Parameters:
1708441846f8SBarry Smith +  tao - the Tao context
1709f55353a2SBarry Smith -  ctx - PetscViewer context
1710a7e14dcfSSatish Balay 
1711a7e14dcfSSatish Balay    Options Database Keys:
1712a7e14dcfSSatish Balay .  -tao_draw_solution
1713a7e14dcfSSatish Balay 
1714a7e14dcfSSatish Balay    Level: advanced
1715a7e14dcfSSatish Balay 
1716a7e14dcfSSatish Balay .seealso: TaoSolutionMonitor(), TaoSetMonitor(), TaoDrawGradientMonitor
1717a7e14dcfSSatish Balay @*/
1718441846f8SBarry Smith PetscErrorCode TaoDrawSolutionMonitor(Tao tao, void *ctx)
1719a7e14dcfSSatish Balay {
1720a7e14dcfSSatish Balay   PetscErrorCode ierr;
1721a7e14dcfSSatish Balay   PetscViewer    viewer = (PetscViewer) ctx;
172247a47007SBarry Smith 
1723a7e14dcfSSatish Balay   PetscFunctionBegin;
1724d393f493SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
1725a7e14dcfSSatish Balay   ierr = VecView(tao->solution, viewer);CHKERRQ(ierr);
1726a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1727a7e14dcfSSatish Balay }
1728a7e14dcfSSatish Balay 
1729a7e14dcfSSatish Balay #undef __FUNCT__
1730a7e14dcfSSatish Balay #define __FUNCT__ "TaoDrawGradientMonitor"
1731a7e14dcfSSatish Balay /*@C
1732a7e14dcfSSatish Balay    TaoDrawGradientMonitor - Plots the gradient at each iteration
1733a7e14dcfSSatish Balay    It can be turned on from the command line using the
1734a7e14dcfSSatish Balay    -tao_draw_gradient option
1735a7e14dcfSSatish Balay 
1736441846f8SBarry Smith    Collective on Tao
1737a7e14dcfSSatish Balay 
1738a7e14dcfSSatish Balay    Input Parameters:
1739441846f8SBarry Smith +  tao - the Tao context
1740f55353a2SBarry Smith -  ctx - PetscViewer context
1741a7e14dcfSSatish Balay 
1742a7e14dcfSSatish Balay    Options Database Keys:
1743a7e14dcfSSatish Balay .  -tao_draw_gradient
1744a7e14dcfSSatish Balay 
1745a7e14dcfSSatish Balay    Level: advanced
1746a7e14dcfSSatish Balay 
1747a7e14dcfSSatish Balay .seealso: TaoGradientMonitor(), TaoSetMonitor(), TaoDrawSolutionMonitor
1748a7e14dcfSSatish Balay @*/
1749441846f8SBarry Smith PetscErrorCode TaoDrawGradientMonitor(Tao tao, void *ctx)
1750a7e14dcfSSatish Balay {
1751a7e14dcfSSatish Balay   PetscErrorCode ierr;
1752a7e14dcfSSatish Balay   PetscViewer    viewer = (PetscViewer)ctx;
175347a47007SBarry Smith 
1754a7e14dcfSSatish Balay   PetscFunctionBegin;
1755d393f493SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
1756a7e14dcfSSatish Balay   ierr = VecView(tao->gradient, viewer);CHKERRQ(ierr);
1757a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1758a7e14dcfSSatish Balay }
1759a7e14dcfSSatish Balay 
1760a7e14dcfSSatish Balay #undef __FUNCT__
1761a7e14dcfSSatish Balay #define __FUNCT__ "TaoDrawStepMonitor"
1762a7e14dcfSSatish Balay /*@C
1763a7e14dcfSSatish Balay    TaoDrawStepMonitor - Plots the step direction at each iteration
1764a7e14dcfSSatish Balay    It can be turned on from the command line using the
1765a7e14dcfSSatish Balay    -tao_draw_step option
1766a7e14dcfSSatish Balay 
1767441846f8SBarry Smith    Collective on Tao
1768a7e14dcfSSatish Balay 
1769a7e14dcfSSatish Balay    Input Parameters:
1770441846f8SBarry Smith +  tao - the Tao context
1771f55353a2SBarry Smith -  ctx - PetscViewer context
1772a7e14dcfSSatish Balay 
1773a7e14dcfSSatish Balay    Options Database Keys:
1774a7e14dcfSSatish Balay .  -tao_draw_step
1775a7e14dcfSSatish Balay 
1776a7e14dcfSSatish Balay    Level: advanced
1777a7e14dcfSSatish Balay 
1778a7e14dcfSSatish Balay .seealso: TaoSetMonitor(), TaoDrawSolutionMonitor
1779a7e14dcfSSatish Balay @*/
1780441846f8SBarry Smith PetscErrorCode TaoDrawStepMonitor(Tao tao, void *ctx)
1781a7e14dcfSSatish Balay {
1782a7e14dcfSSatish Balay   PetscErrorCode ierr;
1783a7e14dcfSSatish Balay   PetscViewer    viewer = (PetscViewer)(ctx);
178447a47007SBarry Smith 
1785a7e14dcfSSatish Balay   PetscFunctionBegin;
1786a7e14dcfSSatish Balay   ierr = VecView(tao->stepdirection, viewer);CHKERRQ(ierr);
1787a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1788a7e14dcfSSatish Balay }
1789a7e14dcfSSatish Balay 
1790a7e14dcfSSatish Balay #undef __FUNCT__
1791a7e14dcfSSatish Balay #define __FUNCT__ "TaoSeparableObjectiveMonitor"
1792a7e14dcfSSatish Balay /*@C
1793a7e14dcfSSatish Balay    TaoSeparableObjectiveMonitor - Views the separable objective function at each iteration
1794a7e14dcfSSatish Balay    It can be turned on from the command line using the
1795a7e14dcfSSatish Balay    -tao_view_separableobjective option
1796a7e14dcfSSatish Balay 
1797441846f8SBarry Smith    Collective on Tao
1798a7e14dcfSSatish Balay 
1799a7e14dcfSSatish Balay    Input Parameters:
1800441846f8SBarry Smith +  tao - the Tao context
180147a47007SBarry Smith -  ctx - PetscViewer context or NULL
1802a7e14dcfSSatish Balay 
1803a7e14dcfSSatish Balay    Options Database Keys:
1804a7e14dcfSSatish Balay .  -tao_view_separableobjective
1805a7e14dcfSSatish Balay 
1806a7e14dcfSSatish Balay    Level: advanced
1807a7e14dcfSSatish Balay 
1808a7e14dcfSSatish Balay .seealso: TaoDefaultSMonitor(), TaoSetMonitor()
1809a7e14dcfSSatish Balay @*/
1810441846f8SBarry Smith PetscErrorCode TaoSeparableObjectiveMonitor(Tao tao, void *ctx)
1811a7e14dcfSSatish Balay {
1812a7e14dcfSSatish Balay   PetscErrorCode ierr;
1813d393f493SBarry Smith   PetscViewer    viewer  = (PetscViewer)ctx;
181447a47007SBarry Smith 
1815a7e14dcfSSatish Balay   PetscFunctionBegin;
1816d393f493SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
1817a7e14dcfSSatish Balay   ierr = VecView(tao->sep_objective,viewer);CHKERRQ(ierr);
1818a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1819a7e14dcfSSatish Balay }
1820a7e14dcfSSatish Balay 
1821a7e14dcfSSatish Balay #undef __FUNCT__
1822a7e14dcfSSatish Balay #define __FUNCT__ "TaoDefaultConvergenceTest"
18237fab98ebSJason Sarich /*@
1824a7e14dcfSSatish Balay    TaoDefaultConvergenceTest - Determines whether the solver should continue iterating
1825a7e14dcfSSatish Balay    or terminate.
1826a7e14dcfSSatish Balay 
1827441846f8SBarry Smith    Collective on Tao
1828a7e14dcfSSatish Balay 
1829a7e14dcfSSatish Balay    Input Parameters:
1830441846f8SBarry Smith +  tao - the Tao context
1831a7e14dcfSSatish Balay -  dummy - unused dummy context
1832a7e14dcfSSatish Balay 
1833a7e14dcfSSatish Balay    Output Parameter:
1834a7e14dcfSSatish Balay .  reason - for terminating
1835a7e14dcfSSatish Balay 
1836a7e14dcfSSatish Balay    Notes:
1837a7e14dcfSSatish Balay    This routine checks the residual in the optimality conditions, the
1838a7e14dcfSSatish Balay    relative residual in the optimity conditions, the number of function
1839a7e14dcfSSatish Balay    evaluations, and the function value to test convergence.  Some
1840a7e14dcfSSatish Balay    solvers may use different convergence routines.
1841a7e14dcfSSatish Balay 
1842a7e14dcfSSatish Balay    Level: developer
1843a7e14dcfSSatish Balay 
1844e4cb33bbSBarry Smith .seealso: TaoSetTolerances(),TaoGetConvergedReason(),TaoSetConvergedReason()
1845a7e14dcfSSatish Balay @*/
1846a7e14dcfSSatish Balay 
1847441846f8SBarry Smith PetscErrorCode TaoDefaultConvergenceTest(Tao tao,void *dummy)
1848a7e14dcfSSatish Balay {
1849a7e14dcfSSatish Balay   PetscInt           niter=tao->niter, nfuncs=PetscMax(tao->nfuncs,tao->nfuncgrads);
1850a7e14dcfSSatish Balay   PetscInt           max_funcs=tao->max_funcs;
1851a7e14dcfSSatish Balay   PetscReal          gnorm=tao->residual, gnorm0=tao->gnorm0;
1852a7e14dcfSSatish Balay   PetscReal          f=tao->fc, steptol=tao->steptol,trradius=tao->step;
1853a7e14dcfSSatish Balay   PetscReal          gatol=tao->gatol,grtol=tao->grtol,gttol=tao->gttol;
1854e52336cbSBarry Smith   PetscReal          catol=tao->catol,crtol=tao->crtol;
1855e52336cbSBarry Smith   PetscReal          fmin=tao->fmin, cnorm=tao->cnorm;
1856e4cb33bbSBarry Smith   TaoConvergedReason reason=tao->reason;
1857a7e14dcfSSatish Balay   PetscErrorCode     ierr;
1858a7e14dcfSSatish Balay 
1859a7e14dcfSSatish Balay   PetscFunctionBegin;
1860441846f8SBarry Smith   PetscValidHeaderSpecific(tao, TAO_CLASSID,1);
1861a7e14dcfSSatish Balay   if (reason != TAO_CONTINUE_ITERATING) {
1862a7e14dcfSSatish Balay     PetscFunctionReturn(0);
1863a7e14dcfSSatish Balay   }
1864a7e14dcfSSatish Balay 
1865a7e14dcfSSatish Balay   if (PetscIsInfOrNanReal(f)) {
1866a7e14dcfSSatish Balay     ierr = PetscInfo(tao,"Failed to converged, function value is Inf or NaN\n");CHKERRQ(ierr);
1867a7e14dcfSSatish Balay     reason = TAO_DIVERGED_NAN;
1868a7e14dcfSSatish Balay   } else if (f <= fmin && cnorm <=catol) {
186947a47007SBarry Smith     ierr = PetscInfo2(tao,"Converged due to function value %g < minimum function value %g\n", (double)f,(double)fmin);CHKERRQ(ierr);
1870a7e14dcfSSatish Balay     reason = TAO_CONVERGED_MINF;
1871a7e14dcfSSatish Balay   } else if (gnorm<= gatol && cnorm <=catol) {
187247a47007SBarry Smith     ierr = PetscInfo2(tao,"Converged due to residual norm ||g(X)||=%g < %g\n",(double)gnorm,(double)gatol);CHKERRQ(ierr);
1873a7e14dcfSSatish Balay     reason = TAO_CONVERGED_GATOL;
1874a7e14dcfSSatish Balay   } else if ( f!=0 && PetscAbsReal(gnorm/f) <= grtol && cnorm <= crtol) {
187547a47007SBarry Smith     ierr = PetscInfo2(tao,"Converged due to residual ||g(X)||/|f(X)| =%g < %g\n",(double)(gnorm/f),(double)grtol);CHKERRQ(ierr);
1876a7e14dcfSSatish Balay     reason = TAO_CONVERGED_GRTOL;
1877e1d16bb9SBarry Smith   } else if (gnorm0 != 0 && ((gttol == 0 && gnorm == 0) || gnorm/gnorm0 < gttol) && cnorm <= crtol) {
187847a47007SBarry Smith     ierr = PetscInfo2(tao,"Converged due to relative residual norm ||g(X)||/||g(X0)|| = %g < %g\n",(double)(gnorm/gnorm0),(double)gttol);CHKERRQ(ierr);
1879a7e14dcfSSatish Balay     reason = TAO_CONVERGED_GTTOL;
1880a7e14dcfSSatish Balay   } else if (nfuncs > max_funcs){
1881a7e14dcfSSatish Balay     ierr = PetscInfo2(tao,"Exceeded maximum number of function evaluations: %D > %D\n", nfuncs,max_funcs);CHKERRQ(ierr);
1882a7e14dcfSSatish Balay     reason = TAO_DIVERGED_MAXFCN;
1883a7e14dcfSSatish Balay   } else if ( tao->lsflag != 0 ){
1884a7e14dcfSSatish Balay     ierr = PetscInfo(tao,"Tao Line Search failure.\n");CHKERRQ(ierr);
1885a7e14dcfSSatish Balay     reason = TAO_DIVERGED_LS_FAILURE;
1886a7e14dcfSSatish Balay   } else if (trradius < steptol && niter > 0){
188747a47007SBarry Smith     ierr = PetscInfo2(tao,"Trust region/step size too small: %g < %g\n", (double)trradius,(double)steptol);CHKERRQ(ierr);
1888a7e14dcfSSatish Balay     reason = TAO_CONVERGED_STEPTOL;
1889a7e14dcfSSatish Balay   } else if (niter > tao->max_it) {
1890a7e14dcfSSatish Balay     ierr = PetscInfo2(tao,"Exceeded maximum number of iterations: %D > %D\n",niter,tao->max_it);CHKERRQ(ierr);
1891a7e14dcfSSatish Balay     reason = TAO_DIVERGED_MAXITS;
1892a7e14dcfSSatish Balay   } else {
1893a7e14dcfSSatish Balay     reason = TAO_CONTINUE_ITERATING;
1894a7e14dcfSSatish Balay   }
1895a7e14dcfSSatish Balay   tao->reason = reason;
1896a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1897a7e14dcfSSatish Balay }
1898a7e14dcfSSatish Balay 
1899a7e14dcfSSatish Balay #undef __FUNCT__
1900a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetOptionsPrefix"
1901a7e14dcfSSatish Balay /*@C
1902a7e14dcfSSatish Balay    TaoSetOptionsPrefix - Sets the prefix used for searching for all
1903a7e14dcfSSatish Balay    TAO options in the database.
1904a7e14dcfSSatish Balay 
1905a7e14dcfSSatish Balay 
1906441846f8SBarry Smith    Logically Collective on Tao
1907a7e14dcfSSatish Balay 
1908a7e14dcfSSatish Balay    Input Parameters:
1909441846f8SBarry Smith +  tao - the Tao context
1910a7e14dcfSSatish Balay -  prefix - the prefix string to prepend to all TAO option requests
1911a7e14dcfSSatish Balay 
1912a7e14dcfSSatish Balay    Notes:
1913a7e14dcfSSatish Balay    A hyphen (-) must NOT be given at the beginning of the prefix name.
1914a7e14dcfSSatish Balay    The first character of all runtime options is AUTOMATICALLY the hyphen.
1915a7e14dcfSSatish Balay 
1916a7e14dcfSSatish Balay    For example, to distinguish between the runtime options for two
1917a7e14dcfSSatish Balay    different TAO solvers, one could call
1918a7e14dcfSSatish Balay .vb
1919a7e14dcfSSatish Balay       TaoSetOptionsPrefix(tao1,"sys1_")
1920a7e14dcfSSatish Balay       TaoSetOptionsPrefix(tao2,"sys2_")
1921a7e14dcfSSatish Balay .ve
1922a7e14dcfSSatish Balay 
1923a7e14dcfSSatish Balay    This would enable use of different options for each system, such as
1924a7e14dcfSSatish Balay .vb
1925a7e14dcfSSatish Balay       -sys1_tao_method blmvm -sys1_tao_gtol 1.e-3
1926a7e14dcfSSatish Balay       -sys2_tao_method lmvm  -sys2_tao_gtol 1.e-4
1927a7e14dcfSSatish Balay .ve
1928a7e14dcfSSatish Balay 
1929a7e14dcfSSatish Balay 
1930a7e14dcfSSatish Balay    Level: advanced
1931a7e14dcfSSatish Balay 
1932a7e14dcfSSatish Balay .seealso: TaoAppendOptionsPrefix(), TaoGetOptionsPrefix()
1933a7e14dcfSSatish Balay @*/
1934a7e14dcfSSatish Balay 
1935441846f8SBarry Smith PetscErrorCode TaoSetOptionsPrefix(Tao tao, const char p[])
1936a7e14dcfSSatish Balay {
1937a7e14dcfSSatish Balay   PetscErrorCode ierr;
193847a47007SBarry Smith 
1939a7e14dcfSSatish Balay   PetscFunctionBegin;
1940a7e14dcfSSatish Balay   ierr = PetscObjectSetOptionsPrefix((PetscObject)tao,p);CHKERRQ(ierr);
1941a7e14dcfSSatish Balay   if (tao->linesearch) {
1942a7e14dcfSSatish Balay     ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,p);CHKERRQ(ierr);
1943a7e14dcfSSatish Balay   }
1944a7e14dcfSSatish Balay   if (tao->ksp) {
1945a7e14dcfSSatish Balay     ierr = KSPSetOptionsPrefix(tao->ksp,p);CHKERRQ(ierr);
1946a7e14dcfSSatish Balay   }
1947a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1948a7e14dcfSSatish Balay }
1949a7e14dcfSSatish Balay 
1950a7e14dcfSSatish Balay #undef __FUNCT__
1951a7e14dcfSSatish Balay #define __FUNCT__ "TaoAppendOptionsPrefix"
1952a7e14dcfSSatish Balay /*@C
1953a7e14dcfSSatish Balay    TaoAppendOptionsPrefix - Appends to the prefix used for searching for all
1954a7e14dcfSSatish Balay    TAO options in the database.
1955a7e14dcfSSatish Balay 
1956a7e14dcfSSatish Balay 
1957441846f8SBarry Smith    Logically Collective on Tao
1958a7e14dcfSSatish Balay 
1959a7e14dcfSSatish Balay    Input Parameters:
1960441846f8SBarry Smith +  tao - the Tao solver context
1961a7e14dcfSSatish Balay -  prefix - the prefix string to prepend to all TAO option requests
1962a7e14dcfSSatish Balay 
1963a7e14dcfSSatish Balay    Notes:
1964a7e14dcfSSatish Balay    A hyphen (-) must NOT be given at the beginning of the prefix name.
1965a7e14dcfSSatish Balay    The first character of all runtime options is AUTOMATICALLY the hyphen.
1966a7e14dcfSSatish Balay 
1967a7e14dcfSSatish Balay 
1968a7e14dcfSSatish Balay    Level: advanced
1969a7e14dcfSSatish Balay 
1970a7e14dcfSSatish Balay .seealso: TaoSetOptionsPrefix(), TaoGetOptionsPrefix()
1971a7e14dcfSSatish Balay @*/
1972441846f8SBarry Smith PetscErrorCode TaoAppendOptionsPrefix(Tao tao, const char p[])
1973a7e14dcfSSatish Balay {
1974a7e14dcfSSatish Balay   PetscErrorCode ierr;
197547a47007SBarry Smith 
1976a7e14dcfSSatish Balay   PetscFunctionBegin;
1977a7e14dcfSSatish Balay   ierr = PetscObjectAppendOptionsPrefix((PetscObject)tao,p);CHKERRQ(ierr);
1978a7e14dcfSSatish Balay   if (tao->linesearch) {
1979a7e14dcfSSatish Balay     ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,p);CHKERRQ(ierr);
1980a7e14dcfSSatish Balay   }
1981a7e14dcfSSatish Balay   if (tao->ksp) {
1982a7e14dcfSSatish Balay     ierr = KSPSetOptionsPrefix(tao->ksp,p);CHKERRQ(ierr);
1983a7e14dcfSSatish Balay   }
1984a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1985a7e14dcfSSatish Balay }
1986a7e14dcfSSatish Balay 
1987a7e14dcfSSatish Balay #undef __FUNCT__
1988a7e14dcfSSatish Balay #define __FUNCT__ "TaoGetOptionsPrefix"
1989a7e14dcfSSatish Balay /*@C
1990a7e14dcfSSatish Balay   TaoGetOptionsPrefix - Gets the prefix used for searching for all
1991a7e14dcfSSatish Balay   TAO options in the database
1992a7e14dcfSSatish Balay 
1993a7e14dcfSSatish Balay   Not Collective
1994a7e14dcfSSatish Balay 
1995a7e14dcfSSatish Balay   Input Parameters:
1996441846f8SBarry Smith . tao - the Tao context
1997a7e14dcfSSatish Balay 
1998a7e14dcfSSatish Balay   Output Parameters:
1999a7e14dcfSSatish Balay . prefix - pointer to the prefix string used is returned
2000a7e14dcfSSatish Balay 
2001a7e14dcfSSatish Balay   Notes: On the fortran side, the user should pass in a string 'prefix' of
2002a7e14dcfSSatish Balay   sufficient length to hold the prefix.
2003a7e14dcfSSatish Balay 
2004a7e14dcfSSatish Balay   Level: advanced
2005a7e14dcfSSatish Balay 
2006a7e14dcfSSatish Balay .seealso: TaoSetOptionsPrefix(), TaoAppendOptionsPrefix()
2007a7e14dcfSSatish Balay @*/
2008441846f8SBarry Smith PetscErrorCode TaoGetOptionsPrefix(Tao tao, const char *p[])
2009a7e14dcfSSatish Balay {
201047a47007SBarry Smith    return PetscObjectGetOptionsPrefix((PetscObject)tao,p);
2011a7e14dcfSSatish Balay }
2012a7e14dcfSSatish Balay 
2013a7e14dcfSSatish Balay #undef __FUNCT__
2014a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetType"
2015a7e14dcfSSatish Balay /*@C
2016a7e14dcfSSatish Balay    TaoSetType - Sets the method for the unconstrained minimization solver.
2017a7e14dcfSSatish Balay 
2018441846f8SBarry Smith    Collective on Tao
2019a7e14dcfSSatish Balay 
2020a7e14dcfSSatish Balay    Input Parameters:
2021441846f8SBarry Smith +  solver - the Tao solver context
2022a7e14dcfSSatish Balay -  type - a known method
2023a7e14dcfSSatish Balay 
2024a7e14dcfSSatish Balay    Options Database Key:
202558417fe7SBarry Smith .  -tao_type <type> - Sets the method; use -help for a list
202658417fe7SBarry Smith    of available methods (for instance, "-tao_type lmvm" or "-tao_type tron")
2027a7e14dcfSSatish Balay 
2028a7e14dcfSSatish Balay    Available methods include:
202958417fe7SBarry Smith +    nls - Newton's method with line search for unconstrained minimization
203058417fe7SBarry Smith .    ntr - Newton's method with trust region for unconstrained minimization
203158417fe7SBarry Smith .    ntl - Newton's method with trust region, line search for unconstrained minimization
203258417fe7SBarry Smith .    lmvm - Limited memory variable metric method for unconstrained minimization
203358417fe7SBarry Smith .    cg - Nonlinear conjugate gradient method for unconstrained minimization
203458417fe7SBarry Smith .    nm - Nelder-Mead algorithm for derivate-free unconstrained minimization
203558417fe7SBarry Smith .    tron - Newton Trust Region method for bound constrained minimization
203658417fe7SBarry Smith .    gpcg - Newton Trust Region method for quadratic bound constrained minimization
203758417fe7SBarry Smith .    blmvm - Limited memory variable metric method for bound constrained minimization
20381eb8069cSJason Sarich -    pounders - Model-based algorithm pounder extended for nonlinear least squares
2039a7e14dcfSSatish Balay 
2040a7e14dcfSSatish Balay   Level: intermediate
2041a7e14dcfSSatish Balay 
20421eb8069cSJason Sarich .seealso: TaoCreate(), TaoGetType(), TaoType
2043a7e14dcfSSatish Balay 
2044a7e14dcfSSatish Balay @*/
2045441846f8SBarry Smith PetscErrorCode TaoSetType(Tao tao, const TaoType type)
2046a7e14dcfSSatish Balay {
2047a7e14dcfSSatish Balay   PetscErrorCode ierr;
2048441846f8SBarry Smith   PetscErrorCode (*create_xxx)(Tao);
2049a7e14dcfSSatish Balay   PetscBool      issame;
2050a7e14dcfSSatish Balay 
2051a7e14dcfSSatish Balay   PetscFunctionBegin;
2052441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2053a7e14dcfSSatish Balay 
2054a7e14dcfSSatish Balay   ierr = PetscObjectTypeCompare((PetscObject)tao,type,&issame);CHKERRQ(ierr);
2055a7e14dcfSSatish Balay   if (issame) PetscFunctionReturn(0);
2056a7e14dcfSSatish Balay 
2057441846f8SBarry Smith   ierr = PetscFunctionListFind(TaoList, type, (void(**)(void))&create_xxx);CHKERRQ(ierr);
2058441846f8SBarry Smith   if (!create_xxx) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested Tao type %s",type);
2059a7e14dcfSSatish Balay 
2060a7e14dcfSSatish Balay   /* Destroy the existing solver information */
2061a7e14dcfSSatish Balay   if (tao->ops->destroy) {
2062a7e14dcfSSatish Balay     ierr = (*tao->ops->destroy)(tao);CHKERRQ(ierr);
2063a7e14dcfSSatish Balay   }
2064a7e14dcfSSatish Balay   ierr = KSPDestroy(&tao->ksp);CHKERRQ(ierr);
2065a7e14dcfSSatish Balay   ierr = TaoLineSearchDestroy(&tao->linesearch);CHKERRQ(ierr);
2066a7e14dcfSSatish Balay   ierr = VecDestroy(&tao->gradient);CHKERRQ(ierr);
2067a7e14dcfSSatish Balay   ierr = VecDestroy(&tao->stepdirection);CHKERRQ(ierr);
2068a7e14dcfSSatish Balay 
2069a7e14dcfSSatish Balay   tao->ops->setup = 0;
2070a7e14dcfSSatish Balay   tao->ops->solve = 0;
2071a7e14dcfSSatish Balay   tao->ops->view  = 0;
2072a7e14dcfSSatish Balay   tao->ops->setfromoptions = 0;
2073a7e14dcfSSatish Balay   tao->ops->destroy = 0;
2074a7e14dcfSSatish Balay 
2075a7e14dcfSSatish Balay   tao->setupcalled = PETSC_FALSE;
2076a7e14dcfSSatish Balay 
2077a7e14dcfSSatish Balay   ierr = (*create_xxx)(tao);CHKERRQ(ierr);
2078a7e14dcfSSatish Balay   ierr = PetscObjectChangeTypeName((PetscObject)tao,type);CHKERRQ(ierr);
2079a7e14dcfSSatish Balay   PetscFunctionReturn(0);
2080a7e14dcfSSatish Balay }
208147a47007SBarry Smith 
2082a7e14dcfSSatish Balay #undef __FUNCT__
2083441846f8SBarry Smith #define __FUNCT__ "TaoRegister"
2084a7e14dcfSSatish Balay /*MC
2085441846f8SBarry Smith    TaoRegister - Adds a method to the TAO package for unconstrained minimization.
2086a7e14dcfSSatish Balay 
2087a7e14dcfSSatish Balay    Synopsis:
2088441846f8SBarry Smith    TaoRegister(char *name_solver,char *path,char *name_Create,int (*routine_Create)(Tao))
2089a7e14dcfSSatish Balay 
2090a7e14dcfSSatish Balay    Not collective
2091a7e14dcfSSatish Balay 
2092a7e14dcfSSatish Balay    Input Parameters:
2093a7e14dcfSSatish Balay +  sname - name of a new user-defined solver
2094a7e14dcfSSatish Balay -  func - routine to Create method context
2095a7e14dcfSSatish Balay 
2096a7e14dcfSSatish Balay    Notes:
2097441846f8SBarry Smith    TaoRegister() may be called multiple times to add several user-defined solvers.
2098a7e14dcfSSatish Balay 
2099a7e14dcfSSatish Balay    Sample usage:
2100a7e14dcfSSatish Balay .vb
2101441846f8SBarry Smith    TaoRegister("my_solver",MySolverCreate);
2102a7e14dcfSSatish Balay .ve
2103a7e14dcfSSatish Balay 
2104a7e14dcfSSatish Balay    Then, your solver can be chosen with the procedural interface via
2105a7e14dcfSSatish Balay $     TaoSetType(tao,"my_solver")
2106a7e14dcfSSatish Balay    or at runtime via the option
210758417fe7SBarry Smith $     -tao_type my_solver
2108a7e14dcfSSatish Balay 
2109a7e14dcfSSatish Balay    Level: advanced
2110a7e14dcfSSatish Balay 
2111441846f8SBarry Smith .seealso: TaoRegisterAll(), TaoRegisterDestroy()
2112a7e14dcfSSatish Balay M*/
2113441846f8SBarry Smith PetscErrorCode TaoRegister(const char sname[], PetscErrorCode (*func)(Tao))
2114a7e14dcfSSatish Balay {
2115a7e14dcfSSatish Balay   PetscErrorCode ierr;
2116a7e14dcfSSatish Balay 
2117a7e14dcfSSatish Balay   PetscFunctionBegin;
2118441846f8SBarry Smith   ierr = PetscFunctionListAdd(&TaoList,sname, (void (*)(void))func);CHKERRQ(ierr);
2119a7e14dcfSSatish Balay   PetscFunctionReturn(0);
2120a7e14dcfSSatish Balay }
2121a7e14dcfSSatish Balay 
2122a7e14dcfSSatish Balay #undef __FUNCT__
2123441846f8SBarry Smith #define __FUNCT__ "TaoRegisterDestroy"
2124a7e14dcfSSatish Balay /*@C
2125441846f8SBarry Smith    TaoRegisterDestroy - Frees the list of minimization solvers that were
2126441846f8SBarry Smith    registered by TaoRegisterDynamic().
2127a7e14dcfSSatish Balay 
2128a7e14dcfSSatish Balay    Not Collective
2129a7e14dcfSSatish Balay 
2130a7e14dcfSSatish Balay    Level: advanced
2131a7e14dcfSSatish Balay 
2132441846f8SBarry Smith .seealso: TaoRegisterAll(), TaoRegister()
2133a7e14dcfSSatish Balay @*/
2134441846f8SBarry Smith PetscErrorCode TaoRegisterDestroy(void)
2135a7e14dcfSSatish Balay {
2136a7e14dcfSSatish Balay   PetscErrorCode ierr;
2137a7e14dcfSSatish Balay   PetscFunctionBegin;
2138441846f8SBarry Smith   ierr = PetscFunctionListDestroy(&TaoList);CHKERRQ(ierr);
2139441846f8SBarry Smith   TaoRegisterAllCalled = PETSC_FALSE;
2140a7e14dcfSSatish Balay   PetscFunctionReturn(0);
2141a7e14dcfSSatish Balay }
214247a47007SBarry Smith 
2143a7e14dcfSSatish Balay #undef __FUNCT__
21448931d482SJason Sarich #define __FUNCT__ "TaoGetIterationNumber"
21458931d482SJason Sarich /*@
21468931d482SJason Sarich    TaoGetIterationNumber - Gets the number of Tao iterations completed
21478931d482SJason Sarich    at this time.
21488931d482SJason Sarich 
21498931d482SJason Sarich    Not Collective
21508931d482SJason Sarich 
21518931d482SJason Sarich    Input Parameter:
21528931d482SJason Sarich .  tao - Tao context
21538931d482SJason Sarich 
21548931d482SJason Sarich    Output Parameter:
21558931d482SJason Sarich .  iter - iteration number
21568931d482SJason Sarich 
21578931d482SJason Sarich    Notes:
21588931d482SJason Sarich    For example, during the computation of iteration 2 this would return 1.
21598931d482SJason Sarich 
21608931d482SJason Sarich 
21618931d482SJason Sarich    Level: intermediate
21628931d482SJason Sarich 
21638931d482SJason Sarich .keywords: Tao, nonlinear, get, iteration, number,
21648931d482SJason Sarich 
21658931d482SJason Sarich .seealso:   TaoGetLinearSolveIterations()
21668931d482SJason Sarich @*/
21678931d482SJason Sarich PetscErrorCode  TaoGetIterationNumber(Tao tao,PetscInt *iter)
21688931d482SJason Sarich {
21698931d482SJason Sarich   PetscFunctionBegin;
21708931d482SJason Sarich   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
21718931d482SJason Sarich   PetscValidIntPointer(iter,2);
21728931d482SJason Sarich   *iter = tao->niter;
21738931d482SJason Sarich   PetscFunctionReturn(0);
21748931d482SJason Sarich }
21758931d482SJason Sarich 
21768931d482SJason Sarich #undef __FUNCT__
21778931d482SJason Sarich #define __FUNCT__ "TaoSetIterationNumber"
21788931d482SJason Sarich /*@
21798931d482SJason Sarich    TaoSetIterationNumber - Sets the current iteration number.
21808931d482SJason Sarich 
21818931d482SJason Sarich    Not Collective
21828931d482SJason Sarich 
21838931d482SJason Sarich    Input Parameter:
21848931d482SJason Sarich .  tao - Tao context
21858931d482SJason Sarich .  iter - iteration number
21868931d482SJason Sarich 
21878931d482SJason Sarich    Level: developer
21888931d482SJason Sarich 
21898931d482SJason Sarich .keywords: Tao, nonlinear, set, iteration, number,
21908931d482SJason Sarich 
21918931d482SJason Sarich .seealso:   TaoGetLinearSolveIterations()
21928931d482SJason Sarich @*/
21938931d482SJason Sarich PetscErrorCode  TaoSetIterationNumber(Tao tao,PetscInt iter)
21948931d482SJason Sarich {
21958931d482SJason Sarich   PetscErrorCode ierr;
21968931d482SJason Sarich 
21978931d482SJason Sarich   PetscFunctionBegin;
21988931d482SJason Sarich   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
21998931d482SJason Sarich   ierr       = PetscObjectSAWsTakeAccess((PetscObject)tao);CHKERRQ(ierr);
22008931d482SJason Sarich   tao->niter = iter;
22018931d482SJason Sarich   ierr       = PetscObjectSAWsGrantAccess((PetscObject)tao);CHKERRQ(ierr);
22028931d482SJason Sarich   PetscFunctionReturn(0);
22038931d482SJason Sarich }
22048931d482SJason Sarich 
22058931d482SJason Sarich #undef __FUNCT__
22068931d482SJason Sarich #define __FUNCT__ "TaoGetTotalIterationNumber"
22078931d482SJason Sarich /*@
22088931d482SJason Sarich    TaoGetTotalIterationNumber - Gets the total number of Tao iterations
22098931d482SJason Sarich    completed. This number keeps accumulating if multiple solves
22108931d482SJason Sarich    are called with the Tao object.
22118931d482SJason Sarich 
22128931d482SJason Sarich    Not Collective
22138931d482SJason Sarich 
22148931d482SJason Sarich    Input Parameter:
22158931d482SJason Sarich .  tao - Tao context
22168931d482SJason Sarich 
22178931d482SJason Sarich    Output Parameter:
22188931d482SJason Sarich .  iter - iteration number
22198931d482SJason Sarich 
22208931d482SJason Sarich    Notes:
22218931d482SJason Sarich    The total iteration count is updated after each solve, if there is a current
22228931d482SJason Sarich    TaoSolve() in progress then those iterations are not yet counted.
22238931d482SJason Sarich 
22248931d482SJason Sarich    Level: intermediate
22258931d482SJason Sarich 
22268931d482SJason Sarich .keywords: Tao, nonlinear, get, iteration, number,
22278931d482SJason Sarich 
22288931d482SJason Sarich .seealso:   TaoGetLinearSolveIterations()
22298931d482SJason Sarich @*/
22308931d482SJason Sarich PetscErrorCode  TaoGetTotalIterationNumber(Tao tao,PetscInt *iter)
22318931d482SJason Sarich {
22328931d482SJason Sarich   PetscFunctionBegin;
22338931d482SJason Sarich   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
22348931d482SJason Sarich   PetscValidIntPointer(iter,2);
22358931d482SJason Sarich   *iter = tao->ntotalits;
22368931d482SJason Sarich   PetscFunctionReturn(0);
22378931d482SJason Sarich }
22388931d482SJason Sarich 
22398931d482SJason Sarich #undef __FUNCT__
22408931d482SJason Sarich #define __FUNCT__ "TaoSetTotalIterationNumber"
22418931d482SJason Sarich /*@
22428931d482SJason Sarich    TaoSetTotalIterationNumber - Sets the current total iteration number.
22438931d482SJason Sarich 
22448931d482SJason Sarich    Not Collective
22458931d482SJason Sarich 
22468931d482SJason Sarich    Input Parameter:
22478931d482SJason Sarich .  tao - Tao context
22488931d482SJason Sarich .  iter - iteration number
22498931d482SJason Sarich 
22508931d482SJason Sarich    Level: developer
22518931d482SJason Sarich 
22528931d482SJason Sarich .keywords: Tao, nonlinear, set, iteration, number,
22538931d482SJason Sarich 
22548931d482SJason Sarich .seealso:   TaoGetLinearSolveIterations()
22558931d482SJason Sarich @*/
22568931d482SJason Sarich PetscErrorCode  TaoSetTotalIterationNumber(Tao tao,PetscInt iter)
22578931d482SJason Sarich {
22588931d482SJason Sarich   PetscErrorCode ierr;
22598931d482SJason Sarich 
22608931d482SJason Sarich   PetscFunctionBegin;
22618931d482SJason Sarich   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
22628931d482SJason Sarich   ierr       = PetscObjectSAWsTakeAccess((PetscObject)tao);CHKERRQ(ierr);
22638931d482SJason Sarich   tao->ntotalits = iter;
22648931d482SJason Sarich   ierr       = PetscObjectSAWsGrantAccess((PetscObject)tao);CHKERRQ(ierr);
22658931d482SJason Sarich   PetscFunctionReturn(0);
22668931d482SJason Sarich }
22678931d482SJason Sarich 
22688931d482SJason Sarich #undef __FUNCT__
2269e4cb33bbSBarry Smith #define __FUNCT__ "TaoSetConvergedReason"
2270a7e14dcfSSatish Balay /*@
2271e4cb33bbSBarry Smith   TaoSetConvergedReason - Sets the termination flag on a Tao object
2272a7e14dcfSSatish Balay 
2273441846f8SBarry Smith   Logically Collective on Tao
2274a7e14dcfSSatish Balay 
2275a7e14dcfSSatish Balay   Input Parameters:
2276441846f8SBarry Smith + tao - the Tao context
2277a7e14dcfSSatish Balay - reason - one of
2278a7e14dcfSSatish Balay $     TAO_CONVERGED_ATOL (2),
2279a7e14dcfSSatish Balay $     TAO_CONVERGED_RTOL (3),
2280a7e14dcfSSatish Balay $     TAO_CONVERGED_STEPTOL (4),
2281a7e14dcfSSatish Balay $     TAO_CONVERGED_MINF (5),
2282a7e14dcfSSatish Balay $     TAO_CONVERGED_USER (6),
2283a7e14dcfSSatish Balay $     TAO_DIVERGED_MAXITS (-2),
2284a7e14dcfSSatish Balay $     TAO_DIVERGED_NAN (-4),
2285a7e14dcfSSatish Balay $     TAO_DIVERGED_MAXFCN (-5),
2286a7e14dcfSSatish Balay $     TAO_DIVERGED_LS_FAILURE (-6),
2287a7e14dcfSSatish Balay $     TAO_DIVERGED_TR_REDUCTION (-7),
2288a7e14dcfSSatish Balay $     TAO_DIVERGED_USER (-8),
2289a7e14dcfSSatish Balay $     TAO_CONTINUE_ITERATING (0)
2290a7e14dcfSSatish Balay 
2291a7e14dcfSSatish Balay    Level: intermediate
2292a7e14dcfSSatish Balay 
2293a7e14dcfSSatish Balay @*/
2294e4cb33bbSBarry Smith PetscErrorCode TaoSetConvergedReason(Tao tao, TaoConvergedReason reason)
2295a7e14dcfSSatish Balay {
2296441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2297a7e14dcfSSatish Balay   PetscFunctionBegin;
2298a7e14dcfSSatish Balay   tao->reason = reason;
2299a7e14dcfSSatish Balay   PetscFunctionReturn(0);
2300a7e14dcfSSatish Balay }
2301a7e14dcfSSatish Balay 
2302a7e14dcfSSatish Balay #undef __FUNCT__
2303e4cb33bbSBarry Smith #define __FUNCT__ "TaoGetConvergedReason"
2304a7e14dcfSSatish Balay /*@
2305e4cb33bbSBarry Smith    TaoGetConvergedReason - Gets the reason the Tao iteration was stopped.
2306a7e14dcfSSatish Balay 
2307a7e14dcfSSatish Balay    Not Collective
2308a7e14dcfSSatish Balay 
2309a7e14dcfSSatish Balay    Input Parameter:
2310441846f8SBarry Smith .  tao - the Tao solver context
2311a7e14dcfSSatish Balay 
2312a7e14dcfSSatish Balay    Output Parameter:
2313a7e14dcfSSatish Balay .  reason - one of
2314a7e14dcfSSatish Balay $  TAO_CONVERGED_GATOL (3)           ||g(X)|| < gatol
2315a7e14dcfSSatish Balay $  TAO_CONVERGED_GRTOL (4)           ||g(X)|| / f(X)  < grtol
2316a7e14dcfSSatish Balay $  TAO_CONVERGED_GTTOL (5)           ||g(X)|| / ||g(X0)|| < gttol
2317a7e14dcfSSatish Balay $  TAO_CONVERGED_STEPTOL (6)         step size small
2318a7e14dcfSSatish Balay $  TAO_CONVERGED_MINF (7)            F < F_min
2319a7e14dcfSSatish Balay $  TAO_CONVERGED_USER (8)            User defined
2320a7e14dcfSSatish Balay $  TAO_DIVERGED_MAXITS (-2)          its > maxits
2321a7e14dcfSSatish Balay $  TAO_DIVERGED_NAN (-4)             Numerical problems
2322a7e14dcfSSatish Balay $  TAO_DIVERGED_MAXFCN (-5)          fevals > max_funcsals
2323a7e14dcfSSatish Balay $  TAO_DIVERGED_LS_FAILURE (-6)      line search failure
2324a7e14dcfSSatish Balay $  TAO_DIVERGED_TR_REDUCTION (-7)    trust region failure
2325a7e14dcfSSatish Balay $  TAO_DIVERGED_USER(-8)             (user defined)
2326a7e14dcfSSatish Balay  $  TAO_CONTINUE_ITERATING (0)
2327a7e14dcfSSatish Balay 
2328a7e14dcfSSatish Balay    where
2329a7e14dcfSSatish Balay +  X - current solution
2330a7e14dcfSSatish Balay .  X0 - initial guess
2331a7e14dcfSSatish Balay .  f(X) - current function value
2332a7e14dcfSSatish Balay .  f(X*) - true solution (estimated)
2333a7e14dcfSSatish Balay .  g(X) - current gradient
2334a7e14dcfSSatish Balay .  its - current iterate number
2335a7e14dcfSSatish Balay .  maxits - maximum number of iterates
2336a7e14dcfSSatish Balay .  fevals - number of function evaluations
2337a7e14dcfSSatish Balay -  max_funcsals - maximum number of function evaluations
2338a7e14dcfSSatish Balay 
2339a7e14dcfSSatish Balay    Level: intermediate
2340a7e14dcfSSatish Balay 
2341a7e14dcfSSatish Balay .seealso: TaoSetConvergenceTest(), TaoSetTolerances()
2342a7e14dcfSSatish Balay 
2343a7e14dcfSSatish Balay @*/
2344e4cb33bbSBarry Smith PetscErrorCode TaoGetConvergedReason(Tao tao, TaoConvergedReason *reason)
2345a7e14dcfSSatish Balay {
2346a7e14dcfSSatish Balay   PetscFunctionBegin;
2347441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2348a7e14dcfSSatish Balay   PetscValidPointer(reason,2);
2349a7e14dcfSSatish Balay   *reason = tao->reason;
2350a7e14dcfSSatish Balay   PetscFunctionReturn(0);
2351a7e14dcfSSatish Balay }
2352a7e14dcfSSatish Balay 
2353a7e14dcfSSatish Balay #undef __FUNCT__
2354a7e14dcfSSatish Balay #define __FUNCT__ "TaoGetSolutionStatus"
2355a7e14dcfSSatish Balay /*@
2356a7e14dcfSSatish Balay   TaoGetSolutionStatus - Get the current iterate, objective value,
2357a7e14dcfSSatish Balay   residual, infeasibility, and termination
2358a7e14dcfSSatish Balay 
2359a7e14dcfSSatish Balay   Not Collective
2360a7e14dcfSSatish Balay 
2361a7e14dcfSSatish Balay    Input Parameters:
2362441846f8SBarry Smith .  tao - the Tao context
2363a7e14dcfSSatish Balay 
2364a7e14dcfSSatish Balay    Output Parameters:
2365a7e14dcfSSatish Balay +  iterate - the current iterate number (>=0)
2366a7e14dcfSSatish Balay .  f - the current function value
236758417fe7SBarry Smith .  gnorm - the square of the gradient norm, duality gap, or other measure indicating distance from optimality.
2368a7e14dcfSSatish Balay .  cnorm - the infeasibility of the current solution with regard to the constraints.
2369a7e14dcfSSatish Balay .  xdiff - the step length or trust region radius of the most recent iterate.
2370a7e14dcfSSatish Balay -  reason - The termination reason, which can equal TAO_CONTINUE_ITERATING
2371a7e14dcfSSatish Balay 
2372a7e14dcfSSatish Balay    Level: intermediate
2373a7e14dcfSSatish Balay 
2374a7e14dcfSSatish Balay    Note:
2375a7e14dcfSSatish Balay    TAO returns the values set by the solvers in the routine TaoMonitor().
2376a7e14dcfSSatish Balay 
2377a7e14dcfSSatish Balay    Note:
237847a47007SBarry Smith    If any of the output arguments are set to NULL, no corresponding value will be returned.
2379a7e14dcfSSatish Balay 
2380e4cb33bbSBarry Smith .seealso: TaoMonitor(), TaoGetConvergedReason()
2381a7e14dcfSSatish Balay @*/
2382e4cb33bbSBarry Smith PetscErrorCode TaoGetSolutionStatus(Tao tao, PetscInt *its, PetscReal *f, PetscReal *gnorm, PetscReal *cnorm, PetscReal *xdiff, TaoConvergedReason *reason)
2383a7e14dcfSSatish Balay {
2384a7e14dcfSSatish Balay   PetscFunctionBegin;
2385a7e14dcfSSatish Balay   if (its) *its=tao->niter;
2386a7e14dcfSSatish Balay   if (f) *f=tao->fc;
2387a7e14dcfSSatish Balay   if (gnorm) *gnorm=tao->residual;
2388a7e14dcfSSatish Balay   if (cnorm) *cnorm=tao->cnorm;
2389a7e14dcfSSatish Balay   if (reason) *reason=tao->reason;
2390a7e14dcfSSatish Balay   if (xdiff) *xdiff=tao->step;
2391a7e14dcfSSatish Balay   PetscFunctionReturn(0);
2392a7e14dcfSSatish Balay }
2393a7e14dcfSSatish Balay 
2394a7e14dcfSSatish Balay #undef __FUNCT__
2395a7e14dcfSSatish Balay #define __FUNCT__ "TaoGetType"
2396a7e14dcfSSatish Balay /*@C
2397441846f8SBarry Smith    TaoGetType - Gets the current Tao algorithm.
2398a7e14dcfSSatish Balay 
2399a7e14dcfSSatish Balay    Not Collective
2400a7e14dcfSSatish Balay 
2401a7e14dcfSSatish Balay    Input Parameter:
2402441846f8SBarry Smith .  tao - the Tao solver context
2403a7e14dcfSSatish Balay 
2404a7e14dcfSSatish Balay    Output Parameter:
2405441846f8SBarry Smith .  type - Tao method
2406a7e14dcfSSatish Balay 
2407a7e14dcfSSatish Balay    Level: intermediate
2408a7e14dcfSSatish Balay 
2409a7e14dcfSSatish Balay @*/
2410441846f8SBarry Smith PetscErrorCode TaoGetType(Tao tao, const TaoType *type)
2411a7e14dcfSSatish Balay {
2412a7e14dcfSSatish Balay   PetscFunctionBegin;
2413441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2414a7e14dcfSSatish Balay   PetscValidPointer(type,2);
2415a7e14dcfSSatish Balay   *type=((PetscObject)tao)->type_name;
2416a7e14dcfSSatish Balay   PetscFunctionReturn(0);
2417a7e14dcfSSatish Balay }
2418a7e14dcfSSatish Balay 
2419a7e14dcfSSatish Balay #undef __FUNCT__
2420a7e14dcfSSatish Balay #define __FUNCT__ "TaoMonitor"
2421a7e14dcfSSatish Balay /*@C
2422a7e14dcfSSatish Balay   TaoMonitor - Monitor the solver and the current solution.  This
2423a7e14dcfSSatish Balay   routine will record the iteration number and residual statistics,
2424a7e14dcfSSatish Balay   call any monitors specified by the user, and calls the convergence-check routine.
2425a7e14dcfSSatish Balay 
2426a7e14dcfSSatish Balay    Input Parameters:
2427441846f8SBarry Smith +  tao - the Tao context
2428a7e14dcfSSatish Balay .  its - the current iterate number (>=0)
2429a7e14dcfSSatish Balay .  f - the current objective function value
243058417fe7SBarry Smith .  res - the gradient norm, square root of the duality gap, or other measure indicating distince from optimality.  This measure will be recorded and
2431a7e14dcfSSatish Balay           used for some termination tests.
2432a7e14dcfSSatish Balay .  cnorm - the infeasibility of the current solution with regard to the constraints.
2433a7e14dcfSSatish Balay -  steplength - multiple of the step direction added to the previous iterate.
2434a7e14dcfSSatish Balay 
2435a7e14dcfSSatish Balay    Output Parameters:
2436a7e14dcfSSatish Balay .  reason - The termination reason, which can equal TAO_CONTINUE_ITERATING
2437a7e14dcfSSatish Balay 
2438a7e14dcfSSatish Balay    Options Database Key:
2439a7e14dcfSSatish Balay .  -tao_monitor - Use the default monitor, which prints statistics to standard output
2440a7e14dcfSSatish Balay 
2441e4cb33bbSBarry Smith .seealso TaoGetConvergedReason(), TaoDefaultMonitor(), TaoSetMonitor()
2442a7e14dcfSSatish Balay 
2443a7e14dcfSSatish Balay    Level: developer
2444a7e14dcfSSatish Balay 
2445a7e14dcfSSatish Balay @*/
2446e4cb33bbSBarry Smith PetscErrorCode TaoMonitor(Tao tao, PetscInt its, PetscReal f, PetscReal res, PetscReal cnorm, PetscReal steplength, TaoConvergedReason *reason)
2447a7e14dcfSSatish Balay {
2448a7e14dcfSSatish Balay   PetscErrorCode ierr;
244947a47007SBarry Smith   PetscInt       i;
245047a47007SBarry Smith 
2451a7e14dcfSSatish Balay   PetscFunctionBegin;
2452441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2453a7e14dcfSSatish Balay   tao->fc = f;
2454a7e14dcfSSatish Balay   tao->residual = res;
2455a7e14dcfSSatish Balay   tao->cnorm = cnorm;
2456a7e14dcfSSatish Balay   tao->step = steplength;
2457e52336cbSBarry Smith   if (!its) {
2458a7e14dcfSSatish Balay     tao->cnorm0 = cnorm; tao->gnorm0 = res;
2459a7e14dcfSSatish Balay   }
2460ae93cb3cSJason Sarich   TaoLogConvergenceHistory(tao,f,res,cnorm,tao->ksp_its);
246147a47007SBarry Smith   if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(res)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
2462a7e14dcfSSatish Balay   if (tao->ops->convergencetest) {
2463a7e14dcfSSatish Balay     ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
2464a7e14dcfSSatish Balay   }
2465a7e14dcfSSatish Balay   for (i=0;i<tao->numbermonitors;i++) {
2466a7e14dcfSSatish Balay     ierr = (*tao->monitor[i])(tao,tao->monitorcontext[i]);CHKERRQ(ierr);
2467a7e14dcfSSatish Balay   }
2468a7e14dcfSSatish Balay   *reason = tao->reason;
2469a7e14dcfSSatish Balay   PetscFunctionReturn(0);
2470a7e14dcfSSatish Balay }
2471a7e14dcfSSatish Balay 
2472a7e14dcfSSatish Balay #undef __FUNCT__
2473ae93cb3cSJason Sarich #define __FUNCT__ "TaoSetConvergenceHistory"
2474a7e14dcfSSatish Balay /*@
2475ae93cb3cSJason Sarich    TaoSetConvergenceHistory - Sets the array used to hold the convergence history.
2476a7e14dcfSSatish Balay 
2477441846f8SBarry Smith    Logically Collective on Tao
2478a7e14dcfSSatish Balay 
2479a7e14dcfSSatish Balay    Input Parameters:
2480441846f8SBarry Smith +  tao - the Tao solver context
2481a7e14dcfSSatish Balay .  obj   - array to hold objective value history
2482a7e14dcfSSatish Balay .  resid - array to hold residual history
2483a7e14dcfSSatish Balay .  cnorm - array to hold constraint violation history
2484be1558d8SJason Sarich .  lits - integer array holds the number of linear iterations for each Tao iteration
2485a7e14dcfSSatish Balay .  na  - size of obj, resid, and cnorm
2486a7e14dcfSSatish Balay -  reset - PetscTrue indicates each new minimization resets the history counter to zero,
2487a7e14dcfSSatish Balay            else it continues storing new values for new minimizations after the old ones
2488a7e14dcfSSatish Balay 
2489a7e14dcfSSatish Balay    Notes:
2490a7e14dcfSSatish Balay    If set, TAO will fill the given arrays with the indicated
2491ae93cb3cSJason Sarich    information at each iteration.  If 'obj','resid','cnorm','lits' are
2492ae93cb3cSJason Sarich    *all* NULL then space (using size na, or 1000 if na is PETSC_DECIDE or
2493ae93cb3cSJason Sarich    PETSC_DEFAULT) is allocated for the history.
2494be1558d8SJason Sarich    If not all are NULL, then only the non-NULL information categories
2495ae93cb3cSJason Sarich    will be stored, the others will be ignored.
2496ae93cb3cSJason Sarich 
2497ae93cb3cSJason Sarich    Any convergence information after iteration number 'na' will not be stored.
2498a7e14dcfSSatish Balay 
2499a7e14dcfSSatish Balay    This routine is useful, e.g., when running a code for purposes
2500a7e14dcfSSatish Balay    of accurate performance monitoring, when no I/O should be done
2501a7e14dcfSSatish Balay    during the section of code that is being timed.
2502a7e14dcfSSatish Balay 
2503a7e14dcfSSatish Balay    Level: intermediate
2504a7e14dcfSSatish Balay 
2505ae93cb3cSJason Sarich .seealso: TaoGetConvergenceHistory()
2506a7e14dcfSSatish Balay 
2507a7e14dcfSSatish Balay @*/
2508ae93cb3cSJason Sarich PetscErrorCode TaoSetConvergenceHistory(Tao tao, PetscReal *obj, PetscReal *resid, PetscReal *cnorm, PetscInt *lits, PetscInt na,PetscBool reset)
2509a7e14dcfSSatish Balay {
2510ae93cb3cSJason Sarich   PetscErrorCode ierr;
2511ae93cb3cSJason Sarich 
2512a7e14dcfSSatish Balay   PetscFunctionBegin;
2513441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2514ae93cb3cSJason Sarich   if (obj) PetscValidScalarPointer(obj,2);
2515ae93cb3cSJason Sarich   if (resid) PetscValidScalarPointer(resid,3);
2516ae93cb3cSJason Sarich   if (cnorm) PetscValidScalarPointer(cnorm,4);
2517ae93cb3cSJason Sarich   if (lits) PetscValidIntPointer(lits,5);
2518ae93cb3cSJason Sarich 
2519ae93cb3cSJason Sarich   if (na == PETSC_DECIDE || na == PETSC_DEFAULT) na = 1000;
2520ae93cb3cSJason Sarich   if (!obj && !resid && !cnorm && !lits) {
2521ae93cb3cSJason Sarich     ierr = PetscCalloc1(na,&obj);CHKERRQ(ierr);
2522ae93cb3cSJason Sarich     ierr = PetscCalloc1(na,&resid);CHKERRQ(ierr);
2523ae93cb3cSJason Sarich     ierr = PetscCalloc1(na,&cnorm);CHKERRQ(ierr);
2524ae93cb3cSJason Sarich     ierr = PetscCalloc1(na,&lits);CHKERRQ(ierr);
2525ae93cb3cSJason Sarich     tao->hist_malloc=PETSC_TRUE;
2526ae93cb3cSJason Sarich   }
2527ae93cb3cSJason Sarich 
2528a7e14dcfSSatish Balay   tao->hist_obj = obj;
2529a7e14dcfSSatish Balay   tao->hist_resid = resid;
2530a7e14dcfSSatish Balay   tao->hist_cnorm = cnorm;
2531ae93cb3cSJason Sarich   tao->hist_lits = lits;
2532a7e14dcfSSatish Balay   tao->hist_max   = na;
2533a7e14dcfSSatish Balay   tao->hist_reset = reset;
2534ae93cb3cSJason Sarich   tao->hist_len = 0;
2535a7e14dcfSSatish Balay   PetscFunctionReturn(0);
2536a7e14dcfSSatish Balay }
2537a7e14dcfSSatish Balay 
2538a7e14dcfSSatish Balay #undef __FUNCT__
2539ae93cb3cSJason Sarich #define __FUNCT__ "TaoGetConvergenceHistory"
2540a7e14dcfSSatish Balay /*@C
2541ae93cb3cSJason Sarich    TaoGetConvergenceHistory - Gets the arrays used to hold the convergence history.
2542a7e14dcfSSatish Balay 
2543441846f8SBarry Smith    Collective on Tao
2544a7e14dcfSSatish Balay 
2545a7e14dcfSSatish Balay    Input Parameter:
2546441846f8SBarry Smith .  tao - the Tao context
2547a7e14dcfSSatish Balay 
2548a7e14dcfSSatish Balay    Output Parameters:
2549a7e14dcfSSatish Balay +  obj   - array used to hold objective value history
2550a7e14dcfSSatish Balay .  resid - array used to hold residual history
2551a7e14dcfSSatish Balay .  cnorm - array used to hold constraint violation history
2552ae93cb3cSJason Sarich .  lits  - integer array used to hold linear solver iteration count
2553ae93cb3cSJason Sarich -  nhist  - size of obj, resid, cnorm, and lits (will be less than or equal to na given in TaoSetHistory)
2554a7e14dcfSSatish Balay 
2555a7e14dcfSSatish Balay    Notes:
2556ae93cb3cSJason Sarich     This routine must be preceded by calls to TaoSetConvergenceHistory()
2557ae93cb3cSJason Sarich     and TaoSolve(), otherwise it returns useless information.
2558ae93cb3cSJason Sarich 
2559a7e14dcfSSatish Balay     The calling sequence for this routine in Fortran is
2560ae93cb3cSJason Sarich $   call TaoGetConvergenceHistory(Tao tao, PetscInt nhist, PetscErrorCode ierr)
2561a7e14dcfSSatish Balay 
2562a7e14dcfSSatish Balay    This routine is useful, e.g., when running a code for purposes
2563a7e14dcfSSatish Balay    of accurate performance monitoring, when no I/O should be done
2564a7e14dcfSSatish Balay    during the section of code that is being timed.
2565a7e14dcfSSatish Balay 
2566a7e14dcfSSatish Balay    Level: advanced
2567a7e14dcfSSatish Balay 
2568ae93cb3cSJason Sarich .seealso: TaoSetConvergenceHistory()
2569a7e14dcfSSatish Balay 
2570a7e14dcfSSatish Balay @*/
2571ae93cb3cSJason Sarich PetscErrorCode TaoGetConvergenceHistory(Tao tao, PetscReal **obj, PetscReal **resid, PetscReal **cnorm, PetscInt **lits, PetscInt *nhist)
2572a7e14dcfSSatish Balay {
2573a7e14dcfSSatish Balay   PetscFunctionBegin;
2574441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2575a7e14dcfSSatish Balay   if (obj)   *obj   = tao->hist_obj;
2576a7e14dcfSSatish Balay   if (cnorm) *cnorm = tao->hist_cnorm;
2577a7e14dcfSSatish Balay   if (resid) *resid = tao->hist_resid;
2578a7e14dcfSSatish Balay   if (nhist) *nhist   = tao->hist_len;
2579a7e14dcfSSatish Balay   PetscFunctionReturn(0);
2580a7e14dcfSSatish Balay }
2581a7e14dcfSSatish Balay 
2582a7e14dcfSSatish Balay #undef __FUNCT__
2583a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetApplicationContext"
2584a7e14dcfSSatish Balay /*@
2585a7e14dcfSSatish Balay    TaoSetApplicationContext - Sets the optional user-defined context for
2586a7e14dcfSSatish Balay    a solver.
2587a7e14dcfSSatish Balay 
2588441846f8SBarry Smith    Logically Collective on Tao
2589a7e14dcfSSatish Balay 
2590a7e14dcfSSatish Balay    Input Parameters:
2591441846f8SBarry Smith +  tao  - the Tao context
2592a7e14dcfSSatish Balay -  usrP - optional user context
2593a7e14dcfSSatish Balay 
2594a7e14dcfSSatish Balay    Level: intermediate
2595a7e14dcfSSatish Balay 
2596a7e14dcfSSatish Balay .seealso: TaoGetApplicationContext(), TaoSetApplicationContext()
2597a7e14dcfSSatish Balay @*/
2598441846f8SBarry Smith PetscErrorCode  TaoSetApplicationContext(Tao tao,void *usrP)
2599a7e14dcfSSatish Balay {
2600a7e14dcfSSatish Balay   PetscFunctionBegin;
2601441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2602a7e14dcfSSatish Balay   tao->user = usrP;
2603a7e14dcfSSatish Balay   PetscFunctionReturn(0);
2604a7e14dcfSSatish Balay }
2605a7e14dcfSSatish Balay 
2606a7e14dcfSSatish Balay #undef __FUNCT__
2607a7e14dcfSSatish Balay #define __FUNCT__ "TaoGetApplicationContext"
2608a7e14dcfSSatish Balay /*@
2609a7e14dcfSSatish Balay    TaoGetApplicationContext - Gets the user-defined context for a
2610a7e14dcfSSatish Balay    TAO solvers.
2611a7e14dcfSSatish Balay 
2612a7e14dcfSSatish Balay    Not Collective
2613a7e14dcfSSatish Balay 
2614a7e14dcfSSatish Balay    Input Parameter:
2615441846f8SBarry Smith .  tao  - Tao context
2616a7e14dcfSSatish Balay 
2617a7e14dcfSSatish Balay    Output Parameter:
2618a7e14dcfSSatish Balay .  usrP - user context
2619a7e14dcfSSatish Balay 
2620a7e14dcfSSatish Balay    Level: intermediate
2621a7e14dcfSSatish Balay 
2622a7e14dcfSSatish Balay .seealso: TaoSetApplicationContext()
2623a7e14dcfSSatish Balay @*/
2624441846f8SBarry Smith PetscErrorCode  TaoGetApplicationContext(Tao tao,void *usrP)
2625a7e14dcfSSatish Balay {
2626a7e14dcfSSatish Balay   PetscFunctionBegin;
2627441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2628a7e14dcfSSatish Balay   *(void**)usrP = tao->user;
2629a7e14dcfSSatish Balay   PetscFunctionReturn(0);
2630a7e14dcfSSatish Balay }
2631a9603a14SPatrick Farrell 
2632a9603a14SPatrick Farrell #undef __FUNCT__
2633a9603a14SPatrick Farrell #define __FUNCT__ "TaoSetGradientNorm"
2634a9603a14SPatrick Farrell /*@
2635a9603a14SPatrick Farrell    TaoSetGradientNorm - Sets the matrix used to define the inner product that measures the size of the gradient.
2636a9603a14SPatrick Farrell 
2637a9603a14SPatrick Farrell    Collective on tao
2638a9603a14SPatrick Farrell 
2639a9603a14SPatrick Farrell    Input Parameters:
2640a9603a14SPatrick Farrell +  tao  - the Tao context
2641a9603a14SPatrick Farrell -  M    - gradient norm
2642a9603a14SPatrick Farrell 
2643a9603a14SPatrick Farrell    Level: beginner
2644a9603a14SPatrick Farrell 
2645a9603a14SPatrick Farrell .seealso: TaoGetGradientNorm(), TaoGradientNorm()
2646a9603a14SPatrick Farrell @*/
2647a9603a14SPatrick Farrell PetscErrorCode  TaoSetGradientNorm(Tao tao, Mat M)
2648a9603a14SPatrick Farrell {
2649a9603a14SPatrick Farrell   PetscErrorCode ierr;
2650a9603a14SPatrick Farrell 
2651a9603a14SPatrick Farrell   PetscFunctionBegin;
2652a9603a14SPatrick Farrell   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2653a9603a14SPatrick Farrell 
2654a9603a14SPatrick Farrell   if (tao->gradient_norm) {
2655a9603a14SPatrick Farrell     ierr = PetscObjectDereference((PetscObject)tao->gradient_norm);CHKERRQ(ierr);
2656a9603a14SPatrick Farrell     ierr = VecDestroy(&tao->gradient_norm_tmp);CHKERRQ(ierr);
2657a9603a14SPatrick Farrell   }
2658a9603a14SPatrick Farrell 
2659a9603a14SPatrick Farrell   ierr = PetscObjectReference((PetscObject)M);CHKERRQ(ierr);
2660a9603a14SPatrick Farrell   tao->gradient_norm = M;
2661a9603a14SPatrick Farrell   ierr = MatCreateVecs(M, NULL, &tao->gradient_norm_tmp);CHKERRQ(ierr);
2662a9603a14SPatrick Farrell   PetscFunctionReturn(0);
2663a9603a14SPatrick Farrell }
2664a9603a14SPatrick Farrell 
2665a9603a14SPatrick Farrell #undef __FUNCT__
2666a9603a14SPatrick Farrell #define __FUNCT__ "TaoGetGradientNorm"
2667a9603a14SPatrick Farrell /*@
2668a9603a14SPatrick Farrell    TaoGetGradientNorm - Returns the matrix used to define the inner product for measuring the size of the gradient.
2669a9603a14SPatrick Farrell 
2670a9603a14SPatrick Farrell    Not Collective
2671a9603a14SPatrick Farrell 
2672a9603a14SPatrick Farrell    Input Parameter:
2673a9603a14SPatrick Farrell .  tao  - Tao context
2674a9603a14SPatrick Farrell 
2675a9603a14SPatrick Farrell    Output Parameter:
2676a9603a14SPatrick Farrell .  M - gradient norm
2677a9603a14SPatrick Farrell 
2678a9603a14SPatrick Farrell    Level: beginner
2679a9603a14SPatrick Farrell 
2680a9603a14SPatrick Farrell .seealso: TaoSetGradientNorm(), TaoGradientNorm()
2681a9603a14SPatrick Farrell @*/
2682a9603a14SPatrick Farrell PetscErrorCode  TaoGetGradientNorm(Tao tao, Mat *M)
2683a9603a14SPatrick Farrell {
2684a9603a14SPatrick Farrell   PetscFunctionBegin;
2685a9603a14SPatrick Farrell   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2686a9603a14SPatrick Farrell   *M = tao->gradient_norm;
2687a9603a14SPatrick Farrell   PetscFunctionReturn(0);
2688a9603a14SPatrick Farrell }
2689a9603a14SPatrick Farrell 
2690a9603a14SPatrick Farrell #undef __FUNCT__
2691a9603a14SPatrick Farrell #define __FUNCT__ "TaoGradientNorm"
2692a9603a14SPatrick Farrell /*c
2693a9603a14SPatrick Farrell    TaoGradientNorm - Compute the norm with respect to the inner product the user has set.
2694a9603a14SPatrick Farrell 
2695a9603a14SPatrick Farrell    Collective on tao
2696a9603a14SPatrick Farrell 
2697a9603a14SPatrick Farrell    Input Parameter:
2698a9603a14SPatrick Farrell .  tao      - the Tao context
2699a9603a14SPatrick Farrell .  gradient - the gradient to be computed
2700a9603a14SPatrick Farrell .  norm     - the norm type
2701a9603a14SPatrick Farrell 
2702a9603a14SPatrick Farrell    Output Parameter:
2703a9603a14SPatrick Farrell .  gnorm    - the gradient norm
2704a9603a14SPatrick Farrell 
2705a9603a14SPatrick Farrell    Level: developer
2706a9603a14SPatrick Farrell 
2707a9603a14SPatrick Farrell .seealso: TaoSetGradientNorm(), TaoGetGradientNorm()
2708a9603a14SPatrick Farrell @*/
2709a9603a14SPatrick Farrell PetscErrorCode  TaoGradientNorm(Tao tao, Vec gradient, NormType type, PetscReal *gnorm)
2710a9603a14SPatrick Farrell {
2711a9603a14SPatrick Farrell   PetscErrorCode ierr;
2712a9603a14SPatrick Farrell 
2713a9603a14SPatrick Farrell   PetscFunctionBegin;
2714a9603a14SPatrick Farrell   PetscValidHeaderSpecific(gradient,VEC_CLASSID,1);
2715a9603a14SPatrick Farrell 
2716a9603a14SPatrick Farrell   if (tao->gradient_norm) {
2717680e2bc7SPatrick Farrell     PetscScalar gnorms;
2718680e2bc7SPatrick Farrell 
2719a9603a14SPatrick Farrell     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.");
2720a9603a14SPatrick Farrell     ierr = MatMult(tao->gradient_norm, gradient, tao->gradient_norm_tmp);CHKERRQ(ierr);
2721680e2bc7SPatrick Farrell     ierr = VecDot(gradient, tao->gradient_norm_tmp, &gnorms);CHKERRQ(ierr);
27227bb79937SPatrick Farrell     *gnorm = PetscRealPart(PetscSqrtScalar(gnorms));
2723a9603a14SPatrick Farrell   } else {
2724a9603a14SPatrick Farrell     ierr = VecNorm(gradient, type, gnorm);CHKERRQ(ierr);
2725a9603a14SPatrick Farrell   }
2726a9603a14SPatrick Farrell   PetscFunctionReturn(0);
2727a9603a14SPatrick Farrell }
2728a9603a14SPatrick Farrell 
2729a9603a14SPatrick Farrell 
2730