Actual source code: bnk.c

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

  5: static const char *BNK_INIT[64]   = {"constant", "direction", "interpolation"};
  6: static const char *BNK_UPDATE[64] = {"step", "reduction", "interpolation"};
  7: static const char *BNK_AS[64]     = {"none", "bertsekas"};

  9: /* Extracts from the full Hessian the part associated with the current bnk->inactive_idx and set the PCLMVM preconditioner */

 11: static PetscErrorCode TaoBNKComputeSubHessian(Tao tao)
 12: {
 13:   TAO_BNK *bnk = (TAO_BNK *)tao->data;

 15:   PetscFunctionBegin;
 16:   PetscCall(MatDestroy(&bnk->Hpre_inactive));
 17:   PetscCall(MatDestroy(&bnk->H_inactive));
 18:   if (bnk->active_idx) {
 19:     PetscCall(MatCreateSubMatrix(tao->hessian, bnk->inactive_idx, bnk->inactive_idx, MAT_INITIAL_MATRIX, &bnk->H_inactive));
 20:     if (tao->hessian == tao->hessian_pre) {
 21:       PetscCall(PetscObjectReference((PetscObject)bnk->H_inactive));
 22:       bnk->Hpre_inactive = bnk->H_inactive;
 23:     } else {
 24:       PetscCall(MatCreateSubMatrix(tao->hessian_pre, bnk->inactive_idx, bnk->inactive_idx, MAT_INITIAL_MATRIX, &bnk->Hpre_inactive));
 25:     }
 26:     if (bnk->bfgs_pre) PetscCall(PCLMVMSetIS(bnk->bfgs_pre, bnk->inactive_idx));
 27:   } else {
 28:     PetscCall(PetscObjectReference((PetscObject)tao->hessian));
 29:     bnk->H_inactive = tao->hessian;
 30:     PetscCall(PetscObjectReference((PetscObject)tao->hessian_pre));
 31:     bnk->Hpre_inactive = tao->hessian_pre;
 32:     if (bnk->bfgs_pre) PetscCall(PCLMVMClearIS(bnk->bfgs_pre));
 33:   }
 34:   PetscFunctionReturn(PETSC_SUCCESS);
 35: }

 37: /* Initializes the KSP solver, the BFGS preconditioner, and the initial trust radius estimation */

 39: PetscErrorCode TaoBNKInitialize(Tao tao, PetscInt initType, PetscBool *needH)
 40: {
 41:   TAO_BNK  *bnk = (TAO_BNK *)tao->data;
 42:   PC        pc;
 43:   PetscReal f_min, ftrial, prered, actred, kappa, sigma, resnorm;
 44:   PetscReal tau, tau_1, tau_2, tau_max, tau_min, max_radius;
 45:   PetscBool is_bfgs, is_jacobi, is_symmetric, sym_set;
 46:   PetscInt  n, N, nDiff;
 47:   PetscInt  i_max = 5;
 48:   PetscInt  j_max = 1;
 49:   PetscInt  i, j;
 50:   PetscBool kspTR;

 52:   PetscFunctionBegin;
 53:   /* Project the current point onto the feasible set */
 54:   PetscCall(TaoComputeVariableBounds(tao));
 55:   PetscCall(TaoSetVariableBounds(bnk->bncg, tao->XL, tao->XU));
 56:   if (tao->bounded) PetscCall(TaoLineSearchSetVariableBounds(tao->linesearch, tao->XL, tao->XU));

 58:   /* Project the initial point onto the feasible region */
 59:   PetscCall(TaoBoundSolution(tao->solution, tao->XL, tao->XU, 0.0, &nDiff, tao->solution));

 61:   /* Check convergence criteria */
 62:   PetscCall(TaoComputeObjectiveAndGradient(tao, tao->solution, &bnk->f, bnk->unprojected_gradient));
 63:   PetscCall(TaoBNKEstimateActiveSet(tao, bnk->as_type));
 64:   PetscCall(VecCopy(bnk->unprojected_gradient, tao->gradient));
 65:   if (bnk->active_idx) PetscCall(VecISSet(tao->gradient, bnk->active_idx, 0.0));
 66:   PetscCall(TaoGradientNorm(tao, tao->gradient, NORM_2, &bnk->gnorm));

 68:   /* Test the initial point for convergence */
 69:   PetscCall(VecFischer(tao->solution, bnk->unprojected_gradient, tao->XL, tao->XU, bnk->W));
 70:   PetscCall(VecNorm(bnk->W, NORM_2, &resnorm));
 71:   PetscCheck(!PetscIsInfOrNanReal(bnk->f) && !PetscIsInfOrNanReal(resnorm), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute function generated infinity or NaN");
 72:   PetscCall(TaoLogConvergenceHistory(tao, bnk->f, resnorm, 0.0, tao->ksp_its));
 73:   PetscCall(TaoMonitor(tao, tao->niter, bnk->f, resnorm, 0.0, 1.0));
 74:   PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
 75:   if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(PETSC_SUCCESS);

 77:   /* Reset KSP stopping reason counters */
 78:   bnk->ksp_atol = 0;
 79:   bnk->ksp_rtol = 0;
 80:   bnk->ksp_dtol = 0;
 81:   bnk->ksp_ctol = 0;
 82:   bnk->ksp_negc = 0;
 83:   bnk->ksp_iter = 0;
 84:   bnk->ksp_othr = 0;

 86:   /* Reset accepted step type counters */
 87:   bnk->tot_cg_its = 0;
 88:   bnk->newt       = 0;
 89:   bnk->bfgs       = 0;
 90:   bnk->sgrad      = 0;
 91:   bnk->grad       = 0;

 93:   /* Initialize the Hessian perturbation */
 94:   bnk->pert = bnk->sval;

 96:   /* Reset initial steplength to zero (this helps BNCG reset its direction internally) */
 97:   PetscCall(VecSet(tao->stepdirection, 0.0));

 99:   /* Allocate the vectors needed for the BFGS approximation */
100:   PetscCall(KSPGetPC(tao->ksp, &pc));
101:   PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCLMVM, &is_bfgs));
102:   PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCJACOBI, &is_jacobi));
103:   if (is_bfgs) {
104:     bnk->bfgs_pre = pc;
105:     PetscCall(PCLMVMGetMatLMVM(bnk->bfgs_pre, &bnk->M));
106:     PetscCall(VecGetLocalSize(tao->solution, &n));
107:     PetscCall(VecGetSize(tao->solution, &N));
108:     PetscCall(MatSetSizes(bnk->M, n, n, N, N));
109:     PetscCall(MatLMVMAllocate(bnk->M, tao->solution, bnk->unprojected_gradient));
110:     PetscCall(MatIsSymmetricKnown(bnk->M, &sym_set, &is_symmetric));
111:     PetscCheck(sym_set && is_symmetric, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM matrix in the LMVM preconditioner must be symmetric.");
112:   } else if (is_jacobi) PetscCall(PCJacobiSetUseAbs(pc, PETSC_TRUE));

114:   /* Prepare the min/max vectors for safeguarding diagonal scales */
115:   PetscCall(VecSet(bnk->Diag_min, bnk->dmin));
116:   PetscCall(VecSet(bnk->Diag_max, bnk->dmax));

118:   /* Initialize trust-region radius.  The initialization is only performed
119:      when we are using Nash, Steihaug-Toint or the Generalized Lanczos method. */
120:   *needH = PETSC_TRUE;
121:   PetscCall(PetscObjectHasFunction((PetscObject)tao->ksp, "KSPCGSetRadius_C", &kspTR));
122:   if (kspTR) {
123:     switch (initType) {
124:     case BNK_INIT_CONSTANT:
125:       /* Use the initial radius specified */
126:       tao->trust = tao->trust0;
127:       break;

129:     case BNK_INIT_INTERPOLATION:
130:       /* Use interpolation based on the initial Hessian */
131:       max_radius = 0.0;
132:       tao->trust = tao->trust0;
133:       for (j = 0; j < j_max; ++j) {
134:         f_min = bnk->f;
135:         sigma = 0.0;

137:         if (*needH) {
138:           /* Compute the Hessian at the new step, and extract the inactive subsystem */
139:           PetscCall((*bnk->computehessian)(tao));
140:           PetscCall(TaoBNKEstimateActiveSet(tao, BNK_AS_NONE));
141:           PetscCall(TaoBNKComputeSubHessian(tao));
142:           *needH = PETSC_FALSE;
143:         }

145:         for (i = 0; i < i_max; ++i) {
146:           /* Take a steepest descent step and snap it to bounds */
147:           PetscCall(VecCopy(tao->solution, bnk->Xold));
148:           PetscCall(VecAXPY(tao->solution, -tao->trust / bnk->gnorm, tao->gradient));
149:           PetscCall(TaoBoundSolution(tao->solution, tao->XL, tao->XU, 0.0, &nDiff, tao->solution));
150:           /* Compute the step we actually accepted */
151:           PetscCall(VecCopy(tao->solution, bnk->W));
152:           PetscCall(VecAXPY(bnk->W, -1.0, bnk->Xold));
153:           /* Compute the objective at the trial */
154:           PetscCall(TaoComputeObjective(tao, tao->solution, &ftrial));
155:           PetscCheck(!PetscIsInfOrNanReal(bnk->f), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute function generated infinity or NaN");
156:           PetscCall(VecCopy(bnk->Xold, tao->solution));
157:           if (PetscIsInfOrNanReal(ftrial)) {
158:             tau = bnk->gamma1_i;
159:           } else {
160:             if (ftrial < f_min) {
161:               f_min = ftrial;
162:               sigma = -tao->trust / bnk->gnorm;
163:             }

165:             /* Compute the predicted and actual reduction */
166:             if (bnk->active_idx) {
167:               PetscCall(VecGetSubVector(bnk->W, bnk->inactive_idx, &bnk->X_inactive));
168:               PetscCall(VecGetSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work));
169:             } else {
170:               bnk->X_inactive    = bnk->W;
171:               bnk->inactive_work = bnk->Xwork;
172:             }
173:             PetscCall(MatMult(bnk->H_inactive, bnk->X_inactive, bnk->inactive_work));
174:             PetscCall(VecDot(bnk->X_inactive, bnk->inactive_work, &prered));
175:             if (bnk->active_idx) {
176:               PetscCall(VecRestoreSubVector(bnk->W, bnk->inactive_idx, &bnk->X_inactive));
177:               PetscCall(VecRestoreSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work));
178:             }
179:             prered = tao->trust * (bnk->gnorm - 0.5 * tao->trust * prered / (bnk->gnorm * bnk->gnorm));
180:             actred = bnk->f - ftrial;
181:             if ((PetscAbsScalar(actred) <= bnk->epsilon) && (PetscAbsScalar(prered) <= bnk->epsilon)) {
182:               kappa = 1.0;
183:             } else {
184:               kappa = actred / prered;
185:             }

187:             tau_1   = bnk->theta_i * bnk->gnorm * tao->trust / (bnk->theta_i * bnk->gnorm * tao->trust + (1.0 - bnk->theta_i) * prered - actred);
188:             tau_2   = bnk->theta_i * bnk->gnorm * tao->trust / (bnk->theta_i * bnk->gnorm * tao->trust - (1.0 + bnk->theta_i) * prered + actred);
189:             tau_min = PetscMin(tau_1, tau_2);
190:             tau_max = PetscMax(tau_1, tau_2);

192:             if (PetscAbsScalar(kappa - (PetscReal)1.0) <= bnk->mu1_i) {
193:               /*  Great agreement */
194:               max_radius = PetscMax(max_radius, tao->trust);

196:               if (tau_max < 1.0) {
197:                 tau = bnk->gamma3_i;
198:               } else if (tau_max > bnk->gamma4_i) {
199:                 tau = bnk->gamma4_i;
200:               } else {
201:                 tau = tau_max;
202:               }
203:             } else if (PetscAbsScalar(kappa - (PetscReal)1.0) <= bnk->mu2_i) {
204:               /*  Good agreement */
205:               max_radius = PetscMax(max_radius, tao->trust);

207:               if (tau_max < bnk->gamma2_i) {
208:                 tau = bnk->gamma2_i;
209:               } else if (tau_max > bnk->gamma3_i) {
210:                 tau = bnk->gamma3_i;
211:               } else {
212:                 tau = tau_max;
213:               }
214:             } else {
215:               /*  Not good agreement */
216:               if (tau_min > 1.0) {
217:                 tau = bnk->gamma2_i;
218:               } else if (tau_max < bnk->gamma1_i) {
219:                 tau = bnk->gamma1_i;
220:               } else if ((tau_min < bnk->gamma1_i) && (tau_max >= 1.0)) {
221:                 tau = bnk->gamma1_i;
222:               } else if ((tau_1 >= bnk->gamma1_i) && (tau_1 < 1.0) && ((tau_2 < bnk->gamma1_i) || (tau_2 >= 1.0))) {
223:                 tau = tau_1;
224:               } else if ((tau_2 >= bnk->gamma1_i) && (tau_2 < 1.0) && ((tau_1 < bnk->gamma1_i) || (tau_2 >= 1.0))) {
225:                 tau = tau_2;
226:               } else {
227:                 tau = tau_max;
228:               }
229:             }
230:           }
231:           tao->trust = tau * tao->trust;
232:         }

234:         if (f_min < bnk->f) {
235:           /* We accidentally found a solution better than the initial, so accept it */
236:           bnk->f = f_min;
237:           PetscCall(VecCopy(tao->solution, bnk->Xold));
238:           PetscCall(VecAXPY(tao->solution, sigma, tao->gradient));
239:           PetscCall(TaoBoundSolution(tao->solution, tao->XL, tao->XU, 0.0, &nDiff, tao->solution));
240:           PetscCall(VecCopy(tao->solution, tao->stepdirection));
241:           PetscCall(VecAXPY(tao->stepdirection, -1.0, bnk->Xold));
242:           PetscCall(TaoComputeGradient(tao, tao->solution, bnk->unprojected_gradient));
243:           PetscCall(TaoBNKEstimateActiveSet(tao, bnk->as_type));
244:           PetscCall(VecCopy(bnk->unprojected_gradient, tao->gradient));
245:           if (bnk->active_idx) PetscCall(VecISSet(tao->gradient, bnk->active_idx, 0.0));
246:           /* Compute gradient at the new iterate and flip switch to compute the Hessian later */
247:           PetscCall(TaoGradientNorm(tao, tao->gradient, NORM_2, &bnk->gnorm));
248:           *needH = PETSC_TRUE;
249:           /* Test the new step for convergence */
250:           PetscCall(VecFischer(tao->solution, bnk->unprojected_gradient, tao->XL, tao->XU, bnk->W));
251:           PetscCall(VecNorm(bnk->W, NORM_2, &resnorm));
252:           PetscCheck(!PetscIsInfOrNanReal(resnorm), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute function generated infinity or NaN");
253:           PetscCall(TaoLogConvergenceHistory(tao, bnk->f, resnorm, 0.0, tao->ksp_its));
254:           PetscCall(TaoMonitor(tao, tao->niter, bnk->f, resnorm, 0.0, 1.0));
255:           PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
256:           if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(PETSC_SUCCESS);
257:           /* active BNCG recycling early because we have a stepdirection computed */
258:           PetscCall(TaoSetRecycleHistory(bnk->bncg, PETSC_TRUE));
259:         }
260:       }
261:       tao->trust = PetscMax(tao->trust, max_radius);

263:       /* Ensure that the trust radius is within the limits */
264:       tao->trust = PetscMax(tao->trust, bnk->min_radius);
265:       tao->trust = PetscMin(tao->trust, bnk->max_radius);
266:       break;

268:     default:
269:       /* Norm of the first direction will initialize radius */
270:       tao->trust = 0.0;
271:       break;
272:     }
273:   }
274:   PetscFunctionReturn(PETSC_SUCCESS);
275: }

277: /* Computes the exact Hessian and extracts its subHessian */

279: PetscErrorCode TaoBNKComputeHessian(Tao tao)
280: {
281:   TAO_BNK *bnk = (TAO_BNK *)tao->data;

283:   PetscFunctionBegin;
284:   /* Compute the Hessian */
285:   PetscCall(TaoComputeHessian(tao, tao->solution, tao->hessian, tao->hessian_pre));
286:   /* Add a correction to the BFGS preconditioner */
287:   if (bnk->M) PetscCall(MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient));
288:   /* Prepare the reduced sub-matrices for the inactive set */
289:   PetscCall(TaoBNKComputeSubHessian(tao));
290:   PetscFunctionReturn(PETSC_SUCCESS);
291: }

293: /* Routine for estimating the active set */

295: PetscErrorCode TaoBNKEstimateActiveSet(Tao tao, PetscInt asType)
296: {
297:   TAO_BNK  *bnk = (TAO_BNK *)tao->data;
298:   PetscBool hessComputed, diagExists, hadactive;

300:   PetscFunctionBegin;
301:   hadactive = bnk->active_idx ? PETSC_TRUE : PETSC_FALSE;
302:   switch (asType) {
303:   case BNK_AS_NONE:
304:     PetscCall(ISDestroy(&bnk->inactive_idx));
305:     PetscCall(VecWhichInactive(tao->XL, tao->solution, bnk->unprojected_gradient, tao->XU, PETSC_TRUE, &bnk->inactive_idx));
306:     PetscCall(ISDestroy(&bnk->active_idx));
307:     PetscCall(ISComplementVec(bnk->inactive_idx, tao->solution, &bnk->active_idx));
308:     break;

310:   case BNK_AS_BERTSEKAS:
311:     /* Compute the trial step vector with which we will estimate the active set at the next iteration */
312:     if (bnk->M) {
313:       /* If the BFGS matrix is available, we will construct a trial step with it */
314:       PetscCall(MatSolve(bnk->M, bnk->unprojected_gradient, bnk->W));
315:     } else {
316:       hessComputed = diagExists = PETSC_FALSE;
317:       if (tao->hessian) PetscCall(MatAssembled(tao->hessian, &hessComputed));
318:       if (hessComputed) PetscCall(MatHasOperation(tao->hessian, MATOP_GET_DIAGONAL, &diagExists));
319:       if (diagExists) {
320:         /* BFGS preconditioner doesn't exist so let's invert the absolute diagonal of the Hessian instead onto the gradient */
321:         PetscCall(MatGetDiagonal(tao->hessian, bnk->Xwork));
322:         PetscCall(VecAbs(bnk->Xwork));
323:         PetscCall(VecMedian(bnk->Diag_min, bnk->Xwork, bnk->Diag_max, bnk->Xwork));
324:         PetscCall(VecReciprocal(bnk->Xwork));
325:         PetscCall(VecPointwiseMult(bnk->W, bnk->Xwork, bnk->unprojected_gradient));
326:       } else {
327:         /* If the Hessian or its diagonal does not exist, we will simply use gradient step */
328:         PetscCall(VecCopy(bnk->unprojected_gradient, bnk->W));
329:       }
330:     }
331:     PetscCall(VecScale(bnk->W, -1.0));
332:     PetscCall(TaoEstimateActiveBounds(tao->solution, tao->XL, tao->XU, bnk->unprojected_gradient, bnk->W, bnk->Xwork, bnk->as_step, &bnk->as_tol, &bnk->active_lower, &bnk->active_upper, &bnk->active_fixed, &bnk->active_idx, &bnk->inactive_idx));
333:     break;

335:   default:
336:     break;
337:   }
338:   bnk->resetksp = (PetscBool)(bnk->active_idx || hadactive); /* inactive Hessian size may have changed, need to reset operators */
339:   PetscFunctionReturn(PETSC_SUCCESS);
340: }

342: /* Routine for bounding the step direction */

344: PetscErrorCode TaoBNKBoundStep(Tao tao, PetscInt asType, Vec step)
345: {
346:   TAO_BNK *bnk = (TAO_BNK *)tao->data;

348:   PetscFunctionBegin;
349:   switch (asType) {
350:   case BNK_AS_NONE:
351:     if (bnk->active_idx) PetscCall(VecISSet(step, bnk->active_idx, 0.0));
352:     break;
353:   case BNK_AS_BERTSEKAS:
354:     PetscCall(TaoBoundStep(tao->solution, tao->XL, tao->XU, bnk->active_lower, bnk->active_upper, bnk->active_fixed, 1.0, step));
355:     break;
356:   default:
357:     break;
358:   }
359:   PetscFunctionReturn(PETSC_SUCCESS);
360: }

362: /* Routine for taking a finite number of BNCG iterations to
363:    accelerate Newton convergence.

365:    In practice, this approach simply trades off Hessian evaluations
366:    for more gradient evaluations.
367: */

369: PetscErrorCode TaoBNKTakeCGSteps(Tao tao, PetscBool *terminate)
370: {
371:   TAO_BNK *bnk = (TAO_BNK *)tao->data;

373:   PetscFunctionBegin;
374:   *terminate = PETSC_FALSE;
375:   if (bnk->max_cg_its > 0) {
376:     /* Copy the current function value (important vectors are already shared) */
377:     bnk->bncg_ctx->f = bnk->f;
378:     /* Take some small finite number of BNCG iterations */
379:     PetscCall(TaoSolve(bnk->bncg));
380:     /* Add the number of gradient and function evaluations to the total     *
381:      * Note: nfuncs are not copied as tao and subsolvers share same TaoTerm */
382:     bnk->tot_cg_its += bnk->bncg->niter;
383:     /* Extract the BNCG function value out and save it into BNK */
384:     bnk->f = bnk->bncg_ctx->f;
385:     if (bnk->bncg->reason == TAO_CONVERGED_GATOL || bnk->bncg->reason == TAO_CONVERGED_GRTOL || bnk->bncg->reason == TAO_CONVERGED_GTTOL || bnk->bncg->reason == TAO_CONVERGED_MINF) {
386:       *terminate = PETSC_TRUE;
387:     } else {
388:       PetscCall(TaoBNKEstimateActiveSet(tao, bnk->as_type));
389:     }
390:   }
391:   PetscFunctionReturn(PETSC_SUCCESS);
392: }

394: /* Routine for computing the Newton step. */

396: PetscErrorCode TaoBNKComputeStep(Tao tao, PetscBool shift, KSPConvergedReason *ksp_reason, PetscInt *step_type)
397: {
398:   TAO_BNK  *bnk         = (TAO_BNK *)tao->data;
399:   PetscInt  bfgsUpdates = 0;
400:   PetscInt  kspits;
401:   PetscBool is_lmvm;
402:   PetscBool kspTR;

404:   PetscFunctionBegin;
405:   /* If there are no inactive variables left, save some computation and return an adjusted zero step
406:      that has (l-x) and (u-x) for lower and upper bounded variables. */
407:   if (!bnk->inactive_idx) {
408:     PetscCall(VecSet(tao->stepdirection, 0.0));
409:     PetscCall(TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection));
410:     PetscFunctionReturn(PETSC_SUCCESS);
411:   }

413:   /* Shift the reduced Hessian matrix */
414:   if (shift && bnk->pert > 0) {
415:     PetscCall(PetscObjectTypeCompare((PetscObject)tao->hessian, MATLMVM, &is_lmvm));
416:     if (is_lmvm) {
417:       PetscCall(MatShift(tao->hessian, bnk->pert));
418:     } else {
419:       PetscCall(MatShift(bnk->H_inactive, bnk->pert));
420:       if (bnk->H_inactive != bnk->Hpre_inactive) PetscCall(MatShift(bnk->Hpre_inactive, bnk->pert));
421:     }
422:   }

424:   /* Solve the Newton system of equations */
425:   tao->ksp_its = 0;
426:   PetscCall(VecSet(tao->stepdirection, 0.0));
427:   if (bnk->resetksp) {
428:     PetscCall(KSPReset(tao->ksp));
429:     PetscCall(KSPResetFromOptions(tao->ksp));
430:     bnk->resetksp = PETSC_FALSE;
431:   }
432:   PetscCall(KSPSetOperators(tao->ksp, bnk->H_inactive, bnk->Hpre_inactive));
433:   PetscCall(VecCopy(bnk->unprojected_gradient, bnk->Gwork));
434:   if (bnk->active_idx) {
435:     PetscCall(VecGetSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive));
436:     PetscCall(VecGetSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive));
437:   } else {
438:     bnk->G_inactive = bnk->unprojected_gradient;
439:     bnk->X_inactive = tao->stepdirection;
440:   }
441:   PetscCall(KSPCGSetRadius(tao->ksp, tao->trust));
442:   PetscCall(KSPSolve(tao->ksp, bnk->G_inactive, bnk->X_inactive));
443:   PetscCall(KSPGetIterationNumber(tao->ksp, &kspits));
444:   tao->ksp_its += kspits;
445:   tao->ksp_tot_its += kspits;
446:   PetscCall(PetscObjectHasFunction((PetscObject)tao->ksp, "KSPCGGetNormD_C", &kspTR));
447:   if (kspTR) {
448:     PetscCall(KSPCGGetNormD(tao->ksp, &bnk->dnorm));

450:     if (0.0 == tao->trust) {
451:       /* Radius was uninitialized; use the norm of the direction */
452:       if (bnk->dnorm > 0.0) {
453:         tao->trust = bnk->dnorm;

455:         /* Modify the radius if it is too large or small */
456:         tao->trust = PetscMax(tao->trust, bnk->min_radius);
457:         tao->trust = PetscMin(tao->trust, bnk->max_radius);
458:       } else {
459:         /* The direction was bad; set radius to default value and re-solve
460:            the trust-region subproblem to get a direction */
461:         tao->trust = tao->trust0;

463:         /* Modify the radius if it is too large or small */
464:         tao->trust = PetscMax(tao->trust, bnk->min_radius);
465:         tao->trust = PetscMin(tao->trust, bnk->max_radius);

467:         PetscCall(KSPCGSetRadius(tao->ksp, tao->trust));
468:         PetscCall(KSPSolve(tao->ksp, bnk->G_inactive, bnk->X_inactive));
469:         PetscCall(KSPGetIterationNumber(tao->ksp, &kspits));
470:         tao->ksp_its += kspits;
471:         tao->ksp_tot_its += kspits;
472:         PetscCall(KSPCGGetNormD(tao->ksp, &bnk->dnorm));

474:         PetscCheck(bnk->dnorm != 0.0, PetscObjectComm((PetscObject)tao), PETSC_ERR_PLIB, "Initial direction zero");
475:       }
476:     }
477:   }
478:   /* Restore sub vectors back */
479:   if (bnk->active_idx) {
480:     PetscCall(VecRestoreSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive));
481:     PetscCall(VecRestoreSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive));
482:   }
483:   /* Make sure the safeguarded fall-back step is zero for actively bounded variables */
484:   PetscCall(VecScale(tao->stepdirection, -1.0));
485:   PetscCall(TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection));

487:   /* Record convergence reasons */
488:   PetscCall(KSPGetConvergedReason(tao->ksp, ksp_reason));
489:   if (KSP_CONVERGED_ATOL == *ksp_reason) {
490:     ++bnk->ksp_atol;
491:   } else if (KSP_CONVERGED_RTOL == *ksp_reason) {
492:     ++bnk->ksp_rtol;
493:   } else if (KSP_CONVERGED_STEP_LENGTH == *ksp_reason) {
494:     ++bnk->ksp_ctol;
495:   } else if (KSP_CONVERGED_NEG_CURVE == *ksp_reason) {
496:     ++bnk->ksp_negc;
497:   } else if (KSP_DIVERGED_DTOL == *ksp_reason) {
498:     ++bnk->ksp_dtol;
499:   } else if (KSP_DIVERGED_ITS == *ksp_reason) {
500:     ++bnk->ksp_iter;
501:   } else {
502:     ++bnk->ksp_othr;
503:   }

505:   /* Make sure the BFGS preconditioner is healthy */
506:   if (bnk->M) {
507:     PetscCall(MatLMVMGetUpdateCount(bnk->M, &bfgsUpdates));
508:     if ((KSP_DIVERGED_INDEFINITE_PC == *ksp_reason) && (bfgsUpdates > 0)) {
509:       /* Preconditioner is numerically indefinite; reset the approximation. */
510:       PetscCall(MatLMVMReset(bnk->M, PETSC_FALSE));
511:       PetscCall(MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient));
512:     }
513:   }
514:   *step_type = BNK_NEWTON;
515:   PetscFunctionReturn(PETSC_SUCCESS);
516: }

518: /* Routine for recomputing the predicted reduction for a given step vector */

520: PetscErrorCode TaoBNKRecomputePred(Tao tao, Vec S, PetscReal *prered)
521: {
522:   TAO_BNK *bnk = (TAO_BNK *)tao->data;

524:   PetscFunctionBegin;
525:   /* Extract subvectors associated with the inactive set */
526:   if (bnk->active_idx) {
527:     PetscCall(VecGetSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive));
528:     PetscCall(VecGetSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work));
529:     PetscCall(VecGetSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive));
530:   } else {
531:     bnk->X_inactive    = tao->stepdirection;
532:     bnk->inactive_work = bnk->Xwork;
533:     bnk->G_inactive    = bnk->Gwork;
534:   }
535:   /* Recompute the predicted decrease based on the quadratic model */
536:   PetscCall(MatMult(bnk->H_inactive, bnk->X_inactive, bnk->inactive_work));
537:   PetscCall(VecAYPX(bnk->inactive_work, -0.5, bnk->G_inactive));
538:   PetscCall(VecDot(bnk->inactive_work, bnk->X_inactive, prered));
539:   /* Restore the sub vectors */
540:   if (bnk->active_idx) {
541:     PetscCall(VecRestoreSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive));
542:     PetscCall(VecRestoreSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work));
543:     PetscCall(VecRestoreSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive));
544:   }
545:   PetscFunctionReturn(PETSC_SUCCESS);
546: }

548: /* Routine for ensuring that the Newton step is a descent direction.

550:    The step direction falls back onto BFGS, scaled gradient and gradient steps
551:    in the event that the Newton step fails the test.
552: */

554: PetscErrorCode TaoBNKSafeguardStep(Tao tao, KSPConvergedReason ksp_reason, PetscInt *stepType)
555: {
556:   TAO_BNK  *bnk = (TAO_BNK *)tao->data;
557:   PetscReal gdx, e_min;
558:   PetscInt  bfgsUpdates;

560:   PetscFunctionBegin;
561:   switch (*stepType) {
562:   case BNK_NEWTON:
563:     PetscCall(VecDot(tao->stepdirection, tao->gradient, &gdx));
564:     if ((gdx >= 0.0) || PetscIsInfOrNanReal(gdx)) {
565:       /* Newton step is not descent or direction produced infinity or NaN
566:         Update the perturbation for next time */
567:       if (bnk->pert <= 0.0) {
568:         PetscBool is_gltr;

570:         /* Initialize the perturbation */
571:         bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
572:         PetscCall(PetscObjectTypeCompare((PetscObject)tao->ksp, KSPGLTR, &is_gltr));
573:         if (is_gltr) {
574:           PetscCall(KSPGLTRGetMinEig(tao->ksp, &e_min));
575:           bnk->pert = PetscMax(bnk->pert, -e_min);
576:         }
577:       } else {
578:         /* Increase the perturbation */
579:         bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
580:       }

582:       if (!bnk->M) {
583:         /* We don't have the bfgs matrix around and updated
584:           Must use gradient direction in this case */
585:         PetscCall(VecCopy(tao->gradient, tao->stepdirection));
586:         *stepType = BNK_GRADIENT;
587:       } else {
588:         /* Attempt to use the BFGS direction */
589:         PetscCall(MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection));

591:         /* Check for success (descent direction)
592:           NOTE: Negative gdx here means not a descent direction because
593:           the fall-back step is missing a negative sign. */
594:         PetscCall(VecDot(tao->gradient, tao->stepdirection, &gdx));
595:         if ((gdx <= 0.0) || PetscIsInfOrNanReal(gdx)) {
596:           /* BFGS direction is not descent or direction produced not a number
597:             We can assert bfgsUpdates > 1 in this case because
598:             the first solve produces the scaled gradient direction,
599:             which is guaranteed to be descent */

601:           /* Use steepest descent direction (scaled) */
602:           PetscCall(MatLMVMReset(bnk->M, PETSC_FALSE));
603:           PetscCall(MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient));
604:           PetscCall(MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection));

606:           *stepType = BNK_SCALED_GRADIENT;
607:         } else {
608:           PetscCall(MatLMVMGetUpdateCount(bnk->M, &bfgsUpdates));
609:           if (1 == bfgsUpdates) {
610:             /* The first BFGS direction is always the scaled gradient */
611:             *stepType = BNK_SCALED_GRADIENT;
612:           } else {
613:             *stepType = BNK_BFGS;
614:           }
615:         }
616:       }
617:       /* Make sure the safeguarded fall-back step is zero for actively bounded variables */
618:       PetscCall(VecScale(tao->stepdirection, -1.0));
619:       PetscCall(TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection));
620:     } else {
621:       /* Computed Newton step is descent */
622:       switch (ksp_reason) {
623:       case KSP_DIVERGED_NANORINF:
624:       case KSP_DIVERGED_BREAKDOWN:
625:       case KSP_DIVERGED_INDEFINITE_MAT:
626:       case KSP_DIVERGED_INDEFINITE_PC:
627:       case KSP_CONVERGED_NEG_CURVE:
628:         /* Matrix or preconditioner is indefinite; increase perturbation */
629:         if (bnk->pert <= 0.0) {
630:           PetscBool is_gltr;

632:           /* Initialize the perturbation */
633:           bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
634:           PetscCall(PetscObjectTypeCompare((PetscObject)tao->ksp, KSPGLTR, &is_gltr));
635:           if (is_gltr) {
636:             PetscCall(KSPGLTRGetMinEig(tao->ksp, &e_min));
637:             bnk->pert = PetscMax(bnk->pert, -e_min);
638:           }
639:         } else {
640:           /* Increase the perturbation */
641:           bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
642:         }
643:         break;

645:       default:
646:         /* Newton step computation is good; decrease perturbation */
647:         bnk->pert = PetscMin(bnk->psfac * bnk->pert, bnk->pmsfac * bnk->gnorm);
648:         if (bnk->pert < bnk->pmin) bnk->pert = 0.0;
649:         break;
650:       }
651:       *stepType = BNK_NEWTON;
652:     }
653:     break;

655:   case BNK_BFGS:
656:     /* Check for success (descent direction) */
657:     PetscCall(VecDot(tao->stepdirection, tao->gradient, &gdx));
658:     if (gdx >= 0 || PetscIsInfOrNanReal(gdx)) {
659:       /* Step is not descent or solve was not successful
660:          Use steepest descent direction (scaled) */
661:       PetscCall(MatLMVMReset(bnk->M, PETSC_FALSE));
662:       PetscCall(MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient));
663:       PetscCall(MatSolve(bnk->M, tao->gradient, tao->stepdirection));
664:       PetscCall(VecScale(tao->stepdirection, -1.0));
665:       PetscCall(TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection));
666:       *stepType = BNK_SCALED_GRADIENT;
667:     } else {
668:       *stepType = BNK_BFGS;
669:     }
670:     break;

672:   case BNK_SCALED_GRADIENT:
673:     break;

675:   default:
676:     break;
677:   }
678:   PetscFunctionReturn(PETSC_SUCCESS);
679: }

681: /* Routine for performing a bound-projected More-Thuente line search.

683:   Includes fallbacks to BFGS, scaled gradient, and unscaled gradient steps if the
684:   Newton step does not produce a valid step length.
685: */

687: PetscErrorCode TaoBNKPerformLineSearch(Tao tao, PetscInt *stepType, PetscReal *steplen, TaoLineSearchConvergedReason *reason)
688: {
689:   TAO_BNK                     *bnk = (TAO_BNK *)tao->data;
690:   TaoLineSearchConvergedReason ls_reason;
691:   PetscReal                    e_min, gdx;
692:   PetscInt                     bfgsUpdates;

694:   PetscFunctionBegin;
695:   /* Perform the linesearch */
696:   PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, steplen, &ls_reason));
697:   PetscCall(TaoAddLineSearchCounts(tao));

699:   while (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER && *stepType != BNK_SCALED_GRADIENT && *stepType != BNK_GRADIENT) {
700:     /* Linesearch failed, revert solution */
701:     bnk->f = bnk->fold;
702:     PetscCall(VecCopy(bnk->Xold, tao->solution));
703:     PetscCall(VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient));

705:     switch (*stepType) {
706:     case BNK_NEWTON:
707:       /* Failed to obtain acceptable iterate with Newton step
708:          Update the perturbation for next time */
709:       if (bnk->pert <= 0.0) {
710:         PetscBool is_gltr;

712:         /* Initialize the perturbation */
713:         bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
714:         PetscCall(PetscObjectTypeCompare((PetscObject)tao->ksp, KSPGLTR, &is_gltr));
715:         if (is_gltr) {
716:           PetscCall(KSPGLTRGetMinEig(tao->ksp, &e_min));
717:           bnk->pert = PetscMax(bnk->pert, -e_min);
718:         }
719:       } else {
720:         /* Increase the perturbation */
721:         bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
722:       }

724:       if (!bnk->M) {
725:         /* We don't have the bfgs matrix around and being updated
726:            Must use gradient direction in this case */
727:         PetscCall(VecCopy(bnk->unprojected_gradient, tao->stepdirection));
728:         *stepType = BNK_GRADIENT;
729:       } else {
730:         /* Attempt to use the BFGS direction */
731:         PetscCall(MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection));
732:         /* Check for success (descent direction)
733:            NOTE: Negative gdx means not a descent direction because the step here is missing a negative sign. */
734:         PetscCall(VecDot(tao->gradient, tao->stepdirection, &gdx));
735:         if ((gdx <= 0.0) || PetscIsInfOrNanReal(gdx)) {
736:           /* BFGS direction is not descent or direction produced not a number
737:              We can assert bfgsUpdates > 1 in this case
738:              Use steepest descent direction (scaled) */
739:           PetscCall(MatLMVMReset(bnk->M, PETSC_FALSE));
740:           PetscCall(MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient));
741:           PetscCall(MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection));

743:           bfgsUpdates = 1;
744:           *stepType   = BNK_SCALED_GRADIENT;
745:         } else {
746:           PetscCall(MatLMVMGetUpdateCount(bnk->M, &bfgsUpdates));
747:           if (1 == bfgsUpdates) {
748:             /* The first BFGS direction is always the scaled gradient */
749:             *stepType = BNK_SCALED_GRADIENT;
750:           } else {
751:             *stepType = BNK_BFGS;
752:           }
753:         }
754:       }
755:       break;

757:     case BNK_BFGS:
758:       /* Can only enter if pc_type == BNK_PC_BFGS
759:          Failed to obtain acceptable iterate with BFGS step
760:          Attempt to use the scaled gradient direction */
761:       PetscCall(MatLMVMReset(bnk->M, PETSC_FALSE));
762:       PetscCall(MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient));
763:       PetscCall(MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection));

765:       bfgsUpdates = 1;
766:       *stepType   = BNK_SCALED_GRADIENT;
767:       break;
768:     }
769:     /* Make sure the safeguarded fall-back step is zero for actively bounded variables */
770:     PetscCall(VecScale(tao->stepdirection, -1.0));
771:     PetscCall(TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection));

773:     /* Perform one last line search with the fall-back step */
774:     PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, steplen, &ls_reason));
775:     PetscCall(TaoAddLineSearchCounts(tao));
776:   }
777:   *reason = ls_reason;
778:   PetscFunctionReturn(PETSC_SUCCESS);
779: }

781: /* Routine for updating the trust radius.

783:   Function features three different update methods:
784:   1) Line-search step length based
785:   2) Predicted decrease on the CG quadratic model
786:   3) Interpolation
787: */

789: PetscErrorCode TaoBNKUpdateTrustRadius(Tao tao, PetscReal prered, PetscReal actred, PetscInt updateType, PetscInt stepType, PetscBool *accept)
790: {
791:   TAO_BNK *bnk = (TAO_BNK *)tao->data;

793:   PetscReal step, kappa;
794:   PetscReal gdx, tau_1, tau_2, tau_min, tau_max;

796:   PetscFunctionBegin;
797:   /* Update trust region radius */
798:   *accept = PETSC_FALSE;
799:   switch (updateType) {
800:   case BNK_UPDATE_STEP:
801:     *accept = PETSC_TRUE; /* always accept here because line search succeeded */
802:     if (stepType == BNK_NEWTON) {
803:       PetscCall(TaoLineSearchGetStepLength(tao->linesearch, &step));
804:       if (step < bnk->nu1) {
805:         /* Very bad step taken; reduce radius */
806:         tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust);
807:       } else if (step < bnk->nu2) {
808:         /* Reasonably bad step taken; reduce radius */
809:         tao->trust = bnk->omega2 * PetscMin(bnk->dnorm, tao->trust);
810:       } else if (step < bnk->nu3) {
811:         /*  Reasonable step was taken; leave radius alone */
812:         if (bnk->omega3 < 1.0) {
813:           tao->trust = bnk->omega3 * PetscMin(bnk->dnorm, tao->trust);
814:         } else if (bnk->omega3 > 1.0) {
815:           tao->trust = PetscMax(bnk->omega3 * bnk->dnorm, tao->trust);
816:         }
817:       } else if (step < bnk->nu4) {
818:         /*  Full step taken; increase the radius */
819:         tao->trust = PetscMax(bnk->omega4 * bnk->dnorm, tao->trust);
820:       } else {
821:         /*  More than full step taken; increase the radius */
822:         tao->trust = PetscMax(bnk->omega5 * bnk->dnorm, tao->trust);
823:       }
824:     } else {
825:       /*  Newton step was not good; reduce the radius */
826:       tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust);
827:     }
828:     break;

830:   case BNK_UPDATE_REDUCTION:
831:     if (stepType == BNK_NEWTON) {
832:       if ((prered < 0.0) || PetscIsInfOrNanReal(prered)) {
833:         /* The predicted reduction has the wrong sign.  This cannot
834:            happen in infinite precision arithmetic.  Step should
835:            be rejected! */
836:         tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
837:       } else {
838:         if (PetscIsInfOrNanReal(actred)) {
839:           tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
840:         } else {
841:           if ((PetscAbsScalar(actred) <= PetscMax(1.0, PetscAbsScalar(bnk->f)) * bnk->epsilon) && (PetscAbsScalar(prered) <= PetscMax(1.0, PetscAbsScalar(bnk->f)) * bnk->epsilon)) {
842:             kappa = 1.0;
843:           } else {
844:             kappa = actred / prered;
845:           }
846:           /* Accept or reject the step and update radius */
847:           if (kappa < bnk->eta1) {
848:             /* Reject the step */
849:             tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
850:           } else {
851:             /* Accept the step */
852:             *accept = PETSC_TRUE;
853:             /* Update the trust region radius only if the computed step is at the trust radius boundary */
854:             if (bnk->dnorm == tao->trust) {
855:               if (kappa < bnk->eta2) {
856:                 /* Marginal bad step */
857:                 tao->trust = bnk->alpha2 * tao->trust;
858:               } else if (kappa < bnk->eta3) {
859:                 /* Reasonable step */
860:                 tao->trust = bnk->alpha3 * tao->trust;
861:               } else if (kappa < bnk->eta4) {
862:                 /* Good step */
863:                 tao->trust = bnk->alpha4 * tao->trust;
864:               } else {
865:                 /* Very good step */
866:                 tao->trust = bnk->alpha5 * tao->trust;
867:               }
868:             }
869:           }
870:         }
871:       }
872:     } else {
873:       /*  Newton step was not good; reduce the radius */
874:       tao->trust = bnk->alpha1 * PetscMin(bnk->dnorm, tao->trust);
875:     }
876:     break;

878:   default:
879:     if (stepType == BNK_NEWTON) {
880:       if (prered < 0.0) {
881:         /*  The predicted reduction has the wrong sign.  This cannot */
882:         /*  happen in infinite precision arithmetic.  Step should */
883:         /*  be rejected! */
884:         tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
885:       } else {
886:         if (PetscIsInfOrNanReal(actred)) {
887:           tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
888:         } else {
889:           if ((PetscAbsScalar(actred) <= bnk->epsilon) && (PetscAbsScalar(prered) <= bnk->epsilon)) {
890:             kappa = 1.0;
891:           } else {
892:             kappa = actred / prered;
893:           }

895:           PetscCall(VecDot(tao->gradient, tao->stepdirection, &gdx));
896:           tau_1   = bnk->theta * gdx / (bnk->theta * gdx - (1.0 - bnk->theta) * prered + actred);
897:           tau_2   = bnk->theta * gdx / (bnk->theta * gdx + (1.0 + bnk->theta) * prered - actred);
898:           tau_min = PetscMin(tau_1, tau_2);
899:           tau_max = PetscMax(tau_1, tau_2);

901:           if (kappa >= 1.0 - bnk->mu1) {
902:             /*  Great agreement */
903:             *accept = PETSC_TRUE;
904:             if (tau_max < 1.0) {
905:               tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm);
906:             } else if (tau_max > bnk->gamma4) {
907:               tao->trust = PetscMax(tao->trust, bnk->gamma4 * bnk->dnorm);
908:             } else {
909:               tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm);
910:             }
911:           } else if (kappa >= 1.0 - bnk->mu2) {
912:             /*  Good agreement */
913:             *accept = PETSC_TRUE;
914:             if (tau_max < bnk->gamma2) {
915:               tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm);
916:             } else if (tau_max > bnk->gamma3) {
917:               tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm);
918:             } else if (tau_max < 1.0) {
919:               tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm);
920:             } else {
921:               tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm);
922:             }
923:           } else {
924:             /*  Not good agreement */
925:             if (tau_min > 1.0) {
926:               tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm);
927:             } else if (tau_max < bnk->gamma1) {
928:               tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
929:             } else if ((tau_min < bnk->gamma1) && (tau_max >= 1.0)) {
930:               tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
931:             } else if ((tau_1 >= bnk->gamma1) && (tau_1 < 1.0) && ((tau_2 < bnk->gamma1) || (tau_2 >= 1.0))) {
932:               tao->trust = tau_1 * PetscMin(tao->trust, bnk->dnorm);
933:             } else if ((tau_2 >= bnk->gamma1) && (tau_2 < 1.0) && ((tau_1 < bnk->gamma1) || (tau_2 >= 1.0))) {
934:               tao->trust = tau_2 * PetscMin(tao->trust, bnk->dnorm);
935:             } else {
936:               tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm);
937:             }
938:           }
939:         }
940:       }
941:     } else {
942:       /*  Newton step was not good; reduce the radius */
943:       tao->trust = bnk->gamma1 * PetscMin(bnk->dnorm, tao->trust);
944:     }
945:     break;
946:   }
947:   /* Make sure the radius does not violate min and max settings */
948:   tao->trust = PetscMin(tao->trust, bnk->max_radius);
949:   tao->trust = PetscMax(tao->trust, bnk->min_radius);
950:   PetscFunctionReturn(PETSC_SUCCESS);
951: }

953: PetscErrorCode TaoBNKAddStepCounts(Tao tao, PetscInt stepType)
954: {
955:   TAO_BNK *bnk = (TAO_BNK *)tao->data;

957:   PetscFunctionBegin;
958:   switch (stepType) {
959:   case BNK_NEWTON:
960:     ++bnk->newt;
961:     break;
962:   case BNK_BFGS:
963:     ++bnk->bfgs;
964:     break;
965:   case BNK_SCALED_GRADIENT:
966:     ++bnk->sgrad;
967:     break;
968:   case BNK_GRADIENT:
969:     ++bnk->grad;
970:     break;
971:   default:
972:     break;
973:   }
974:   PetscFunctionReturn(PETSC_SUCCESS);
975: }

977: PetscErrorCode TaoSetUp_BNK(Tao tao)
978: {
979:   TAO_BNK *bnk = (TAO_BNK *)tao->data;

981:   PetscFunctionBegin;
982:   if (!tao->gradient) PetscCall(VecDuplicate(tao->solution, &tao->gradient));
983:   if (!tao->stepdirection) PetscCall(VecDuplicate(tao->solution, &tao->stepdirection));
984:   if (!bnk->W) PetscCall(VecDuplicate(tao->solution, &bnk->W));
985:   if (!bnk->Xold) PetscCall(VecDuplicate(tao->solution, &bnk->Xold));
986:   if (!bnk->Gold) PetscCall(VecDuplicate(tao->solution, &bnk->Gold));
987:   if (!bnk->Xwork) PetscCall(VecDuplicate(tao->solution, &bnk->Xwork));
988:   if (!bnk->Gwork) PetscCall(VecDuplicate(tao->solution, &bnk->Gwork));
989:   if (!bnk->unprojected_gradient) PetscCall(VecDuplicate(tao->solution, &bnk->unprojected_gradient));
990:   if (!bnk->unprojected_gradient_old) PetscCall(VecDuplicate(tao->solution, &bnk->unprojected_gradient_old));
991:   if (!bnk->Diag_min) PetscCall(VecDuplicate(tao->solution, &bnk->Diag_min));
992:   if (!bnk->Diag_max) PetscCall(VecDuplicate(tao->solution, &bnk->Diag_max));
993:   PetscCall(TaoSetSolution(bnk->bncg, tao->solution));
994:   if (bnk->max_cg_its > 0) {
995:     /* Ensure that the important common vectors are shared between BNK and embedded BNCG */
996:     bnk->bncg_ctx = (TAO_BNCG *)bnk->bncg->data;
997:     PetscCall(PetscObjectReference((PetscObject)bnk->unprojected_gradient_old));
998:     PetscCall(VecDestroy(&bnk->bncg_ctx->unprojected_gradient_old));
999:     bnk->bncg_ctx->unprojected_gradient_old = bnk->unprojected_gradient_old;
1000:     PetscCall(PetscObjectReference((PetscObject)bnk->unprojected_gradient));
1001:     PetscCall(VecDestroy(&bnk->bncg_ctx->unprojected_gradient));
1002:     bnk->bncg_ctx->unprojected_gradient = bnk->unprojected_gradient;
1003:     PetscCall(PetscObjectReference((PetscObject)bnk->Gold));
1004:     PetscCall(VecDestroy(&bnk->bncg_ctx->G_old));
1005:     bnk->bncg_ctx->G_old = bnk->Gold;
1006:     PetscCall(PetscObjectReference((PetscObject)tao->gradient));
1007:     PetscCall(VecDestroy(&bnk->bncg->gradient));
1008:     bnk->bncg->gradient = tao->gradient;
1009:     PetscCall(PetscObjectReference((PetscObject)tao->stepdirection));
1010:     PetscCall(VecDestroy(&bnk->bncg->stepdirection));
1011:     bnk->bncg->stepdirection = tao->stepdirection;
1012:     /* Copy over some settings from BNK into BNCG */
1013:     PetscCall(TaoSetMaximumIterations(bnk->bncg, bnk->max_cg_its));
1014:     PetscCall(TaoSetTolerances(bnk->bncg, tao->gatol, tao->grtol, tao->gttol));
1015:     PetscCall(TaoSetFunctionLowerBound(bnk->bncg, tao->fmin));
1016:     PetscCall(TaoSetConvergenceTest(bnk->bncg, tao->ops->convergencetest, tao->cnvP));
1017:     {
1018:       TaoTerm   term;
1019:       Vec       params;
1020:       PetscReal scale;
1021:       Mat       map;

1023:       // Note: tao->objective_term.term and bnk->bncg->objective_term.term will point to same address
1024:       PetscCall(TaoGetTerm(tao, &scale, &term, &params, &map));
1025:       PetscCall(TaoAddTerm(bnk->bncg, NULL, scale, term, params, map));
1026:     }
1027:     PetscCall(PetscObjectCopyFortranFunctionPointers((PetscObject)tao, (PetscObject)bnk->bncg));
1028:   }
1029:   bnk->X_inactive    = NULL;
1030:   bnk->G_inactive    = NULL;
1031:   bnk->inactive_work = NULL;
1032:   bnk->active_work   = NULL;
1033:   bnk->inactive_idx  = NULL;
1034:   bnk->active_idx    = NULL;
1035:   bnk->active_lower  = NULL;
1036:   bnk->active_upper  = NULL;
1037:   bnk->active_fixed  = NULL;
1038:   bnk->M             = NULL;
1039:   bnk->H_inactive    = NULL;
1040:   bnk->Hpre_inactive = NULL;
1041:   PetscFunctionReturn(PETSC_SUCCESS);
1042: }

1044: PetscErrorCode TaoDestroy_BNK(Tao tao)
1045: {
1046:   TAO_BNK *bnk = (TAO_BNK *)tao->data;

1048:   PetscFunctionBegin;
1049:   PetscCall(VecDestroy(&bnk->W));
1050:   PetscCall(VecDestroy(&bnk->Xold));
1051:   PetscCall(VecDestroy(&bnk->Gold));
1052:   PetscCall(VecDestroy(&bnk->Xwork));
1053:   PetscCall(VecDestroy(&bnk->Gwork));
1054:   PetscCall(VecDestroy(&bnk->unprojected_gradient));
1055:   PetscCall(VecDestroy(&bnk->unprojected_gradient_old));
1056:   PetscCall(VecDestroy(&bnk->Diag_min));
1057:   PetscCall(VecDestroy(&bnk->Diag_max));
1058:   PetscCall(ISDestroy(&bnk->active_lower));
1059:   PetscCall(ISDestroy(&bnk->active_upper));
1060:   PetscCall(ISDestroy(&bnk->active_fixed));
1061:   PetscCall(ISDestroy(&bnk->active_idx));
1062:   PetscCall(ISDestroy(&bnk->inactive_idx));
1063:   PetscCall(MatDestroy(&bnk->Hpre_inactive));
1064:   PetscCall(MatDestroy(&bnk->H_inactive));
1065:   PetscCall(TaoDestroy(&bnk->bncg));
1066:   PetscCall(KSPDestroy(&tao->ksp));
1067:   PetscCall(PetscFree(tao->data));
1068:   PetscFunctionReturn(PETSC_SUCCESS);
1069: }

1071: PetscErrorCode TaoSetFromOptions_BNK(Tao tao, PetscOptionItems PetscOptionsObject)
1072: {
1073:   TAO_BNK *bnk = (TAO_BNK *)tao->data;

1075:   PetscFunctionBegin;
1076:   PetscOptionsHeadBegin(PetscOptionsObject, "Newton-Krylov method for bound constrained optimization");
1077:   PetscCall(PetscOptionsEList("-tao_bnk_init_type", "radius initialization type", "", BNK_INIT, BNK_INIT_TYPES, BNK_INIT[bnk->init_type], &bnk->init_type, NULL));
1078:   PetscCall(PetscOptionsEList("-tao_bnk_update_type", "radius update type", "", BNK_UPDATE, BNK_UPDATE_TYPES, BNK_UPDATE[bnk->update_type], &bnk->update_type, NULL));
1079:   PetscCall(PetscOptionsEList("-tao_bnk_as_type", "active set estimation method", "", BNK_AS, BNK_AS_TYPES, BNK_AS[bnk->as_type], &bnk->as_type, NULL));
1080:   PetscCall(PetscOptionsReal("-tao_bnk_sval", "(developer) Hessian perturbation starting value", "", bnk->sval, &bnk->sval, NULL));
1081:   PetscCall(PetscOptionsReal("-tao_bnk_imin", "(developer) minimum initial Hessian perturbation", "", bnk->imin, &bnk->imin, NULL));
1082:   PetscCall(PetscOptionsReal("-tao_bnk_imax", "(developer) maximum initial Hessian perturbation", "", bnk->imax, &bnk->imax, NULL));
1083:   PetscCall(PetscOptionsReal("-tao_bnk_imfac", "(developer) initial merit factor for Hessian perturbation", "", bnk->imfac, &bnk->imfac, NULL));
1084:   PetscCall(PetscOptionsReal("-tao_bnk_pmin", "(developer) minimum Hessian perturbation", "", bnk->pmin, &bnk->pmin, NULL));
1085:   PetscCall(PetscOptionsReal("-tao_bnk_pmax", "(developer) maximum Hessian perturbation", "", bnk->pmax, &bnk->pmax, NULL));
1086:   PetscCall(PetscOptionsReal("-tao_bnk_pgfac", "(developer) Hessian perturbation growth factor", "", bnk->pgfac, &bnk->pgfac, NULL));
1087:   PetscCall(PetscOptionsReal("-tao_bnk_psfac", "(developer) Hessian perturbation shrink factor", "", bnk->psfac, &bnk->psfac, NULL));
1088:   PetscCall(PetscOptionsReal("-tao_bnk_pmgfac", "(developer) merit growth factor for Hessian perturbation", "", bnk->pmgfac, &bnk->pmgfac, NULL));
1089:   PetscCall(PetscOptionsReal("-tao_bnk_pmsfac", "(developer) merit shrink factor for Hessian perturbation", "", bnk->pmsfac, &bnk->pmsfac, NULL));
1090:   PetscCall(PetscOptionsReal("-tao_bnk_eta1", "(developer) threshold for rejecting step (-tao_bnk_update_type reduction)", "", bnk->eta1, &bnk->eta1, NULL));
1091:   PetscCall(PetscOptionsReal("-tao_bnk_eta2", "(developer) threshold for accepting marginal step (-tao_bnk_update_type reduction)", "", bnk->eta2, &bnk->eta2, NULL));
1092:   PetscCall(PetscOptionsReal("-tao_bnk_eta3", "(developer) threshold for accepting reasonable step (-tao_bnk_update_type reduction)", "", bnk->eta3, &bnk->eta3, NULL));
1093:   PetscCall(PetscOptionsReal("-tao_bnk_eta4", "(developer) threshold for accepting good step (-tao_bnk_update_type reduction)", "", bnk->eta4, &bnk->eta4, NULL));
1094:   PetscCall(PetscOptionsReal("-tao_bnk_alpha1", "(developer) radius reduction factor for rejected step (-tao_bnk_update_type reduction)", "", bnk->alpha1, &bnk->alpha1, NULL));
1095:   PetscCall(PetscOptionsReal("-tao_bnk_alpha2", "(developer) radius reduction factor for marginally accepted bad step (-tao_bnk_update_type reduction)", "", bnk->alpha2, &bnk->alpha2, NULL));
1096:   PetscCall(PetscOptionsReal("-tao_bnk_alpha3", "(developer) radius increase factor for reasonable accepted step (-tao_bnk_update_type reduction)", "", bnk->alpha3, &bnk->alpha3, NULL));
1097:   PetscCall(PetscOptionsReal("-tao_bnk_alpha4", "(developer) radius increase factor for good accepted step (-tao_bnk_update_type reduction)", "", bnk->alpha4, &bnk->alpha4, NULL));
1098:   PetscCall(PetscOptionsReal("-tao_bnk_alpha5", "(developer) radius increase factor for very good accepted step (-tao_bnk_update_type reduction)", "", bnk->alpha5, &bnk->alpha5, NULL));
1099:   PetscCall(PetscOptionsReal("-tao_bnk_nu1", "(developer) threshold for small line-search step length (-tao_bnk_update_type step)", "", bnk->nu1, &bnk->nu1, NULL));
1100:   PetscCall(PetscOptionsReal("-tao_bnk_nu2", "(developer) threshold for reasonable line-search step length (-tao_bnk_update_type step)", "", bnk->nu2, &bnk->nu2, NULL));
1101:   PetscCall(PetscOptionsReal("-tao_bnk_nu3", "(developer) threshold for large line-search step length (-tao_bnk_update_type step)", "", bnk->nu3, &bnk->nu3, NULL));
1102:   PetscCall(PetscOptionsReal("-tao_bnk_nu4", "(developer) threshold for very large line-search step length (-tao_bnk_update_type step)", "", bnk->nu4, &bnk->nu4, NULL));
1103:   PetscCall(PetscOptionsReal("-tao_bnk_omega1", "(developer) radius reduction factor for very small line-search step length (-tao_bnk_update_type step)", "", bnk->omega1, &bnk->omega1, NULL));
1104:   PetscCall(PetscOptionsReal("-tao_bnk_omega2", "(developer) radius reduction factor for small line-search step length (-tao_bnk_update_type step)", "", bnk->omega2, &bnk->omega2, NULL));
1105:   PetscCall(PetscOptionsReal("-tao_bnk_omega3", "(developer) radius factor for decent line-search step length (-tao_bnk_update_type step)", "", bnk->omega3, &bnk->omega3, NULL));
1106:   PetscCall(PetscOptionsReal("-tao_bnk_omega4", "(developer) radius increase factor for large line-search step length (-tao_bnk_update_type step)", "", bnk->omega4, &bnk->omega4, NULL));
1107:   PetscCall(PetscOptionsReal("-tao_bnk_omega5", "(developer) radius increase factor for very large line-search step length (-tao_bnk_update_type step)", "", bnk->omega5, &bnk->omega5, NULL));
1108:   PetscCall(PetscOptionsReal("-tao_bnk_mu1_i", "(developer) threshold for accepting very good step (-tao_bnk_init_type interpolation)", "", bnk->mu1_i, &bnk->mu1_i, NULL));
1109:   PetscCall(PetscOptionsReal("-tao_bnk_mu2_i", "(developer) threshold for accepting good step (-tao_bnk_init_type interpolation)", "", bnk->mu2_i, &bnk->mu2_i, NULL));
1110:   PetscCall(PetscOptionsReal("-tao_bnk_gamma1_i", "(developer) radius reduction factor for rejected very bad step (-tao_bnk_init_type interpolation)", "", bnk->gamma1_i, &bnk->gamma1_i, NULL));
1111:   PetscCall(PetscOptionsReal("-tao_bnk_gamma2_i", "(developer) radius reduction factor for rejected bad step (-tao_bnk_init_type interpolation)", "", bnk->gamma2_i, &bnk->gamma2_i, NULL));
1112:   PetscCall(PetscOptionsReal("-tao_bnk_gamma3_i", "(developer) radius increase factor for accepted good step (-tao_bnk_init_type interpolation)", "", bnk->gamma3_i, &bnk->gamma3_i, NULL));
1113:   PetscCall(PetscOptionsReal("-tao_bnk_gamma4_i", "(developer) radius increase factor for accepted very good step (-tao_bnk_init_type interpolation)", "", bnk->gamma4_i, &bnk->gamma4_i, NULL));
1114:   PetscCall(PetscOptionsReal("-tao_bnk_theta_i", "(developer) trust region interpolation factor (-tao_bnk_init_type interpolation)", "", bnk->theta_i, &bnk->theta_i, NULL));
1115:   PetscCall(PetscOptionsReal("-tao_bnk_mu1", "(developer) threshold for accepting very good step (-tao_bnk_update_type interpolation)", "", bnk->mu1, &bnk->mu1, NULL));
1116:   PetscCall(PetscOptionsReal("-tao_bnk_mu2", "(developer) threshold for accepting good step (-tao_bnk_update_type interpolation)", "", bnk->mu2, &bnk->mu2, NULL));
1117:   PetscCall(PetscOptionsReal("-tao_bnk_gamma1", "(developer) radius reduction factor for rejected very bad step (-tao_bnk_update_type interpolation)", "", bnk->gamma1, &bnk->gamma1, NULL));
1118:   PetscCall(PetscOptionsReal("-tao_bnk_gamma2", "(developer) radius reduction factor for rejected bad step (-tao_bnk_update_type interpolation)", "", bnk->gamma2, &bnk->gamma2, NULL));
1119:   PetscCall(PetscOptionsReal("-tao_bnk_gamma3", "(developer) radius increase factor for accepted good step (-tao_bnk_update_type interpolation)", "", bnk->gamma3, &bnk->gamma3, NULL));
1120:   PetscCall(PetscOptionsReal("-tao_bnk_gamma4", "(developer) radius increase factor for accepted very good step (-tao_bnk_update_type interpolation)", "", bnk->gamma4, &bnk->gamma4, NULL));
1121:   PetscCall(PetscOptionsReal("-tao_bnk_theta", "(developer) trust region interpolation factor (-tao_bnk_update_type interpolation)", "", bnk->theta, &bnk->theta, NULL));
1122:   PetscCall(PetscOptionsReal("-tao_bnk_min_radius", "(developer) lower bound on initial radius", "", bnk->min_radius, &bnk->min_radius, NULL));
1123:   PetscCall(PetscOptionsReal("-tao_bnk_max_radius", "(developer) upper bound on radius", "", bnk->max_radius, &bnk->max_radius, NULL));
1124:   PetscCall(PetscOptionsReal("-tao_bnk_epsilon", "(developer) tolerance used when computing actual and predicted reduction", "", bnk->epsilon, &bnk->epsilon, NULL));
1125:   PetscCall(PetscOptionsReal("-tao_bnk_as_tol", "(developer) initial tolerance used when estimating actively bounded variables", "", bnk->as_tol, &bnk->as_tol, NULL));
1126:   PetscCall(PetscOptionsReal("-tao_bnk_as_step", "(developer) step length used when estimating actively bounded variables", "", bnk->as_step, &bnk->as_step, NULL));
1127:   PetscCall(PetscOptionsInt("-tao_bnk_max_cg_its", "number of BNCG iterations to take for each Newton step", "", bnk->max_cg_its, &bnk->max_cg_its, NULL));
1128:   PetscOptionsHeadEnd();

1130:   PetscCall(TaoSetOptionsPrefix(bnk->bncg, ((PetscObject)tao)->prefix));
1131:   PetscCall(TaoAppendOptionsPrefix(bnk->bncg, "tao_bnk_cg_"));
1132:   PetscCall(TaoSetFromOptions(bnk->bncg));

1134:   PetscCall(KSPSetOptionsPrefix(tao->ksp, ((PetscObject)tao)->prefix));
1135:   PetscCall(KSPAppendOptionsPrefix(tao->ksp, "tao_bnk_"));
1136:   PetscCall(KSPSetFromOptions(tao->ksp));
1137:   PetscFunctionReturn(PETSC_SUCCESS);
1138: }

1140: PetscErrorCode TaoView_BNK(Tao tao, PetscViewer viewer)
1141: {
1142:   TAO_BNK  *bnk = (TAO_BNK *)tao->data;
1143:   PetscInt  nrejects;
1144:   PetscBool isascii;

1146:   PetscFunctionBegin;
1147:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
1148:   if (isascii) {
1149:     PetscCall(PetscViewerASCIIPushTab(viewer));
1150:     PetscCall(TaoView(bnk->bncg, viewer));
1151:     if (bnk->M) {
1152:       PetscCall(MatLMVMGetRejectCount(bnk->M, &nrejects));
1153:       PetscCall(PetscViewerASCIIPrintf(viewer, "Rejected BFGS updates: %" PetscInt_FMT "\n", nrejects));
1154:     }
1155:     PetscCall(PetscViewerASCIIPrintf(viewer, "CG steps: %" PetscInt_FMT "\n", bnk->tot_cg_its));
1156:     PetscCall(PetscViewerASCIIPrintf(viewer, "Newton steps: %" PetscInt_FMT "\n", bnk->newt));
1157:     if (bnk->M) PetscCall(PetscViewerASCIIPrintf(viewer, "BFGS steps: %" PetscInt_FMT "\n", bnk->bfgs));
1158:     PetscCall(PetscViewerASCIIPrintf(viewer, "Scaled gradient steps: %" PetscInt_FMT "\n", bnk->sgrad));
1159:     PetscCall(PetscViewerASCIIPrintf(viewer, "Gradient steps: %" PetscInt_FMT "\n", bnk->grad));
1160:     PetscCall(PetscViewerASCIIPrintf(viewer, "KSP termination reasons:\n"));
1161:     PetscCall(PetscViewerASCIIPrintf(viewer, "  atol: %" PetscInt_FMT "\n", bnk->ksp_atol));
1162:     PetscCall(PetscViewerASCIIPrintf(viewer, "  rtol: %" PetscInt_FMT "\n", bnk->ksp_rtol));
1163:     PetscCall(PetscViewerASCIIPrintf(viewer, "  ctol: %" PetscInt_FMT "\n", bnk->ksp_ctol));
1164:     PetscCall(PetscViewerASCIIPrintf(viewer, "  negc: %" PetscInt_FMT "\n", bnk->ksp_negc));
1165:     PetscCall(PetscViewerASCIIPrintf(viewer, "  dtol: %" PetscInt_FMT "\n", bnk->ksp_dtol));
1166:     PetscCall(PetscViewerASCIIPrintf(viewer, "  iter: %" PetscInt_FMT "\n", bnk->ksp_iter));
1167:     PetscCall(PetscViewerASCIIPrintf(viewer, "  othr: %" PetscInt_FMT "\n", bnk->ksp_othr));
1168:     PetscCall(PetscViewerASCIIPopTab(viewer));
1169:   }
1170:   PetscFunctionReturn(PETSC_SUCCESS);
1171: }

1173: /*MC
1174:   TAOBNK - Shared base-type for Bounded Newton-Krylov type algorithms.
1175:   At each iteration, the BNK methods solve the symmetric
1176:   system of equations to obtain the step direction $d_k$:
1177:              $ H_k d_k = -g_k $
1178:   for free variables only. The step can be globalized either through
1179:   trust-region methods, or a line search, or a heuristic mixture of both.

1181:     Options Database Keys:
1182: + -tao_bnk_max_cg_its  - maximum number of bounded conjugate-gradient iterations taken in each Newton loop
1183: . -tao_bnk_init_type   - trust radius initialization method ("constant", "direction", "interpolation")
1184: . -tao_bnk_update_type - trust radius update method ("step", "direction", "interpolation")
1185: . -tao_bnk_as_type     - active-set estimation method ("none", "bertsekas")
1186: . -tao_bnk_as_tol      - (developer) initial tolerance used in estimating bounded active variables (-as_type bertsekas)
1187: . -tao_bnk_as_step     - (developer) trial step length used in estimating bounded active variables (-as_type bertsekas)
1188: . -tao_bnk_sval        - (developer) Hessian perturbation starting value
1189: . -tao_bnk_imin        - (developer) minimum initial Hessian perturbation
1190: . -tao_bnk_imax        - (developer) maximum initial Hessian perturbation
1191: . -tao_bnk_pmin        - (developer) minimum Hessian perturbation
1192: . -tao_bnk_pmax        - (developer) aximum Hessian perturbation
1193: . -tao_bnk_pgfac       - (developer) Hessian perturbation growth factor
1194: . -tao_bnk_psfac       - (developer) Hessian perturbation shrink factor
1195: . -tao_bnk_imfac       - (developer) initial merit factor for Hessian perturbation
1196: . -tao_bnk_pmgfac      - (developer) merit growth factor for Hessian perturbation
1197: . -tao_bnk_pmsfac      - (developer) merit shrink factor for Hessian perturbation
1198: . -tao_bnk_eta1        - (developer) threshold for rejecting step (-update_type reduction)
1199: . -tao_bnk_eta2        - (developer) threshold for accepting marginal step (-update_type reduction)
1200: . -tao_bnk_eta3        - (developer) threshold for accepting reasonable step (-update_type reduction)
1201: . -tao_bnk_eta4        - (developer) threshold for accepting good step (-update_type reduction)
1202: . -tao_bnk_alpha1      - (developer) radius reduction factor for rejected step (-update_type reduction)
1203: . -tao_bnk_alpha2      - (developer) radius reduction factor for marginally accepted bad step (-update_type reduction)
1204: . -tao_bnk_alpha3      - (developer) radius increase factor for reasonable accepted step (-update_type reduction)
1205: . -tao_bnk_alpha4      - (developer) radius increase factor for good accepted step (-update_type reduction)
1206: . -tao_bnk_alpha5      - (developer) radius increase factor for very good accepted step (-update_type reduction)
1207: . -tao_bnk_epsilon     - (developer) tolerance for small pred/actual ratios that trigger automatic step acceptance (-update_type reduction)
1208: . -tao_bnk_mu1         - (developer) threshold for accepting very good step (-update_type interpolation)
1209: . -tao_bnk_mu2         - (developer) threshold for accepting good step (-update_type interpolation)
1210: . -tao_bnk_gamma1      - (developer) radius reduction factor for rejected very bad step (-update_type interpolation)
1211: . -tao_bnk_gamma2      - (developer) radius reduction factor for rejected bad step (-update_type interpolation)
1212: . -tao_bnk_gamma3      - (developer) radius increase factor for accepted good step (-update_type interpolation)
1213: . -tao_bnk_gamma4      - (developer) radius increase factor for accepted very good step (-update_type interpolation)
1214: . -tao_bnk_theta       - (developer) trust region interpolation factor (-update_type interpolation)
1215: . -tao_bnk_nu1         - (developer) threshold for small line-search step length (-update_type step)
1216: . -tao_bnk_nu2         - (developer) threshold for reasonable line-search step length (-update_type step)
1217: . -tao_bnk_nu3         - (developer) threshold for large line-search step length (-update_type step)
1218: . -tao_bnk_nu4         - (developer) threshold for very large line-search step length (-update_type step)
1219: . -tao_bnk_omega1      - (developer) radius reduction factor for very small line-search step length (-update_type step)
1220: . -tao_bnk_omega2      - (developer) radius reduction factor for small line-search step length (-update_type step)
1221: . -tao_bnk_omega3      - (developer) radius factor for decent line-search step length (-update_type step)
1222: . -tao_bnk_omega4      - (developer) radius increase factor for large line-search step length (-update_type step)
1223: . -tao_bnk_omega5      - (developer) radius increase factor for very large line-search step length (-update_type step)
1224: . -tao_bnk_mu1_i       - (developer) threshold for accepting very good step (-init_type interpolation)
1225: . -tao_bnk_mu2_i       - (developer) threshold for accepting good step (-init_type interpolation)
1226: . -tao_bnk_gamma1_i    - (developer) radius reduction factor for rejected very bad step (-init_type interpolation)
1227: . -tao_bnk_gamma2_i    - (developer) radius reduction factor for rejected bad step (-init_type interpolation)
1228: . -tao_bnk_gamma3_i    - (developer) radius increase factor for accepted good step (-init_type interpolation)
1229: . -tao_bnk_gamma4_i    - (developer) radius increase factor for accepted very good step (-init_type interpolation)
1230: - -tao_bnk_theta_i     - (developer) trust region interpolation factor (-init_type interpolation)

1232:   Level: beginner

1234:   The various algorithmic factors can only be supplied via the options database

1236: .seealso: `Tao`, `TAONLS`, `TAONTL`, `TAONM`, `TaoType`, `TaoCreate()`
1237: M*/

1239: PetscErrorCode TaoCreate_BNK(Tao tao)
1240: {
1241:   TAO_BNK *bnk;
1242:   PC       pc;

1244:   PetscFunctionBegin;
1245:   PetscCall(PetscNew(&bnk));

1247:   tao->ops->setup            = TaoSetUp_BNK;
1248:   tao->ops->view             = TaoView_BNK;
1249:   tao->ops->setfromoptions   = TaoSetFromOptions_BNK;
1250:   tao->ops->destroy          = TaoDestroy_BNK;
1251:   tao->uses_gradient         = PETSC_TRUE;
1252:   tao->uses_hessian_matrices = PETSC_TRUE;

1254:   /*  Override default settings (unless already changed) */
1255:   PetscCall(TaoParametersInitialize(tao));
1256:   PetscObjectParameterSetDefault(tao, max_it, 50);
1257:   PetscObjectParameterSetDefault(tao, trust0, 100.0);

1259:   tao->data = (void *)bnk;

1261:   /*  Hessian shifting parameters */
1262:   bnk->computehessian = TaoBNKComputeHessian;
1263:   bnk->computestep    = TaoBNKComputeStep;

1265:   bnk->sval  = 0.0;
1266:   bnk->imin  = 1.0e-4;
1267:   bnk->imax  = 1.0e+2;
1268:   bnk->imfac = 1.0e-1;

1270:   bnk->pmin   = 1.0e-12;
1271:   bnk->pmax   = 1.0e+2;
1272:   bnk->pgfac  = 1.0e+1;
1273:   bnk->psfac  = 4.0e-1;
1274:   bnk->pmgfac = 1.0e-1;
1275:   bnk->pmsfac = 1.0e-1;

1277:   /*  Default values for trust-region radius update based on steplength */
1278:   bnk->nu1 = 0.25;
1279:   bnk->nu2 = 0.50;
1280:   bnk->nu3 = 1.00;
1281:   bnk->nu4 = 1.25;

1283:   bnk->omega1 = 0.25;
1284:   bnk->omega2 = 0.50;
1285:   bnk->omega3 = 1.00;
1286:   bnk->omega4 = 2.00;
1287:   bnk->omega5 = 4.00;

1289:   /*  Default values for trust-region radius update based on reduction */
1290:   bnk->eta1 = 1.0e-4;
1291:   bnk->eta2 = 0.25;
1292:   bnk->eta3 = 0.50;
1293:   bnk->eta4 = 0.90;

1295:   bnk->alpha1 = 0.25;
1296:   bnk->alpha2 = 0.50;
1297:   bnk->alpha3 = 1.00;
1298:   bnk->alpha4 = 2.00;
1299:   bnk->alpha5 = 4.00;

1301:   /*  Default values for trust-region radius update based on interpolation */
1302:   bnk->mu1 = 0.10;
1303:   bnk->mu2 = 0.50;

1305:   bnk->gamma1 = 0.25;
1306:   bnk->gamma2 = 0.50;
1307:   bnk->gamma3 = 2.00;
1308:   bnk->gamma4 = 4.00;

1310:   bnk->theta = 0.05;

1312:   /*  Default values for trust region initialization based on interpolation */
1313:   bnk->mu1_i = 0.35;
1314:   bnk->mu2_i = 0.50;

1316:   bnk->gamma1_i = 0.0625;
1317:   bnk->gamma2_i = 0.5;
1318:   bnk->gamma3_i = 2.0;
1319:   bnk->gamma4_i = 5.0;

1321:   bnk->theta_i = 0.25;

1323:   /*  Remaining parameters */
1324:   bnk->max_cg_its = 0;
1325:   bnk->min_radius = 1.0e-10;
1326:   bnk->max_radius = 1.0e10;
1327:   bnk->epsilon    = PetscPowReal(PETSC_MACHINE_EPSILON, 2.0 / 3.0);
1328:   bnk->as_tol     = 1.0e-3;
1329:   bnk->as_step    = 1.0e-3;
1330:   bnk->dmin       = 1.0e-6;
1331:   bnk->dmax       = 1.0e6;

1333:   bnk->M           = NULL;
1334:   bnk->bfgs_pre    = NULL;
1335:   bnk->init_type   = BNK_INIT_INTERPOLATION;
1336:   bnk->update_type = BNK_UPDATE_REDUCTION;
1337:   bnk->as_type     = BNK_AS_BERTSEKAS;

1339:   /* Create the embedded BNCG solver */
1340:   PetscCall(TaoCreate(PetscObjectComm((PetscObject)tao), &bnk->bncg));
1341:   PetscCall(PetscObjectIncrementTabLevel((PetscObject)bnk->bncg, (PetscObject)tao, 1));
1342:   PetscCall(TaoSetType(bnk->bncg, TAOBNCG));

1344:   /* Create the line search */
1345:   PetscCall(TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch));
1346:   PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1));
1347:   PetscCall(TaoLineSearchSetType(tao->linesearch, TAOLINESEARCHMT));
1348:   PetscCall(TaoLineSearchUseTaoRoutines(tao->linesearch, tao));

1350:   /*  Set linear solver to default for symmetric matrices */
1351:   PetscCall(KSPCreate(((PetscObject)tao)->comm, &tao->ksp));
1352:   PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1));
1353:   PetscCall(KSPSetType(tao->ksp, KSPSTCG));
1354:   PetscCall(KSPGetPC(tao->ksp, &pc));
1355:   PetscCall(PCSetType(pc, PCLMVM));
1356:   PetscFunctionReturn(PETSC_SUCCESS);
1357: }