Actual source code: bnls.c

petsc-3.10.5 2019-03-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 line search approach for 
  6:  solving 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_DIRECTION)
 15:  niter = 0
 16:  step_accepted = true
 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

 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
 36:     
 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 p > 0
 54:       Hr_k += p*
 55:     end
 56:     if pc_type == BNK_PC_BFGS && scale_type == BNK_SCALE_PHESS
 57:       D = VecMedian(1e-6, abs(diag(Hr_k)), 1e6)
 58:       scale BFGS with VecReciprocal(D)
 59:     end
 60:     solve Hr_k dr_k = -gr_k
 61:     set d_k to (l - x) for variables in L, (u - x) for variables in U, and 0 for variables in F
 62:     
 63:     if dot(d_k, pg_k)) >= 0 || norm(d_k) == NaN || norm(d_k) == Inf
 64:       dr_k = -BFGS*gr_k for variables in I
 65:       if dot(d_k, pg_k)) >= 0 || norm(d_k) == NaN || norm(d_k) == Inf
 66:         reset the BFGS preconditioner
 67:         calculate scale delta and apply it to BFGS
 68:         dr_k = -BFGS*gr_k for variables in I
 69:         if dot(d_k, pg_k)) >= 0 || norm(d_k) == NaN || norm(d_k) == Inf
 70:           dr_k = -gr_k for variables in I
 71:         end
 72:       end
 73:     end
 74:     
 75:     x_{k+1}, f_{k+1}, g_{k+1}, ls_failed = TaoBNKPerformLineSearch()
 76:     if ls_failed
 77:       f_{k+1} = f_k
 78:       x_{k+1} = x_k
 79:       g_{k+1} = g_k
 80:       pg_{k+1} = pg_k
 81:       terminate
 82:     else
 83:       pg_{k+1} = project(g_{k+1})
 84:       count the accepted step type (Newton, BFGS, scaled grad or grad)
 85:     end 
 86:     
 87:     check convergence at pg_{k+1}
 88:  end
 89: */

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

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

108:   /* Have not converged; continue with Newton method */
109:   while (tao->reason == TAO_CONTINUE_ITERATING) {
110:     ++tao->niter;
111: 
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:     }
123: 
124:     /* Use the common BNK kernel to compute the safeguarded Newton step (for inactive variables only) */
125:     (*bnk->computestep)(tao, shift, &ksp_reason, &stepType);
126:     TaoBNKSafeguardStep(tao, ksp_reason, &stepType);

128:     /* Store current solution before it changes */
129:     bnk->fold = bnk->f;
130:     VecCopy(tao->solution, bnk->Xold);
131:     VecCopy(tao->gradient, bnk->Gold);
132:     VecCopy(bnk->unprojected_gradient, bnk->unprojected_gradient_old);
133: 
134:     /* Trigger the line search */
135:     TaoBNKPerformLineSearch(tao, &stepType, &steplen, &ls_reason);

137:     if (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER) {
138:       /* Failed to find an improving point */
139:       needH = PETSC_FALSE;
140:       bnk->f = bnk->fold;
141:       VecCopy(bnk->Xold, tao->solution);
142:       VecCopy(bnk->Gold, tao->gradient);
143:       VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient);
144:       steplen = 0.0;
145:       tao->reason = TAO_DIVERGED_LS_FAILURE;
146:     } else {
147:       /* new iterate so we need to recompute the Hessian */
148:       needH = PETSC_TRUE;
149:       /* compute the projected gradient */
150:       TaoBNKEstimateActiveSet(tao, bnk->as_type);
151:       VecCopy(bnk->unprojected_gradient, tao->gradient);
152:       VecISSet(tao->gradient, bnk->active_idx, 0.0);
153:       TaoGradientNorm(tao, tao->gradient, NORM_2, &bnk->gnorm);
154:       /* update the trust radius based on the step length */
155:       TaoBNKUpdateTrustRadius(tao, 0.0, 0.0, BNK_UPDATE_STEP, stepType, &stepAccepted);
156:       /* count the accepted step type */
157:       TaoBNKAddStepCounts(tao, stepType);
158:       /* active BNCG recycling for next iteration */
159:       TaoBNCGSetRecycleFlag(bnk->bncg, PETSC_TRUE);
160:     }

162:     /*  Check for termination */
163:     VecFischer(tao->solution, bnk->unprojected_gradient, tao->XL, tao->XU, bnk->W);
164:     VecNorm(bnk->W, NORM_2, &resnorm);
165:     if (PetscIsInfOrNanReal(resnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
166:     TaoLogConvergenceHistory(tao, bnk->f, resnorm, 0.0, tao->ksp_its);
167:     TaoMonitor(tao, tao->niter, bnk->f, resnorm, 0.0, steplen);
168:     (*tao->ops->convergencetest)(tao, tao->cnvP);
169:   }
170:   return(0);
171: }

173: /*------------------------------------------------------------*/

175: PETSC_EXTERN PetscErrorCode TaoCreate_BNLS(Tao tao)
176: {
177:   TAO_BNK        *bnk;
179: 
181:   TaoCreate_BNK(tao);
182:   tao->ops->solve = TaoSolve_BNLS;
183: 
184:   bnk = (TAO_BNK *)tao->data;
185:   bnk->init_type = BNK_INIT_DIRECTION;
186:   bnk->update_type = BNK_UPDATE_STEP; /* trust region updates based on line search step length */
187:   return(0);
188: }