Actual source code: tron.c
petsc-dev 2014-02-02
1: #include <../src/tao/bound/impls/tron/tron.h>
2: #include <petsc-private/kspimpl.h>
3: #include <petsc-private/matimpl.h>
4: #include <../src/tao/matrix/submatfree.h>
7: /* TRON Routines */
8: static PetscErrorCode TronGradientProjections(Tao,TAO_TRON*);
9: /*------------------------------------------------------------*/
12: static PetscErrorCode TaoDestroy_TRON(Tao tao)
13: {
14: TAO_TRON *tron = (TAO_TRON *)tao->data;
18: VecDestroy(&tron->X_New);
19: VecDestroy(&tron->G_New);
20: VecDestroy(&tron->Work);
21: VecDestroy(&tron->DXFree);
22: VecDestroy(&tron->R);
23: VecDestroy(&tron->diag);
24: VecScatterDestroy(&tron->scatter);
25: ISDestroy(&tron->Free_Local);
26: MatDestroy(&tron->H_sub);
27: MatDestroy(&tron->Hpre_sub);
28: PetscFree(tao->data);
29: return(0);
30: }
32: /*------------------------------------------------------------*/
35: static PetscErrorCode TaoSetFromOptions_TRON(Tao tao)
36: {
37: TAO_TRON *tron = (TAO_TRON *)tao->data;
39: PetscBool flg;
42: PetscOptionsHead("Newton Trust Region Method for bound constrained optimization");
43: PetscOptionsInt("-tron_maxgpits","maximum number of gradient projections per TRON iterate","TaoSetMaxGPIts",tron->maxgpits,&tron->maxgpits,&flg);
44: PetscOptionsTail();
45: TaoLineSearchSetFromOptions(tao->linesearch);
46: KSPSetFromOptions(tao->ksp);
47: return(0);
48: }
50: /*------------------------------------------------------------*/
53: static PetscErrorCode TaoView_TRON(Tao tao, PetscViewer viewer)
54: {
55: TAO_TRON *tron = (TAO_TRON *)tao->data;
56: PetscBool isascii;
57: PetscErrorCode ierr;
60: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
61: if (isascii) {
62: PetscViewerASCIIPushTab(viewer);
63: PetscViewerASCIIPrintf(viewer,"Total PG its: %D,",tron->total_gp_its);
64: PetscViewerASCIIPrintf(viewer,"PG tolerance: %g \n",(double)tron->pg_ftol);
65: PetscViewerASCIIPopTab(viewer);
66: }
67: return(0);
68: }
71: /* ---------------------------------------------------------- */
74: static PetscErrorCode TaoSetup_TRON(Tao tao)
75: {
77: TAO_TRON *tron = (TAO_TRON *)tao->data;
81: /* Allocate some arrays */
82: VecDuplicate(tao->solution, &tron->diag);
83: VecDuplicate(tao->solution, &tron->X_New);
84: VecDuplicate(tao->solution, &tron->G_New);
85: VecDuplicate(tao->solution, &tron->Work);
86: VecDuplicate(tao->solution, &tao->gradient);
87: VecDuplicate(tao->solution, &tao->stepdirection);
88: if (!tao->XL) {
89: VecDuplicate(tao->solution, &tao->XL);
90: VecSet(tao->XL, PETSC_NINFINITY);
91: }
92: if (!tao->XU) {
93: VecDuplicate(tao->solution, &tao->XU);
94: VecSet(tao->XU, PETSC_INFINITY);
95: }
96: return(0);
97: }
103: static PetscErrorCode TaoSolve_TRON(Tao tao)
104: {
105: TAO_TRON *tron = (TAO_TRON *)tao->data;
106: PetscErrorCode ierr;
107: PetscInt iter=0,its;
108: TaoTerminationReason reason = TAO_CONTINUE_ITERATING;
109: TaoLineSearchTerminationReason ls_reason = TAOLINESEARCH_CONTINUE_ITERATING;
110: PetscReal prered,actred,delta,f,f_new,rhok,gdx,xdiff,stepsize;
113: tron->pgstepsize=1.0;
114: tao->trust = tao->trust0;
115: /* Project the current point onto the feasible set */
116: TaoComputeVariableBounds(tao);
117: VecMedian(tao->XL,tao->solution,tao->XU,tao->solution);
118: TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU);
120: TaoComputeObjectiveAndGradient(tao,tao->solution,&tron->f,tao->gradient);
121: ISDestroy(&tron->Free_Local);
123: VecWhichBetween(tao->XL,tao->solution,tao->XU,&tron->Free_Local);
125: /* Project the gradient and calculate the norm */
126: VecBoundGradientProjection(tao->gradient,tao->solution, tao->XL, tao->XU, tao->gradient);
127: VecNorm(tao->gradient,NORM_2,&tron->gnorm);
129: if (PetscIsInfOrNanReal(tron->f) || PetscIsInfOrNanReal(tron->gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf pr NaN");
130: if (tao->trust <= 0) {
131: tao->trust=PetscMax(tron->gnorm*tron->gnorm,1.0);
132: }
134: tron->stepsize=tao->trust;
135: TaoMonitor(tao, iter, tron->f, tron->gnorm, 0.0, tron->stepsize, &reason);
136: while (reason==TAO_CONTINUE_ITERATING){
138: TronGradientProjections(tao,tron);
139: f=tron->f; delta=tao->trust;
140: tron->n_free_last = tron->n_free;
141: TaoComputeHessian(tao,tao->solution,&tao->hessian, &tao->hessian_pre, &tron->matflag);
143: ISGetSize(tron->Free_Local, &tron->n_free);
145: /* If no free variables */
146: if (tron->n_free == 0) {
147: actred=0;
148: PetscInfo(tao,"No free variables in tron iteration.");
149: break;
151: }
152: /* use free_local to mask/submat gradient, hessian, stepdirection */
153: VecGetSubVec(tao->gradient,tron->Free_Local,tao->subset_type,0.0,&tron->R);
154: VecGetSubVec(tao->gradient,tron->Free_Local,tao->subset_type,0.0,&tron->DXFree);
155: VecSet(tron->DXFree,0.0);
156: VecScale(tron->R, -1.0);
157: MatGetSubMat(tao->hessian, tron->Free_Local, tron->diag, tao->subset_type, &tron->H_sub);
158: if (tao->hessian == tao->hessian_pre) {
159: MatDestroy(&tron->Hpre_sub);
160: PetscObjectReference((PetscObject)(tron->H_sub));
161: tron->Hpre_sub = tron->H_sub;
162: } else {
163: MatGetSubMat(tao->hessian_pre, tron->Free_Local, tron->diag, tao->subset_type,&tron->Hpre_sub);
164: }
165: KSPReset(tao->ksp);
166: KSPSetOperators(tao->ksp, tron->H_sub, tron->Hpre_sub, tron->matflag);
167: while (1) {
169: /* Approximately solve the reduced linear system */
170: KSPSTCGSetRadius(tao->ksp,delta);
171: KSPSolve(tao->ksp, tron->R, tron->DXFree);
172: KSPGetIterationNumber(tao->ksp,&its);
173: tao->ksp_its+=its;
174: VecSet(tao->stepdirection,0.0);
176: /* Add dxfree matrix to compute step direction vector */
177: VecISAXPY(tao->stepdirection,tron->Free_Local,1.0,tron->DXFree);
178: if (0) {
179: PetscReal rhs,stepnorm;
180: VecNorm(tron->R,NORM_2,&rhs);
181: VecNorm(tron->DXFree,NORM_2,&stepnorm);
182: PetscPrintf(PETSC_COMM_WORLD,"|rhs|=%g\t|s|=%g\n",(double)rhs,(double)stepnorm);
183: }
186: VecDot(tao->gradient, tao->stepdirection, &gdx);
187: PetscInfo1(tao,"Expected decrease in function value: %14.12e\n",(double)gdx);
189: VecCopy(tao->solution, tron->X_New);
190: VecCopy(tao->gradient, tron->G_New);
192: stepsize=1.0;f_new=f;
194: TaoLineSearchSetInitialStepLength(tao->linesearch,1.0);
195: TaoLineSearchApply(tao->linesearch, tron->X_New, &f_new, tron->G_New, tao->stepdirection,&stepsize,&ls_reason);
196: TaoAddLineSearchCounts(tao);
198: MatMult(tao->hessian, tao->stepdirection, tron->Work);
199: VecAYPX(tron->Work, 0.5, tao->gradient);
200: VecDot(tao->stepdirection, tron->Work, &prered);
201: actred = f_new - f;
202: if (actred<0) {
203: rhok=PetscAbs(-actred/prered);
204: } else {
205: rhok=0.0;
206: }
208: /* Compare actual improvement to the quadratic model */
209: if (rhok > tron->eta1) { /* Accept the point */
210: /* d = x_new - x */
211: VecCopy(tron->X_New, tao->stepdirection);
212: VecAXPY(tao->stepdirection, -1.0, tao->solution);
214: VecNorm(tao->stepdirection, NORM_2, &xdiff);
215: xdiff *= stepsize;
217: /* Adjust trust region size */
218: if (rhok < tron->eta2 ){
219: delta = PetscMin(xdiff,delta)*tron->sigma1;
220: } else if (rhok > tron->eta4 ){
221: delta= PetscMin(xdiff,delta)*tron->sigma3;
222: } else if (rhok > tron->eta3 ){
223: delta=PetscMin(xdiff,delta)*tron->sigma2;
224: }
225: VecBoundGradientProjection(tron->G_New,tron->X_New, tao->XL, tao->XU, tao->gradient);
226: if (tron->Free_Local) {
227: ISDestroy(&tron->Free_Local);
228: }
229: VecWhichBetween(tao->XL, tron->X_New, tao->XU, &tron->Free_Local);
230: f=f_new;
231: VecNorm(tao->gradient,NORM_2,&tron->gnorm);
232: VecCopy(tron->X_New, tao->solution);
233: VecCopy(tron->G_New, tao->gradient);
234: break;
235: }
236: else if (delta <= 1e-30) {
237: break;
238: }
239: else {
240: delta /= 4.0;
241: }
242: } /* end linear solve loop */
245: tron->f=f; tron->actred=actred; tao->trust=delta;
246: iter++;
247: TaoMonitor(tao, iter, tron->f, tron->gnorm, 0.0, delta, &reason);
248: } /* END MAIN LOOP */
250: return(0);
251: }
256: static PetscErrorCode TronGradientProjections(Tao tao,TAO_TRON *tron)
257: {
258: PetscErrorCode ierr;
259: PetscInt i;
260: TaoLineSearchTerminationReason ls_reason;
261: PetscReal actred=-1.0,actred_max=0.0;
262: PetscReal f_new;
263: /*
264: The gradient and function value passed into and out of this
265: routine should be current and correct.
267: The free, active, and binding variables should be already identified
268: */
270: if (tron->Free_Local) {
271: ISDestroy(&tron->Free_Local);
272: }
273: VecWhichBetween(tao->XL,tao->solution,tao->XU,&tron->Free_Local);
275: for (i=0;i<tron->maxgpits;i++){
277: if ( -actred <= (tron->pg_ftol)*actred_max) break;
279: tron->gp_iterates++; tron->total_gp_its++;
280: f_new=tron->f;
282: VecCopy(tao->gradient, tao->stepdirection);
283: VecScale(tao->stepdirection, -1.0);
284: TaoLineSearchSetInitialStepLength(tao->linesearch,tron->pgstepsize);
285: TaoLineSearchApply(tao->linesearch, tao->solution, &f_new, tao->gradient, tao->stepdirection,
286: &tron->pgstepsize, &ls_reason);
287: TaoAddLineSearchCounts(tao);
290: /* Update the iterate */
291: actred = f_new - tron->f;
292: actred_max = PetscMax(actred_max,-(f_new - tron->f));
293: tron->f = f_new;
294: if (tron->Free_Local) {
295: ISDestroy(&tron->Free_Local);
296: }
297: VecWhichBetween(tao->XL,tao->solution,tao->XU,&tron->Free_Local);
298: }
300: return(0);
301: }
305: static PetscErrorCode TaoComputeDual_TRON(Tao tao, Vec DXL, Vec DXU) {
307: TAO_TRON *tron = (TAO_TRON *)tao->data;
314: if (!tron->Work || !tao->gradient) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Dual variables don't exist yet or no longer exist.\n");
316: VecBoundGradientProjection(tao->gradient,tao->solution,tao->XL,tao->XU,tron->Work);
317: VecCopy(tron->Work,DXL);
318: VecAXPY(DXL,-1.0,tao->gradient);
319: VecSet(DXU,0.0);
320: VecPointwiseMax(DXL,DXL,DXU);
322: VecCopy(tao->gradient,DXU);
323: VecAXPY(DXU,-1.0,tron->Work);
324: VecSet(tron->Work,0.0);
325: VecPointwiseMin(DXU,tron->Work,DXU);
326: return(0);
327: }
329: /*------------------------------------------------------------*/
330: EXTERN_C_BEGIN
333: PetscErrorCode TaoCreate_TRON(Tao tao)
334: {
335: TAO_TRON *tron;
337: const char *morethuente_type = TAOLINESEARCH_MT;
340: tao->ops->setup = TaoSetup_TRON;
341: tao->ops->solve = TaoSolve_TRON;
342: tao->ops->view = TaoView_TRON;
343: tao->ops->setfromoptions = TaoSetFromOptions_TRON;
344: tao->ops->destroy = TaoDestroy_TRON;
345: tao->ops->computedual = TaoComputeDual_TRON;
347: PetscNewLog(tao,&tron);
349: tao->max_it = 50;
350: #if defined(PETSC_USE_REAL_SINGLE)
351: tao->fatol = 1e-5;
352: tao->frtol = 1e-5;
353: tao->steptol = 1e-6;
354: #else
355: tao->fatol = 1e-10;
356: tao->frtol = 1e-10;
357: tao->steptol = 1e-12;
358: #endif
359: tao->data = (void*)tron;
360: tao->trust0 = 1.0;
362: /* Initialize pointers and variables */
363: tron->n = 0;
364: tron->maxgpits = 3;
365: tron->pg_ftol = 0.001;
367: tron->eta1 = 1.0e-4;
368: tron->eta2 = 0.25;
369: tron->eta3 = 0.50;
370: tron->eta4 = 0.90;
372: tron->sigma1 = 0.5;
373: tron->sigma2 = 2.0;
374: tron->sigma3 = 4.0;
376: tron->gp_iterates = 0; /* Cumulative number */
377: tron->total_gp_its = 0;
378: tron->n_free = 0;
380: tron->DXFree=NULL;
381: tron->R=NULL;
382: tron->X_New=NULL;
383: tron->G_New=NULL;
384: tron->Work=NULL;
385: tron->Free_Local=NULL;
386: tron->H_sub=NULL;
387: tron->Hpre_sub=NULL;
388: tao->subset_type = TAO_SUBSET_SUBVEC;
390: TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch);
391: TaoLineSearchSetType(tao->linesearch,morethuente_type);
392: TaoLineSearchUseTaoRoutines(tao->linesearch,tao);
394: KSPCreate(((PetscObject)tao)->comm, &tao->ksp);
395: KSPSetType(tao->ksp,KSPSTCG);
396: return(0);
397: }
398: EXTERN_C_END