Actual source code: bntr.c

petsc-3.11.4 2019-09-28
Report Typos and Errors
  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.
  7:  
  8:  ------------------------------------------------------------
  9:  
 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
 17:  
 18:  while niter <= max_it
 19:     niter += 1
 20:     
 21:     if needH
 22:       If max_cg_steps > 0
 23:         x_k, g_k, pg_k = TaoSolve(BNCG)
 24:       end
 25:     
 26:       H_k = TaoComputeHessian(x_k)
 27:       if pc_type == BNK_PC_BFGS
 28:         add correction to BFGS approx
 29:         if scale_type == BNK_SCALE_AHESS
 30:           D = VecMedian(1e-6, abs(diag(H_k)), 1e6)
 31:           scale BFGS with VecReciprocal(D)
 32:         end
 33:       end
 34:       needH = False
 35:     end

 37:     if pc_type = BNK_PC_BFGS
 38:       B_k = BFGS
 39:     else
 40:       B_k = VecMedian(1e-6, abs(diag(H_k)), 1e6)
 41:       B_k = VecReciprocal(B_k)
 42:     end
 43:     w = x_k - VecMedian(x_k - 0.001*B_k*g_k)
 44:     eps = min(eps, norm2(w))
 45:     determine the active and inactive index sets such that
 46:       L = {i : (x_k)_i <= l_i + eps && (g_k)_i > 0}
 47:       U = {i : (x_k)_i >= u_i - eps && (g_k)_i < 0}
 48:       F = {i : l_i = (x_k)_i = u_i}
 49:       A = {L + U + F}
 50:       IA = {i : i not in A}
 51:       
 52:     generate the reduced system Hr_k dr_k = -gr_k for variables in IA
 53:     if pc_type == BNK_PC_BFGS && scale_type == BNK_SCALE_PHESS
 54:       D = VecMedian(1e-6, abs(diag(Hr_k)), 1e6)
 55:       scale BFGS with VecReciprocal(D)
 56:     end
 57:     
 58:     while !stepAccepted
 59:       solve Hr_k dr_k = -gr_k
 60:       set d_k to (l - x) for variables in L, (u - x) for variables in U, and 0 for variables in F
 61:       
 62:       x_{k+1} = VecMedian(x_k + d_k)
 63:       s = x_{k+1} - x_k
 64:       prered = dot(s, 0.5*gr_k - Hr_k*s)
 65:       f_{k+1} = TaoComputeObjective(x_{k+1})
 66:       actred = f_k - f_{k+1}

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

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

 97:   PetscReal                    oldTrust, prered, actred, steplen, resnorm;
 98:   PetscBool                    cgTerminate, needH = PETSC_TRUE, stepAccepted, shift = PETSC_FALSE;
 99:   PetscInt                     stepType, nDiff;
100: 
102:   /* Initialize the preconditioner, KSP solver and trust radius/line search */
103:   tao->reason = TAO_CONTINUE_ITERATING;
104:   TaoBNKInitialize(tao, bnk->init_type, &needH);
105:   if (tao->reason != TAO_CONTINUE_ITERATING) return(0);

107:   /* Have not converged; continue with Newton method */
108:   while (tao->reason == TAO_CONTINUE_ITERATING) {
109:     /* Call general purpose update function */
110:     if (tao->ops->update) {
111:       (*tao->ops->update)(tao, tao->niter, tao->user_update);
112:     }
113:     ++tao->niter;
114: 
115:     if (needH && bnk->inactive_idx) {
116:       /* Take BNCG steps (if enabled) to trade-off Hessian evaluations for more gradient evaluations */
117:       TaoBNKTakeCGSteps(tao, &cgTerminate);
118:       if (cgTerminate) {
119:         tao->reason = bnk->bncg->reason;
120:         return(0);
121:       }
122:       /* Compute the hessian and update the BFGS preconditioner at the new iterate */
123:       (*bnk->computehessian)(tao);
124:       needH = PETSC_FALSE;
125:     }
126: 
127:     /* Store current solution before it changes */
128:     bnk->fold = bnk->f;
129:     VecCopy(tao->solution, bnk->Xold);
130:     VecCopy(tao->gradient, bnk->Gold);
131:     VecCopy(bnk->unprojected_gradient, bnk->unprojected_gradient_old);
132: 
133:     /* Enter into trust region loops */
134:     stepAccepted = PETSC_FALSE;
135:     while (!stepAccepted && tao->reason == TAO_CONTINUE_ITERATING) {
136:       tao->ksp_its=0;
137: 
138:       /* Use the common BNK kernel to compute the Newton step (for inactive variables only) */
139:       (*bnk->computestep)(tao, shift, &ksp_reason, &stepType);

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

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

158:       /* Compute the actual reduction and update the trust radius */
159:       TaoComputeObjective(tao, tao->solution, &bnk->f);
160:       if (PetscIsInfOrNanReal(bnk->f)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
161:       actred = bnk->fold - bnk->f;
162:       oldTrust = tao->trust;
163:       TaoBNKUpdateTrustRadius(tao, prered, actred, bnk->update_type, stepType, &stepAccepted);

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

189:       /*  Check for termination */
190:       VecFischer(tao->solution, bnk->unprojected_gradient, tao->XL, tao->XU, bnk->W);
191:       VecNorm(bnk->W, NORM_2, &resnorm);
192:       if (PetscIsInfOrNanReal(resnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
193:       TaoLogConvergenceHistory(tao, bnk->f, resnorm, 0.0, tao->ksp_its);
194:       TaoMonitor(tao, tao->niter, bnk->f, resnorm, 0.0, steplen);
195:       (*tao->ops->convergencetest)(tao, tao->cnvP);
196:     }
197:   }
198:   return(0);
199: }

201: /*------------------------------------------------------------*/

203: static PetscErrorCode TaoSetFromOptions_BNTR(PetscOptionItems *PetscOptionsObject,Tao tao)
204: {
205:   TAO_BNK        *bnk = (TAO_BNK *)tao->data;

209:   TaoSetFromOptions_BNK(PetscOptionsObject, tao);
210:   if (bnk->update_type == BNK_UPDATE_STEP) bnk->update_type = BNK_UPDATE_REDUCTION;
211:   if (!bnk->is_nash && !bnk->is_stcg && !bnk->is_gltr) SETERRQ(PETSC_COMM_SELF,1,"Must use a trust-region CG method for KSP (KSPNASH, KSPSTCG, KSPGLTR)");
212:   return(0);
213: }

215: /*------------------------------------------------------------*/
216: /*MC
217:   TAOBNTR - Bounded Newton Trust Region for nonlinear minimization with bound constraints.

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

225:   Level: beginner
226: M*/
227: PETSC_EXTERN PetscErrorCode TaoCreate_BNTR(Tao tao)
228: {
229:   TAO_BNK        *bnk;
231: 
233:   TaoCreate_BNK(tao);
234:   tao->ops->solve=TaoSolve_BNTR;
235:   tao->ops->setfromoptions=TaoSetFromOptions_BNTR;
236: 
237:   bnk = (TAO_BNK *)tao->data;
238:   bnk->update_type = BNK_UPDATE_REDUCTION; /* trust region updates based on predicted/actual reduction */
239:   return(0);
240: }