Actual source code: bntr.c

  1: #include <../src/tao/bound/impls/bnk/bnk.h>
  2: #include <petscksp.h>

  4: /*
  5:  Implements Newton's Method with a trust region approach for solving
  6:  bound constrained minimization problems.

  8:  ------------------------------------------------------------

 10:  x_0 = VecMedian(x_0)
 11:  f_0, g_0= TaoComputeObjectiveAndGradient(x_0)
 12:  pg_0 = project(g_0)
 13:  check convergence at pg_0
 14:  needH = TaoBNKInitialize(default:BNK_INIT_INTERPOLATION)
 15:  niter = 0
 16:  step_accepted = false

 18:  while niter <= max_it

 20:     if needH
 21:       If max_cg_steps > 0
 22:         x_k, g_k, pg_k = TaoSolve(BNCG)
 23:       end

 25:       H_k = TaoComputeHessian(x_k)
 26:       if pc_type == BNK_PC_BFGS
 27:         add correction to BFGS approx
 28:         if scale_type == BNK_SCALE_AHESS
 29:           D = VecMedian(1e-6, abs(diag(H_k)), 1e6)
 30:           scale BFGS with VecReciprocal(D)
 31:         end
 32:       end
 33:       needH = False
 34:     end

 36:     if pc_type = BNK_PC_BFGS
 37:       B_k = BFGS
 38:     else
 39:       B_k = VecMedian(1e-6, abs(diag(H_k)), 1e6)
 40:       B_k = VecReciprocal(B_k)
 41:     end
 42:     w = x_k - VecMedian(x_k - 0.001*B_k*g_k)
 43:     eps = min(eps, norm2(w))
 44:     determine the active and inactive index sets such that
 45:       L = {i : (x_k)_i <= l_i + eps && (g_k)_i > 0}
 46:       U = {i : (x_k)_i >= u_i - eps && (g_k)_i < 0}
 47:       F = {i : l_i = (x_k)_i = u_i}
 48:       A = {L + U + F}
 49:       IA = {i : i not in A}

 51:     generate the reduced system Hr_k dr_k = -gr_k for variables in IA
 52:     if pc_type == BNK_PC_BFGS && scale_type == BNK_SCALE_PHESS
 53:       D = VecMedian(1e-6, abs(diag(Hr_k)), 1e6)
 54:       scale BFGS with VecReciprocal(D)
 55:     end

 57:     while !stepAccepted
 58:       solve Hr_k dr_k = -gr_k
 59:       set d_k to (l - x) for variables in L, (u - x) for variables in U, and 0 for variables in F

 61:       x_{k+1} = VecMedian(x_k + d_k)
 62:       s = x_{k+1} - x_k
 63:       prered = dot(s, 0.5*gr_k - Hr_k*s)
 64:       f_{k+1} = TaoComputeObjective(x_{k+1})
 65:       actred = f_k - f_{k+1}

 67:       oldTrust = trust
 68:       step_accepted, trust = TaoBNKUpdateTrustRadius(default: BNK_UPDATE_REDUCTION)
 69:       if step_accepted
 70:         g_{k+1} = TaoComputeGradient(x_{k+1})
 71:         pg_{k+1} = project(g_{k+1})
 72:         count the accepted Newton step
 73:         needH = True
 74:       else
 75:         f_{k+1} = f_k
 76:         x_{k+1} = x_k
 77:         g_{k+1} = g_k
 78:         pg_{k+1} = pg_k
 79:         if trust == oldTrust
 80:           terminate because we cannot shrink the radius any further
 81:         end
 82:       end

 84:       check convergence at pg_{k+1}
 85:     end
 86:     niter += 1

 88:  end
 89: */

 91: PetscErrorCode TaoSolve_BNTR(Tao tao)
 92: {
 93:   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
 94:   KSPConvergedReason           ksp_reason;

 96:   PetscReal                    oldTrust, prered, actred, steplen, resnorm;
 97:   PetscBool                    cgTerminate, needH = PETSC_TRUE, stepAccepted, shift = PETSC_FALSE;
 98:   PetscInt                     stepType, nDiff;

100:   /* Initialize the preconditioner, KSP solver and trust radius/line search */
101:   tao->reason = TAO_CONTINUE_ITERATING;
102:   TaoBNKInitialize(tao, bnk->init_type, &needH);
103:   if (tao->reason != TAO_CONTINUE_ITERATING) return 0;

105:   /* Have not converged; continue with Newton method */
106:   while (tao->reason == TAO_CONTINUE_ITERATING) {
107:     /* Call general purpose update function */
108:     if (tao->ops->update) {
109:       (*tao->ops->update)(tao, tao->niter, tao->user_update);
110:     }

112:     if (needH && bnk->inactive_idx) {
113:       /* Take BNCG steps (if enabled) to trade-off Hessian evaluations for more gradient evaluations */
114:       TaoBNKTakeCGSteps(tao, &cgTerminate);
115:       if (cgTerminate) {
116:         tao->reason = bnk->bncg->reason;
117:         return 0;
118:       }
119:       /* Compute the hessian and update the BFGS preconditioner at the new iterate */
120:       (*bnk->computehessian)(tao);
121:       needH = PETSC_FALSE;
122:     }

124:     /* Store current solution before it changes */
125:     bnk->fold = bnk->f;
126:     VecCopy(tao->solution, bnk->Xold);
127:     VecCopy(tao->gradient, bnk->Gold);
128:     VecCopy(bnk->unprojected_gradient, bnk->unprojected_gradient_old);

130:     /* Enter into trust region loops */
131:     stepAccepted = PETSC_FALSE;
132:     while (!stepAccepted && tao->reason == TAO_CONTINUE_ITERATING) {
133:       tao->ksp_its=0;

135:       /* Use the common BNK kernel to compute the Newton step (for inactive variables only) */
136:       (*bnk->computestep)(tao, shift, &ksp_reason, &stepType);

138:       /* Temporarily accept the step and project it into the bounds */
139:       VecAXPY(tao->solution, 1.0, tao->stepdirection);
140:       TaoBoundSolution(tao->solution, tao->XL,tao->XU, 0.0, &nDiff, tao->solution);

142:       /* Check if the projection changed the step direction */
143:       if (nDiff > 0) {
144:         /* Projection changed the step, so we have to recompute the step and
145:            the predicted reduction. Leave the trust radius unchanged. */
146:         VecCopy(tao->solution, tao->stepdirection);
147:         VecAXPY(tao->stepdirection, -1.0, bnk->Xold);
148:         TaoBNKRecomputePred(tao, tao->stepdirection, &prered);
149:       } else {
150:         /* Step did not change, so we can just recover the pre-computed prediction */
151:         KSPCGGetObjFcn(tao->ksp, &prered);
152:       }
153:       prered = -prered;

155:       /* Compute the actual reduction and update the trust radius */
156:       TaoComputeObjective(tao, tao->solution, &bnk->f);
158:       actred = bnk->fold - bnk->f;
159:       oldTrust = tao->trust;
160:       TaoBNKUpdateTrustRadius(tao, prered, actred, bnk->update_type, stepType, &stepAccepted);

162:       if (stepAccepted) {
163:         /* Step is good, evaluate the gradient and flip the need-Hessian switch */
164:         steplen = 1.0;
165:         needH = PETSC_TRUE;
166:         ++bnk->newt;
167:         TaoComputeGradient(tao, tao->solution, bnk->unprojected_gradient);
168:         TaoBNKEstimateActiveSet(tao, bnk->as_type);
169:         VecCopy(bnk->unprojected_gradient, tao->gradient);
170:         VecISSet(tao->gradient, bnk->active_idx, 0.0);
171:         TaoGradientNorm(tao, tao->gradient, NORM_2, &bnk->gnorm);
172:       } else {
173:         /* Step is bad, revert old solution and re-solve with new radius*/
174:         steplen = 0.0;
175:         needH = PETSC_FALSE;
176:         bnk->f = bnk->fold;
177:         VecCopy(bnk->Xold, tao->solution);
178:         VecCopy(bnk->Gold, tao->gradient);
179:         VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient);
180:         if (oldTrust == tao->trust) {
181:           /* Can't change the radius anymore so just terminate */
182:           tao->reason = TAO_DIVERGED_TR_REDUCTION;
183:         }
184:       }

186:       /*  Check for termination */
187:       VecFischer(tao->solution, bnk->unprojected_gradient, tao->XL, tao->XU, bnk->W);
188:       VecNorm(bnk->W, NORM_2, &resnorm);
190:       ++tao->niter; /* temporarily increment the iteration number */
191:       TaoLogConvergenceHistory(tao, bnk->f, resnorm, 0.0, tao->ksp_its);
192:       TaoMonitor(tao, tao->niter, bnk->f, resnorm, 0.0, steplen);
193:       (*tao->ops->convergencetest)(tao, tao->cnvP);
194:       --tao->niter; /* decrement */
195:     }
196:     ++tao->niter; /* finished iterating */
197:   }
198:   return 0;
199: }

201: /*------------------------------------------------------------*/
202: static PetscErrorCode TaoSetUp_BNTR(Tao tao)
203: {
204:   KSP               ksp;
205:   PetscVoidFunction valid;

207:   TaoSetUp_BNK(tao);
208:   TaoGetKSP(tao,&ksp);
209:   PetscObjectQueryFunction((PetscObject)ksp,"KSPCGSetRadius_C",&valid);
211:   return 0;
212: }

214: /*------------------------------------------------------------*/

216: static PetscErrorCode TaoSetFromOptions_BNTR(PetscOptionItems *PetscOptionsObject,Tao tao)
217: {
218:   TAO_BNK        *bnk = (TAO_BNK *)tao->data;

220:   TaoSetFromOptions_BNK(PetscOptionsObject, tao);
221:   if (bnk->update_type == BNK_UPDATE_STEP) bnk->update_type = BNK_UPDATE_REDUCTION;
222:   return 0;
223: }

225: /*------------------------------------------------------------*/
226: /*MC
227:   TAOBNTR - Bounded Newton Trust Region for nonlinear minimization with bound constraints.

229:   Options Database Keys:
230: + -tao_bnk_max_cg_its - maximum number of bounded conjugate-gradient iterations taken in each Newton loop
231: . -tao_bnk_init_type - trust radius initialization method ("constant", "direction", "interpolation")
232: . -tao_bnk_update_type - trust radius update method ("step", "direction", "interpolation")
233: - -tao_bnk_as_type - active-set estimation method ("none", "bertsekas")

235:   Level: beginner
236: M*/
237: PETSC_EXTERN PetscErrorCode TaoCreate_BNTR(Tao tao)
238: {
239:   TAO_BNK        *bnk;

241:   TaoCreate_BNK(tao);
242:   tao->ops->solve=TaoSolve_BNTR;
243:   tao->ops->setup=TaoSetUp_BNTR;
244:   tao->ops->setfromoptions=TaoSetFromOptions_BNTR;

246:   bnk = (TAO_BNK *)tao->data;
247:   bnk->update_type = BNK_UPDATE_REDUCTION; /* trust region updates based on predicted/actual reduction */
248:   return 0;
249: }