xref: /petsc/src/tao/bound/impls/tron/tron.c (revision f4f49eeac7efa77fffa46b7ff95a3ed169f659ed)
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 /*------------------------------------------------------------*/
7d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoDestroy_TRON(Tao tao)
8d71ae5a4SJacob Faibussowitsch {
9a7e14dcfSSatish Balay   TAO_TRON *tron = (TAO_TRON *)tao->data;
10a7e14dcfSSatish Balay 
11a7e14dcfSSatish Balay   PetscFunctionBegin;
129566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&tron->X_New));
139566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&tron->G_New));
149566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&tron->Work));
159566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&tron->DXFree));
169566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&tron->R));
179566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&tron->diag));
189566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&tron->scatter));
199566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&tron->Free_Local));
209566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&tron->H_sub));
219566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&tron->Hpre_sub));
22a958fbfcSStefano Zampini   PetscCall(KSPDestroy(&tao->ksp));
239566063dSJacob Faibussowitsch   PetscCall(PetscFree(tao->data));
243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
25a7e14dcfSSatish Balay }
26a7e14dcfSSatish Balay 
27a7e14dcfSSatish Balay /*------------------------------------------------------------*/
28d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSetFromOptions_TRON(Tao tao, PetscOptionItems *PetscOptionsObject)
29d71ae5a4SJacob Faibussowitsch {
30a7e14dcfSSatish Balay   TAO_TRON *tron = (TAO_TRON *)tao->data;
31a7e14dcfSSatish Balay   PetscBool flg;
32a7e14dcfSSatish Balay 
33a7e14dcfSSatish Balay   PetscFunctionBegin;
34d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Newton Trust Region Method for bound constrained optimization");
359566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-tao_tron_maxgpits", "maximum number of gradient projections per TRON iterate", "TaoSetMaxGPIts", tron->maxgpits, &tron->maxgpits, &flg));
36d0609cedSBarry Smith   PetscOptionsHeadEnd();
379566063dSJacob Faibussowitsch   PetscCall(KSPSetFromOptions(tao->ksp));
383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
39a7e14dcfSSatish Balay }
40a7e14dcfSSatish Balay 
41a7e14dcfSSatish Balay /*------------------------------------------------------------*/
42d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoView_TRON(Tao tao, PetscViewer viewer)
43d71ae5a4SJacob Faibussowitsch {
44a7e14dcfSSatish Balay   TAO_TRON *tron = (TAO_TRON *)tao->data;
45a7e14dcfSSatish Balay   PetscBool isascii;
46a7e14dcfSSatish Balay 
47a7e14dcfSSatish Balay   PetscFunctionBegin;
489566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
49a7e14dcfSSatish Balay   if (isascii) {
5063a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Total PG its: %" PetscInt_FMT ",", tron->total_gp_its));
519566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "PG tolerance: %g \n", (double)tron->pg_ftol));
52a7e14dcfSSatish Balay   }
533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
54a7e14dcfSSatish Balay }
55a7e14dcfSSatish Balay 
56a7e14dcfSSatish Balay /* ---------------------------------------------------------- */
57d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSetup_TRON(Tao tao)
58d71ae5a4SJacob Faibussowitsch {
59a7e14dcfSSatish Balay   TAO_TRON *tron = (TAO_TRON *)tao->data;
60a7e14dcfSSatish Balay 
61a7e14dcfSSatish Balay   PetscFunctionBegin;
62a7e14dcfSSatish Balay   /* Allocate some arrays */
639566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tao->solution, &tron->diag));
649566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tao->solution, &tron->X_New));
659566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tao->solution, &tron->G_New));
669566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tao->solution, &tron->Work));
679566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tao->solution, &tao->gradient));
689566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(tao->solution, &tao->stepdirection));
693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
70a7e14dcfSSatish Balay }
71a7e14dcfSSatish Balay 
72d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSolve_TRON(Tao tao)
73d71ae5a4SJacob Faibussowitsch {
74a7e14dcfSSatish Balay   TAO_TRON                    *tron = (TAO_TRON *)tao->data;
758931d482SJason Sarich   PetscInt                     its;
76e4cb33bbSBarry Smith   TaoLineSearchConvergedReason ls_reason = TAOLINESEARCH_CONTINUE_ITERATING;
77a7e14dcfSSatish Balay   PetscReal                    prered, actred, delta, f, f_new, rhok, gdx, xdiff, stepsize;
7847a47007SBarry Smith 
79a7e14dcfSSatish Balay   PetscFunctionBegin;
80a7e14dcfSSatish Balay   tron->pgstepsize = 1.0;
81a7e14dcfSSatish Balay   tao->trust       = tao->trust0;
82a7e14dcfSSatish Balay   /*   Project the current point onto the feasible set */
839566063dSJacob Faibussowitsch   PetscCall(TaoComputeVariableBounds(tao));
849566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetVariableBounds(tao->linesearch, tao->XL, tao->XU));
85a7e14dcfSSatish Balay 
86ce902467SBarry Smith   /* Project the initial point onto the feasible region */
879566063dSJacob Faibussowitsch   PetscCall(VecMedian(tao->XL, tao->solution, tao->XU, tao->solution));
8853506e15SBarry Smith 
89ce902467SBarry Smith   /* Compute the objective function and gradient */
909566063dSJacob Faibussowitsch   PetscCall(TaoComputeObjectiveAndGradient(tao, tao->solution, &tron->f, tao->gradient));
919566063dSJacob Faibussowitsch   PetscCall(VecNorm(tao->gradient, NORM_2, &tron->gnorm));
923c859ba3SBarry Smith   PetscCheck(!PetscIsInfOrNanReal(tron->f) && !PetscIsInfOrNanReal(tron->gnorm), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute function generated Inf or NaN");
93a7e14dcfSSatish Balay 
94a7e14dcfSSatish Balay   /* Project the gradient and calculate the norm */
959566063dSJacob Faibussowitsch   PetscCall(VecBoundGradientProjection(tao->gradient, tao->solution, tao->XL, tao->XU, tao->gradient));
969566063dSJacob Faibussowitsch   PetscCall(VecNorm(tao->gradient, NORM_2, &tron->gnorm));
97a7e14dcfSSatish Balay 
98ce902467SBarry Smith   /* Initialize trust region radius */
99ce902467SBarry Smith   tao->trust = tao->trust0;
100ad540459SPierre Jolivet   if (tao->trust <= 0) tao->trust = PetscMax(tron->gnorm * tron->gnorm, 1.0);
101a7e14dcfSSatish Balay 
102ce902467SBarry Smith   /* Initialize step sizes for the line searches */
103ce902467SBarry Smith   tron->pgstepsize = 1.0;
104a7e14dcfSSatish Balay   tron->stepsize   = tao->trust;
105ce902467SBarry Smith 
1063ecd9318SAlp Dener   tao->reason = TAO_CONTINUE_ITERATING;
1079566063dSJacob Faibussowitsch   PetscCall(TaoLogConvergenceHistory(tao, tron->f, tron->gnorm, 0.0, tao->ksp_its));
1089566063dSJacob Faibussowitsch   PetscCall(TaoMonitor(tao, tao->niter, tron->f, tron->gnorm, 0.0, tron->stepsize));
109dbbe0bcdSBarry Smith   PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
1103ecd9318SAlp Dener   while (tao->reason == TAO_CONTINUE_ITERATING) {
111e1e80dc8SAlp Dener     /* Call general purpose update function */
112dbbe0bcdSBarry Smith     PetscTryTypeMethod(tao, update, tao->niter, tao->user_update);
113ce902467SBarry Smith 
114ce902467SBarry Smith     /* Perform projected gradient iterations */
1159566063dSJacob Faibussowitsch     PetscCall(TronGradientProjections(tao, tron));
116ce902467SBarry Smith 
1179566063dSJacob Faibussowitsch     PetscCall(VecBoundGradientProjection(tao->gradient, tao->solution, tao->XL, tao->XU, tao->gradient));
1189566063dSJacob Faibussowitsch     PetscCall(VecNorm(tao->gradient, NORM_2, &tron->gnorm));
119ce902467SBarry Smith 
120ce902467SBarry Smith     tao->ksp_its      = 0;
1219371c9d4SSatish Balay     f                 = tron->f;
1229371c9d4SSatish Balay     delta             = tao->trust;
123a7e14dcfSSatish Balay     tron->n_free_last = tron->n_free;
1249566063dSJacob Faibussowitsch     PetscCall(TaoComputeHessian(tao, tao->solution, tao->hessian, tao->hessian_pre));
125a7e14dcfSSatish Balay 
126baa3a814SAlp Dener     /* Generate index set (IS) of which bound constraints are active */
1279566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&tron->Free_Local));
1289566063dSJacob Faibussowitsch     PetscCall(VecWhichInactive(tao->XL, tao->solution, tao->gradient, tao->XU, PETSC_TRUE, &tron->Free_Local));
1299566063dSJacob Faibussowitsch     PetscCall(ISGetSize(tron->Free_Local, &tron->n_free));
130a7e14dcfSSatish Balay 
131a7e14dcfSSatish Balay     /* If no free variables */
132a7e14dcfSSatish Balay     if (tron->n_free == 0) {
1339566063dSJacob Faibussowitsch       PetscCall(VecNorm(tao->gradient, NORM_2, &tron->gnorm));
1349566063dSJacob Faibussowitsch       PetscCall(TaoLogConvergenceHistory(tao, tron->f, tron->gnorm, 0.0, tao->ksp_its));
1359566063dSJacob Faibussowitsch       PetscCall(TaoMonitor(tao, tao->niter, tron->f, tron->gnorm, 0.0, delta));
136dbbe0bcdSBarry Smith       PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
137ad540459SPierre Jolivet       if (!tao->reason) tao->reason = TAO_CONVERGED_STEPTOL;
138a7e14dcfSSatish Balay       break;
139a7e14dcfSSatish Balay     }
140a7e14dcfSSatish Balay     /* use free_local to mask/submat gradient, hessian, stepdirection */
1419566063dSJacob Faibussowitsch     PetscCall(TaoVecGetSubVec(tao->gradient, tron->Free_Local, tao->subset_type, 0.0, &tron->R));
1429566063dSJacob Faibussowitsch     PetscCall(TaoVecGetSubVec(tao->gradient, tron->Free_Local, tao->subset_type, 0.0, &tron->DXFree));
1439566063dSJacob Faibussowitsch     PetscCall(VecSet(tron->DXFree, 0.0));
1449566063dSJacob Faibussowitsch     PetscCall(VecScale(tron->R, -1.0));
1459566063dSJacob Faibussowitsch     PetscCall(TaoMatGetSubMat(tao->hessian, tron->Free_Local, tron->diag, tao->subset_type, &tron->H_sub));
146a7e14dcfSSatish Balay     if (tao->hessian == tao->hessian_pre) {
1479566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&tron->Hpre_sub));
148*f4f49eeaSPierre Jolivet       PetscCall(PetscObjectReference((PetscObject)tron->H_sub));
149a7e14dcfSSatish Balay       tron->Hpre_sub = tron->H_sub;
150a7e14dcfSSatish Balay     } else {
1519566063dSJacob Faibussowitsch       PetscCall(TaoMatGetSubMat(tao->hessian_pre, tron->Free_Local, tron->diag, tao->subset_type, &tron->Hpre_sub));
152a7e14dcfSSatish Balay     }
1539566063dSJacob Faibussowitsch     PetscCall(KSPReset(tao->ksp));
1549566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(tao->ksp, tron->H_sub, tron->Hpre_sub));
155a7e14dcfSSatish Balay     while (1) {
156a7e14dcfSSatish Balay       /* Approximately solve the reduced linear system */
1579566063dSJacob Faibussowitsch       PetscCall(KSPCGSetRadius(tao->ksp, delta));
158d8ec8e84SJason Sarich 
1599566063dSJacob Faibussowitsch       PetscCall(KSPSolve(tao->ksp, tron->R, tron->DXFree));
1609566063dSJacob Faibussowitsch       PetscCall(KSPGetIterationNumber(tao->ksp, &its));
161a7e14dcfSSatish Balay       tao->ksp_its += its;
162ae93cb3cSJason Sarich       tao->ksp_tot_its += its;
1639566063dSJacob Faibussowitsch       PetscCall(VecSet(tao->stepdirection, 0.0));
164a7e14dcfSSatish Balay 
165a7e14dcfSSatish Balay       /* Add dxfree matrix to compute step direction vector */
1669566063dSJacob Faibussowitsch       PetscCall(VecISAXPY(tao->stepdirection, tron->Free_Local, 1.0, tron->DXFree));
167a7e14dcfSSatish Balay 
1689566063dSJacob Faibussowitsch       PetscCall(VecDot(tao->gradient, tao->stepdirection, &gdx));
1699566063dSJacob Faibussowitsch       PetscCall(VecCopy(tao->solution, tron->X_New));
1709566063dSJacob Faibussowitsch       PetscCall(VecCopy(tao->gradient, tron->G_New));
171a7e14dcfSSatish Balay 
1729371c9d4SSatish Balay       stepsize = 1.0;
1739371c9d4SSatish Balay       f_new    = f;
174a7e14dcfSSatish Balay 
1759566063dSJacob Faibussowitsch       PetscCall(TaoLineSearchSetInitialStepLength(tao->linesearch, 1.0));
1769566063dSJacob Faibussowitsch       PetscCall(TaoLineSearchApply(tao->linesearch, tron->X_New, &f_new, tron->G_New, tao->stepdirection, &stepsize, &ls_reason));
1779566063dSJacob Faibussowitsch       PetscCall(TaoAddLineSearchCounts(tao));
178a7e14dcfSSatish Balay 
1799566063dSJacob Faibussowitsch       PetscCall(MatMult(tao->hessian, tao->stepdirection, tron->Work));
1809566063dSJacob Faibussowitsch       PetscCall(VecAYPX(tron->Work, 0.5, tao->gradient));
1819566063dSJacob Faibussowitsch       PetscCall(VecDot(tao->stepdirection, tron->Work, &prered));
182a7e14dcfSSatish Balay       actred = f_new - f;
183ce902467SBarry Smith       if ((PetscAbsScalar(actred) <= 1e-6) && (PetscAbsScalar(prered) <= 1e-6)) {
184ce902467SBarry Smith         rhok = 1.0;
185ce902467SBarry Smith       } else if (actred < 0) {
186a7e14dcfSSatish Balay         rhok = PetscAbs(-actred / prered);
187a7e14dcfSSatish Balay       } else {
188a7e14dcfSSatish Balay         rhok = 0.0;
189a7e14dcfSSatish Balay       }
190a7e14dcfSSatish Balay 
191a7e14dcfSSatish Balay       /* Compare actual improvement to the quadratic model */
192a7e14dcfSSatish Balay       if (rhok > tron->eta1) { /* Accept the point */
193a7e14dcfSSatish Balay         /* d = x_new - x */
1949566063dSJacob Faibussowitsch         PetscCall(VecCopy(tron->X_New, tao->stepdirection));
1959566063dSJacob Faibussowitsch         PetscCall(VecAXPY(tao->stepdirection, -1.0, tao->solution));
196a7e14dcfSSatish Balay 
1979566063dSJacob Faibussowitsch         PetscCall(VecNorm(tao->stepdirection, NORM_2, &xdiff));
198a7e14dcfSSatish Balay         xdiff *= stepsize;
199a7e14dcfSSatish Balay 
200a7e14dcfSSatish Balay         /* Adjust trust region size */
201a7e14dcfSSatish Balay         if (rhok < tron->eta2) {
202a7e14dcfSSatish Balay           delta = PetscMin(xdiff, delta) * tron->sigma1;
203a7e14dcfSSatish Balay         } else if (rhok > tron->eta4) {
204a7e14dcfSSatish Balay           delta = PetscMin(xdiff, delta) * tron->sigma3;
205a7e14dcfSSatish Balay         } else if (rhok > tron->eta3) {
206a7e14dcfSSatish Balay           delta = PetscMin(xdiff, delta) * tron->sigma2;
207a7e14dcfSSatish Balay         }
2089566063dSJacob Faibussowitsch         PetscCall(VecBoundGradientProjection(tron->G_New, tron->X_New, tao->XL, tao->XU, tao->gradient));
2099566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&tron->Free_Local));
2109566063dSJacob Faibussowitsch         PetscCall(VecWhichInactive(tao->XL, tron->X_New, tao->gradient, tao->XU, PETSC_TRUE, &tron->Free_Local));
211a7e14dcfSSatish Balay         f = f_new;
2129566063dSJacob Faibussowitsch         PetscCall(VecNorm(tao->gradient, NORM_2, &tron->gnorm));
2139566063dSJacob Faibussowitsch         PetscCall(VecCopy(tron->X_New, tao->solution));
2149566063dSJacob Faibussowitsch         PetscCall(VecCopy(tron->G_New, tao->gradient));
215a7e14dcfSSatish Balay         break;
2169371c9d4SSatish Balay       } else if (delta <= 1e-30) {
217a7e14dcfSSatish Balay         break;
2189371c9d4SSatish Balay       } else {
219a7e14dcfSSatish Balay         delta /= 4.0;
220a7e14dcfSSatish Balay       }
221a7e14dcfSSatish Balay     } /* end linear solve loop */
222a7e14dcfSSatish Balay 
2239371c9d4SSatish Balay     tron->f      = f;
2249371c9d4SSatish Balay     tron->actred = actred;
2259371c9d4SSatish Balay     tao->trust   = delta;
2268931d482SJason Sarich     tao->niter++;
2279566063dSJacob Faibussowitsch     PetscCall(TaoLogConvergenceHistory(tao, tron->f, tron->gnorm, 0.0, tao->ksp_its));
2289566063dSJacob Faibussowitsch     PetscCall(TaoMonitor(tao, tao->niter, tron->f, tron->gnorm, 0.0, stepsize));
229dbbe0bcdSBarry Smith     PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
230a7e14dcfSSatish Balay   } /* END MAIN LOOP  */
2313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
232a7e14dcfSSatish Balay }
233a7e14dcfSSatish Balay 
234d71ae5a4SJacob Faibussowitsch static PetscErrorCode TronGradientProjections(Tao tao, TAO_TRON *tron)
235d71ae5a4SJacob Faibussowitsch {
236a7e14dcfSSatish Balay   PetscInt                     i;
237e4cb33bbSBarry Smith   TaoLineSearchConvergedReason ls_reason;
238a7e14dcfSSatish Balay   PetscReal                    actred = -1.0, actred_max = 0.0;
239a7e14dcfSSatish Balay   PetscReal                    f_new;
240a7e14dcfSSatish Balay   /*
241a7e14dcfSSatish Balay      The gradient and function value passed into and out of this
242a7e14dcfSSatish Balay      routine should be current and correct.
243a7e14dcfSSatish Balay 
244a7e14dcfSSatish Balay      The free, active, and binding variables should be already identified
245a7e14dcfSSatish Balay   */
246a7e14dcfSSatish Balay 
2474d86920dSPierre Jolivet   PetscFunctionBegin;
248ce902467SBarry Smith   for (i = 0; i < tron->maxgpits; ++i) {
249a7e14dcfSSatish Balay     if (-actred <= (tron->pg_ftol) * actred_max) break;
250a7e14dcfSSatish Balay 
251ce902467SBarry Smith     ++tron->gp_iterates;
252ce902467SBarry Smith     ++tron->total_gp_its;
253a7e14dcfSSatish Balay     f_new = tron->f;
254a7e14dcfSSatish Balay 
2559566063dSJacob Faibussowitsch     PetscCall(VecCopy(tao->gradient, tao->stepdirection));
2569566063dSJacob Faibussowitsch     PetscCall(VecScale(tao->stepdirection, -1.0));
2579566063dSJacob Faibussowitsch     PetscCall(TaoLineSearchSetInitialStepLength(tao->linesearch, tron->pgstepsize));
258d0609cedSBarry Smith     PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &f_new, tao->gradient, tao->stepdirection, &tron->pgstepsize, &ls_reason));
2599566063dSJacob Faibussowitsch     PetscCall(TaoAddLineSearchCounts(tao));
260a7e14dcfSSatish Balay 
2619566063dSJacob Faibussowitsch     PetscCall(VecBoundGradientProjection(tao->gradient, tao->solution, tao->XL, tao->XU, tao->gradient));
2629566063dSJacob Faibussowitsch     PetscCall(VecNorm(tao->gradient, NORM_2, &tron->gnorm));
263ce902467SBarry Smith 
264a7e14dcfSSatish Balay     /* Update the iterate */
265a7e14dcfSSatish Balay     actred     = f_new - tron->f;
266a7e14dcfSSatish Balay     actred_max = PetscMax(actred_max, -(f_new - tron->f));
267a7e14dcfSSatish Balay     tron->f    = f_new;
268a7e14dcfSSatish Balay   }
2693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
270a7e14dcfSSatish Balay }
271a7e14dcfSSatish Balay 
272d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoComputeDual_TRON(Tao tao, Vec DXL, Vec DXU)
273d71ae5a4SJacob Faibussowitsch {
274a7e14dcfSSatish Balay   TAO_TRON *tron = (TAO_TRON *)tao->data;
275a7e14dcfSSatish Balay 
276a7e14dcfSSatish Balay   PetscFunctionBegin;
277441846f8SBarry Smith   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
278a7e14dcfSSatish Balay   PetscValidHeaderSpecific(DXL, VEC_CLASSID, 2);
279a7e14dcfSSatish Balay   PetscValidHeaderSpecific(DXU, VEC_CLASSID, 3);
2803c859ba3SBarry Smith   PetscCheck(tron->Work && tao->gradient, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Dual variables don't exist yet or no longer exist.");
281a7e14dcfSSatish Balay 
2829566063dSJacob Faibussowitsch   PetscCall(VecBoundGradientProjection(tao->gradient, tao->solution, tao->XL, tao->XU, tron->Work));
2839566063dSJacob Faibussowitsch   PetscCall(VecCopy(tron->Work, DXL));
2849566063dSJacob Faibussowitsch   PetscCall(VecAXPY(DXL, -1.0, tao->gradient));
2859566063dSJacob Faibussowitsch   PetscCall(VecSet(DXU, 0.0));
2869566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMax(DXL, DXL, DXU));
287a7e14dcfSSatish Balay 
2889566063dSJacob Faibussowitsch   PetscCall(VecCopy(tao->gradient, DXU));
2899566063dSJacob Faibussowitsch   PetscCall(VecAXPY(DXU, -1.0, tron->Work));
2909566063dSJacob Faibussowitsch   PetscCall(VecSet(tron->Work, 0.0));
2919566063dSJacob Faibussowitsch   PetscCall(VecPointwiseMin(DXU, tron->Work, DXU));
2923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
293a7e14dcfSSatish Balay }
294a7e14dcfSSatish Balay 
295a7e14dcfSSatish Balay /*------------------------------------------------------------*/
2961522df2eSJason Sarich /*MC
2971522df2eSJason Sarich   TAOTRON - The TRON algorithm is an active-set Newton trust region method
2981522df2eSJason Sarich   for bound-constrained minimization.
2991522df2eSJason Sarich 
3001522df2eSJason Sarich   Options Database Keys:
3011522df2eSJason Sarich + -tao_tron_maxgpits - maximum number of gradient projections per TRON iterate
3021522df2eSJason Sarich - -tao_subset_type - "subvec","mask","matrix-free", strategies for handling active-sets
3031522df2eSJason Sarich 
3041eb8069cSJason Sarich   Level: beginner
3051522df2eSJason Sarich M*/
306d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TaoCreate_TRON(Tao tao)
307d71ae5a4SJacob Faibussowitsch {
308a7e14dcfSSatish Balay   TAO_TRON   *tron;
3098caf6e8cSBarry Smith   const char *morethuente_type = TAOLINESEARCHMT;
310a7e14dcfSSatish Balay 
31147a47007SBarry Smith   PetscFunctionBegin;
312a7e14dcfSSatish Balay   tao->ops->setup          = TaoSetup_TRON;
313a7e14dcfSSatish Balay   tao->ops->solve          = TaoSolve_TRON;
314a7e14dcfSSatish Balay   tao->ops->view           = TaoView_TRON;
315a7e14dcfSSatish Balay   tao->ops->setfromoptions = TaoSetFromOptions_TRON;
316a7e14dcfSSatish Balay   tao->ops->destroy        = TaoDestroy_TRON;
317a7e14dcfSSatish Balay   tao->ops->computedual    = TaoComputeDual_TRON;
318a7e14dcfSSatish Balay 
3194dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&tron));
3206f4723b1SBarry Smith   tao->data = (void *)tron;
3216552cf8aSJason Sarich 
3226552cf8aSJason Sarich   /* Override default settings (unless already changed) */
3236552cf8aSJason Sarich   if (!tao->max_it_changed) tao->max_it = 50;
3246552cf8aSJason Sarich   if (!tao->trust0_changed) tao->trust0 = 1.0;
3254f1535bcSTodd Munson   if (!tao->steptol_changed) tao->steptol = 0.0;
326a7e14dcfSSatish Balay 
327a7e14dcfSSatish Balay   /* Initialize pointers and variables */
328a7e14dcfSSatish Balay   tron->n        = 0;
329a7e14dcfSSatish Balay   tron->maxgpits = 3;
330a7e14dcfSSatish Balay   tron->pg_ftol  = 0.001;
331a7e14dcfSSatish Balay 
332a7e14dcfSSatish Balay   tron->eta1 = 1.0e-4;
333a7e14dcfSSatish Balay   tron->eta2 = 0.25;
334a7e14dcfSSatish Balay   tron->eta3 = 0.50;
335a7e14dcfSSatish Balay   tron->eta4 = 0.90;
336a7e14dcfSSatish Balay 
337a7e14dcfSSatish Balay   tron->sigma1 = 0.5;
338a7e14dcfSSatish Balay   tron->sigma2 = 2.0;
339a7e14dcfSSatish Balay   tron->sigma3 = 4.0;
340a7e14dcfSSatish Balay 
341a7e14dcfSSatish Balay   tron->gp_iterates  = 0; /* Cumulative number */
342a7e14dcfSSatish Balay   tron->total_gp_its = 0;
343a7e14dcfSSatish Balay   tron->n_free       = 0;
344a7e14dcfSSatish Balay 
3456c23d075SBarry Smith   tron->DXFree     = NULL;
3466c23d075SBarry Smith   tron->R          = NULL;
3476c23d075SBarry Smith   tron->X_New      = NULL;
3486c23d075SBarry Smith   tron->G_New      = NULL;
3496c23d075SBarry Smith   tron->Work       = NULL;
3506c23d075SBarry Smith   tron->Free_Local = NULL;
3516c23d075SBarry Smith   tron->H_sub      = NULL;
3526c23d075SBarry Smith   tron->Hpre_sub   = NULL;
353a7e14dcfSSatish Balay   tao->subset_type = TAO_SUBSET_SUBVEC;
354a7e14dcfSSatish Balay 
3559566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch));
3569566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1));
3579566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetType(tao->linesearch, morethuente_type));
3589566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchUseTaoRoutines(tao->linesearch, tao));
3599566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetOptionsPrefix(tao->linesearch, tao->hdr.prefix));
360a7e14dcfSSatish Balay 
3619566063dSJacob Faibussowitsch   PetscCall(KSPCreate(((PetscObject)tao)->comm, &tao->ksp));
3629566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1));
3639566063dSJacob Faibussowitsch   PetscCall(KSPSetOptionsPrefix(tao->ksp, tao->hdr.prefix));
3649566063dSJacob Faibussowitsch   PetscCall(KSPSetType(tao->ksp, KSPSTCG));
3653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
366a7e14dcfSSatish Balay }
367