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