xref: /petsc/src/tao/interface/taosolver.c (revision 1690c2ae071c7584458d4e437df7b47bc4686b3c)
1af0996ceSBarry Smith #include <petsc/private/taoimpl.h> /*I "petsctao.h" I*/
20f0abf79SStefano Zampini #include <petsc/private/snesimpl.h>
3a7e14dcfSSatish Balay 
4441846f8SBarry Smith PetscBool         TaoRegisterAllCalled = PETSC_FALSE;
5441846f8SBarry Smith PetscFunctionList TaoList              = NULL;
6a7e14dcfSSatish Balay 
7441846f8SBarry Smith PetscClassId TAO_CLASSID;
80ebee16dSLisandro Dalcin 
90ebee16dSLisandro Dalcin PetscLogEvent TAO_Solve;
100ebee16dSLisandro Dalcin PetscLogEvent TAO_ObjectiveEval;
110ebee16dSLisandro Dalcin PetscLogEvent TAO_GradientEval;
120ebee16dSLisandro Dalcin PetscLogEvent TAO_ObjGradEval;
130ebee16dSLisandro Dalcin PetscLogEvent TAO_HessianEval;
140ebee16dSLisandro Dalcin PetscLogEvent TAO_JacobianEval;
150ebee16dSLisandro Dalcin PetscLogEvent TAO_ConstraintsEval;
16a7e14dcfSSatish Balay 
1783c8fe1dSLisandro Dalcin const char *TaoSubSetTypes[] = {"subvec", "mask", "matrixfree", "TaoSubSetType", "TAO_SUBSET_", NULL};
18a7e14dcfSSatish Balay 
19e882e171SHong Zhang struct _n_TaoMonitorDrawCtx {
20e882e171SHong Zhang   PetscViewer viewer;
21e882e171SHong Zhang   PetscInt    howoften; /* when > 0 uses iteration % howoften, when negative only final solution plotted */
22e882e171SHong Zhang };
23e882e171SHong Zhang 
24f1e99121SPierre Jolivet static PetscErrorCode KSPPreSolve_TAOEW_Private(KSP ksp, Vec b, Vec x, void *ctx)
25d71ae5a4SJacob Faibussowitsch {
26f1e99121SPierre Jolivet   Tao  tao          = (Tao)ctx;
270f0abf79SStefano Zampini   SNES snes_ewdummy = tao->snes_ewdummy;
280f0abf79SStefano Zampini 
290f0abf79SStefano Zampini   PetscFunctionBegin;
303ba16761SJacob Faibussowitsch   if (!snes_ewdummy) PetscFunctionReturn(PETSC_SUCCESS);
310f0abf79SStefano Zampini   /* populate snes_ewdummy struct values used in KSPPreSolve_SNESEW */
320f0abf79SStefano Zampini   snes_ewdummy->vec_func = b;
330f0abf79SStefano Zampini   snes_ewdummy->rtol     = tao->gttol;
340f0abf79SStefano Zampini   snes_ewdummy->iter     = tao->niter;
350f0abf79SStefano Zampini   PetscCall(VecNorm(b, NORM_2, &snes_ewdummy->norm));
360f0abf79SStefano Zampini   PetscCall(KSPPreSolve_SNESEW(ksp, b, x, snes_ewdummy));
370f0abf79SStefano Zampini   snes_ewdummy->vec_func = NULL;
383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
390f0abf79SStefano Zampini }
400f0abf79SStefano Zampini 
41f1e99121SPierre Jolivet static PetscErrorCode KSPPostSolve_TAOEW_Private(KSP ksp, Vec b, Vec x, void *ctx)
42d71ae5a4SJacob Faibussowitsch {
43f1e99121SPierre Jolivet   Tao  tao          = (Tao)ctx;
440f0abf79SStefano Zampini   SNES snes_ewdummy = tao->snes_ewdummy;
450f0abf79SStefano Zampini 
460f0abf79SStefano Zampini   PetscFunctionBegin;
473ba16761SJacob Faibussowitsch   if (!snes_ewdummy) PetscFunctionReturn(PETSC_SUCCESS);
480f0abf79SStefano Zampini   PetscCall(KSPPostSolve_SNESEW(ksp, b, x, snes_ewdummy));
493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
500f0abf79SStefano Zampini }
510f0abf79SStefano Zampini 
52d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSetUpEW_Private(Tao tao)
53d71ae5a4SJacob Faibussowitsch {
540f0abf79SStefano Zampini   SNESKSPEW  *kctx;
550f0abf79SStefano Zampini   const char *ewprefix;
560f0abf79SStefano Zampini 
570f0abf79SStefano Zampini   PetscFunctionBegin;
583ba16761SJacob Faibussowitsch   if (!tao->ksp) PetscFunctionReturn(PETSC_SUCCESS);
590f0abf79SStefano Zampini   if (tao->ksp_ewconv) {
600f0abf79SStefano Zampini     if (!tao->snes_ewdummy) PetscCall(SNESCreate(PetscObjectComm((PetscObject)tao), &tao->snes_ewdummy));
610f0abf79SStefano Zampini     tao->snes_ewdummy->ksp_ewconv = PETSC_TRUE;
62f1e99121SPierre Jolivet     PetscCall(KSPSetPreSolve(tao->ksp, KSPPreSolve_TAOEW_Private, tao));
63f1e99121SPierre Jolivet     PetscCall(KSPSetPostSolve(tao->ksp, KSPPostSolve_TAOEW_Private, tao));
640f0abf79SStefano Zampini 
650f0abf79SStefano Zampini     PetscCall(KSPGetOptionsPrefix(tao->ksp, &ewprefix));
660f0abf79SStefano Zampini     kctx = (SNESKSPEW *)tao->snes_ewdummy->kspconvctx;
67a4598233SStefano Zampini     PetscCall(SNESEWSetFromOptions_Private(kctx, PETSC_FALSE, PetscObjectComm((PetscObject)tao), ewprefix));
680f0abf79SStefano Zampini   } else PetscCall(SNESDestroy(&tao->snes_ewdummy));
693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
700f0abf79SStefano Zampini }
710f0abf79SStefano Zampini 
72a7e14dcfSSatish Balay /*@
73606f75f6SBarry Smith   TaoParametersInitialize - Sets all the parameters in `tao` to their default value (when `TaoCreate()` was called) if they
74606f75f6SBarry Smith   currently contain default values. Default values are the parameter values when the object's type is set.
75606f75f6SBarry Smith 
76606f75f6SBarry Smith   Collective
77606f75f6SBarry Smith 
78606f75f6SBarry Smith   Input Parameter:
79606f75f6SBarry Smith . tao - the `Tao` object
80606f75f6SBarry Smith 
81606f75f6SBarry Smith   Level: developer
82606f75f6SBarry Smith 
83606f75f6SBarry Smith   Developer Note:
84606f75f6SBarry Smith   This is called by all the `TaoCreate_XXX()` routines.
85606f75f6SBarry Smith 
86606f75f6SBarry Smith .seealso: [](ch_snes), `Tao`, `TaoSolve()`, `TaoDestroy()`,
87606f75f6SBarry Smith           `PetscObjectParameterSetDefault()`
88606f75f6SBarry Smith @*/
89606f75f6SBarry Smith PetscErrorCode TaoParametersInitialize(Tao tao)
90606f75f6SBarry Smith {
91606f75f6SBarry Smith   PetscObjectParameterSetDefault(tao, max_it, 10000);
92606f75f6SBarry Smith   PetscObjectParameterSetDefault(tao, max_funcs, PETSC_UNLIMITED);
93606f75f6SBarry Smith   PetscObjectParameterSetDefault(tao, gatol, PetscDefined(USE_REAL_SINGLE) ? 1e-5 : 1e-8);
94606f75f6SBarry Smith   PetscObjectParameterSetDefault(tao, grtol, PetscDefined(USE_REAL_SINGLE) ? 1e-5 : 1e-8);
95606f75f6SBarry Smith   PetscObjectParameterSetDefault(tao, crtol, PetscDefined(USE_REAL_SINGLE) ? 1e-5 : 1e-8);
96606f75f6SBarry Smith   PetscObjectParameterSetDefault(tao, catol, PetscDefined(USE_REAL_SINGLE) ? 1e-5 : 1e-8);
97606f75f6SBarry Smith   PetscObjectParameterSetDefault(tao, gttol, 0.0);
98606f75f6SBarry Smith   PetscObjectParameterSetDefault(tao, steptol, 0.0);
99606f75f6SBarry Smith   PetscObjectParameterSetDefault(tao, fmin, PETSC_NINFINITY);
100606f75f6SBarry Smith   PetscObjectParameterSetDefault(tao, trust0, PETSC_INFINITY);
101606f75f6SBarry Smith   return PETSC_SUCCESS;
102606f75f6SBarry Smith }
103606f75f6SBarry Smith 
104606f75f6SBarry Smith /*@
10565ba42b6SBarry Smith   TaoCreate - Creates a Tao solver
106a7e14dcfSSatish Balay 
107d083f849SBarry Smith   Collective
108a7e14dcfSSatish Balay 
109a7e14dcfSSatish Balay   Input Parameter:
110a7e14dcfSSatish Balay . comm - MPI communicator
111a7e14dcfSSatish Balay 
112a7e14dcfSSatish Balay   Output Parameter:
11347450a7bSBarry Smith . newtao - the new `Tao` context
114a7e14dcfSSatish Balay 
11547450a7bSBarry Smith   Options Database Key:
11665ba42b6SBarry Smith . -tao_type - select which method Tao should use
117a7e14dcfSSatish Balay 
118a7e14dcfSSatish Balay   Level: beginner
119a7e14dcfSSatish Balay 
1200241b274SPierre Jolivet .seealso: [](ch_tao), `Tao`, `TaoSolve()`, `TaoDestroy()`, `TaoSetFromOptions()`, `TaoSetType()`
121a7e14dcfSSatish Balay @*/
122d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoCreate(MPI_Comm comm, Tao *newtao)
123d71ae5a4SJacob Faibussowitsch {
124441846f8SBarry Smith   Tao tao;
125a7e14dcfSSatish Balay 
126a7e14dcfSSatish Balay   PetscFunctionBegin;
1274f572ea9SToby Isaac   PetscAssertPointer(newtao, 2);
1289566063dSJacob Faibussowitsch   PetscCall(TaoInitializePackage());
1299566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchInitializePackage());
130a7e14dcfSSatish Balay 
131606f75f6SBarry Smith   PetscCall(PetscHeaderCreate(tao, TAO_CLASSID, "Tao", "Optimization solver", "Tao", comm, TaoDestroy, TaoView));
13294511df7SStefano Zampini   tao->ops->convergencetest = TaoDefaultConvergenceTest;
133606f75f6SBarry Smith 
134a7e14dcfSSatish Balay   tao->hist_reset = PETSC_TRUE;
1359566063dSJacob Faibussowitsch   PetscCall(TaoResetStatistics(tao));
136a7e14dcfSSatish Balay   *newtao = tao;
1373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
138a7e14dcfSSatish Balay }
139a7e14dcfSSatish Balay 
140a7e14dcfSSatish Balay /*@
141a7e14dcfSSatish Balay   TaoSolve - Solves an optimization problem min F(x) s.t. l <= x <= u
142a7e14dcfSSatish Balay 
143c3339decSBarry Smith   Collective
144a7e14dcfSSatish Balay 
14547450a7bSBarry Smith   Input Parameter:
14647450a7bSBarry Smith . tao - the `Tao` context
147a7e14dcfSSatish Balay 
14867be906fSBarry Smith   Level: beginner
14967be906fSBarry Smith 
150a7e14dcfSSatish Balay   Notes:
15147450a7bSBarry Smith   The user must set up the `Tao` object  with calls to `TaoSetSolution()`, `TaoSetObjective()`, `TaoSetGradient()`, and (if using 2nd order method) `TaoSetHessian()`.
152a7e14dcfSSatish Balay 
15367be906fSBarry Smith   You should call `TaoGetConvergedReason()` or run with `-tao_converged_reason` to determine if the optimization algorithm actually succeeded or
154a35d58b8SBarry Smith   why it failed.
155a35d58b8SBarry Smith 
1561cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoCreate()`, `TaoSetObjective()`, `TaoSetGradient()`, `TaoSetHessian()`, `TaoGetConvergedReason()`, `TaoSetUp()`
157a7e14dcfSSatish Balay  @*/
158d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoSolve(Tao tao)
159d71ae5a4SJacob Faibussowitsch {
160e2379f4fSBarry Smith   static PetscBool set = PETSC_FALSE;
16147a47007SBarry Smith 
162a7e14dcfSSatish Balay   PetscFunctionBegin;
163441846f8SBarry Smith   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
164d0609cedSBarry Smith   PetscCall(PetscCitationsRegister("@TechReport{tao-user-ref,\n"
165e2379f4fSBarry Smith                                    "title   = {Toolkit for Advanced Optimization (TAO) Users Manual},\n"
166e2379f4fSBarry Smith                                    "author  = {Todd Munson and Jason Sarich and Stefan Wild and Steve Benson and Lois Curfman McInnes},\n"
167e2379f4fSBarry Smith                                    "Institution = {Argonne National Laboratory},\n"
168e2379f4fSBarry Smith                                    "Year   = 2014,\n"
169e2379f4fSBarry Smith                                    "Number = {ANL/MCS-TM-322 - Revision 3.5},\n"
1709371c9d4SSatish Balay                                    "url    = {https://www.mcs.anl.gov/research/projects/tao/}\n}\n",
1719371c9d4SSatish Balay                                    &set));
172494bef23SAlp Dener   tao->header_printed = PETSC_FALSE;
1739566063dSJacob Faibussowitsch   PetscCall(TaoSetUp(tao));
1749566063dSJacob Faibussowitsch   PetscCall(TaoResetStatistics(tao));
1751baa6e33SBarry Smith   if (tao->linesearch) PetscCall(TaoLineSearchReset(tao->linesearch));
176a7e14dcfSSatish Balay 
1779566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(TAO_Solve, tao, 0, 0, 0));
178dbbe0bcdSBarry Smith   PetscTryTypeMethod(tao, solve);
1799566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(TAO_Solve, tao, 0, 0, 0));
180a7e14dcfSSatish Balay 
1819566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(tao->solution, (PetscObject)tao, "-tao_view_solution"));
1827cf37e64SBarry Smith 
1838931d482SJason Sarich   tao->ntotalits += tao->niter;
184a7e14dcfSSatish Balay 
185a7e14dcfSSatish Balay   if (tao->printreason) {
186f642d338SBarry Smith     PetscViewer viewer = PETSC_VIEWER_STDOUT_(((PetscObject)tao)->comm);
187f642d338SBarry Smith     PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)tao)->tablevel));
188a7e14dcfSSatish Balay     if (tao->reason > 0) {
189f642d338SBarry Smith       PetscCall(PetscViewerASCIIPrintf(viewer, "  TAO %s solve converged due to %s iterations %" PetscInt_FMT "\n", ((PetscObject)tao)->prefix ? ((PetscObject)tao)->prefix : "", TaoConvergedReasons[tao->reason], tao->niter));
190a7e14dcfSSatish Balay     } else {
191f642d338SBarry Smith       PetscCall(PetscViewerASCIIPrintf(viewer, "  TAO %s solve did not converge due to %s iteration %" PetscInt_FMT "\n", ((PetscObject)tao)->prefix ? ((PetscObject)tao)->prefix : "", TaoConvergedReasons[tao->reason], tao->niter));
192a7e14dcfSSatish Balay     }
193f642d338SBarry Smith     PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)tao)->tablevel));
194a7e14dcfSSatish Balay   }
195f642d338SBarry Smith   PetscCall(TaoViewFromOptions(tao, NULL, "-tao_view"));
1963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
197a7e14dcfSSatish Balay }
198a7e14dcfSSatish Balay 
199a7e14dcfSSatish Balay /*@
200a7e14dcfSSatish Balay   TaoSetUp - Sets up the internal data structures for the later use
201a7e14dcfSSatish Balay   of a Tao solver
202a7e14dcfSSatish Balay 
203c3339decSBarry Smith   Collective
204a7e14dcfSSatish Balay 
20547450a7bSBarry Smith   Input Parameter:
20647450a7bSBarry Smith . tao - the `Tao` context
207a7e14dcfSSatish Balay 
20867be906fSBarry Smith   Level: advanced
20967be906fSBarry Smith 
21047450a7bSBarry Smith   Note:
21165ba42b6SBarry Smith   The user will not need to explicitly call `TaoSetUp()`, as it will
21265ba42b6SBarry Smith   automatically be called in `TaoSolve()`.  However, if the user
21365ba42b6SBarry Smith   desires to call it explicitly, it should come after `TaoCreate()`
21465ba42b6SBarry Smith   and any TaoSetSomething() routines, but before `TaoSolve()`.
215a7e14dcfSSatish Balay 
2161cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoCreate()`, `TaoSolve()`
217a7e14dcfSSatish Balay @*/
218d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoSetUp(Tao tao)
219d71ae5a4SJacob Faibussowitsch {
220a7e14dcfSSatish Balay   PetscFunctionBegin;
221441846f8SBarry Smith   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
2223ba16761SJacob Faibussowitsch   if (tao->setupcalled) PetscFunctionReturn(PETSC_SUCCESS);
2230f0abf79SStefano Zampini   PetscCall(TaoSetUpEW_Private(tao));
2243c859ba3SBarry Smith   PetscCheck(tao->solution, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "Must call TaoSetSolution");
225dbbe0bcdSBarry Smith   PetscTryTypeMethod(tao, setup);
226a7e14dcfSSatish Balay   tao->setupcalled = PETSC_TRUE;
2273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
228a7e14dcfSSatish Balay }
229a7e14dcfSSatish Balay 
2300764c050SBarry Smith /*@
23147450a7bSBarry Smith   TaoDestroy - Destroys the `Tao` context that was created with `TaoCreate()`
232a7e14dcfSSatish Balay 
233c3339decSBarry Smith   Collective
234a7e14dcfSSatish Balay 
235a7e14dcfSSatish Balay   Input Parameter:
23647450a7bSBarry Smith . tao - the `Tao` context
237a7e14dcfSSatish Balay 
238a7e14dcfSSatish Balay   Level: beginner
239a7e14dcfSSatish Balay 
2401cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoCreate()`, `TaoSolve()`
241a7e14dcfSSatish Balay @*/
242d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoDestroy(Tao *tao)
243d71ae5a4SJacob Faibussowitsch {
244a7e14dcfSSatish Balay   PetscFunctionBegin;
2453ba16761SJacob Faibussowitsch   if (!*tao) PetscFunctionReturn(PETSC_SUCCESS);
246441846f8SBarry Smith   PetscValidHeaderSpecific(*tao, TAO_CLASSID, 1);
2479371c9d4SSatish Balay   if (--((PetscObject)*tao)->refct > 0) {
2489371c9d4SSatish Balay     *tao = NULL;
2493ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
250a7e14dcfSSatish Balay   }
2519371c9d4SSatish Balay 
252213acdd3SPierre Jolivet   PetscTryTypeMethod(*tao, destroy);
2539566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&(*tao)->ksp));
2540f0abf79SStefano Zampini   PetscCall(SNESDestroy(&(*tao)->snes_ewdummy));
2559566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchDestroy(&(*tao)->linesearch));
256a7e14dcfSSatish Balay 
257a7e14dcfSSatish Balay   if ((*tao)->ops->convergencedestroy) {
2589566063dSJacob Faibussowitsch     PetscCall((*(*tao)->ops->convergencedestroy)((*tao)->cnvP));
25948a46eb9SPierre Jolivet     if ((*tao)->jacobian_state_inv) PetscCall(MatDestroy(&(*tao)->jacobian_state_inv));
260a7e14dcfSSatish Balay   }
2619566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&(*tao)->solution));
2629566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&(*tao)->gradient));
2639566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&(*tao)->ls_res));
264a7e14dcfSSatish Balay 
265a9603a14SPatrick Farrell   if ((*tao)->gradient_norm) {
2669566063dSJacob Faibussowitsch     PetscCall(PetscObjectDereference((PetscObject)(*tao)->gradient_norm));
2679566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&(*tao)->gradient_norm_tmp));
268a9603a14SPatrick Farrell   }
269a9603a14SPatrick Farrell 
2709566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&(*tao)->XL));
2719566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&(*tao)->XU));
2729566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&(*tao)->IL));
2739566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&(*tao)->IU));
2749566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&(*tao)->DE));
2759566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&(*tao)->DI));
2769566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&(*tao)->constraints));
2779566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&(*tao)->constraints_equality));
2789566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&(*tao)->constraints_inequality));
2799566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&(*tao)->stepdirection));
2809566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&(*tao)->hessian_pre));
2819566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&(*tao)->hessian));
2829566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&(*tao)->ls_jac));
2839566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&(*tao)->ls_jac_pre));
2849566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&(*tao)->jacobian_pre));
2859566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&(*tao)->jacobian));
2869566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&(*tao)->jacobian_state_pre));
2879566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&(*tao)->jacobian_state));
2889566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&(*tao)->jacobian_state_inv));
2899566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&(*tao)->jacobian_design));
2909566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&(*tao)->jacobian_equality));
2919566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&(*tao)->jacobian_equality_pre));
2929566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&(*tao)->jacobian_inequality));
2939566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&(*tao)->jacobian_inequality_pre));
2949566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&(*tao)->state_is));
2959566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&(*tao)->design_is));
2969566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&(*tao)->res_weights_v));
297b2dc45ddSPierre Jolivet   PetscCall(TaoMonitorCancel(*tao));
29848a46eb9SPierre Jolivet   if ((*tao)->hist_malloc) PetscCall(PetscFree4((*tao)->hist_obj, (*tao)->hist_resid, (*tao)->hist_cnorm, (*tao)->hist_lits));
299737f463aSAlp Dener   if ((*tao)->res_weights_n) {
3009566063dSJacob Faibussowitsch     PetscCall(PetscFree((*tao)->res_weights_rows));
3019566063dSJacob Faibussowitsch     PetscCall(PetscFree((*tao)->res_weights_cols));
3029566063dSJacob Faibussowitsch     PetscCall(PetscFree((*tao)->res_weights_w));
303125ddc8aSJason Sarich   }
3049566063dSJacob Faibussowitsch   PetscCall(PetscHeaderDestroy(tao));
3053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
306a7e14dcfSSatish Balay }
307a7e14dcfSSatish Balay 
308a7e14dcfSSatish Balay /*@
3091d27aa22SBarry Smith   TaoKSPSetUseEW - Sets `SNES` to use Eisenstat-Walker method {cite}`ew96`for computing relative tolerance for linear solvers.
3100f0abf79SStefano Zampini 
311c3339decSBarry Smith   Logically Collective
3120f0abf79SStefano Zampini 
3130f0abf79SStefano Zampini   Input Parameters:
3140f0abf79SStefano Zampini + tao  - Tao context
31565ba42b6SBarry Smith - flag - `PETSC_TRUE` or `PETSC_FALSE`
3160f0abf79SStefano Zampini 
31767be906fSBarry Smith   Level: advanced
31867be906fSBarry Smith 
31947450a7bSBarry Smith   Note:
32065ba42b6SBarry Smith   See `SNESKSPSetUseEW()` for customization details.
3210f0abf79SStefano Zampini 
3221cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `SNESKSPSetUseEW()`
3230f0abf79SStefano Zampini @*/
324d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoKSPSetUseEW(Tao tao, PetscBool flag)
325d71ae5a4SJacob Faibussowitsch {
3260f0abf79SStefano Zampini   PetscFunctionBegin;
3270f0abf79SStefano Zampini   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
3280f0abf79SStefano Zampini   PetscValidLogicalCollectiveBool(tao, flag, 2);
3290f0abf79SStefano Zampini   tao->ksp_ewconv = flag;
3303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3310f0abf79SStefano Zampini }
3320f0abf79SStefano Zampini 
3330f0abf79SStefano Zampini /*@
33447450a7bSBarry Smith   TaoSetFromOptions - Sets various Tao parameters from the options database
335a7e14dcfSSatish Balay 
336c3339decSBarry Smith   Collective
337a7e14dcfSSatish Balay 
33801d2d390SJose E. Roman   Input Parameter:
33947450a7bSBarry Smith . tao - the `Tao` solver context
340a7e14dcfSSatish Balay 
34147450a7bSBarry Smith   Options Database Keys:
34265ba42b6SBarry Smith + -tao_type <type>             - The algorithm that Tao uses (lmvm, nls, etc.)
343a7e14dcfSSatish Balay . -tao_gatol <gatol>           - absolute error tolerance for ||gradient||
344a7e14dcfSSatish Balay . -tao_grtol <grtol>           - relative error tolerance for ||gradient||
345a7e14dcfSSatish Balay . -tao_gttol <gttol>           - reduction of ||gradient|| relative to initial gradient
346a7e14dcfSSatish Balay . -tao_max_it <max>            - sets maximum number of iterations
347a7e14dcfSSatish Balay . -tao_max_funcs <max>         - sets maximum number of function evaluations
348a7e14dcfSSatish Balay . -tao_fmin <fmin>             - stop if function value reaches fmin
349a7e14dcfSSatish Balay . -tao_steptol <tol>           - stop if trust region radius less than <tol>
350a7e14dcfSSatish Balay . -tao_trust0 <t>              - initial trust region radius
35110978b7dSBarry Smith . -tao_view_solution           - view the solution at the end of the optimization process
35247450a7bSBarry Smith . -tao_monitor                 - prints function value and residual norm at each iteration
35310978b7dSBarry Smith . -tao_monitor_short           - same as `-tao_monitor`, but truncates very small values
35410978b7dSBarry Smith . -tao_monitor_constraint_norm - prints objective value, gradient, and constraint norm at each iteration
35510978b7dSBarry Smith . -tao_monitor_globalization   - prints information about the globalization at each iteration
35610978b7dSBarry Smith . -tao_monitor_solution        - prints solution vector at each iteration
35710978b7dSBarry Smith . -tao_monitor_ls_residual     - prints least-squares residual vector at each iteration
35810978b7dSBarry Smith . -tao_monitor_step            - prints step vector at each iteration
35910978b7dSBarry Smith . -tao_monitor_gradient        - prints gradient vector at each iteration
36010978b7dSBarry Smith . -tao_monitor_solution_draw   - graphically view solution vector at each iteration
36110978b7dSBarry Smith . -tao_monitor_step_draw       - graphically view step vector at each iteration
36210978b7dSBarry Smith . -tao_monitor_gradient_draw   - graphically view gradient at each iteration
36310978b7dSBarry Smith . -tao_monitor_cancel          - cancels all monitors (except those set with command line)
364a7e14dcfSSatish Balay . -tao_fd_gradient             - use gradient computed with finite differences
365f4c1ad5cSStefano Zampini . -tao_fd_hessian              - use hessian computed with finite differences
36610978b7dSBarry Smith . -tao_mf_hessian              - use matrix-free Hessian computed with finite differences
367441846f8SBarry Smith . -tao_view                    - prints information about the Tao after solving
36865ba42b6SBarry Smith - -tao_converged_reason        - prints the reason Tao stopped iterating
369a7e14dcfSSatish Balay 
37067be906fSBarry Smith   Level: beginner
37167be906fSBarry Smith 
37267be906fSBarry Smith   Note:
37347450a7bSBarry Smith   To see all options, run your program with the `-help` option or consult the
37465ba42b6SBarry Smith   user's manual. Should be called after `TaoCreate()` but before `TaoSolve()`
375a7e14dcfSSatish Balay 
3761cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoCreate()`, `TaoSolve()`
377a7e14dcfSSatish Balay @*/
378d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoSetFromOptions(Tao tao)
379d71ae5a4SJacob Faibussowitsch {
380b625d6c7SJed Brown   TaoType     default_type = TAOLMVM;
381a7e14dcfSSatish Balay   char        type[256], monfilename[PETSC_MAX_PATH_LEN];
382a7e14dcfSSatish Balay   PetscViewer monviewer;
383e876fe75SStephan Köhler   PetscBool   flg, found;
384a7e14dcfSSatish Balay   MPI_Comm    comm;
385606f75f6SBarry Smith   PetscReal   catol, crtol, gatol, grtol, gttol;
386a7e14dcfSSatish Balay 
387a7e14dcfSSatish Balay   PetscFunctionBegin;
388441846f8SBarry Smith   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
3899566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)tao, &comm));
390125ddc8aSJason Sarich 
39160b70c5cSStefano Zampini   if (((PetscObject)tao)->type_name) default_type = ((PetscObject)tao)->type_name;
39260b70c5cSStefano Zampini 
393d0609cedSBarry Smith   PetscObjectOptionsBegin((PetscObject)tao);
394a7e14dcfSSatish Balay   /* Check for type from options */
3959566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFList("-tao_type", "Tao Solver type", "TaoSetType", TaoList, default_type, type, 256, &flg));
396a7e14dcfSSatish Balay   if (flg) {
3979566063dSJacob Faibussowitsch     PetscCall(TaoSetType(tao, type));
398a7e14dcfSSatish Balay   } else if (!((PetscObject)tao)->type_name) {
3999566063dSJacob Faibussowitsch     PetscCall(TaoSetType(tao, default_type));
400a7e14dcfSSatish Balay   }
401a7e14dcfSSatish Balay 
40260b70c5cSStefano Zampini   /* Tao solvers do not set the prefix, set it here if not yet done
40360b70c5cSStefano Zampini      We do it after SetType since solver may have been changed */
40460b70c5cSStefano Zampini   if (tao->linesearch) {
40560b70c5cSStefano Zampini     const char *prefix;
40660b70c5cSStefano Zampini     PetscCall(TaoLineSearchGetOptionsPrefix(tao->linesearch, &prefix));
407f4f49eeaSPierre Jolivet     if (!prefix) PetscCall(TaoLineSearchSetOptionsPrefix(tao->linesearch, ((PetscObject)tao)->prefix));
40860b70c5cSStefano Zampini   }
40960b70c5cSStefano Zampini 
410606f75f6SBarry Smith   catol = tao->catol;
411606f75f6SBarry Smith   crtol = tao->crtol;
412606f75f6SBarry Smith   PetscCall(PetscOptionsReal("-tao_catol", "Stop if constraints violations within", "TaoSetConstraintTolerances", tao->catol, &catol, NULL));
413606f75f6SBarry Smith   PetscCall(PetscOptionsReal("-tao_crtol", "Stop if relative constraint violations within", "TaoSetConstraintTolerances", tao->crtol, &crtol, NULL));
414606f75f6SBarry Smith   PetscCall(TaoSetConstraintTolerances(tao, catol, crtol));
415606f75f6SBarry Smith 
416606f75f6SBarry Smith   gatol = tao->gatol;
417606f75f6SBarry Smith   grtol = tao->grtol;
418606f75f6SBarry Smith   gttol = tao->gttol;
419606f75f6SBarry Smith   PetscCall(PetscOptionsReal("-tao_gatol", "Stop if norm of gradient less than", "TaoSetTolerances", tao->gatol, &gatol, NULL));
420606f75f6SBarry Smith   PetscCall(PetscOptionsReal("-tao_grtol", "Stop if norm of gradient divided by the function value is less than", "TaoSetTolerances", tao->grtol, &grtol, NULL));
421606f75f6SBarry Smith   PetscCall(PetscOptionsReal("-tao_gttol", "Stop if the norm of the gradient is less than the norm of the initial gradient times tol", "TaoSetTolerances", tao->gttol, &gttol, NULL));
422606f75f6SBarry Smith   PetscCall(TaoSetTolerances(tao, gatol, grtol, gttol));
423606f75f6SBarry Smith 
4249566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-tao_max_it", "Stop if iteration number exceeds", "TaoSetMaximumIterations", tao->max_it, &tao->max_it, &flg));
425606f75f6SBarry Smith   if (flg) PetscCall(TaoSetMaximumIterations(tao, tao->max_it));
426606f75f6SBarry Smith 
4279566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-tao_max_funcs", "Stop if number of function evaluations exceeds", "TaoSetMaximumFunctionEvaluations", tao->max_funcs, &tao->max_funcs, &flg));
428606f75f6SBarry Smith   if (flg) PetscCall(TaoSetMaximumFunctionEvaluations(tao, tao->max_funcs));
429606f75f6SBarry Smith 
430606f75f6SBarry Smith   PetscCall(PetscOptionsReal("-tao_fmin", "Stop if function less than", "TaoSetFunctionLowerBound", tao->fmin, &tao->fmin, NULL));
431606f75f6SBarry Smith   PetscCall(PetscOptionsBoundedReal("-tao_steptol", "Stop if step size or trust region radius less than", "", tao->steptol, &tao->steptol, NULL, 0));
432606f75f6SBarry Smith   PetscCall(PetscOptionsReal("-tao_trust0", "Initial trust region radius", "TaoSetInitialTrustRegionRadius", tao->trust0, &tao->trust0, &flg));
433606f75f6SBarry Smith   if (flg) PetscCall(TaoSetInitialTrustRegionRadius(tao, tao->trust0));
43410978b7dSBarry Smith 
43510978b7dSBarry Smith   PetscCall(PetscOptionsDeprecated("-tao_solution_monitor", "-tao_monitor_solution", "3.21", NULL));
43610978b7dSBarry Smith   PetscCall(PetscOptionsDeprecated("-tao_gradient_monitor", "-tao_monitor_gradient", "3.21", NULL));
43710978b7dSBarry Smith   PetscCall(PetscOptionsDeprecated("-tao_stepdirection_monitor", "-tao_monitor_step", "3.21", NULL));
43810978b7dSBarry Smith   PetscCall(PetscOptionsDeprecated("-tao_residual_monitor", "-tao_monitor_residual", "3.21", NULL));
43910978b7dSBarry Smith   PetscCall(PetscOptionsDeprecated("-tao_smonitor", "-tao_monitor_short", "3.21", NULL));
44010978b7dSBarry Smith   PetscCall(PetscOptionsDeprecated("-tao_cmonitor", "-tao_monitor_constraint_norm", "3.21", NULL));
44110978b7dSBarry Smith   PetscCall(PetscOptionsDeprecated("-tao_gmonitor", "-tao_monitor_globalization", "3.21", NULL));
44210978b7dSBarry Smith   PetscCall(PetscOptionsDeprecated("-tao_draw_solution", "-tao_monitor_solution_draw", "3.21", NULL));
44310978b7dSBarry Smith   PetscCall(PetscOptionsDeprecated("-tao_draw_gradient", "-tao_monitor_gradient_draw", "3.21", NULL));
44410978b7dSBarry Smith   PetscCall(PetscOptionsDeprecated("-tao_draw_step", "-tao_monitor_step_draw", "3.21", NULL));
44510978b7dSBarry Smith 
44610978b7dSBarry Smith   PetscCall(PetscOptionsString("-tao_monitor_solution", "View solution vector after each iteration", "TaoMonitorSet", "stdout", monfilename, sizeof(monfilename), &flg));
447a7e14dcfSSatish Balay   if (flg) {
4489566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIOpen(comm, monfilename, &monviewer));
44910978b7dSBarry Smith     PetscCall(TaoMonitorSet(tao, TaoMonitorSolution, monviewer, (PetscErrorCode(*)(void **))PetscViewerDestroy));
450a7e14dcfSSatish Balay   }
451a7e14dcfSSatish Balay 
45265ba42b6SBarry Smith   PetscCall(PetscOptionsBool("-tao_converged_reason", "Print reason for Tao converged", "TaoSolve", tao->printreason, &tao->printreason, NULL));
45310978b7dSBarry Smith   PetscCall(PetscOptionsString("-tao_monitor_gradient", "View gradient vector for each iteration", "TaoMonitorSet", "stdout", monfilename, sizeof(monfilename), &flg));
454a7e14dcfSSatish Balay   if (flg) {
4559566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIOpen(comm, monfilename, &monviewer));
45610978b7dSBarry Smith     PetscCall(TaoMonitorSet(tao, TaoMonitorGradient, monviewer, (PetscErrorCode(*)(void **))PetscViewerDestroy));
457a7e14dcfSSatish Balay   }
458a7e14dcfSSatish Balay 
45910978b7dSBarry Smith   PetscCall(PetscOptionsString("-tao_monitor_step", "View step vector after each iteration", "TaoMonitorSet", "stdout", monfilename, sizeof(monfilename), &flg));
460a7e14dcfSSatish Balay   if (flg) {
4619566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIOpen(comm, monfilename, &monviewer));
46210978b7dSBarry Smith     PetscCall(TaoMonitorSet(tao, TaoMonitorStep, monviewer, (PetscErrorCode(*)(void **))PetscViewerDestroy));
463a7e14dcfSSatish Balay   }
464a7e14dcfSSatish Balay 
46510978b7dSBarry Smith   PetscCall(PetscOptionsString("-tao_monitor_residual", "View least-squares residual vector after each iteration", "TaoMonitorSet", "stdout", monfilename, sizeof(monfilename), &flg));
466a7e14dcfSSatish Balay   if (flg) {
4679566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIOpen(comm, monfilename, &monviewer));
46810978b7dSBarry Smith     PetscCall(TaoMonitorSet(tao, TaoMonitorResidual, monviewer, (PetscErrorCode(*)(void **))PetscViewerDestroy));
469a7e14dcfSSatish Balay   }
470a7e14dcfSSatish Balay 
47110978b7dSBarry Smith   PetscCall(PetscOptionsString("-tao_monitor", "Use the default convergence monitor", "TaoMonitorSet", "stdout", monfilename, sizeof(monfilename), &flg));
472a7e14dcfSSatish Balay   if (flg) {
4739566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIOpen(comm, monfilename, &monviewer));
47410978b7dSBarry Smith     PetscCall(TaoMonitorSet(tao, TaoMonitorDefault, monviewer, (PetscErrorCode(*)(void **))PetscViewerDestroy));
475a7e14dcfSSatish Balay   }
476a7e14dcfSSatish Balay 
47710978b7dSBarry Smith   PetscCall(PetscOptionsString("-tao_monitor_globalization", "Use the convergence monitor with extra globalization info", "TaoMonitorSet", "stdout", monfilename, sizeof(monfilename), &flg));
4788d5ead36SAlp Dener   if (flg) {
4799566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIOpen(comm, monfilename, &monviewer));
48010978b7dSBarry Smith     PetscCall(TaoMonitorSet(tao, TaoMonitorGlobalization, monviewer, (PetscErrorCode(*)(void **))PetscViewerDestroy));
4818d5ead36SAlp Dener   }
4828d5ead36SAlp Dener 
48310978b7dSBarry Smith   PetscCall(PetscOptionsString("-tao_monitor_short", "Use the short convergence monitor", "TaoMonitorSet", "stdout", monfilename, sizeof(monfilename), &flg));
484a7e14dcfSSatish Balay   if (flg) {
4859566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIOpen(comm, monfilename, &monviewer));
48610978b7dSBarry Smith     PetscCall(TaoMonitorSet(tao, TaoMonitorDefaultShort, monviewer, (PetscErrorCode(*)(void **))PetscViewerDestroy));
487a7e14dcfSSatish Balay   }
488a7e14dcfSSatish Balay 
48910978b7dSBarry Smith   PetscCall(PetscOptionsString("-tao_monitor_constraint_norm", "Use the default convergence monitor with constraint norm", "TaoMonitorSet", "stdout", monfilename, sizeof(monfilename), &flg));
490a7e14dcfSSatish Balay   if (flg) {
4919566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIOpen(comm, monfilename, &monviewer));
49210978b7dSBarry Smith     PetscCall(TaoMonitorSet(tao, TaoMonitorConstraintNorm, monviewer, (PetscErrorCode(*)(void **))PetscViewerDestroy));
493a7e14dcfSSatish Balay   }
494a7e14dcfSSatish Balay 
4958afaa268SBarry Smith   flg = PETSC_FALSE;
496b2dc45ddSPierre Jolivet   PetscCall(PetscOptionsDeprecated("-tao_cancelmonitors", "-tao_monitor_cancel", "3.21", NULL));
497b2dc45ddSPierre Jolivet   PetscCall(PetscOptionsBool("-tao_monitor_cancel", "cancel all monitors and call any registered destroy routines", "TaoMonitorCancel", flg, &flg, NULL));
498b2dc45ddSPierre Jolivet   if (flg) PetscCall(TaoMonitorCancel(tao));
499a7e14dcfSSatish Balay 
5008afaa268SBarry Smith   flg = PETSC_FALSE;
50110978b7dSBarry Smith   PetscCall(PetscOptionsBool("-tao_monitor_solution_draw", "Plot solution vector at each iteration", "TaoMonitorSet", flg, &flg, NULL));
502a7e14dcfSSatish Balay   if (flg) {
503e882e171SHong Zhang     TaoMonitorDrawCtx drawctx;
504e882e171SHong Zhang     PetscInt          howoften = 1;
5059566063dSJacob Faibussowitsch     PetscCall(TaoMonitorDrawCtxCreate(PetscObjectComm((PetscObject)tao), NULL, NULL, PETSC_DECIDE, PETSC_DECIDE, 300, 300, howoften, &drawctx));
50610978b7dSBarry Smith     PetscCall(TaoMonitorSet(tao, TaoMonitorSolutionDraw, drawctx, (PetscErrorCode(*)(void **))TaoMonitorDrawCtxDestroy));
507a7e14dcfSSatish Balay   }
508a7e14dcfSSatish Balay 
5098afaa268SBarry Smith   flg = PETSC_FALSE;
51010978b7dSBarry Smith   PetscCall(PetscOptionsBool("-tao_monitor_step_draw", "Plots step at each iteration", "TaoMonitorSet", flg, &flg, NULL));
51110978b7dSBarry Smith   if (flg) PetscCall(TaoMonitorSet(tao, TaoMonitorStepDraw, NULL, NULL));
512a7e14dcfSSatish Balay 
5138afaa268SBarry Smith   flg = PETSC_FALSE;
51410978b7dSBarry Smith   PetscCall(PetscOptionsBool("-tao_monitor_gradient_draw", "plots gradient at each iteration", "TaoMonitorSet", flg, &flg, NULL));
515a7e14dcfSSatish Balay   if (flg) {
516e882e171SHong Zhang     TaoMonitorDrawCtx drawctx;
517e882e171SHong Zhang     PetscInt          howoften = 1;
5189566063dSJacob Faibussowitsch     PetscCall(TaoMonitorDrawCtxCreate(PetscObjectComm((PetscObject)tao), NULL, NULL, PETSC_DECIDE, PETSC_DECIDE, 300, 300, howoften, &drawctx));
51910978b7dSBarry Smith     PetscCall(TaoMonitorSet(tao, TaoMonitorGradientDraw, drawctx, (PetscErrorCode(*)(void **))TaoMonitorDrawCtxDestroy));
520a7e14dcfSSatish Balay   }
5218afaa268SBarry Smith   flg = PETSC_FALSE;
5229566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-tao_fd_gradient", "compute gradient using finite differences", "TaoDefaultComputeGradient", flg, &flg, NULL));
5231baa6e33SBarry Smith   if (flg) PetscCall(TaoSetGradient(tao, NULL, TaoDefaultComputeGradient, NULL));
524f4c1ad5cSStefano Zampini   flg = PETSC_FALSE;
52510978b7dSBarry Smith   PetscCall(PetscOptionsBool("-tao_fd_hessian", "compute Hessian using finite differences", "TaoDefaultComputeHessian", flg, &flg, NULL));
526f4c1ad5cSStefano Zampini   if (flg) {
527f4c1ad5cSStefano Zampini     Mat H;
528f4c1ad5cSStefano Zampini 
5299566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)tao), &H));
5309566063dSJacob Faibussowitsch     PetscCall(MatSetType(H, MATAIJ));
5319566063dSJacob Faibussowitsch     PetscCall(TaoSetHessian(tao, H, H, TaoDefaultComputeHessian, NULL));
5329566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&H));
533f4c1ad5cSStefano Zampini   }
534f4c1ad5cSStefano Zampini   flg = PETSC_FALSE;
53510978b7dSBarry Smith   PetscCall(PetscOptionsBool("-tao_mf_hessian", "compute matrix-free Hessian using finite differences", "TaoDefaultComputeHessianMFFD", flg, &flg, NULL));
536f4c1ad5cSStefano Zampini   if (flg) {
537f4c1ad5cSStefano Zampini     Mat H;
538f4c1ad5cSStefano Zampini 
5399566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)tao), &H));
5409566063dSJacob Faibussowitsch     PetscCall(TaoSetHessian(tao, H, H, TaoDefaultComputeHessianMFFD, NULL));
5419566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&H));
542f4c1ad5cSStefano Zampini   }
543e876fe75SStephan Köhler   PetscCall(PetscOptionsBool("-tao_recycle_history", "enable recycling/re-using information from the previous TaoSolve() call for some algorithms", "TaoSetRecycleHistory", flg, &flg, &found));
544e876fe75SStephan Köhler   if (found) PetscCall(TaoSetRecycleHistory(tao, flg));
5459566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-tao_subset_type", "subset type", "", TaoSubSetTypes, (PetscEnum)tao->subset_type, (PetscEnum *)&tao->subset_type, NULL));
546a7e14dcfSSatish Balay 
5470f0abf79SStefano Zampini   if (tao->ksp) {
5480f0abf79SStefano Zampini     PetscCall(PetscOptionsBool("-tao_ksp_ew", "Use Eisentat-Walker linear system convergence test", "TaoKSPSetUseEW", tao->ksp_ewconv, &tao->ksp_ewconv, NULL));
5490f0abf79SStefano Zampini     PetscCall(TaoKSPSetUseEW(tao, tao->ksp_ewconv));
5500f0abf79SStefano Zampini   }
5510f0abf79SStefano Zampini 
552dbbe0bcdSBarry Smith   PetscTryTypeMethod(tao, setfromoptions, PetscOptionsObject);
55360b70c5cSStefano Zampini 
55460b70c5cSStefano Zampini   /* process any options handlers added with PetscObjectAddOptionsHandler() */
55560b70c5cSStefano Zampini   PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)tao, PetscOptionsObject));
556d0609cedSBarry Smith   PetscOptionsEnd();
55760b70c5cSStefano Zampini 
55860b70c5cSStefano Zampini   if (tao->linesearch) PetscCall(TaoLineSearchSetFromOptions(tao->linesearch));
5593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
560a7e14dcfSSatish Balay }
561a7e14dcfSSatish Balay 
562ffeef943SBarry Smith /*@
56347450a7bSBarry Smith   TaoViewFromOptions - View a `Tao` object based on values in the options database
564fe2efc57SMark 
565c3339decSBarry Smith   Collective
566fe2efc57SMark 
567fe2efc57SMark   Input Parameters:
56847450a7bSBarry Smith + A    - the  `Tao` context
56947450a7bSBarry Smith . obj  - Optional object that provides the prefix for the options database
570736c3998SJose E. Roman - name - command line option
571fe2efc57SMark 
572fe2efc57SMark   Level: intermediate
57367be906fSBarry Smith 
5741cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoView`, `PetscObjectViewFromOptions()`, `TaoCreate()`
575fe2efc57SMark @*/
576d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoViewFromOptions(Tao A, PetscObject obj, const char name[])
577d71ae5a4SJacob Faibussowitsch {
578fe2efc57SMark   PetscFunctionBegin;
579fe2efc57SMark   PetscValidHeaderSpecific(A, TAO_CLASSID, 1);
5809566063dSJacob Faibussowitsch   PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name));
5813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
582fe2efc57SMark }
583fe2efc57SMark 
584ffeef943SBarry Smith /*@
58547450a7bSBarry Smith   TaoView - Prints information about the `Tao` object
586a7e14dcfSSatish Balay 
587c3339decSBarry Smith   Collective
588a7e14dcfSSatish Balay 
589a7e14dcfSSatish Balay   Input Parameters:
59047450a7bSBarry Smith + tao    - the `Tao` context
591a7e14dcfSSatish Balay - viewer - visualization context
592a7e14dcfSSatish Balay 
593a7e14dcfSSatish Balay   Options Database Key:
59465ba42b6SBarry Smith . -tao_view - Calls `TaoView()` at the end of `TaoSolve()`
595a7e14dcfSSatish Balay 
59667be906fSBarry Smith   Level: beginner
59767be906fSBarry Smith 
598a7e14dcfSSatish Balay   Notes:
599a7e14dcfSSatish Balay   The available visualization contexts include
60065ba42b6SBarry Smith +     `PETSC_VIEWER_STDOUT_SELF` - standard output (default)
60165ba42b6SBarry Smith -     `PETSC_VIEWER_STDOUT_WORLD` - synchronized standard
602a7e14dcfSSatish Balay   output where only the first processor opens
603a7e14dcfSSatish Balay   the file.  All other processors send their
604a7e14dcfSSatish Balay   data to the first processor to print.
605a7e14dcfSSatish Balay 
6061cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `PetscViewerASCIIOpen()`
607a7e14dcfSSatish Balay @*/
608d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoView(Tao tao, PetscViewer viewer)
609d71ae5a4SJacob Faibussowitsch {
610a7e14dcfSSatish Balay   PetscBool isascii, isstring;
611b625d6c7SJed Brown   TaoType   type;
61247a47007SBarry Smith 
613a7e14dcfSSatish Balay   PetscFunctionBegin;
614441846f8SBarry Smith   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
61548a46eb9SPierre Jolivet   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(((PetscObject)tao)->comm, &viewer));
616a7e14dcfSSatish Balay   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
617a7e14dcfSSatish Balay   PetscCheckSameComm(tao, 1, viewer, 2);
618a7e14dcfSSatish Balay 
6199566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
6209566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring));
621a7e14dcfSSatish Balay   if (isascii) {
6229566063dSJacob Faibussowitsch     PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)tao, viewer));
623a7e14dcfSSatish Balay 
6249566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
625606f75f6SBarry Smith     PetscTryTypeMethod(tao, view, viewer);
626606f75f6SBarry Smith     if (tao->linesearch) PetscCall(TaoLineSearchView(tao->linesearch, viewer));
627a7e14dcfSSatish Balay     if (tao->ksp) {
6289566063dSJacob Faibussowitsch       PetscCall(KSPView(tao->ksp, viewer));
62963a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "total KSP iterations: %" PetscInt_FMT "\n", tao->ksp_tot_its));
630a7e14dcfSSatish Balay     }
63119b5c101SAlp Dener 
63248a46eb9SPierre Jolivet     if (tao->XL || tao->XU) PetscCall(PetscViewerASCIIPrintf(viewer, "Active Set subset type: %s\n", TaoSubSetTypes[tao->subset_type]));
633a7e14dcfSSatish Balay 
6349566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "convergence tolerances: gatol=%g,", (double)tao->gatol));
635606f75f6SBarry Smith     PetscCall(PetscViewerASCIIPrintf(viewer, " grtol=%g,", (double)tao->grtol));
6369566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, " steptol=%g,", (double)tao->steptol));
6379566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, " gttol=%g\n", (double)tao->gttol));
6389566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Residual in Function/Gradient:=%g\n", (double)tao->residual));
639a7e14dcfSSatish Balay 
6406246e8cdSAlp Dener     if (tao->constrained) {
6419566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "convergence tolerances:"));
6429566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, " catol=%g,", (double)tao->catol));
6439566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, " crtol=%g\n", (double)tao->crtol));
6449566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "Residual in Constraints:=%g\n", (double)tao->cnorm));
645a7e14dcfSSatish Balay     }
646a7e14dcfSSatish Balay 
647a7e14dcfSSatish Balay     if (tao->trust < tao->steptol) {
6489566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "convergence tolerances: steptol=%g\n", (double)tao->steptol));
6499566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "Final trust region radius:=%g\n", (double)tao->trust));
650a7e14dcfSSatish Balay     }
651a7e14dcfSSatish Balay 
65248a46eb9SPierre Jolivet     if (tao->fmin > -1.e25) PetscCall(PetscViewerASCIIPrintf(viewer, "convergence tolerances: function minimum=%g\n", (double)tao->fmin));
6539566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Objective value=%g\n", (double)tao->fc));
654a7e14dcfSSatish Balay 
65563a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "total number of iterations=%" PetscInt_FMT ",          ", tao->niter));
65663a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "              (max: %" PetscInt_FMT ")\n", tao->max_it));
657a7e14dcfSSatish Balay 
658a7e14dcfSSatish Balay     if (tao->nfuncs > 0) {
65963a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "total number of function evaluations=%" PetscInt_FMT ",", tao->nfuncs));
660606f75f6SBarry Smith       if (tao->max_funcs == PETSC_UNLIMITED) PetscCall(PetscViewerASCIIPrintf(viewer, "                (max: unlimited)\n"));
661606f75f6SBarry Smith       else PetscCall(PetscViewerASCIIPrintf(viewer, "                (max: %" PetscInt_FMT ")\n", tao->max_funcs));
662a7e14dcfSSatish Balay     }
663a7e14dcfSSatish Balay     if (tao->ngrads > 0) {
66463a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "total number of gradient evaluations=%" PetscInt_FMT ",", tao->ngrads));
665606f75f6SBarry Smith       if (tao->max_funcs == PETSC_UNLIMITED) PetscCall(PetscViewerASCIIPrintf(viewer, "                (max: unlimited)\n"));
666606f75f6SBarry Smith       else PetscCall(PetscViewerASCIIPrintf(viewer, "                (max: %" PetscInt_FMT ")\n", tao->max_funcs));
667a7e14dcfSSatish Balay     }
668a7e14dcfSSatish Balay     if (tao->nfuncgrads > 0) {
66963a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "total number of function/gradient evaluations=%" PetscInt_FMT ",", tao->nfuncgrads));
670606f75f6SBarry Smith       if (tao->max_funcs == PETSC_UNLIMITED) PetscCall(PetscViewerASCIIPrintf(viewer, "    (max: unlimited)\n"));
671606f75f6SBarry Smith       else PetscCall(PetscViewerASCIIPrintf(viewer, "    (max: %" PetscInt_FMT ")\n", tao->max_funcs));
672a7e14dcfSSatish Balay     }
67348a46eb9SPierre Jolivet     if (tao->nhess > 0) PetscCall(PetscViewerASCIIPrintf(viewer, "total number of Hessian evaluations=%" PetscInt_FMT "\n", tao->nhess));
67448a46eb9SPierre Jolivet     if (tao->nconstraints > 0) PetscCall(PetscViewerASCIIPrintf(viewer, "total number of constraint function evaluations=%" PetscInt_FMT "\n", tao->nconstraints));
67548a46eb9SPierre Jolivet     if (tao->njac > 0) PetscCall(PetscViewerASCIIPrintf(viewer, "total number of Jacobian evaluations=%" PetscInt_FMT "\n", tao->njac));
676a7e14dcfSSatish Balay 
677a7e14dcfSSatish Balay     if (tao->reason > 0) {
6789566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "Solution converged: "));
679a7e14dcfSSatish Balay       switch (tao->reason) {
680d71ae5a4SJacob Faibussowitsch       case TAO_CONVERGED_GATOL:
681d71ae5a4SJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, " ||g(X)|| <= gatol\n"));
682d71ae5a4SJacob Faibussowitsch         break;
683d71ae5a4SJacob Faibussowitsch       case TAO_CONVERGED_GRTOL:
684d71ae5a4SJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, " ||g(X)||/|f(X)| <= grtol\n"));
685d71ae5a4SJacob Faibussowitsch         break;
686d71ae5a4SJacob Faibussowitsch       case TAO_CONVERGED_GTTOL:
687d71ae5a4SJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, " ||g(X)||/||g(X0)|| <= gttol\n"));
688d71ae5a4SJacob Faibussowitsch         break;
689d71ae5a4SJacob Faibussowitsch       case TAO_CONVERGED_STEPTOL:
690d71ae5a4SJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, " Steptol -- step size small\n"));
691d71ae5a4SJacob Faibussowitsch         break;
692d71ae5a4SJacob Faibussowitsch       case TAO_CONVERGED_MINF:
693d71ae5a4SJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, " Minf --  f < fmin\n"));
694d71ae5a4SJacob Faibussowitsch         break;
695d71ae5a4SJacob Faibussowitsch       case TAO_CONVERGED_USER:
696d71ae5a4SJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, " User Terminated\n"));
697d71ae5a4SJacob Faibussowitsch         break;
698d71ae5a4SJacob Faibussowitsch       default:
699606f75f6SBarry Smith         PetscCall(PetscViewerASCIIPrintf(viewer, " %d\n", tao->reason));
700d71ae5a4SJacob Faibussowitsch         break;
701a7e14dcfSSatish Balay       }
702606f75f6SBarry Smith     } else if (tao->reason == TAO_CONTINUE_ITERATING) {
703606f75f6SBarry Smith       PetscCall(PetscViewerASCIIPrintf(viewer, "Solver never run\n"));
704a7e14dcfSSatish Balay     } else {
705606f75f6SBarry Smith       PetscCall(PetscViewerASCIIPrintf(viewer, "Solver failed: "));
706a7e14dcfSSatish Balay       switch (tao->reason) {
707d71ae5a4SJacob Faibussowitsch       case TAO_DIVERGED_MAXITS:
708d71ae5a4SJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, " Maximum Iterations\n"));
709d71ae5a4SJacob Faibussowitsch         break;
710d71ae5a4SJacob Faibussowitsch       case TAO_DIVERGED_NAN:
711d71ae5a4SJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, " NAN or Inf encountered\n"));
712d71ae5a4SJacob Faibussowitsch         break;
713d71ae5a4SJacob Faibussowitsch       case TAO_DIVERGED_MAXFCN:
714d71ae5a4SJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, " Maximum Function Evaluations\n"));
715d71ae5a4SJacob Faibussowitsch         break;
716d71ae5a4SJacob Faibussowitsch       case TAO_DIVERGED_LS_FAILURE:
717d71ae5a4SJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, " Line Search Failure\n"));
718d71ae5a4SJacob Faibussowitsch         break;
719d71ae5a4SJacob Faibussowitsch       case TAO_DIVERGED_TR_REDUCTION:
720d71ae5a4SJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, " Trust Region too small\n"));
721d71ae5a4SJacob Faibussowitsch         break;
722d71ae5a4SJacob Faibussowitsch       case TAO_DIVERGED_USER:
723d71ae5a4SJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, " User Terminated\n"));
724d71ae5a4SJacob Faibussowitsch         break;
725d71ae5a4SJacob Faibussowitsch       default:
726606f75f6SBarry Smith         PetscCall(PetscViewerASCIIPrintf(viewer, " %d\n", tao->reason));
727d71ae5a4SJacob Faibussowitsch         break;
728a7e14dcfSSatish Balay       }
729a7e14dcfSSatish Balay     }
7309566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
731a7e14dcfSSatish Balay   } else if (isstring) {
7329566063dSJacob Faibussowitsch     PetscCall(TaoGetType(tao, &type));
7339566063dSJacob Faibussowitsch     PetscCall(PetscViewerStringSPrintf(viewer, " %-3.3s", type));
734a7e14dcfSSatish Balay   }
7353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
736a7e14dcfSSatish Balay }
737a7e14dcfSSatish Balay 
738a7e14dcfSSatish Balay /*@
739414d97d3SAlp Dener   TaoSetRecycleHistory - Sets the boolean flag to enable/disable re-using
74065ba42b6SBarry Smith   iterate information from the previous `TaoSolve()`. This feature is disabled by
741414d97d3SAlp Dener   default.
742414d97d3SAlp Dener 
74367be906fSBarry Smith   Logically Collective
74467be906fSBarry Smith 
74567be906fSBarry Smith   Input Parameters:
74647450a7bSBarry Smith + tao     - the `Tao` context
74767be906fSBarry Smith - recycle - boolean flag
74867be906fSBarry Smith 
74947450a7bSBarry Smith   Options Database Key:
75067be906fSBarry Smith . -tao_recycle_history <true,false> - reuse the history
75167be906fSBarry Smith 
75267be906fSBarry Smith   Level: intermediate
75367be906fSBarry Smith 
75467be906fSBarry Smith   Notes:
75565ba42b6SBarry Smith   For conjugate gradient methods (`TAOBNCG`), this re-uses the latest search direction
75665ba42b6SBarry Smith   from the previous `TaoSolve()` call when computing the first search direction in a
757414d97d3SAlp Dener   new solution. By default, CG methods set the first search direction to the
758414d97d3SAlp Dener   negative gradient.
759414d97d3SAlp Dener 
76065ba42b6SBarry Smith   For quasi-Newton family of methods (`TAOBQNLS`, `TAOBQNKLS`, `TAOBQNKTR`, `TAOBQNKTL`), this re-uses
76165ba42b6SBarry Smith   the accumulated quasi-Newton Hessian approximation from the previous `TaoSolve()`
762414d97d3SAlp Dener   call. By default, QN family of methods reset the initial Hessian approximation to
763414d97d3SAlp Dener   the identity matrix.
764414d97d3SAlp Dener 
765414d97d3SAlp Dener   For any other algorithm, this setting has no effect.
766414d97d3SAlp Dener 
7671cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoGetRecycleHistory()`, `TAOBNCG`, `TAOBQNLS`, `TAOBQNKLS`, `TAOBQNKTR`, `TAOBQNKTL`
768414d97d3SAlp Dener @*/
769d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoSetRecycleHistory(Tao tao, PetscBool recycle)
770d71ae5a4SJacob Faibussowitsch {
771414d97d3SAlp Dener   PetscFunctionBegin;
772414d97d3SAlp Dener   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
77394511df7SStefano Zampini   PetscValidLogicalCollectiveBool(tao, recycle, 2);
774414d97d3SAlp Dener   tao->recycle = recycle;
7753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
776414d97d3SAlp Dener }
777414d97d3SAlp Dener 
778414d97d3SAlp Dener /*@
779414d97d3SAlp Dener   TaoGetRecycleHistory - Retrieve the boolean flag for re-using iterate information
78065ba42b6SBarry Smith   from the previous `TaoSolve()`. This feature is disabled by default.
781414d97d3SAlp Dener 
78267be906fSBarry Smith   Logically Collective
783414d97d3SAlp Dener 
78447450a7bSBarry Smith   Input Parameter:
78547450a7bSBarry Smith . tao - the `Tao` context
786414d97d3SAlp Dener 
78747450a7bSBarry Smith   Output Parameter:
788147403d9SBarry Smith . recycle - boolean flag
789414d97d3SAlp Dener 
790414d97d3SAlp Dener   Level: intermediate
791414d97d3SAlp Dener 
7921cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoSetRecycleHistory()`, `TAOBNCG`, `TAOBQNLS`, `TAOBQNKLS`, `TAOBQNKTR`, `TAOBQNKTL`
793414d97d3SAlp Dener @*/
794d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoGetRecycleHistory(Tao tao, PetscBool *recycle)
795d71ae5a4SJacob Faibussowitsch {
796414d97d3SAlp Dener   PetscFunctionBegin;
797414d97d3SAlp Dener   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
7984f572ea9SToby Isaac   PetscAssertPointer(recycle, 2);
799414d97d3SAlp Dener   *recycle = tao->recycle;
8003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
801414d97d3SAlp Dener }
802414d97d3SAlp Dener 
803414d97d3SAlp Dener /*@
80447450a7bSBarry Smith   TaoSetTolerances - Sets parameters used in `TaoSolve()` convergence tests
805a7e14dcfSSatish Balay 
80667be906fSBarry Smith   Logically Collective
807a7e14dcfSSatish Balay 
808a7e14dcfSSatish Balay   Input Parameters:
80947450a7bSBarry Smith + tao   - the `Tao` context
810a7e14dcfSSatish Balay . gatol - stop if norm of gradient is less than this
811a7e14dcfSSatish Balay . grtol - stop if relative norm of gradient is less than this
812a7e14dcfSSatish Balay - gttol - stop if norm of gradient is reduced by this factor
813a7e14dcfSSatish Balay 
814a7e14dcfSSatish Balay   Options Database Keys:
815e52336cbSBarry Smith + -tao_gatol <gatol> - Sets gatol
816a7e14dcfSSatish Balay . -tao_grtol <grtol> - Sets grtol
817a7e14dcfSSatish Balay - -tao_gttol <gttol> - Sets gttol
818a7e14dcfSSatish Balay 
8193b242c63SJacob Faibussowitsch   Stopping Criteria\:
82067be906fSBarry Smith .vb
82167be906fSBarry Smith   ||g(X)||                            <= gatol
82267be906fSBarry Smith   ||g(X)|| / |f(X)|                   <= grtol
82367be906fSBarry Smith   ||g(X)|| / ||g(X0)||                <= gttol
82467be906fSBarry Smith .ve
825a7e14dcfSSatish Balay 
826a7e14dcfSSatish Balay   Level: beginner
827a7e14dcfSSatish Balay 
828606f75f6SBarry Smith   Notes:
829606f75f6SBarry Smith   Use `PETSC_CURRENT` to leave one or more tolerances unchanged.
830606f75f6SBarry Smith 
831606f75f6SBarry Smith   Use `PETSC_DETERMINE` to set one or more tolerances to their values when the `tao`object's type was set
832606f75f6SBarry Smith 
833606f75f6SBarry Smith   Fortran Note:
834606f75f6SBarry Smith   Use `PETSC_CURRENT_REAL` or `PETSC_DETERMINE_REAL`
835a7e14dcfSSatish Balay 
8361cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoConvergedReason`, `TaoGetTolerances()`
837a7e14dcfSSatish Balay @*/
838d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoSetTolerances(Tao tao, PetscReal gatol, PetscReal grtol, PetscReal gttol)
839d71ae5a4SJacob Faibussowitsch {
840a7e14dcfSSatish Balay   PetscFunctionBegin;
841441846f8SBarry Smith   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
84294511df7SStefano Zampini   PetscValidLogicalCollectiveReal(tao, gatol, 2);
84394511df7SStefano Zampini   PetscValidLogicalCollectiveReal(tao, grtol, 3);
84494511df7SStefano Zampini   PetscValidLogicalCollectiveReal(tao, gttol, 4);
845a7e14dcfSSatish Balay 
846606f75f6SBarry Smith   if (gatol == (PetscReal)PETSC_DETERMINE) {
847606f75f6SBarry Smith     tao->gatol = tao->default_gatol;
848606f75f6SBarry Smith   } else if (gatol != (PetscReal)PETSC_CURRENT) {
849606f75f6SBarry Smith     PetscCheck(gatol >= 0, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_OUTOFRANGE, "Negative gatol not allowed");
850606f75f6SBarry Smith     tao->gatol = gatol;
851a7e14dcfSSatish Balay   }
852a7e14dcfSSatish Balay 
853606f75f6SBarry Smith   if (grtol == (PetscReal)PETSC_DETERMINE) {
854606f75f6SBarry Smith     tao->grtol = tao->default_grtol;
855606f75f6SBarry Smith   } else if (grtol != (PetscReal)PETSC_CURRENT) {
856606f75f6SBarry Smith     PetscCheck(grtol >= 0, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_OUTOFRANGE, "Negative grtol not allowed");
857606f75f6SBarry Smith     tao->grtol = grtol;
858a7e14dcfSSatish Balay   }
859a7e14dcfSSatish Balay 
860606f75f6SBarry Smith   if (gttol == (PetscReal)PETSC_DETERMINE) {
861606f75f6SBarry Smith     tao->gttol = tao->default_gttol;
862606f75f6SBarry Smith   } else if (gttol != (PetscReal)PETSC_CURRENT) {
863606f75f6SBarry Smith     PetscCheck(gttol >= 0, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_OUTOFRANGE, "Negative gttol not allowed");
864606f75f6SBarry Smith     tao->gttol = gttol;
865a7e14dcfSSatish Balay   }
8663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
867a7e14dcfSSatish Balay }
868a7e14dcfSSatish Balay 
869a7e14dcfSSatish Balay /*@
87047450a7bSBarry Smith   TaoSetConstraintTolerances - Sets constraint tolerance parameters used in `TaoSolve()` convergence tests
871a7e14dcfSSatish Balay 
87267be906fSBarry Smith   Logically Collective
873a7e14dcfSSatish Balay 
874a7e14dcfSSatish Balay   Input Parameters:
87547450a7bSBarry Smith + tao   - the `Tao` context
876606f75f6SBarry Smith . catol - absolute constraint tolerance, constraint norm must be less than `catol` for used for `gatol` convergence criteria
877606f75f6SBarry Smith - crtol - relative constraint tolerance, constraint norm must be less than `crtol` for used for `gatol`, `gttol` convergence criteria
878a7e14dcfSSatish Balay 
879a7e14dcfSSatish Balay   Options Database Keys:
880a7e14dcfSSatish Balay + -tao_catol <catol> - Sets catol
881a7e14dcfSSatish Balay - -tao_crtol <crtol> - Sets crtol
882a7e14dcfSSatish Balay 
883a7e14dcfSSatish Balay   Level: intermediate
884a7e14dcfSSatish Balay 
88567be906fSBarry Smith   Notes:
886606f75f6SBarry Smith   Use `PETSC_CURRENT` to leave one or tolerance unchanged.
887606f75f6SBarry Smith 
888606f75f6SBarry Smith   Use `PETSC_DETERMINE` to set one or more tolerances to their values when the `tao` object's type was set
889606f75f6SBarry Smith 
890606f75f6SBarry Smith   Fortran Note:
891606f75f6SBarry Smith   Use `PETSC_CURRENT_REAL` or `PETSC_DETERMINE_REAL`
892a7e14dcfSSatish Balay 
8931cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoConvergedReason`, `TaoGetTolerances()`, `TaoGetConstraintTolerances()`, `TaoSetTolerances()`
894a7e14dcfSSatish Balay @*/
895d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoSetConstraintTolerances(Tao tao, PetscReal catol, PetscReal crtol)
896d71ae5a4SJacob Faibussowitsch {
897a7e14dcfSSatish Balay   PetscFunctionBegin;
898441846f8SBarry Smith   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
89994511df7SStefano Zampini   PetscValidLogicalCollectiveReal(tao, catol, 2);
90094511df7SStefano Zampini   PetscValidLogicalCollectiveReal(tao, crtol, 3);
901a7e14dcfSSatish Balay 
902606f75f6SBarry Smith   if (catol == (PetscReal)PETSC_DETERMINE) {
903606f75f6SBarry Smith     tao->catol = tao->default_catol;
904606f75f6SBarry Smith   } else if (catol != (PetscReal)PETSC_CURRENT) {
905606f75f6SBarry Smith     PetscCheck(catol >= 0, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_OUTOFRANGE, "Negative catol not allowed");
906606f75f6SBarry Smith     tao->catol = catol;
907a7e14dcfSSatish Balay   }
908a7e14dcfSSatish Balay 
909606f75f6SBarry Smith   if (crtol == (PetscReal)PETSC_DETERMINE) {
910606f75f6SBarry Smith     tao->crtol = tao->default_crtol;
911606f75f6SBarry Smith   } else if (crtol != (PetscReal)PETSC_CURRENT) {
912606f75f6SBarry Smith     PetscCheck(crtol >= 0, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_OUTOFRANGE, "Negative crtol not allowed");
913606f75f6SBarry Smith     tao->crtol = crtol;
914a7e14dcfSSatish Balay   }
9153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
916a7e14dcfSSatish Balay }
917a7e14dcfSSatish Balay 
91858417fe7SBarry Smith /*@
91947450a7bSBarry Smith   TaoGetConstraintTolerances - Gets constraint tolerance parameters used in `TaoSolve()` convergence tests
92058417fe7SBarry Smith 
92167be906fSBarry Smith   Not Collective
92258417fe7SBarry Smith 
92358417fe7SBarry Smith   Input Parameter:
92447450a7bSBarry Smith . tao - the `Tao` context
92558417fe7SBarry Smith 
926d8d19677SJose E. Roman   Output Parameters:
927606f75f6SBarry Smith + catol - absolute constraint tolerance, constraint norm must be less than `catol` for used for `gatol` convergence criteria
928606f75f6SBarry Smith - crtol - relative constraint tolerance, constraint norm must be less than `crtol` for used for `gatol`, `gttol` convergence criteria
92958417fe7SBarry Smith 
93058417fe7SBarry Smith   Level: intermediate
93158417fe7SBarry Smith 
9321cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoConvergedReasons`,`TaoGetTolerances()`, `TaoSetTolerances()`, `TaoSetConstraintTolerances()`
93358417fe7SBarry Smith @*/
934d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoGetConstraintTolerances(Tao tao, PetscReal *catol, PetscReal *crtol)
935d71ae5a4SJacob Faibussowitsch {
93658417fe7SBarry Smith   PetscFunctionBegin;
93758417fe7SBarry Smith   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
93858417fe7SBarry Smith   if (catol) *catol = tao->catol;
93958417fe7SBarry Smith   if (crtol) *crtol = tao->crtol;
9403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
94158417fe7SBarry Smith }
94258417fe7SBarry Smith 
943a7e14dcfSSatish Balay /*@
944a7e14dcfSSatish Balay   TaoSetFunctionLowerBound - Sets a bound on the solution objective value.
945a7e14dcfSSatish Balay   When an approximate solution with an objective value below this number
946a7e14dcfSSatish Balay   has been found, the solver will terminate.
947a7e14dcfSSatish Balay 
948c3339decSBarry Smith   Logically Collective
949a7e14dcfSSatish Balay 
950a7e14dcfSSatish Balay   Input Parameters:
951441846f8SBarry Smith + tao  - the Tao solver context
952a7e14dcfSSatish Balay - fmin - the tolerance
953a7e14dcfSSatish Balay 
95447450a7bSBarry Smith   Options Database Key:
955a7e14dcfSSatish Balay . -tao_fmin <fmin> - sets the minimum function value
956a7e14dcfSSatish Balay 
957a7e14dcfSSatish Balay   Level: intermediate
958a7e14dcfSSatish Balay 
9591cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoConvergedReason`, `TaoSetTolerances()`
960a7e14dcfSSatish Balay @*/
961d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoSetFunctionLowerBound(Tao tao, PetscReal fmin)
962d71ae5a4SJacob Faibussowitsch {
963a7e14dcfSSatish Balay   PetscFunctionBegin;
964441846f8SBarry Smith   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
96594511df7SStefano Zampini   PetscValidLogicalCollectiveReal(tao, fmin, 2);
966a7e14dcfSSatish Balay   tao->fmin = fmin;
9673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
968a7e14dcfSSatish Balay }
969a7e14dcfSSatish Balay 
970a7e14dcfSSatish Balay /*@
9718e96d397SJason Sarich   TaoGetFunctionLowerBound - Gets the bound on the solution objective value.
972a7e14dcfSSatish Balay   When an approximate solution with an objective value below this number
973a7e14dcfSSatish Balay   has been found, the solver will terminate.
974a7e14dcfSSatish Balay 
97567be906fSBarry Smith   Not Collective
976a7e14dcfSSatish Balay 
97747450a7bSBarry Smith   Input Parameter:
97847450a7bSBarry Smith . tao - the `Tao` solver context
979a7e14dcfSSatish Balay 
98047450a7bSBarry Smith   Output Parameter:
981a7e14dcfSSatish Balay . fmin - the minimum function value
982a7e14dcfSSatish Balay 
983a7e14dcfSSatish Balay   Level: intermediate
984a7e14dcfSSatish Balay 
9851cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoConvergedReason`, `TaoSetFunctionLowerBound()`
986a7e14dcfSSatish Balay @*/
987d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoGetFunctionLowerBound(Tao tao, PetscReal *fmin)
988d71ae5a4SJacob Faibussowitsch {
989a7e14dcfSSatish Balay   PetscFunctionBegin;
990441846f8SBarry Smith   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
9914f572ea9SToby Isaac   PetscAssertPointer(fmin, 2);
992a7e14dcfSSatish Balay   *fmin = tao->fmin;
9933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
994a7e14dcfSSatish Balay }
995a7e14dcfSSatish Balay 
996a7e14dcfSSatish Balay /*@
99747450a7bSBarry Smith   TaoSetMaximumFunctionEvaluations - Sets a maximum number of function evaluations allowed for a `TaoSolve()`.
998a7e14dcfSSatish Balay 
999c3339decSBarry Smith   Logically Collective
1000a7e14dcfSSatish Balay 
1001a7e14dcfSSatish Balay   Input Parameters:
100247450a7bSBarry Smith + tao  - the `Tao` solver context
1003606f75f6SBarry Smith - nfcn - the maximum number of function evaluations (>=0), use `PETSC_UNLIMITED` to have no bound
1004a7e14dcfSSatish Balay 
100547450a7bSBarry Smith   Options Database Key:
1006a7e14dcfSSatish Balay . -tao_max_funcs <nfcn> - sets the maximum number of function evaluations
1007a7e14dcfSSatish Balay 
1008a7e14dcfSSatish Balay   Level: intermediate
1009a7e14dcfSSatish Balay 
1010606f75f6SBarry Smith   Note:
1011606f75f6SBarry Smith   Use `PETSC_DETERMINE` to use the default maximum number of function evaluations that was set when the object type was set.
1012606f75f6SBarry Smith 
1013606f75f6SBarry Smith   Developer Note:
1014606f75f6SBarry Smith   Deprecated support for an unlimited number of function evaluations by passing a negative value.
1015606f75f6SBarry Smith 
10161cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoSetTolerances()`, `TaoSetMaximumIterations()`
1017a7e14dcfSSatish Balay @*/
1018d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoSetMaximumFunctionEvaluations(Tao tao, PetscInt nfcn)
1019d71ae5a4SJacob Faibussowitsch {
1020a7e14dcfSSatish Balay   PetscFunctionBegin;
1021441846f8SBarry Smith   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
102294511df7SStefano Zampini   PetscValidLogicalCollectiveInt(tao, nfcn, 2);
1023606f75f6SBarry Smith   if (nfcn == PETSC_DETERMINE) {
1024606f75f6SBarry Smith     tao->max_funcs = tao->default_max_funcs;
1025606f75f6SBarry Smith   } else if (nfcn == PETSC_UNLIMITED || nfcn < 0) {
1026606f75f6SBarry Smith     tao->max_funcs = PETSC_UNLIMITED;
10279371c9d4SSatish Balay   } else {
1028606f75f6SBarry Smith     PetscCheck(nfcn >= 0, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_OUTOFRANGE, "Maxium number of function evaluations  must be positive");
1029606f75f6SBarry Smith     tao->max_funcs = nfcn;
10309371c9d4SSatish Balay   }
10313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1032a7e14dcfSSatish Balay }
1033a7e14dcfSSatish Balay 
1034a7e14dcfSSatish Balay /*@
103547450a7bSBarry Smith   TaoGetMaximumFunctionEvaluations - Gets a maximum number of function evaluations allowed for a `TaoSolve()`
1036a7e14dcfSSatish Balay 
1037c3339decSBarry Smith   Logically Collective
1038a7e14dcfSSatish Balay 
103947450a7bSBarry Smith   Input Parameter:
104047450a7bSBarry Smith . tao - the `Tao` solver context
1041a7e14dcfSSatish Balay 
104247450a7bSBarry Smith   Output Parameter:
1043a7e14dcfSSatish Balay . nfcn - the maximum number of function evaluations
1044a7e14dcfSSatish Balay 
1045a7e14dcfSSatish Balay   Level: intermediate
1046a7e14dcfSSatish Balay 
10471cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoSetMaximumFunctionEvaluations()`, `TaoGetMaximumIterations()`
1048a7e14dcfSSatish Balay @*/
1049d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoGetMaximumFunctionEvaluations(Tao tao, PetscInt *nfcn)
1050d71ae5a4SJacob Faibussowitsch {
1051a7e14dcfSSatish Balay   PetscFunctionBegin;
1052441846f8SBarry Smith   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
10534f572ea9SToby Isaac   PetscAssertPointer(nfcn, 2);
1054a7e14dcfSSatish Balay   *nfcn = tao->max_funcs;
10553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1056a7e14dcfSSatish Balay }
1057a7e14dcfSSatish Balay 
1058770232b9SCe Qin /*@
105947450a7bSBarry Smith   TaoGetCurrentFunctionEvaluations - Get current number of function evaluations used by a `Tao` object
1060770232b9SCe Qin 
1061770232b9SCe Qin   Not Collective
1062770232b9SCe Qin 
106347450a7bSBarry Smith   Input Parameter:
106447450a7bSBarry Smith . tao - the `Tao` solver context
1065770232b9SCe Qin 
106647450a7bSBarry Smith   Output Parameter:
106794511df7SStefano Zampini . nfuncs - the current number of function evaluations (maximum between gradient and function evaluations)
1068770232b9SCe Qin 
1069770232b9SCe Qin   Level: intermediate
1070770232b9SCe Qin 
10711cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoSetMaximumFunctionEvaluations()`, `TaoGetMaximumFunctionEvaluations()`, `TaoGetMaximumIterations()`
1072770232b9SCe Qin @*/
1073d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoGetCurrentFunctionEvaluations(Tao tao, PetscInt *nfuncs)
1074d71ae5a4SJacob Faibussowitsch {
1075770232b9SCe Qin   PetscFunctionBegin;
1076770232b9SCe Qin   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
10774f572ea9SToby Isaac   PetscAssertPointer(nfuncs, 2);
1078770232b9SCe Qin   *nfuncs = PetscMax(tao->nfuncs, tao->nfuncgrads);
10793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1080770232b9SCe Qin }
1081770232b9SCe Qin 
1082a7e14dcfSSatish Balay /*@
108347450a7bSBarry Smith   TaoSetMaximumIterations - Sets a maximum number of iterates to be used in `TaoSolve()`
1084a7e14dcfSSatish Balay 
1085c3339decSBarry Smith   Logically Collective
1086a7e14dcfSSatish Balay 
1087a7e14dcfSSatish Balay   Input Parameters:
108847450a7bSBarry Smith + tao    - the `Tao` solver context
1089606f75f6SBarry Smith - maxits - the maximum number of iterates (>=0), use `PETSC_UNLIMITED` to have no bound
1090a7e14dcfSSatish Balay 
109147450a7bSBarry Smith   Options Database Key:
1092a7e14dcfSSatish Balay . -tao_max_it <its> - sets the maximum number of iterations
1093a7e14dcfSSatish Balay 
1094a7e14dcfSSatish Balay   Level: intermediate
1095a7e14dcfSSatish Balay 
1096606f75f6SBarry Smith   Note:
1097606f75f6SBarry Smith   Use `PETSC_DETERMINE` to use the default maximum number of iterations that was set when the object's type was set.
1098606f75f6SBarry Smith 
1099606f75f6SBarry Smith   Developer Note:
1100606f75f6SBarry Smith   DeprAlso accepts the deprecated negative values to indicate no limit
1101606f75f6SBarry Smith 
11021cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoSetTolerances()`, `TaoSetMaximumFunctionEvaluations()`
1103a7e14dcfSSatish Balay @*/
1104d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoSetMaximumIterations(Tao tao, PetscInt maxits)
1105d71ae5a4SJacob Faibussowitsch {
1106a7e14dcfSSatish Balay   PetscFunctionBegin;
1107441846f8SBarry Smith   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
110894511df7SStefano Zampini   PetscValidLogicalCollectiveInt(tao, maxits, 2);
1109606f75f6SBarry Smith   if (maxits == PETSC_DETERMINE) {
1110606f75f6SBarry Smith     tao->max_it = tao->default_max_it;
1111606f75f6SBarry Smith   } else if (maxits == PETSC_UNLIMITED) {
1112*1690c2aeSBarry Smith     tao->max_it = PETSC_INT_MAX;
1113606f75f6SBarry Smith   } else {
1114606f75f6SBarry Smith     PetscCheck(maxits > 0, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_OUTOFRANGE, "Maxium number of iterations must be positive");
1115606f75f6SBarry Smith     tao->max_it = maxits;
1116606f75f6SBarry Smith   }
11173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1118a7e14dcfSSatish Balay }
1119a7e14dcfSSatish Balay 
1120a7e14dcfSSatish Balay /*@
112165ba42b6SBarry Smith   TaoGetMaximumIterations - Gets a maximum number of iterates that will be used
1122a7e14dcfSSatish Balay 
1123a7e14dcfSSatish Balay   Not Collective
1124a7e14dcfSSatish Balay 
112547450a7bSBarry Smith   Input Parameter:
112647450a7bSBarry Smith . tao - the `Tao` solver context
1127a7e14dcfSSatish Balay 
112847450a7bSBarry Smith   Output Parameter:
1129a7e14dcfSSatish Balay . maxits - the maximum number of iterates
1130a7e14dcfSSatish Balay 
1131a7e14dcfSSatish Balay   Level: intermediate
1132a7e14dcfSSatish Balay 
11331cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoSetMaximumIterations()`, `TaoGetMaximumFunctionEvaluations()`
1134a7e14dcfSSatish Balay @*/
1135d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoGetMaximumIterations(Tao tao, PetscInt *maxits)
1136d71ae5a4SJacob Faibussowitsch {
1137a7e14dcfSSatish Balay   PetscFunctionBegin;
1138441846f8SBarry Smith   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
11394f572ea9SToby Isaac   PetscAssertPointer(maxits, 2);
1140a7e14dcfSSatish Balay   *maxits = tao->max_it;
11413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1142a7e14dcfSSatish Balay }
1143a7e14dcfSSatish Balay 
1144a7e14dcfSSatish Balay /*@
1145a7e14dcfSSatish Balay   TaoSetInitialTrustRegionRadius - Sets the initial trust region radius.
1146a7e14dcfSSatish Balay 
114767be906fSBarry Smith   Logically Collective
1148a7e14dcfSSatish Balay 
1149d8d19677SJose E. Roman   Input Parameters:
115047450a7bSBarry Smith + tao    - a `Tao` optimization solver
1151a7e14dcfSSatish Balay - radius - the trust region radius
1152a7e14dcfSSatish Balay 
1153a7e14dcfSSatish Balay   Options Database Key:
1154a7e14dcfSSatish Balay . -tao_trust0 <t0> - sets initial trust region radius
1155a7e14dcfSSatish Balay 
115647450a7bSBarry Smith   Level: intermediate
115747450a7bSBarry Smith 
1158606f75f6SBarry Smith   Note:
1159606f75f6SBarry Smith   Use `PETSC_DETERMINE` to use the default radius that was set when the object's type was set.
1160606f75f6SBarry Smith 
11611cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoGetTrustRegionRadius()`, `TaoSetTrustRegionTolerance()`, `TAONTR`
1162a7e14dcfSSatish Balay @*/
1163d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoSetInitialTrustRegionRadius(Tao tao, PetscReal radius)
1164d71ae5a4SJacob Faibussowitsch {
1165a7e14dcfSSatish Balay   PetscFunctionBegin;
1166441846f8SBarry Smith   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
116794511df7SStefano Zampini   PetscValidLogicalCollectiveReal(tao, radius, 2);
1168606f75f6SBarry Smith   if (radius == PETSC_DETERMINE) {
1169606f75f6SBarry Smith     tao->trust0 = tao->default_trust0;
1170606f75f6SBarry Smith   } else {
1171606f75f6SBarry Smith     PetscCheck(radius > 0, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_OUTOFRANGE, "Radius must be positive");
1172606f75f6SBarry Smith     tao->trust0 = radius;
1173606f75f6SBarry Smith   }
11743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1175a7e14dcfSSatish Balay }
1176a7e14dcfSSatish Balay 
1177a7e14dcfSSatish Balay /*@
117865ba42b6SBarry Smith   TaoGetInitialTrustRegionRadius - Gets the initial trust region radius.
1179a7e14dcfSSatish Balay 
1180a7e14dcfSSatish Balay   Not Collective
1181a7e14dcfSSatish Balay 
1182a7e14dcfSSatish Balay   Input Parameter:
118347450a7bSBarry Smith . tao - a `Tao` optimization solver
1184a7e14dcfSSatish Balay 
1185a7e14dcfSSatish Balay   Output Parameter:
1186a7e14dcfSSatish Balay . radius - the trust region radius
1187a7e14dcfSSatish Balay 
1188a7e14dcfSSatish Balay   Level: intermediate
1189a7e14dcfSSatish Balay 
11901cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoSetInitialTrustRegionRadius()`, `TaoGetCurrentTrustRegionRadius()`, `TAONTR`
1191a7e14dcfSSatish Balay @*/
1192d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoGetInitialTrustRegionRadius(Tao tao, PetscReal *radius)
1193d71ae5a4SJacob Faibussowitsch {
1194a7e14dcfSSatish Balay   PetscFunctionBegin;
1195441846f8SBarry Smith   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
11964f572ea9SToby Isaac   PetscAssertPointer(radius, 2);
1197a7e14dcfSSatish Balay   *radius = tao->trust0;
11983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1199a7e14dcfSSatish Balay }
1200a7e14dcfSSatish Balay 
1201a7e14dcfSSatish Balay /*@
1202a7e14dcfSSatish Balay   TaoGetCurrentTrustRegionRadius - Gets the current trust region radius.
1203a7e14dcfSSatish Balay 
1204a7e14dcfSSatish Balay   Not Collective
1205a7e14dcfSSatish Balay 
1206a7e14dcfSSatish Balay   Input Parameter:
120747450a7bSBarry Smith . tao - a `Tao` optimization solver
1208a7e14dcfSSatish Balay 
1209a7e14dcfSSatish Balay   Output Parameter:
1210a7e14dcfSSatish Balay . radius - the trust region radius
1211a7e14dcfSSatish Balay 
1212a7e14dcfSSatish Balay   Level: intermediate
1213a7e14dcfSSatish Balay 
12141cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoSetInitialTrustRegionRadius()`, `TaoGetInitialTrustRegionRadius()`, `TAONTR`
1215a7e14dcfSSatish Balay @*/
1216d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoGetCurrentTrustRegionRadius(Tao tao, PetscReal *radius)
1217d71ae5a4SJacob Faibussowitsch {
1218a7e14dcfSSatish Balay   PetscFunctionBegin;
1219441846f8SBarry Smith   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
12204f572ea9SToby Isaac   PetscAssertPointer(radius, 2);
1221a7e14dcfSSatish Balay   *radius = tao->trust;
12223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1223a7e14dcfSSatish Balay }
1224a7e14dcfSSatish Balay 
1225a7e14dcfSSatish Balay /*@
122647450a7bSBarry Smith   TaoGetTolerances - gets the current values of some tolerances used for the convergence testing of `TaoSolve()`
1227a7e14dcfSSatish Balay 
1228a7e14dcfSSatish Balay   Not Collective
1229a7e14dcfSSatish Balay 
1230f899ff85SJose E. Roman   Input Parameter:
123147450a7bSBarry Smith . tao - the `Tao` context
1232a7e14dcfSSatish Balay 
1233a7e14dcfSSatish Balay   Output Parameters:
1234e52336cbSBarry Smith + gatol - stop if norm of gradient is less than this
1235a7e14dcfSSatish Balay . grtol - stop if relative norm of gradient is less than this
1236a7e14dcfSSatish Balay - gttol - stop if norm of gradient is reduced by a this factor
1237a7e14dcfSSatish Balay 
1238a7e14dcfSSatish Balay   Level: intermediate
1239a1cb98faSBarry Smith 
1240a1cb98faSBarry Smith   Note:
124147450a7bSBarry Smith   `NULL` can be used as an argument if not all tolerances values are needed
1242a1cb98faSBarry Smith 
12431cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoSetTolerances()`
1244a7e14dcfSSatish Balay @*/
1245d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoGetTolerances(Tao tao, PetscReal *gatol, PetscReal *grtol, PetscReal *gttol)
1246d71ae5a4SJacob Faibussowitsch {
1247a7e14dcfSSatish Balay   PetscFunctionBegin;
1248441846f8SBarry Smith   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
1249a7e14dcfSSatish Balay   if (gatol) *gatol = tao->gatol;
1250a7e14dcfSSatish Balay   if (grtol) *grtol = tao->grtol;
1251a7e14dcfSSatish Balay   if (gttol) *gttol = tao->gttol;
12523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1253a7e14dcfSSatish Balay }
1254a7e14dcfSSatish Balay 
1255a7e14dcfSSatish Balay /*@
1256a7e14dcfSSatish Balay   TaoGetKSP - Gets the linear solver used by the optimization solver.
1257a7e14dcfSSatish Balay 
1258a7e14dcfSSatish Balay   Not Collective
1259a7e14dcfSSatish Balay 
126047450a7bSBarry Smith   Input Parameter:
126147450a7bSBarry Smith . tao - the `Tao` solver
1262a7e14dcfSSatish Balay 
126347450a7bSBarry Smith   Output Parameter:
126447450a7bSBarry Smith . ksp - the `KSP` linear solver used in the optimization solver
1265a7e14dcfSSatish Balay 
1266a7e14dcfSSatish Balay   Level: intermediate
1267a7e14dcfSSatish Balay 
12681cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `KSP`
1269a7e14dcfSSatish Balay @*/
1270d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoGetKSP(Tao tao, KSP *ksp)
1271d71ae5a4SJacob Faibussowitsch {
1272a7e14dcfSSatish Balay   PetscFunctionBegin;
127394511df7SStefano Zampini   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
12744f572ea9SToby Isaac   PetscAssertPointer(ksp, 2);
1275a7e14dcfSSatish Balay   *ksp = tao->ksp;
12763ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1277a7e14dcfSSatish Balay }
1278a7e14dcfSSatish Balay 
1279025e9500SJason Sarich /*@
1280025e9500SJason Sarich   TaoGetLinearSolveIterations - Gets the total number of linear iterations
128147450a7bSBarry Smith   used by the `Tao` solver
1282025e9500SJason Sarich 
1283025e9500SJason Sarich   Not Collective
1284025e9500SJason Sarich 
1285025e9500SJason Sarich   Input Parameter:
128647450a7bSBarry Smith . tao - the `Tao` context
1287025e9500SJason Sarich 
1288025e9500SJason Sarich   Output Parameter:
1289025e9500SJason Sarich . lits - number of linear iterations
1290025e9500SJason Sarich 
1291025e9500SJason Sarich   Level: intermediate
1292025e9500SJason Sarich 
129347450a7bSBarry Smith   Note:
129447450a7bSBarry Smith   This counter is reset to zero for each successive call to `TaoSolve()`
129547450a7bSBarry Smith 
12961cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoGetKSP()`
1297025e9500SJason Sarich @*/
1298d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoGetLinearSolveIterations(Tao tao, PetscInt *lits)
1299d71ae5a4SJacob Faibussowitsch {
1300025e9500SJason Sarich   PetscFunctionBegin;
1301025e9500SJason Sarich   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
13024f572ea9SToby Isaac   PetscAssertPointer(lits, 2);
13032d9aa51bSJason Sarich   *lits = tao->ksp_tot_its;
13043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1305025e9500SJason Sarich }
1306025e9500SJason Sarich 
1307a7e14dcfSSatish Balay /*@
1308a7e14dcfSSatish Balay   TaoGetLineSearch - Gets the line search used by the optimization solver.
1309a7e14dcfSSatish Balay 
1310a7e14dcfSSatish Balay   Not Collective
1311a7e14dcfSSatish Balay 
131247450a7bSBarry Smith   Input Parameter:
131347450a7bSBarry Smith . tao - the `Tao` solver
1314a7e14dcfSSatish Balay 
131547450a7bSBarry Smith   Output Parameter:
1316a7e14dcfSSatish Balay . ls - the line search used in the optimization solver
1317a7e14dcfSSatish Balay 
1318a7e14dcfSSatish Balay   Level: intermediate
1319a7e14dcfSSatish Balay 
13201cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoLineSearch`, `TaoLineSearchType`
1321a7e14dcfSSatish Balay @*/
1322d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoGetLineSearch(Tao tao, TaoLineSearch *ls)
1323d71ae5a4SJacob Faibussowitsch {
1324a7e14dcfSSatish Balay   PetscFunctionBegin;
132594511df7SStefano Zampini   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
13264f572ea9SToby Isaac   PetscAssertPointer(ls, 2);
1327a7e14dcfSSatish Balay   *ls = tao->linesearch;
13283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1329a7e14dcfSSatish Balay }
1330a7e14dcfSSatish Balay 
1331a7e14dcfSSatish Balay /*@
1332a7e14dcfSSatish Balay   TaoAddLineSearchCounts - Adds the number of function evaluations spent
1333a7e14dcfSSatish Balay   in the line search to the running total.
1334a7e14dcfSSatish Balay 
1335a7e14dcfSSatish Balay   Input Parameters:
13362fe279fdSBarry Smith . tao - the `Tao` solver
1337a7e14dcfSSatish Balay 
1338a7e14dcfSSatish Balay   Level: developer
1339a7e14dcfSSatish Balay 
13401cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoGetLineSearch()`, `TaoLineSearchApply()`
1341a7e14dcfSSatish Balay @*/
1342d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoAddLineSearchCounts(Tao tao)
1343d71ae5a4SJacob Faibussowitsch {
1344a7e14dcfSSatish Balay   PetscBool flg;
1345a7e14dcfSSatish Balay   PetscInt  nfeval, ngeval, nfgeval;
134647a47007SBarry Smith 
1347a7e14dcfSSatish Balay   PetscFunctionBegin;
1348441846f8SBarry Smith   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
1349a7e14dcfSSatish Balay   if (tao->linesearch) {
13509566063dSJacob Faibussowitsch     PetscCall(TaoLineSearchIsUsingTaoRoutines(tao->linesearch, &flg));
135194ae4db5SBarry Smith     if (!flg) {
13529566063dSJacob Faibussowitsch       PetscCall(TaoLineSearchGetNumberFunctionEvaluations(tao->linesearch, &nfeval, &ngeval, &nfgeval));
1353a7e14dcfSSatish Balay       tao->nfuncs += nfeval;
1354a7e14dcfSSatish Balay       tao->ngrads += ngeval;
1355a7e14dcfSSatish Balay       tao->nfuncgrads += nfgeval;
1356a7e14dcfSSatish Balay     }
1357a7e14dcfSSatish Balay   }
13583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1359a7e14dcfSSatish Balay }
1360a7e14dcfSSatish Balay 
1361a7e14dcfSSatish Balay /*@
136247450a7bSBarry Smith   TaoGetSolution - Returns the vector with the current solution from the `Tao` object
1363a7e14dcfSSatish Balay 
1364a7e14dcfSSatish Balay   Not Collective
1365a7e14dcfSSatish Balay 
1366a7e14dcfSSatish Balay   Input Parameter:
136747450a7bSBarry Smith . tao - the `Tao` context
1368a7e14dcfSSatish Balay 
1369a7e14dcfSSatish Balay   Output Parameter:
1370a7e14dcfSSatish Balay . X - the current solution
1371a7e14dcfSSatish Balay 
1372a7e14dcfSSatish Balay   Level: intermediate
1373a7e14dcfSSatish Balay 
137465ba42b6SBarry Smith   Note:
137565ba42b6SBarry Smith   The returned vector will be the same object that was passed into `TaoSetSolution()`
137665ba42b6SBarry Smith 
13771cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoSetSolution()`, `TaoSolve()`
1378a7e14dcfSSatish Balay @*/
1379d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoGetSolution(Tao tao, Vec *X)
1380d71ae5a4SJacob Faibussowitsch {
1381a7e14dcfSSatish Balay   PetscFunctionBegin;
1382441846f8SBarry Smith   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
13834f572ea9SToby Isaac   PetscAssertPointer(X, 2);
1384a7e14dcfSSatish Balay   *X = tao->solution;
13853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1386a7e14dcfSSatish Balay }
1387a7e14dcfSSatish Balay 
1388a7e14dcfSSatish Balay /*@
138947450a7bSBarry Smith   TaoResetStatistics - Initialize the statistics collected by the `Tao` object.
1390a7e14dcfSSatish Balay   These statistics include the iteration number, residual norms, and convergence status.
1391a7e14dcfSSatish Balay   This routine gets called before solving each optimization problem.
1392a7e14dcfSSatish Balay 
1393c3339decSBarry Smith   Collective
1394a7e14dcfSSatish Balay 
139547450a7bSBarry Smith   Input Parameter:
1396e056e8ceSJacob Faibussowitsch . tao - the `Tao` context
1397a7e14dcfSSatish Balay 
1398a7e14dcfSSatish Balay   Level: developer
1399a7e14dcfSSatish Balay 
14001cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoCreate()`, `TaoSolve()`
1401a7e14dcfSSatish Balay @*/
1402d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoResetStatistics(Tao tao)
1403d71ae5a4SJacob Faibussowitsch {
1404a7e14dcfSSatish Balay   PetscFunctionBegin;
1405441846f8SBarry Smith   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
1406a7e14dcfSSatish Balay   tao->niter        = 0;
1407a7e14dcfSSatish Balay   tao->nfuncs       = 0;
1408a7e14dcfSSatish Balay   tao->nfuncgrads   = 0;
1409a7e14dcfSSatish Balay   tao->ngrads       = 0;
1410a7e14dcfSSatish Balay   tao->nhess        = 0;
1411a7e14dcfSSatish Balay   tao->njac         = 0;
1412a7e14dcfSSatish Balay   tao->nconstraints = 0;
1413a7e14dcfSSatish Balay   tao->ksp_its      = 0;
1414ae93cb3cSJason Sarich   tao->ksp_tot_its  = 0;
1415a7e14dcfSSatish Balay   tao->reason       = TAO_CONTINUE_ITERATING;
1416a7e14dcfSSatish Balay   tao->residual     = 0.0;
1417a7e14dcfSSatish Balay   tao->cnorm        = 0.0;
1418a7e14dcfSSatish Balay   tao->step         = 0.0;
1419a7e14dcfSSatish Balay   tao->lsflag       = PETSC_FALSE;
1420a7e14dcfSSatish Balay   if (tao->hist_reset) tao->hist_len = 0;
14213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1422a7e14dcfSSatish Balay }
1423a7e14dcfSSatish Balay 
1424a7e14dcfSSatish Balay /*@C
1425e1e80dc8SAlp Dener   TaoSetUpdate - Sets the general-purpose update function called
14262fe279fdSBarry Smith   at the beginning of every iteration of the optimization algorithm. Called after the new solution and the gradient
1427e1e80dc8SAlp Dener   is determined, but before the Hessian is computed (if applicable).
1428e1e80dc8SAlp Dener 
1429c3339decSBarry Smith   Logically Collective
1430e1e80dc8SAlp Dener 
1431e1e80dc8SAlp Dener   Input Parameters:
1432270bebe6SStefano Zampini + tao  - The `Tao` solver
1433270bebe6SStefano Zampini . func - The function
1434270bebe6SStefano Zampini - ctx  - The update function context
1435e1e80dc8SAlp Dener 
143620f4b53cSBarry Smith   Calling sequence of `func`:
1437270bebe6SStefano Zampini + tao - The optimizer context
1438270bebe6SStefano Zampini . it  - The current iteration index
1439270bebe6SStefano Zampini - ctx - The update context
1440e1e80dc8SAlp Dener 
1441e1e80dc8SAlp Dener   Level: advanced
1442e1e80dc8SAlp Dener 
1443270bebe6SStefano Zampini   Notes:
1444270bebe6SStefano Zampini   Users can modify the gradient direction or any other vector associated to the specific solver used.
1445270bebe6SStefano Zampini   The objective function value is always recomputed after a call to the update hook.
1446270bebe6SStefano Zampini 
14471cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoSolve()`
1448e1e80dc8SAlp Dener @*/
1449270bebe6SStefano Zampini PetscErrorCode TaoSetUpdate(Tao tao, PetscErrorCode (*func)(Tao tao, PetscInt it, void *ctx), void *ctx)
1450d71ae5a4SJacob Faibussowitsch {
1451e1e80dc8SAlp Dener   PetscFunctionBegin;
1452e1e80dc8SAlp Dener   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
1453e1e80dc8SAlp Dener   tao->ops->update = func;
1454e1e80dc8SAlp Dener   tao->user_update = ctx;
14553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1456e1e80dc8SAlp Dener }
1457e1e80dc8SAlp Dener 
1458e1e80dc8SAlp Dener /*@C
1459a7e14dcfSSatish Balay   TaoSetConvergenceTest - Sets the function that is to be used to test
1460a7e14dcfSSatish Balay   for convergence of the iterative minimization solution.  The new convergence
146165ba42b6SBarry Smith   testing routine will replace Tao's default convergence test.
1462a7e14dcfSSatish Balay 
1463c3339decSBarry Smith   Logically Collective
1464a7e14dcfSSatish Balay 
1465a7e14dcfSSatish Balay   Input Parameters:
146647450a7bSBarry Smith + tao  - the `Tao` object
1467a7e14dcfSSatish Balay . conv - the routine to test for convergence
1468a7e14dcfSSatish Balay - ctx  - [optional] context for private data for the convergence routine
146947450a7bSBarry Smith         (may be `NULL`)
1470a7e14dcfSSatish Balay 
147120f4b53cSBarry Smith   Calling sequence of `conv`:
147247450a7bSBarry Smith + tao - the `Tao` object
1473a7e14dcfSSatish Balay - ctx - [optional] convergence context
1474a7e14dcfSSatish Balay 
147547450a7bSBarry Smith   Level: advanced
147647450a7bSBarry Smith 
147765ba42b6SBarry Smith   Note:
147865ba42b6SBarry Smith   The new convergence testing routine should call `TaoSetConvergedReason()`.
1479a7e14dcfSSatish Balay 
148010978b7dSBarry Smith .seealso: [](ch_tao), `Tao`, `TaoSolve()`, `TaoSetConvergedReason()`, `TaoGetSolutionStatus()`, `TaoGetTolerances()`, `TaoMonitorSet()`
1481a7e14dcfSSatish Balay @*/
1482d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoSetConvergenceTest(Tao tao, PetscErrorCode (*conv)(Tao, void *), void *ctx)
1483d71ae5a4SJacob Faibussowitsch {
1484a7e14dcfSSatish Balay   PetscFunctionBegin;
1485441846f8SBarry Smith   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
148694511df7SStefano Zampini   tao->ops->convergencetest = conv;
148794511df7SStefano Zampini   tao->cnvP                 = ctx;
14883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1489a7e14dcfSSatish Balay }
1490a7e14dcfSSatish Balay 
1491a7e14dcfSSatish Balay /*@C
149210978b7dSBarry Smith   TaoMonitorSet - Sets an additional function that is to be used at every
1493a7e14dcfSSatish Balay   iteration of the solver to display the iteration's
1494a7e14dcfSSatish Balay   progress.
1495a7e14dcfSSatish Balay 
1496c3339decSBarry Smith   Logically Collective
1497a7e14dcfSSatish Balay 
1498a7e14dcfSSatish Balay   Input Parameters:
149947450a7bSBarry Smith + tao  - the `Tao` solver context
15002fe279fdSBarry Smith . func - monitoring routine
15012fe279fdSBarry Smith . ctx  - [optional] user-defined context for private data for the monitor routine (may be `NULL`)
15022fe279fdSBarry Smith - dest - [optional] function to destroy the context when the `Tao` is destroyed
1503a7e14dcfSSatish Balay 
15042fe279fdSBarry Smith   Calling sequence of `func`:
150547450a7bSBarry Smith + tao - the `Tao` solver context
15062fe279fdSBarry Smith - ctx - [optional] monitoring context
15072fe279fdSBarry Smith 
15082fe279fdSBarry Smith   Calling sequence of `dest`:
15092fe279fdSBarry Smith . ctx - monitoring context
1510a7e14dcfSSatish Balay 
151147450a7bSBarry Smith   Level: intermediate
151247450a7bSBarry Smith 
151387497f52SBarry Smith   Notes:
151410978b7dSBarry Smith   See `TaoSetFromOptions()` for a monitoring options.
151510978b7dSBarry Smith 
1516a7e14dcfSSatish Balay   Several different monitoring routines may be set by calling
151710978b7dSBarry Smith   `TaoMonitorSet()` multiple times; all will be called in the
1518a7e14dcfSSatish Balay   order in which they were set.
1519a7e14dcfSSatish Balay 
1520e056e8ceSJacob Faibussowitsch   Fortran Notes:
152195452b02SPatrick Sanan   Only one monitor function may be set
1522a7e14dcfSSatish Balay 
1523b2dc45ddSPierre Jolivet .seealso: [](ch_tao), `Tao`, `TaoSolve()`, `TaoMonitorDefault()`, `TaoMonitorCancel()`, `TaoSetDestroyRoutine()`, `TaoView()`
1524a7e14dcfSSatish Balay @*/
152510978b7dSBarry Smith PetscErrorCode TaoMonitorSet(Tao tao, PetscErrorCode (*func)(Tao, void *), void *ctx, PetscErrorCode (*dest)(void **))
1526d71ae5a4SJacob Faibussowitsch {
152708d19d1fSJason Sarich   PetscInt  i;
15288732526dSAlp Dener   PetscBool identical;
152908d19d1fSJason Sarich 
1530a7e14dcfSSatish Balay   PetscFunctionBegin;
1531441846f8SBarry Smith   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
15323c859ba3SBarry Smith   PetscCheck(tao->numbermonitors < MAXTAOMONITORS, PetscObjectComm((PetscObject)tao), PETSC_ERR_SUP, "Cannot attach another monitor -- max=%d", MAXTAOMONITORS);
153308d19d1fSJason Sarich 
153408d19d1fSJason Sarich   for (i = 0; i < tao->numbermonitors; i++) {
15359566063dSJacob Faibussowitsch     PetscCall(PetscMonitorCompare((PetscErrorCode(*)(void))func, ctx, dest, (PetscErrorCode(*)(void))tao->monitor[i], tao->monitorcontext[i], tao->monitordestroy[i], &identical));
15363ba16761SJacob Faibussowitsch     if (identical) PetscFunctionReturn(PETSC_SUCCESS);
153708d19d1fSJason Sarich   }
1538a7e14dcfSSatish Balay   tao->monitor[tao->numbermonitors]        = func;
15398732526dSAlp Dener   tao->monitorcontext[tao->numbermonitors] = (void *)ctx;
1540a7e14dcfSSatish Balay   tao->monitordestroy[tao->numbermonitors] = dest;
1541a7e14dcfSSatish Balay   ++tao->numbermonitors;
15423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1543a7e14dcfSSatish Balay }
1544a7e14dcfSSatish Balay 
1545a7e14dcfSSatish Balay /*@
1546b2dc45ddSPierre Jolivet   TaoMonitorCancel - Clears all the monitor functions for a `Tao` object.
1547a7e14dcfSSatish Balay 
1548c3339decSBarry Smith   Logically Collective
1549a7e14dcfSSatish Balay 
155047450a7bSBarry Smith   Input Parameter:
155147450a7bSBarry Smith . tao - the `Tao` solver context
1552a7e14dcfSSatish Balay 
15533c7db156SBarry Smith   Options Database Key:
1554b2dc45ddSPierre Jolivet . -tao_monitor_cancel - cancels all monitors that have been hardwired
155510978b7dSBarry Smith     into a code by calls to `TaoMonitorSet()`, but does not cancel those
1556a7e14dcfSSatish Balay     set via the options database
1557a7e14dcfSSatish Balay 
1558a7e14dcfSSatish Balay   Level: advanced
1559a7e14dcfSSatish Balay 
156047450a7bSBarry Smith   Note:
156147450a7bSBarry Smith   There is no way to clear one specific monitor from a `Tao` object.
156247450a7bSBarry Smith 
156310978b7dSBarry Smith .seealso: [](ch_tao), `Tao`, `TaoMonitorDefault()`, `TaoMonitorSet()`
1564a7e14dcfSSatish Balay @*/
1565b2dc45ddSPierre Jolivet PetscErrorCode TaoMonitorCancel(Tao tao)
1566d71ae5a4SJacob Faibussowitsch {
1567a7e14dcfSSatish Balay   PetscInt i;
156847a47007SBarry Smith 
1569a7e14dcfSSatish Balay   PetscFunctionBegin;
1570441846f8SBarry Smith   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
1571a7e14dcfSSatish Balay   for (i = 0; i < tao->numbermonitors; i++) {
157248a46eb9SPierre Jolivet     if (tao->monitordestroy[i]) PetscCall((*tao->monitordestroy[i])(&tao->monitorcontext[i]));
1573a7e14dcfSSatish Balay   }
1574a7e14dcfSSatish Balay   tao->numbermonitors = 0;
15753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1576a7e14dcfSSatish Balay }
1577a7e14dcfSSatish Balay 
15787fab98ebSJason Sarich /*@
157947450a7bSBarry Smith   TaoMonitorDefault - Default routine for monitoring progress of `TaoSolve()`
1580a7e14dcfSSatish Balay 
1581c3339decSBarry Smith   Collective
1582a7e14dcfSSatish Balay 
1583a7e14dcfSSatish Balay   Input Parameters:
158447450a7bSBarry Smith + tao - the `Tao` context
158547450a7bSBarry Smith - ctx - `PetscViewer` context or `NULL`
1586a7e14dcfSSatish Balay 
158747450a7bSBarry Smith   Options Database Key:
1588147403d9SBarry Smith . -tao_monitor - turn on default monitoring
1589a7e14dcfSSatish Balay 
1590a7e14dcfSSatish Balay   Level: advanced
1591a7e14dcfSSatish Balay 
159247450a7bSBarry Smith   Note:
159347450a7bSBarry Smith   This monitor prints the function value and gradient
159447450a7bSBarry Smith   norm at each iteration.
159547450a7bSBarry Smith 
159610978b7dSBarry Smith .seealso: [](ch_tao), `Tao`, `TaoMonitorDefaultShort()`, `TaoMonitorSet()`
1597a7e14dcfSSatish Balay @*/
1598d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoMonitorDefault(Tao tao, void *ctx)
1599d71ae5a4SJacob Faibussowitsch {
160063b15415SAlp Dener   PetscInt    its, tabs;
1601a7e14dcfSSatish Balay   PetscReal   fct, gnorm;
16028163d661SBarry Smith   PetscViewer viewer = (PetscViewer)ctx;
160347a47007SBarry Smith 
1604a7e14dcfSSatish Balay   PetscFunctionBegin;
160594511df7SStefano Zampini   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
16068163d661SBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
1607a7e14dcfSSatish Balay   its   = tao->niter;
1608a7e14dcfSSatish Balay   fct   = tao->fc;
1609a7e14dcfSSatish Balay   gnorm = tao->residual;
16109566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIGetTab(viewer, &tabs));
16119566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIISetTab(viewer, ((PetscObject)tao)->tablevel));
1612494bef23SAlp Dener   if (its == 0 && ((PetscObject)tao)->prefix && !tao->header_printed) {
16139566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Iteration information for %s solve.\n", ((PetscObject)tao)->prefix));
1614494bef23SAlp Dener     tao->header_printed = PETSC_TRUE;
161563b15415SAlp Dener   }
161663a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT " TAO,", its));
16179566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "  Function value: %g,", (double)fct));
161808d19d1fSJason Sarich   if (gnorm >= PETSC_INFINITY) {
16199566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Residual: Inf \n"));
162008d19d1fSJason Sarich   } else {
16219566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Residual: %g \n", (double)gnorm));
162208d19d1fSJason Sarich   }
16239566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIISetTab(viewer, tabs));
16243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1625a7e14dcfSSatish Balay }
1626a7e14dcfSSatish Balay 
16277fab98ebSJason Sarich /*@
162810978b7dSBarry Smith   TaoMonitorGlobalization - Default routine for monitoring progress of `TaoSolve()` with extra detail on the globalization method.
16298d5ead36SAlp Dener 
1630c3339decSBarry Smith   Collective
16318d5ead36SAlp Dener 
16328d5ead36SAlp Dener   Input Parameters:
163347450a7bSBarry Smith + tao - the `Tao` context
163447450a7bSBarry Smith - ctx - `PetscViewer` context or `NULL`
16358d5ead36SAlp Dener 
163647450a7bSBarry Smith   Options Database Key:
163710978b7dSBarry Smith . -tao_monitor_globalization - turn on monitoring with globalization information
16388d5ead36SAlp Dener 
16398d5ead36SAlp Dener   Level: advanced
16408d5ead36SAlp Dener 
164147450a7bSBarry Smith   Note:
164247450a7bSBarry Smith   This monitor prints the function value and gradient norm at each
164347450a7bSBarry Smith   iteration, as well as the step size and trust radius. Note that the
164447450a7bSBarry Smith   step size and trust radius may be the same for some algorithms.
164547450a7bSBarry Smith 
164610978b7dSBarry Smith .seealso: [](ch_tao), `Tao`, `TaoMonitorDefaultShort()`, `TaoMonitorSet()`
16478d5ead36SAlp Dener @*/
164810978b7dSBarry Smith PetscErrorCode TaoMonitorGlobalization(Tao tao, void *ctx)
1649d71ae5a4SJacob Faibussowitsch {
16508d5ead36SAlp Dener   PetscInt    its, tabs;
16518d5ead36SAlp Dener   PetscReal   fct, gnorm, stp, tr;
16528d5ead36SAlp Dener   PetscViewer viewer = (PetscViewer)ctx;
16538d5ead36SAlp Dener 
16548d5ead36SAlp Dener   PetscFunctionBegin;
165594511df7SStefano Zampini   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
16568d5ead36SAlp Dener   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
16578d5ead36SAlp Dener   its   = tao->niter;
16588d5ead36SAlp Dener   fct   = tao->fc;
16598d5ead36SAlp Dener   gnorm = tao->residual;
16608d5ead36SAlp Dener   stp   = tao->step;
16618d5ead36SAlp Dener   tr    = tao->trust;
16629566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIGetTab(viewer, &tabs));
16639566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIISetTab(viewer, ((PetscObject)tao)->tablevel));
1664494bef23SAlp Dener   if (its == 0 && ((PetscObject)tao)->prefix && !tao->header_printed) {
16659566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Iteration information for %s solve.\n", ((PetscObject)tao)->prefix));
1666494bef23SAlp Dener     tao->header_printed = PETSC_TRUE;
16678d5ead36SAlp Dener   }
166863a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "%3" PetscInt_FMT " TAO,", its));
16699566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "  Function value: %g,", (double)fct));
16708d5ead36SAlp Dener   if (gnorm >= PETSC_INFINITY) {
16719566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Residual: Inf,"));
16728d5ead36SAlp Dener   } else {
16739566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  Residual: %g,", (double)gnorm));
16748d5ead36SAlp Dener   }
16759566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "  Step: %g,  Trust: %g\n", (double)stp, (double)tr));
16769566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIISetTab(viewer, tabs));
16773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16788d5ead36SAlp Dener }
16798d5ead36SAlp Dener 
16808d5ead36SAlp Dener /*@
168110978b7dSBarry Smith   TaoMonitorDefaultShort - Routine for monitoring progress of `TaoSolve()` that displays fewer digits than `TaoMonitorDefault()`
1682a7e14dcfSSatish Balay 
1683c3339decSBarry Smith   Collective
1684a7e14dcfSSatish Balay 
1685a7e14dcfSSatish Balay   Input Parameters:
168647450a7bSBarry Smith + tao - the `Tao` context
168747450a7bSBarry Smith - ctx - `PetscViewer` context of type `PETSCVIEWERASCII`
1688a7e14dcfSSatish Balay 
168947450a7bSBarry Smith   Options Database Key:
169010978b7dSBarry Smith . -tao_monitor_short - turn on default short monitoring
1691a7e14dcfSSatish Balay 
1692a7e14dcfSSatish Balay   Level: advanced
1693a7e14dcfSSatish Balay 
169447450a7bSBarry Smith   Note:
169547450a7bSBarry Smith   Same as `TaoMonitorDefault()` except
169647450a7bSBarry Smith   it prints fewer digits of the residual as the residual gets smaller.
169747450a7bSBarry Smith   This is because the later digits are meaningless and are often
169847450a7bSBarry Smith   different on different machines; by using this routine different
169947450a7bSBarry Smith   machines will usually generate the same output.
170047450a7bSBarry Smith 
170110978b7dSBarry Smith .seealso: [](ch_tao), `Tao`, `TaoMonitorDefault()`, `TaoMonitorSet()`
1702a7e14dcfSSatish Balay @*/
170310978b7dSBarry Smith PetscErrorCode TaoMonitorDefaultShort(Tao tao, void *ctx)
1704d71ae5a4SJacob Faibussowitsch {
17052c9b8b83SAlp Dener   PetscInt    its, tabs;
1706a7e14dcfSSatish Balay   PetscReal   fct, gnorm;
17074d4332d5SBarry Smith   PetscViewer viewer = (PetscViewer)ctx;
170847a47007SBarry Smith 
1709a7e14dcfSSatish Balay   PetscFunctionBegin;
171094511df7SStefano Zampini   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
17114d4332d5SBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
1712a7e14dcfSSatish Balay   its   = tao->niter;
1713a7e14dcfSSatish Balay   fct   = tao->fc;
1714a7e14dcfSSatish Balay   gnorm = tao->residual;
17159566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIGetTab(viewer, &tabs));
17169566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIISetTab(viewer, ((PetscObject)tao)->tablevel));
171763a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "iter = %3" PetscInt_FMT ",", its));
17189566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, " Function value %g,", (double)fct));
1719d393f493SBarry Smith   if (gnorm >= PETSC_INFINITY) {
17209566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, " Residual: Inf \n"));
172108d19d1fSJason Sarich   } else if (gnorm > 1.e-6) {
17229566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, " Residual: %g \n", (double)gnorm));
1723a7e14dcfSSatish Balay   } else if (gnorm > 1.e-11) {
17249566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, " Residual: < 1.0e-6 \n"));
1725a7e14dcfSSatish Balay   } else {
17269566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, " Residual: < 1.0e-11 \n"));
1727a7e14dcfSSatish Balay   }
17289566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIISetTab(viewer, tabs));
17293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1730a7e14dcfSSatish Balay }
1731a7e14dcfSSatish Balay 
17327fab98ebSJason Sarich /*@
173310978b7dSBarry Smith   TaoMonitorConstraintNorm - same as `TaoMonitorDefault()` except
173447450a7bSBarry Smith   it prints the norm of the constraint function.
1735a7e14dcfSSatish Balay 
1736c3339decSBarry Smith   Collective
1737a7e14dcfSSatish Balay 
1738a7e14dcfSSatish Balay   Input Parameters:
173947450a7bSBarry Smith + tao - the `Tao` context
174047450a7bSBarry Smith - ctx - `PetscViewer` context or `NULL`
1741a7e14dcfSSatish Balay 
174247450a7bSBarry Smith   Options Database Key:
174310978b7dSBarry Smith . -tao_monitor_constraint_norm - monitor the constraints
1744a7e14dcfSSatish Balay 
1745a7e14dcfSSatish Balay   Level: advanced
1746a7e14dcfSSatish Balay 
174710978b7dSBarry Smith .seealso: [](ch_tao), `Tao`, `TaoMonitorDefault()`, `TaoMonitorSet()`
1748a7e14dcfSSatish Balay @*/
174910978b7dSBarry Smith PetscErrorCode TaoMonitorConstraintNorm(Tao tao, void *ctx)
1750d71ae5a4SJacob Faibussowitsch {
17512c9b8b83SAlp Dener   PetscInt    its, tabs;
1752a7e14dcfSSatish Balay   PetscReal   fct, gnorm;
1753d393f493SBarry Smith   PetscViewer viewer = (PetscViewer)ctx;
1754a7e14dcfSSatish Balay 
1755a7e14dcfSSatish Balay   PetscFunctionBegin;
175694511df7SStefano Zampini   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
1757d393f493SBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
1758a7e14dcfSSatish Balay   its   = tao->niter;
1759a7e14dcfSSatish Balay   fct   = tao->fc;
1760a7e14dcfSSatish Balay   gnorm = tao->residual;
17619566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIGetTab(viewer, &tabs));
17629566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIISetTab(viewer, ((PetscObject)tao)->tablevel));
176363a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "iter = %" PetscInt_FMT ",", its));
17649566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, " Function value: %g,", (double)fct));
17659566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "  Residual: %g ", (double)gnorm));
17669566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "  Constraint: %g \n", (double)tao->cnorm));
17679566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIISetTab(viewer, tabs));
17683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1769a7e14dcfSSatish Balay }
1770a7e14dcfSSatish Balay 
1771a7e14dcfSSatish Balay /*@C
177210978b7dSBarry Smith   TaoMonitorSolution - Views the solution at each iteration of `TaoSolve()`
1773a7e14dcfSSatish Balay 
1774c3339decSBarry Smith   Collective
1775a7e14dcfSSatish Balay 
1776a7e14dcfSSatish Balay   Input Parameters:
177747450a7bSBarry Smith + tao - the `Tao` context
177847450a7bSBarry Smith - ctx - `PetscViewer` context or `NULL`
1779a7e14dcfSSatish Balay 
178047450a7bSBarry Smith   Options Database Key:
178110978b7dSBarry Smith . -tao_monitor_solution - view the solution
1782a7e14dcfSSatish Balay 
1783a7e14dcfSSatish Balay   Level: advanced
1784a7e14dcfSSatish Balay 
178510978b7dSBarry Smith .seealso: [](ch_tao), `Tao`, `TaoMonitorDefaultShort()`, `TaoMonitorSet()`
1786a7e14dcfSSatish Balay @*/
178710978b7dSBarry Smith PetscErrorCode TaoMonitorSolution(Tao tao, void *ctx)
1788d71ae5a4SJacob Faibussowitsch {
1789feb237baSPierre Jolivet   PetscViewer viewer = (PetscViewer)ctx;
179047a47007SBarry Smith 
1791a7e14dcfSSatish Balay   PetscFunctionBegin;
179294511df7SStefano Zampini   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
1793d393f493SBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
17949566063dSJacob Faibussowitsch   PetscCall(VecView(tao->solution, viewer));
17953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1796a7e14dcfSSatish Balay }
1797a7e14dcfSSatish Balay 
1798a7e14dcfSSatish Balay /*@C
179910978b7dSBarry Smith   TaoMonitorGradient - Views the gradient at each iteration of `TaoSolve()`
1800a7e14dcfSSatish Balay 
1801c3339decSBarry Smith   Collective
1802a7e14dcfSSatish Balay 
1803a7e14dcfSSatish Balay   Input Parameters:
180447450a7bSBarry Smith + tao - the `Tao` context
180547450a7bSBarry Smith - ctx - `PetscViewer` context or `NULL`
1806a7e14dcfSSatish Balay 
180747450a7bSBarry Smith   Options Database Key:
180810978b7dSBarry Smith . -tao_monitor_gradient - view the gradient at each iteration
1809a7e14dcfSSatish Balay 
1810a7e14dcfSSatish Balay   Level: advanced
1811a7e14dcfSSatish Balay 
181210978b7dSBarry Smith .seealso: [](ch_tao), `Tao`, `TaoMonitorDefaultShort()`, `TaoMonitorSet()`
1813a7e14dcfSSatish Balay @*/
181410978b7dSBarry Smith PetscErrorCode TaoMonitorGradient(Tao tao, void *ctx)
1815d71ae5a4SJacob Faibussowitsch {
1816d393f493SBarry Smith   PetscViewer viewer = (PetscViewer)ctx;
181747a47007SBarry Smith 
1818a7e14dcfSSatish Balay   PetscFunctionBegin;
181994511df7SStefano Zampini   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
1820d393f493SBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
18219566063dSJacob Faibussowitsch   PetscCall(VecView(tao->gradient, viewer));
18223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1823a7e14dcfSSatish Balay }
1824a7e14dcfSSatish Balay 
1825a7e14dcfSSatish Balay /*@C
182610978b7dSBarry Smith   TaoMonitorStep - Views the step-direction at each iteration of `TaoSolve()`
1827a7e14dcfSSatish Balay 
1828c3339decSBarry Smith   Collective
1829a7e14dcfSSatish Balay 
1830a7e14dcfSSatish Balay   Input Parameters:
183147450a7bSBarry Smith + tao - the `Tao` context
183247450a7bSBarry Smith - ctx - `PetscViewer` context or `NULL`
1833a7e14dcfSSatish Balay 
183447450a7bSBarry Smith   Options Database Key:
183510978b7dSBarry Smith . -tao_monitor_step - view the step vector at each iteration
1836a7e14dcfSSatish Balay 
1837a7e14dcfSSatish Balay   Level: advanced
1838a7e14dcfSSatish Balay 
183910978b7dSBarry Smith .seealso: [](ch_tao), `Tao`, `TaoMonitorDefaultShort()`, `TaoMonitorSet()`
1840a7e14dcfSSatish Balay @*/
184110978b7dSBarry Smith PetscErrorCode TaoMonitorStep(Tao tao, void *ctx)
1842d71ae5a4SJacob Faibussowitsch {
1843d393f493SBarry Smith   PetscViewer viewer = (PetscViewer)ctx;
1844d393f493SBarry Smith 
1845a7e14dcfSSatish Balay   PetscFunctionBegin;
184694511df7SStefano Zampini   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
1847d393f493SBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
18489566063dSJacob Faibussowitsch   PetscCall(VecView(tao->stepdirection, viewer));
18493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1850a7e14dcfSSatish Balay }
1851a7e14dcfSSatish Balay 
1852a7e14dcfSSatish Balay /*@C
185310978b7dSBarry Smith   TaoMonitorSolutionDraw - Plots the solution at each iteration of `TaoSolve()`
1854a7e14dcfSSatish Balay 
1855c3339decSBarry Smith   Collective
1856a7e14dcfSSatish Balay 
1857a7e14dcfSSatish Balay   Input Parameters:
185847450a7bSBarry Smith + tao - the `Tao` context
185965ba42b6SBarry Smith - ctx - `TaoMonitorDraw` context
1860a7e14dcfSSatish Balay 
186147450a7bSBarry Smith   Options Database Key:
186210978b7dSBarry Smith . -tao_monitor_solution_draw - draw the solution at each iteration
1863a7e14dcfSSatish Balay 
1864a7e14dcfSSatish Balay   Level: advanced
1865a7e14dcfSSatish Balay 
1866a81c3919SBarry Smith   Note:
1867a81c3919SBarry Smith   The context created by `TaoMonitorDrawCtxCreate()`, along with `TaoMonitorSolutionDraw()`, and `TaoMonitorDrawCtxDestroy()`
1868a81c3919SBarry Smith   are passed to `TaoMonitorSet()` to monitor the solution graphically.
1869a81c3919SBarry Smith 
1870a81c3919SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoMonitorSolution()`, `TaoMonitorSet()`, `TaoMonitorGradientDraw()`, `TaoMonitorDrawCtxCreate()`,
1871a81c3919SBarry Smith           `TaoMonitorDrawCtxDestroy()`
1872a7e14dcfSSatish Balay @*/
187310978b7dSBarry Smith PetscErrorCode TaoMonitorSolutionDraw(Tao tao, void *ctx)
1874d71ae5a4SJacob Faibussowitsch {
1875e882e171SHong Zhang   TaoMonitorDrawCtx ictx = (TaoMonitorDrawCtx)ctx;
187647a47007SBarry Smith 
1877a7e14dcfSSatish Balay   PetscFunctionBegin;
187894511df7SStefano Zampini   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
18793ba16761SJacob Faibussowitsch   if (!(((ictx->howoften > 0) && (!(tao->niter % ictx->howoften))) || ((ictx->howoften == -1) && tao->reason))) PetscFunctionReturn(PETSC_SUCCESS);
18809566063dSJacob Faibussowitsch   PetscCall(VecView(tao->solution, ictx->viewer));
18813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1882a7e14dcfSSatish Balay }
1883a7e14dcfSSatish Balay 
1884a7e14dcfSSatish Balay /*@C
188510978b7dSBarry Smith   TaoMonitorGradientDraw - Plots the gradient at each iteration of `TaoSolve()`
1886a7e14dcfSSatish Balay 
1887c3339decSBarry Smith   Collective
1888a7e14dcfSSatish Balay 
1889a7e14dcfSSatish Balay   Input Parameters:
189047450a7bSBarry Smith + tao - the `Tao` context
189165ba42b6SBarry Smith - ctx - `PetscViewer` context
1892a7e14dcfSSatish Balay 
189347450a7bSBarry Smith   Options Database Key:
189410978b7dSBarry Smith . -tao_monitor_gradient_draw - draw the gradient at each iteration
1895a7e14dcfSSatish Balay 
1896a7e14dcfSSatish Balay   Level: advanced
1897a7e14dcfSSatish Balay 
189810978b7dSBarry Smith .seealso: [](ch_tao), `Tao`, `TaoMonitorGradient()`, `TaoMonitorSet()`, `TaoMonitorSolutionDraw()`
1899a7e14dcfSSatish Balay @*/
190010978b7dSBarry Smith PetscErrorCode TaoMonitorGradientDraw(Tao tao, void *ctx)
1901d71ae5a4SJacob Faibussowitsch {
1902e882e171SHong Zhang   TaoMonitorDrawCtx ictx = (TaoMonitorDrawCtx)ctx;
190347a47007SBarry Smith 
1904a7e14dcfSSatish Balay   PetscFunctionBegin;
190594511df7SStefano Zampini   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
19063ba16761SJacob Faibussowitsch   if (!(((ictx->howoften > 0) && (!(tao->niter % ictx->howoften))) || ((ictx->howoften == -1) && tao->reason))) PetscFunctionReturn(PETSC_SUCCESS);
19079566063dSJacob Faibussowitsch   PetscCall(VecView(tao->gradient, ictx->viewer));
19083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1909a7e14dcfSSatish Balay }
1910a7e14dcfSSatish Balay 
1911a7e14dcfSSatish Balay /*@C
191210978b7dSBarry Smith   TaoMonitorStepDraw - Plots the step direction at each iteration of `TaoSolve()`
1913a7e14dcfSSatish Balay 
1914c3339decSBarry Smith   Collective
1915a7e14dcfSSatish Balay 
1916a7e14dcfSSatish Balay   Input Parameters:
191747450a7bSBarry Smith + tao - the `Tao` context
191847450a7bSBarry Smith - ctx - the `PetscViewer` context
1919a7e14dcfSSatish Balay 
192047450a7bSBarry Smith   Options Database Key:
192110978b7dSBarry Smith . -tao_monitor_step_draw - draw the step direction at each iteration
1922a7e14dcfSSatish Balay 
1923a7e14dcfSSatish Balay   Level: advanced
1924a7e14dcfSSatish Balay 
192510978b7dSBarry Smith .seealso: [](ch_tao), `Tao`, `TaoMonitorSet()`, `TaoMonitorSolutionDraw`
1926a7e14dcfSSatish Balay @*/
192710978b7dSBarry Smith PetscErrorCode TaoMonitorStepDraw(Tao tao, void *ctx)
1928d71ae5a4SJacob Faibussowitsch {
192994511df7SStefano Zampini   PetscViewer viewer = (PetscViewer)ctx;
193047a47007SBarry Smith 
1931a7e14dcfSSatish Balay   PetscFunctionBegin;
193294511df7SStefano Zampini   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
193394511df7SStefano Zampini   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
19349566063dSJacob Faibussowitsch   PetscCall(VecView(tao->stepdirection, viewer));
19353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1936a7e14dcfSSatish Balay }
1937a7e14dcfSSatish Balay 
1938a7e14dcfSSatish Balay /*@C
193910978b7dSBarry Smith   TaoMonitorResidual - Views the least-squares residual at each iteration of `TaoSolve()`
1940a7e14dcfSSatish Balay 
1941c3339decSBarry Smith   Collective
1942a7e14dcfSSatish Balay 
1943a7e14dcfSSatish Balay   Input Parameters:
194447450a7bSBarry Smith + tao - the `Tao` context
194547450a7bSBarry Smith - ctx - the `PetscViewer` context or `NULL`
1946a7e14dcfSSatish Balay 
194747450a7bSBarry Smith   Options Database Key:
194810978b7dSBarry Smith . -tao_monitor_ls_residual - view the residual at each iteration
1949a7e14dcfSSatish Balay 
1950a7e14dcfSSatish Balay   Level: advanced
1951a7e14dcfSSatish Balay 
195210978b7dSBarry Smith .seealso: [](ch_tao), `Tao`, `TaoMonitorDefaultShort()`, `TaoMonitorSet()`
1953a7e14dcfSSatish Balay @*/
195410978b7dSBarry Smith PetscErrorCode TaoMonitorResidual(Tao tao, void *ctx)
1955d71ae5a4SJacob Faibussowitsch {
1956d393f493SBarry Smith   PetscViewer viewer = (PetscViewer)ctx;
195747a47007SBarry Smith 
1958a7e14dcfSSatish Balay   PetscFunctionBegin;
195994511df7SStefano Zampini   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
1960d393f493SBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
19619566063dSJacob Faibussowitsch   PetscCall(VecView(tao->ls_res, viewer));
19623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1963a7e14dcfSSatish Balay }
1964a7e14dcfSSatish Balay 
19657fab98ebSJason Sarich /*@
1966a7e14dcfSSatish Balay   TaoDefaultConvergenceTest - Determines whether the solver should continue iterating
1967a7e14dcfSSatish Balay   or terminate.
1968a7e14dcfSSatish Balay 
1969c3339decSBarry Smith   Collective
1970a7e14dcfSSatish Balay 
1971a7e14dcfSSatish Balay   Input Parameters:
197247450a7bSBarry Smith + tao   - the `Tao` context
1973a7e14dcfSSatish Balay - dummy - unused dummy context
1974a7e14dcfSSatish Balay 
197547450a7bSBarry Smith   Level: developer
197647450a7bSBarry Smith 
1977a7e14dcfSSatish Balay   Notes:
1978a7e14dcfSSatish Balay   This routine checks the residual in the optimality conditions, the
1979a7e14dcfSSatish Balay   relative residual in the optimity conditions, the number of function
1980a7e14dcfSSatish Balay   evaluations, and the function value to test convergence.  Some
1981a7e14dcfSSatish Balay   solvers may use different convergence routines.
1982a7e14dcfSSatish Balay 
19831cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoSetTolerances()`, `TaoGetConvergedReason()`, `TaoSetConvergedReason()`
1984a7e14dcfSSatish Balay @*/
1985d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoDefaultConvergenceTest(Tao tao, void *dummy)
1986d71ae5a4SJacob Faibussowitsch {
1987a7e14dcfSSatish Balay   PetscInt           niter = tao->niter, nfuncs = PetscMax(tao->nfuncs, tao->nfuncgrads);
1988a7e14dcfSSatish Balay   PetscInt           max_funcs = tao->max_funcs;
1989a7e14dcfSSatish Balay   PetscReal          gnorm = tao->residual, gnorm0 = tao->gnorm0;
1990a7e14dcfSSatish Balay   PetscReal          f = tao->fc, steptol = tao->steptol, trradius = tao->step;
1991a7e14dcfSSatish Balay   PetscReal          gatol = tao->gatol, grtol = tao->grtol, gttol = tao->gttol;
1992e52336cbSBarry Smith   PetscReal          catol = tao->catol, crtol = tao->crtol;
1993e52336cbSBarry Smith   PetscReal          fmin = tao->fmin, cnorm = tao->cnorm;
1994e4cb33bbSBarry Smith   TaoConvergedReason reason = tao->reason;
1995a7e14dcfSSatish Balay 
1996a7e14dcfSSatish Balay   PetscFunctionBegin;
1997441846f8SBarry Smith   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
19983ba16761SJacob Faibussowitsch   if (reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(PETSC_SUCCESS);
1999a7e14dcfSSatish Balay 
2000a7e14dcfSSatish Balay   if (PetscIsInfOrNanReal(f)) {
20019566063dSJacob Faibussowitsch     PetscCall(PetscInfo(tao, "Failed to converged, function value is Inf or NaN\n"));
2002a7e14dcfSSatish Balay     reason = TAO_DIVERGED_NAN;
2003a7e14dcfSSatish Balay   } else if (f <= fmin && cnorm <= catol) {
20049566063dSJacob Faibussowitsch     PetscCall(PetscInfo(tao, "Converged due to function value %g < minimum function value %g\n", (double)f, (double)fmin));
2005a7e14dcfSSatish Balay     reason = TAO_CONVERGED_MINF;
2006a7e14dcfSSatish Balay   } else if (gnorm <= gatol && cnorm <= catol) {
20079566063dSJacob Faibussowitsch     PetscCall(PetscInfo(tao, "Converged due to residual norm ||g(X)||=%g < %g\n", (double)gnorm, (double)gatol));
2008a7e14dcfSSatish Balay     reason = TAO_CONVERGED_GATOL;
2009a7e14dcfSSatish Balay   } else if (f != 0 && PetscAbsReal(gnorm / f) <= grtol && cnorm <= crtol) {
20109566063dSJacob Faibussowitsch     PetscCall(PetscInfo(tao, "Converged due to residual ||g(X)||/|f(X)| =%g < %g\n", (double)(gnorm / f), (double)grtol));
2011a7e14dcfSSatish Balay     reason = TAO_CONVERGED_GRTOL;
2012e1d16bb9SBarry Smith   } else if (gnorm0 != 0 && ((gttol == 0 && gnorm == 0) || gnorm / gnorm0 < gttol) && cnorm <= crtol) {
20139566063dSJacob Faibussowitsch     PetscCall(PetscInfo(tao, "Converged due to relative residual norm ||g(X)||/||g(X0)|| = %g < %g\n", (double)(gnorm / gnorm0), (double)gttol));
2014a7e14dcfSSatish Balay     reason = TAO_CONVERGED_GTTOL;
2015606f75f6SBarry Smith   } else if (max_funcs != PETSC_UNLIMITED && nfuncs > max_funcs) {
20169566063dSJacob Faibussowitsch     PetscCall(PetscInfo(tao, "Exceeded maximum number of function evaluations: %" PetscInt_FMT " > %" PetscInt_FMT "\n", nfuncs, max_funcs));
2017a7e14dcfSSatish Balay     reason = TAO_DIVERGED_MAXFCN;
2018a7e14dcfSSatish Balay   } else if (tao->lsflag != 0) {
20199566063dSJacob Faibussowitsch     PetscCall(PetscInfo(tao, "Tao Line Search failure.\n"));
2020a7e14dcfSSatish Balay     reason = TAO_DIVERGED_LS_FAILURE;
2021a7e14dcfSSatish Balay   } else if (trradius < steptol && niter > 0) {
20229566063dSJacob Faibussowitsch     PetscCall(PetscInfo(tao, "Trust region/step size too small: %g < %g\n", (double)trradius, (double)steptol));
2023a7e14dcfSSatish Balay     reason = TAO_CONVERGED_STEPTOL;
2024e031d6f5SAlp Dener   } else if (niter >= tao->max_it) {
202563a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(tao, "Exceeded maximum number of iterations: %" PetscInt_FMT " > %" PetscInt_FMT "\n", niter, tao->max_it));
2026a7e14dcfSSatish Balay     reason = TAO_DIVERGED_MAXITS;
2027a7e14dcfSSatish Balay   } else {
2028a7e14dcfSSatish Balay     reason = TAO_CONTINUE_ITERATING;
2029a7e14dcfSSatish Balay   }
2030a7e14dcfSSatish Balay   tao->reason = reason;
20313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2032a7e14dcfSSatish Balay }
2033a7e14dcfSSatish Balay 
2034cc4c1da9SBarry Smith /*@
2035a7e14dcfSSatish Balay   TaoSetOptionsPrefix - Sets the prefix used for searching for all
203665ba42b6SBarry Smith   Tao options in the database.
2037a7e14dcfSSatish Balay 
2038c3339decSBarry Smith   Logically Collective
2039a7e14dcfSSatish Balay 
2040a7e14dcfSSatish Balay   Input Parameters:
204147450a7bSBarry Smith + tao - the `Tao` context
2042e056e8ceSJacob Faibussowitsch - p   - the prefix string to prepend to all Tao option requests
2043a7e14dcfSSatish Balay 
20442fe279fdSBarry Smith   Level: advanced
20452fe279fdSBarry Smith 
2046a7e14dcfSSatish Balay   Notes:
2047a7e14dcfSSatish Balay   A hyphen (-) must NOT be given at the beginning of the prefix name.
2048a7e14dcfSSatish Balay   The first character of all runtime options is AUTOMATICALLY the hyphen.
2049a7e14dcfSSatish Balay 
2050a7e14dcfSSatish Balay   For example, to distinguish between the runtime options for two
205165ba42b6SBarry Smith   different Tao solvers, one could call
2052a7e14dcfSSatish Balay .vb
2053a7e14dcfSSatish Balay       TaoSetOptionsPrefix(tao1,"sys1_")
2054a7e14dcfSSatish Balay       TaoSetOptionsPrefix(tao2,"sys2_")
2055a7e14dcfSSatish Balay .ve
2056a7e14dcfSSatish Balay 
2057a7e14dcfSSatish Balay   This would enable use of different options for each system, such as
2058a7e14dcfSSatish Balay .vb
20599fa2b5dcSStefano Zampini       -sys1_tao_method blmvm -sys1_tao_grtol 1.e-3
20609fa2b5dcSStefano Zampini       -sys2_tao_method lmvm  -sys2_tao_grtol 1.e-4
2061a7e14dcfSSatish Balay .ve
2062a7e14dcfSSatish Balay 
20631cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoSetFromOptions()`, `TaoAppendOptionsPrefix()`, `TaoGetOptionsPrefix()`
2064a7e14dcfSSatish Balay @*/
2065d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoSetOptionsPrefix(Tao tao, const char p[])
2066d71ae5a4SJacob Faibussowitsch {
2067a7e14dcfSSatish Balay   PetscFunctionBegin;
206894511df7SStefano Zampini   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
20699566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)tao, p));
20701baa6e33SBarry Smith   if (tao->linesearch) PetscCall(TaoLineSearchSetOptionsPrefix(tao->linesearch, p));
20711baa6e33SBarry Smith   if (tao->ksp) PetscCall(KSPSetOptionsPrefix(tao->ksp, p));
20723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2073a7e14dcfSSatish Balay }
2074a7e14dcfSSatish Balay 
2075cc4c1da9SBarry Smith /*@
207647450a7bSBarry Smith   TaoAppendOptionsPrefix - Appends to the prefix used for searching for all Tao options in the database.
2077a7e14dcfSSatish Balay 
2078c3339decSBarry Smith   Logically Collective
2079a7e14dcfSSatish Balay 
2080a7e14dcfSSatish Balay   Input Parameters:
208147450a7bSBarry Smith + tao - the `Tao` solver context
2082e056e8ceSJacob Faibussowitsch - p   - the prefix string to prepend to all `Tao` option requests
2083a7e14dcfSSatish Balay 
20842fe279fdSBarry Smith   Level: advanced
20852fe279fdSBarry Smith 
208665ba42b6SBarry Smith   Note:
2087a7e14dcfSSatish Balay   A hyphen (-) must NOT be given at the beginning of the prefix name.
208865ba42b6SBarry Smith   The first character of all runtime options is automatically the hyphen.
2089a7e14dcfSSatish Balay 
20901cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoSetFromOptions()`, `TaoSetOptionsPrefix()`, `TaoGetOptionsPrefix()`
2091a7e14dcfSSatish Balay @*/
2092d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoAppendOptionsPrefix(Tao tao, const char p[])
2093d71ae5a4SJacob Faibussowitsch {
2094a7e14dcfSSatish Balay   PetscFunctionBegin;
209594511df7SStefano Zampini   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
20969566063dSJacob Faibussowitsch   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)tao, p));
20971baa6e33SBarry Smith   if (tao->linesearch) PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)tao->linesearch, p));
20981baa6e33SBarry Smith   if (tao->ksp) PetscCall(KSPAppendOptionsPrefix(tao->ksp, p));
20993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2100a7e14dcfSSatish Balay }
2101a7e14dcfSSatish Balay 
2102cc4c1da9SBarry Smith /*@
2103a7e14dcfSSatish Balay   TaoGetOptionsPrefix - Gets the prefix used for searching for all
210465ba42b6SBarry Smith   Tao options in the database
2105a7e14dcfSSatish Balay 
2106a7e14dcfSSatish Balay   Not Collective
2107a7e14dcfSSatish Balay 
210847450a7bSBarry Smith   Input Parameter:
210947450a7bSBarry Smith . tao - the `Tao` context
2110a7e14dcfSSatish Balay 
211147450a7bSBarry Smith   Output Parameter:
2112e056e8ceSJacob Faibussowitsch . p - pointer to the prefix string used is returned
2113a7e14dcfSSatish Balay 
2114e056e8ceSJacob Faibussowitsch   Fortran Notes:
211547450a7bSBarry Smith   Pass in a string 'prefix' of sufficient length to hold the prefix.
2116a7e14dcfSSatish Balay 
2117a7e14dcfSSatish Balay   Level: advanced
2118a7e14dcfSSatish Balay 
21191cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoSetFromOptions()`, `TaoSetOptionsPrefix()`, `TaoAppendOptionsPrefix()`
2120a7e14dcfSSatish Balay @*/
2121d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoGetOptionsPrefix(Tao tao, const char *p[])
2122d71ae5a4SJacob Faibussowitsch {
212394511df7SStefano Zampini   PetscFunctionBegin;
212494511df7SStefano Zampini   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
21259566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)tao, p));
21263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2127a7e14dcfSSatish Balay }
2128a7e14dcfSSatish Balay 
2129cc4c1da9SBarry Smith /*@
213047450a7bSBarry Smith   TaoSetType - Sets the `TaoType` for the minimization solver.
2131a7e14dcfSSatish Balay 
2132c3339decSBarry Smith   Collective
2133a7e14dcfSSatish Balay 
2134a7e14dcfSSatish Balay   Input Parameters:
2135e056e8ceSJacob Faibussowitsch + tao  - the `Tao` solver context
2136a7e14dcfSSatish Balay - type - a known method
2137a7e14dcfSSatish Balay 
2138a7e14dcfSSatish Balay   Options Database Key:
213958417fe7SBarry Smith . -tao_type <type> - Sets the method; use -help for a list
214058417fe7SBarry Smith    of available methods (for instance, "-tao_type lmvm" or "-tao_type tron")
2141a7e14dcfSSatish Balay 
2142a7e14dcfSSatish Balay   Level: intermediate
2143a7e14dcfSSatish Balay 
21441cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoCreate()`, `TaoGetType()`, `TaoType`
2145a7e14dcfSSatish Balay @*/
2146d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoSetType(Tao tao, TaoType type)
2147d71ae5a4SJacob Faibussowitsch {
2148441846f8SBarry Smith   PetscErrorCode (*create_xxx)(Tao);
2149a7e14dcfSSatish Balay   PetscBool issame;
2150a7e14dcfSSatish Balay 
2151a7e14dcfSSatish Balay   PetscFunctionBegin;
2152441846f8SBarry Smith   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
2153a7e14dcfSSatish Balay 
21549566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)tao, type, &issame));
21553ba16761SJacob Faibussowitsch   if (issame) PetscFunctionReturn(PETSC_SUCCESS);
2156a7e14dcfSSatish Balay 
21579566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListFind(TaoList, type, (void (**)(void)) & create_xxx));
21583c859ba3SBarry Smith   PetscCheck(create_xxx, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unable to find requested Tao type %s", type);
2159a7e14dcfSSatish Balay 
2160a7e14dcfSSatish Balay   /* Destroy the existing solver information */
2161dbbe0bcdSBarry Smith   PetscTryTypeMethod(tao, destroy);
21621baa6e33SBarry Smith   PetscCall(KSPDestroy(&tao->ksp));
21639566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchDestroy(&tao->linesearch));
216483c8fe1dSLisandro Dalcin   tao->ops->setup          = NULL;
216583c8fe1dSLisandro Dalcin   tao->ops->solve          = NULL;
216683c8fe1dSLisandro Dalcin   tao->ops->view           = NULL;
216783c8fe1dSLisandro Dalcin   tao->ops->setfromoptions = NULL;
216883c8fe1dSLisandro Dalcin   tao->ops->destroy        = NULL;
2169a7e14dcfSSatish Balay 
2170a7e14dcfSSatish Balay   tao->setupcalled = PETSC_FALSE;
2171a7e14dcfSSatish Balay 
21729566063dSJacob Faibussowitsch   PetscCall((*create_xxx)(tao));
21739566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)tao, type));
21743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2175a7e14dcfSSatish Balay }
217647a47007SBarry Smith 
217767be906fSBarry Smith /*@C
217847450a7bSBarry Smith   TaoRegister - Adds a method to the Tao package for minimization.
2179a7e14dcfSSatish Balay 
2180cc4c1da9SBarry Smith   Not Collective, No Fortran Support
2181a7e14dcfSSatish Balay 
2182a7e14dcfSSatish Balay   Input Parameters:
2183a7e14dcfSSatish Balay + sname - name of a new user-defined solver
2184a7e14dcfSSatish Balay - func  - routine to Create method context
2185a7e14dcfSSatish Balay 
2186e056e8ceSJacob Faibussowitsch   Example Usage:
2187a7e14dcfSSatish Balay .vb
2188441846f8SBarry Smith    TaoRegister("my_solver", MySolverCreate);
2189a7e14dcfSSatish Balay .ve
2190a7e14dcfSSatish Balay 
2191a7e14dcfSSatish Balay   Then, your solver can be chosen with the procedural interface via
2192a7e14dcfSSatish Balay $     TaoSetType(tao, "my_solver")
2193a7e14dcfSSatish Balay   or at runtime via the option
219458417fe7SBarry Smith $     -tao_type my_solver
2195a7e14dcfSSatish Balay 
2196a7e14dcfSSatish Balay   Level: advanced
2197a7e14dcfSSatish Balay 
219867be906fSBarry Smith   Note:
219967be906fSBarry Smith   `TaoRegister()` may be called multiple times to add several user-defined solvers.
220067be906fSBarry Smith 
22011cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoSetType()`, `TaoRegisterAll()`, `TaoRegisterDestroy()`
220267be906fSBarry Smith @*/
2203d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoRegister(const char sname[], PetscErrorCode (*func)(Tao))
2204d71ae5a4SJacob Faibussowitsch {
2205a7e14dcfSSatish Balay   PetscFunctionBegin;
22069566063dSJacob Faibussowitsch   PetscCall(TaoInitializePackage());
22079566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&TaoList, sname, (void (*)(void))func));
22083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2209a7e14dcfSSatish Balay }
2210a7e14dcfSSatish Balay 
2211a7e14dcfSSatish Balay /*@C
2212441846f8SBarry Smith   TaoRegisterDestroy - Frees the list of minimization solvers that were
221347450a7bSBarry Smith   registered by `TaoRegister()`.
2214a7e14dcfSSatish Balay 
2215a7e14dcfSSatish Balay   Not Collective
2216a7e14dcfSSatish Balay 
2217a7e14dcfSSatish Balay   Level: advanced
2218a7e14dcfSSatish Balay 
22191cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoRegisterAll()`, `TaoRegister()`
2220a7e14dcfSSatish Balay @*/
2221d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoRegisterDestroy(void)
2222d71ae5a4SJacob Faibussowitsch {
2223a7e14dcfSSatish Balay   PetscFunctionBegin;
22249566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListDestroy(&TaoList));
2225441846f8SBarry Smith   TaoRegisterAllCalled = PETSC_FALSE;
22263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2227a7e14dcfSSatish Balay }
222847a47007SBarry Smith 
22298931d482SJason Sarich /*@
223047450a7bSBarry Smith   TaoGetIterationNumber - Gets the number of `TaoSolve()` iterations completed
22318931d482SJason Sarich   at this time.
22328931d482SJason Sarich 
22338931d482SJason Sarich   Not Collective
22348931d482SJason Sarich 
22358931d482SJason Sarich   Input Parameter:
223647450a7bSBarry Smith . tao - the `Tao` context
22378931d482SJason Sarich 
22388931d482SJason Sarich   Output Parameter:
22398931d482SJason Sarich . iter - iteration number
22408931d482SJason Sarich 
22418931d482SJason Sarich   Notes:
22428931d482SJason Sarich   For example, during the computation of iteration 2 this would return 1.
22438931d482SJason Sarich 
22448931d482SJason Sarich   Level: intermediate
22458931d482SJason Sarich 
22461cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoGetLinearSolveIterations()`, `TaoGetResidualNorm()`, `TaoGetObjective()`
22478931d482SJason Sarich @*/
2248d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoGetIterationNumber(Tao tao, PetscInt *iter)
2249d71ae5a4SJacob Faibussowitsch {
22508931d482SJason Sarich   PetscFunctionBegin;
22518931d482SJason Sarich   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
22524f572ea9SToby Isaac   PetscAssertPointer(iter, 2);
22538931d482SJason Sarich   *iter = tao->niter;
22543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
22558931d482SJason Sarich }
22568931d482SJason Sarich 
22578931d482SJason Sarich /*@
225847450a7bSBarry Smith   TaoGetResidualNorm - Gets the current value of the norm of the residual (gradient)
225979f5d8caSBarry Smith   at this time.
226079f5d8caSBarry Smith 
226179f5d8caSBarry Smith   Not Collective
226279f5d8caSBarry Smith 
226379f5d8caSBarry Smith   Input Parameter:
226447450a7bSBarry Smith . tao - the `Tao` context
226579f5d8caSBarry Smith 
226679f5d8caSBarry Smith   Output Parameter:
226779f5d8caSBarry Smith . value - the current value
226879f5d8caSBarry Smith 
226979f5d8caSBarry Smith   Level: intermediate
227079f5d8caSBarry Smith 
2271e056e8ceSJacob Faibussowitsch   Developer Notes:
227267be906fSBarry Smith   This is the 2-norm of the residual, we cannot use `TaoGetGradientNorm()` because that has
227347450a7bSBarry Smith   a different meaning. For some reason `Tao` sometimes calls the gradient the residual.
227479f5d8caSBarry Smith 
22751cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoGetLinearSolveIterations()`, `TaoGetIterationNumber()`, `TaoGetObjective()`
227679f5d8caSBarry Smith @*/
2277d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoGetResidualNorm(Tao tao, PetscReal *value)
2278d71ae5a4SJacob Faibussowitsch {
227979f5d8caSBarry Smith   PetscFunctionBegin;
228079f5d8caSBarry Smith   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
22814f572ea9SToby Isaac   PetscAssertPointer(value, 2);
228279f5d8caSBarry Smith   *value = tao->residual;
22833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
228479f5d8caSBarry Smith }
228579f5d8caSBarry Smith 
228679f5d8caSBarry Smith /*@
22878931d482SJason Sarich   TaoSetIterationNumber - Sets the current iteration number.
22888931d482SJason Sarich 
2289c3339decSBarry Smith   Logically Collective
22908931d482SJason Sarich 
2291d8d19677SJose E. Roman   Input Parameters:
229247450a7bSBarry Smith + tao  - the `Tao` context
2293a2b725a8SWilliam Gropp - iter - iteration number
22948931d482SJason Sarich 
22958931d482SJason Sarich   Level: developer
22968931d482SJason Sarich 
22971cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoGetLinearSolveIterations()`
22988931d482SJason Sarich @*/
2299d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoSetIterationNumber(Tao tao, PetscInt iter)
2300d71ae5a4SJacob Faibussowitsch {
23018931d482SJason Sarich   PetscFunctionBegin;
23028931d482SJason Sarich   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
230394511df7SStefano Zampini   PetscValidLogicalCollectiveInt(tao, iter, 2);
23049566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)tao));
23058931d482SJason Sarich   tao->niter = iter;
23069566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)tao));
23073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
23088931d482SJason Sarich }
23098931d482SJason Sarich 
23108931d482SJason Sarich /*@
231147450a7bSBarry Smith   TaoGetTotalIterationNumber - Gets the total number of `TaoSolve()` iterations
23128931d482SJason Sarich   completed. This number keeps accumulating if multiple solves
231347450a7bSBarry Smith   are called with the `Tao` object.
23148931d482SJason Sarich 
23158931d482SJason Sarich   Not Collective
23168931d482SJason Sarich 
23178931d482SJason Sarich   Input Parameter:
231847450a7bSBarry Smith . tao - the `Tao` context
23198931d482SJason Sarich 
23208931d482SJason Sarich   Output Parameter:
232147450a7bSBarry Smith . iter - number of iterations
23228931d482SJason Sarich 
23238931d482SJason Sarich   Level: intermediate
23248931d482SJason Sarich 
232520f4b53cSBarry Smith   Note:
232667be906fSBarry Smith   The total iteration count is updated after each solve, if there is a current
232747450a7bSBarry Smith   `TaoSolve()` in progress then those iterations are not included in the count
232867be906fSBarry Smith 
23291cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoGetLinearSolveIterations()`
23308931d482SJason Sarich @*/
2331d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoGetTotalIterationNumber(Tao tao, PetscInt *iter)
2332d71ae5a4SJacob Faibussowitsch {
23338931d482SJason Sarich   PetscFunctionBegin;
23348931d482SJason Sarich   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
23354f572ea9SToby Isaac   PetscAssertPointer(iter, 2);
23368931d482SJason Sarich   *iter = tao->ntotalits;
23373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
23388931d482SJason Sarich }
23398931d482SJason Sarich 
23408931d482SJason Sarich /*@
23418931d482SJason Sarich   TaoSetTotalIterationNumber - Sets the current total iteration number.
23428931d482SJason Sarich 
2343c3339decSBarry Smith   Logically Collective
23448931d482SJason Sarich 
2345d8d19677SJose E. Roman   Input Parameters:
234647450a7bSBarry Smith + tao  - the `Tao` context
234747450a7bSBarry Smith - iter - the iteration number
23488931d482SJason Sarich 
23498931d482SJason Sarich   Level: developer
23508931d482SJason Sarich 
23511cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoGetLinearSolveIterations()`
23528931d482SJason Sarich @*/
2353d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoSetTotalIterationNumber(Tao tao, PetscInt iter)
2354d71ae5a4SJacob Faibussowitsch {
23558931d482SJason Sarich   PetscFunctionBegin;
23568931d482SJason Sarich   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
235794511df7SStefano Zampini   PetscValidLogicalCollectiveInt(tao, iter, 2);
23589566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)tao));
23598931d482SJason Sarich   tao->ntotalits = iter;
23609566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)tao));
23613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
23628931d482SJason Sarich }
23638931d482SJason Sarich 
2364a7e14dcfSSatish Balay /*@
236547450a7bSBarry Smith   TaoSetConvergedReason - Sets the termination flag on a `Tao` object
2366a7e14dcfSSatish Balay 
2367c3339decSBarry Smith   Logically Collective
2368a7e14dcfSSatish Balay 
2369a7e14dcfSSatish Balay   Input Parameters:
237047450a7bSBarry Smith + tao    - the `Tao` context
237147450a7bSBarry Smith - reason - the `TaoConvergedReason`
2372a7e14dcfSSatish Balay 
2373a7e14dcfSSatish Balay   Level: intermediate
2374a7e14dcfSSatish Balay 
23751cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoConvergedReason`
2376a7e14dcfSSatish Balay @*/
2377d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoSetConvergedReason(Tao tao, TaoConvergedReason reason)
2378d71ae5a4SJacob Faibussowitsch {
2379a7e14dcfSSatish Balay   PetscFunctionBegin;
238094511df7SStefano Zampini   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
238194511df7SStefano Zampini   PetscValidLogicalCollectiveEnum(tao, reason, 2);
2382a7e14dcfSSatish Balay   tao->reason = reason;
23833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2384a7e14dcfSSatish Balay }
2385a7e14dcfSSatish Balay 
2386a7e14dcfSSatish Balay /*@
238747450a7bSBarry Smith   TaoGetConvergedReason - Gets the reason the `TaoSolve()` was stopped.
2388a7e14dcfSSatish Balay 
2389a7e14dcfSSatish Balay   Not Collective
2390a7e14dcfSSatish Balay 
2391a7e14dcfSSatish Balay   Input Parameter:
239247450a7bSBarry Smith . tao - the `Tao` solver context
2393a7e14dcfSSatish Balay 
2394a7e14dcfSSatish Balay   Output Parameter:
239547450a7bSBarry Smith . reason - value of `TaoConvergedReason`
2396a7e14dcfSSatish Balay 
2397a7e14dcfSSatish Balay   Level: intermediate
2398a7e14dcfSSatish Balay 
23991cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoConvergedReason`, `TaoSetConvergenceTest()`, `TaoSetTolerances()`
2400a7e14dcfSSatish Balay @*/
2401d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoGetConvergedReason(Tao tao, TaoConvergedReason *reason)
2402d71ae5a4SJacob Faibussowitsch {
2403a7e14dcfSSatish Balay   PetscFunctionBegin;
2404441846f8SBarry Smith   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
24054f572ea9SToby Isaac   PetscAssertPointer(reason, 2);
2406a7e14dcfSSatish Balay   *reason = tao->reason;
24073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2408a7e14dcfSSatish Balay }
2409a7e14dcfSSatish Balay 
2410a7e14dcfSSatish Balay /*@
2411a7e14dcfSSatish Balay   TaoGetSolutionStatus - Get the current iterate, objective value,
241247450a7bSBarry Smith   residual, infeasibility, and termination from a `Tao` object
2413a7e14dcfSSatish Balay 
2414a7e14dcfSSatish Balay   Not Collective
2415a7e14dcfSSatish Balay 
2416f899ff85SJose E. Roman   Input Parameter:
241747450a7bSBarry Smith . tao - the `Tao` context
2418a7e14dcfSSatish Balay 
2419a7e14dcfSSatish Balay   Output Parameters:
2420e056e8ceSJacob Faibussowitsch + its    - the current iterate number (>=0)
2421a7e14dcfSSatish Balay . f      - the current function value
242258417fe7SBarry Smith . gnorm  - the square of the gradient norm, duality gap, or other measure indicating distance from optimality.
2423a7e14dcfSSatish Balay . cnorm  - the infeasibility of the current solution with regard to the constraints.
2424a7e14dcfSSatish Balay . xdiff  - the step length or trust region radius of the most recent iterate.
242565ba42b6SBarry Smith - reason - The termination reason, which can equal `TAO_CONTINUE_ITERATING`
2426a7e14dcfSSatish Balay 
2427a7e14dcfSSatish Balay   Level: intermediate
2428a7e14dcfSSatish Balay 
242965ba42b6SBarry Smith   Notes:
243065ba42b6SBarry Smith   Tao returns the values set by the solvers in the routine `TaoMonitor()`.
2431a7e14dcfSSatish Balay 
243267be906fSBarry Smith   If any of the output arguments are set to `NULL`, no corresponding value will be returned.
2433a7e14dcfSSatish Balay 
24341cc06b55SBarry Smith .seealso: [](ch_tao), `TaoMonitor()`, `TaoGetConvergedReason()`
2435a7e14dcfSSatish Balay @*/
2436d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoGetSolutionStatus(Tao tao, PetscInt *its, PetscReal *f, PetscReal *gnorm, PetscReal *cnorm, PetscReal *xdiff, TaoConvergedReason *reason)
2437d71ae5a4SJacob Faibussowitsch {
2438a7e14dcfSSatish Balay   PetscFunctionBegin;
243994511df7SStefano Zampini   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
2440a7e14dcfSSatish Balay   if (its) *its = tao->niter;
2441a7e14dcfSSatish Balay   if (f) *f = tao->fc;
2442a7e14dcfSSatish Balay   if (gnorm) *gnorm = tao->residual;
2443a7e14dcfSSatish Balay   if (cnorm) *cnorm = tao->cnorm;
2444a7e14dcfSSatish Balay   if (reason) *reason = tao->reason;
2445a7e14dcfSSatish Balay   if (xdiff) *xdiff = tao->step;
24463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2447a7e14dcfSSatish Balay }
2448a7e14dcfSSatish Balay 
2449cc4c1da9SBarry Smith /*@
245047450a7bSBarry Smith   TaoGetType - Gets the current `TaoType` being used in the `Tao` object
2451a7e14dcfSSatish Balay 
2452a7e14dcfSSatish Balay   Not Collective
2453a7e14dcfSSatish Balay 
2454a7e14dcfSSatish Balay   Input Parameter:
245547450a7bSBarry Smith . tao - the `Tao` solver context
2456a7e14dcfSSatish Balay 
2457a7e14dcfSSatish Balay   Output Parameter:
245847450a7bSBarry Smith . type - the `TaoType`
2459a7e14dcfSSatish Balay 
2460a7e14dcfSSatish Balay   Level: intermediate
2461a7e14dcfSSatish Balay 
24621cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoType`, `TaoSetType()`
2463a7e14dcfSSatish Balay @*/
2464d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoGetType(Tao tao, TaoType *type)
2465d71ae5a4SJacob Faibussowitsch {
2466a7e14dcfSSatish Balay   PetscFunctionBegin;
2467441846f8SBarry Smith   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
24684f572ea9SToby Isaac   PetscAssertPointer(type, 2);
2469a7e14dcfSSatish Balay   *type = ((PetscObject)tao)->type_name;
24703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2471a7e14dcfSSatish Balay }
2472a7e14dcfSSatish Balay 
2473a7e14dcfSSatish Balay /*@C
2474a7e14dcfSSatish Balay   TaoMonitor - Monitor the solver and the current solution.  This
2475a7e14dcfSSatish Balay   routine will record the iteration number and residual statistics,
2476670c031eSStefano Zampini   and call any monitors specified by the user.
2477a7e14dcfSSatish Balay 
2478a7e14dcfSSatish Balay   Input Parameters:
247947450a7bSBarry Smith + tao        - the `Tao` context
2480a7e14dcfSSatish Balay . its        - the current iterate number (>=0)
2481a7e14dcfSSatish Balay . f          - the current objective function value
24821cb3f120SPierre Jolivet . res        - the gradient norm, square root of the duality gap, or other measure indicating distance from optimality.  This measure will be recorded and
2483a7e14dcfSSatish Balay           used for some termination tests.
2484a7e14dcfSSatish Balay . cnorm      - the infeasibility of the current solution with regard to the constraints.
2485a7e14dcfSSatish Balay - steplength - multiple of the step direction added to the previous iterate.
2486a7e14dcfSSatish Balay 
2487a7e14dcfSSatish Balay   Options Database Key:
2488a7e14dcfSSatish Balay . -tao_monitor - Use the default monitor, which prints statistics to standard output
2489a7e14dcfSSatish Balay 
2490a7e14dcfSSatish Balay   Level: developer
2491a7e14dcfSSatish Balay 
249210978b7dSBarry Smith .seealso: [](ch_tao), `Tao`, `TaoGetConvergedReason()`, `TaoMonitorDefault()`, `TaoMonitorSet()`
2493a7e14dcfSSatish Balay @*/
2494d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoMonitor(Tao tao, PetscInt its, PetscReal f, PetscReal res, PetscReal cnorm, PetscReal steplength)
2495d71ae5a4SJacob Faibussowitsch {
249647a47007SBarry Smith   PetscInt i;
249747a47007SBarry Smith 
2498a7e14dcfSSatish Balay   PetscFunctionBegin;
2499441846f8SBarry Smith   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
2500a7e14dcfSSatish Balay   tao->fc       = f;
2501a7e14dcfSSatish Balay   tao->residual = res;
2502a7e14dcfSSatish Balay   tao->cnorm    = cnorm;
2503a7e14dcfSSatish Balay   tao->step     = steplength;
2504e52336cbSBarry Smith   if (!its) {
250594511df7SStefano Zampini     tao->cnorm0 = cnorm;
250694511df7SStefano Zampini     tao->gnorm0 = res;
2507a7e14dcfSSatish Balay   }
2508bfdcf5a2SStefano Zampini   PetscCall(VecLockReadPush(tao->solution));
250948a46eb9SPierre Jolivet   for (i = 0; i < tao->numbermonitors; i++) PetscCall((*tao->monitor[i])(tao, tao->monitorcontext[i]));
2510bfdcf5a2SStefano Zampini   PetscCall(VecLockReadPop(tao->solution));
25113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2512a7e14dcfSSatish Balay }
2513a7e14dcfSSatish Balay 
2514a7e14dcfSSatish Balay /*@
2515ae93cb3cSJason Sarich   TaoSetConvergenceHistory - Sets the array used to hold the convergence history.
2516a7e14dcfSSatish Balay 
2517c3339decSBarry Smith   Logically Collective
2518a7e14dcfSSatish Balay 
2519a7e14dcfSSatish Balay   Input Parameters:
252047450a7bSBarry Smith + tao   - the `Tao` solver context
2521a7e14dcfSSatish Balay . obj   - array to hold objective value history
2522a7e14dcfSSatish Balay . resid - array to hold residual history
2523a7e14dcfSSatish Balay . cnorm - array to hold constraint violation history
2524be1558d8SJason Sarich . lits  - integer array holds the number of linear iterations for each Tao iteration
252567be906fSBarry Smith . na    - size of `obj`, `resid`, and `cnorm`
252665ba42b6SBarry Smith - reset - `PETSC_TRUE` indicates each new minimization resets the history counter to zero,
2527a7e14dcfSSatish Balay            else it continues storing new values for new minimizations after the old ones
2528a7e14dcfSSatish Balay 
252967be906fSBarry Smith   Level: intermediate
253067be906fSBarry Smith 
2531a7e14dcfSSatish Balay   Notes:
253247450a7bSBarry Smith   If set, `Tao` will fill the given arrays with the indicated
2533ae93cb3cSJason Sarich   information at each iteration.  If 'obj','resid','cnorm','lits' are
2534606f75f6SBarry Smith   *all* `NULL` then space (using size `na`, or 1000 if `na` is `PETSC_DECIDE`) is allocated for the history.
253567be906fSBarry Smith   If not all are `NULL`, then only the non-`NULL` information categories
2536ae93cb3cSJason Sarich   will be stored, the others will be ignored.
2537ae93cb3cSJason Sarich 
2538ae93cb3cSJason Sarich   Any convergence information after iteration number 'na' will not be stored.
2539a7e14dcfSSatish Balay 
2540a7e14dcfSSatish Balay   This routine is useful, e.g., when running a code for purposes
2541a7e14dcfSSatish Balay   of accurate performance monitoring, when no I/O should be done
2542a7e14dcfSSatish Balay   during the section of code that is being timed.
2543a7e14dcfSSatish Balay 
25441cc06b55SBarry Smith .seealso: [](ch_tao), `TaoGetConvergenceHistory()`
2545a7e14dcfSSatish Balay @*/
2546d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoSetConvergenceHistory(Tao tao, PetscReal obj[], PetscReal resid[], PetscReal cnorm[], PetscInt lits[], PetscInt na, PetscBool reset)
2547d71ae5a4SJacob Faibussowitsch {
2548a7e14dcfSSatish Balay   PetscFunctionBegin;
2549441846f8SBarry Smith   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
25504f572ea9SToby Isaac   if (obj) PetscAssertPointer(obj, 2);
25514f572ea9SToby Isaac   if (resid) PetscAssertPointer(resid, 3);
25524f572ea9SToby Isaac   if (cnorm) PetscAssertPointer(cnorm, 4);
25534f572ea9SToby Isaac   if (lits) PetscAssertPointer(lits, 5);
2554ae93cb3cSJason Sarich 
2555606f75f6SBarry Smith   if (na == PETSC_DECIDE || na == PETSC_CURRENT) na = 1000;
2556ae93cb3cSJason Sarich   if (!obj && !resid && !cnorm && !lits) {
25579566063dSJacob Faibussowitsch     PetscCall(PetscCalloc4(na, &obj, na, &resid, na, &cnorm, na, &lits));
2558ae93cb3cSJason Sarich     tao->hist_malloc = PETSC_TRUE;
2559ae93cb3cSJason Sarich   }
2560ae93cb3cSJason Sarich 
2561a7e14dcfSSatish Balay   tao->hist_obj   = obj;
2562a7e14dcfSSatish Balay   tao->hist_resid = resid;
2563a7e14dcfSSatish Balay   tao->hist_cnorm = cnorm;
2564ae93cb3cSJason Sarich   tao->hist_lits  = lits;
2565a7e14dcfSSatish Balay   tao->hist_max   = na;
2566a7e14dcfSSatish Balay   tao->hist_reset = reset;
2567ae93cb3cSJason Sarich   tao->hist_len   = 0;
25683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2569a7e14dcfSSatish Balay }
2570a7e14dcfSSatish Balay 
2571a7e14dcfSSatish Balay /*@C
257265ba42b6SBarry Smith   TaoGetConvergenceHistory - Gets the arrays used that hold the convergence history.
2573a7e14dcfSSatish Balay 
2574c3339decSBarry Smith   Collective
2575a7e14dcfSSatish Balay 
2576a7e14dcfSSatish Balay   Input Parameter:
257747450a7bSBarry Smith . tao - the `Tao` context
2578a7e14dcfSSatish Balay 
2579a7e14dcfSSatish Balay   Output Parameters:
2580a7e14dcfSSatish Balay + obj   - array used to hold objective value history
2581a7e14dcfSSatish Balay . resid - array used to hold residual history
2582a7e14dcfSSatish Balay . cnorm - array used to hold constraint violation history
2583ae93cb3cSJason Sarich . lits  - integer array used to hold linear solver iteration count
258467be906fSBarry Smith - nhist - size of `obj`, `resid`, `cnorm`, and `lits`
258567be906fSBarry Smith 
258667be906fSBarry Smith   Level: advanced
2587a7e14dcfSSatish Balay 
2588a7e14dcfSSatish Balay   Notes:
258965ba42b6SBarry Smith   This routine must be preceded by calls to `TaoSetConvergenceHistory()`
259065ba42b6SBarry Smith   and `TaoSolve()`, otherwise it returns useless information.
2591ae93cb3cSJason Sarich 
2592a7e14dcfSSatish Balay   This routine is useful, e.g., when running a code for purposes
2593a7e14dcfSSatish Balay   of accurate performance monitoring, when no I/O should be done
2594a7e14dcfSSatish Balay   during the section of code that is being timed.
2595a7e14dcfSSatish Balay 
2596e056e8ceSJacob Faibussowitsch   Fortran Notes:
259767be906fSBarry Smith   The calling sequence is
259867be906fSBarry Smith .vb
259967be906fSBarry Smith    call TaoGetConvergenceHistory(Tao tao, PetscInt nhist, PetscErrorCode ierr)
260067be906fSBarry Smith .ve
2601a3b724e8SBarry Smith   In other words this gets the current number of entries in the history. Access the history through the array you passed to `TaoSetConvergenceHistory()`
2602a7e14dcfSSatish Balay 
26031cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoSolve()`, `TaoSetConvergenceHistory()`
2604a7e14dcfSSatish Balay @*/
2605d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoGetConvergenceHistory(Tao tao, PetscReal **obj, PetscReal **resid, PetscReal **cnorm, PetscInt **lits, PetscInt *nhist)
2606d71ae5a4SJacob Faibussowitsch {
2607a7e14dcfSSatish Balay   PetscFunctionBegin;
2608441846f8SBarry Smith   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
2609a7e14dcfSSatish Balay   if (obj) *obj = tao->hist_obj;
2610a7e14dcfSSatish Balay   if (cnorm) *cnorm = tao->hist_cnorm;
2611a7e14dcfSSatish Balay   if (resid) *resid = tao->hist_resid;
26125b494b05SBarry Smith   if (lits) *lits = tao->hist_lits;
2613a7e14dcfSSatish Balay   if (nhist) *nhist = tao->hist_len;
26143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2615a7e14dcfSSatish Balay }
2616a7e14dcfSSatish Balay 
2617a7e14dcfSSatish Balay /*@
261847450a7bSBarry Smith   TaoSetApplicationContext - Sets the optional user-defined context for a `Tao` solver.
2619a7e14dcfSSatish Balay 
2620c3339decSBarry Smith   Logically Collective
2621a7e14dcfSSatish Balay 
2622a7e14dcfSSatish Balay   Input Parameters:
262347450a7bSBarry Smith + tao  - the `Tao` context
2624a7e14dcfSSatish Balay - usrP - optional user context
2625a7e14dcfSSatish Balay 
2626a7e14dcfSSatish Balay   Level: intermediate
2627a7e14dcfSSatish Balay 
262842747ad1SJacob Faibussowitsch .seealso: [](ch_tao), `Tao`, `TaoGetApplicationContext()`
2629a7e14dcfSSatish Balay @*/
2630d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoSetApplicationContext(Tao tao, void *usrP)
2631d71ae5a4SJacob Faibussowitsch {
2632a7e14dcfSSatish Balay   PetscFunctionBegin;
2633441846f8SBarry Smith   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
2634a7e14dcfSSatish Balay   tao->user = usrP;
26353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2636a7e14dcfSSatish Balay }
2637a7e14dcfSSatish Balay 
2638a7e14dcfSSatish Balay /*@
263947450a7bSBarry Smith   TaoGetApplicationContext - Gets the user-defined context for a `Tao` solver
2640a7e14dcfSSatish Balay 
2641a7e14dcfSSatish Balay   Not Collective
2642a7e14dcfSSatish Balay 
2643a7e14dcfSSatish Balay   Input Parameter:
264447450a7bSBarry Smith . tao - the `Tao` context
2645a7e14dcfSSatish Balay 
2646a7e14dcfSSatish Balay   Output Parameter:
2647a7e14dcfSSatish Balay . usrP - user context
2648a7e14dcfSSatish Balay 
2649a7e14dcfSSatish Balay   Level: intermediate
2650a7e14dcfSSatish Balay 
26511cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoSetApplicationContext()`
2652a7e14dcfSSatish Balay @*/
2653d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoGetApplicationContext(Tao tao, void *usrP)
2654d71ae5a4SJacob Faibussowitsch {
2655a7e14dcfSSatish Balay   PetscFunctionBegin;
2656441846f8SBarry Smith   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
26574f572ea9SToby Isaac   PetscAssertPointer(usrP, 2);
2658a7e14dcfSSatish Balay   *(void **)usrP = tao->user;
26593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2660a7e14dcfSSatish Balay }
2661a9603a14SPatrick Farrell 
2662a9603a14SPatrick Farrell /*@
2663a81c3919SBarry Smith   TaoSetGradientNorm - Sets the matrix used to define the norm that measures the size of the gradient in some of the `Tao` algorithms
2664a9603a14SPatrick Farrell 
2665c3339decSBarry Smith   Collective
2666a9603a14SPatrick Farrell 
2667a9603a14SPatrick Farrell   Input Parameters:
266847450a7bSBarry Smith + tao - the `Tao` context
266965ba42b6SBarry Smith - M   - matrix that defines the norm
2670a9603a14SPatrick Farrell 
2671a9603a14SPatrick Farrell   Level: beginner
2672a9603a14SPatrick Farrell 
26731cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoGetGradientNorm()`, `TaoGradientNorm()`
2674a9603a14SPatrick Farrell @*/
2675d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoSetGradientNorm(Tao tao, Mat M)
2676d71ae5a4SJacob Faibussowitsch {
2677a9603a14SPatrick Farrell   PetscFunctionBegin;
2678a9603a14SPatrick Farrell   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
267994511df7SStefano Zampini   PetscValidHeaderSpecific(M, MAT_CLASSID, 2);
26809566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)M));
26819566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&tao->gradient_norm));
26829566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&tao->gradient_norm_tmp));
2683a9603a14SPatrick Farrell   tao->gradient_norm = M;
26849566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(M, NULL, &tao->gradient_norm_tmp));
26853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2686a9603a14SPatrick Farrell }
2687a9603a14SPatrick Farrell 
2688a9603a14SPatrick Farrell /*@
2689a81c3919SBarry Smith   TaoGetGradientNorm - Returns the matrix used to define the norm used for measuring the size of the gradient in some of the `Tao` algorithms
2690a9603a14SPatrick Farrell 
2691a9603a14SPatrick Farrell   Not Collective
2692a9603a14SPatrick Farrell 
2693a9603a14SPatrick Farrell   Input Parameter:
269447450a7bSBarry Smith . tao - the `Tao` context
2695a9603a14SPatrick Farrell 
2696a9603a14SPatrick Farrell   Output Parameter:
2697a9603a14SPatrick Farrell . M - gradient norm
2698a9603a14SPatrick Farrell 
2699a9603a14SPatrick Farrell   Level: beginner
2700a9603a14SPatrick Farrell 
27011cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoSetGradientNorm()`, `TaoGradientNorm()`
2702a9603a14SPatrick Farrell @*/
2703d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoGetGradientNorm(Tao tao, Mat *M)
2704d71ae5a4SJacob Faibussowitsch {
2705a9603a14SPatrick Farrell   PetscFunctionBegin;
2706a9603a14SPatrick Farrell   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
27074f572ea9SToby Isaac   PetscAssertPointer(M, 2);
2708a9603a14SPatrick Farrell   *M = tao->gradient_norm;
27093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2710a9603a14SPatrick Farrell }
2711a9603a14SPatrick Farrell 
2712cc4c1da9SBarry Smith /*@
271347450a7bSBarry Smith   TaoGradientNorm - Compute the norm using the `NormType`, the user has selected
2714a9603a14SPatrick Farrell 
2715c3339decSBarry Smith   Collective
2716a9603a14SPatrick Farrell 
2717d8d19677SJose E. Roman   Input Parameters:
271847450a7bSBarry Smith + tao      - the `Tao` context
2719a81c3919SBarry Smith . gradient - the gradient
27203b242c63SJacob Faibussowitsch - type     - the norm type
2721a9603a14SPatrick Farrell 
2722a9603a14SPatrick Farrell   Output Parameter:
2723a9603a14SPatrick Farrell . gnorm - the gradient norm
2724a9603a14SPatrick Farrell 
272547450a7bSBarry Smith   Level: advanced
2726a9603a14SPatrick Farrell 
2727a81c3919SBarry Smith   Note:
2728a81c3919SBarry Smith   If `TaoSetGradientNorm()` has been set and `type` is `NORM_2` then the norm provided with `TaoSetGradientNorm()` is used.
2729a81c3919SBarry Smith 
2730a81c3919SBarry Smith   Developer Notes:
2731a81c3919SBarry Smith   Should be named `TaoComputeGradientNorm()`.
2732a81c3919SBarry Smith 
2733a81c3919SBarry Smith   The usage is a bit confusing, with `TaoSetGradientNorm()` plus `NORM_2` resulting in the computation of the user provided
2734a81c3919SBarry Smith   norm, perhaps a refactorization is in order.
2735a81c3919SBarry Smith 
27361cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoSetGradientNorm()`, `TaoGetGradientNorm()`
2737a9603a14SPatrick Farrell @*/
2738d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoGradientNorm(Tao tao, Vec gradient, NormType type, PetscReal *gnorm)
2739d71ae5a4SJacob Faibussowitsch {
2740a9603a14SPatrick Farrell   PetscFunctionBegin;
27418854b543SStefano Zampini   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
27428854b543SStefano Zampini   PetscValidHeaderSpecific(gradient, VEC_CLASSID, 2);
27438854b543SStefano Zampini   PetscValidLogicalCollectiveEnum(tao, type, 3);
27444f572ea9SToby Isaac   PetscAssertPointer(gnorm, 4);
2745a9603a14SPatrick Farrell   if (tao->gradient_norm) {
2746680e2bc7SPatrick Farrell     PetscScalar gnorms;
2747680e2bc7SPatrick Farrell 
27483c859ba3SBarry Smith     PetscCheck(type == NORM_2, PetscObjectComm((PetscObject)gradient), PETSC_ERR_ARG_WRONG, "Norm type must be NORM_2 if an inner product for the gradient norm is set.");
27499566063dSJacob Faibussowitsch     PetscCall(MatMult(tao->gradient_norm, gradient, tao->gradient_norm_tmp));
27509566063dSJacob Faibussowitsch     PetscCall(VecDot(gradient, tao->gradient_norm_tmp, &gnorms));
27517bb79937SPatrick Farrell     *gnorm = PetscRealPart(PetscSqrtScalar(gnorms));
2752a9603a14SPatrick Farrell   } else {
27539566063dSJacob Faibussowitsch     PetscCall(VecNorm(gradient, type, gnorm));
2754a9603a14SPatrick Farrell   }
27553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2756a9603a14SPatrick Farrell }
2757a9603a14SPatrick Farrell 
2758e882e171SHong Zhang /*@C
2759a81c3919SBarry Smith   TaoMonitorDrawCtxCreate - Creates the monitor context for `TaoMonitorSolutionDraw()`
2760a9603a14SPatrick Farrell 
2761c3339decSBarry Smith   Collective
2762e882e171SHong Zhang 
27632fe279fdSBarry Smith   Input Parameters:
27642fe279fdSBarry Smith + comm     - the communicator to share the context
27652fe279fdSBarry Smith . host     - the name of the X Windows host that will display the monitor
27662fe279fdSBarry Smith . label    - the label to put at the top of the display window
27672fe279fdSBarry Smith . x        - the horizontal coordinate of the lower left corner of the window to open
27682fe279fdSBarry Smith . y        - the vertical coordinate of the lower left corner of the window to open
27692fe279fdSBarry Smith . m        - the width of the window
27702fe279fdSBarry Smith . n        - the height of the window
27712fe279fdSBarry Smith - howoften - how many `Tao` iterations between displaying the monitor information
27722fe279fdSBarry Smith 
2773d5b43468SJose E. Roman   Output Parameter:
2774e882e171SHong Zhang . ctx - the monitor context
2775e882e171SHong Zhang 
2776a81c3919SBarry Smith   Options Database Keys:
2777a81c3919SBarry Smith + -tao_monitor_solution_draw - use `TaoMonitorSolutionDraw()` to monitor the solution
2778a81c3919SBarry Smith - -tao_draw_solution_initial - show initial guess as well as current solution
2779e882e171SHong Zhang 
2780e882e171SHong Zhang   Level: intermediate
2781e882e171SHong Zhang 
2782a81c3919SBarry Smith   Note:
2783a81c3919SBarry Smith   The context this creates, along with `TaoMonitorSolutionDraw()`, and `TaoMonitorDrawCtxDestroy()`
2784a81c3919SBarry Smith   are passed to `TaoMonitorSet()`.
2785a81c3919SBarry Smith 
27861cc06b55SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoMonitorSet()`, `TaoMonitorDefault()`, `VecView()`, `TaoMonitorDrawCtx()`
2787e882e171SHong Zhang @*/
2788d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoMonitorDrawCtxCreate(MPI_Comm comm, const char host[], const char label[], int x, int y, int m, int n, PetscInt howoften, TaoMonitorDrawCtx *ctx)
2789d71ae5a4SJacob Faibussowitsch {
2790e882e171SHong Zhang   PetscFunctionBegin;
27919566063dSJacob Faibussowitsch   PetscCall(PetscNew(ctx));
27929566063dSJacob Faibussowitsch   PetscCall(PetscViewerDrawOpen(comm, host, label, x, y, m, n, &(*ctx)->viewer));
27939566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetFromOptions((*ctx)->viewer));
2794e882e171SHong Zhang   (*ctx)->howoften = howoften;
27953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2796e882e171SHong Zhang }
2797e882e171SHong Zhang 
2798e882e171SHong Zhang /*@C
2799a81c3919SBarry Smith   TaoMonitorDrawCtxDestroy - Destroys the monitor context for `TaoMonitorSolutionDraw()`
2800e882e171SHong Zhang 
2801c3339decSBarry Smith   Collective
2802e882e171SHong Zhang 
280347450a7bSBarry Smith   Input Parameter:
2804e056e8ceSJacob Faibussowitsch . ictx - the monitor context
2805e882e171SHong Zhang 
2806e882e171SHong Zhang   Level: intermediate
2807e882e171SHong Zhang 
2808a81c3919SBarry Smith   Note:
2809a81c3919SBarry Smith   This is passed to `TaoMonitorSet()` as the final argument, along with `TaoMonitorSolutionDraw()`, and the context
2810a81c3919SBarry Smith   obtained with `TaoMonitorDrawCtxCreate()`.
2811a81c3919SBarry Smith 
2812a81c3919SBarry Smith .seealso: [](ch_tao), `Tao`, `TaoMonitorSet()`, `TaoMonitorDefault()`, `VecView()`, `TaoMonitorSolutionDraw()`
2813e882e171SHong Zhang @*/
2814d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoMonitorDrawCtxDestroy(TaoMonitorDrawCtx *ictx)
2815d71ae5a4SJacob Faibussowitsch {
2816e882e171SHong Zhang   PetscFunctionBegin;
28179566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&(*ictx)->viewer));
28189566063dSJacob Faibussowitsch   PetscCall(PetscFree(*ictx));
28193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2820e882e171SHong Zhang }
2821