xref: /petsc/src/tao/interface/taosolver.c (revision e4cb33bb7dbdbae9285fba102465ca0f1dcb3977)
1441846f8SBarry Smith #define TAO_DLL
2a7e14dcfSSatish Balay 
3ba92ff59SBarry 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 
58441846f8SBarry Smith   ierr = PetscHeaderCreate(tao,_p_Tao, struct _TaoOps, TAO_CLASSID,"Tao",0,0,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;
8547a47007SBarry Smith   tao->stepdirection=NULL;
8647a47007SBarry Smith   tao->XL = NULL;
8747a47007SBarry Smith   tao->XU = NULL;
8847a47007SBarry Smith   tao->IL = NULL;
8947a47007SBarry Smith   tao->IU = NULL;
9047a47007SBarry Smith   tao->DI = NULL;
9147a47007SBarry Smith   tao->DE = NULL;
9247a47007SBarry Smith   tao->hessian = NULL;
9347a47007SBarry Smith   tao->hessian_pre = NULL;
9447a47007SBarry Smith   tao->jacobian = NULL;
9547a47007SBarry Smith   tao->jacobian_pre = NULL;
9647a47007SBarry Smith   tao->jacobian_state = NULL;
9747a47007SBarry Smith   tao->jacobian_state_pre = NULL;
9847a47007SBarry Smith   tao->jacobian_state_inv = NULL;
9947a47007SBarry Smith   tao->jacobian_design = NULL;
10047a47007SBarry Smith   tao->jacobian_design_pre = NULL;
10147a47007SBarry Smith   tao->jacobian_equality = NULL;
10247a47007SBarry Smith   tao->jacobian_equality_pre = NULL;
10347a47007SBarry Smith   tao->jacobian_inequality = NULL;
10447a47007SBarry Smith   tao->jacobian_inequality_pre = NULL;
10547a47007SBarry Smith   tao->state_is = NULL;
10647a47007SBarry Smith   tao->design_is = NULL;
107a7e14dcfSSatish Balay 
108a7e14dcfSSatish Balay   tao->max_it     = 10000;
109a7e14dcfSSatish Balay   tao->max_funcs   = 10000;
1106f4723b1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE)
1116f4723b1SBarry Smith   tao->fatol       = 1e-5;
1126f4723b1SBarry Smith   tao->frtol       = 1e-5;
1136f4723b1SBarry Smith   tao->gatol       = 1e-5;
1146f4723b1SBarry Smith   tao->grtol       = 1e-5;
1156f4723b1SBarry Smith #else
116a7e14dcfSSatish Balay   tao->fatol       = 1e-8;
117a7e14dcfSSatish Balay   tao->frtol       = 1e-8;
118a7e14dcfSSatish Balay   tao->gatol       = 1e-8;
119a7e14dcfSSatish Balay   tao->grtol       = 1e-8;
1206f4723b1SBarry Smith #endif
121a7e14dcfSSatish Balay   tao->gttol       = 0.0;
122a7e14dcfSSatish Balay   tao->catol       = 0.0;
123a7e14dcfSSatish Balay   tao->crtol       = 0.0;
124a7e14dcfSSatish Balay   tao->xtol        = 0.0;
125a7e14dcfSSatish Balay   tao->steptol     = 0.0;
126e270355aSBarry Smith   tao->trust0      = PETSC_INFINITY;
1276f4723b1SBarry Smith   tao->fmin        = PETSC_NINFINITY;
128a7e14dcfSSatish Balay   tao->hist_reset = PETSC_TRUE;
129a7e14dcfSSatish Balay   tao->hist_max = 0;
130a7e14dcfSSatish Balay   tao->hist_len = 0;
13147a47007SBarry Smith   tao->hist_obj = NULL;
13247a47007SBarry Smith   tao->hist_resid = NULL;
13347a47007SBarry Smith   tao->hist_cnorm = NULL;
134a7e14dcfSSatish Balay 
135a7e14dcfSSatish Balay   tao->numbermonitors=0;
136a7e14dcfSSatish Balay   tao->viewsolution=PETSC_FALSE;
137a7e14dcfSSatish Balay   tao->viewhessian=PETSC_FALSE;
138a7e14dcfSSatish Balay   tao->viewgradient=PETSC_FALSE;
139a7e14dcfSSatish Balay   tao->viewjacobian=PETSC_FALSE;
140a7e14dcfSSatish Balay   tao->viewconstraints = PETSC_FALSE;
141a7e14dcfSSatish Balay   tao->viewtao = PETSC_FALSE;
142a7e14dcfSSatish Balay 
143a7e14dcfSSatish Balay   ierr = TaoResetStatistics(tao);CHKERRQ(ierr);
144a7e14dcfSSatish Balay   *newtao = tao;
145a7e14dcfSSatish Balay   PetscFunctionReturn(0);
146a7e14dcfSSatish Balay }
147a7e14dcfSSatish Balay 
148a7e14dcfSSatish Balay #undef __FUNCT__
149a7e14dcfSSatish Balay #define __FUNCT__ "TaoSolve"
150a7e14dcfSSatish Balay /*@
151a7e14dcfSSatish Balay   TaoSolve - Solves an optimization problem min F(x) s.t. l <= x <= u
152a7e14dcfSSatish Balay 
153441846f8SBarry Smith   Collective on Tao
154a7e14dcfSSatish Balay 
155a7e14dcfSSatish Balay   Input Parameters:
156441846f8SBarry Smith . tao - the Tao context
157a7e14dcfSSatish Balay 
158a7e14dcfSSatish Balay   Notes:
159441846f8SBarry Smith   The user must set up the Tao with calls to TaoSetInitialVector(),
160a7e14dcfSSatish Balay   TaoSetObjectiveRoutine(),
161a7e14dcfSSatish Balay   TaoSetGradientRoutine(), and (if using 2nd order method) TaoSetHessianRoutine().
162a7e14dcfSSatish Balay 
163a7e14dcfSSatish Balay   Level: beginner
164a7e14dcfSSatish Balay 
165a7e14dcfSSatish Balay .seealso: TaoCreate(), TaoSetObjectiveRoutine(), TaoSetGradientRoutine(), TaoSetHessianRoutine()
166a7e14dcfSSatish Balay  @*/
167441846f8SBarry Smith PetscErrorCode TaoSolve(Tao tao)
168a7e14dcfSSatish Balay {
169a7e14dcfSSatish Balay   PetscErrorCode   ierr;
170a7e14dcfSSatish Balay   char             filename[PETSC_MAX_PATH_LEN];
171a7e14dcfSSatish Balay   PetscBool        flg;
172a7e14dcfSSatish Balay   PetscViewer      viewer;
173e2379f4fSBarry Smith   static PetscBool set = PETSC_FALSE;
17447a47007SBarry Smith 
175a7e14dcfSSatish Balay   PetscFunctionBegin;
176441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
177e2379f4fSBarry Smith   ierr = PetscCitationsRegister("@TechReport{tao-user-ref,\n"
178e2379f4fSBarry Smith                                 "title   = {Toolkit for Advanced Optimization (TAO) Users Manual},\n"
179e2379f4fSBarry Smith                                 "author  = {Todd Munson and Jason Sarich and Stefan Wild and Steve Benson and Lois Curfman McInnes},\n"
180e2379f4fSBarry Smith                                 "Institution = {Argonne National Laboratory},\n"
181e2379f4fSBarry Smith                                 "Year   = 2014,\n"
182e2379f4fSBarry Smith                                 "Number = {ANL/MCS-TM-322 - Revision 3.5},\n"
183e2379f4fSBarry Smith                                 "url    = {http://www.mcs.anl.gov/tao}\n}\n",&set);
184e2379f4fSBarry Smith 
185a7e14dcfSSatish Balay   ierr = TaoSetUp(tao);CHKERRQ(ierr);
186a7e14dcfSSatish Balay   ierr = TaoResetStatistics(tao);CHKERRQ(ierr);
187a7e14dcfSSatish Balay   if (tao->linesearch) {
188a7e14dcfSSatish Balay     ierr = TaoLineSearchReset(tao->linesearch);CHKERRQ(ierr);
189a7e14dcfSSatish Balay   }
190a7e14dcfSSatish Balay 
191441846f8SBarry Smith   ierr = PetscLogEventBegin(Tao_Solve,tao,0,0,0);CHKERRQ(ierr);
192a7e14dcfSSatish Balay   if (tao->ops->solve){ ierr = (*tao->ops->solve)(tao);CHKERRQ(ierr); }
193441846f8SBarry Smith   ierr = PetscLogEventEnd(Tao_Solve,tao,0,0,0);CHKERRQ(ierr);
194a7e14dcfSSatish Balay 
195a7e14dcfSSatish Balay   ierr = PetscOptionsGetString(((PetscObject)tao)->prefix,"-tao_view",filename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
196a7e14dcfSSatish Balay   if (flg && !PetscPreLoadingOn) {
197a7e14dcfSSatish Balay     ierr = PetscViewerASCIIOpen(((PetscObject)tao)->comm,filename,&viewer);CHKERRQ(ierr);
198a7e14dcfSSatish Balay     ierr = TaoView(tao,viewer);CHKERRQ(ierr);
199a7e14dcfSSatish Balay     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
200a7e14dcfSSatish Balay   }
201a7e14dcfSSatish Balay 
202a7e14dcfSSatish Balay   if (tao->printreason) {
203a7e14dcfSSatish Balay     if (tao->reason > 0) {
204*e4cb33bbSBarry Smith       ierr = PetscPrintf(((PetscObject)tao)->comm,"TAO solve converged due to %s\n",TaoConvergedReasons[tao->reason]);CHKERRQ(ierr);
205a7e14dcfSSatish Balay     } else {
206*e4cb33bbSBarry Smith       ierr = PetscPrintf(((PetscObject)tao)->comm,"TAO solve did not converge due to %s\n",TaoConvergedReasons[tao->reason]);CHKERRQ(ierr);
207a7e14dcfSSatish Balay     }
208a7e14dcfSSatish Balay   }
209a7e14dcfSSatish Balay   PetscFunctionReturn(0);
210a7e14dcfSSatish Balay }
211a7e14dcfSSatish Balay 
212a7e14dcfSSatish Balay #undef __FUNCT__
213a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetUp"
214a7e14dcfSSatish Balay /*@
215a7e14dcfSSatish Balay   TaoSetUp - Sets up the internal data structures for the later use
216a7e14dcfSSatish Balay   of a Tao solver
217a7e14dcfSSatish Balay 
218a7e14dcfSSatish Balay   Collective on tao
219a7e14dcfSSatish Balay 
220a7e14dcfSSatish Balay   Input Parameters:
221a7e14dcfSSatish Balay . tao - the TAO context
222a7e14dcfSSatish Balay 
223a7e14dcfSSatish Balay   Notes:
224a7e14dcfSSatish Balay   The user will not need to explicitly call TaoSetUp(), as it will
225a7e14dcfSSatish Balay   automatically be called in TaoSolve().  However, if the user
226a7e14dcfSSatish Balay   desires to call it explicitly, it should come after TaoCreate()
227a7e14dcfSSatish Balay   and any TaoSetSomething() routines, but before TaoSolve().
228a7e14dcfSSatish Balay 
229a7e14dcfSSatish Balay   Level: advanced
230a7e14dcfSSatish Balay 
231a7e14dcfSSatish Balay .seealso: TaoCreate(), TaoSolve()
232a7e14dcfSSatish Balay @*/
233441846f8SBarry Smith PetscErrorCode TaoSetUp(Tao tao)
234a7e14dcfSSatish Balay {
235a7e14dcfSSatish Balay   PetscErrorCode ierr;
23647a47007SBarry Smith 
237a7e14dcfSSatish Balay   PetscFunctionBegin;
238441846f8SBarry Smith   PetscValidHeaderSpecific(tao, TAO_CLASSID,1);
239a7e14dcfSSatish Balay   if (tao->setupcalled) PetscFunctionReturn(0);
240a7e14dcfSSatish Balay 
24147a47007SBarry Smith   if (!tao->solution) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TaoSetInitialVector");
242a7e14dcfSSatish Balay   if (tao->ops->setup) {
243a7e14dcfSSatish Balay     ierr = (*tao->ops->setup)(tao);CHKERRQ(ierr);
244a7e14dcfSSatish Balay   }
245a7e14dcfSSatish Balay   tao->setupcalled = PETSC_TRUE;
246a7e14dcfSSatish Balay   PetscFunctionReturn(0);
247a7e14dcfSSatish Balay }
248a7e14dcfSSatish Balay 
249a7e14dcfSSatish Balay #undef __FUNCT__
250a7e14dcfSSatish Balay #define __FUNCT__ "TaoDestroy"
251a7e14dcfSSatish Balay /*@
252a7e14dcfSSatish Balay   TaoDestroy - Destroys the TAO context that was created with
253a7e14dcfSSatish Balay   TaoCreate()
254a7e14dcfSSatish Balay 
255441846f8SBarry Smith   Collective on Tao
256a7e14dcfSSatish Balay 
257a7e14dcfSSatish Balay   Input Parameter:
258441846f8SBarry Smith . tao - the Tao context
259a7e14dcfSSatish Balay 
260a7e14dcfSSatish Balay   Level: beginner
261a7e14dcfSSatish Balay 
262a7e14dcfSSatish Balay .seealso: TaoCreate(), TaoSolve()
263a7e14dcfSSatish Balay @*/
264441846f8SBarry Smith PetscErrorCode TaoDestroy(Tao *tao)
265a7e14dcfSSatish Balay {
266a7e14dcfSSatish Balay   PetscErrorCode ierr;
26747a47007SBarry Smith 
268a7e14dcfSSatish Balay   PetscFunctionBegin;
269a7e14dcfSSatish Balay   if (!*tao) PetscFunctionReturn(0);
270441846f8SBarry Smith   PetscValidHeaderSpecific(*tao,TAO_CLASSID,1);
271a7e14dcfSSatish Balay   if (--((PetscObject)*tao)->refct > 0) {*tao=0;PetscFunctionReturn(0);}
272a7e14dcfSSatish Balay 
273a7e14dcfSSatish Balay   if ((*tao)->ops->destroy) {
274a7e14dcfSSatish Balay     ierr = (*((*tao))->ops->destroy)(*tao);CHKERRQ(ierr);
275a7e14dcfSSatish Balay   }
276a7e14dcfSSatish Balay   ierr = KSPDestroy(&(*tao)->ksp);CHKERRQ(ierr);
277a7e14dcfSSatish Balay   ierr = TaoLineSearchDestroy(&(*tao)->linesearch);CHKERRQ(ierr);
278a7e14dcfSSatish Balay 
279a7e14dcfSSatish Balay   if ((*tao)->ops->convergencedestroy) {
280a7e14dcfSSatish Balay     ierr = (*(*tao)->ops->convergencedestroy)((*tao)->cnvP);CHKERRQ(ierr);
281a7e14dcfSSatish Balay     if ((*tao)->jacobian_state_inv) {
282a7e14dcfSSatish Balay       ierr = MatDestroy(&(*tao)->jacobian_state_inv);CHKERRQ(ierr);
283a7e14dcfSSatish Balay     }
284a7e14dcfSSatish Balay   }
285a7e14dcfSSatish Balay   ierr = VecDestroy(&(*tao)->solution);CHKERRQ(ierr);
286a7e14dcfSSatish Balay   ierr = VecDestroy(&(*tao)->gradient);CHKERRQ(ierr);
287a7e14dcfSSatish Balay 
288a7e14dcfSSatish Balay   ierr = VecDestroy(&(*tao)->XL);CHKERRQ(ierr);
289a7e14dcfSSatish Balay   ierr = VecDestroy(&(*tao)->XU);CHKERRQ(ierr);
290a7e14dcfSSatish Balay   ierr = VecDestroy(&(*tao)->IL);CHKERRQ(ierr);
291a7e14dcfSSatish Balay   ierr = VecDestroy(&(*tao)->IU);CHKERRQ(ierr);
292a7e14dcfSSatish Balay   ierr = VecDestroy(&(*tao)->DE);CHKERRQ(ierr);
293a7e14dcfSSatish Balay   ierr = VecDestroy(&(*tao)->DI);CHKERRQ(ierr);
294a7e14dcfSSatish Balay   ierr = VecDestroy(&(*tao)->constraints_equality);CHKERRQ(ierr);
295a7e14dcfSSatish Balay   ierr = VecDestroy(&(*tao)->constraints_inequality);CHKERRQ(ierr);
296a7e14dcfSSatish Balay   ierr = VecDestroy(&(*tao)->stepdirection);CHKERRQ(ierr);
297a7e14dcfSSatish Balay   ierr = MatDestroy(&(*tao)->hessian_pre);CHKERRQ(ierr);
298a7e14dcfSSatish Balay   ierr = MatDestroy(&(*tao)->hessian);CHKERRQ(ierr);
299a7e14dcfSSatish Balay   ierr = MatDestroy(&(*tao)->jacobian_pre);CHKERRQ(ierr);
300a7e14dcfSSatish Balay   ierr = MatDestroy(&(*tao)->jacobian);CHKERRQ(ierr);
301a7e14dcfSSatish Balay   ierr = MatDestroy(&(*tao)->jacobian_state_pre);CHKERRQ(ierr);
302a7e14dcfSSatish Balay   ierr = MatDestroy(&(*tao)->jacobian_state);CHKERRQ(ierr);
303a7e14dcfSSatish Balay   ierr = MatDestroy(&(*tao)->jacobian_state_inv);CHKERRQ(ierr);
304067f330cSJason Sarich   ierr = MatDestroy(&(*tao)->jacobian_design);CHKERRQ(ierr);
305a7e14dcfSSatish Balay   ierr = MatDestroy(&(*tao)->jacobian_equality);CHKERRQ(ierr);
306a7e14dcfSSatish Balay   ierr = MatDestroy(&(*tao)->jacobian_equality_pre);CHKERRQ(ierr);
307a7e14dcfSSatish Balay   ierr = MatDestroy(&(*tao)->jacobian_inequality);CHKERRQ(ierr);
308a7e14dcfSSatish Balay   ierr = MatDestroy(&(*tao)->jacobian_inequality_pre);CHKERRQ(ierr);
309a7e14dcfSSatish Balay   ierr = ISDestroy(&(*tao)->state_is);CHKERRQ(ierr);
310a7e14dcfSSatish Balay   ierr = ISDestroy(&(*tao)->design_is);CHKERRQ(ierr);
311a7e14dcfSSatish Balay   ierr = TaoCancelMonitors(*tao);CHKERRQ(ierr);
312a7e14dcfSSatish Balay   ierr = PetscHeaderDestroy(tao);CHKERRQ(ierr);
313a7e14dcfSSatish Balay   PetscFunctionReturn(0);
314a7e14dcfSSatish Balay }
315a7e14dcfSSatish Balay 
316a7e14dcfSSatish Balay #undef __FUNCT__
317a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetFromOptions"
318a7e14dcfSSatish Balay /*@
319441846f8SBarry Smith   TaoSetFromOptions - Sets various Tao parameters from user
320a7e14dcfSSatish Balay   options.
321a7e14dcfSSatish Balay 
322441846f8SBarry Smith   Collective on Tao
323a7e14dcfSSatish Balay 
324a7e14dcfSSatish Balay   Input Paremeter:
325441846f8SBarry Smith . tao - the Tao solver context
326a7e14dcfSSatish Balay 
327a7e14dcfSSatish Balay   options Database Keys:
32858417fe7SBarry Smith + -tao_type <type> - The algorithm that TAO uses (lmvm, nls, etc.)
329a7e14dcfSSatish Balay . -tao_fatol <fatol> - absolute error tolerance in function value
330a7e14dcfSSatish Balay . -tao_frtol <frtol> - relative error tolerance in function value
331a7e14dcfSSatish Balay . -tao_gatol <gatol> - absolute error tolerance for ||gradient||
332a7e14dcfSSatish Balay . -tao_grtol <grtol> - relative error tolerance for ||gradient||
333a7e14dcfSSatish Balay . -tao_gttol <gttol> - reduction of ||gradient|| relative to initial gradient
334a7e14dcfSSatish Balay . -tao_max_it <max> - sets maximum number of iterations
335a7e14dcfSSatish Balay . -tao_max_funcs <max> - sets maximum number of function evaluations
336a7e14dcfSSatish Balay . -tao_fmin <fmin> - stop if function value reaches fmin
337a7e14dcfSSatish Balay . -tao_steptol <tol> - stop if trust region radius less than <tol>
338a7e14dcfSSatish Balay . -tao_trust0 <t> - initial trust region radius
339a7e14dcfSSatish Balay . -tao_monitor - prints function value and residual at each iteration
340a7e14dcfSSatish Balay . -tao_smonitor - same as tao_monitor, but truncates very small values
341a7e14dcfSSatish Balay . -tao_cmonitor - prints function value, residual, and constraint norm at each iteration
342a7e14dcfSSatish Balay . -tao_view_solution - prints solution vector at each iteration
343a7e14dcfSSatish Balay . -tao_view_separableobjective - prints separable objective vector at each iteration
344a7e14dcfSSatish Balay . -tao_view_step - prints step direction vector at each iteration
345a7e14dcfSSatish Balay . -tao_view_gradient - prints gradient vector at each iteration
346a7e14dcfSSatish Balay . -tao_draw_solution - graphically view solution vector at each iteration
347a7e14dcfSSatish Balay . -tao_draw_step - graphically view step vector at each iteration
348a7e14dcfSSatish Balay . -tao_draw_gradient - graphically view gradient at each iteration
349a7e14dcfSSatish Balay . -tao_fd_gradient - use gradient computed with finite differences
350a7e14dcfSSatish Balay . -tao_cancelmonitors - cancels all monitors (except those set with command line)
351441846f8SBarry Smith . -tao_view - prints information about the Tao after solving
352a7e14dcfSSatish Balay - -tao_converged_reason - prints the reason TAO stopped iterating
353a7e14dcfSSatish Balay 
354a7e14dcfSSatish Balay   Notes:
355a7e14dcfSSatish Balay   To see all options, run your program with the -help option or consult the
356a7e14dcfSSatish Balay   user's manual. Should be called after TaoCreate() but before TaoSolve()
357a7e14dcfSSatish Balay 
358a7e14dcfSSatish Balay   Level: beginner
359a7e14dcfSSatish Balay @*/
360441846f8SBarry Smith PetscErrorCode TaoSetFromOptions(Tao tao)
361a7e14dcfSSatish Balay {
362a7e14dcfSSatish Balay   PetscErrorCode ierr;
36358417fe7SBarry Smith   const TaoType  default_type = TAOLMVM;
364a7e14dcfSSatish Balay   const char     *prefix;
365a7e14dcfSSatish Balay   char           type[256], monfilename[PETSC_MAX_PATH_LEN];
366a7e14dcfSSatish Balay   PetscViewer    monviewer;
367a7e14dcfSSatish Balay   PetscBool      flg;
368a7e14dcfSSatish Balay   MPI_Comm       comm;
369a7e14dcfSSatish Balay 
370a7e14dcfSSatish Balay   PetscFunctionBegin;
371441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
372a7e14dcfSSatish Balay   ierr = PetscObjectGetComm((PetscObject)tao,&comm);CHKERRQ(ierr);
373a7e14dcfSSatish Balay   ierr = TaoGetOptionsPrefix(tao,&prefix);
374a7e14dcfSSatish Balay   /* So no warnings are given about unused options */
375a7e14dcfSSatish Balay   ierr = PetscOptionsHasName(prefix,"-tao_ls_type",&flg);
376a7e14dcfSSatish Balay 
377a7e14dcfSSatish Balay   ierr = PetscObjectOptionsBegin((PetscObject)tao);CHKERRQ(ierr);
378a7e14dcfSSatish Balay   {
379441846f8SBarry Smith     if (!TaoRegisterAllCalled) {
380441846f8SBarry Smith       ierr = TaoRegisterAll();CHKERRQ(ierr);
381a7e14dcfSSatish Balay     }
382a7e14dcfSSatish Balay     if (((PetscObject)tao)->type_name) {
383a7e14dcfSSatish Balay       default_type = ((PetscObject)tao)->type_name;
384a7e14dcfSSatish Balay     }
385a7e14dcfSSatish Balay     /* Check for type from options */
386441846f8SBarry Smith     ierr = PetscOptionsFList("-tao_type","Tao Solver type","TaoSetType",TaoList,default_type,type,256,&flg);CHKERRQ(ierr);
387a7e14dcfSSatish Balay     if (flg) {
388a7e14dcfSSatish Balay       ierr = TaoSetType(tao,type);CHKERRQ(ierr);
389a7e14dcfSSatish Balay     } else if (!((PetscObject)tao)->type_name) {
390a7e14dcfSSatish Balay       ierr = TaoSetType(tao,default_type);
391a7e14dcfSSatish Balay     }
392a7e14dcfSSatish Balay 
393441846f8SBarry Smith     ierr = PetscOptionsBool("-tao_view","view Tao info after each minimization has completed","TaoView",PETSC_FALSE,&tao->viewtao,&flg);CHKERRQ(ierr);
394a7e14dcfSSatish Balay     ierr = PetscOptionsReal("-tao_fatol","Stop if solution within","TaoSetTolerances",tao->fatol,&tao->fatol,&flg);CHKERRQ(ierr);
395a7e14dcfSSatish Balay     ierr = PetscOptionsReal("-tao_frtol","Stop if relative solution within","TaoSetTolerances",tao->frtol,&tao->frtol,&flg);CHKERRQ(ierr);
396a7e14dcfSSatish Balay     ierr = PetscOptionsReal("-tao_catol","Stop if constraints violations within","TaoSetConstraintTolerances",tao->catol,&tao->catol,&flg);CHKERRQ(ierr);
397a7e14dcfSSatish Balay     ierr = PetscOptionsReal("-tao_crtol","Stop if relative contraint violations within","TaoSetConstraintTolerances",tao->crtol,&tao->crtol,&flg);CHKERRQ(ierr);
398a7e14dcfSSatish Balay     ierr = PetscOptionsReal("-tao_gatol","Stop if norm of gradient less than","TaoSetTolerances",tao->gatol,&tao->gatol,&flg);CHKERRQ(ierr);
399a7e14dcfSSatish Balay     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);
400a7e14dcfSSatish Balay     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);
4016f4723b1SBarry Smith     ierr = PetscOptionsInt("-tao_max_it","Stop if iteration number exceeds","TaoSetMaximumIterations",tao->max_it,&tao->max_it,&flg);CHKERRQ(ierr);
402a7e14dcfSSatish Balay     ierr = PetscOptionsInt("-tao_max_funcs","Stop if number of function evaluations exceeds","TaoSetMaximumFunctionEvaluations",tao->max_funcs,&tao->max_funcs,&flg);CHKERRQ(ierr);
403a7e14dcfSSatish Balay     ierr = PetscOptionsReal("-tao_fmin","Stop if function less than","TaoSetFunctionLowerBound",tao->fmin,&tao->fmin,&flg);CHKERRQ(ierr);
404a7e14dcfSSatish Balay     ierr = PetscOptionsReal("-tao_steptol","Stop if step size or trust region radius less than","",tao->steptol,&tao->steptol,&flg);CHKERRQ(ierr);
405a7e14dcfSSatish Balay     ierr = PetscOptionsReal("-tao_trust0","Initial trust region radius","TaoSetTrustRegionRadius",tao->trust0,&tao->trust0,&flg);CHKERRQ(ierr);
406a7e14dcfSSatish Balay 
407a7e14dcfSSatish Balay     ierr = PetscOptionsString("-tao_view_solution","view solution vector after each evaluation","TaoSetMonitor","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
408a7e14dcfSSatish Balay     if (flg) {
409a7e14dcfSSatish Balay       ierr = PetscViewerASCIIOpen(comm,monfilename,&monviewer);CHKERRQ(ierr);
410a7e14dcfSSatish Balay       ierr = TaoSetMonitor(tao,TaoSolutionMonitor,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
411a7e14dcfSSatish Balay     }
412a7e14dcfSSatish Balay 
413*e4cb33bbSBarry Smith     ierr = PetscOptionsBool("-tao_converged_reason","Print reason for TAO converged","TaoSolve",flg,&flg,NULL);CHKERRQ(ierr);
414a7e14dcfSSatish Balay     if (flg) {
415a7e14dcfSSatish Balay       tao->printreason = PETSC_TRUE;
416a7e14dcfSSatish Balay     }
417a7e14dcfSSatish Balay     ierr = PetscOptionsString("-tao_view_gradient","view gradient vector after each evaluation","TaoSetMonitor","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
418a7e14dcfSSatish Balay     if (flg) {
419a7e14dcfSSatish Balay       ierr = PetscViewerASCIIOpen(comm,monfilename,&monviewer);CHKERRQ(ierr);
420a7e14dcfSSatish Balay       ierr = TaoSetMonitor(tao,TaoGradientMonitor,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
421a7e14dcfSSatish Balay     }
422a7e14dcfSSatish Balay 
423a7e14dcfSSatish Balay     ierr = PetscOptionsString("-tao_view_stepdirection","view step direction vector after each iteration","TaoSetMonitor","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
424a7e14dcfSSatish Balay     if (flg) {
425a7e14dcfSSatish Balay       ierr = PetscViewerASCIIOpen(comm,monfilename,&monviewer);CHKERRQ(ierr);
426a7e14dcfSSatish Balay       ierr = TaoSetMonitor(tao,TaoStepDirectionMonitor,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
427a7e14dcfSSatish Balay     }
428a7e14dcfSSatish Balay 
429a7e14dcfSSatish Balay     ierr = PetscOptionsString("-tao_view_separableobjective","view separable objective vector after each evaluation","TaoSetMonitor","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
430a7e14dcfSSatish Balay     if (flg) {
431a7e14dcfSSatish Balay       ierr = PetscViewerASCIIOpen(comm,monfilename,&monviewer);CHKERRQ(ierr);
432a7e14dcfSSatish Balay       ierr = TaoSetMonitor(tao,TaoSeparableObjectiveMonitor,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
433a7e14dcfSSatish Balay     }
434a7e14dcfSSatish Balay 
435a7e14dcfSSatish Balay     ierr = PetscOptionsString("-tao_monitor","Use the default convergence monitor","TaoSetMonitor","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
436a7e14dcfSSatish Balay     if (flg) {
437a7e14dcfSSatish Balay       ierr = PetscViewerASCIIOpen(comm,monfilename,&monviewer);CHKERRQ(ierr);
438a7e14dcfSSatish Balay       ierr = TaoSetMonitor(tao,TaoDefaultMonitor,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
439a7e14dcfSSatish Balay     }
440a7e14dcfSSatish Balay 
441a7e14dcfSSatish Balay     ierr = PetscOptionsString("-tao_smonitor","Use the short convergence monitor","TaoSetMonitor","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
442a7e14dcfSSatish Balay     if (flg) {
443a7e14dcfSSatish Balay       ierr = PetscViewerASCIIOpen(comm,monfilename,&monviewer);CHKERRQ(ierr);
444a7e14dcfSSatish Balay       ierr = TaoSetMonitor(tao,TaoDefaultSMonitor,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
445a7e14dcfSSatish Balay     }
446a7e14dcfSSatish Balay 
447a7e14dcfSSatish Balay     ierr = PetscOptionsString("-tao_cmonitor","Use the default convergence monitor with constraint norm","TaoSetMonitor","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
448a7e14dcfSSatish Balay     if (flg) {
449a7e14dcfSSatish Balay       ierr = PetscViewerASCIIOpen(comm,monfilename,&monviewer);CHKERRQ(ierr);
450a7e14dcfSSatish Balay       ierr = TaoSetMonitor(tao,TaoDefaultCMonitor,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
451a7e14dcfSSatish Balay     }
452a7e14dcfSSatish Balay 
453a7e14dcfSSatish Balay 
45447a47007SBarry Smith     ierr = PetscOptionsBool("-tao_cancelmonitors","cancel all monitors and call any registered destroy routines","TaoCancelMonitors",PETSC_FALSE,&flg,NULL);CHKERRQ(ierr);
455a7e14dcfSSatish Balay     if (flg) {ierr = TaoCancelMonitors(tao);CHKERRQ(ierr);}
456a7e14dcfSSatish Balay 
45747a47007SBarry Smith     ierr = PetscOptionsBool("-tao_draw_solution","Plot solution vector at each iteration","TaoSetMonitor",PETSC_FALSE,&flg,NULL);CHKERRQ(ierr);
458a7e14dcfSSatish Balay     if (flg) {
45947a47007SBarry Smith       ierr = TaoSetMonitor(tao,TaoDrawSolutionMonitor,NULL,NULL);CHKERRQ(ierr);
460a7e14dcfSSatish Balay     }
461a7e14dcfSSatish Balay 
46247a47007SBarry Smith     ierr = PetscOptionsBool("-tao_draw_step","plots step direction at each iteration","TaoSetMonitor",PETSC_FALSE,&flg,NULL);CHKERRQ(ierr);
463a7e14dcfSSatish Balay     if (flg) {
46447a47007SBarry Smith       ierr = TaoSetMonitor(tao,TaoDrawStepMonitor,NULL,NULL);CHKERRQ(ierr);
465a7e14dcfSSatish Balay     }
466a7e14dcfSSatish Balay 
46747a47007SBarry Smith     ierr = PetscOptionsBool("-tao_draw_gradient","plots gradient at each iteration","TaoSetMonitor",PETSC_FALSE,&flg,NULL);CHKERRQ(ierr);
468a7e14dcfSSatish Balay     if (flg) {
46947a47007SBarry Smith       ierr = TaoSetMonitor(tao,TaoDrawGradientMonitor,NULL,NULL);CHKERRQ(ierr);
470a7e14dcfSSatish Balay     }
47147a47007SBarry Smith     ierr = PetscOptionsBool("-tao_fd_gradient","compute gradient using finite differences","TaoDefaultComputeGradient",PETSC_FALSE,&flg,NULL);CHKERRQ(ierr);
472a7e14dcfSSatish Balay     if (flg) {
47347a47007SBarry Smith       ierr = TaoSetGradientRoutine(tao,TaoDefaultComputeGradient,NULL);CHKERRQ(ierr);
474a7e14dcfSSatish Balay     }
475e1cc520bSBarry Smith     ierr = PetscOptionsEnum("-tao_subset_type","subset type", "", TaoSubSetTypes,(PetscEnum)tao->subset_type, (PetscEnum*)&tao->subset_type, 0);CHKERRQ(ierr);
476a7e14dcfSSatish Balay 
477a7e14dcfSSatish Balay     if (tao->ops->setfromoptions) {
478a7e14dcfSSatish Balay       ierr = (*tao->ops->setfromoptions)(tao);CHKERRQ(ierr);
479a7e14dcfSSatish Balay     }
480a7e14dcfSSatish Balay   }
481a7e14dcfSSatish Balay   ierr = PetscOptionsEnd();CHKERRQ(ierr);
482a7e14dcfSSatish Balay   PetscFunctionReturn(0);
483a7e14dcfSSatish Balay }
484a7e14dcfSSatish Balay 
485a7e14dcfSSatish Balay #undef __FUNCT__
486a7e14dcfSSatish Balay #define __FUNCT__ "TaoView"
487a7e14dcfSSatish Balay /*@C
488441846f8SBarry Smith   TaoView - Prints information about the Tao
489a7e14dcfSSatish Balay 
490441846f8SBarry Smith   Collective on Tao
491a7e14dcfSSatish Balay 
492a7e14dcfSSatish Balay   InputParameters:
493441846f8SBarry Smith + tao - the Tao context
494a7e14dcfSSatish Balay - viewer - visualization context
495a7e14dcfSSatish Balay 
496a7e14dcfSSatish Balay   Options Database Key:
497a7e14dcfSSatish Balay . -tao_view - Calls TaoView() at the end of TaoSolve()
498a7e14dcfSSatish Balay 
499a7e14dcfSSatish Balay   Notes:
500a7e14dcfSSatish Balay   The available visualization contexts include
501a7e14dcfSSatish Balay +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
502a7e14dcfSSatish Balay -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
503a7e14dcfSSatish Balay          output where only the first processor opens
504a7e14dcfSSatish Balay          the file.  All other processors send their
505a7e14dcfSSatish Balay          data to the first processor to print.
506a7e14dcfSSatish Balay 
507a7e14dcfSSatish Balay   Level: beginner
508a7e14dcfSSatish Balay 
509a7e14dcfSSatish Balay .seealso: PetscViewerASCIIOpen()
510a7e14dcfSSatish Balay @*/
511441846f8SBarry Smith PetscErrorCode TaoView(Tao tao, PetscViewer viewer)
512a7e14dcfSSatish Balay {
513a7e14dcfSSatish Balay   PetscErrorCode      ierr;
514a7e14dcfSSatish Balay   PetscBool           isascii,isstring;
515441846f8SBarry Smith   const TaoType type;
51647a47007SBarry Smith 
517a7e14dcfSSatish Balay   PetscFunctionBegin;
518441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
519a7e14dcfSSatish Balay   if (!viewer) {
520a7e14dcfSSatish Balay     ierr = PetscViewerASCIIGetStdout(((PetscObject)tao)->comm,&viewer);CHKERRQ(ierr);
521a7e14dcfSSatish Balay   }
522a7e14dcfSSatish Balay   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
523a7e14dcfSSatish Balay   PetscCheckSameComm(tao,1,viewer,2);
524a7e14dcfSSatish Balay 
525a7e14dcfSSatish Balay   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
526a7e14dcfSSatish Balay   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
527a7e14dcfSSatish Balay   if (isascii) {
528a7e14dcfSSatish Balay     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)tao,viewer);CHKERRQ(ierr);
529a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
530a7e14dcfSSatish Balay 
531a7e14dcfSSatish Balay     if (tao->ops->view) {
532a7e14dcfSSatish Balay       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
533a7e14dcfSSatish Balay       ierr = (*tao->ops->view)(tao,viewer);CHKERRQ(ierr);
534a7e14dcfSSatish Balay       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
535a7e14dcfSSatish Balay     }
536a7e14dcfSSatish Balay     if (tao->linesearch) {
537a7e14dcfSSatish Balay       ierr = PetscObjectPrintClassNamePrefixType((PetscObject)(tao->linesearch),viewer);CHKERRQ(ierr);
538a7e14dcfSSatish Balay     }
539a7e14dcfSSatish Balay     if (tao->ksp) {
540a7e14dcfSSatish Balay       ierr = PetscObjectPrintClassNamePrefixType((PetscObject)(tao->ksp),viewer);CHKERRQ(ierr);
541a7e14dcfSSatish Balay       ierr = PetscViewerASCIIPrintf(viewer,"total KSP iterations: %D\n",tao->ksp_its);CHKERRQ(ierr);
542a7e14dcfSSatish Balay     }
543a7e14dcfSSatish Balay     if (tao->XL || tao->XU) {
544e1cc520bSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"Active Set subset type: %s\n",TaoSubSetTypes[tao->subset_type]);
545a7e14dcfSSatish Balay     }
546a7e14dcfSSatish Balay 
54747a47007SBarry Smith     ierr=PetscViewerASCIIPrintf(viewer,"convergence tolerances: fatol=%g,",(double)tao->fatol);CHKERRQ(ierr);
54847a47007SBarry Smith     ierr=PetscViewerASCIIPrintf(viewer," frtol=%g\n",(double)tao->frtol);CHKERRQ(ierr);
549a7e14dcfSSatish Balay 
55047a47007SBarry Smith     ierr=PetscViewerASCIIPrintf(viewer,"convergence tolerances: gatol=%g,",(double)tao->gatol);CHKERRQ(ierr);
55147a47007SBarry Smith     ierr=PetscViewerASCIIPrintf(viewer," steptol=%g,",(double)tao->steptol);CHKERRQ(ierr);
55247a47007SBarry Smith     ierr=PetscViewerASCIIPrintf(viewer," gttol=%g\n",(double)tao->gttol);CHKERRQ(ierr);
553a7e14dcfSSatish Balay 
55447a47007SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"Residual in Function/Gradient:=%g\n",(double)tao->residual);CHKERRQ(ierr);
555a7e14dcfSSatish Balay 
556a7e14dcfSSatish Balay     if (tao->cnorm>0 || tao->catol>0 || tao->crtol>0){
557a7e14dcfSSatish Balay       ierr=PetscViewerASCIIPrintf(viewer,"convergence tolerances:");CHKERRQ(ierr);
55847a47007SBarry Smith       ierr=PetscViewerASCIIPrintf(viewer," catol=%g,",(double)tao->catol);CHKERRQ(ierr);
55947a47007SBarry Smith       ierr=PetscViewerASCIIPrintf(viewer," crtol=%g\n",(double)tao->crtol);CHKERRQ(ierr);
56047a47007SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"Residual in Constraints:=%g\n",(double)tao->cnorm);CHKERRQ(ierr);
561a7e14dcfSSatish Balay     }
562a7e14dcfSSatish Balay 
563a7e14dcfSSatish Balay     if (tao->trust < tao->steptol){
56447a47007SBarry Smith       ierr=PetscViewerASCIIPrintf(viewer,"convergence tolerances: steptol=%g\n",(double)tao->steptol);CHKERRQ(ierr);
56547a47007SBarry Smith       ierr=PetscViewerASCIIPrintf(viewer,"Final trust region radius:=%g\n",(double)tao->trust);CHKERRQ(ierr);
566a7e14dcfSSatish Balay     }
567a7e14dcfSSatish Balay 
568a7e14dcfSSatish Balay     if (tao->fmin>-1.e25){
56947a47007SBarry Smith       ierr=PetscViewerASCIIPrintf(viewer,"convergence tolerances: function minimum=%g\n",(double)tao->fmin);CHKERRQ(ierr);
570a7e14dcfSSatish Balay     }
57147a47007SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"Objective value=%g\n",(double)tao->fc);CHKERRQ(ierr);
572a7e14dcfSSatish Balay 
57347a47007SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"total number of iterations=%D,          ",tao->niter);CHKERRQ(ierr);
574a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"              (max: %D)\n",tao->max_it);CHKERRQ(ierr);
575a7e14dcfSSatish Balay 
576a7e14dcfSSatish Balay     if (tao->nfuncs>0){
57747a47007SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"total number of function evaluations=%D,",tao->nfuncs);CHKERRQ(ierr);
57847a47007SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"                max: %D\n",tao->max_funcs);CHKERRQ(ierr);
579a7e14dcfSSatish Balay     }
580a7e14dcfSSatish Balay     if (tao->ngrads>0){
58147a47007SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"total number of gradient evaluations=%D,",tao->ngrads);CHKERRQ(ierr);
58247a47007SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"                max: %D\n",tao->max_funcs);CHKERRQ(ierr);
583a7e14dcfSSatish Balay     }
584a7e14dcfSSatish Balay     if (tao->nfuncgrads>0){
58547a47007SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"total number of function/gradient evaluations=%D,",tao->nfuncgrads);CHKERRQ(ierr);
58647a47007SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    (max: %D)\n",tao->max_funcs);CHKERRQ(ierr);
587a7e14dcfSSatish Balay     }
588a7e14dcfSSatish Balay     if (tao->nhess>0){
58947a47007SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"total number of Hessian evaluations=%D\n",tao->nhess);CHKERRQ(ierr);
590a7e14dcfSSatish Balay     }
591a7e14dcfSSatish Balay     /*  if (tao->linear_its>0){
59247a47007SBarry Smith      ierr = PetscViewerASCIIPrintf(viewer,"  total Krylov method iterations=%D\n",tao->linear_its);CHKERRQ(ierr);
593a7e14dcfSSatish Balay      }*/
594a7e14dcfSSatish Balay     if (tao->nconstraints>0){
59547a47007SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"total number of constraint function evaluations=%D\n",tao->nconstraints);CHKERRQ(ierr);
596a7e14dcfSSatish Balay     }
597a7e14dcfSSatish Balay     if (tao->njac>0){
59847a47007SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"total number of Jacobian evaluations=%D\n",tao->njac);CHKERRQ(ierr);
599a7e14dcfSSatish Balay     }
600a7e14dcfSSatish Balay 
601a7e14dcfSSatish Balay     if (tao->reason>0){
602a7e14dcfSSatish Balay       ierr = PetscViewerASCIIPrintf(viewer,    "Solution converged: ");CHKERRQ(ierr);
603a7e14dcfSSatish Balay       switch (tao->reason) {
604a7e14dcfSSatish Balay       case TAO_CONVERGED_FATOL:
605a7e14dcfSSatish Balay         ierr = PetscViewerASCIIPrintf(viewer,"estimated f(x)-f(X*) <= fatol\n");CHKERRQ(ierr);
606a7e14dcfSSatish Balay         break;
607a7e14dcfSSatish Balay       case TAO_CONVERGED_FRTOL:
608a7e14dcfSSatish Balay         ierr = PetscViewerASCIIPrintf(viewer,"estimated |f(x)-f(X*)|/|f(X*)| <= frtol\n");CHKERRQ(ierr);
609a7e14dcfSSatish Balay         break;
610a7e14dcfSSatish Balay       case TAO_CONVERGED_GATOL:
611a7e14dcfSSatish Balay         ierr = PetscViewerASCIIPrintf(viewer," ||g(X)|| <= gatol\n");CHKERRQ(ierr);
612a7e14dcfSSatish Balay         break;
613a7e14dcfSSatish Balay       case TAO_CONVERGED_GRTOL:
614a7e14dcfSSatish Balay         ierr = PetscViewerASCIIPrintf(viewer," ||g(X)||/|f(X)| <= grtol\n");CHKERRQ(ierr);
615a7e14dcfSSatish Balay         break;
616a7e14dcfSSatish Balay       case TAO_CONVERGED_GTTOL:
617a7e14dcfSSatish Balay         ierr = PetscViewerASCIIPrintf(viewer," ||g(X)||/||g(X0)|| <= gttol\n");CHKERRQ(ierr);
618a7e14dcfSSatish Balay         break;
619a7e14dcfSSatish Balay       case TAO_CONVERGED_STEPTOL:
620a7e14dcfSSatish Balay         ierr = PetscViewerASCIIPrintf(viewer," Steptol -- step size small\n");CHKERRQ(ierr);
621a7e14dcfSSatish Balay         break;
622a7e14dcfSSatish Balay       case TAO_CONVERGED_MINF:
623a7e14dcfSSatish Balay         ierr = PetscViewerASCIIPrintf(viewer," Minf --  f < fmin\n");CHKERRQ(ierr);
624a7e14dcfSSatish Balay         break;
625a7e14dcfSSatish Balay       case TAO_CONVERGED_USER:
626a7e14dcfSSatish Balay         ierr = PetscViewerASCIIPrintf(viewer," User Terminated\n");CHKERRQ(ierr);
627a7e14dcfSSatish Balay         break;
628a7e14dcfSSatish Balay       default:
629a7e14dcfSSatish Balay         ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
630a7e14dcfSSatish Balay         break;
631a7e14dcfSSatish Balay       }
632a7e14dcfSSatish Balay 
633a7e14dcfSSatish Balay     } else {
634a7e14dcfSSatish Balay       ierr = PetscViewerASCIIPrintf(viewer,"Solver terminated: %D",tao->reason);CHKERRQ(ierr);
635a7e14dcfSSatish Balay       switch (tao->reason) {
636a7e14dcfSSatish Balay       case TAO_DIVERGED_MAXITS:
63747a47007SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer," Maximum Iterations\n");CHKERRQ(ierr);
638a7e14dcfSSatish Balay         break;
639a7e14dcfSSatish Balay       case TAO_DIVERGED_NAN:
64047a47007SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer," NAN or Inf encountered\n");CHKERRQ(ierr);
641a7e14dcfSSatish Balay         break;
642a7e14dcfSSatish Balay       case TAO_DIVERGED_MAXFCN:
64347a47007SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer," Maximum Function Evaluations\n");CHKERRQ(ierr);
644a7e14dcfSSatish Balay         break;
645a7e14dcfSSatish Balay       case TAO_DIVERGED_LS_FAILURE:
64647a47007SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer," Line Search Failure\n");CHKERRQ(ierr);
647a7e14dcfSSatish Balay         break;
648a7e14dcfSSatish Balay       case TAO_DIVERGED_TR_REDUCTION:
64947a47007SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer," Trust Region too small\n");CHKERRQ(ierr);
650a7e14dcfSSatish Balay         break;
651a7e14dcfSSatish Balay       case TAO_DIVERGED_USER:
65247a47007SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer," User Terminated\n");CHKERRQ(ierr);
653a7e14dcfSSatish Balay         break;
654a7e14dcfSSatish Balay       default:
655a7e14dcfSSatish Balay         ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
656a7e14dcfSSatish Balay         break;
657a7e14dcfSSatish Balay       }
658a7e14dcfSSatish Balay     }
659a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
660a7e14dcfSSatish Balay   } else if (isstring) {
661a7e14dcfSSatish Balay     ierr = TaoGetType(tao,&type);CHKERRQ(ierr);
662a7e14dcfSSatish Balay     ierr = PetscViewerStringSPrintf(viewer," %-3.3s",type);CHKERRQ(ierr);
663a7e14dcfSSatish Balay   }
664a7e14dcfSSatish Balay   PetscFunctionReturn(0);
665a7e14dcfSSatish Balay }
666a7e14dcfSSatish Balay 
667a7e14dcfSSatish Balay #undef __FUNCT__
668a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetTolerances"
669a7e14dcfSSatish Balay /*@
670a7e14dcfSSatish Balay   TaoSetTolerances - Sets parameters used in TAO convergence tests
671a7e14dcfSSatish Balay 
672441846f8SBarry Smith   Logically collective on Tao
673a7e14dcfSSatish Balay 
674a7e14dcfSSatish Balay   Input Parameters:
675441846f8SBarry Smith + tao - the Tao context
676a7e14dcfSSatish Balay . fatol - absolute convergence tolerance
677a7e14dcfSSatish Balay . frtol - relative convergence tolerance
678a7e14dcfSSatish Balay . gatol - stop if norm of gradient is less than this
679a7e14dcfSSatish Balay . grtol - stop if relative norm of gradient is less than this
680a7e14dcfSSatish Balay - gttol - stop if norm of gradient is reduced by this factor
681a7e14dcfSSatish Balay 
682a7e14dcfSSatish Balay   Options Database Keys:
683a7e14dcfSSatish Balay + -tao_fatol <fatol> - Sets fatol
684a7e14dcfSSatish Balay . -tao_frtol <frtol> - Sets frtol
685a7e14dcfSSatish Balay . -tao_gatol <gatol> - Sets gatol
686a7e14dcfSSatish Balay . -tao_grtol <grtol> - Sets grtol
687a7e14dcfSSatish Balay - -tao_gttol <gttol> - Sets gttol
688a7e14dcfSSatish Balay 
689a7e14dcfSSatish Balay   Stopping Criteria:
690a7e14dcfSSatish Balay $ f(X) - f(X*) (estimated)            <= fatol
691a7e14dcfSSatish Balay $ |f(X) - f(X*)| (estimated) / |f(X)| <= frtol
692a7e14dcfSSatish Balay $ ||g(X)||                            <= gatol
693a7e14dcfSSatish Balay $ ||g(X)|| / |f(X)|                   <= grtol
694a7e14dcfSSatish Balay $ ||g(X)|| / ||g(X0)||                <= gttol
695a7e14dcfSSatish Balay 
696a7e14dcfSSatish Balay   Notes:
697a7e14dcfSSatish Balay   Use PETSC_DEFAULT to leave one or more tolerances unchanged.
698a7e14dcfSSatish Balay 
699a7e14dcfSSatish Balay   Level: beginner
700a7e14dcfSSatish Balay 
701a7e14dcfSSatish Balay .seealso: TaoGetTolerances()
702a7e14dcfSSatish Balay 
703a7e14dcfSSatish Balay @*/
704441846f8SBarry Smith PetscErrorCode TaoSetTolerances(Tao tao, PetscReal fatol, PetscReal frtol, PetscReal gatol, PetscReal grtol, PetscReal gttol)
705a7e14dcfSSatish Balay {
706a7e14dcfSSatish Balay   PetscErrorCode ierr;
70747a47007SBarry Smith 
708a7e14dcfSSatish Balay   PetscFunctionBegin;
709441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
710a7e14dcfSSatish Balay 
711a7e14dcfSSatish Balay   if (fatol != PETSC_DEFAULT) {
712a7e14dcfSSatish Balay     if (fatol<0) {
71347a47007SBarry Smith       ierr = PetscInfo(tao,"Tried to set negative fatol -- ignored.");CHKERRQ(ierr);
714a7e14dcfSSatish Balay     } else {
715a7e14dcfSSatish Balay       tao->fatol = PetscMax(0,fatol);
716a7e14dcfSSatish Balay     }
717a7e14dcfSSatish Balay   }
718a7e14dcfSSatish Balay   if (frtol != PETSC_DEFAULT) {
719a7e14dcfSSatish Balay     if (frtol<0) {
72047a47007SBarry Smith       ierr = PetscInfo(tao,"Tried to set negative frtol -- ignored.");CHKERRQ(ierr);
721a7e14dcfSSatish Balay     } else {
722a7e14dcfSSatish Balay       tao->frtol = PetscMax(0,frtol);
723a7e14dcfSSatish Balay     }
724a7e14dcfSSatish Balay   }
725a7e14dcfSSatish Balay 
726a7e14dcfSSatish Balay   if (gatol != PETSC_DEFAULT) {
727a7e14dcfSSatish Balay     if (gatol<0) {
72847a47007SBarry Smith       ierr = PetscInfo(tao,"Tried to set negative gatol -- ignored.");CHKERRQ(ierr);
729a7e14dcfSSatish Balay     } else {
730a7e14dcfSSatish Balay       tao->gatol = PetscMax(0,gatol);
731a7e14dcfSSatish Balay     }
732a7e14dcfSSatish Balay   }
733a7e14dcfSSatish Balay 
734a7e14dcfSSatish Balay   if (grtol != PETSC_DEFAULT) {
735a7e14dcfSSatish Balay     if (grtol<0) {
73647a47007SBarry Smith       ierr = PetscInfo(tao,"Tried to set negative grtol -- ignored.");CHKERRQ(ierr);
737a7e14dcfSSatish Balay     } else {
738a7e14dcfSSatish Balay       tao->grtol = PetscMax(0,grtol);
739a7e14dcfSSatish Balay     }
740a7e14dcfSSatish Balay   }
741a7e14dcfSSatish Balay 
742a7e14dcfSSatish Balay   if (gttol != PETSC_DEFAULT) {
743a7e14dcfSSatish Balay     if (gttol<0) {
74447a47007SBarry Smith       ierr = PetscInfo(tao,"Tried to set negative gttol -- ignored.");CHKERRQ(ierr);
745a7e14dcfSSatish Balay     } else {
746a7e14dcfSSatish Balay       tao->gttol = PetscMax(0,gttol);
747a7e14dcfSSatish Balay     }
748a7e14dcfSSatish Balay   }
749a7e14dcfSSatish Balay   PetscFunctionReturn(0);
750a7e14dcfSSatish Balay }
751a7e14dcfSSatish Balay 
752a7e14dcfSSatish Balay #undef __FUNCT__
753a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetConstraintTolerances"
754a7e14dcfSSatish Balay /*@
75558417fe7SBarry Smith   TaoSetConstraintTolerances - Sets constraint tolerance parameters used in TAO  convergence tests
756a7e14dcfSSatish Balay 
757441846f8SBarry Smith   Logically collective on Tao
758a7e14dcfSSatish Balay 
759a7e14dcfSSatish Balay   Input Parameters:
760441846f8SBarry Smith + tao - the Tao context
761a7e14dcfSSatish Balay . catol - absolute constraint tolerance, constraint norm must be less than catol for used for fatol, gatol convergence criteria
762a7e14dcfSSatish Balay - crtol - relative contraint tolerance, constraint norm must be less than crtol for used for fatol, gatol, gttol convergence criteria
763a7e14dcfSSatish Balay 
764a7e14dcfSSatish Balay   Options Database Keys:
765a7e14dcfSSatish Balay + -tao_catol <catol> - Sets catol
766a7e14dcfSSatish Balay - -tao_crtol <crtol> - Sets crtol
767a7e14dcfSSatish Balay 
768a7e14dcfSSatish Balay   Level: intermediate
769a7e14dcfSSatish Balay 
77058417fe7SBarry Smith .seealso: TaoGetTolerances(), TaoGetConstraintTolerances(), TaoSetTolerances()
771a7e14dcfSSatish Balay 
772a7e14dcfSSatish Balay @*/
773441846f8SBarry Smith PetscErrorCode TaoSetConstraintTolerances(Tao tao, PetscReal catol, PetscReal crtol)
774a7e14dcfSSatish Balay {
775a7e14dcfSSatish Balay   PetscErrorCode ierr;
77647a47007SBarry Smith 
777a7e14dcfSSatish Balay   PetscFunctionBegin;
778441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
779a7e14dcfSSatish Balay 
780a7e14dcfSSatish Balay   if (catol != PETSC_DEFAULT) {
781a7e14dcfSSatish Balay     if (catol<0) {
78247a47007SBarry Smith       ierr = PetscInfo(tao,"Tried to set negative catol -- ignored.");CHKERRQ(ierr);
783a7e14dcfSSatish Balay     } else {
784a7e14dcfSSatish Balay       tao->catol = PetscMax(0,catol);
785a7e14dcfSSatish Balay     }
786a7e14dcfSSatish Balay   }
787a7e14dcfSSatish Balay 
788a7e14dcfSSatish Balay   if (crtol != PETSC_DEFAULT) {
789a7e14dcfSSatish Balay     if (crtol<0) {
79047a47007SBarry Smith       ierr = PetscInfo(tao,"Tried to set negative crtol -- ignored.");CHKERRQ(ierr);
791a7e14dcfSSatish Balay     } else {
792a7e14dcfSSatish Balay       tao->crtol = PetscMax(0,crtol);
793a7e14dcfSSatish Balay     }
794a7e14dcfSSatish Balay   }
795a7e14dcfSSatish Balay   PetscFunctionReturn(0);
796a7e14dcfSSatish Balay }
797a7e14dcfSSatish Balay 
798a7e14dcfSSatish Balay #undef __FUNCT__
79958417fe7SBarry Smith #define __FUNCT__ "TaoGetConstraintTolerances"
80058417fe7SBarry Smith /*@
80158417fe7SBarry Smith   TaoGetConstraintTolerances - Gets constraint tolerance parameters used in TAO  convergence tests
80258417fe7SBarry Smith 
80358417fe7SBarry Smith   Not ollective
80458417fe7SBarry Smith 
80558417fe7SBarry Smith   Input Parameter:
80658417fe7SBarry Smith . tao - the Tao context
80758417fe7SBarry Smith 
80858417fe7SBarry Smith   Output Parameter:
80958417fe7SBarry Smith + catol - absolute constraint tolerance, constraint norm must be less than catol for used for fatol, gatol convergence criteria
81058417fe7SBarry Smith - crtol - relative contraint tolerance, constraint norm must be less than crtol for used for fatol, gatol, gttol convergence criteria
81158417fe7SBarry Smith 
81258417fe7SBarry Smith   Level: intermediate
81358417fe7SBarry Smith 
81458417fe7SBarry Smith .seealso: TaoGetTolerances(), TaoSetTolerances(), TaoSetConstraintTolerances()
81558417fe7SBarry Smith 
81658417fe7SBarry Smith @*/
81758417fe7SBarry Smith PetscErrorCode TaoGetConstraintTolerances(Tao tao, PetscReal *catol, PetscReal *crtol)
81858417fe7SBarry Smith {
81958417fe7SBarry Smith   PetscFunctionBegin;
82058417fe7SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
82158417fe7SBarry Smith   if (catol) *catol = tao->catol;
82258417fe7SBarry Smith   if (crtol) *crtol = tao->crtol;
82358417fe7SBarry Smith   PetscFunctionReturn(0);
82458417fe7SBarry Smith }
82558417fe7SBarry Smith 
82658417fe7SBarry Smith #undef __FUNCT__
827a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetFunctionLowerBound"
828a7e14dcfSSatish Balay /*@
829a7e14dcfSSatish Balay    TaoSetFunctionLowerBound - Sets a bound on the solution objective value.
830a7e14dcfSSatish Balay    When an approximate solution with an objective value below this number
831a7e14dcfSSatish Balay    has been found, the solver will terminate.
832a7e14dcfSSatish Balay 
833441846f8SBarry Smith    Logically Collective on Tao
834a7e14dcfSSatish Balay 
835a7e14dcfSSatish Balay    Input Parameters:
836441846f8SBarry Smith +  tao - the Tao solver context
837a7e14dcfSSatish Balay -  fmin - the tolerance
838a7e14dcfSSatish Balay 
839a7e14dcfSSatish Balay    Options Database Keys:
840a7e14dcfSSatish Balay .    -tao_fmin <fmin> - sets the minimum function value
841a7e14dcfSSatish Balay 
842a7e14dcfSSatish Balay    Level: intermediate
843a7e14dcfSSatish Balay 
844a7e14dcfSSatish Balay .seealso: TaoSetTolerances()
845a7e14dcfSSatish Balay @*/
846441846f8SBarry Smith PetscErrorCode TaoSetFunctionLowerBound(Tao tao,PetscReal fmin)
847a7e14dcfSSatish Balay {
848a7e14dcfSSatish Balay   PetscFunctionBegin;
849441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
850a7e14dcfSSatish Balay   tao->fmin = fmin;
851a7e14dcfSSatish Balay   PetscFunctionReturn(0);
852a7e14dcfSSatish Balay }
853a7e14dcfSSatish Balay 
854a7e14dcfSSatish Balay #undef __FUNCT__
855a7e14dcfSSatish Balay #define __FUNCT__ "TaoGetFunctionLowerBound"
856a7e14dcfSSatish Balay /*@
857a7e14dcfSSatish Balay    TaoGetFunctionLowerBound - Sets a bound on the solution objective value.
858a7e14dcfSSatish Balay    When an approximate solution with an objective value below this number
859a7e14dcfSSatish Balay    has been found, the solver will terminate.
860a7e14dcfSSatish Balay 
861441846f8SBarry Smith    Not collective on Tao
862a7e14dcfSSatish Balay 
863a7e14dcfSSatish Balay    Input Parameters:
864441846f8SBarry Smith .  tao - the Tao solver context
865a7e14dcfSSatish Balay 
866a7e14dcfSSatish Balay    OutputParameters:
867a7e14dcfSSatish Balay .  fmin - the minimum function value
868a7e14dcfSSatish Balay 
869a7e14dcfSSatish Balay    Level: intermediate
870a7e14dcfSSatish Balay 
871a7e14dcfSSatish Balay .seealso: TaoSetFunctionLowerBound()
872a7e14dcfSSatish Balay @*/
873441846f8SBarry Smith PetscErrorCode TaoGetFunctionLowerBound(Tao tao,PetscReal *fmin)
874a7e14dcfSSatish Balay {
875a7e14dcfSSatish Balay   PetscFunctionBegin;
876441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
877a7e14dcfSSatish Balay   *fmin = tao->fmin;
878a7e14dcfSSatish Balay   PetscFunctionReturn(0);
879a7e14dcfSSatish Balay }
880a7e14dcfSSatish Balay 
881a7e14dcfSSatish Balay #undef __FUNCT__
882a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetMaximumFunctionEvaluations"
883a7e14dcfSSatish Balay /*@
884a7e14dcfSSatish Balay    TaoSetMaximumFunctionEvaluations - Sets a maximum number of
885a7e14dcfSSatish Balay    function evaluations.
886a7e14dcfSSatish Balay 
887441846f8SBarry Smith    Logically Collective on Tao
888a7e14dcfSSatish Balay 
889a7e14dcfSSatish Balay    Input Parameters:
890441846f8SBarry Smith +  tao - the Tao solver context
891a7e14dcfSSatish Balay -  nfcn - the maximum number of function evaluations (>=0)
892a7e14dcfSSatish Balay 
893a7e14dcfSSatish Balay    Options Database Keys:
894a7e14dcfSSatish Balay .    -tao_max_funcs <nfcn> - sets the maximum number of function evaluations
895a7e14dcfSSatish Balay 
896a7e14dcfSSatish Balay    Level: intermediate
897a7e14dcfSSatish Balay 
898a7e14dcfSSatish Balay .seealso: TaoSetTolerances(), TaoSetMaximumIterations()
899a7e14dcfSSatish Balay @*/
900a7e14dcfSSatish Balay 
901441846f8SBarry Smith PetscErrorCode TaoSetMaximumFunctionEvaluations(Tao tao,PetscInt nfcn)
902a7e14dcfSSatish Balay {
903a7e14dcfSSatish Balay   PetscFunctionBegin;
904441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
90547a47007SBarry Smith   tao->max_funcs = PetscMax(0,nfcn);
906a7e14dcfSSatish Balay   PetscFunctionReturn(0);
907a7e14dcfSSatish Balay }
908a7e14dcfSSatish Balay 
909a7e14dcfSSatish Balay #undef __FUNCT__
910a7e14dcfSSatish Balay #define __FUNCT__ "TaoGetMaximumFunctionEvaluations"
911a7e14dcfSSatish Balay /*@
912a7e14dcfSSatish Balay    TaoGetMaximumFunctionEvaluations - Sets a maximum number of
913a7e14dcfSSatish Balay    function evaluations.
914a7e14dcfSSatish Balay 
915a7e14dcfSSatish Balay    Not Collective
916a7e14dcfSSatish Balay 
917a7e14dcfSSatish Balay    Input Parameters:
918441846f8SBarry Smith .  tao - the Tao solver context
919a7e14dcfSSatish Balay 
920a7e14dcfSSatish Balay    Output Parameters:
921a7e14dcfSSatish Balay .  nfcn - the maximum number of function evaluations
922a7e14dcfSSatish Balay 
923a7e14dcfSSatish Balay    Level: intermediate
924a7e14dcfSSatish Balay 
925a7e14dcfSSatish Balay .seealso: TaoSetMaximumFunctionEvaluations(), TaoGetMaximumIterations()
926a7e14dcfSSatish Balay @*/
927a7e14dcfSSatish Balay 
928441846f8SBarry Smith PetscErrorCode TaoGetMaximumFunctionEvaluations(Tao tao,PetscInt *nfcn)
929a7e14dcfSSatish Balay {
930a7e14dcfSSatish Balay   PetscFunctionBegin;
931441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
932a7e14dcfSSatish Balay   *nfcn = tao->max_funcs;
933a7e14dcfSSatish Balay   PetscFunctionReturn(0);
934a7e14dcfSSatish Balay }
935a7e14dcfSSatish Balay 
936a7e14dcfSSatish Balay #undef __FUNCT__
937a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetMaximumIterations"
938a7e14dcfSSatish Balay /*@
939a7e14dcfSSatish Balay    TaoSetMaximumIterations - Sets a maximum number of iterates.
940a7e14dcfSSatish Balay 
941441846f8SBarry Smith    Logically Collective on Tao
942a7e14dcfSSatish Balay 
943a7e14dcfSSatish Balay    Input Parameters:
944441846f8SBarry Smith +  tao - the Tao solver context
945a7e14dcfSSatish Balay -  maxits - the maximum number of iterates (>=0)
946a7e14dcfSSatish Balay 
947a7e14dcfSSatish Balay    Options Database Keys:
948a7e14dcfSSatish Balay .    -tao_max_it <its> - sets the maximum number of iterations
949a7e14dcfSSatish Balay 
950a7e14dcfSSatish Balay    Level: intermediate
951a7e14dcfSSatish Balay 
952a7e14dcfSSatish Balay .seealso: TaoSetTolerances(), TaoSetMaximumFunctionEvaluations()
953a7e14dcfSSatish Balay @*/
954441846f8SBarry Smith PetscErrorCode TaoSetMaximumIterations(Tao tao,PetscInt maxits)
955a7e14dcfSSatish Balay {
956a7e14dcfSSatish Balay   PetscFunctionBegin;
957441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
95847a47007SBarry Smith   tao->max_it = PetscMax(0,maxits);
959a7e14dcfSSatish Balay   PetscFunctionReturn(0);
960a7e14dcfSSatish Balay }
961a7e14dcfSSatish Balay 
962a7e14dcfSSatish Balay #undef __FUNCT__
963a7e14dcfSSatish Balay #define __FUNCT__ "TaoGetMaximumIterations"
964a7e14dcfSSatish Balay /*@
965a7e14dcfSSatish Balay    TaoGetMaximumIterations - Sets a maximum number of iterates.
966a7e14dcfSSatish Balay 
967a7e14dcfSSatish Balay    Not Collective
968a7e14dcfSSatish Balay 
969a7e14dcfSSatish Balay    Input Parameters:
970441846f8SBarry Smith .  tao - the Tao solver context
971a7e14dcfSSatish Balay 
972a7e14dcfSSatish Balay    Output Parameters:
973a7e14dcfSSatish Balay .  maxits - the maximum number of iterates
974a7e14dcfSSatish Balay 
975a7e14dcfSSatish Balay    Level: intermediate
976a7e14dcfSSatish Balay 
977a7e14dcfSSatish Balay .seealso: TaoSetMaximumIterations(), TaoGetMaximumFunctionEvaluations()
978a7e14dcfSSatish Balay @*/
979441846f8SBarry Smith PetscErrorCode TaoGetMaximumIterations(Tao tao,PetscInt *maxits)
980a7e14dcfSSatish Balay {
981a7e14dcfSSatish Balay   PetscFunctionBegin;
982441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
983a7e14dcfSSatish Balay   *maxits = tao->max_it;
984a7e14dcfSSatish Balay   PetscFunctionReturn(0);
985a7e14dcfSSatish Balay }
986a7e14dcfSSatish Balay 
987a7e14dcfSSatish Balay #undef __FUNCT__
988a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetInitialTrustRegionRadius"
989a7e14dcfSSatish Balay /*@
990a7e14dcfSSatish Balay    TaoSetInitialTrustRegionRadius - Sets the initial trust region radius.
991a7e14dcfSSatish Balay 
992441846f8SBarry Smith    Logically collective on Tao
993a7e14dcfSSatish Balay 
994a7e14dcfSSatish Balay    Input Parameter:
995a7e14dcfSSatish Balay +  tao - a TAO optimization solver
996a7e14dcfSSatish Balay -  radius - the trust region radius
997a7e14dcfSSatish Balay 
998a7e14dcfSSatish Balay    Level: intermediate
999a7e14dcfSSatish Balay 
1000a7e14dcfSSatish Balay    Options Database Key:
1001a7e14dcfSSatish Balay .  -tao_trust0 <t0> - sets initial trust region radius
1002a7e14dcfSSatish Balay 
1003a7e14dcfSSatish Balay .seealso: TaoGetTrustRegionRadius(), TaoSetTrustRegionTolerance()
1004a7e14dcfSSatish Balay @*/
1005441846f8SBarry Smith PetscErrorCode TaoSetInitialTrustRegionRadius(Tao tao, PetscReal radius)
1006a7e14dcfSSatish Balay {
1007a7e14dcfSSatish Balay   PetscFunctionBegin;
1008441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1009a7e14dcfSSatish Balay   tao->trust0 = PetscMax(0.0,radius);
1010a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1011a7e14dcfSSatish Balay }
1012a7e14dcfSSatish Balay 
1013a7e14dcfSSatish Balay #undef __FUNCT__
1014a7e14dcfSSatish Balay #define __FUNCT__ "TaoGetInitialTrustRegionRadius"
1015a7e14dcfSSatish Balay /*@
1016a7e14dcfSSatish Balay    TaoGetInitialTrustRegionRadius - Sets the initial trust region radius.
1017a7e14dcfSSatish Balay 
1018a7e14dcfSSatish Balay    Not Collective
1019a7e14dcfSSatish Balay 
1020a7e14dcfSSatish Balay    Input Parameter:
1021a7e14dcfSSatish Balay .  tao - a TAO optimization solver
1022a7e14dcfSSatish Balay 
1023a7e14dcfSSatish Balay    Output Parameter:
1024a7e14dcfSSatish Balay .  radius - the trust region radius
1025a7e14dcfSSatish Balay 
1026a7e14dcfSSatish Balay    Level: intermediate
1027a7e14dcfSSatish Balay 
1028a7e14dcfSSatish Balay .seealso: TaoSetInitialTrustRegionRadius(), TaoGetCurrentTrustRegionRadius()
1029a7e14dcfSSatish Balay @*/
1030441846f8SBarry Smith PetscErrorCode TaoGetInitialTrustRegionRadius(Tao tao, PetscReal *radius)
1031a7e14dcfSSatish Balay {
1032a7e14dcfSSatish Balay   PetscFunctionBegin;
1033441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1034a7e14dcfSSatish Balay   *radius = tao->trust0;
1035a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1036a7e14dcfSSatish Balay }
1037a7e14dcfSSatish Balay 
1038a7e14dcfSSatish Balay #undef __FUNCT__
1039a7e14dcfSSatish Balay #define __FUNCT__ "TaoGetCurrentTrustRegionRadius"
1040a7e14dcfSSatish Balay /*@
1041a7e14dcfSSatish Balay    TaoGetCurrentTrustRegionRadius - Gets the current trust region radius.
1042a7e14dcfSSatish Balay 
1043a7e14dcfSSatish Balay    Not Collective
1044a7e14dcfSSatish Balay 
1045a7e14dcfSSatish Balay    Input Parameter:
1046a7e14dcfSSatish Balay .  tao - a TAO optimization solver
1047a7e14dcfSSatish Balay 
1048a7e14dcfSSatish Balay    Output Parameter:
1049a7e14dcfSSatish Balay .  radius - the trust region radius
1050a7e14dcfSSatish Balay 
1051a7e14dcfSSatish Balay    Level: intermediate
1052a7e14dcfSSatish Balay 
1053a7e14dcfSSatish Balay .seealso: TaoSetInitialTrustRegionRadius(), TaoGetInitialTrustRegionRadius()
1054a7e14dcfSSatish Balay @*/
1055441846f8SBarry Smith PetscErrorCode TaoGetCurrentTrustRegionRadius(Tao tao, PetscReal *radius)
1056a7e14dcfSSatish Balay {
1057a7e14dcfSSatish Balay   PetscFunctionBegin;
1058441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1059a7e14dcfSSatish Balay   *radius = tao->trust;
1060a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1061a7e14dcfSSatish Balay }
1062a7e14dcfSSatish Balay 
1063a7e14dcfSSatish Balay #undef __FUNCT__
1064a7e14dcfSSatish Balay #define __FUNCT__ "TaoGetTolerances"
1065a7e14dcfSSatish Balay /*@
1066a7e14dcfSSatish Balay   TaoGetTolerances - gets the current values of tolerances
1067a7e14dcfSSatish Balay 
1068a7e14dcfSSatish Balay   Not Collective
1069a7e14dcfSSatish Balay 
1070a7e14dcfSSatish Balay   Input Parameters:
1071441846f8SBarry Smith . tao - the Tao context
1072a7e14dcfSSatish Balay 
1073a7e14dcfSSatish Balay   Output Parameters:
1074a7e14dcfSSatish Balay + fatol - absolute convergence tolerance
1075a7e14dcfSSatish Balay . frtol - relative convergence tolerance
1076a7e14dcfSSatish Balay . gatol - stop if norm of gradient is less than this
1077a7e14dcfSSatish Balay . grtol - stop if relative norm of gradient is less than this
1078a7e14dcfSSatish Balay - gttol - stop if norm of gradient is reduced by a this factor
1079a7e14dcfSSatish Balay 
108047a47007SBarry Smith   Note: NULL can be used as an argument if not all tolerances values are needed
1081a7e14dcfSSatish Balay 
1082a7e14dcfSSatish Balay .seealso TaoSetTolerances()
1083a7e14dcfSSatish Balay 
1084a7e14dcfSSatish Balay   Level: intermediate
1085a7e14dcfSSatish Balay @*/
1086441846f8SBarry Smith PetscErrorCode TaoGetTolerances(Tao tao, PetscReal *fatol, PetscReal *frtol, PetscReal *gatol, PetscReal *grtol, PetscReal *gttol)
1087a7e14dcfSSatish Balay {
1088a7e14dcfSSatish Balay   PetscFunctionBegin;
1089441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1090a7e14dcfSSatish Balay   if (fatol) *fatol=tao->fatol;
1091a7e14dcfSSatish Balay   if (frtol) *frtol=tao->frtol;
1092a7e14dcfSSatish Balay   if (gatol) *gatol=tao->gatol;
1093a7e14dcfSSatish Balay   if (grtol) *grtol=tao->grtol;
1094a7e14dcfSSatish Balay   if (gttol) *gttol=tao->gttol;
1095a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1096a7e14dcfSSatish Balay }
1097a7e14dcfSSatish Balay 
1098a7e14dcfSSatish Balay #undef __FUNCT__
1099a7e14dcfSSatish Balay #define __FUNCT__ "TaoGetKSP"
1100a7e14dcfSSatish Balay /*@
1101a7e14dcfSSatish Balay   TaoGetKSP - Gets the linear solver used by the optimization solver.
1102a7e14dcfSSatish Balay   Application writers should use TaoGetKSP if they need direct access
1103a7e14dcfSSatish Balay   to the PETSc KSP object.
1104a7e14dcfSSatish Balay 
1105a7e14dcfSSatish Balay   Not Collective
1106a7e14dcfSSatish Balay 
1107a7e14dcfSSatish Balay    Input Parameters:
1108a7e14dcfSSatish Balay .  tao - the TAO solver
1109a7e14dcfSSatish Balay 
1110a7e14dcfSSatish Balay    Output Parameters:
1111a7e14dcfSSatish Balay .  ksp - the KSP linear solver used in the optimization solver
1112a7e14dcfSSatish Balay 
1113a7e14dcfSSatish Balay    Level: intermediate
1114a7e14dcfSSatish Balay 
1115a7e14dcfSSatish Balay @*/
1116441846f8SBarry Smith PetscErrorCode TaoGetKSP(Tao tao, KSP *ksp)
111747a47007SBarry Smith {
1118a7e14dcfSSatish Balay   PetscFunctionBegin;
1119a7e14dcfSSatish Balay   *ksp = tao->ksp;
1120a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1121a7e14dcfSSatish Balay }
1122a7e14dcfSSatish Balay 
1123a7e14dcfSSatish Balay #undef __FUNCT__
1124a7e14dcfSSatish Balay #define __FUNCT__ "TaoGetLineSearch"
1125a7e14dcfSSatish Balay /*@
1126a7e14dcfSSatish Balay   TaoGetLineSearch - Gets the line search used by the optimization solver.
1127a7e14dcfSSatish Balay   Application writers should use TaoGetLineSearch if they need direct access
1128a7e14dcfSSatish Balay   to the TaoLineSearch object.
1129a7e14dcfSSatish Balay 
1130a7e14dcfSSatish Balay   Not Collective
1131a7e14dcfSSatish Balay 
1132a7e14dcfSSatish Balay    Input Parameters:
1133a7e14dcfSSatish Balay .  tao - the TAO solver
1134a7e14dcfSSatish Balay 
1135a7e14dcfSSatish Balay    Output Parameters:
1136a7e14dcfSSatish Balay .  ls - the line search used in the optimization solver
1137a7e14dcfSSatish Balay 
1138a7e14dcfSSatish Balay    Level: intermediate
1139a7e14dcfSSatish Balay 
1140a7e14dcfSSatish Balay @*/
1141441846f8SBarry Smith PetscErrorCode TaoGetLineSearch(Tao tao, TaoLineSearch *ls)
114247a47007SBarry Smith {
1143a7e14dcfSSatish Balay   PetscFunctionBegin;
1144a7e14dcfSSatish Balay   *ls = tao->linesearch;
1145a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1146a7e14dcfSSatish Balay }
1147a7e14dcfSSatish Balay 
1148a7e14dcfSSatish Balay #undef __FUNCT__
1149a7e14dcfSSatish Balay #define __FUNCT__ "TaoAddLineSearchCounts"
1150a7e14dcfSSatish Balay /*@
1151a7e14dcfSSatish Balay   TaoAddLineSearchCounts - Adds the number of function evaluations spent
1152a7e14dcfSSatish Balay   in the line search to the running total.
1153a7e14dcfSSatish Balay 
1154a7e14dcfSSatish Balay    Input Parameters:
1155a7e14dcfSSatish Balay +  tao - the TAO solver
1156a7e14dcfSSatish Balay -  ls - the line search used in the optimization solver
1157a7e14dcfSSatish Balay 
1158a7e14dcfSSatish Balay    Level: developer
1159a7e14dcfSSatish Balay 
1160a7e14dcfSSatish Balay .seealso: TaoLineSearchApply()
1161a7e14dcfSSatish Balay @*/
1162441846f8SBarry Smith PetscErrorCode TaoAddLineSearchCounts(Tao tao)
116347a47007SBarry Smith {
1164a7e14dcfSSatish Balay   PetscErrorCode ierr;
1165a7e14dcfSSatish Balay   PetscBool      flg;
1166a7e14dcfSSatish Balay   PetscInt       nfeval,ngeval,nfgeval;
116747a47007SBarry Smith 
1168a7e14dcfSSatish Balay   PetscFunctionBegin;
1169441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1170a7e14dcfSSatish Balay   if (tao->linesearch) {
1171441846f8SBarry Smith     ierr = TaoLineSearchIsUsingTaoRoutines(tao->linesearch,&flg);
1172a7e14dcfSSatish Balay     if (flg == PETSC_FALSE) {
117347a47007SBarry Smith       ierr = TaoLineSearchGetNumberFunctionEvaluations(tao->linesearch,&nfeval,&ngeval,&nfgeval);CHKERRQ(ierr);
1174a7e14dcfSSatish Balay       tao->nfuncs+=nfeval;
1175a7e14dcfSSatish Balay       tao->ngrads+=ngeval;
1176a7e14dcfSSatish Balay       tao->nfuncgrads+=nfgeval;
1177a7e14dcfSSatish Balay     }
1178a7e14dcfSSatish Balay   }
1179a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1180a7e14dcfSSatish Balay }
1181a7e14dcfSSatish Balay 
1182a7e14dcfSSatish Balay #undef __FUNCT__
1183a7e14dcfSSatish Balay #define __FUNCT__ "TaoGetSolutionVector"
1184a7e14dcfSSatish Balay /*@
1185a7e14dcfSSatish Balay   TaoGetSolutionVector - Returns the vector with the current TAO solution
1186a7e14dcfSSatish Balay 
1187a7e14dcfSSatish Balay   Not Collective
1188a7e14dcfSSatish Balay 
1189a7e14dcfSSatish Balay   Input Parameter:
1190441846f8SBarry Smith . tao - the Tao context
1191a7e14dcfSSatish Balay 
1192a7e14dcfSSatish Balay   Output Parameter:
1193a7e14dcfSSatish Balay . X - the current solution
1194a7e14dcfSSatish Balay 
1195a7e14dcfSSatish Balay   Level: intermediate
1196a7e14dcfSSatish Balay 
1197a7e14dcfSSatish Balay   Note:  The returned vector will be the same object that was passed into TaoSetInitialVector()
1198a7e14dcfSSatish Balay @*/
1199441846f8SBarry Smith PetscErrorCode TaoGetSolutionVector(Tao tao, Vec *X)
1200a7e14dcfSSatish Balay {
1201a7e14dcfSSatish Balay   PetscFunctionBegin;
1202441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1203a7e14dcfSSatish Balay   *X = tao->solution;
1204a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1205a7e14dcfSSatish Balay }
1206a7e14dcfSSatish Balay 
1207a7e14dcfSSatish Balay #undef __FUNCT__
1208a7e14dcfSSatish Balay #define __FUNCT__ "TaoGetGradientVector"
1209a7e14dcfSSatish Balay /*@
1210a7e14dcfSSatish Balay   TaoGetGradientVector - Returns the vector with the current TAO gradient
1211a7e14dcfSSatish Balay 
1212a7e14dcfSSatish Balay   Not Collective
1213a7e14dcfSSatish Balay 
1214a7e14dcfSSatish Balay   Input Parameter:
1215441846f8SBarry Smith . tao - the Tao context
1216a7e14dcfSSatish Balay 
1217a7e14dcfSSatish Balay   Output Parameter:
1218a7e14dcfSSatish Balay . G - the current solution
1219a7e14dcfSSatish Balay 
1220a7e14dcfSSatish Balay   Level: intermediate
1221a7e14dcfSSatish Balay @*/
1222441846f8SBarry Smith PetscErrorCode TaoGetGradientVector(Tao tao, Vec *G)
1223a7e14dcfSSatish Balay {
1224a7e14dcfSSatish Balay   PetscFunctionBegin;
1225441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1226a7e14dcfSSatish Balay   *G = tao->gradient;
1227a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1228a7e14dcfSSatish Balay }
1229a7e14dcfSSatish Balay 
1230a7e14dcfSSatish Balay #undef __FUNCT__
1231a7e14dcfSSatish Balay #define __FUNCT__ "TaoResetStatistics"
1232a7e14dcfSSatish Balay /*@
1233a7e14dcfSSatish Balay    TaoResetStatistics - Initialize the statistics used by TAO for all of the solvers.
1234a7e14dcfSSatish Balay    These statistics include the iteration number, residual norms, and convergence status.
1235a7e14dcfSSatish Balay    This routine gets called before solving each optimization problem.
1236a7e14dcfSSatish Balay 
1237441846f8SBarry Smith    Collective on Tao
1238a7e14dcfSSatish Balay 
1239a7e14dcfSSatish Balay    Input Parameters:
1240441846f8SBarry Smith .  solver - the Tao context
1241a7e14dcfSSatish Balay 
1242a7e14dcfSSatish Balay    Level: developer
1243a7e14dcfSSatish Balay 
1244a7e14dcfSSatish Balay .seealso: TaoCreate(), TaoSolve()
1245a7e14dcfSSatish Balay @*/
1246441846f8SBarry Smith PetscErrorCode TaoResetStatistics(Tao tao)
1247a7e14dcfSSatish Balay {
1248a7e14dcfSSatish Balay   PetscFunctionBegin;
1249441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1250a7e14dcfSSatish Balay   tao->niter        = 0;
1251a7e14dcfSSatish Balay   tao->nfuncs       = 0;
1252a7e14dcfSSatish Balay   tao->nfuncgrads   = 0;
1253a7e14dcfSSatish Balay   tao->ngrads       = 0;
1254a7e14dcfSSatish Balay   tao->nhess        = 0;
1255a7e14dcfSSatish Balay   tao->njac         = 0;
1256a7e14dcfSSatish Balay   tao->nconstraints = 0;
1257a7e14dcfSSatish Balay   tao->ksp_its      = 0;
1258a7e14dcfSSatish Balay   tao->reason       = TAO_CONTINUE_ITERATING;
1259a7e14dcfSSatish Balay   tao->residual     = 0.0;
1260a7e14dcfSSatish Balay   tao->cnorm        = 0.0;
1261a7e14dcfSSatish Balay   tao->step         = 0.0;
1262a7e14dcfSSatish Balay   tao->lsflag       = PETSC_FALSE;
1263a7e14dcfSSatish Balay   if (tao->hist_reset) tao->hist_len=0;
1264a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1265a7e14dcfSSatish Balay }
1266a7e14dcfSSatish Balay 
1267a7e14dcfSSatish Balay #undef __FUNCT__
1268a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetConvergenceTest"
1269a7e14dcfSSatish Balay /*@C
1270a7e14dcfSSatish Balay   TaoSetConvergenceTest - Sets the function that is to be used to test
1271a7e14dcfSSatish Balay   for convergence o fthe iterative minimization solution.  The new convergence
1272a7e14dcfSSatish Balay   testing routine will replace TAO's default convergence test.
1273a7e14dcfSSatish Balay 
1274441846f8SBarry Smith   Logically Collective on Tao
1275a7e14dcfSSatish Balay 
1276a7e14dcfSSatish Balay   Input Parameters:
1277441846f8SBarry Smith + tao - the Tao object
1278a7e14dcfSSatish Balay . conv - the routine to test for convergence
1279a7e14dcfSSatish Balay - ctx - [optional] context for private data for the convergence routine
128047a47007SBarry Smith         (may be NULL)
1281a7e14dcfSSatish Balay 
1282a7e14dcfSSatish Balay   Calling sequence of conv:
1283441846f8SBarry Smith $   PetscErrorCode conv(Tao tao, void *ctx)
1284a7e14dcfSSatish Balay 
1285441846f8SBarry Smith + tao - the Tao object
1286a7e14dcfSSatish Balay - ctx - [optional] convergence context
1287a7e14dcfSSatish Balay 
1288*e4cb33bbSBarry Smith   Note: The new convergence testing routine should call TaoSetConvergedReason().
1289a7e14dcfSSatish Balay 
1290a7e14dcfSSatish Balay   Level: advanced
1291a7e14dcfSSatish Balay 
1292*e4cb33bbSBarry Smith .seealso: TaoSetConvergedReason(), TaoGetSolutionStatus(), TaoGetTolerances(), TaoSetMonitor
1293a7e14dcfSSatish Balay 
1294a7e14dcfSSatish Balay @*/
1295441846f8SBarry Smith PetscErrorCode TaoSetConvergenceTest(Tao tao, PetscErrorCode (*conv)(Tao,void*), void *ctx)
1296a7e14dcfSSatish Balay {
1297a7e14dcfSSatish Balay   PetscFunctionBegin;
1298441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1299a7e14dcfSSatish Balay   (tao)->ops->convergencetest = conv;
1300a7e14dcfSSatish Balay   (tao)->cnvP = ctx;
1301a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1302a7e14dcfSSatish Balay }
1303a7e14dcfSSatish Balay 
1304a7e14dcfSSatish Balay #undef __FUNCT__
1305a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetMonitor"
1306a7e14dcfSSatish Balay /*@C
1307a7e14dcfSSatish Balay    TaoSetMonitor - Sets an ADDITIONAL function that is to be used at every
1308a7e14dcfSSatish Balay    iteration of the solver to display the iteration's
1309a7e14dcfSSatish Balay    progress.
1310a7e14dcfSSatish Balay 
1311441846f8SBarry Smith    Logically Collective on Tao
1312a7e14dcfSSatish Balay 
1313a7e14dcfSSatish Balay    Input Parameters:
1314441846f8SBarry Smith +  tao - the Tao solver context
1315a7e14dcfSSatish Balay .  mymonitor - monitoring routine
1316a7e14dcfSSatish Balay -  mctx - [optional] user-defined context for private data for the
131747a47007SBarry Smith           monitor routine (may be NULL)
1318a7e14dcfSSatish Balay 
1319a7e14dcfSSatish Balay    Calling sequence of mymonitor:
1320441846f8SBarry Smith $     int mymonitor(Tao tao,void *mctx)
1321a7e14dcfSSatish Balay 
1322441846f8SBarry Smith +    tao - the Tao solver context
1323a7e14dcfSSatish Balay -    mctx - [optional] monitoring context
1324a7e14dcfSSatish Balay 
1325a7e14dcfSSatish Balay 
1326a7e14dcfSSatish Balay    Options Database Keys:
1327a7e14dcfSSatish Balay +    -tao_monitor        - sets TaoDefaultMonitor()
1328a7e14dcfSSatish Balay .    -tao_smonitor       - sets short monitor
1329a7e14dcfSSatish Balay .    -tao_cmonitor       - same as smonitor plus constraint norm
1330a7e14dcfSSatish Balay .    -tao_view_solution   - view solution at each iteration
1331a7e14dcfSSatish Balay .    -tao_view_gradient   - view gradient at each iteration
1332a7e14dcfSSatish Balay .    -tao_view_separableobjective - view separable objective function at each iteration
1333a7e14dcfSSatish 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.
1334a7e14dcfSSatish Balay 
1335a7e14dcfSSatish Balay 
1336a7e14dcfSSatish Balay    Notes:
1337a7e14dcfSSatish Balay    Several different monitoring routines may be set by calling
1338a7e14dcfSSatish Balay    TaoSetMonitor() multiple times; all will be called in the
1339a7e14dcfSSatish Balay    order in which they were set.
1340a7e14dcfSSatish Balay 
1341a7e14dcfSSatish Balay    Fortran Notes: Only one monitor function may be set
1342a7e14dcfSSatish Balay 
1343a7e14dcfSSatish Balay    Level: intermediate
1344a7e14dcfSSatish Balay 
1345a7e14dcfSSatish Balay .seealso: TaoDefaultMonitor(), TaoCancelMonitors(),  TaoSetDestroyRoutine()
1346a7e14dcfSSatish Balay @*/
1347441846f8SBarry Smith PetscErrorCode TaoSetMonitor(Tao tao, PetscErrorCode (*func)(Tao, void*), void *ctx,PetscErrorCode (*dest)(void**))
1348a7e14dcfSSatish Balay {
1349a7e14dcfSSatish Balay   PetscFunctionBegin;
1350441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
135147a47007SBarry Smith   if (tao->numbermonitors >= MAXTAOMONITORS) SETERRQ1(PETSC_COMM_SELF,1,"Cannot attach another monitor -- max=",MAXTAOMONITORS);
1352a7e14dcfSSatish Balay   tao->monitor[tao->numbermonitors] = func;
1353a7e14dcfSSatish Balay   tao->monitorcontext[tao->numbermonitors] = ctx;
1354a7e14dcfSSatish Balay   tao->monitordestroy[tao->numbermonitors] = dest;
1355a7e14dcfSSatish Balay   ++tao->numbermonitors;
1356a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1357a7e14dcfSSatish Balay }
1358a7e14dcfSSatish Balay 
1359a7e14dcfSSatish Balay #undef __FUNCT__
1360a7e14dcfSSatish Balay #define __FUNCT__ "TaoCancelMonitors"
1361a7e14dcfSSatish Balay /*@
1362441846f8SBarry Smith    TaoCancelMonitors - Clears all the monitor functions for a Tao object.
1363a7e14dcfSSatish Balay 
1364441846f8SBarry Smith    Logically Collective on Tao
1365a7e14dcfSSatish Balay 
1366a7e14dcfSSatish Balay    Input Parameters:
1367441846f8SBarry Smith .  tao - the Tao solver context
1368a7e14dcfSSatish Balay 
1369a7e14dcfSSatish Balay    Options Database:
1370a7e14dcfSSatish Balay .  -tao_cancelmonitors - cancels all monitors that have been hardwired
1371a7e14dcfSSatish Balay     into a code by calls to TaoSetMonitor(), but does not cancel those
1372a7e14dcfSSatish Balay     set via the options database
1373a7e14dcfSSatish Balay 
1374a7e14dcfSSatish Balay    Notes:
1375441846f8SBarry Smith    There is no way to clear one specific monitor from a Tao object.
1376a7e14dcfSSatish Balay 
1377a7e14dcfSSatish Balay    Level: advanced
1378a7e14dcfSSatish Balay 
1379a7e14dcfSSatish Balay .seealso: TaoDefaultMonitor(), TaoSetMonitor()
1380a7e14dcfSSatish Balay @*/
1381441846f8SBarry Smith PetscErrorCode TaoCancelMonitors(Tao tao)
1382a7e14dcfSSatish Balay {
1383a7e14dcfSSatish Balay   PetscInt       i;
1384a7e14dcfSSatish Balay   PetscErrorCode ierr;
138547a47007SBarry Smith 
1386a7e14dcfSSatish Balay   PetscFunctionBegin;
1387441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
1388a7e14dcfSSatish Balay   for (i=0;i<tao->numbermonitors;i++) {
1389a7e14dcfSSatish Balay     if (tao->monitordestroy[i]) {
1390a7e14dcfSSatish Balay       ierr = (*tao->monitordestroy[i])(&tao->monitorcontext[i]);CHKERRQ(ierr);
1391a7e14dcfSSatish Balay     }
1392a7e14dcfSSatish Balay   }
1393a7e14dcfSSatish Balay   tao->numbermonitors=0;
1394a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1395a7e14dcfSSatish Balay }
1396a7e14dcfSSatish Balay 
1397a7e14dcfSSatish Balay #undef __FUNCT__
1398a7e14dcfSSatish Balay #define __FUNCT__ "TaoDefaultMonitor"
1399a7e14dcfSSatish Balay /*@C
1400a7e14dcfSSatish Balay    TaoDefaultMonitor - Default routine for monitoring progress of the
1401441846f8SBarry Smith    Tao solvers (default).  This monitor prints the function value and gradient
1402a7e14dcfSSatish Balay    norm at each iteration.  It can be turned on from the command line using the
1403a7e14dcfSSatish Balay    -tao_monitor option
1404a7e14dcfSSatish Balay 
1405441846f8SBarry Smith    Collective on Tao
1406a7e14dcfSSatish Balay 
1407a7e14dcfSSatish Balay    Input Parameters:
1408441846f8SBarry Smith +  tao - the Tao context
140947a47007SBarry Smith -  ctx - PetscViewer context or NULL
1410a7e14dcfSSatish Balay 
1411a7e14dcfSSatish Balay    Options Database Keys:
1412a7e14dcfSSatish Balay .  -tao_monitor
1413a7e14dcfSSatish Balay 
1414a7e14dcfSSatish Balay    Level: advanced
1415a7e14dcfSSatish Balay 
1416a7e14dcfSSatish Balay .seealso: TaoDefaultSMonitor(), TaoSetMonitor()
1417a7e14dcfSSatish Balay @*/
1418441846f8SBarry Smith PetscErrorCode TaoDefaultMonitor(Tao tao, void *ctx)
1419a7e14dcfSSatish Balay {
1420a7e14dcfSSatish Balay   PetscErrorCode ierr;
1421a7e14dcfSSatish Balay   PetscInt       its;
1422a7e14dcfSSatish Balay   PetscReal      fct,gnorm;
1423a7e14dcfSSatish Balay   PetscViewer    viewer;
142447a47007SBarry Smith 
1425a7e14dcfSSatish Balay   PetscFunctionBegin;
1426a7e14dcfSSatish Balay   if (ctx) {
1427a7e14dcfSSatish Balay     viewer = (PetscViewer)ctx;
1428a7e14dcfSSatish Balay   } else {
1429a7e14dcfSSatish Balay     viewer = PETSC_VIEWER_STDOUT_(((PetscObject)tao)->comm);
1430a7e14dcfSSatish Balay   }
1431a7e14dcfSSatish Balay   its=tao->niter;
1432a7e14dcfSSatish Balay   fct=tao->fc;
1433a7e14dcfSSatish Balay   gnorm=tao->residual;
1434a7e14dcfSSatish Balay   ierr=PetscViewerASCIIPrintf(viewer,"iter = %3D,",its);CHKERRQ(ierr);
143547a47007SBarry Smith   ierr=PetscViewerASCIIPrintf(viewer," Function value: %g,",(double)fct);CHKERRQ(ierr);
143647a47007SBarry Smith   ierr=PetscViewerASCIIPrintf(viewer,"  Residual: %g \n",(double)gnorm);CHKERRQ(ierr);
1437a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1438a7e14dcfSSatish Balay }
1439a7e14dcfSSatish Balay 
1440a7e14dcfSSatish Balay #undef __FUNCT__
1441a7e14dcfSSatish Balay #define __FUNCT__ "TaoDefaultSMonitor"
1442a7e14dcfSSatish Balay /*@C
1443a7e14dcfSSatish Balay    TaoDefaultSMonitor - Default routine for monitoring progress of the
1444a7e14dcfSSatish Balay    solver. Same as TaoDefaultMonitor() except
1445a7e14dcfSSatish Balay    it prints fewer digits of the residual as the residual gets smaller.
1446a7e14dcfSSatish Balay    This is because the later digits are meaningless and are often
1447a7e14dcfSSatish Balay    different on different machines; by using this routine different
1448a7e14dcfSSatish Balay    machines will usually generate the same output. It can be turned on
1449a7e14dcfSSatish Balay    by using the -tao_smonitor option
1450a7e14dcfSSatish Balay 
1451441846f8SBarry Smith    Collective on Tao
1452a7e14dcfSSatish Balay 
1453a7e14dcfSSatish Balay    Input Parameters:
1454441846f8SBarry Smith +  tao - the Tao context
145547a47007SBarry Smith -  ctx - PetscViewer context or NULL
1456a7e14dcfSSatish Balay 
1457a7e14dcfSSatish Balay    Options Database Keys:
1458a7e14dcfSSatish Balay .  -tao_smonitor
1459a7e14dcfSSatish Balay 
1460a7e14dcfSSatish Balay    Level: advanced
1461a7e14dcfSSatish Balay 
1462a7e14dcfSSatish Balay .seealso: TaoDefaultMonitor(), TaoSetMonitor()
1463a7e14dcfSSatish Balay @*/
1464441846f8SBarry Smith PetscErrorCode TaoDefaultSMonitor(Tao tao, void *ctx)
1465a7e14dcfSSatish Balay {
1466a7e14dcfSSatish Balay   PetscErrorCode ierr;
1467a7e14dcfSSatish Balay   PetscInt       its;
1468a7e14dcfSSatish Balay   PetscReal      fct,gnorm;
1469a7e14dcfSSatish Balay   PetscViewer    viewer;
147047a47007SBarry Smith 
1471a7e14dcfSSatish Balay   PetscFunctionBegin;
1472a7e14dcfSSatish Balay   if (ctx) {
1473a7e14dcfSSatish Balay     viewer = (PetscViewer)ctx;
1474a7e14dcfSSatish Balay   } else {
1475a7e14dcfSSatish Balay     viewer = PETSC_VIEWER_STDOUT_(((PetscObject)tao)->comm);
1476a7e14dcfSSatish Balay   }
1477a7e14dcfSSatish Balay   its=tao->niter;
1478a7e14dcfSSatish Balay   fct=tao->fc;
1479a7e14dcfSSatish Balay   gnorm=tao->residual;
1480a7e14dcfSSatish Balay   ierr=PetscViewerASCIIPrintf(viewer,"iter = %3D,",its);CHKERRQ(ierr);
148147a47007SBarry Smith   ierr=PetscViewerASCIIPrintf(viewer," Function value %g,",(double)fct);CHKERRQ(ierr);
1482a7e14dcfSSatish Balay   if (gnorm > 1.e-6) {
148347a47007SBarry Smith     ierr=PetscViewerASCIIPrintf(viewer," Residual: %g \n",(double)gnorm);CHKERRQ(ierr);
1484a7e14dcfSSatish Balay   } else if (gnorm > 1.e-11) {
1485a7e14dcfSSatish Balay     ierr=PetscViewerASCIIPrintf(viewer," Residual: < 1.0e-6 \n");CHKERRQ(ierr);
1486a7e14dcfSSatish Balay   } else {
1487a7e14dcfSSatish Balay     ierr=PetscViewerASCIIPrintf(viewer," Residual: < 1.0e-11 \n");CHKERRQ(ierr);
1488a7e14dcfSSatish Balay   }
1489a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1490a7e14dcfSSatish Balay }
1491a7e14dcfSSatish Balay 
1492a7e14dcfSSatish Balay #undef __FUNCT__
1493a7e14dcfSSatish Balay #define __FUNCT__ "TaoDefaultCMonitor"
1494a7e14dcfSSatish Balay /*@C
1495a7e14dcfSSatish Balay    TaoDefaultCMonitor - same as TaoDefaultMonitor() except
1496a7e14dcfSSatish Balay    it prints the norm of the constraints function. It can be turned on
1497a7e14dcfSSatish Balay    from the command line using the -tao_cmonitor option
1498a7e14dcfSSatish Balay 
1499441846f8SBarry Smith    Collective on Tao
1500a7e14dcfSSatish Balay 
1501a7e14dcfSSatish Balay    Input Parameters:
1502441846f8SBarry Smith +  tao - the Tao context
150347a47007SBarry Smith -  ctx - PetscViewer context or NULL
1504a7e14dcfSSatish Balay 
1505a7e14dcfSSatish Balay    Options Database Keys:
1506a7e14dcfSSatish Balay .  -tao_cmonitor
1507a7e14dcfSSatish Balay 
1508a7e14dcfSSatish Balay    Level: advanced
1509a7e14dcfSSatish Balay 
1510a7e14dcfSSatish Balay .seealso: TaoDefaultMonitor(), TaoSetMonitor()
1511a7e14dcfSSatish Balay @*/
1512441846f8SBarry Smith PetscErrorCode TaoDefaultCMonitor(Tao tao, void *ctx)
1513a7e14dcfSSatish Balay {
1514a7e14dcfSSatish Balay   PetscErrorCode ierr;
1515a7e14dcfSSatish Balay   PetscInt       its;
1516a7e14dcfSSatish Balay   PetscReal      fct,gnorm;
1517a7e14dcfSSatish Balay   PetscViewer    viewer;
1518a7e14dcfSSatish Balay 
1519a7e14dcfSSatish Balay   PetscFunctionBegin;
1520a7e14dcfSSatish Balay   if (ctx) {
1521a7e14dcfSSatish Balay     viewer = (PetscViewer)ctx;
1522a7e14dcfSSatish Balay   } else {
1523a7e14dcfSSatish Balay     viewer = PETSC_VIEWER_STDOUT_(((PetscObject)tao)->comm);
1524a7e14dcfSSatish Balay   }
1525a7e14dcfSSatish Balay   its=tao->niter;
1526a7e14dcfSSatish Balay   fct=tao->fc;
1527a7e14dcfSSatish Balay   gnorm=tao->residual;
1528a7e14dcfSSatish Balay   ierr=PetscViewerASCIIPrintf(viewer,"iter = %D,",its);CHKERRQ(ierr);
152947a47007SBarry Smith   ierr=PetscViewerASCIIPrintf(viewer," Function value: %g,",(double)fct);CHKERRQ(ierr);
153047a47007SBarry Smith   ierr=PetscViewerASCIIPrintf(viewer,"  Residual: %g ",(double)gnorm);CHKERRQ(ierr);
153147a47007SBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"  Constraint: %g \n",(double)tao->cnorm);CHKERRQ(ierr);
1532a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1533a7e14dcfSSatish Balay }
1534a7e14dcfSSatish Balay 
1535a7e14dcfSSatish Balay #undef __FUNCT__
1536a7e14dcfSSatish Balay #define __FUNCT__ "TaoSolutionMonitor"
1537a7e14dcfSSatish Balay /*@C
1538a7e14dcfSSatish Balay    TaoSolutionMonitor - Views the solution at each iteration
1539a7e14dcfSSatish Balay    It can be turned on from the command line using the
1540a7e14dcfSSatish Balay    -tao_view_solution option
1541a7e14dcfSSatish Balay 
1542441846f8SBarry Smith    Collective on Tao
1543a7e14dcfSSatish Balay 
1544a7e14dcfSSatish Balay    Input Parameters:
1545441846f8SBarry Smith +  tao - the Tao context
154647a47007SBarry Smith -  ctx - PetscViewer context or NULL
1547a7e14dcfSSatish Balay 
1548a7e14dcfSSatish Balay    Options Database Keys:
1549a7e14dcfSSatish Balay .  -tao_view_solution
1550a7e14dcfSSatish Balay 
1551a7e14dcfSSatish Balay    Level: advanced
1552a7e14dcfSSatish Balay 
1553a7e14dcfSSatish Balay .seealso: TaoDefaultSMonitor(), TaoSetMonitor()
1554a7e14dcfSSatish Balay @*/
1555441846f8SBarry Smith PetscErrorCode TaoSolutionMonitor(Tao tao, void *ctx)
1556a7e14dcfSSatish Balay {
1557a7e14dcfSSatish Balay   PetscErrorCode ierr;
1558a7e14dcfSSatish Balay   PetscViewer viewer;
155947a47007SBarry Smith 
1560a7e14dcfSSatish Balay   PetscFunctionBegin;
1561a7e14dcfSSatish Balay   if (ctx) {
1562a7e14dcfSSatish Balay     viewer = (PetscViewer)ctx;
1563a7e14dcfSSatish Balay   } else {
1564a7e14dcfSSatish Balay     viewer = PETSC_VIEWER_STDOUT_(((PetscObject)tao)->comm);
1565a7e14dcfSSatish Balay   }
1566a7e14dcfSSatish Balay   ierr = VecView(tao->solution, viewer);CHKERRQ(ierr);
1567a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1568a7e14dcfSSatish Balay }
1569a7e14dcfSSatish Balay 
1570a7e14dcfSSatish Balay #undef __FUNCT__
1571a7e14dcfSSatish Balay #define __FUNCT__ "TaoGradientMonitor"
1572a7e14dcfSSatish Balay /*@C
1573a7e14dcfSSatish Balay    TaoGradientMonitor - Views the gradient at each iteration
1574a7e14dcfSSatish Balay    It can be turned on from the command line using the
1575a7e14dcfSSatish Balay    -tao_view_gradient option
1576a7e14dcfSSatish Balay 
1577441846f8SBarry Smith    Collective on Tao
1578a7e14dcfSSatish Balay 
1579a7e14dcfSSatish Balay    Input Parameters:
1580441846f8SBarry Smith +  tao - the Tao context
158147a47007SBarry Smith -  ctx - PetscViewer context or NULL
1582a7e14dcfSSatish Balay 
1583a7e14dcfSSatish Balay    Options Database Keys:
1584a7e14dcfSSatish Balay .  -tao_view_gradient
1585a7e14dcfSSatish Balay 
1586a7e14dcfSSatish Balay    Level: advanced
1587a7e14dcfSSatish Balay 
1588a7e14dcfSSatish Balay .seealso: TaoDefaultSMonitor(), TaoSetMonitor()
1589a7e14dcfSSatish Balay @*/
1590441846f8SBarry Smith PetscErrorCode TaoGradientMonitor(Tao tao, void *ctx)
1591a7e14dcfSSatish Balay {
1592a7e14dcfSSatish Balay   PetscErrorCode ierr;
1593a7e14dcfSSatish Balay   PetscViewer viewer;
159447a47007SBarry Smith 
1595a7e14dcfSSatish Balay   PetscFunctionBegin;
1596a7e14dcfSSatish Balay   if (ctx) {
1597a7e14dcfSSatish Balay     viewer = (PetscViewer)ctx;
1598a7e14dcfSSatish Balay   } else {
1599a7e14dcfSSatish Balay     viewer = PETSC_VIEWER_STDOUT_(((PetscObject)tao)->comm);
1600a7e14dcfSSatish Balay   }
1601a7e14dcfSSatish Balay   ierr = VecView(tao->gradient, viewer);CHKERRQ(ierr);
1602a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1603a7e14dcfSSatish Balay }
1604a7e14dcfSSatish Balay 
1605a7e14dcfSSatish Balay #undef __FUNCT__
1606a7e14dcfSSatish Balay #define __FUNCT__ "TaoStepDirectionMonitor"
1607a7e14dcfSSatish Balay /*@C
1608a7e14dcfSSatish Balay    TaoStepDirectionMonitor - Views the gradient at each iteration
1609a7e14dcfSSatish Balay    It can be turned on from the command line using the
1610a7e14dcfSSatish Balay    -tao_view_gradient 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_gradient
1620a7e14dcfSSatish Balay 
1621a7e14dcfSSatish Balay    Level: advanced
1622a7e14dcfSSatish Balay 
1623a7e14dcfSSatish Balay .seealso: TaoDefaultSMonitor(), TaoSetMonitor()
1624a7e14dcfSSatish Balay @*/
1625441846f8SBarry Smith PetscErrorCode TaoStepDirectionMonitor(Tao tao, void *ctx)
1626a7e14dcfSSatish Balay {
1627a7e14dcfSSatish Balay   PetscErrorCode ierr;
1628a7e14dcfSSatish Balay   PetscViewer viewer;
1629a7e14dcfSSatish Balay   PetscFunctionBegin;
1630a7e14dcfSSatish Balay   if (ctx) {
1631a7e14dcfSSatish Balay     viewer = (PetscViewer)ctx;
1632a7e14dcfSSatish Balay   } else {
1633a7e14dcfSSatish Balay     viewer = PETSC_VIEWER_STDOUT_(((PetscObject)tao)->comm);
1634a7e14dcfSSatish Balay   }
1635a7e14dcfSSatish Balay   ierr = VecView(tao->stepdirection, viewer);CHKERRQ(ierr);
1636a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1637a7e14dcfSSatish Balay }
1638a7e14dcfSSatish Balay 
1639a7e14dcfSSatish Balay #undef __FUNCT__
1640a7e14dcfSSatish Balay #define __FUNCT__ "TaoDrawSolutionMonitor"
1641a7e14dcfSSatish Balay /*@C
1642a7e14dcfSSatish Balay    TaoDrawSolutionMonitor - Plots the solution at each iteration
1643a7e14dcfSSatish Balay    It can be turned on from the command line using the
1644a7e14dcfSSatish Balay    -tao_draw_solution option
1645a7e14dcfSSatish Balay 
1646441846f8SBarry Smith    Collective on Tao
1647a7e14dcfSSatish Balay 
1648a7e14dcfSSatish Balay    Input Parameters:
1649441846f8SBarry Smith +  tao - the Tao context
165047a47007SBarry Smith -  ctx - PetscViewer context or NULL
1651a7e14dcfSSatish Balay 
1652a7e14dcfSSatish Balay    Options Database Keys:
1653a7e14dcfSSatish Balay .  -tao_draw_solution
1654a7e14dcfSSatish Balay 
1655a7e14dcfSSatish Balay    Level: advanced
1656a7e14dcfSSatish Balay 
1657a7e14dcfSSatish Balay .seealso: TaoSolutionMonitor(), TaoSetMonitor(), TaoDrawGradientMonitor
1658a7e14dcfSSatish Balay @*/
1659441846f8SBarry Smith PetscErrorCode TaoDrawSolutionMonitor(Tao tao, void *ctx)
1660a7e14dcfSSatish Balay {
1661a7e14dcfSSatish Balay   PetscErrorCode ierr;
1662a7e14dcfSSatish Balay   PetscViewer    viewer = (PetscViewer) ctx;
1663a7e14dcfSSatish Balay   MPI_Comm       comm;
166447a47007SBarry Smith 
1665a7e14dcfSSatish Balay   PetscFunctionBegin;
1666a7e14dcfSSatish Balay   if (!viewer) {
1667a7e14dcfSSatish Balay     ierr = PetscObjectGetComm((PetscObject)tao,&comm);CHKERRQ(ierr);
1668a7e14dcfSSatish Balay     viewer = PETSC_VIEWER_DRAW_(comm);
1669a7e14dcfSSatish Balay   }
1670a7e14dcfSSatish Balay   ierr = VecView(tao->solution, viewer);CHKERRQ(ierr);
1671a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1672a7e14dcfSSatish Balay }
1673a7e14dcfSSatish Balay 
1674a7e14dcfSSatish Balay #undef __FUNCT__
1675a7e14dcfSSatish Balay #define __FUNCT__ "TaoDrawGradientMonitor"
1676a7e14dcfSSatish Balay /*@C
1677a7e14dcfSSatish Balay    TaoDrawGradientMonitor - Plots the gradient at each iteration
1678a7e14dcfSSatish Balay    It can be turned on from the command line using the
1679a7e14dcfSSatish Balay    -tao_draw_gradient option
1680a7e14dcfSSatish Balay 
1681441846f8SBarry Smith    Collective on Tao
1682a7e14dcfSSatish Balay 
1683a7e14dcfSSatish Balay    Input Parameters:
1684441846f8SBarry Smith +  tao - the Tao context
168547a47007SBarry Smith -  ctx - PetscViewer context or NULL
1686a7e14dcfSSatish Balay 
1687a7e14dcfSSatish Balay    Options Database Keys:
1688a7e14dcfSSatish Balay .  -tao_draw_gradient
1689a7e14dcfSSatish Balay 
1690a7e14dcfSSatish Balay    Level: advanced
1691a7e14dcfSSatish Balay 
1692a7e14dcfSSatish Balay .seealso: TaoGradientMonitor(), TaoSetMonitor(), TaoDrawSolutionMonitor
1693a7e14dcfSSatish Balay @*/
1694441846f8SBarry Smith PetscErrorCode TaoDrawGradientMonitor(Tao tao, void *ctx)
1695a7e14dcfSSatish Balay {
1696a7e14dcfSSatish Balay   PetscErrorCode ierr;
1697a7e14dcfSSatish Balay   PetscViewer    viewer = (PetscViewer)ctx;
1698a7e14dcfSSatish Balay   MPI_Comm       comm;
169947a47007SBarry Smith 
1700a7e14dcfSSatish Balay   PetscFunctionBegin;
1701a7e14dcfSSatish Balay   if (!viewer) {
1702a7e14dcfSSatish Balay     ierr = PetscObjectGetComm((PetscObject)tao,&comm);CHKERRQ(ierr);
1703a7e14dcfSSatish Balay     viewer = PETSC_VIEWER_DRAW_(comm);
1704a7e14dcfSSatish Balay   }
1705a7e14dcfSSatish Balay   ierr = VecView(tao->gradient, viewer);CHKERRQ(ierr);
1706a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1707a7e14dcfSSatish Balay }
1708a7e14dcfSSatish Balay 
1709a7e14dcfSSatish Balay #undef __FUNCT__
1710a7e14dcfSSatish Balay #define __FUNCT__ "TaoDrawStepMonitor"
1711a7e14dcfSSatish Balay /*@C
1712a7e14dcfSSatish Balay    TaoDrawStepMonitor - Plots the step direction at each iteration
1713a7e14dcfSSatish Balay    It can be turned on from the command line using the
1714a7e14dcfSSatish Balay    -tao_draw_step option
1715a7e14dcfSSatish Balay 
1716441846f8SBarry Smith    Collective on Tao
1717a7e14dcfSSatish Balay 
1718a7e14dcfSSatish Balay    Input Parameters:
1719441846f8SBarry Smith +  tao - the Tao context
172047a47007SBarry Smith -  ctx - PetscViewer context or NULL
1721a7e14dcfSSatish Balay 
1722a7e14dcfSSatish Balay    Options Database Keys:
1723a7e14dcfSSatish Balay .  -tao_draw_step
1724a7e14dcfSSatish Balay 
1725a7e14dcfSSatish Balay    Level: advanced
1726a7e14dcfSSatish Balay 
1727a7e14dcfSSatish Balay .seealso: TaoSetMonitor(), TaoDrawSolutionMonitor
1728a7e14dcfSSatish Balay @*/
1729441846f8SBarry Smith PetscErrorCode TaoDrawStepMonitor(Tao tao, void *ctx)
1730a7e14dcfSSatish Balay {
1731a7e14dcfSSatish Balay   PetscErrorCode ierr;
1732a7e14dcfSSatish Balay   PetscViewer    viewer = (PetscViewer)(ctx);
1733a7e14dcfSSatish Balay   MPI_Comm       comm;
173447a47007SBarry Smith 
1735a7e14dcfSSatish Balay   PetscFunctionBegin;
1736a7e14dcfSSatish Balay   if (!viewer) {
1737a7e14dcfSSatish Balay     ierr = PetscObjectGetComm((PetscObject)tao,&comm);CHKERRQ(ierr);
1738a7e14dcfSSatish Balay     viewer = PETSC_VIEWER_DRAW_(comm);
1739a7e14dcfSSatish Balay   }
1740a7e14dcfSSatish Balay   ierr = VecView(tao->stepdirection, viewer);CHKERRQ(ierr);
1741a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1742a7e14dcfSSatish Balay }
1743a7e14dcfSSatish Balay 
1744a7e14dcfSSatish Balay #undef __FUNCT__
1745a7e14dcfSSatish Balay #define __FUNCT__ "TaoSeparableObjectiveMonitor"
1746a7e14dcfSSatish Balay /*@C
1747a7e14dcfSSatish Balay    TaoSeparableObjectiveMonitor - Views the separable objective function at each iteration
1748a7e14dcfSSatish Balay    It can be turned on from the command line using the
1749a7e14dcfSSatish Balay    -tao_view_separableobjective option
1750a7e14dcfSSatish Balay 
1751441846f8SBarry Smith    Collective on Tao
1752a7e14dcfSSatish Balay 
1753a7e14dcfSSatish Balay    Input Parameters:
1754441846f8SBarry Smith +  tao - the Tao context
175547a47007SBarry Smith -  ctx - PetscViewer context or NULL
1756a7e14dcfSSatish Balay 
1757a7e14dcfSSatish Balay    Options Database Keys:
1758a7e14dcfSSatish Balay .  -tao_view_separableobjective
1759a7e14dcfSSatish Balay 
1760a7e14dcfSSatish Balay    Level: advanced
1761a7e14dcfSSatish Balay 
1762a7e14dcfSSatish Balay .seealso: TaoDefaultSMonitor(), TaoSetMonitor()
1763a7e14dcfSSatish Balay @*/
1764441846f8SBarry Smith PetscErrorCode TaoSeparableObjectiveMonitor(Tao tao, void *ctx)
1765a7e14dcfSSatish Balay {
1766a7e14dcfSSatish Balay   PetscErrorCode ierr;
1767a7e14dcfSSatish Balay   PetscViewer    viewer;
176847a47007SBarry Smith 
1769a7e14dcfSSatish Balay   PetscFunctionBegin;
1770a7e14dcfSSatish Balay   if (ctx) {
1771a7e14dcfSSatish Balay     viewer = (PetscViewer)ctx;
1772a7e14dcfSSatish Balay   } else {
1773a7e14dcfSSatish Balay     viewer = PETSC_VIEWER_STDOUT_(((PetscObject)tao)->comm);
1774a7e14dcfSSatish Balay   }
1775a7e14dcfSSatish Balay   ierr = VecView(tao->sep_objective,viewer);CHKERRQ(ierr);
1776a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1777a7e14dcfSSatish Balay }
1778a7e14dcfSSatish Balay 
1779a7e14dcfSSatish Balay #undef __FUNCT__
1780a7e14dcfSSatish Balay #define __FUNCT__ "TaoDefaultConvergenceTest"
1781a7e14dcfSSatish Balay /*@C
1782a7e14dcfSSatish Balay    TaoDefaultConvergenceTest - Determines whether the solver should continue iterating
1783a7e14dcfSSatish Balay    or terminate.
1784a7e14dcfSSatish Balay 
1785441846f8SBarry Smith    Collective on Tao
1786a7e14dcfSSatish Balay 
1787a7e14dcfSSatish Balay    Input Parameters:
1788441846f8SBarry Smith +  tao - the Tao context
1789a7e14dcfSSatish Balay -  dummy - unused dummy context
1790a7e14dcfSSatish Balay 
1791a7e14dcfSSatish Balay    Output Parameter:
1792a7e14dcfSSatish Balay .  reason - for terminating
1793a7e14dcfSSatish Balay 
1794a7e14dcfSSatish Balay    Notes:
1795a7e14dcfSSatish Balay    This routine checks the residual in the optimality conditions, the
1796a7e14dcfSSatish Balay    relative residual in the optimity conditions, the number of function
1797a7e14dcfSSatish Balay    evaluations, and the function value to test convergence.  Some
1798a7e14dcfSSatish Balay    solvers may use different convergence routines.
1799a7e14dcfSSatish Balay 
1800a7e14dcfSSatish Balay    Level: developer
1801a7e14dcfSSatish Balay 
1802*e4cb33bbSBarry Smith .seealso: TaoSetTolerances(),TaoGetConvergedReason(),TaoSetConvergedReason()
1803a7e14dcfSSatish Balay @*/
1804a7e14dcfSSatish Balay 
1805441846f8SBarry Smith PetscErrorCode TaoDefaultConvergenceTest(Tao tao,void *dummy)
1806a7e14dcfSSatish Balay {
1807a7e14dcfSSatish Balay   PetscInt           niter=tao->niter, nfuncs=PetscMax(tao->nfuncs,tao->nfuncgrads);
1808a7e14dcfSSatish Balay   PetscInt           max_funcs=tao->max_funcs;
1809a7e14dcfSSatish Balay   PetscReal          gnorm=tao->residual, gnorm0=tao->gnorm0;
1810a7e14dcfSSatish Balay   PetscReal          f=tao->fc, steptol=tao->steptol,trradius=tao->step;
1811a7e14dcfSSatish Balay   PetscReal          gatol=tao->gatol,grtol=tao->grtol,gttol=tao->gttol;
1812a7e14dcfSSatish Balay   PetscReal          fatol=tao->fatol,frtol=tao->frtol,catol=tao->catol,crtol=tao->crtol;
1813a7e14dcfSSatish Balay   PetscReal          fmin=tao->fmin, cnorm=tao->cnorm, cnorm0=tao->cnorm0;
1814a7e14dcfSSatish Balay   PetscReal          gnorm2;
1815*e4cb33bbSBarry Smith   TaoConvergedReason reason=tao->reason;
1816a7e14dcfSSatish Balay   PetscErrorCode     ierr;
1817a7e14dcfSSatish Balay 
1818a7e14dcfSSatish Balay   PetscFunctionBegin;
1819441846f8SBarry Smith   PetscValidHeaderSpecific(tao, TAO_CLASSID,1);
1820a7e14dcfSSatish Balay   if (reason != TAO_CONTINUE_ITERATING) {
1821a7e14dcfSSatish Balay     PetscFunctionReturn(0);
1822a7e14dcfSSatish Balay   }
1823a7e14dcfSSatish Balay   gnorm2=gnorm*gnorm;
1824a7e14dcfSSatish Balay 
1825a7e14dcfSSatish Balay   if (PetscIsInfOrNanReal(f)) {
1826a7e14dcfSSatish Balay     ierr = PetscInfo(tao,"Failed to converged, function value is Inf or NaN\n");CHKERRQ(ierr);
1827a7e14dcfSSatish Balay     reason = TAO_DIVERGED_NAN;
1828a7e14dcfSSatish Balay   } else if (f <= fmin && cnorm <=catol) {
182947a47007SBarry Smith     ierr = PetscInfo2(tao,"Converged due to function value %g < minimum function value %g\n", (double)f,(double)fmin);CHKERRQ(ierr);
1830a7e14dcfSSatish Balay     reason = TAO_CONVERGED_MINF;
1831a7e14dcfSSatish Balay   } else if (gnorm2 <= fatol && cnorm <=catol) {
183247a47007SBarry Smith     ierr = PetscInfo2(tao,"Converged due to estimated f(X) - f(X*) = %g < %g\n",(double)gnorm2,(double)fatol);CHKERRQ(ierr);
1833a7e14dcfSSatish Balay     reason = TAO_CONVERGED_FATOL;
1834a7e14dcfSSatish Balay   } else if (f != 0 && gnorm2 / PetscAbsReal(f)<= frtol && cnorm/PetscMax(cnorm0,1.0) <= crtol) {
183547a47007SBarry Smith     ierr = PetscInfo2(tao,"Converged due to estimated |f(X)-f(X*)|/f(X) = %g < %g\n",(double)(gnorm2/PetscAbsReal(f)),(double)frtol);CHKERRQ(ierr);
1836a7e14dcfSSatish Balay     reason = TAO_CONVERGED_FRTOL;
1837a7e14dcfSSatish Balay   } else if (gnorm<= gatol && cnorm <=catol) {
183847a47007SBarry Smith     ierr = PetscInfo2(tao,"Converged due to residual norm ||g(X)||=%g < %g\n",(double)gnorm,(double)gatol);CHKERRQ(ierr);
1839a7e14dcfSSatish Balay     reason = TAO_CONVERGED_GATOL;
1840a7e14dcfSSatish Balay   } else if ( f!=0 && PetscAbsReal(gnorm/f) <= grtol && cnorm <= crtol) {
184147a47007SBarry Smith     ierr = PetscInfo2(tao,"Converged due to residual ||g(X)||/|f(X)| =%g < %g\n",(double)(gnorm/f),(double)grtol);CHKERRQ(ierr);
1842a7e14dcfSSatish Balay     reason = TAO_CONVERGED_GRTOL;
1843a7e14dcfSSatish Balay   } else if (gnorm0 != 0 && gnorm/gnorm0 <= gttol && cnorm <= crtol) {
184447a47007SBarry Smith     ierr = PetscInfo2(tao,"Converged due to relative residual norm ||g(X)||/||g(X0)|| = %g < %g\n",(double)(gnorm/gnorm0),(double)gttol);CHKERRQ(ierr);
1845a7e14dcfSSatish Balay     reason = TAO_CONVERGED_GTTOL;
1846a7e14dcfSSatish Balay   } else if (nfuncs > max_funcs){
1847a7e14dcfSSatish Balay     ierr = PetscInfo2(tao,"Exceeded maximum number of function evaluations: %D > %D\n", nfuncs,max_funcs);CHKERRQ(ierr);
1848a7e14dcfSSatish Balay     reason = TAO_DIVERGED_MAXFCN;
1849a7e14dcfSSatish Balay   } else if ( tao->lsflag != 0 ){
1850a7e14dcfSSatish Balay     ierr = PetscInfo(tao,"Tao Line Search failure.\n");CHKERRQ(ierr);
1851a7e14dcfSSatish Balay     reason = TAO_DIVERGED_LS_FAILURE;
1852a7e14dcfSSatish Balay   } else if (trradius < steptol && niter > 0){
185347a47007SBarry Smith     ierr = PetscInfo2(tao,"Trust region/step size too small: %g < %g\n", (double)trradius,(double)steptol);CHKERRQ(ierr);
1854a7e14dcfSSatish Balay     reason = TAO_CONVERGED_STEPTOL;
1855a7e14dcfSSatish Balay   } else if (niter > tao->max_it) {
1856a7e14dcfSSatish Balay     ierr = PetscInfo2(tao,"Exceeded maximum number of iterations: %D > %D\n",niter,tao->max_it);CHKERRQ(ierr);
1857a7e14dcfSSatish Balay     reason = TAO_DIVERGED_MAXITS;
1858a7e14dcfSSatish Balay   } else {
1859a7e14dcfSSatish Balay     reason = TAO_CONTINUE_ITERATING;
1860a7e14dcfSSatish Balay   }
1861a7e14dcfSSatish Balay   tao->reason = reason;
1862a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1863a7e14dcfSSatish Balay }
1864a7e14dcfSSatish Balay 
1865a7e14dcfSSatish Balay #undef __FUNCT__
1866a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetOptionsPrefix"
1867a7e14dcfSSatish Balay /*@C
1868a7e14dcfSSatish Balay    TaoSetOptionsPrefix - Sets the prefix used for searching for all
1869a7e14dcfSSatish Balay    TAO options in the database.
1870a7e14dcfSSatish Balay 
1871a7e14dcfSSatish Balay 
1872441846f8SBarry Smith    Logically Collective on Tao
1873a7e14dcfSSatish Balay 
1874a7e14dcfSSatish Balay    Input Parameters:
1875441846f8SBarry Smith +  tao - the Tao context
1876a7e14dcfSSatish Balay -  prefix - the prefix string to prepend to all TAO option requests
1877a7e14dcfSSatish Balay 
1878a7e14dcfSSatish Balay    Notes:
1879a7e14dcfSSatish Balay    A hyphen (-) must NOT be given at the beginning of the prefix name.
1880a7e14dcfSSatish Balay    The first character of all runtime options is AUTOMATICALLY the hyphen.
1881a7e14dcfSSatish Balay 
1882a7e14dcfSSatish Balay    For example, to distinguish between the runtime options for two
1883a7e14dcfSSatish Balay    different TAO solvers, one could call
1884a7e14dcfSSatish Balay .vb
1885a7e14dcfSSatish Balay       TaoSetOptionsPrefix(tao1,"sys1_")
1886a7e14dcfSSatish Balay       TaoSetOptionsPrefix(tao2,"sys2_")
1887a7e14dcfSSatish Balay .ve
1888a7e14dcfSSatish Balay 
1889a7e14dcfSSatish Balay    This would enable use of different options for each system, such as
1890a7e14dcfSSatish Balay .vb
1891a7e14dcfSSatish Balay       -sys1_tao_method blmvm -sys1_tao_gtol 1.e-3
1892a7e14dcfSSatish Balay       -sys2_tao_method lmvm  -sys2_tao_gtol 1.e-4
1893a7e14dcfSSatish Balay .ve
1894a7e14dcfSSatish Balay 
1895a7e14dcfSSatish Balay 
1896a7e14dcfSSatish Balay    Level: advanced
1897a7e14dcfSSatish Balay 
1898a7e14dcfSSatish Balay .seealso: TaoAppendOptionsPrefix(), TaoGetOptionsPrefix()
1899a7e14dcfSSatish Balay @*/
1900a7e14dcfSSatish Balay 
1901441846f8SBarry Smith PetscErrorCode TaoSetOptionsPrefix(Tao tao, const char p[])
1902a7e14dcfSSatish Balay {
1903a7e14dcfSSatish Balay   PetscErrorCode ierr;
190447a47007SBarry Smith 
1905a7e14dcfSSatish Balay   PetscFunctionBegin;
1906a7e14dcfSSatish Balay   ierr = PetscObjectSetOptionsPrefix((PetscObject)tao,p);CHKERRQ(ierr);
1907a7e14dcfSSatish Balay   if (tao->linesearch) {
1908a7e14dcfSSatish Balay     ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,p);CHKERRQ(ierr);
1909a7e14dcfSSatish Balay   }
1910a7e14dcfSSatish Balay   if (tao->ksp) {
1911a7e14dcfSSatish Balay     ierr = KSPSetOptionsPrefix(tao->ksp,p);CHKERRQ(ierr);
1912a7e14dcfSSatish Balay   }
1913a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1914a7e14dcfSSatish Balay }
1915a7e14dcfSSatish Balay 
1916a7e14dcfSSatish Balay #undef __FUNCT__
1917a7e14dcfSSatish Balay #define __FUNCT__ "TaoAppendOptionsPrefix"
1918a7e14dcfSSatish Balay /*@C
1919a7e14dcfSSatish Balay    TaoAppendOptionsPrefix - Appends to the prefix used for searching for all
1920a7e14dcfSSatish Balay    TAO options in the database.
1921a7e14dcfSSatish Balay 
1922a7e14dcfSSatish Balay 
1923441846f8SBarry Smith    Logically Collective on Tao
1924a7e14dcfSSatish Balay 
1925a7e14dcfSSatish Balay    Input Parameters:
1926441846f8SBarry Smith +  tao - the Tao solver context
1927a7e14dcfSSatish Balay -  prefix - the prefix string to prepend to all TAO option requests
1928a7e14dcfSSatish Balay 
1929a7e14dcfSSatish Balay    Notes:
1930a7e14dcfSSatish Balay    A hyphen (-) must NOT be given at the beginning of the prefix name.
1931a7e14dcfSSatish Balay    The first character of all runtime options is AUTOMATICALLY the hyphen.
1932a7e14dcfSSatish Balay 
1933a7e14dcfSSatish Balay 
1934a7e14dcfSSatish Balay    Level: advanced
1935a7e14dcfSSatish Balay 
1936a7e14dcfSSatish Balay .seealso: TaoSetOptionsPrefix(), TaoGetOptionsPrefix()
1937a7e14dcfSSatish Balay @*/
1938441846f8SBarry Smith PetscErrorCode TaoAppendOptionsPrefix(Tao tao, const char p[])
1939a7e14dcfSSatish Balay {
1940a7e14dcfSSatish Balay   PetscErrorCode ierr;
194147a47007SBarry Smith 
1942a7e14dcfSSatish Balay   PetscFunctionBegin;
1943a7e14dcfSSatish Balay   ierr = PetscObjectAppendOptionsPrefix((PetscObject)tao,p);CHKERRQ(ierr);
1944a7e14dcfSSatish Balay   if (tao->linesearch) {
1945a7e14dcfSSatish Balay     ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,p);CHKERRQ(ierr);
1946a7e14dcfSSatish Balay   }
1947a7e14dcfSSatish Balay   if (tao->ksp) {
1948a7e14dcfSSatish Balay     ierr = KSPSetOptionsPrefix(tao->ksp,p);CHKERRQ(ierr);
1949a7e14dcfSSatish Balay   }
1950a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1951a7e14dcfSSatish Balay }
1952a7e14dcfSSatish Balay 
1953a7e14dcfSSatish Balay #undef __FUNCT__
1954a7e14dcfSSatish Balay #define __FUNCT__ "TaoGetOptionsPrefix"
1955a7e14dcfSSatish Balay /*@C
1956a7e14dcfSSatish Balay   TaoGetOptionsPrefix - Gets the prefix used for searching for all
1957a7e14dcfSSatish Balay   TAO options in the database
1958a7e14dcfSSatish Balay 
1959a7e14dcfSSatish Balay   Not Collective
1960a7e14dcfSSatish Balay 
1961a7e14dcfSSatish Balay   Input Parameters:
1962441846f8SBarry Smith . tao - the Tao context
1963a7e14dcfSSatish Balay 
1964a7e14dcfSSatish Balay   Output Parameters:
1965a7e14dcfSSatish Balay . prefix - pointer to the prefix string used is returned
1966a7e14dcfSSatish Balay 
1967a7e14dcfSSatish Balay   Notes: On the fortran side, the user should pass in a string 'prefix' of
1968a7e14dcfSSatish Balay   sufficient length to hold the prefix.
1969a7e14dcfSSatish Balay 
1970a7e14dcfSSatish Balay   Level: advanced
1971a7e14dcfSSatish Balay 
1972a7e14dcfSSatish Balay .seealso: TaoSetOptionsPrefix(), TaoAppendOptionsPrefix()
1973a7e14dcfSSatish Balay @*/
1974441846f8SBarry Smith PetscErrorCode TaoGetOptionsPrefix(Tao tao, const char *p[])
1975a7e14dcfSSatish Balay {
197647a47007SBarry Smith    return PetscObjectGetOptionsPrefix((PetscObject)tao,p);
1977a7e14dcfSSatish Balay }
1978a7e14dcfSSatish Balay 
1979a7e14dcfSSatish Balay #undef __FUNCT__
1980a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetType"
1981a7e14dcfSSatish Balay /*@C
1982a7e14dcfSSatish Balay    TaoSetType - Sets the method for the unconstrained minimization solver.
1983a7e14dcfSSatish Balay 
1984441846f8SBarry Smith    Collective on Tao
1985a7e14dcfSSatish Balay 
1986a7e14dcfSSatish Balay    Input Parameters:
1987441846f8SBarry Smith +  solver - the Tao solver context
1988a7e14dcfSSatish Balay -  type - a known method
1989a7e14dcfSSatish Balay 
1990a7e14dcfSSatish Balay    Options Database Key:
199158417fe7SBarry Smith .  -tao_type <type> - Sets the method; use -help for a list
199258417fe7SBarry Smith    of available methods (for instance, "-tao_type lmvm" or "-tao_type tron")
1993a7e14dcfSSatish Balay 
1994a7e14dcfSSatish Balay    Available methods include:
199558417fe7SBarry Smith +    nls - Newton's method with line search for unconstrained minimization
199658417fe7SBarry Smith .    ntr - Newton's method with trust region for unconstrained minimization
199758417fe7SBarry Smith .    ntl - Newton's method with trust region, line search for unconstrained minimization
199858417fe7SBarry Smith .    lmvm - Limited memory variable metric method for unconstrained minimization
199958417fe7SBarry Smith .    cg - Nonlinear conjugate gradient method for unconstrained minimization
200058417fe7SBarry Smith .    nm - Nelder-Mead algorithm for derivate-free unconstrained minimization
200158417fe7SBarry Smith .    tron - Newton Trust Region method for bound constrained minimization
200258417fe7SBarry Smith .    gpcg - Newton Trust Region method for quadratic bound constrained minimization
200358417fe7SBarry Smith .    blmvm - Limited memory variable metric method for bound constrained minimization
200458417fe7SBarry Smith +    pounders - Model-based algorithm pounder extended for nonlinear least squares
2005a7e14dcfSSatish Balay 
2006a7e14dcfSSatish Balay   Level: intermediate
2007a7e14dcfSSatish Balay 
2008a7e14dcfSSatish Balay .seealso: TaoCreate(), TaoGetType()
2009a7e14dcfSSatish Balay 
2010a7e14dcfSSatish Balay @*/
2011441846f8SBarry Smith PetscErrorCode TaoSetType(Tao tao, const TaoType type)
2012a7e14dcfSSatish Balay {
2013a7e14dcfSSatish Balay   PetscErrorCode ierr;
2014441846f8SBarry Smith   PetscErrorCode (*create_xxx)(Tao);
2015a7e14dcfSSatish Balay   PetscBool      issame;
2016a7e14dcfSSatish Balay 
2017a7e14dcfSSatish Balay   PetscFunctionBegin;
2018441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2019a7e14dcfSSatish Balay 
2020a7e14dcfSSatish Balay   ierr = PetscObjectTypeCompare((PetscObject)tao,type,&issame);CHKERRQ(ierr);
2021a7e14dcfSSatish Balay   if (issame) PetscFunctionReturn(0);
2022a7e14dcfSSatish Balay 
2023441846f8SBarry Smith   ierr = PetscFunctionListFind(TaoList, type, (void(**)(void))&create_xxx);CHKERRQ(ierr);
2024441846f8SBarry Smith   if (!create_xxx) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested Tao type %s",type);
2025a7e14dcfSSatish Balay 
2026a7e14dcfSSatish Balay   /* Destroy the existing solver information */
2027a7e14dcfSSatish Balay   if (tao->ops->destroy) {
2028a7e14dcfSSatish Balay     ierr = (*tao->ops->destroy)(tao);CHKERRQ(ierr);
2029a7e14dcfSSatish Balay   }
2030a7e14dcfSSatish Balay   ierr = KSPDestroy(&tao->ksp);CHKERRQ(ierr);
2031a7e14dcfSSatish Balay   ierr = TaoLineSearchDestroy(&tao->linesearch);CHKERRQ(ierr);
2032a7e14dcfSSatish Balay   ierr = VecDestroy(&tao->gradient);CHKERRQ(ierr);
2033a7e14dcfSSatish Balay   ierr = VecDestroy(&tao->stepdirection);CHKERRQ(ierr);
2034a7e14dcfSSatish Balay 
2035a7e14dcfSSatish Balay   tao->ops->setup = 0;
2036a7e14dcfSSatish Balay   tao->ops->solve = 0;
2037a7e14dcfSSatish Balay   tao->ops->view  = 0;
2038a7e14dcfSSatish Balay   tao->ops->setfromoptions = 0;
2039a7e14dcfSSatish Balay   tao->ops->destroy = 0;
2040a7e14dcfSSatish Balay 
2041a7e14dcfSSatish Balay   tao->setupcalled = PETSC_FALSE;
2042a7e14dcfSSatish Balay 
2043a7e14dcfSSatish Balay   ierr = (*create_xxx)(tao);CHKERRQ(ierr);
2044a7e14dcfSSatish Balay   ierr = PetscObjectChangeTypeName((PetscObject)tao,type);CHKERRQ(ierr);
2045a7e14dcfSSatish Balay   PetscFunctionReturn(0);
2046a7e14dcfSSatish Balay }
204747a47007SBarry Smith 
2048a7e14dcfSSatish Balay #undef __FUNCT__
2049441846f8SBarry Smith #define __FUNCT__ "TaoRegister"
2050a7e14dcfSSatish Balay /*MC
2051441846f8SBarry Smith    TaoRegister - Adds a method to the TAO package for unconstrained minimization.
2052a7e14dcfSSatish Balay 
2053a7e14dcfSSatish Balay    Synopsis:
2054441846f8SBarry Smith    TaoRegister(char *name_solver,char *path,char *name_Create,int (*routine_Create)(Tao))
2055a7e14dcfSSatish Balay 
2056a7e14dcfSSatish Balay    Not collective
2057a7e14dcfSSatish Balay 
2058a7e14dcfSSatish Balay    Input Parameters:
2059a7e14dcfSSatish Balay +  sname - name of a new user-defined solver
2060a7e14dcfSSatish Balay -  func - routine to Create method context
2061a7e14dcfSSatish Balay 
2062a7e14dcfSSatish Balay    Notes:
2063441846f8SBarry Smith    TaoRegister() may be called multiple times to add several user-defined solvers.
2064a7e14dcfSSatish Balay 
2065a7e14dcfSSatish Balay    Sample usage:
2066a7e14dcfSSatish Balay .vb
2067441846f8SBarry Smith    TaoRegister("my_solver",MySolverCreate);
2068a7e14dcfSSatish Balay .ve
2069a7e14dcfSSatish Balay 
2070a7e14dcfSSatish Balay    Then, your solver can be chosen with the procedural interface via
2071a7e14dcfSSatish Balay $     TaoSetType(tao,"my_solver")
2072a7e14dcfSSatish Balay    or at runtime via the option
207358417fe7SBarry Smith $     -tao_type my_solver
2074a7e14dcfSSatish Balay 
2075a7e14dcfSSatish Balay    Level: advanced
2076a7e14dcfSSatish Balay 
2077441846f8SBarry Smith .seealso: TaoRegisterAll(), TaoRegisterDestroy()
2078a7e14dcfSSatish Balay M*/
2079441846f8SBarry Smith PetscErrorCode TaoRegister(const char sname[], PetscErrorCode (*func)(Tao))
2080a7e14dcfSSatish Balay {
2081a7e14dcfSSatish Balay   PetscErrorCode ierr;
2082a7e14dcfSSatish Balay 
2083a7e14dcfSSatish Balay   PetscFunctionBegin;
2084441846f8SBarry Smith   ierr = PetscFunctionListAdd(&TaoList,sname, (void (*)(void))func);CHKERRQ(ierr);
2085a7e14dcfSSatish Balay   PetscFunctionReturn(0);
2086a7e14dcfSSatish Balay }
2087a7e14dcfSSatish Balay 
2088a7e14dcfSSatish Balay #undef __FUNCT__
2089441846f8SBarry Smith #define __FUNCT__ "TaoRegisterDestroy"
2090a7e14dcfSSatish Balay /*@C
2091441846f8SBarry Smith    TaoRegisterDestroy - Frees the list of minimization solvers that were
2092441846f8SBarry Smith    registered by TaoRegisterDynamic().
2093a7e14dcfSSatish Balay 
2094a7e14dcfSSatish Balay    Not Collective
2095a7e14dcfSSatish Balay 
2096a7e14dcfSSatish Balay    Level: advanced
2097a7e14dcfSSatish Balay 
2098441846f8SBarry Smith .seealso: TaoRegisterAll(), TaoRegister()
2099a7e14dcfSSatish Balay @*/
2100441846f8SBarry Smith PetscErrorCode TaoRegisterDestroy(void)
2101a7e14dcfSSatish Balay {
2102a7e14dcfSSatish Balay   PetscErrorCode ierr;
2103a7e14dcfSSatish Balay   PetscFunctionBegin;
2104441846f8SBarry Smith   ierr = PetscFunctionListDestroy(&TaoList);CHKERRQ(ierr);
2105441846f8SBarry Smith   TaoRegisterAllCalled = PETSC_FALSE;
2106a7e14dcfSSatish Balay   PetscFunctionReturn(0);
2107a7e14dcfSSatish Balay }
210847a47007SBarry Smith 
2109a7e14dcfSSatish Balay #undef __FUNCT__
2110*e4cb33bbSBarry Smith #define __FUNCT__ "TaoSetConvergedReason"
2111a7e14dcfSSatish Balay /*@
2112*e4cb33bbSBarry Smith   TaoSetConvergedReason - Sets the termination flag on a Tao object
2113a7e14dcfSSatish Balay 
2114441846f8SBarry Smith   Logically Collective on Tao
2115a7e14dcfSSatish Balay 
2116a7e14dcfSSatish Balay   Input Parameters:
2117441846f8SBarry Smith + tao - the Tao context
2118a7e14dcfSSatish Balay - reason - one of
2119a7e14dcfSSatish Balay $     TAO_CONVERGED_ATOL (2),
2120a7e14dcfSSatish Balay $     TAO_CONVERGED_RTOL (3),
2121a7e14dcfSSatish Balay $     TAO_CONVERGED_STEPTOL (4),
2122a7e14dcfSSatish Balay $     TAO_CONVERGED_MINF (5),
2123a7e14dcfSSatish Balay $     TAO_CONVERGED_USER (6),
2124a7e14dcfSSatish Balay $     TAO_DIVERGED_MAXITS (-2),
2125a7e14dcfSSatish Balay $     TAO_DIVERGED_NAN (-4),
2126a7e14dcfSSatish Balay $     TAO_DIVERGED_MAXFCN (-5),
2127a7e14dcfSSatish Balay $     TAO_DIVERGED_LS_FAILURE (-6),
2128a7e14dcfSSatish Balay $     TAO_DIVERGED_TR_REDUCTION (-7),
2129a7e14dcfSSatish Balay $     TAO_DIVERGED_USER (-8),
2130a7e14dcfSSatish Balay $     TAO_CONTINUE_ITERATING (0)
2131a7e14dcfSSatish Balay 
2132a7e14dcfSSatish Balay    Level: intermediate
2133a7e14dcfSSatish Balay 
2134a7e14dcfSSatish Balay @*/
2135*e4cb33bbSBarry Smith PetscErrorCode TaoSetConvergedReason(Tao tao, TaoConvergedReason reason)
2136a7e14dcfSSatish Balay {
2137441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2138a7e14dcfSSatish Balay   PetscFunctionBegin;
2139a7e14dcfSSatish Balay   tao->reason = reason;
2140a7e14dcfSSatish Balay   PetscFunctionReturn(0);
2141a7e14dcfSSatish Balay }
2142a7e14dcfSSatish Balay 
2143a7e14dcfSSatish Balay #undef __FUNCT__
2144*e4cb33bbSBarry Smith #define __FUNCT__ "TaoGetConvergedReason"
2145a7e14dcfSSatish Balay /*@
2146*e4cb33bbSBarry Smith    TaoGetConvergedReason - Gets the reason the Tao iteration was stopped.
2147a7e14dcfSSatish Balay 
2148a7e14dcfSSatish Balay    Not Collective
2149a7e14dcfSSatish Balay 
2150a7e14dcfSSatish Balay    Input Parameter:
2151441846f8SBarry Smith .  tao - the Tao solver context
2152a7e14dcfSSatish Balay 
2153a7e14dcfSSatish Balay    Output Parameter:
2154a7e14dcfSSatish Balay .  reason - one of
2155a7e14dcfSSatish Balay 
2156a7e14dcfSSatish Balay $  TAO_CONVERGED_FATOL (1)           f(X)-f(X*) <= fatol
2157a7e14dcfSSatish Balay $  TAO_CONVERGED_FRTOL (2)           |f(X) - f(X*)|/|f(X)| < frtol
2158a7e14dcfSSatish Balay $  TAO_CONVERGED_GATOL (3)           ||g(X)|| < gatol
2159a7e14dcfSSatish Balay $  TAO_CONVERGED_GRTOL (4)           ||g(X)|| / f(X)  < grtol
2160a7e14dcfSSatish Balay $  TAO_CONVERGED_GTTOL (5)           ||g(X)|| / ||g(X0)|| < gttol
2161a7e14dcfSSatish Balay $  TAO_CONVERGED_STEPTOL (6)         step size small
2162a7e14dcfSSatish Balay $  TAO_CONVERGED_MINF (7)            F < F_min
2163a7e14dcfSSatish Balay $  TAO_CONVERGED_USER (8)            User defined
2164a7e14dcfSSatish Balay 
2165a7e14dcfSSatish Balay $  TAO_DIVERGED_MAXITS (-2)          its > maxits
2166a7e14dcfSSatish Balay $  TAO_DIVERGED_NAN (-4)             Numerical problems
2167a7e14dcfSSatish Balay $  TAO_DIVERGED_MAXFCN (-5)          fevals > max_funcsals
2168a7e14dcfSSatish Balay $  TAO_DIVERGED_LS_FAILURE (-6)      line search failure
2169a7e14dcfSSatish Balay $  TAO_DIVERGED_TR_REDUCTION (-7)    trust region failure
2170a7e14dcfSSatish Balay $  TAO_DIVERGED_USER(-8)             (user defined)
2171a7e14dcfSSatish Balay 
2172a7e14dcfSSatish Balay $  TAO_CONTINUE_ITERATING (0)
2173a7e14dcfSSatish Balay 
2174a7e14dcfSSatish Balay    where
2175a7e14dcfSSatish Balay +  X - current solution
2176a7e14dcfSSatish Balay .  X0 - initial guess
2177a7e14dcfSSatish Balay .  f(X) - current function value
2178a7e14dcfSSatish Balay .  f(X*) - true solution (estimated)
2179a7e14dcfSSatish Balay .  g(X) - current gradient
2180a7e14dcfSSatish Balay .  its - current iterate number
2181a7e14dcfSSatish Balay .  maxits - maximum number of iterates
2182a7e14dcfSSatish Balay .  fevals - number of function evaluations
2183a7e14dcfSSatish Balay -  max_funcsals - maximum number of function evaluations
2184a7e14dcfSSatish Balay 
2185a7e14dcfSSatish Balay    Level: intermediate
2186a7e14dcfSSatish Balay 
2187a7e14dcfSSatish Balay .seealso: TaoSetConvergenceTest(), TaoSetTolerances()
2188a7e14dcfSSatish Balay 
2189a7e14dcfSSatish Balay @*/
2190*e4cb33bbSBarry Smith PetscErrorCode TaoGetConvergedReason(Tao tao, TaoConvergedReason *reason)
2191a7e14dcfSSatish Balay {
2192a7e14dcfSSatish Balay   PetscFunctionBegin;
2193441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2194a7e14dcfSSatish Balay   PetscValidPointer(reason,2);
2195a7e14dcfSSatish Balay   *reason = tao->reason;
2196a7e14dcfSSatish Balay   PetscFunctionReturn(0);
2197a7e14dcfSSatish Balay }
2198a7e14dcfSSatish Balay 
2199a7e14dcfSSatish Balay #undef __FUNCT__
2200a7e14dcfSSatish Balay #define __FUNCT__ "TaoGetSolutionStatus"
2201a7e14dcfSSatish Balay /*@
2202a7e14dcfSSatish Balay   TaoGetSolutionStatus - Get the current iterate, objective value,
2203a7e14dcfSSatish Balay   residual, infeasibility, and termination
2204a7e14dcfSSatish Balay 
2205a7e14dcfSSatish Balay   Not Collective
2206a7e14dcfSSatish Balay 
2207a7e14dcfSSatish Balay    Input Parameters:
2208441846f8SBarry Smith .  tao - the Tao context
2209a7e14dcfSSatish Balay 
2210a7e14dcfSSatish Balay    Output Parameters:
2211a7e14dcfSSatish Balay +  iterate - the current iterate number (>=0)
2212a7e14dcfSSatish Balay .  f - the current function value
221358417fe7SBarry Smith .  gnorm - the square of the gradient norm, duality gap, or other measure indicating distance from optimality.
2214a7e14dcfSSatish Balay .  cnorm - the infeasibility of the current solution with regard to the constraints.
2215a7e14dcfSSatish Balay .  xdiff - the step length or trust region radius of the most recent iterate.
2216a7e14dcfSSatish Balay -  reason - The termination reason, which can equal TAO_CONTINUE_ITERATING
2217a7e14dcfSSatish Balay 
2218a7e14dcfSSatish Balay    Level: intermediate
2219a7e14dcfSSatish Balay 
2220a7e14dcfSSatish Balay    Note:
2221a7e14dcfSSatish Balay    TAO returns the values set by the solvers in the routine TaoMonitor().
2222a7e14dcfSSatish Balay 
2223a7e14dcfSSatish Balay    Note:
222447a47007SBarry Smith    If any of the output arguments are set to NULL, no corresponding value will be returned.
2225a7e14dcfSSatish Balay 
2226*e4cb33bbSBarry Smith .seealso: TaoMonitor(), TaoGetConvergedReason()
2227a7e14dcfSSatish Balay @*/
2228*e4cb33bbSBarry Smith PetscErrorCode TaoGetSolutionStatus(Tao tao, PetscInt *its, PetscReal *f, PetscReal *gnorm, PetscReal *cnorm, PetscReal *xdiff, TaoConvergedReason *reason)
2229a7e14dcfSSatish Balay {
2230a7e14dcfSSatish Balay   PetscFunctionBegin;
2231a7e14dcfSSatish Balay   if (its) *its=tao->niter;
2232a7e14dcfSSatish Balay   if (f) *f=tao->fc;
2233a7e14dcfSSatish Balay   if (gnorm) *gnorm=tao->residual;
2234a7e14dcfSSatish Balay   if (cnorm) *cnorm=tao->cnorm;
2235a7e14dcfSSatish Balay   if (reason) *reason=tao->reason;
2236a7e14dcfSSatish Balay   if (xdiff) *xdiff=tao->step;
2237a7e14dcfSSatish Balay   PetscFunctionReturn(0);
2238a7e14dcfSSatish Balay }
2239a7e14dcfSSatish Balay 
2240a7e14dcfSSatish Balay #undef __FUNCT__
2241a7e14dcfSSatish Balay #define __FUNCT__ "TaoGetType"
2242a7e14dcfSSatish Balay /*@C
2243441846f8SBarry Smith    TaoGetType - Gets the current Tao algorithm.
2244a7e14dcfSSatish Balay 
2245a7e14dcfSSatish Balay    Not Collective
2246a7e14dcfSSatish Balay 
2247a7e14dcfSSatish Balay    Input Parameter:
2248441846f8SBarry Smith .  tao - the Tao solver context
2249a7e14dcfSSatish Balay 
2250a7e14dcfSSatish Balay    Output Parameter:
2251441846f8SBarry Smith .  type - Tao method
2252a7e14dcfSSatish Balay 
2253a7e14dcfSSatish Balay    Level: intermediate
2254a7e14dcfSSatish Balay 
2255a7e14dcfSSatish Balay @*/
2256441846f8SBarry Smith PetscErrorCode TaoGetType(Tao tao, const TaoType *type)
2257a7e14dcfSSatish Balay {
2258a7e14dcfSSatish Balay   PetscFunctionBegin;
2259441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2260a7e14dcfSSatish Balay   PetscValidPointer(type,2);
2261a7e14dcfSSatish Balay   *type=((PetscObject)tao)->type_name;
2262a7e14dcfSSatish Balay   PetscFunctionReturn(0);
2263a7e14dcfSSatish Balay }
2264a7e14dcfSSatish Balay 
2265a7e14dcfSSatish Balay #undef __FUNCT__
2266a7e14dcfSSatish Balay #define __FUNCT__ "TaoMonitor"
2267a7e14dcfSSatish Balay /*@C
2268a7e14dcfSSatish Balay   TaoMonitor - Monitor the solver and the current solution.  This
2269a7e14dcfSSatish Balay   routine will record the iteration number and residual statistics,
2270a7e14dcfSSatish Balay   call any monitors specified by the user, and calls the convergence-check routine.
2271a7e14dcfSSatish Balay 
2272a7e14dcfSSatish Balay    Input Parameters:
2273441846f8SBarry Smith +  tao - the Tao context
2274a7e14dcfSSatish Balay .  its - the current iterate number (>=0)
2275a7e14dcfSSatish Balay .  f - the current objective function value
227658417fe7SBarry Smith .  res - the gradient norm, square root of the duality gap, or other measure indicating distince from optimality.  This measure will be recorded and
2277a7e14dcfSSatish Balay           used for some termination tests.
2278a7e14dcfSSatish Balay .  cnorm - the infeasibility of the current solution with regard to the constraints.
2279a7e14dcfSSatish Balay -  steplength - multiple of the step direction added to the previous iterate.
2280a7e14dcfSSatish Balay 
2281a7e14dcfSSatish Balay    Output Parameters:
2282a7e14dcfSSatish Balay .  reason - The termination reason, which can equal TAO_CONTINUE_ITERATING
2283a7e14dcfSSatish Balay 
2284a7e14dcfSSatish Balay    Options Database Key:
2285a7e14dcfSSatish Balay .  -tao_monitor - Use the default monitor, which prints statistics to standard output
2286a7e14dcfSSatish Balay 
2287*e4cb33bbSBarry Smith .seealso TaoGetConvergedReason(), TaoDefaultMonitor(), TaoSetMonitor()
2288a7e14dcfSSatish Balay 
2289a7e14dcfSSatish Balay    Level: developer
2290a7e14dcfSSatish Balay 
2291a7e14dcfSSatish Balay @*/
2292*e4cb33bbSBarry Smith PetscErrorCode TaoMonitor(Tao tao, PetscInt its, PetscReal f, PetscReal res, PetscReal cnorm, PetscReal steplength, TaoConvergedReason *reason)
2293a7e14dcfSSatish Balay {
2294a7e14dcfSSatish Balay   PetscErrorCode ierr;
229547a47007SBarry Smith   PetscInt       i;
229647a47007SBarry Smith 
2297a7e14dcfSSatish Balay   PetscFunctionBegin;
2298441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2299a7e14dcfSSatish Balay   tao->fc = f;
2300a7e14dcfSSatish Balay   tao->residual = res;
2301a7e14dcfSSatish Balay   tao->cnorm = cnorm;
2302a7e14dcfSSatish Balay   tao->step = steplength;
2303a7e14dcfSSatish Balay   tao->niter=its;
2304a7e14dcfSSatish Balay   if (its == 0) {
2305a7e14dcfSSatish Balay     tao->cnorm0 = cnorm; tao->gnorm0 = res;
2306a7e14dcfSSatish Balay   }
2307a7e14dcfSSatish Balay   TaoLogHistory(tao,f,res,cnorm);
230847a47007SBarry Smith   if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(res)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
2309a7e14dcfSSatish Balay   if (tao->ops->convergencetest) {
2310a7e14dcfSSatish Balay     ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
2311a7e14dcfSSatish Balay   }
2312a7e14dcfSSatish Balay   for (i=0;i<tao->numbermonitors;i++) {
2313a7e14dcfSSatish Balay     ierr = (*tao->monitor[i])(tao,tao->monitorcontext[i]);CHKERRQ(ierr);
2314a7e14dcfSSatish Balay   }
2315a7e14dcfSSatish Balay   *reason = tao->reason;
2316a7e14dcfSSatish Balay   PetscFunctionReturn(0);
2317a7e14dcfSSatish Balay }
2318a7e14dcfSSatish Balay 
2319a7e14dcfSSatish Balay #undef __FUNCT__
2320a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetHistory"
2321a7e14dcfSSatish Balay /*@
2322a7e14dcfSSatish Balay    TaoSetHistory - Sets the array used to hold the convergence history.
2323a7e14dcfSSatish Balay 
2324441846f8SBarry Smith    Logically Collective on Tao
2325a7e14dcfSSatish Balay 
2326a7e14dcfSSatish Balay    Input Parameters:
2327441846f8SBarry Smith +  tao - the Tao solver context
2328a7e14dcfSSatish Balay .  obj   - array to hold objective value history
2329a7e14dcfSSatish Balay .  resid - array to hold residual history
2330a7e14dcfSSatish Balay .  cnorm - array to hold constraint violation history
2331a7e14dcfSSatish Balay .  na  - size of obj, resid, and cnorm
2332a7e14dcfSSatish Balay -  reset - PetscTrue indicates each new minimization resets the history counter to zero,
2333a7e14dcfSSatish Balay            else it continues storing new values for new minimizations after the old ones
2334a7e14dcfSSatish Balay 
2335a7e14dcfSSatish Balay    Notes:
2336a7e14dcfSSatish Balay    If set, TAO will fill the given arrays with the indicated
2337a7e14dcfSSatish Balay    information at each iteration.  If no information is desired
233847a47007SBarry Smith    for a given array, then NULL may be used.
2339a7e14dcfSSatish Balay 
2340a7e14dcfSSatish Balay    This routine is useful, e.g., when running a code for purposes
2341a7e14dcfSSatish Balay    of accurate performance monitoring, when no I/O should be done
2342a7e14dcfSSatish Balay    during the section of code that is being timed.
2343a7e14dcfSSatish Balay 
2344a7e14dcfSSatish Balay    Level: intermediate
2345a7e14dcfSSatish Balay 
2346a7e14dcfSSatish Balay .seealso: TaoGetHistory()
2347a7e14dcfSSatish Balay 
2348a7e14dcfSSatish Balay @*/
2349441846f8SBarry Smith PetscErrorCode TaoSetHistory(Tao tao, PetscReal *obj, PetscReal *resid, PetscReal *cnorm, PetscInt na,PetscBool reset)
2350a7e14dcfSSatish Balay {
2351a7e14dcfSSatish Balay   PetscFunctionBegin;
2352441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2353a7e14dcfSSatish Balay   tao->hist_obj = obj;
2354a7e14dcfSSatish Balay   tao->hist_resid = resid;
2355a7e14dcfSSatish Balay   tao->hist_cnorm = cnorm;
2356a7e14dcfSSatish Balay   tao->hist_max   = na;
2357a7e14dcfSSatish Balay   tao->hist_reset = reset;
2358a7e14dcfSSatish Balay   PetscFunctionReturn(0);
2359a7e14dcfSSatish Balay }
2360a7e14dcfSSatish Balay 
2361a7e14dcfSSatish Balay #undef __FUNCT__
2362a7e14dcfSSatish Balay #define __FUNCT__ "TaoGetHistory"
2363a7e14dcfSSatish Balay /*@C
2364a7e14dcfSSatish Balay    TaoGetHistory - Gets the array used to hold the convergence history.
2365a7e14dcfSSatish Balay 
2366441846f8SBarry Smith    Collective on Tao
2367a7e14dcfSSatish Balay 
2368a7e14dcfSSatish Balay    Input Parameter:
2369441846f8SBarry Smith .  tao - the Tao context
2370a7e14dcfSSatish Balay 
2371a7e14dcfSSatish Balay    Output Parameters:
2372a7e14dcfSSatish Balay +  obj   - array used to hold objective value history
2373a7e14dcfSSatish Balay .  resid - array used to hold residual history
2374a7e14dcfSSatish Balay .  cnorm - array used to hold constraint violation history
2375a7e14dcfSSatish Balay -  nhist  - size of obj, resid, and cnorm (will be less than or equal to na given in TaoSetHistory)
2376a7e14dcfSSatish Balay 
2377a7e14dcfSSatish Balay    Notes:
2378a7e14dcfSSatish Balay     The calling sequence for this routine in Fortran is
2379441846f8SBarry Smith $   call TaoGetHistory(Tao tao, integer nhist, integer info)
2380a7e14dcfSSatish Balay 
2381a7e14dcfSSatish Balay    This routine is useful, e.g., when running a code for purposes
2382a7e14dcfSSatish Balay    of accurate performance monitoring, when no I/O should be done
2383a7e14dcfSSatish Balay    during the section of code that is being timed.
2384a7e14dcfSSatish Balay 
2385a7e14dcfSSatish Balay    Level: advanced
2386a7e14dcfSSatish Balay 
2387a7e14dcfSSatish Balay .seealso: TaoSetHistory()
2388a7e14dcfSSatish Balay 
2389a7e14dcfSSatish Balay @*/
2390441846f8SBarry Smith PetscErrorCode TaoGetHistory(Tao tao, PetscReal **obj, PetscReal **resid, PetscReal **cnorm, PetscInt *nhist)
2391a7e14dcfSSatish Balay {
2392a7e14dcfSSatish Balay   PetscFunctionBegin;
2393441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2394a7e14dcfSSatish Balay   if (obj)   *obj   = tao->hist_obj;
2395a7e14dcfSSatish Balay   if (cnorm) *cnorm = tao->hist_cnorm;
2396a7e14dcfSSatish Balay   if (resid) *resid = tao->hist_resid;
2397a7e14dcfSSatish Balay   if (nhist) *nhist   = tao->hist_len;
2398a7e14dcfSSatish Balay   PetscFunctionReturn(0);
2399a7e14dcfSSatish Balay }
2400a7e14dcfSSatish Balay 
2401a7e14dcfSSatish Balay #undef __FUNCT__
2402a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetApplicationContext"
2403a7e14dcfSSatish Balay /*@
2404a7e14dcfSSatish Balay    TaoSetApplicationContext - Sets the optional user-defined context for
2405a7e14dcfSSatish Balay    a solver.
2406a7e14dcfSSatish Balay 
2407441846f8SBarry Smith    Logically Collective on Tao
2408a7e14dcfSSatish Balay 
2409a7e14dcfSSatish Balay    Input Parameters:
2410441846f8SBarry Smith +  tao  - the Tao context
2411a7e14dcfSSatish Balay -  usrP - optional user context
2412a7e14dcfSSatish Balay 
2413a7e14dcfSSatish Balay    Level: intermediate
2414a7e14dcfSSatish Balay 
2415a7e14dcfSSatish Balay .seealso: TaoGetApplicationContext(), TaoSetApplicationContext()
2416a7e14dcfSSatish Balay @*/
2417441846f8SBarry Smith PetscErrorCode  TaoSetApplicationContext(Tao tao,void *usrP)
2418a7e14dcfSSatish Balay {
2419a7e14dcfSSatish Balay   PetscFunctionBegin;
2420441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2421a7e14dcfSSatish Balay   tao->user = usrP;
2422a7e14dcfSSatish Balay   PetscFunctionReturn(0);
2423a7e14dcfSSatish Balay }
2424a7e14dcfSSatish Balay 
2425a7e14dcfSSatish Balay #undef __FUNCT__
2426a7e14dcfSSatish Balay #define __FUNCT__ "TaoGetApplicationContext"
2427a7e14dcfSSatish Balay /*@
2428a7e14dcfSSatish Balay    TaoGetApplicationContext - Gets the user-defined context for a
2429a7e14dcfSSatish Balay    TAO solvers.
2430a7e14dcfSSatish Balay 
2431a7e14dcfSSatish Balay    Not Collective
2432a7e14dcfSSatish Balay 
2433a7e14dcfSSatish Balay    Input Parameter:
2434441846f8SBarry Smith .  tao  - Tao context
2435a7e14dcfSSatish Balay 
2436a7e14dcfSSatish Balay    Output Parameter:
2437a7e14dcfSSatish Balay .  usrP - user context
2438a7e14dcfSSatish Balay 
2439a7e14dcfSSatish Balay    Level: intermediate
2440a7e14dcfSSatish Balay 
2441a7e14dcfSSatish Balay .seealso: TaoSetApplicationContext()
2442a7e14dcfSSatish Balay @*/
2443441846f8SBarry Smith PetscErrorCode  TaoGetApplicationContext(Tao tao,void *usrP)
2444a7e14dcfSSatish Balay {
2445a7e14dcfSSatish Balay   PetscFunctionBegin;
2446441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
2447a7e14dcfSSatish Balay   *(void**)usrP = tao->user;
2448a7e14dcfSSatish Balay   PetscFunctionReturn(0);
2449a7e14dcfSSatish Balay }
2450