xref: /petsc/src/tao/bound/impls/tron/tron.c (revision 8fcddce65efd55a8fe3f87d4c08c15577ce4cbef)
1aaa7dc30SBarry Smith #include <../src/tao/bound/impls/tron/tron.h>
2aaa7dc30SBarry Smith #include <../src/tao/matrix/submatfree.h>
3a7e14dcfSSatish Balay 
4a7e14dcfSSatish Balay 
5a7e14dcfSSatish Balay /* TRON Routines */
6441846f8SBarry Smith static PetscErrorCode TronGradientProjections(Tao,TAO_TRON*);
7a7e14dcfSSatish Balay /*------------------------------------------------------------*/
8441846f8SBarry Smith static PetscErrorCode TaoDestroy_TRON(Tao tao)
9a7e14dcfSSatish Balay {
10a7e14dcfSSatish Balay   TAO_TRON       *tron = (TAO_TRON *)tao->data;
11a7e14dcfSSatish Balay   PetscErrorCode ierr;
12a7e14dcfSSatish Balay 
13a7e14dcfSSatish Balay   PetscFunctionBegin;
14a7e14dcfSSatish Balay   ierr = VecDestroy(&tron->X_New);CHKERRQ(ierr);
15a7e14dcfSSatish Balay   ierr = VecDestroy(&tron->G_New);CHKERRQ(ierr);
16a7e14dcfSSatish Balay   ierr = VecDestroy(&tron->Work);CHKERRQ(ierr);
17a7e14dcfSSatish Balay   ierr = VecDestroy(&tron->DXFree);CHKERRQ(ierr);
18a7e14dcfSSatish Balay   ierr = VecDestroy(&tron->R);CHKERRQ(ierr);
19a7e14dcfSSatish Balay   ierr = VecDestroy(&tron->diag);CHKERRQ(ierr);
20a7e14dcfSSatish Balay   ierr = VecScatterDestroy(&tron->scatter);CHKERRQ(ierr);
21a7e14dcfSSatish Balay   ierr = ISDestroy(&tron->Free_Local);CHKERRQ(ierr);
22a7e14dcfSSatish Balay   ierr = MatDestroy(&tron->H_sub);CHKERRQ(ierr);
23a7e14dcfSSatish Balay   ierr = MatDestroy(&tron->Hpre_sub);CHKERRQ(ierr);
24302440fdSBarry Smith   ierr = PetscFree(tao->data);CHKERRQ(ierr);
25a7e14dcfSSatish Balay   PetscFunctionReturn(0);
26a7e14dcfSSatish Balay }
27a7e14dcfSSatish Balay 
28a7e14dcfSSatish Balay /*------------------------------------------------------------*/
294416b707SBarry Smith static PetscErrorCode TaoSetFromOptions_TRON(PetscOptionItems *PetscOptionsObject,Tao tao)
30a7e14dcfSSatish Balay {
31a7e14dcfSSatish Balay   TAO_TRON       *tron = (TAO_TRON *)tao->data;
32a7e14dcfSSatish Balay   PetscErrorCode ierr;
33a7e14dcfSSatish Balay   PetscBool      flg;
34a7e14dcfSSatish Balay 
35a7e14dcfSSatish Balay   PetscFunctionBegin;
361a1499c8SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Newton Trust Region Method for bound constrained optimization");CHKERRQ(ierr);
371522df2eSJason Sarich   ierr = PetscOptionsInt("-tao_tron_maxgpits","maximum number of gradient projections per TRON iterate","TaoSetMaxGPIts",tron->maxgpits,&tron->maxgpits,&flg);CHKERRQ(ierr);
38a7e14dcfSSatish Balay   ierr = PetscOptionsTail();CHKERRQ(ierr);
39a7e14dcfSSatish Balay   ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr);
40a7e14dcfSSatish Balay   ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
41a7e14dcfSSatish Balay   PetscFunctionReturn(0);
42a7e14dcfSSatish Balay }
43a7e14dcfSSatish Balay 
44a7e14dcfSSatish Balay /*------------------------------------------------------------*/
45441846f8SBarry Smith static PetscErrorCode TaoView_TRON(Tao tao, PetscViewer viewer)
46a7e14dcfSSatish Balay {
47a7e14dcfSSatish Balay   TAO_TRON         *tron = (TAO_TRON *)tao->data;
48a7e14dcfSSatish Balay   PetscBool        isascii;
49a7e14dcfSSatish Balay   PetscErrorCode   ierr;
50a7e14dcfSSatish Balay 
51a7e14dcfSSatish Balay   PetscFunctionBegin;
52a7e14dcfSSatish Balay   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
53a7e14dcfSSatish Balay   if (isascii) {
54a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"Total PG its: %D,",tron->total_gp_its);CHKERRQ(ierr);
5547a47007SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"PG tolerance: %g \n",(double)tron->pg_ftol);CHKERRQ(ierr);
56a7e14dcfSSatish Balay   }
57a7e14dcfSSatish Balay   PetscFunctionReturn(0);
58a7e14dcfSSatish Balay }
59a7e14dcfSSatish Balay 
60a7e14dcfSSatish Balay /* ---------------------------------------------------------- */
61441846f8SBarry Smith static PetscErrorCode TaoSetup_TRON(Tao tao)
62a7e14dcfSSatish Balay {
63a7e14dcfSSatish Balay   PetscErrorCode ierr;
64a7e14dcfSSatish Balay   TAO_TRON       *tron = (TAO_TRON *)tao->data;
65a7e14dcfSSatish Balay 
66a7e14dcfSSatish Balay   PetscFunctionBegin;
67a7e14dcfSSatish Balay 
68a7e14dcfSSatish Balay   /* Allocate some arrays */
69a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution, &tron->diag);CHKERRQ(ierr);
70a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution, &tron->X_New);CHKERRQ(ierr);
71a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution, &tron->G_New);CHKERRQ(ierr);
72a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution, &tron->Work);CHKERRQ(ierr);
73a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution, &tao->gradient);CHKERRQ(ierr);
74a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution, &tao->stepdirection);CHKERRQ(ierr);
75a7e14dcfSSatish Balay   if (!tao->XL) {
76a7e14dcfSSatish Balay     ierr = VecDuplicate(tao->solution, &tao->XL);CHKERRQ(ierr);
77e270355aSBarry Smith     ierr = VecSet(tao->XL, PETSC_NINFINITY);CHKERRQ(ierr);
78a7e14dcfSSatish Balay   }
79a7e14dcfSSatish Balay   if (!tao->XU) {
80a7e14dcfSSatish Balay     ierr = VecDuplicate(tao->solution, &tao->XU);CHKERRQ(ierr);
81e270355aSBarry Smith     ierr = VecSet(tao->XU, PETSC_INFINITY);CHKERRQ(ierr);
82a7e14dcfSSatish Balay   }
83a7e14dcfSSatish Balay   PetscFunctionReturn(0);
84a7e14dcfSSatish Balay }
85a7e14dcfSSatish Balay 
86441846f8SBarry Smith static PetscErrorCode TaoSolve_TRON(Tao tao)
8747a47007SBarry Smith {
88a7e14dcfSSatish Balay   TAO_TRON                     *tron = (TAO_TRON *)tao->data;
89a7e14dcfSSatish Balay   PetscErrorCode               ierr;
908931d482SJason Sarich   PetscInt                     its;
91e4cb33bbSBarry Smith   TaoLineSearchConvergedReason ls_reason = TAOLINESEARCH_CONTINUE_ITERATING;
92a7e14dcfSSatish Balay   PetscReal                    prered,actred,delta,f,f_new,rhok,gdx,xdiff,stepsize;
9347a47007SBarry Smith 
94a7e14dcfSSatish Balay   PetscFunctionBegin;
95a7e14dcfSSatish Balay   tron->pgstepsize = 1.0;
96a7e14dcfSSatish Balay   tao->trust = tao->trust0;
97a7e14dcfSSatish Balay   /*   Project the current point onto the feasible set */
98a7e14dcfSSatish Balay   ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr);
99a7e14dcfSSatish Balay   ierr = TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU);CHKERRQ(ierr);
100a7e14dcfSSatish Balay 
101ce902467SBarry Smith   /* Project the initial point onto the feasible region */
102ce902467SBarry Smith   ierr = VecMedian(tao->XL,tao->solution,tao->XU,tao->solution);CHKERRQ(ierr);
10353506e15SBarry Smith 
104ce902467SBarry Smith   /* Compute the objective function and gradient */
105ce902467SBarry Smith   ierr = TaoComputeObjectiveAndGradient(tao,tao->solution,&tron->f,tao->gradient);CHKERRQ(ierr);
106ce902467SBarry Smith   ierr = VecNorm(tao->gradient,NORM_2,&tron->gnorm);CHKERRQ(ierr);
107ce902467SBarry Smith   if (PetscIsInfOrNanReal(tron->f) || PetscIsInfOrNanReal(tron->gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
108a7e14dcfSSatish Balay 
109a7e14dcfSSatish Balay   /* Project the gradient and calculate the norm */
110a7e14dcfSSatish Balay   ierr = VecBoundGradientProjection(tao->gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr);
111a7e14dcfSSatish Balay   ierr = VecNorm(tao->gradient,NORM_2,&tron->gnorm);CHKERRQ(ierr);
112a7e14dcfSSatish Balay 
113ce902467SBarry Smith   /* Initialize trust region radius */
114ce902467SBarry Smith   tao->trust=tao->trust0;
115a7e14dcfSSatish Balay   if (tao->trust <= 0) {
116a7e14dcfSSatish Balay     tao->trust=PetscMax(tron->gnorm*tron->gnorm,1.0);
117a7e14dcfSSatish Balay   }
118a7e14dcfSSatish Balay 
119ce902467SBarry Smith   /* Initialize step sizes for the line searches */
120ce902467SBarry Smith   tron->pgstepsize=1.0;
121a7e14dcfSSatish Balay   tron->stepsize=tao->trust;
122ce902467SBarry Smith 
1233ecd9318SAlp Dener   tao->reason = TAO_CONTINUE_ITERATING;
1243ecd9318SAlp Dener   ierr = TaoLogConvergenceHistory(tao,tron->f,tron->gnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
1253ecd9318SAlp Dener   ierr = TaoMonitor(tao,tao->niter,tron->f,tron->gnorm,0.0,tron->stepsize);CHKERRQ(ierr);
1263ecd9318SAlp Dener   ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
1273ecd9318SAlp Dener   while (tao->reason==TAO_CONTINUE_ITERATING){
128e1e80dc8SAlp Dener     /* Call general purpose update function */
129e1e80dc8SAlp Dener     if (tao->ops->update) {
130*8fcddce6SStefano Zampini       ierr = (*tao->ops->update)(tao, tao->niter, tao->user_update);CHKERRQ(ierr);
131e1e80dc8SAlp Dener     }
132ce902467SBarry Smith 
133ce902467SBarry Smith     /* Perform projected gradient iterations */
134a7e14dcfSSatish Balay     ierr = TronGradientProjections(tao,tron);CHKERRQ(ierr);
135ce902467SBarry Smith 
136ce902467SBarry Smith     ierr = VecBoundGradientProjection(tao->gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr);
137ce902467SBarry Smith     ierr = VecNorm(tao->gradient,NORM_2,&tron->gnorm);CHKERRQ(ierr);
138ce902467SBarry Smith 
139ce902467SBarry Smith     tao->ksp_its=0;
140a7e14dcfSSatish Balay     f=tron->f; delta=tao->trust;
141a7e14dcfSSatish Balay     tron->n_free_last = tron->n_free;
142ffad9901SBarry Smith     ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr);
143a7e14dcfSSatish Balay 
144baa3a814SAlp Dener     /* Generate index set (IS) of which bound constraints are active */
145ce902467SBarry Smith     ierr = ISDestroy(&tron->Free_Local);CHKERRQ(ierr);
146ce902467SBarry Smith     ierr = VecWhichInactive(tao->XL,tao->solution,tao->gradient,tao->XU,PETSC_TRUE,&tron->Free_Local);CHKERRQ(ierr);
147a7e14dcfSSatish Balay     ierr = ISGetSize(tron->Free_Local, &tron->n_free);CHKERRQ(ierr);
148a7e14dcfSSatish Balay 
149a7e14dcfSSatish Balay     /* If no free variables */
150a7e14dcfSSatish Balay     if (tron->n_free == 0) {
15155723ba5SJason Sarich       ierr = VecNorm(tao->gradient,NORM_2,&tron->gnorm);CHKERRQ(ierr);
1523ecd9318SAlp Dener       ierr = TaoLogConvergenceHistory(tao,tron->f,tron->gnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
1533ecd9318SAlp Dener       ierr = TaoMonitor(tao,tao->niter,tron->f,tron->gnorm,0.0,delta);CHKERRQ(ierr);
1543ecd9318SAlp Dener       ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
1553ecd9318SAlp Dener       if (!tao->reason) {
1563ecd9318SAlp Dener         tao->reason = TAO_CONVERGED_STEPTOL;
15755723ba5SJason Sarich       }
158a7e14dcfSSatish Balay       break;
159a7e14dcfSSatish Balay     }
160a7e14dcfSSatish Balay     /* use free_local to mask/submat gradient, hessian, stepdirection */
161b98f30f2SJason Sarich     ierr = TaoVecGetSubVec(tao->gradient,tron->Free_Local,tao->subset_type,0.0,&tron->R);CHKERRQ(ierr);
162b98f30f2SJason Sarich     ierr = TaoVecGetSubVec(tao->gradient,tron->Free_Local,tao->subset_type,0.0,&tron->DXFree);CHKERRQ(ierr);
163a7e14dcfSSatish Balay     ierr = VecSet(tron->DXFree,0.0);CHKERRQ(ierr);
164a7e14dcfSSatish Balay     ierr = VecScale(tron->R, -1.0);CHKERRQ(ierr);
165b98f30f2SJason Sarich     ierr = TaoMatGetSubMat(tao->hessian, tron->Free_Local, tron->diag, tao->subset_type, &tron->H_sub);CHKERRQ(ierr);
166a7e14dcfSSatish Balay     if (tao->hessian == tao->hessian_pre) {
167a7e14dcfSSatish Balay       ierr = MatDestroy(&tron->Hpre_sub);CHKERRQ(ierr);
168a7e14dcfSSatish Balay       ierr = PetscObjectReference((PetscObject)(tron->H_sub));CHKERRQ(ierr);
169a7e14dcfSSatish Balay       tron->Hpre_sub = tron->H_sub;
170a7e14dcfSSatish Balay     } else {
171b98f30f2SJason Sarich       ierr = TaoMatGetSubMat(tao->hessian_pre, tron->Free_Local, tron->diag, tao->subset_type,&tron->Hpre_sub);CHKERRQ(ierr);
172a7e14dcfSSatish Balay     }
173a7e14dcfSSatish Balay     ierr = KSPReset(tao->ksp);CHKERRQ(ierr);
17423ee1639SBarry Smith     ierr = KSPSetOperators(tao->ksp, tron->H_sub, tron->Hpre_sub);CHKERRQ(ierr);
175a7e14dcfSSatish Balay     while (1) {
176a7e14dcfSSatish Balay 
177a7e14dcfSSatish Balay       /* Approximately solve the reduced linear system */
1788f1e51d3STodd Munson       ierr = KSPCGSetRadius(tao->ksp,delta);CHKERRQ(ierr);
179d8ec8e84SJason Sarich 
180a7e14dcfSSatish Balay       ierr = KSPSolve(tao->ksp, tron->R, tron->DXFree);CHKERRQ(ierr);
181a7e14dcfSSatish Balay       ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr);
182a7e14dcfSSatish Balay       tao->ksp_its+=its;
183ae93cb3cSJason Sarich       tao->ksp_tot_its+=its;
184a7e14dcfSSatish Balay       ierr = VecSet(tao->stepdirection,0.0);CHKERRQ(ierr);
185a7e14dcfSSatish Balay 
186a7e14dcfSSatish Balay       /* Add dxfree matrix to compute step direction vector */
1874473680cSBarry Smith       ierr = VecISAXPY(tao->stepdirection,tron->Free_Local,1.0,tron->DXFree);CHKERRQ(ierr);
188a7e14dcfSSatish Balay 
189a7e14dcfSSatish Balay       ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr);
190a7e14dcfSSatish Balay       ierr = VecCopy(tao->solution, tron->X_New);CHKERRQ(ierr);
191a7e14dcfSSatish Balay       ierr = VecCopy(tao->gradient, tron->G_New);CHKERRQ(ierr);
192a7e14dcfSSatish Balay 
193a7e14dcfSSatish Balay       stepsize=1.0;f_new=f;
194a7e14dcfSSatish Balay 
195a7e14dcfSSatish Balay       ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,1.0);CHKERRQ(ierr);
19647a47007SBarry Smith       ierr = TaoLineSearchApply(tao->linesearch, tron->X_New, &f_new, tron->G_New, tao->stepdirection,&stepsize,&ls_reason);CHKERRQ(ierr);CHKERRQ(ierr);
197a7e14dcfSSatish Balay       ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
198a7e14dcfSSatish Balay 
199a7e14dcfSSatish Balay       ierr = MatMult(tao->hessian, tao->stepdirection, tron->Work);CHKERRQ(ierr);
200a7e14dcfSSatish Balay       ierr = VecAYPX(tron->Work, 0.5, tao->gradient);CHKERRQ(ierr);
201a7e14dcfSSatish Balay       ierr = VecDot(tao->stepdirection, tron->Work, &prered);CHKERRQ(ierr);
202a7e14dcfSSatish Balay       actred = f_new - f;
203ce902467SBarry Smith       if ((PetscAbsScalar(actred) <= 1e-6) && (PetscAbsScalar(prered) <= 1e-6)) {
204ce902467SBarry Smith         rhok = 1.0;
205ce902467SBarry Smith       } else if (actred<0) {
206a7e14dcfSSatish Balay         rhok=PetscAbs(-actred/prered);
207a7e14dcfSSatish Balay       } else {
208a7e14dcfSSatish Balay         rhok=0.0;
209a7e14dcfSSatish Balay       }
210a7e14dcfSSatish Balay 
211a7e14dcfSSatish Balay       /* Compare actual improvement to the quadratic model */
212a7e14dcfSSatish Balay       if (rhok > tron->eta1) { /* Accept the point */
213a7e14dcfSSatish Balay         /* d = x_new - x */
214a7e14dcfSSatish Balay         ierr = VecCopy(tron->X_New, tao->stepdirection);CHKERRQ(ierr);
215a7e14dcfSSatish Balay         ierr = VecAXPY(tao->stepdirection, -1.0, tao->solution);CHKERRQ(ierr);
216a7e14dcfSSatish Balay 
217a7e14dcfSSatish Balay         ierr = VecNorm(tao->stepdirection, NORM_2, &xdiff);CHKERRQ(ierr);
218a7e14dcfSSatish Balay         xdiff *= stepsize;
219a7e14dcfSSatish Balay 
220a7e14dcfSSatish Balay         /* Adjust trust region size */
221a7e14dcfSSatish Balay         if (rhok < tron->eta2 ){
222a7e14dcfSSatish Balay           delta = PetscMin(xdiff,delta)*tron->sigma1;
223a7e14dcfSSatish Balay         } else if (rhok > tron->eta4 ){
224a7e14dcfSSatish Balay           delta= PetscMin(xdiff,delta)*tron->sigma3;
225a7e14dcfSSatish Balay         } else if (rhok > tron->eta3 ){
226a7e14dcfSSatish Balay           delta=PetscMin(xdiff,delta)*tron->sigma2;
227a7e14dcfSSatish Balay         }
228a7e14dcfSSatish Balay         ierr = VecBoundGradientProjection(tron->G_New,tron->X_New, tao->XL, tao->XU, tao->gradient);CHKERRQ(ierr);
229a7e14dcfSSatish Balay         ierr = ISDestroy(&tron->Free_Local);CHKERRQ(ierr);
230ce902467SBarry Smith         ierr = VecWhichInactive(tao->XL,tron->X_New,tao->gradient,tao->XU,PETSC_TRUE,&tron->Free_Local);CHKERRQ(ierr);
231a7e14dcfSSatish Balay         f=f_new;
232a7e14dcfSSatish Balay         ierr = VecNorm(tao->gradient,NORM_2,&tron->gnorm);CHKERRQ(ierr);
233a7e14dcfSSatish Balay         ierr = VecCopy(tron->X_New, tao->solution);CHKERRQ(ierr);
234a7e14dcfSSatish Balay         ierr = VecCopy(tron->G_New, tao->gradient);CHKERRQ(ierr);
235a7e14dcfSSatish Balay         break;
236a7e14dcfSSatish Balay       }
237a7e14dcfSSatish Balay       else if (delta <= 1e-30) {
238a7e14dcfSSatish Balay         break;
239a7e14dcfSSatish Balay       }
240a7e14dcfSSatish Balay       else {
241a7e14dcfSSatish Balay         delta /= 4.0;
242a7e14dcfSSatish Balay       }
243a7e14dcfSSatish Balay     } /* end linear solve loop */
244a7e14dcfSSatish Balay 
245a7e14dcfSSatish Balay     tron->f=f; tron->actred=actred; tao->trust=delta;
2468931d482SJason Sarich     tao->niter++;
2473ecd9318SAlp Dener     ierr = TaoLogConvergenceHistory(tao,tron->f,tron->gnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
248770b7498SAlp Dener     ierr = TaoMonitor(tao,tao->niter,tron->f,tron->gnorm,0.0,stepsize);CHKERRQ(ierr);
2493ecd9318SAlp Dener     ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
250a7e14dcfSSatish Balay   }  /* END MAIN LOOP  */
251a7e14dcfSSatish Balay   PetscFunctionReturn(0);
252a7e14dcfSSatish Balay }
253a7e14dcfSSatish Balay 
254441846f8SBarry Smith static PetscErrorCode TronGradientProjections(Tao tao,TAO_TRON *tron)
255a7e14dcfSSatish Balay {
256a7e14dcfSSatish Balay   PetscErrorCode               ierr;
257a7e14dcfSSatish Balay   PetscInt                     i;
258e4cb33bbSBarry Smith   TaoLineSearchConvergedReason ls_reason;
259a7e14dcfSSatish Balay   PetscReal                    actred=-1.0,actred_max=0.0;
260a7e14dcfSSatish Balay   PetscReal                    f_new;
261a7e14dcfSSatish Balay   /*
262a7e14dcfSSatish Balay      The gradient and function value passed into and out of this
263a7e14dcfSSatish Balay      routine should be current and correct.
264a7e14dcfSSatish Balay 
265a7e14dcfSSatish Balay      The free, active, and binding variables should be already identified
266a7e14dcfSSatish Balay   */
267a7e14dcfSSatish Balay   PetscFunctionBegin;
268a7e14dcfSSatish Balay 
269ce902467SBarry Smith   for (i=0;i<tron->maxgpits;++i){
270a7e14dcfSSatish Balay 
271a7e14dcfSSatish Balay     if (-actred <= (tron->pg_ftol)*actred_max) break;
272a7e14dcfSSatish Balay 
273ce902467SBarry Smith     ++tron->gp_iterates;
274ce902467SBarry Smith     ++tron->total_gp_its;
275a7e14dcfSSatish Balay     f_new=tron->f;
276a7e14dcfSSatish Balay 
277a7e14dcfSSatish Balay     ierr = VecCopy(tao->gradient,tao->stepdirection);CHKERRQ(ierr);
278a7e14dcfSSatish Balay     ierr = VecScale(tao->stepdirection,-1.0);CHKERRQ(ierr);
279a7e14dcfSSatish Balay     ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,tron->pgstepsize);CHKERRQ(ierr);
280a7e14dcfSSatish Balay     ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f_new, tao->gradient, tao->stepdirection,
281a7e14dcfSSatish Balay                               &tron->pgstepsize, &ls_reason);CHKERRQ(ierr);
282a7e14dcfSSatish Balay     ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
283a7e14dcfSSatish Balay 
284ce902467SBarry Smith     ierr = VecBoundGradientProjection(tao->gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr);
285ce902467SBarry Smith     ierr = VecNorm(tao->gradient,NORM_2,&tron->gnorm);CHKERRQ(ierr);
286ce902467SBarry Smith 
287a7e14dcfSSatish Balay     /* Update the iterate */
288a7e14dcfSSatish Balay     actred = f_new - tron->f;
289a7e14dcfSSatish Balay     actred_max = PetscMax(actred_max,-(f_new - tron->f));
290a7e14dcfSSatish Balay     tron->f = f_new;
291a7e14dcfSSatish Balay   }
292a7e14dcfSSatish Balay   PetscFunctionReturn(0);
293a7e14dcfSSatish Balay }
294a7e14dcfSSatish Balay 
295441846f8SBarry Smith static PetscErrorCode TaoComputeDual_TRON(Tao tao, Vec DXL, Vec DXU) {
296a7e14dcfSSatish Balay 
297a7e14dcfSSatish Balay   TAO_TRON       *tron = (TAO_TRON *)tao->data;
298a7e14dcfSSatish Balay   PetscErrorCode ierr;
299a7e14dcfSSatish Balay 
300a7e14dcfSSatish Balay   PetscFunctionBegin;
301441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
302a7e14dcfSSatish Balay   PetscValidHeaderSpecific(DXL,VEC_CLASSID,2);
303a7e14dcfSSatish Balay   PetscValidHeaderSpecific(DXU,VEC_CLASSID,3);
30447a47007SBarry Smith   if (!tron->Work || !tao->gradient) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Dual variables don't exist yet or no longer exist.\n");
305a7e14dcfSSatish Balay 
306a7e14dcfSSatish Balay   ierr = VecBoundGradientProjection(tao->gradient,tao->solution,tao->XL,tao->XU,tron->Work);CHKERRQ(ierr);
307a7e14dcfSSatish Balay   ierr = VecCopy(tron->Work,DXL);CHKERRQ(ierr);
308a7e14dcfSSatish Balay   ierr = VecAXPY(DXL,-1.0,tao->gradient);CHKERRQ(ierr);
309a7e14dcfSSatish Balay   ierr = VecSet(DXU,0.0);CHKERRQ(ierr);
310a7e14dcfSSatish Balay   ierr = VecPointwiseMax(DXL,DXL,DXU);CHKERRQ(ierr);
311a7e14dcfSSatish Balay 
312a7e14dcfSSatish Balay   ierr = VecCopy(tao->gradient,DXU);CHKERRQ(ierr);
313a7e14dcfSSatish Balay   ierr = VecAXPY(DXU,-1.0,tron->Work);CHKERRQ(ierr);
314a7e14dcfSSatish Balay   ierr = VecSet(tron->Work,0.0);CHKERRQ(ierr);
315a7e14dcfSSatish Balay   ierr = VecPointwiseMin(DXU,tron->Work,DXU);CHKERRQ(ierr);
316a7e14dcfSSatish Balay   PetscFunctionReturn(0);
317a7e14dcfSSatish Balay }
318a7e14dcfSSatish Balay 
319a7e14dcfSSatish Balay /*------------------------------------------------------------*/
3201522df2eSJason Sarich /*MC
3211522df2eSJason Sarich   TAOTRON - The TRON algorithm is an active-set Newton trust region method
3221522df2eSJason Sarich   for bound-constrained minimization.
3231522df2eSJason Sarich 
3241522df2eSJason Sarich   Options Database Keys:
3251522df2eSJason Sarich + -tao_tron_maxgpits - maximum number of gradient projections per TRON iterate
3261522df2eSJason Sarich - -tao_subset_type - "subvec","mask","matrix-free", strategies for handling active-sets
3271522df2eSJason Sarich 
3281eb8069cSJason Sarich   Level: beginner
3291522df2eSJason Sarich M*/
330728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate_TRON(Tao tao)
331a7e14dcfSSatish Balay {
332a7e14dcfSSatish Balay   TAO_TRON       *tron;
333a7e14dcfSSatish Balay   PetscErrorCode ierr;
3348caf6e8cSBarry Smith   const char     *morethuente_type = TAOLINESEARCHMT;
335a7e14dcfSSatish Balay 
33647a47007SBarry Smith   PetscFunctionBegin;
337a7e14dcfSSatish Balay   tao->ops->setup          = TaoSetup_TRON;
338a7e14dcfSSatish Balay   tao->ops->solve          = TaoSolve_TRON;
339a7e14dcfSSatish Balay   tao->ops->view           = TaoView_TRON;
340a7e14dcfSSatish Balay   tao->ops->setfromoptions = TaoSetFromOptions_TRON;
341a7e14dcfSSatish Balay   tao->ops->destroy        = TaoDestroy_TRON;
342a7e14dcfSSatish Balay   tao->ops->computedual    = TaoComputeDual_TRON;
343a7e14dcfSSatish Balay 
3443c9e27cfSGeoffrey Irving   ierr = PetscNewLog(tao,&tron);CHKERRQ(ierr);
3456f4723b1SBarry Smith   tao->data = (void*)tron;
3466552cf8aSJason Sarich 
3476552cf8aSJason Sarich   /* Override default settings (unless already changed) */
3486552cf8aSJason Sarich   if (!tao->max_it_changed) tao->max_it = 50;
3496552cf8aSJason Sarich   if (!tao->trust0_changed) tao->trust0 = 1.0;
3504f1535bcSTodd Munson   if (!tao->steptol_changed) tao->steptol = 0.0;
351a7e14dcfSSatish Balay 
352a7e14dcfSSatish Balay   /* Initialize pointers and variables */
353a7e14dcfSSatish Balay   tron->n            = 0;
354a7e14dcfSSatish Balay   tron->maxgpits     = 3;
355a7e14dcfSSatish Balay   tron->pg_ftol      = 0.001;
356a7e14dcfSSatish Balay 
357a7e14dcfSSatish Balay   tron->eta1         = 1.0e-4;
358a7e14dcfSSatish Balay   tron->eta2         = 0.25;
359a7e14dcfSSatish Balay   tron->eta3         = 0.50;
360a7e14dcfSSatish Balay   tron->eta4         = 0.90;
361a7e14dcfSSatish Balay 
362a7e14dcfSSatish Balay   tron->sigma1       = 0.5;
363a7e14dcfSSatish Balay   tron->sigma2       = 2.0;
364a7e14dcfSSatish Balay   tron->sigma3       = 4.0;
365a7e14dcfSSatish Balay 
366a7e14dcfSSatish Balay   tron->gp_iterates  = 0; /* Cumulative number */
367a7e14dcfSSatish Balay   tron->total_gp_its = 0;
368a7e14dcfSSatish Balay   tron->n_free       = 0;
369a7e14dcfSSatish Balay 
3706c23d075SBarry Smith   tron->DXFree=NULL;
3716c23d075SBarry Smith   tron->R=NULL;
3726c23d075SBarry Smith   tron->X_New=NULL;
3736c23d075SBarry Smith   tron->G_New=NULL;
3746c23d075SBarry Smith   tron->Work=NULL;
3756c23d075SBarry Smith   tron->Free_Local=NULL;
3766c23d075SBarry Smith   tron->H_sub=NULL;
3776c23d075SBarry Smith   tron->Hpre_sub=NULL;
378a7e14dcfSSatish Balay   tao->subset_type = TAO_SUBSET_SUBVEC;
379a7e14dcfSSatish Balay 
380a7e14dcfSSatish Balay   ierr = TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch);CHKERRQ(ierr);
38163b15415SAlp Dener   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr);
382a7e14dcfSSatish Balay   ierr = TaoLineSearchSetType(tao->linesearch,morethuente_type);CHKERRQ(ierr);
383441846f8SBarry Smith   ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr);
3845d527766SPatrick Farrell   ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr);
385a7e14dcfSSatish Balay 
386a7e14dcfSSatish Balay   ierr = KSPCreate(((PetscObject)tao)->comm, &tao->ksp);CHKERRQ(ierr);
38763b15415SAlp Dener   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1);CHKERRQ(ierr);
3885d527766SPatrick Farrell   ierr = KSPSetOptionsPrefix(tao->ksp, tao->hdr.prefix);CHKERRQ(ierr);
3898f1e51d3STodd Munson   ierr = KSPSetType(tao->ksp,KSPCGSTCG);CHKERRQ(ierr);
390a7e14dcfSSatish Balay   PetscFunctionReturn(0);
391a7e14dcfSSatish Balay }
392