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 *);
TaoDestroy_TRON(Tao tao)6d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoDestroy_TRON(Tao tao)
7d71ae5a4SJacob Faibussowitsch {
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));
233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
24a7e14dcfSSatish Balay }
25a7e14dcfSSatish Balay
TaoSetFromOptions_TRON(Tao tao,PetscOptionItems PetscOptionsObject)26ce78bad3SBarry Smith static PetscErrorCode TaoSetFromOptions_TRON(Tao tao, PetscOptionItems PetscOptionsObject)
27d71ae5a4SJacob Faibussowitsch {
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));
363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
37a7e14dcfSSatish Balay }
38a7e14dcfSSatish Balay
TaoView_TRON(Tao tao,PetscViewer viewer)39d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoView_TRON(Tao tao, PetscViewer viewer)
40d71ae5a4SJacob Faibussowitsch {
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 }
503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
51a7e14dcfSSatish Balay }
52a7e14dcfSSatish Balay
TaoSetup_TRON(Tao tao)53d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSetup_TRON(Tao tao)
54d71ae5a4SJacob Faibussowitsch {
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));
653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
66a7e14dcfSSatish Balay }
67a7e14dcfSSatish Balay
TaoSolve_TRON(Tao tao)68d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSolve_TRON(Tao tao)
69d71ae5a4SJacob Faibussowitsch {
70a7e14dcfSSatish Balay TAO_TRON *tron = (TAO_TRON *)tao->data;
718931d482SJason Sarich PetscInt its;
72e4cb33bbSBarry Smith TaoLineSearchConvergedReason ls_reason = TAOLINESEARCH_CONTINUE_ITERATING;
73a7e14dcfSSatish Balay PetscReal prered, actred, delta, f, f_new, rhok, gdx, xdiff, stepsize;
7447a47007SBarry Smith
75a7e14dcfSSatish Balay PetscFunctionBegin;
76a7e14dcfSSatish Balay tron->pgstepsize = 1.0;
77a7e14dcfSSatish Balay tao->trust = tao->trust0;
78a7e14dcfSSatish Balay /* Project the current point onto the feasible set */
799566063dSJacob Faibussowitsch PetscCall(TaoComputeVariableBounds(tao));
809566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetVariableBounds(tao->linesearch, tao->XL, tao->XU));
81a7e14dcfSSatish Balay
82ce902467SBarry Smith /* Project the initial point onto the feasible region */
839566063dSJacob Faibussowitsch PetscCall(VecMedian(tao->XL, tao->solution, tao->XU, tao->solution));
8453506e15SBarry Smith
85ce902467SBarry Smith /* Compute the objective function and gradient */
869566063dSJacob Faibussowitsch PetscCall(TaoComputeObjectiveAndGradient(tao, tao->solution, &tron->f, tao->gradient));
879566063dSJacob Faibussowitsch PetscCall(VecNorm(tao->gradient, NORM_2, &tron->gnorm));
88*76c63389SBarry Smith PetscCheck(!PetscIsInfOrNanReal(tron->f) && !PetscIsInfOrNanReal(tron->gnorm), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute function generated infinity or NaN");
89a7e14dcfSSatish Balay
90a7e14dcfSSatish Balay /* Project the gradient and calculate the norm */
919566063dSJacob Faibussowitsch PetscCall(VecBoundGradientProjection(tao->gradient, tao->solution, tao->XL, tao->XU, tao->gradient));
929566063dSJacob Faibussowitsch PetscCall(VecNorm(tao->gradient, NORM_2, &tron->gnorm));
93a7e14dcfSSatish Balay
94ce902467SBarry Smith /* Initialize trust region radius */
95ce902467SBarry Smith tao->trust = tao->trust0;
96ad540459SPierre Jolivet if (tao->trust <= 0) tao->trust = PetscMax(tron->gnorm * tron->gnorm, 1.0);
97a7e14dcfSSatish Balay
98ce902467SBarry Smith /* Initialize step sizes for the line searches */
99ce902467SBarry Smith tron->pgstepsize = 1.0;
100a7e14dcfSSatish Balay tron->stepsize = tao->trust;
101ce902467SBarry Smith
1023ecd9318SAlp Dener tao->reason = TAO_CONTINUE_ITERATING;
1039566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao, tron->f, tron->gnorm, 0.0, tao->ksp_its));
1049566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao, tao->niter, tron->f, tron->gnorm, 0.0, tron->stepsize));
105dbbe0bcdSBarry Smith PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
1063ecd9318SAlp Dener while (tao->reason == TAO_CONTINUE_ITERATING) {
107e1e80dc8SAlp Dener /* Call general purpose update function */
108270bebe6SStefano Zampini if (tao->ops->update) {
109270bebe6SStefano Zampini PetscUseTypeMethod(tao, update, tao->niter, tao->user_update);
110270bebe6SStefano Zampini PetscCall(TaoComputeObjective(tao, tao->solution, &tron->f));
111270bebe6SStefano Zampini }
112ce902467SBarry Smith
113ce902467SBarry Smith /* Perform projected gradient iterations */
1149566063dSJacob Faibussowitsch PetscCall(TronGradientProjections(tao, tron));
115ce902467SBarry Smith
1169566063dSJacob Faibussowitsch PetscCall(VecBoundGradientProjection(tao->gradient, tao->solution, tao->XL, tao->XU, tao->gradient));
1179566063dSJacob Faibussowitsch PetscCall(VecNorm(tao->gradient, NORM_2, &tron->gnorm));
118ce902467SBarry Smith
119ce902467SBarry Smith tao->ksp_its = 0;
1209371c9d4SSatish Balay f = tron->f;
1219371c9d4SSatish Balay delta = tao->trust;
122a7e14dcfSSatish Balay tron->n_free_last = tron->n_free;
1239566063dSJacob Faibussowitsch PetscCall(TaoComputeHessian(tao, tao->solution, tao->hessian, tao->hessian_pre));
124a7e14dcfSSatish Balay
125baa3a814SAlp Dener /* Generate index set (IS) of which bound constraints are active */
1269566063dSJacob Faibussowitsch PetscCall(ISDestroy(&tron->Free_Local));
1279566063dSJacob Faibussowitsch PetscCall(VecWhichInactive(tao->XL, tao->solution, tao->gradient, tao->XU, PETSC_TRUE, &tron->Free_Local));
1289566063dSJacob Faibussowitsch PetscCall(ISGetSize(tron->Free_Local, &tron->n_free));
129a7e14dcfSSatish Balay
130a7e14dcfSSatish Balay /* If no free variables */
131a7e14dcfSSatish Balay if (tron->n_free == 0) {
1329566063dSJacob Faibussowitsch PetscCall(VecNorm(tao->gradient, NORM_2, &tron->gnorm));
1339566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao, tron->f, tron->gnorm, 0.0, tao->ksp_its));
1349566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao, tao->niter, tron->f, tron->gnorm, 0.0, delta));
135dbbe0bcdSBarry Smith PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
136ad540459SPierre Jolivet if (!tao->reason) tao->reason = TAO_CONVERGED_STEPTOL;
137a7e14dcfSSatish Balay break;
138a7e14dcfSSatish Balay }
139a7e14dcfSSatish Balay /* use free_local to mask/submat gradient, hessian, stepdirection */
1409566063dSJacob Faibussowitsch PetscCall(TaoVecGetSubVec(tao->gradient, tron->Free_Local, tao->subset_type, 0.0, &tron->R));
1419566063dSJacob Faibussowitsch PetscCall(TaoVecGetSubVec(tao->gradient, tron->Free_Local, tao->subset_type, 0.0, &tron->DXFree));
1429566063dSJacob Faibussowitsch PetscCall(VecSet(tron->DXFree, 0.0));
1439566063dSJacob Faibussowitsch PetscCall(VecScale(tron->R, -1.0));
1449566063dSJacob Faibussowitsch PetscCall(TaoMatGetSubMat(tao->hessian, tron->Free_Local, tron->diag, tao->subset_type, &tron->H_sub));
145a7e14dcfSSatish Balay if (tao->hessian == tao->hessian_pre) {
1469566063dSJacob Faibussowitsch PetscCall(MatDestroy(&tron->Hpre_sub));
147f4f49eeaSPierre Jolivet PetscCall(PetscObjectReference((PetscObject)tron->H_sub));
148a7e14dcfSSatish Balay tron->Hpre_sub = tron->H_sub;
149a7e14dcfSSatish Balay } else {
1509566063dSJacob Faibussowitsch PetscCall(TaoMatGetSubMat(tao->hessian_pre, tron->Free_Local, tron->diag, tao->subset_type, &tron->Hpre_sub));
151a7e14dcfSSatish Balay }
1529566063dSJacob Faibussowitsch PetscCall(KSPReset(tao->ksp));
1539566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(tao->ksp, tron->H_sub, tron->Hpre_sub));
154a7e14dcfSSatish Balay while (1) {
155a7e14dcfSSatish Balay /* Approximately solve the reduced linear system */
1569566063dSJacob Faibussowitsch PetscCall(KSPCGSetRadius(tao->ksp, delta));
157d8ec8e84SJason Sarich
1589566063dSJacob Faibussowitsch PetscCall(KSPSolve(tao->ksp, tron->R, tron->DXFree));
1599566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(tao->ksp, &its));
160a7e14dcfSSatish Balay tao->ksp_its += its;
161ae93cb3cSJason Sarich tao->ksp_tot_its += its;
1629566063dSJacob Faibussowitsch PetscCall(VecSet(tao->stepdirection, 0.0));
163a7e14dcfSSatish Balay
164a7e14dcfSSatish Balay /* Add dxfree matrix to compute step direction vector */
1659566063dSJacob Faibussowitsch PetscCall(VecISAXPY(tao->stepdirection, tron->Free_Local, 1.0, tron->DXFree));
166a7e14dcfSSatish Balay
1679566063dSJacob Faibussowitsch PetscCall(VecDot(tao->gradient, tao->stepdirection, &gdx));
1689566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->solution, tron->X_New));
1699566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->gradient, tron->G_New));
170a7e14dcfSSatish Balay
1719371c9d4SSatish Balay stepsize = 1.0;
1729371c9d4SSatish Balay f_new = f;
173a7e14dcfSSatish Balay
1749566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetInitialStepLength(tao->linesearch, 1.0));
1759566063dSJacob Faibussowitsch PetscCall(TaoLineSearchApply(tao->linesearch, tron->X_New, &f_new, tron->G_New, tao->stepdirection, &stepsize, &ls_reason));
1769566063dSJacob Faibussowitsch PetscCall(TaoAddLineSearchCounts(tao));
177a7e14dcfSSatish Balay
1789566063dSJacob Faibussowitsch PetscCall(MatMult(tao->hessian, tao->stepdirection, tron->Work));
1799566063dSJacob Faibussowitsch PetscCall(VecAYPX(tron->Work, 0.5, tao->gradient));
1809566063dSJacob Faibussowitsch PetscCall(VecDot(tao->stepdirection, tron->Work, &prered));
181a7e14dcfSSatish Balay actred = f_new - f;
182ce902467SBarry Smith if ((PetscAbsScalar(actred) <= 1e-6) && (PetscAbsScalar(prered) <= 1e-6)) {
183ce902467SBarry Smith rhok = 1.0;
184ce902467SBarry Smith } else if (actred < 0) {
185a7e14dcfSSatish Balay rhok = PetscAbs(-actred / prered);
186a7e14dcfSSatish Balay } else {
187a7e14dcfSSatish Balay rhok = 0.0;
188a7e14dcfSSatish Balay }
189a7e14dcfSSatish Balay
190a7e14dcfSSatish Balay /* Compare actual improvement to the quadratic model */
191a7e14dcfSSatish Balay if (rhok > tron->eta1) { /* Accept the point */
192a7e14dcfSSatish Balay /* d = x_new - x */
1939566063dSJacob Faibussowitsch PetscCall(VecCopy(tron->X_New, tao->stepdirection));
1949566063dSJacob Faibussowitsch PetscCall(VecAXPY(tao->stepdirection, -1.0, tao->solution));
195a7e14dcfSSatish Balay
1969566063dSJacob Faibussowitsch PetscCall(VecNorm(tao->stepdirection, NORM_2, &xdiff));
197a7e14dcfSSatish Balay xdiff *= stepsize;
198a7e14dcfSSatish Balay
199a7e14dcfSSatish Balay /* Adjust trust region size */
200a7e14dcfSSatish Balay if (rhok < tron->eta2) {
201a7e14dcfSSatish Balay delta = PetscMin(xdiff, delta) * tron->sigma1;
202a7e14dcfSSatish Balay } else if (rhok > tron->eta4) {
203a7e14dcfSSatish Balay delta = PetscMin(xdiff, delta) * tron->sigma3;
204a7e14dcfSSatish Balay } else if (rhok > tron->eta3) {
205a7e14dcfSSatish Balay delta = PetscMin(xdiff, delta) * tron->sigma2;
206a7e14dcfSSatish Balay }
2079566063dSJacob Faibussowitsch PetscCall(VecBoundGradientProjection(tron->G_New, tron->X_New, tao->XL, tao->XU, tao->gradient));
2089566063dSJacob Faibussowitsch PetscCall(ISDestroy(&tron->Free_Local));
2099566063dSJacob Faibussowitsch PetscCall(VecWhichInactive(tao->XL, tron->X_New, tao->gradient, tao->XU, PETSC_TRUE, &tron->Free_Local));
210a7e14dcfSSatish Balay f = f_new;
2119566063dSJacob Faibussowitsch PetscCall(VecNorm(tao->gradient, NORM_2, &tron->gnorm));
2129566063dSJacob Faibussowitsch PetscCall(VecCopy(tron->X_New, tao->solution));
2139566063dSJacob Faibussowitsch PetscCall(VecCopy(tron->G_New, tao->gradient));
214a7e14dcfSSatish Balay break;
2159371c9d4SSatish Balay } else if (delta <= 1e-30) {
216a7e14dcfSSatish Balay break;
2179371c9d4SSatish Balay } else {
218a7e14dcfSSatish Balay delta /= 4.0;
219a7e14dcfSSatish Balay }
220a7e14dcfSSatish Balay } /* end linear solve loop */
221a7e14dcfSSatish Balay
2229371c9d4SSatish Balay tron->f = f;
2239371c9d4SSatish Balay tron->actred = actred;
2249371c9d4SSatish Balay tao->trust = delta;
2258931d482SJason Sarich tao->niter++;
2269566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao, tron->f, tron->gnorm, 0.0, tao->ksp_its));
2279566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao, tao->niter, tron->f, tron->gnorm, 0.0, stepsize));
228dbbe0bcdSBarry Smith PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
229a7e14dcfSSatish Balay } /* END MAIN LOOP */
2303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
231a7e14dcfSSatish Balay }
232a7e14dcfSSatish Balay
TronGradientProjections(Tao tao,TAO_TRON * tron)233d71ae5a4SJacob Faibussowitsch static PetscErrorCode TronGradientProjections(Tao tao, TAO_TRON *tron)
234d71ae5a4SJacob Faibussowitsch {
235a7e14dcfSSatish Balay PetscInt i;
236e4cb33bbSBarry Smith TaoLineSearchConvergedReason ls_reason;
237a7e14dcfSSatish Balay PetscReal actred = -1.0, actred_max = 0.0;
238a7e14dcfSSatish Balay PetscReal f_new;
239a7e14dcfSSatish Balay /*
240a7e14dcfSSatish Balay The gradient and function value passed into and out of this
241a7e14dcfSSatish Balay routine should be current and correct.
242a7e14dcfSSatish Balay
243a7e14dcfSSatish Balay The free, active, and binding variables should be already identified
244a7e14dcfSSatish Balay */
245a7e14dcfSSatish Balay
2464d86920dSPierre Jolivet PetscFunctionBegin;
247ce902467SBarry Smith for (i = 0; i < tron->maxgpits; ++i) {
248a7e14dcfSSatish Balay if (-actred <= (tron->pg_ftol) * actred_max) break;
249a7e14dcfSSatish Balay
250ce902467SBarry Smith ++tron->gp_iterates;
251ce902467SBarry Smith ++tron->total_gp_its;
252a7e14dcfSSatish Balay f_new = tron->f;
253a7e14dcfSSatish Balay
2549566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->gradient, tao->stepdirection));
2559566063dSJacob Faibussowitsch PetscCall(VecScale(tao->stepdirection, -1.0));
2569566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetInitialStepLength(tao->linesearch, tron->pgstepsize));
257d0609cedSBarry Smith PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &f_new, tao->gradient, tao->stepdirection, &tron->pgstepsize, &ls_reason));
2589566063dSJacob Faibussowitsch PetscCall(TaoAddLineSearchCounts(tao));
259a7e14dcfSSatish Balay
2609566063dSJacob Faibussowitsch PetscCall(VecBoundGradientProjection(tao->gradient, tao->solution, tao->XL, tao->XU, tao->gradient));
2619566063dSJacob Faibussowitsch PetscCall(VecNorm(tao->gradient, NORM_2, &tron->gnorm));
262ce902467SBarry Smith
263a7e14dcfSSatish Balay /* Update the iterate */
264a7e14dcfSSatish Balay actred = f_new - tron->f;
265a7e14dcfSSatish Balay actred_max = PetscMax(actred_max, -(f_new - tron->f));
266a7e14dcfSSatish Balay tron->f = f_new;
267a7e14dcfSSatish Balay }
2683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
269a7e14dcfSSatish Balay }
270a7e14dcfSSatish Balay
TaoComputeDual_TRON(Tao tao,Vec DXL,Vec DXU)271d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoComputeDual_TRON(Tao tao, Vec DXL, Vec DXU)
272d71ae5a4SJacob Faibussowitsch {
273a7e14dcfSSatish Balay TAO_TRON *tron = (TAO_TRON *)tao->data;
274a7e14dcfSSatish Balay
275a7e14dcfSSatish Balay PetscFunctionBegin;
276441846f8SBarry Smith PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
277a7e14dcfSSatish Balay PetscValidHeaderSpecific(DXL, VEC_CLASSID, 2);
278a7e14dcfSSatish Balay PetscValidHeaderSpecific(DXU, VEC_CLASSID, 3);
2793c859ba3SBarry Smith PetscCheck(tron->Work && tao->gradient, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Dual variables don't exist yet or no longer exist.");
280a7e14dcfSSatish Balay
2819566063dSJacob Faibussowitsch PetscCall(VecBoundGradientProjection(tao->gradient, tao->solution, tao->XL, tao->XU, tron->Work));
2829566063dSJacob Faibussowitsch PetscCall(VecCopy(tron->Work, DXL));
2839566063dSJacob Faibussowitsch PetscCall(VecAXPY(DXL, -1.0, tao->gradient));
2849566063dSJacob Faibussowitsch PetscCall(VecSet(DXU, 0.0));
2859566063dSJacob Faibussowitsch PetscCall(VecPointwiseMax(DXL, DXL, DXU));
286a7e14dcfSSatish Balay
2879566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->gradient, DXU));
2889566063dSJacob Faibussowitsch PetscCall(VecAXPY(DXU, -1.0, tron->Work));
2899566063dSJacob Faibussowitsch PetscCall(VecSet(tron->Work, 0.0));
2909566063dSJacob Faibussowitsch PetscCall(VecPointwiseMin(DXU, tron->Work, DXU));
2913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
292a7e14dcfSSatish Balay }
293a7e14dcfSSatish Balay
2941522df2eSJason Sarich /*MC
2951522df2eSJason Sarich TAOTRON - The TRON algorithm is an active-set Newton trust region method
2961522df2eSJason Sarich for bound-constrained minimization.
2971522df2eSJason Sarich
2981522df2eSJason Sarich Options Database Keys:
2991522df2eSJason Sarich + -tao_tron_maxgpits - maximum number of gradient projections per TRON iterate
3001522df2eSJason Sarich - -tao_subset_type - "subvec","mask","matrix-free", strategies for handling active-sets
3011522df2eSJason Sarich
3021eb8069cSJason Sarich Level: beginner
3031522df2eSJason Sarich M*/
TaoCreate_TRON(Tao tao)304d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TaoCreate_TRON(Tao tao)
305d71ae5a4SJacob Faibussowitsch {
306a7e14dcfSSatish Balay TAO_TRON *tron;
3078caf6e8cSBarry Smith const char *morethuente_type = TAOLINESEARCHMT;
308a7e14dcfSSatish Balay
30947a47007SBarry Smith PetscFunctionBegin;
310a7e14dcfSSatish Balay tao->ops->setup = TaoSetup_TRON;
311a7e14dcfSSatish Balay tao->ops->solve = TaoSolve_TRON;
312a7e14dcfSSatish Balay tao->ops->view = TaoView_TRON;
313a7e14dcfSSatish Balay tao->ops->setfromoptions = TaoSetFromOptions_TRON;
314a7e14dcfSSatish Balay tao->ops->destroy = TaoDestroy_TRON;
315a7e14dcfSSatish Balay tao->ops->computedual = TaoComputeDual_TRON;
316a7e14dcfSSatish Balay
3174dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&tron));
3186f4723b1SBarry Smith tao->data = (void *)tron;
3196552cf8aSJason Sarich
3206552cf8aSJason Sarich /* Override default settings (unless already changed) */
321606f75f6SBarry Smith PetscCall(TaoParametersInitialize(tao));
322606f75f6SBarry Smith PetscObjectParameterSetDefault(tao, max_it, 50);
323606f75f6SBarry Smith PetscObjectParameterSetDefault(tao, trust0, 1.0);
324606f75f6SBarry Smith PetscObjectParameterSetDefault(tao, steptol, 0.0);
325a7e14dcfSSatish Balay
326a7e14dcfSSatish Balay /* Initialize pointers and variables */
327a7e14dcfSSatish Balay tron->n = 0;
328a7e14dcfSSatish Balay tron->maxgpits = 3;
329a7e14dcfSSatish Balay tron->pg_ftol = 0.001;
330a7e14dcfSSatish Balay
331a7e14dcfSSatish Balay tron->eta1 = 1.0e-4;
332a7e14dcfSSatish Balay tron->eta2 = 0.25;
333a7e14dcfSSatish Balay tron->eta3 = 0.50;
334a7e14dcfSSatish Balay tron->eta4 = 0.90;
335a7e14dcfSSatish Balay
336a7e14dcfSSatish Balay tron->sigma1 = 0.5;
337a7e14dcfSSatish Balay tron->sigma2 = 2.0;
338a7e14dcfSSatish Balay tron->sigma3 = 4.0;
339a7e14dcfSSatish Balay
340a7e14dcfSSatish Balay tron->gp_iterates = 0; /* Cumulative number */
341a7e14dcfSSatish Balay tron->total_gp_its = 0;
342a7e14dcfSSatish Balay tron->n_free = 0;
343a7e14dcfSSatish Balay
3446c23d075SBarry Smith tron->DXFree = NULL;
3456c23d075SBarry Smith tron->R = NULL;
3466c23d075SBarry Smith tron->X_New = NULL;
3476c23d075SBarry Smith tron->G_New = NULL;
3486c23d075SBarry Smith tron->Work = NULL;
3496c23d075SBarry Smith tron->Free_Local = NULL;
3506c23d075SBarry Smith tron->H_sub = NULL;
3516c23d075SBarry Smith tron->Hpre_sub = NULL;
352a7e14dcfSSatish Balay tao->subset_type = TAO_SUBSET_SUBVEC;
353a7e14dcfSSatish Balay
3549566063dSJacob Faibussowitsch PetscCall(TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch));
3559566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1));
3569566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetType(tao->linesearch, morethuente_type));
3579566063dSJacob Faibussowitsch PetscCall(TaoLineSearchUseTaoRoutines(tao->linesearch, tao));
3589566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetOptionsPrefix(tao->linesearch, tao->hdr.prefix));
359a7e14dcfSSatish Balay
3609566063dSJacob Faibussowitsch PetscCall(KSPCreate(((PetscObject)tao)->comm, &tao->ksp));
3619566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1));
3629566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(tao->ksp, tao->hdr.prefix));
3639566063dSJacob Faibussowitsch PetscCall(KSPSetType(tao->ksp, KSPSTCG));
3643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
365a7e14dcfSSatish Balay }
366