xref: /petsc/src/tao/bound/impls/tron/tron.c (revision 4dfa11a44d5adf2389f1d3acbc8f3c1116dc6c3a)
1aaa7dc30SBarry Smith #include <../src/tao/bound/impls/tron/tron.h>
2aaa7dc30SBarry Smith #include <../src/tao/matrix/submatfree.h>
3a7e14dcfSSatish Balay 
4a7e14dcfSSatish Balay /* TRON Routines */
5441846f8SBarry Smith static PetscErrorCode TronGradientProjections(Tao, TAO_TRON *);
6a7e14dcfSSatish Balay /*------------------------------------------------------------*/
79371c9d4SSatish Balay static PetscErrorCode TaoDestroy_TRON(Tao tao) {
8a7e14dcfSSatish Balay   TAO_TRON *tron = (TAO_TRON *)tao->data;
9a7e14dcfSSatish Balay 
10a7e14dcfSSatish Balay   PetscFunctionBegin;
119566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&tron->X_New));
129566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&tron->G_New));
139566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&tron->Work));
149566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&tron->DXFree));
159566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&tron->R));
169566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&tron->diag));
179566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&tron->scatter));
189566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&tron->Free_Local));
199566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&tron->H_sub));
209566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&tron->Hpre_sub));
21a958fbfcSStefano Zampini   PetscCall(KSPDestroy(&tao->ksp));
229566063dSJacob Faibussowitsch   PetscCall(PetscFree(tao->data));
23a7e14dcfSSatish Balay   PetscFunctionReturn(0);
24a7e14dcfSSatish Balay }
25a7e14dcfSSatish Balay 
26a7e14dcfSSatish Balay /*------------------------------------------------------------*/
279371c9d4SSatish Balay static PetscErrorCode TaoSetFromOptions_TRON(Tao tao, PetscOptionItems *PetscOptionsObject) {
28a7e14dcfSSatish Balay   TAO_TRON *tron = (TAO_TRON *)tao->data;
29a7e14dcfSSatish Balay   PetscBool flg;
30a7e14dcfSSatish Balay 
31a7e14dcfSSatish Balay   PetscFunctionBegin;
32d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Newton Trust Region Method for bound constrained optimization");
339566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-tao_tron_maxgpits", "maximum number of gradient projections per TRON iterate", "TaoSetMaxGPIts", tron->maxgpits, &tron->maxgpits, &flg));
34d0609cedSBarry Smith   PetscOptionsHeadEnd();
359566063dSJacob Faibussowitsch   PetscCall(KSPSetFromOptions(tao->ksp));
36a7e14dcfSSatish Balay   PetscFunctionReturn(0);
37a7e14dcfSSatish Balay }
38a7e14dcfSSatish Balay 
39a7e14dcfSSatish Balay /*------------------------------------------------------------*/
409371c9d4SSatish Balay static PetscErrorCode TaoView_TRON(Tao tao, PetscViewer viewer) {
41a7e14dcfSSatish Balay   TAO_TRON *tron = (TAO_TRON *)tao->data;
42a7e14dcfSSatish Balay   PetscBool isascii;
43a7e14dcfSSatish Balay 
44a7e14dcfSSatish Balay   PetscFunctionBegin;
459566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
46a7e14dcfSSatish Balay   if (isascii) {
4763a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Total PG its: %" PetscInt_FMT ",", tron->total_gp_its));
489566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "PG tolerance: %g \n", (double)tron->pg_ftol));
49a7e14dcfSSatish Balay   }
50a7e14dcfSSatish Balay   PetscFunctionReturn(0);
51a7e14dcfSSatish Balay }
52a7e14dcfSSatish Balay 
53a7e14dcfSSatish Balay /* ---------------------------------------------------------- */
549371c9d4SSatish Balay static PetscErrorCode TaoSetup_TRON(Tao tao) {
55a7e14dcfSSatish Balay   TAO_TRON *tron = (TAO_TRON *)tao->data;
56a7e14dcfSSatish Balay 
57a7e14dcfSSatish Balay   PetscFunctionBegin;
58a7e14dcfSSatish Balay   /* Allocate some arrays */
599566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tao->solution, &tron->diag));
609566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tao->solution, &tron->X_New));
619566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tao->solution, &tron->G_New));
629566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tao->solution, &tron->Work));
639566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tao->solution, &tao->gradient));
649566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tao->solution, &tao->stepdirection));
65a7e14dcfSSatish Balay   PetscFunctionReturn(0);
66a7e14dcfSSatish Balay }
67a7e14dcfSSatish Balay 
689371c9d4SSatish Balay static PetscErrorCode TaoSolve_TRON(Tao tao) {
69a7e14dcfSSatish Balay   TAO_TRON                    *tron = (TAO_TRON *)tao->data;
708931d482SJason Sarich   PetscInt                     its;
71e4cb33bbSBarry Smith   TaoLineSearchConvergedReason ls_reason = TAOLINESEARCH_CONTINUE_ITERATING;
72a7e14dcfSSatish Balay   PetscReal                    prered, actred, delta, f, f_new, rhok, gdx, xdiff, stepsize;
7347a47007SBarry Smith 
74a7e14dcfSSatish Balay   PetscFunctionBegin;
75a7e14dcfSSatish Balay   tron->pgstepsize = 1.0;
76a7e14dcfSSatish Balay   tao->trust       = tao->trust0;
77a7e14dcfSSatish Balay   /*   Project the current point onto the feasible set */
789566063dSJacob Faibussowitsch   PetscCall(TaoComputeVariableBounds(tao));
799566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetVariableBounds(tao->linesearch, tao->XL, tao->XU));
80a7e14dcfSSatish Balay 
81ce902467SBarry Smith   /* Project the initial point onto the feasible region */
829566063dSJacob Faibussowitsch   PetscCall(VecMedian(tao->XL, tao->solution, tao->XU, tao->solution));
8353506e15SBarry Smith 
84ce902467SBarry Smith   /* Compute the objective function and gradient */
859566063dSJacob Faibussowitsch   PetscCall(TaoComputeObjectiveAndGradient(tao, tao->solution, &tron->f, tao->gradient));
869566063dSJacob Faibussowitsch   PetscCall(VecNorm(tao->gradient, NORM_2, &tron->gnorm));
873c859ba3SBarry Smith   PetscCheck(!PetscIsInfOrNanReal(tron->f) && !PetscIsInfOrNanReal(tron->gnorm), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute function generated Inf or NaN");
88a7e14dcfSSatish Balay 
89a7e14dcfSSatish Balay   /* Project the gradient and calculate the norm */
909566063dSJacob Faibussowitsch   PetscCall(VecBoundGradientProjection(tao->gradient, tao->solution, tao->XL, tao->XU, tao->gradient));
919566063dSJacob Faibussowitsch   PetscCall(VecNorm(tao->gradient, NORM_2, &tron->gnorm));
92a7e14dcfSSatish Balay 
93ce902467SBarry Smith   /* Initialize trust region radius */
94ce902467SBarry Smith   tao->trust = tao->trust0;
95ad540459SPierre Jolivet   if (tao->trust <= 0) tao->trust = PetscMax(tron->gnorm * tron->gnorm, 1.0);
96a7e14dcfSSatish Balay 
97ce902467SBarry Smith   /* Initialize step sizes for the line searches */
98ce902467SBarry Smith   tron->pgstepsize = 1.0;
99a7e14dcfSSatish Balay   tron->stepsize   = tao->trust;
100ce902467SBarry Smith 
1013ecd9318SAlp Dener   tao->reason = TAO_CONTINUE_ITERATING;
1029566063dSJacob Faibussowitsch   PetscCall(TaoLogConvergenceHistory(tao, tron->f, tron->gnorm, 0.0, tao->ksp_its));
1039566063dSJacob Faibussowitsch   PetscCall(TaoMonitor(tao, tao->niter, tron->f, tron->gnorm, 0.0, tron->stepsize));
104dbbe0bcdSBarry Smith   PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
1053ecd9318SAlp Dener   while (tao->reason == TAO_CONTINUE_ITERATING) {
106e1e80dc8SAlp Dener     /* Call general purpose update function */
107dbbe0bcdSBarry Smith     PetscTryTypeMethod(tao, update, tao->niter, tao->user_update);
108ce902467SBarry Smith 
109ce902467SBarry Smith     /* Perform projected gradient iterations */
1109566063dSJacob Faibussowitsch     PetscCall(TronGradientProjections(tao, tron));
111ce902467SBarry Smith 
1129566063dSJacob Faibussowitsch     PetscCall(VecBoundGradientProjection(tao->gradient, tao->solution, tao->XL, tao->XU, tao->gradient));
1139566063dSJacob Faibussowitsch     PetscCall(VecNorm(tao->gradient, NORM_2, &tron->gnorm));
114ce902467SBarry Smith 
115ce902467SBarry Smith     tao->ksp_its      = 0;
1169371c9d4SSatish Balay     f                 = tron->f;
1179371c9d4SSatish Balay     delta             = tao->trust;
118a7e14dcfSSatish Balay     tron->n_free_last = tron->n_free;
1199566063dSJacob Faibussowitsch     PetscCall(TaoComputeHessian(tao, tao->solution, tao->hessian, tao->hessian_pre));
120a7e14dcfSSatish Balay 
121baa3a814SAlp Dener     /* Generate index set (IS) of which bound constraints are active */
1229566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&tron->Free_Local));
1239566063dSJacob Faibussowitsch     PetscCall(VecWhichInactive(tao->XL, tao->solution, tao->gradient, tao->XU, PETSC_TRUE, &tron->Free_Local));
1249566063dSJacob Faibussowitsch     PetscCall(ISGetSize(tron->Free_Local, &tron->n_free));
125a7e14dcfSSatish Balay 
126a7e14dcfSSatish Balay     /* If no free variables */
127a7e14dcfSSatish Balay     if (tron->n_free == 0) {
1289566063dSJacob Faibussowitsch       PetscCall(VecNorm(tao->gradient, NORM_2, &tron->gnorm));
1299566063dSJacob Faibussowitsch       PetscCall(TaoLogConvergenceHistory(tao, tron->f, tron->gnorm, 0.0, tao->ksp_its));
1309566063dSJacob Faibussowitsch       PetscCall(TaoMonitor(tao, tao->niter, tron->f, tron->gnorm, 0.0, delta));
131dbbe0bcdSBarry Smith       PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
132ad540459SPierre Jolivet       if (!tao->reason) tao->reason = TAO_CONVERGED_STEPTOL;
133a7e14dcfSSatish Balay       break;
134a7e14dcfSSatish Balay     }
135a7e14dcfSSatish Balay     /* use free_local to mask/submat gradient, hessian, stepdirection */
1369566063dSJacob Faibussowitsch     PetscCall(TaoVecGetSubVec(tao->gradient, tron->Free_Local, tao->subset_type, 0.0, &tron->R));
1379566063dSJacob Faibussowitsch     PetscCall(TaoVecGetSubVec(tao->gradient, tron->Free_Local, tao->subset_type, 0.0, &tron->DXFree));
1389566063dSJacob Faibussowitsch     PetscCall(VecSet(tron->DXFree, 0.0));
1399566063dSJacob Faibussowitsch     PetscCall(VecScale(tron->R, -1.0));
1409566063dSJacob Faibussowitsch     PetscCall(TaoMatGetSubMat(tao->hessian, tron->Free_Local, tron->diag, tao->subset_type, &tron->H_sub));
141a7e14dcfSSatish Balay     if (tao->hessian == tao->hessian_pre) {
1429566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&tron->Hpre_sub));
1439566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)(tron->H_sub)));
144a7e14dcfSSatish Balay       tron->Hpre_sub = tron->H_sub;
145a7e14dcfSSatish Balay     } else {
1469566063dSJacob Faibussowitsch       PetscCall(TaoMatGetSubMat(tao->hessian_pre, tron->Free_Local, tron->diag, tao->subset_type, &tron->Hpre_sub));
147a7e14dcfSSatish Balay     }
1489566063dSJacob Faibussowitsch     PetscCall(KSPReset(tao->ksp));
1499566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(tao->ksp, tron->H_sub, tron->Hpre_sub));
150a7e14dcfSSatish Balay     while (1) {
151a7e14dcfSSatish Balay       /* Approximately solve the reduced linear system */
1529566063dSJacob Faibussowitsch       PetscCall(KSPCGSetRadius(tao->ksp, delta));
153d8ec8e84SJason Sarich 
1549566063dSJacob Faibussowitsch       PetscCall(KSPSolve(tao->ksp, tron->R, tron->DXFree));
1559566063dSJacob Faibussowitsch       PetscCall(KSPGetIterationNumber(tao->ksp, &its));
156a7e14dcfSSatish Balay       tao->ksp_its += its;
157ae93cb3cSJason Sarich       tao->ksp_tot_its += its;
1589566063dSJacob Faibussowitsch       PetscCall(VecSet(tao->stepdirection, 0.0));
159a7e14dcfSSatish Balay 
160a7e14dcfSSatish Balay       /* Add dxfree matrix to compute step direction vector */
1619566063dSJacob Faibussowitsch       PetscCall(VecISAXPY(tao->stepdirection, tron->Free_Local, 1.0, tron->DXFree));
162a7e14dcfSSatish Balay 
1639566063dSJacob Faibussowitsch       PetscCall(VecDot(tao->gradient, tao->stepdirection, &gdx));
1649566063dSJacob Faibussowitsch       PetscCall(VecCopy(tao->solution, tron->X_New));
1659566063dSJacob Faibussowitsch       PetscCall(VecCopy(tao->gradient, tron->G_New));
166a7e14dcfSSatish Balay 
1679371c9d4SSatish Balay       stepsize = 1.0;
1689371c9d4SSatish Balay       f_new    = f;
169a7e14dcfSSatish Balay 
1709566063dSJacob Faibussowitsch       PetscCall(TaoLineSearchSetInitialStepLength(tao->linesearch, 1.0));
1719566063dSJacob Faibussowitsch       PetscCall(TaoLineSearchApply(tao->linesearch, tron->X_New, &f_new, tron->G_New, tao->stepdirection, &stepsize, &ls_reason));
1729566063dSJacob Faibussowitsch       PetscCall(TaoAddLineSearchCounts(tao));
173a7e14dcfSSatish Balay 
1749566063dSJacob Faibussowitsch       PetscCall(MatMult(tao->hessian, tao->stepdirection, tron->Work));
1759566063dSJacob Faibussowitsch       PetscCall(VecAYPX(tron->Work, 0.5, tao->gradient));
1769566063dSJacob Faibussowitsch       PetscCall(VecDot(tao->stepdirection, tron->Work, &prered));
177a7e14dcfSSatish Balay       actred = f_new - f;
178ce902467SBarry Smith       if ((PetscAbsScalar(actred) <= 1e-6) && (PetscAbsScalar(prered) <= 1e-6)) {
179ce902467SBarry Smith         rhok = 1.0;
180ce902467SBarry Smith       } else if (actred < 0) {
181a7e14dcfSSatish Balay         rhok = PetscAbs(-actred / prered);
182a7e14dcfSSatish Balay       } else {
183a7e14dcfSSatish Balay         rhok = 0.0;
184a7e14dcfSSatish Balay       }
185a7e14dcfSSatish Balay 
186a7e14dcfSSatish Balay       /* Compare actual improvement to the quadratic model */
187a7e14dcfSSatish Balay       if (rhok > tron->eta1) { /* Accept the point */
188a7e14dcfSSatish Balay         /* d = x_new - x */
1899566063dSJacob Faibussowitsch         PetscCall(VecCopy(tron->X_New, tao->stepdirection));
1909566063dSJacob Faibussowitsch         PetscCall(VecAXPY(tao->stepdirection, -1.0, tao->solution));
191a7e14dcfSSatish Balay 
1929566063dSJacob Faibussowitsch         PetscCall(VecNorm(tao->stepdirection, NORM_2, &xdiff));
193a7e14dcfSSatish Balay         xdiff *= stepsize;
194a7e14dcfSSatish Balay 
195a7e14dcfSSatish Balay         /* Adjust trust region size */
196a7e14dcfSSatish Balay         if (rhok < tron->eta2) {
197a7e14dcfSSatish Balay           delta = PetscMin(xdiff, delta) * tron->sigma1;
198a7e14dcfSSatish Balay         } else if (rhok > tron->eta4) {
199a7e14dcfSSatish Balay           delta = PetscMin(xdiff, delta) * tron->sigma3;
200a7e14dcfSSatish Balay         } else if (rhok > tron->eta3) {
201a7e14dcfSSatish Balay           delta = PetscMin(xdiff, delta) * tron->sigma2;
202a7e14dcfSSatish Balay         }
2039566063dSJacob Faibussowitsch         PetscCall(VecBoundGradientProjection(tron->G_New, tron->X_New, tao->XL, tao->XU, tao->gradient));
2049566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&tron->Free_Local));
2059566063dSJacob Faibussowitsch         PetscCall(VecWhichInactive(tao->XL, tron->X_New, tao->gradient, tao->XU, PETSC_TRUE, &tron->Free_Local));
206a7e14dcfSSatish Balay         f = f_new;
2079566063dSJacob Faibussowitsch         PetscCall(VecNorm(tao->gradient, NORM_2, &tron->gnorm));
2089566063dSJacob Faibussowitsch         PetscCall(VecCopy(tron->X_New, tao->solution));
2099566063dSJacob Faibussowitsch         PetscCall(VecCopy(tron->G_New, tao->gradient));
210a7e14dcfSSatish Balay         break;
2119371c9d4SSatish Balay       } else if (delta <= 1e-30) {
212a7e14dcfSSatish Balay         break;
2139371c9d4SSatish Balay       } else {
214a7e14dcfSSatish Balay         delta /= 4.0;
215a7e14dcfSSatish Balay       }
216a7e14dcfSSatish Balay     } /* end linear solve loop */
217a7e14dcfSSatish Balay 
2189371c9d4SSatish Balay     tron->f      = f;
2199371c9d4SSatish Balay     tron->actred = actred;
2209371c9d4SSatish Balay     tao->trust   = delta;
2218931d482SJason Sarich     tao->niter++;
2229566063dSJacob Faibussowitsch     PetscCall(TaoLogConvergenceHistory(tao, tron->f, tron->gnorm, 0.0, tao->ksp_its));
2239566063dSJacob Faibussowitsch     PetscCall(TaoMonitor(tao, tao->niter, tron->f, tron->gnorm, 0.0, stepsize));
224dbbe0bcdSBarry Smith     PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
225a7e14dcfSSatish Balay   } /* END MAIN LOOP  */
226a7e14dcfSSatish Balay   PetscFunctionReturn(0);
227a7e14dcfSSatish Balay }
228a7e14dcfSSatish Balay 
2299371c9d4SSatish Balay static PetscErrorCode TronGradientProjections(Tao tao, TAO_TRON *tron) {
230a7e14dcfSSatish Balay   PetscInt                     i;
231e4cb33bbSBarry Smith   TaoLineSearchConvergedReason ls_reason;
232a7e14dcfSSatish Balay   PetscReal                    actred = -1.0, actred_max = 0.0;
233a7e14dcfSSatish Balay   PetscReal                    f_new;
234a7e14dcfSSatish Balay   /*
235a7e14dcfSSatish Balay      The gradient and function value passed into and out of this
236a7e14dcfSSatish Balay      routine should be current and correct.
237a7e14dcfSSatish Balay 
238a7e14dcfSSatish Balay      The free, active, and binding variables should be already identified
239a7e14dcfSSatish Balay   */
240a7e14dcfSSatish Balay   PetscFunctionBegin;
241a7e14dcfSSatish Balay 
242ce902467SBarry Smith   for (i = 0; i < tron->maxgpits; ++i) {
243a7e14dcfSSatish Balay     if (-actred <= (tron->pg_ftol) * actred_max) break;
244a7e14dcfSSatish Balay 
245ce902467SBarry Smith     ++tron->gp_iterates;
246ce902467SBarry Smith     ++tron->total_gp_its;
247a7e14dcfSSatish Balay     f_new = tron->f;
248a7e14dcfSSatish Balay 
2499566063dSJacob Faibussowitsch     PetscCall(VecCopy(tao->gradient, tao->stepdirection));
2509566063dSJacob Faibussowitsch     PetscCall(VecScale(tao->stepdirection, -1.0));
2519566063dSJacob Faibussowitsch     PetscCall(TaoLineSearchSetInitialStepLength(tao->linesearch, tron->pgstepsize));
252d0609cedSBarry Smith     PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &f_new, tao->gradient, tao->stepdirection, &tron->pgstepsize, &ls_reason));
2539566063dSJacob Faibussowitsch     PetscCall(TaoAddLineSearchCounts(tao));
254a7e14dcfSSatish Balay 
2559566063dSJacob Faibussowitsch     PetscCall(VecBoundGradientProjection(tao->gradient, tao->solution, tao->XL, tao->XU, tao->gradient));
2569566063dSJacob Faibussowitsch     PetscCall(VecNorm(tao->gradient, NORM_2, &tron->gnorm));
257ce902467SBarry Smith 
258a7e14dcfSSatish Balay     /* Update the iterate */
259a7e14dcfSSatish Balay     actred     = f_new - tron->f;
260a7e14dcfSSatish Balay     actred_max = PetscMax(actred_max, -(f_new - tron->f));
261a7e14dcfSSatish Balay     tron->f    = f_new;
262a7e14dcfSSatish Balay   }
263a7e14dcfSSatish Balay   PetscFunctionReturn(0);
264a7e14dcfSSatish Balay }
265a7e14dcfSSatish Balay 
2669371c9d4SSatish Balay static PetscErrorCode TaoComputeDual_TRON(Tao tao, Vec DXL, Vec DXU) {
267a7e14dcfSSatish Balay   TAO_TRON *tron = (TAO_TRON *)tao->data;
268a7e14dcfSSatish Balay 
269a7e14dcfSSatish Balay   PetscFunctionBegin;
270441846f8SBarry Smith   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
271a7e14dcfSSatish Balay   PetscValidHeaderSpecific(DXL, VEC_CLASSID, 2);
272a7e14dcfSSatish Balay   PetscValidHeaderSpecific(DXU, VEC_CLASSID, 3);
2733c859ba3SBarry Smith   PetscCheck(tron->Work && tao->gradient, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Dual variables don't exist yet or no longer exist.");
274a7e14dcfSSatish Balay 
2759566063dSJacob Faibussowitsch   PetscCall(VecBoundGradientProjection(tao->gradient, tao->solution, tao->XL, tao->XU, tron->Work));
2769566063dSJacob Faibussowitsch   PetscCall(VecCopy(tron->Work, DXL));
2779566063dSJacob Faibussowitsch   PetscCall(VecAXPY(DXL, -1.0, tao->gradient));
2789566063dSJacob Faibussowitsch   PetscCall(VecSet(DXU, 0.0));
2799566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMax(DXL, DXL, DXU));
280a7e14dcfSSatish Balay 
2819566063dSJacob Faibussowitsch   PetscCall(VecCopy(tao->gradient, DXU));
2829566063dSJacob Faibussowitsch   PetscCall(VecAXPY(DXU, -1.0, tron->Work));
2839566063dSJacob Faibussowitsch   PetscCall(VecSet(tron->Work, 0.0));
2849566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMin(DXU, tron->Work, DXU));
285a7e14dcfSSatish Balay   PetscFunctionReturn(0);
286a7e14dcfSSatish Balay }
287a7e14dcfSSatish Balay 
288a7e14dcfSSatish Balay /*------------------------------------------------------------*/
2891522df2eSJason Sarich /*MC
2901522df2eSJason Sarich   TAOTRON - The TRON algorithm is an active-set Newton trust region method
2911522df2eSJason Sarich   for bound-constrained minimization.
2921522df2eSJason Sarich 
2931522df2eSJason Sarich   Options Database Keys:
2941522df2eSJason Sarich + -tao_tron_maxgpits - maximum number of gradient projections per TRON iterate
2951522df2eSJason Sarich - -tao_subset_type - "subvec","mask","matrix-free", strategies for handling active-sets
2961522df2eSJason Sarich 
2971eb8069cSJason Sarich   Level: beginner
2981522df2eSJason Sarich M*/
2999371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode TaoCreate_TRON(Tao tao) {
300a7e14dcfSSatish Balay   TAO_TRON   *tron;
3018caf6e8cSBarry Smith   const char *morethuente_type = TAOLINESEARCHMT;
302a7e14dcfSSatish Balay 
30347a47007SBarry Smith   PetscFunctionBegin;
304a7e14dcfSSatish Balay   tao->ops->setup          = TaoSetup_TRON;
305a7e14dcfSSatish Balay   tao->ops->solve          = TaoSolve_TRON;
306a7e14dcfSSatish Balay   tao->ops->view           = TaoView_TRON;
307a7e14dcfSSatish Balay   tao->ops->setfromoptions = TaoSetFromOptions_TRON;
308a7e14dcfSSatish Balay   tao->ops->destroy        = TaoDestroy_TRON;
309a7e14dcfSSatish Balay   tao->ops->computedual    = TaoComputeDual_TRON;
310a7e14dcfSSatish Balay 
311*4dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&tron));
3126f4723b1SBarry Smith   tao->data = (void *)tron;
3136552cf8aSJason Sarich 
3146552cf8aSJason Sarich   /* Override default settings (unless already changed) */
3156552cf8aSJason Sarich   if (!tao->max_it_changed) tao->max_it = 50;
3166552cf8aSJason Sarich   if (!tao->trust0_changed) tao->trust0 = 1.0;
3174f1535bcSTodd Munson   if (!tao->steptol_changed) tao->steptol = 0.0;
318a7e14dcfSSatish Balay 
319a7e14dcfSSatish Balay   /* Initialize pointers and variables */
320a7e14dcfSSatish Balay   tron->n        = 0;
321a7e14dcfSSatish Balay   tron->maxgpits = 3;
322a7e14dcfSSatish Balay   tron->pg_ftol  = 0.001;
323a7e14dcfSSatish Balay 
324a7e14dcfSSatish Balay   tron->eta1 = 1.0e-4;
325a7e14dcfSSatish Balay   tron->eta2 = 0.25;
326a7e14dcfSSatish Balay   tron->eta3 = 0.50;
327a7e14dcfSSatish Balay   tron->eta4 = 0.90;
328a7e14dcfSSatish Balay 
329a7e14dcfSSatish Balay   tron->sigma1 = 0.5;
330a7e14dcfSSatish Balay   tron->sigma2 = 2.0;
331a7e14dcfSSatish Balay   tron->sigma3 = 4.0;
332a7e14dcfSSatish Balay 
333a7e14dcfSSatish Balay   tron->gp_iterates  = 0; /* Cumulative number */
334a7e14dcfSSatish Balay   tron->total_gp_its = 0;
335a7e14dcfSSatish Balay   tron->n_free       = 0;
336a7e14dcfSSatish Balay 
3376c23d075SBarry Smith   tron->DXFree     = NULL;
3386c23d075SBarry Smith   tron->R          = NULL;
3396c23d075SBarry Smith   tron->X_New      = NULL;
3406c23d075SBarry Smith   tron->G_New      = NULL;
3416c23d075SBarry Smith   tron->Work       = NULL;
3426c23d075SBarry Smith   tron->Free_Local = NULL;
3436c23d075SBarry Smith   tron->H_sub      = NULL;
3446c23d075SBarry Smith   tron->Hpre_sub   = NULL;
345a7e14dcfSSatish Balay   tao->subset_type = TAO_SUBSET_SUBVEC;
346a7e14dcfSSatish Balay 
3479566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch));
3489566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1));
3499566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetType(tao->linesearch, morethuente_type));
3509566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchUseTaoRoutines(tao->linesearch, tao));
3519566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetOptionsPrefix(tao->linesearch, tao->hdr.prefix));
352a7e14dcfSSatish Balay 
3539566063dSJacob Faibussowitsch   PetscCall(KSPCreate(((PetscObject)tao)->comm, &tao->ksp));
3549566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1));
3559566063dSJacob Faibussowitsch   PetscCall(KSPSetOptionsPrefix(tao->ksp, tao->hdr.prefix));
3569566063dSJacob Faibussowitsch   PetscCall(KSPSetType(tao->ksp, KSPSTCG));
357a7e14dcfSSatish Balay   PetscFunctionReturn(0);
358a7e14dcfSSatish Balay }
359