Actual source code: bnk.c

petsc-3.12.5 2020-03-29
Report Typos and Errors
  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: /*------------------------------------------------------------*/

 11: /* Routine for initializing the KSP solver, the BFGS preconditioner, and the initial trust radius estimation */

 13: PetscErrorCode TaoBNKInitialize(Tao tao, PetscInt initType, PetscBool *needH)
 14: {
 15:   PetscErrorCode               ierr;
 16:   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
 17:   PC                           pc;

 19:   PetscReal                    f_min, ftrial, prered, actred, kappa, sigma, resnorm;
 20:   PetscReal                    tau, tau_1, tau_2, tau_max, tau_min, max_radius;
 21:   PetscBool                    is_bfgs, is_jacobi, is_symmetric, sym_set;
 22:   PetscInt                     n, N, nDiff;
 23:   PetscInt                     i_max = 5;
 24:   PetscInt                     j_max = 1;
 25:   PetscInt                     i, j;

 28:   /* Project the current point onto the feasible set */
 29:   TaoComputeVariableBounds(tao);
 30:   TaoSetVariableBounds(bnk->bncg, tao->XL, tao->XU);
 31:   if (tao->bounded) {
 32:     TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU);
 33:   }

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

 38:   /* Check convergence criteria */
 39:   TaoComputeObjectiveAndGradient(tao, tao->solution, &bnk->f, bnk->unprojected_gradient);
 40:   TaoBNKEstimateActiveSet(tao, bnk->as_type);
 41:   VecCopy(bnk->unprojected_gradient, tao->gradient);
 42:   VecISSet(tao->gradient, bnk->active_idx, 0.0);
 43:   TaoGradientNorm(tao, tao->gradient, NORM_2, &bnk->gnorm);

 45:   /* Test the initial point for convergence */
 46:   VecFischer(tao->solution, bnk->unprojected_gradient, tao->XL, tao->XU, bnk->W);
 47:   VecNorm(bnk->W, NORM_2, &resnorm);
 48:   if (PetscIsInfOrNanReal(bnk->f) || PetscIsInfOrNanReal(resnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
 49:   TaoLogConvergenceHistory(tao,bnk->f,resnorm,0.0,tao->ksp_its);
 50:   TaoMonitor(tao,tao->niter,bnk->f,resnorm,0.0,1.0);
 51:   (*tao->ops->convergencetest)(tao,tao->cnvP);
 52:   if (tao->reason != TAO_CONTINUE_ITERATING) return(0);
 53: 
 54:   /* Reset KSP stopping reason counters */
 55:   bnk->ksp_atol = 0;
 56:   bnk->ksp_rtol = 0;
 57:   bnk->ksp_dtol = 0;
 58:   bnk->ksp_ctol = 0;
 59:   bnk->ksp_negc = 0;
 60:   bnk->ksp_iter = 0;
 61:   bnk->ksp_othr = 0;
 62: 
 63:   /* Reset accepted step type counters */
 64:   bnk->tot_cg_its = 0;
 65:   bnk->newt = 0;
 66:   bnk->bfgs = 0;
 67:   bnk->sgrad = 0;
 68:   bnk->grad = 0;
 69: 
 70:   /* Initialize the Hessian perturbation */
 71:   bnk->pert = bnk->sval;
 72: 
 73:   /* Reset initial steplength to zero (this helps BNCG reset its direction internally) */
 74:   VecSet(tao->stepdirection, 0.0);

 76:   /* Allocate the vectors needed for the BFGS approximation */
 77:   KSPGetPC(tao->ksp, &pc);
 78:   PetscObjectTypeCompare((PetscObject)pc, PCLMVM, &is_bfgs);
 79:   PetscObjectTypeCompare((PetscObject)pc, PCJACOBI, &is_jacobi);
 80:   if (is_bfgs) {
 81:     bnk->bfgs_pre = pc;
 82:     PCLMVMGetMatLMVM(bnk->bfgs_pre, &bnk->M);
 83:     VecGetLocalSize(tao->solution, &n);
 84:     VecGetSize(tao->solution, &N);
 85:     MatSetSizes(bnk->M, n, n, N, N);
 86:     MatLMVMAllocate(bnk->M, tao->solution, bnk->unprojected_gradient);
 87:     MatIsSymmetricKnown(bnk->M, &sym_set, &is_symmetric);
 88:     if (!sym_set || !is_symmetric) SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM matrix in the LMVM preconditioner must be symmetric.");
 89:   } else if (is_jacobi) {
 90:     PCJacobiSetUseAbs(pc,PETSC_TRUE);
 91:   }
 92: 
 93:   /* Prepare the min/max vectors for safeguarding diagonal scales */
 94:   VecSet(bnk->Diag_min, bnk->dmin);
 95:   VecSet(bnk->Diag_max, bnk->dmax);

 97:   /* Initialize trust-region radius.  The initialization is only performed
 98:      when we are using Nash, Steihaug-Toint or the Generalized Lanczos method. */
 99:   *needH = PETSC_TRUE;
100:   if (bnk->is_nash || bnk->is_stcg || bnk->is_gltr) {
101:     switch(initType) {
102:     case BNK_INIT_CONSTANT:
103:       /* Use the initial radius specified */
104:       tao->trust = tao->trust0;
105:       break;

107:     case BNK_INIT_INTERPOLATION:
108:       /* Use interpolation based on the initial Hessian */
109:       max_radius = 0.0;
110:       tao->trust = tao->trust0;
111:       for (j = 0; j < j_max; ++j) {
112:         f_min = bnk->f;
113:         sigma = 0.0;

115:         if (*needH) {
116:           /* Compute the Hessian at the new step, and extract the inactive subsystem */
117:           bnk->computehessian(tao);
118:           TaoBNKEstimateActiveSet(tao, BNK_AS_NONE);
119:           MatDestroy(&bnk->H_inactive);
120:           if (bnk->active_idx) {
121:             MatCreateSubMatrix(tao->hessian, bnk->inactive_idx, bnk->inactive_idx, MAT_INITIAL_MATRIX, &bnk->H_inactive);
122:           } else {
123:             MatDuplicate(tao->hessian, MAT_COPY_VALUES, &bnk->H_inactive);
124:           }
125:           *needH = PETSC_FALSE;
126:         }

128:         for (i = 0; i < i_max; ++i) {
129:           /* Take a steepest descent step and snap it to bounds */
130:           VecCopy(tao->solution, bnk->Xold);
131:           VecAXPY(tao->solution, -tao->trust/bnk->gnorm, tao->gradient);
132:           TaoBoundSolution(tao->solution, tao->XL,tao->XU, 0.0, &nDiff, tao->solution);
133:           /* Compute the step we actually accepted */
134:           VecCopy(tao->solution, bnk->W);
135:           VecAXPY(bnk->W, -1.0, bnk->Xold);
136:           /* Compute the objective at the trial */
137:           TaoComputeObjective(tao, tao->solution, &ftrial);
138:           if (PetscIsInfOrNanReal(bnk->f)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
139:           VecCopy(bnk->Xold, tao->solution);
140:           if (PetscIsInfOrNanReal(ftrial)) {
141:             tau = bnk->gamma1_i;
142:           } else {
143:             if (ftrial < f_min) {
144:               f_min = ftrial;
145:               sigma = -tao->trust / bnk->gnorm;
146:             }
147: 
148:             /* Compute the predicted and actual reduction */
149:             if (bnk->active_idx) {
150:               VecGetSubVector(bnk->W, bnk->inactive_idx, &bnk->X_inactive);
151:               VecGetSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work);
152:             } else {
153:               bnk->X_inactive = bnk->W;
154:               bnk->inactive_work = bnk->Xwork;
155:             }
156:             MatMult(bnk->H_inactive, bnk->X_inactive, bnk->inactive_work);
157:             VecDot(bnk->X_inactive, bnk->inactive_work, &prered);
158:             if (bnk->active_idx) {
159:               VecRestoreSubVector(bnk->W, bnk->inactive_idx, &bnk->X_inactive);
160:               VecRestoreSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work);
161:             }
162:             prered = tao->trust * (bnk->gnorm - 0.5 * tao->trust * prered / (bnk->gnorm * bnk->gnorm));
163:             actred = bnk->f - ftrial;
164:             if ((PetscAbsScalar(actred) <= bnk->epsilon) && (PetscAbsScalar(prered) <= bnk->epsilon)) {
165:               kappa = 1.0;
166:             } else {
167:               kappa = actred / prered;
168:             }

170:             tau_1 = bnk->theta_i * bnk->gnorm * tao->trust / (bnk->theta_i * bnk->gnorm * tao->trust + (1.0 - bnk->theta_i) * prered - actred);
171:             tau_2 = bnk->theta_i * bnk->gnorm * tao->trust / (bnk->theta_i * bnk->gnorm * tao->trust - (1.0 + bnk->theta_i) * prered + actred);
172:             tau_min = PetscMin(tau_1, tau_2);
173:             tau_max = PetscMax(tau_1, tau_2);

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

179:               if (tau_max < 1.0) {
180:                 tau = bnk->gamma3_i;
181:               } else if (tau_max > bnk->gamma4_i) {
182:                 tau = bnk->gamma4_i;
183:               } else {
184:                 tau = tau_max;
185:               }
186:             } else if (PetscAbsScalar(kappa - (PetscReal)1.0) <= bnk->mu2_i) {
187:               /*  Good agreement */
188:               max_radius = PetscMax(max_radius, tao->trust);

190:               if (tau_max < bnk->gamma2_i) {
191:                 tau = bnk->gamma2_i;
192:               } else if (tau_max > bnk->gamma3_i) {
193:                 tau = bnk->gamma3_i;
194:               } else {
195:                 tau = tau_max;
196:               }
197:             } else {
198:               /*  Not good agreement */
199:               if (tau_min > 1.0) {
200:                 tau = bnk->gamma2_i;
201:               } else if (tau_max < bnk->gamma1_i) {
202:                 tau = bnk->gamma1_i;
203:               } else if ((tau_min < bnk->gamma1_i) && (tau_max >= 1.0)) {
204:                 tau = bnk->gamma1_i;
205:               } else if ((tau_1 >= bnk->gamma1_i) && (tau_1 < 1.0) && ((tau_2 < bnk->gamma1_i) || (tau_2 >= 1.0))) {
206:                 tau = tau_1;
207:               } else if ((tau_2 >= bnk->gamma1_i) && (tau_2 < 1.0) && ((tau_1 < bnk->gamma1_i) || (tau_2 >= 1.0))) {
208:                 tau = tau_2;
209:               } else {
210:                 tau = tau_max;
211:               }
212:             }
213:           }
214:           tao->trust = tau * tao->trust;
215:         }

217:         if (f_min < bnk->f) {
218:           /* We accidentally found a solution better than the initial, so accept it */
219:           bnk->f = f_min;
220:           VecCopy(tao->solution, bnk->Xold);
221:           VecAXPY(tao->solution,sigma,tao->gradient);
222:           TaoBoundSolution(tao->solution, tao->XL,tao->XU, 0.0, &nDiff, tao->solution);
223:           VecCopy(tao->solution, tao->stepdirection);
224:           VecAXPY(tao->stepdirection, -1.0, bnk->Xold);
225:           TaoComputeGradient(tao,tao->solution,bnk->unprojected_gradient);
226:           TaoBNKEstimateActiveSet(tao, bnk->as_type);
227:           VecCopy(bnk->unprojected_gradient, tao->gradient);
228:           VecISSet(tao->gradient, bnk->active_idx, 0.0);
229:           /* Compute gradient at the new iterate and flip switch to compute the Hessian later */
230:           TaoGradientNorm(tao, tao->gradient, NORM_2, &bnk->gnorm);
231:           *needH = PETSC_TRUE;
232:           /* Test the new step for convergence */
233:           VecFischer(tao->solution, bnk->unprojected_gradient, tao->XL, tao->XU, bnk->W);
234:           VecNorm(bnk->W, NORM_2, &resnorm);
235:           if (PetscIsInfOrNanReal(resnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
236:           TaoLogConvergenceHistory(tao,bnk->f,resnorm,0.0,tao->ksp_its);
237:           TaoMonitor(tao,tao->niter,bnk->f,resnorm,0.0,1.0);
238:           (*tao->ops->convergencetest)(tao,tao->cnvP);
239:           if (tao->reason != TAO_CONTINUE_ITERATING) return(0);
240:           /* active BNCG recycling early because we have a stepdirection computed */
241:           TaoBNCGSetRecycleFlag(bnk->bncg, PETSC_TRUE);
242:         }
243:       }
244:       tao->trust = PetscMax(tao->trust, max_radius);
245: 
246:       /* Ensure that the trust radius is within the limits */
247:       tao->trust = PetscMax(tao->trust, bnk->min_radius);
248:       tao->trust = PetscMin(tao->trust, bnk->max_radius);
249:       break;

251:     default:
252:       /* Norm of the first direction will initialize radius */
253:       tao->trust = 0.0;
254:       break;
255:     }
256:   }
257:   return(0);
258: }

260: /*------------------------------------------------------------*/

262: /* Routine for computing the exact Hessian and preparing the preconditioner at the new iterate */

264: PetscErrorCode TaoBNKComputeHessian(Tao tao)
265: {
266:   PetscErrorCode               ierr;
267:   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;

270:   /* Compute the Hessian */
271:   TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);
272:   /* Add a correction to the BFGS preconditioner */
273:   if (bnk->M) {
274:     MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);
275:   }
276:   /* Prepare the reduced sub-matrices for the inactive set */
277:   if (bnk->Hpre_inactive) {
278:     MatDestroy(&bnk->Hpre_inactive);
279:   }
280:   if (bnk->H_inactive) {
281:     MatDestroy(&bnk->H_inactive);
282:   }
283:   if (bnk->active_idx) {
284:     MatCreateSubMatrix(tao->hessian, bnk->inactive_idx, bnk->inactive_idx, MAT_INITIAL_MATRIX, &bnk->H_inactive);
285:     if (tao->hessian == tao->hessian_pre) {
286:       PetscObjectReference((PetscObject)bnk->H_inactive);
287:       bnk->Hpre_inactive = bnk->H_inactive;
288:     } else {
289:       MatCreateSubMatrix(tao->hessian_pre, bnk->inactive_idx, bnk->inactive_idx, MAT_INITIAL_MATRIX, &bnk->Hpre_inactive);
290:     }
291:     if (bnk->bfgs_pre) {
292:       PCLMVMSetIS(bnk->bfgs_pre, bnk->inactive_idx);
293:     }
294:   } else {
295:     MatDuplicate(tao->hessian, MAT_COPY_VALUES, &bnk->H_inactive);
296:     if (tao->hessian == tao->hessian_pre) {
297:       PetscObjectReference((PetscObject)bnk->H_inactive);
298:       bnk->Hpre_inactive = bnk->H_inactive;
299:     } else {
300:       MatDuplicate(tao->hessian_pre, MAT_COPY_VALUES, &bnk->Hpre_inactive);
301:     }
302:     if (bnk->bfgs_pre) {
303:       PCLMVMClearIS(bnk->bfgs_pre);
304:     }
305:   }
306:   return(0);
307: }

309: /*------------------------------------------------------------*/

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

313: PetscErrorCode TaoBNKEstimateActiveSet(Tao tao, PetscInt asType)
314: {
315:   PetscErrorCode               ierr;
316:   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
317:   PetscBool                    hessComputed, diagExists;

320:   switch (asType) {
321:   case BNK_AS_NONE:
322:     ISDestroy(&bnk->inactive_idx);
323:     VecWhichInactive(tao->XL, tao->solution, bnk->unprojected_gradient, tao->XU, PETSC_TRUE, &bnk->inactive_idx);
324:     ISDestroy(&bnk->active_idx);
325:     ISComplementVec(bnk->inactive_idx, tao->solution, &bnk->active_idx);
326:     break;

328:   case BNK_AS_BERTSEKAS:
329:     /* Compute the trial step vector with which we will estimate the active set at the next iteration */
330:     if (bnk->M) {
331:       /* If the BFGS preconditioner matrix is available, we will construct a trial step with it */
332:       MatSolve(bnk->M, bnk->unprojected_gradient, bnk->W);
333:     } else {
334:       hessComputed = diagExists = PETSC_FALSE;
335:       if (tao->hessian) {
336:         MatAssembled(tao->hessian, &hessComputed);
337:       }
338:       if (hessComputed) {
339:         MatHasOperation(tao->hessian, MATOP_GET_DIAGONAL, &diagExists);
340:       }
341:       if (diagExists) {
342:         /* BFGS preconditioner doesn't exist so let's invert the absolute diagonal of the Hessian instead onto the gradient */
343:         MatGetDiagonal(tao->hessian, bnk->Xwork);
344:         VecAbs(bnk->Xwork);
345:         VecMedian(bnk->Diag_min, bnk->Xwork, bnk->Diag_max, bnk->Xwork);
346:         VecReciprocal(bnk->Xwork);
347:         VecPointwiseMult(bnk->W, bnk->Xwork, bnk->unprojected_gradient);
348:       } else {
349:         /* If the Hessian or its diagonal does not exist, we will simply use gradient step */
350:         VecCopy(bnk->unprojected_gradient, bnk->W);
351:       }
352:     }
353:     VecScale(bnk->W, -1.0);
354:     TaoEstimateActiveBounds(tao->solution, tao->XL, tao->XU, bnk->unprojected_gradient, bnk->W, bnk->Xwork, bnk->as_step, &bnk->as_tol,
355:                                    &bnk->active_lower, &bnk->active_upper, &bnk->active_fixed, &bnk->active_idx, &bnk->inactive_idx);
356:     break;
357: 
358:   default:
359:     break;
360:   }
361:   return(0);
362: }

364: /*------------------------------------------------------------*/

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

368: PetscErrorCode TaoBNKBoundStep(Tao tao, PetscInt asType, Vec step)
369: {
370:   PetscErrorCode               ierr;
371:   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
372: 
374:   switch (asType) {
375:   case BNK_AS_NONE:
376:     VecISSet(step, bnk->active_idx, 0.0);
377:     break;

379:   case BNK_AS_BERTSEKAS:
380:     TaoBoundStep(tao->solution, tao->XL, tao->XU, bnk->active_lower, bnk->active_upper, bnk->active_fixed, 1.0, step);
381:     break;

383:   default:
384:     break;
385:   }
386:   return(0);
387: }

389: /*------------------------------------------------------------*/

391: /* Routine for taking a finite number of BNCG iterations to 
392:    accelerate Newton convergence.
393:    
394:    In practice, this approach simply trades off Hessian evaluations 
395:    for more gradient evaluations.
396: */

398: PetscErrorCode TaoBNKTakeCGSteps(Tao tao, PetscBool *terminate)
399: {
400:   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
401:   PetscErrorCode               ierr;
402: 
404:   *terminate = PETSC_FALSE;
405:   if (bnk->max_cg_its > 0) {
406:     /* Copy the current function value (important vectors are already shared) */
407:     bnk->bncg_ctx->f = bnk->f;
408:     /* Take some small finite number of BNCG iterations */
409:     TaoSolve(bnk->bncg);
410:     /* Add the number of gradient and function evaluations to the total */
411:     tao->nfuncs += bnk->bncg->nfuncs;
412:     tao->nfuncgrads += bnk->bncg->nfuncgrads;
413:     tao->ngrads += bnk->bncg->ngrads;
414:     tao->nhess += bnk->bncg->nhess;
415:     bnk->tot_cg_its += bnk->bncg->niter;
416:     /* Extract the BNCG function value out and save it into BNK */
417:     bnk->f = bnk->bncg_ctx->f;
418:     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) {
419:       *terminate = PETSC_TRUE;
420:     } else {
421:       TaoBNKEstimateActiveSet(tao, bnk->as_type);
422:     }
423:   }
424:   return(0);
425: }

427: /*------------------------------------------------------------*/

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

431: PetscErrorCode TaoBNKComputeStep(Tao tao, PetscBool shift, KSPConvergedReason *ksp_reason, PetscInt *step_type)
432: {
433:   PetscErrorCode               ierr;
434:   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
435:   PetscInt                     bfgsUpdates = 0;
436:   PetscInt                     kspits;
437:   PetscBool                    is_lmvm;
438: 
440:   /* If there are no inactive variables left, save some computation and return an adjusted zero step
441:      that has (l-x) and (u-x) for lower and upper bounded variables. */
442:   if (!bnk->inactive_idx) {
443:     VecSet(tao->stepdirection, 0.0);
444:     TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection);
445:     return(0);
446:   }
447: 
448:   /* Shift the reduced Hessian matrix */
449:   if ((shift) && (bnk->pert > 0)) {
450:     PetscObjectTypeCompare((PetscObject)tao->hessian, MATLMVM, &is_lmvm);
451:     if (is_lmvm) {
452:       MatShift(tao->hessian, bnk->pert);
453:     } else {
454:       MatShift(bnk->H_inactive, bnk->pert);
455:       if (bnk->H_inactive != bnk->Hpre_inactive) {
456:         MatShift(bnk->Hpre_inactive, bnk->pert);
457:       }
458:     }
459:   }
460: 
461:   /* Solve the Newton system of equations */
462:   tao->ksp_its = 0;
463:   VecSet(tao->stepdirection, 0.0);
464:   KSPReset(tao->ksp);
465:   KSPSetOperators(tao->ksp,bnk->H_inactive,bnk->Hpre_inactive);
466:   VecCopy(bnk->unprojected_gradient, bnk->Gwork);
467:   if (bnk->active_idx) {
468:     VecGetSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive);
469:     VecGetSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive);
470:   } else {
471:     bnk->G_inactive = bnk->unprojected_gradient;
472:     bnk->X_inactive = tao->stepdirection;
473:   }
474:   if (bnk->is_nash || bnk->is_stcg || bnk->is_gltr) {
475:     KSPCGSetRadius(tao->ksp,tao->trust);
476:     KSPSolve(tao->ksp, bnk->G_inactive, bnk->X_inactive);
477:     KSPGetIterationNumber(tao->ksp,&kspits);
478:     tao->ksp_its+=kspits;
479:     tao->ksp_tot_its+=kspits;
480:     KSPCGGetNormD(tao->ksp,&bnk->dnorm);

482:     if (0.0 == tao->trust) {
483:       /* Radius was uninitialized; use the norm of the direction */
484:       if (bnk->dnorm > 0.0) {
485:         tao->trust = bnk->dnorm;

487:         /* Modify the radius if it is too large or small */
488:         tao->trust = PetscMax(tao->trust, bnk->min_radius);
489:         tao->trust = PetscMin(tao->trust, bnk->max_radius);
490:       } else {
491:         /* The direction was bad; set radius to default value and re-solve
492:            the trust-region subproblem to get a direction */
493:         tao->trust = tao->trust0;

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

499:         KSPCGSetRadius(tao->ksp,tao->trust);
500:         KSPSolve(tao->ksp, bnk->G_inactive, bnk->X_inactive);
501:         KSPGetIterationNumber(tao->ksp,&kspits);
502:         tao->ksp_its+=kspits;
503:         tao->ksp_tot_its+=kspits;
504:         KSPCGGetNormD(tao->ksp,&bnk->dnorm);

506:         if (bnk->dnorm == 0.0) SETERRQ(PETSC_COMM_SELF,1, "Initial direction zero");
507:       }
508:     }
509:   } else {
510:     KSPSolve(tao->ksp, bnk->G_inactive, bnk->X_inactive);
511:     KSPGetIterationNumber(tao->ksp, &kspits);
512:     tao->ksp_its += kspits;
513:     tao->ksp_tot_its+=kspits;
514:   }
515:   /* Restore sub vectors back */
516:   if (bnk->active_idx) {
517:     VecRestoreSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive);
518:     VecRestoreSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive);
519:   }
520:   /* Make sure the safeguarded fall-back step is zero for actively bounded variables */
521:   VecScale(tao->stepdirection, -1.0);
522:   TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection);
523: 
524:   /* Record convergence reasons */
525:   KSPGetConvergedReason(tao->ksp, ksp_reason);
526:   if (KSP_CONVERGED_ATOL == *ksp_reason) {
527:     ++bnk->ksp_atol;
528:   } else if (KSP_CONVERGED_RTOL == *ksp_reason) {
529:     ++bnk->ksp_rtol;
530:   } else if (KSP_CONVERGED_CG_CONSTRAINED == *ksp_reason) {
531:     ++bnk->ksp_ctol;
532:   } else if (KSP_CONVERGED_CG_NEG_CURVE == *ksp_reason) {
533:     ++bnk->ksp_negc;
534:   } else if (KSP_DIVERGED_DTOL == *ksp_reason) {
535:     ++bnk->ksp_dtol;
536:   } else if (KSP_DIVERGED_ITS == *ksp_reason) {
537:     ++bnk->ksp_iter;
538:   } else {
539:     ++bnk->ksp_othr;
540:   }
541: 
542:   /* Make sure the BFGS preconditioner is healthy */
543:   if (bnk->M) {
544:     MatLMVMGetUpdateCount(bnk->M, &bfgsUpdates);
545:     if ((KSP_DIVERGED_INDEFINITE_PC == *ksp_reason) && (bfgsUpdates > 0)) {
546:       /* Preconditioner is numerically indefinite; reset the approximation. */
547:       MatLMVMReset(bnk->M, PETSC_FALSE);
548:       MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);
549:     }
550:   }
551:   *step_type = BNK_NEWTON;
552:   return(0);
553: }

555: /*------------------------------------------------------------*/

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

559: PetscErrorCode TaoBNKRecomputePred(Tao tao, Vec S, PetscReal *prered)
560: {
561:   PetscErrorCode               ierr;
562:   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
563: 
565:   /* Extract subvectors associated with the inactive set */
566:   if (bnk->active_idx){
567:     VecGetSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive);
568:     VecGetSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work);
569:     VecGetSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive);
570:   } else {
571:     bnk->X_inactive = tao->stepdirection;
572:     bnk->inactive_work = bnk->Xwork;
573:     bnk->G_inactive = bnk->Gwork;
574:   }
575:   /* Recompute the predicted decrease based on the quadratic model */
576:   MatMult(bnk->H_inactive, bnk->X_inactive, bnk->inactive_work);
577:   VecAYPX(bnk->inactive_work, -0.5, bnk->G_inactive);
578:   VecDot(bnk->inactive_work, bnk->X_inactive, prered);
579:   /* Restore the sub vectors */
580:   if (bnk->active_idx){
581:     VecRestoreSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive);
582:     VecRestoreSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work);
583:     VecRestoreSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive);
584:   }
585:   return(0);
586: }

588: /*------------------------------------------------------------*/

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

592:    The step direction falls back onto BFGS, scaled gradient and gradient steps 
593:    in the event that the Newton step fails the test.
594: */

596: PetscErrorCode TaoBNKSafeguardStep(Tao tao, KSPConvergedReason ksp_reason, PetscInt *stepType)
597: {
598:   PetscErrorCode               ierr;
599:   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
600: 
601:   PetscReal                    gdx, e_min;
602:   PetscInt                     bfgsUpdates;
603: 
605:   switch (*stepType) {
606:   case BNK_NEWTON:
607:     VecDot(tao->stepdirection, tao->gradient, &gdx);
608:     if ((gdx >= 0.0) || PetscIsInfOrNanReal(gdx)) {
609:       /* Newton step is not descent or direction produced Inf or NaN
610:         Update the perturbation for next time */
611:       if (bnk->pert <= 0.0) {
612:         /* Initialize the perturbation */
613:         bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
614:         if (bnk->is_gltr) {
615:           KSPGLTRGetMinEig(tao->ksp,&e_min);
616:           bnk->pert = PetscMax(bnk->pert, -e_min);
617:         }
618:       } else {
619:         /* Increase the perturbation */
620:         bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
621:       }

623:       if (!bnk->M) {
624:         /* We don't have the bfgs matrix around and updated
625:           Must use gradient direction in this case */
626:         VecCopy(tao->gradient, tao->stepdirection);
627:         *stepType = BNK_GRADIENT;
628:       } else {
629:         /* Attempt to use the BFGS direction */
630:         MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);

632:         /* Check for success (descent direction) 
633:           NOTE: Negative gdx here means not a descent direction because 
634:           the fall-back step is missing a negative sign. */
635:         VecDot(tao->gradient, tao->stepdirection, &gdx);
636:         if ((gdx <= 0.0) || PetscIsInfOrNanReal(gdx)) {
637:           /* BFGS direction is not descent or direction produced not a number
638:             We can assert bfgsUpdates > 1 in this case because
639:             the first solve produces the scaled gradient direction,
640:             which is guaranteed to be descent */

642:           /* Use steepest descent direction (scaled) */
643:           MatLMVMReset(bnk->M, PETSC_FALSE);
644:           MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);
645:           MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);

647:           *stepType = BNK_SCALED_GRADIENT;
648:         } else {
649:           MatLMVMGetUpdateCount(bnk->M, &bfgsUpdates);
650:           if (1 == bfgsUpdates) {
651:             /* The first BFGS direction is always the scaled gradient */
652:             *stepType = BNK_SCALED_GRADIENT;
653:           } else {
654:             *stepType = BNK_BFGS;
655:           }
656:         }
657:       }
658:       /* Make sure the safeguarded fall-back step is zero for actively bounded variables */
659:       VecScale(tao->stepdirection, -1.0);
660:       TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection);
661:     } else {
662:       /* Computed Newton step is descent */
663:       switch (ksp_reason) {
664:       case KSP_DIVERGED_NANORINF:
665:       case KSP_DIVERGED_BREAKDOWN:
666:       case KSP_DIVERGED_INDEFINITE_MAT:
667:       case KSP_DIVERGED_INDEFINITE_PC:
668:       case KSP_CONVERGED_CG_NEG_CURVE:
669:         /* Matrix or preconditioner is indefinite; increase perturbation */
670:         if (bnk->pert <= 0.0) {
671:           /* Initialize the perturbation */
672:           bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
673:           if (bnk->is_gltr) {
674:             KSPGLTRGetMinEig(tao->ksp, &e_min);
675:             bnk->pert = PetscMax(bnk->pert, -e_min);
676:           }
677:         } else {
678:           /* Increase the perturbation */
679:           bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
680:         }
681:         break;

683:       default:
684:         /* Newton step computation is good; decrease perturbation */
685:         bnk->pert = PetscMin(bnk->psfac * bnk->pert, bnk->pmsfac * bnk->gnorm);
686:         if (bnk->pert < bnk->pmin) {
687:           bnk->pert = 0.0;
688:         }
689:         break;
690:       }
691:       *stepType = BNK_NEWTON;
692:     }
693:     break;
694: 
695:   case BNK_BFGS:
696:     /* Check for success (descent direction) */
697:     VecDot(tao->stepdirection, tao->gradient, &gdx);
698:     if (gdx >= 0 || PetscIsInfOrNanReal(gdx)) {
699:       /* Step is not descent or solve was not successful
700:          Use steepest descent direction (scaled) */
701:       MatLMVMReset(bnk->M, PETSC_FALSE);
702:       MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);
703:       MatSolve(bnk->M, tao->gradient, tao->stepdirection);
704:       VecScale(tao->stepdirection,-1.0);
705:       TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection);
706:       *stepType = BNK_SCALED_GRADIENT;
707:     } else {
708:       *stepType = BNK_BFGS;
709:     }
710:     break;
711: 
712:   case BNK_SCALED_GRADIENT:
713:     break;
714: 
715:   default:
716:     break;
717:   }
718: 
719:   return(0);
720: }

722: /*------------------------------------------------------------*/

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

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

730: PetscErrorCode TaoBNKPerformLineSearch(Tao tao, PetscInt *stepType, PetscReal *steplen, TaoLineSearchConvergedReason *reason)
731: {
732:   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
734:   TaoLineSearchConvergedReason ls_reason;
735: 
736:   PetscReal      e_min, gdx;
737:   PetscInt       bfgsUpdates;
738: 
740:   /* Perform the linesearch */
741:   TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, steplen, &ls_reason);
742:   TaoAddLineSearchCounts(tao);

744:   while (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER && *stepType != BNK_SCALED_GRADIENT && *stepType != BNK_GRADIENT) {
745:     /* Linesearch failed, revert solution */
746:     bnk->f = bnk->fold;
747:     VecCopy(bnk->Xold, tao->solution);
748:     VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient);

750:     switch(*stepType) {
751:     case BNK_NEWTON:
752:       /* Failed to obtain acceptable iterate with Newton step
753:          Update the perturbation for next time */
754:       if (bnk->pert <= 0.0) {
755:         /* Initialize the perturbation */
756:         bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
757:         if (bnk->is_gltr) {
758:           KSPGLTRGetMinEig(tao->ksp,&e_min);
759:           bnk->pert = PetscMax(bnk->pert, -e_min);
760:         }
761:       } else {
762:         /* Increase the perturbation */
763:         bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
764:       }

766:       if (!bnk->M) {
767:         /* We don't have the bfgs matrix around and being updated
768:            Must use gradient direction in this case */
769:         VecCopy(bnk->unprojected_gradient, tao->stepdirection);
770:         *stepType = BNK_GRADIENT;
771:       } else {
772:         /* Attempt to use the BFGS direction */
773:         MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);
774:         /* Check for success (descent direction) 
775:            NOTE: Negative gdx means not a descent direction because the step here is missing a negative sign. */
776:         VecDot(tao->gradient, tao->stepdirection, &gdx);
777:         if ((gdx <= 0.0) || PetscIsInfOrNanReal(gdx)) {
778:           /* BFGS direction is not descent or direction produced not a number
779:              We can assert bfgsUpdates > 1 in this case
780:              Use steepest descent direction (scaled) */
781:           MatLMVMReset(bnk->M, PETSC_FALSE);
782:           MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);
783:           MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);

785:           bfgsUpdates = 1;
786:           *stepType = BNK_SCALED_GRADIENT;
787:         } else {
788:           MatLMVMGetUpdateCount(bnk->M, &bfgsUpdates);
789:           if (1 == bfgsUpdates) {
790:             /* The first BFGS direction is always the scaled gradient */
791:             *stepType = BNK_SCALED_GRADIENT;
792:           } else {
793:             *stepType = BNK_BFGS;
794:           }
795:         }
796:       }
797:       break;

799:     case BNK_BFGS:
800:       /* Can only enter if pc_type == BNK_PC_BFGS
801:          Failed to obtain acceptable iterate with BFGS step
802:          Attempt to use the scaled gradient direction */
803:       MatLMVMReset(bnk->M, PETSC_FALSE);
804:       MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);
805:       MatSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);

807:       bfgsUpdates = 1;
808:       *stepType = BNK_SCALED_GRADIENT;
809:       break;
810:     }
811:     /* Make sure the safeguarded fall-back step is zero for actively bounded variables */
812:     VecScale(tao->stepdirection, -1.0);
813:     TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection);
814: 
815:     /* Perform one last line search with the fall-back step */
816:     TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, steplen, &ls_reason);
817:     TaoAddLineSearchCounts(tao);
818:   }
819:   *reason = ls_reason;
820:   return(0);
821: }

823: /*------------------------------------------------------------*/

825: /* Routine for updating the trust radius. 

827:   Function features three different update methods: 
828:   1) Line-search step length based
829:   2) Predicted decrease on the CG quadratic model
830:   3) Interpolation
831: */

833: PetscErrorCode TaoBNKUpdateTrustRadius(Tao tao, PetscReal prered, PetscReal actred, PetscInt updateType, PetscInt stepType, PetscBool *accept)
834: {
835:   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
837: 
838:   PetscReal      step, kappa;
839:   PetscReal      gdx, tau_1, tau_2, tau_min, tau_max;

842:   /* Update trust region radius */
843:   *accept = PETSC_FALSE;
844:   switch(updateType) {
845:   case BNK_UPDATE_STEP:
846:     *accept = PETSC_TRUE; /* always accept here because line search succeeded */
847:     if (stepType == BNK_NEWTON) {
848:       TaoLineSearchGetStepLength(tao->linesearch, &step);
849:       if (step < bnk->nu1) {
850:         /* Very bad step taken; reduce radius */
851:         tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust);
852:       } else if (step < bnk->nu2) {
853:         /* Reasonably bad step taken; reduce radius */
854:         tao->trust = bnk->omega2 * PetscMin(bnk->dnorm, tao->trust);
855:       } else if (step < bnk->nu3) {
856:         /*  Reasonable step was taken; leave radius alone */
857:         if (bnk->omega3 < 1.0) {
858:           tao->trust = bnk->omega3 * PetscMin(bnk->dnorm, tao->trust);
859:         } else if (bnk->omega3 > 1.0) {
860:           tao->trust = PetscMax(bnk->omega3 * bnk->dnorm, tao->trust);
861:         }
862:       } else if (step < bnk->nu4) {
863:         /*  Full step taken; increase the radius */
864:         tao->trust = PetscMax(bnk->omega4 * bnk->dnorm, tao->trust);
865:       } else {
866:         /*  More than full step taken; increase the radius */
867:         tao->trust = PetscMax(bnk->omega5 * bnk->dnorm, tao->trust);
868:       }
869:     } else {
870:       /*  Newton step was not good; reduce the radius */
871:       tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust);
872:     }
873:     break;

875:   case BNK_UPDATE_REDUCTION:
876:     if (stepType == BNK_NEWTON) {
877:       if ((prered < 0.0) || PetscIsInfOrNanReal(prered)) {
878:         /* The predicted reduction has the wrong sign.  This cannot
879:            happen in infinite precision arithmetic.  Step should
880:            be rejected! */
881:         tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
882:       } else {
883:         if (PetscIsInfOrNanReal(actred)) {
884:           tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
885:         } else {
886:           if ((PetscAbsScalar(actred) <= PetscMax(1.0, PetscAbsScalar(bnk->f))*bnk->epsilon) && (PetscAbsScalar(prered) <= PetscMax(1.0, PetscAbsScalar(bnk->f))*bnk->epsilon)) {
887:             kappa = 1.0;
888:           } else {
889:             kappa = actred / prered;
890:           }
891:           /* Accept or reject the step and update radius */
892:           if (kappa < bnk->eta1) {
893:             /* Reject the step */
894:             tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
895:           } else {
896:             /* Accept the step */
897:             *accept = PETSC_TRUE;
898:             /* Update the trust region radius only if the computed step is at the trust radius boundary */
899:             if (bnk->dnorm == tao->trust) {
900:               if (kappa < bnk->eta2) {
901:                 /* Marginal bad step */
902:                 tao->trust = bnk->alpha2 * tao->trust;
903:               } else if (kappa < bnk->eta3) {
904:                 /* Reasonable step */
905:                 tao->trust = bnk->alpha3 * tao->trust;
906:               } else if (kappa < bnk->eta4) {
907:                 /* Good step */
908:                 tao->trust = bnk->alpha4 * tao->trust;
909:               } else {
910:                 /* Very good step */
911:                 tao->trust = bnk->alpha5 * tao->trust;
912:               }
913:             }
914:           }
915:         }
916:       }
917:     } else {
918:       /*  Newton step was not good; reduce the radius */
919:       tao->trust = bnk->alpha1 * PetscMin(bnk->dnorm, tao->trust);
920:     }
921:     break;

923:   default:
924:     if (stepType == BNK_NEWTON) {
925:       if (prered < 0.0) {
926:         /*  The predicted reduction has the wrong sign.  This cannot */
927:         /*  happen in infinite precision arithmetic.  Step should */
928:         /*  be rejected! */
929:         tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
930:       } else {
931:         if (PetscIsInfOrNanReal(actred)) {
932:           tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
933:         } else {
934:           if ((PetscAbsScalar(actred) <= bnk->epsilon) && (PetscAbsScalar(prered) <= bnk->epsilon)) {
935:             kappa = 1.0;
936:           } else {
937:             kappa = actred / prered;
938:           }
939: 
940:           VecDot(tao->gradient, tao->stepdirection, &gdx);
941:           tau_1 = bnk->theta * gdx / (bnk->theta * gdx - (1.0 - bnk->theta) * prered + actred);
942:           tau_2 = bnk->theta * gdx / (bnk->theta * gdx + (1.0 + bnk->theta) * prered - actred);
943:           tau_min = PetscMin(tau_1, tau_2);
944:           tau_max = PetscMax(tau_1, tau_2);

946:           if (kappa >= 1.0 - bnk->mu1) {
947:             /*  Great agreement */
948:             *accept = PETSC_TRUE;
949:             if (tau_max < 1.0) {
950:               tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm);
951:             } else if (tau_max > bnk->gamma4) {
952:               tao->trust = PetscMax(tao->trust, bnk->gamma4 * bnk->dnorm);
953:             } else {
954:               tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm);
955:             }
956:           } else if (kappa >= 1.0 - bnk->mu2) {
957:             /*  Good agreement */
958:             *accept = PETSC_TRUE;
959:             if (tau_max < bnk->gamma2) {
960:               tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm);
961:             } else if (tau_max > bnk->gamma3) {
962:               tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm);
963:             } else if (tau_max < 1.0) {
964:               tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm);
965:             } else {
966:               tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm);
967:             }
968:           } else {
969:             /*  Not good agreement */
970:             if (tau_min > 1.0) {
971:               tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm);
972:             } else if (tau_max < bnk->gamma1) {
973:               tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
974:             } else if ((tau_min < bnk->gamma1) && (tau_max >= 1.0)) {
975:               tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
976:             } else if ((tau_1 >= bnk->gamma1) && (tau_1 < 1.0) && ((tau_2 < bnk->gamma1) || (tau_2 >= 1.0))) {
977:               tao->trust = tau_1 * PetscMin(tao->trust, bnk->dnorm);
978:             } else if ((tau_2 >= bnk->gamma1) && (tau_2 < 1.0) && ((tau_1 < bnk->gamma1) || (tau_2 >= 1.0))) {
979:               tao->trust = tau_2 * PetscMin(tao->trust, bnk->dnorm);
980:             } else {
981:               tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm);
982:             }
983:           }
984:         }
985:       }
986:     } else {
987:       /*  Newton step was not good; reduce the radius */
988:       tao->trust = bnk->gamma1 * PetscMin(bnk->dnorm, tao->trust);
989:     }
990:     break;
991:   }
992:   /* Make sure the radius does not violate min and max settings */
993:   tao->trust = PetscMin(tao->trust, bnk->max_radius);
994:   tao->trust = PetscMax(tao->trust, bnk->min_radius);
995:   return(0);
996: }

998: /* ---------------------------------------------------------- */

1000: PetscErrorCode TaoBNKAddStepCounts(Tao tao, PetscInt stepType)
1001: {
1002:   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
1003: 
1005:   switch (stepType) {
1006:   case BNK_NEWTON:
1007:     ++bnk->newt;
1008:     break;
1009:   case BNK_BFGS:
1010:     ++bnk->bfgs;
1011:     break;
1012:   case BNK_SCALED_GRADIENT:
1013:     ++bnk->sgrad;
1014:     break;
1015:   case BNK_GRADIENT:
1016:     ++bnk->grad;
1017:     break;
1018:   default:
1019:     break;
1020:   }
1021:   return(0);
1022: }

1024: /* ---------------------------------------------------------- */

1026: PetscErrorCode TaoSetUp_BNK(Tao tao)
1027: {
1028:   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
1030:   PetscInt       i;

1033:   if (!tao->gradient) {
1034:     VecDuplicate(tao->solution,&tao->gradient);
1035:   }
1036:   if (!tao->stepdirection) {
1037:     VecDuplicate(tao->solution,&tao->stepdirection);
1038:   }
1039:   if (!bnk->W) {
1040:     VecDuplicate(tao->solution,&bnk->W);
1041:   }
1042:   if (!bnk->Xold) {
1043:     VecDuplicate(tao->solution,&bnk->Xold);
1044:   }
1045:   if (!bnk->Gold) {
1046:     VecDuplicate(tao->solution,&bnk->Gold);
1047:   }
1048:   if (!bnk->Xwork) {
1049:     VecDuplicate(tao->solution,&bnk->Xwork);
1050:   }
1051:   if (!bnk->Gwork) {
1052:     VecDuplicate(tao->solution,&bnk->Gwork);
1053:   }
1054:   if (!bnk->unprojected_gradient) {
1055:     VecDuplicate(tao->solution,&bnk->unprojected_gradient);
1056:   }
1057:   if (!bnk->unprojected_gradient_old) {
1058:     VecDuplicate(tao->solution,&bnk->unprojected_gradient_old);
1059:   }
1060:   if (!bnk->Diag_min) {
1061:     VecDuplicate(tao->solution,&bnk->Diag_min);
1062:   }
1063:   if (!bnk->Diag_max) {
1064:     VecDuplicate(tao->solution,&bnk->Diag_max);
1065:   }
1066:   if (bnk->max_cg_its > 0) {
1067:     /* Ensure that the important common vectors are shared between BNK and embedded BNCG */
1068:     bnk->bncg_ctx = (TAO_BNCG *)bnk->bncg->data;
1069:     PetscObjectReference((PetscObject)(bnk->unprojected_gradient_old));
1070:     VecDestroy(&bnk->bncg_ctx->unprojected_gradient_old);
1071:     bnk->bncg_ctx->unprojected_gradient_old = bnk->unprojected_gradient_old;
1072:     PetscObjectReference((PetscObject)(bnk->unprojected_gradient));
1073:     VecDestroy(&bnk->bncg_ctx->unprojected_gradient);
1074:     bnk->bncg_ctx->unprojected_gradient = bnk->unprojected_gradient;
1075:     PetscObjectReference((PetscObject)(bnk->Gold));
1076:     VecDestroy(&bnk->bncg_ctx->G_old);
1077:     bnk->bncg_ctx->G_old = bnk->Gold;
1078:     PetscObjectReference((PetscObject)(tao->gradient));
1079:     VecDestroy(&bnk->bncg->gradient);
1080:     bnk->bncg->gradient = tao->gradient;
1081:     PetscObjectReference((PetscObject)(tao->stepdirection));
1082:     VecDestroy(&bnk->bncg->stepdirection);
1083:     bnk->bncg->stepdirection = tao->stepdirection;
1084:     TaoSetInitialVector(bnk->bncg, tao->solution);
1085:     /* Copy over some settings from BNK into BNCG */
1086:     TaoSetMaximumIterations(bnk->bncg, bnk->max_cg_its);
1087:     TaoSetTolerances(bnk->bncg, tao->gatol, tao->grtol, tao->gttol);
1088:     TaoSetFunctionLowerBound(bnk->bncg, tao->fmin);
1089:     TaoSetConvergenceTest(bnk->bncg, tao->ops->convergencetest, tao->cnvP);
1090:     TaoSetObjectiveRoutine(bnk->bncg, tao->ops->computeobjective, tao->user_objP);
1091:     TaoSetGradientRoutine(bnk->bncg, tao->ops->computegradient, tao->user_gradP);
1092:     TaoSetObjectiveAndGradientRoutine(bnk->bncg, tao->ops->computeobjectiveandgradient, tao->user_objgradP);
1093:     PetscObjectCopyFortranFunctionPointers((PetscObject)tao, (PetscObject)(bnk->bncg));
1094:     for (i=0; i<tao->numbermonitors; ++i) {
1095:       TaoSetMonitor(bnk->bncg, tao->monitor[i], tao->monitorcontext[i], tao->monitordestroy[i]);
1096:       PetscObjectReference((PetscObject)(tao->monitorcontext[i]));
1097:     }
1098:   }
1099:   bnk->X_inactive = 0;
1100:   bnk->G_inactive = 0;
1101:   bnk->inactive_work = 0;
1102:   bnk->active_work = 0;
1103:   bnk->inactive_idx = 0;
1104:   bnk->active_idx = 0;
1105:   bnk->active_lower = 0;
1106:   bnk->active_upper = 0;
1107:   bnk->active_fixed = 0;
1108:   bnk->M = 0;
1109:   bnk->H_inactive = 0;
1110:   bnk->Hpre_inactive = 0;
1111:   return(0);
1112: }

1114: /*------------------------------------------------------------*/

1116: PetscErrorCode TaoDestroy_BNK(Tao tao)
1117: {
1118:   TAO_BNK        *bnk = (TAO_BNK *)tao->data;

1122:   if (tao->setupcalled) {
1123:     VecDestroy(&bnk->W);
1124:     VecDestroy(&bnk->Xold);
1125:     VecDestroy(&bnk->Gold);
1126:     VecDestroy(&bnk->Xwork);
1127:     VecDestroy(&bnk->Gwork);
1128:     VecDestroy(&bnk->unprojected_gradient);
1129:     VecDestroy(&bnk->unprojected_gradient_old);
1130:     VecDestroy(&bnk->Diag_min);
1131:     VecDestroy(&bnk->Diag_max);
1132:   }
1133:   ISDestroy(&bnk->active_lower);
1134:   ISDestroy(&bnk->active_upper);
1135:   ISDestroy(&bnk->active_fixed);
1136:   ISDestroy(&bnk->active_idx);
1137:   ISDestroy(&bnk->inactive_idx);
1138:   MatDestroy(&bnk->Hpre_inactive);
1139:   MatDestroy(&bnk->H_inactive);
1140:   TaoDestroy(&bnk->bncg);
1141:   PetscFree(tao->data);
1142:   return(0);
1143: }

1145: /*------------------------------------------------------------*/

1147: PetscErrorCode TaoSetFromOptions_BNK(PetscOptionItems *PetscOptionsObject,Tao tao)
1148: {
1149:   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
1151:   KSPType        ksp_type;

1154:   PetscOptionsHead(PetscOptionsObject,"Newton-Krylov method for bound constrained optimization");
1155:   PetscOptionsEList("-tao_bnk_init_type", "radius initialization type", "", BNK_INIT, BNK_INIT_TYPES, BNK_INIT[bnk->init_type], &bnk->init_type, 0);
1156:   PetscOptionsEList("-tao_bnk_update_type", "radius update type", "", BNK_UPDATE, BNK_UPDATE_TYPES, BNK_UPDATE[bnk->update_type], &bnk->update_type, 0);
1157:   PetscOptionsEList("-tao_bnk_as_type", "active set estimation method", "", BNK_AS, BNK_AS_TYPES, BNK_AS[bnk->as_type], &bnk->as_type, 0);
1158:   PetscOptionsReal("-tao_bnk_sval", "(developer) Hessian perturbation starting value", "", bnk->sval, &bnk->sval,NULL);
1159:   PetscOptionsReal("-tao_bnk_imin", "(developer) minimum initial Hessian perturbation", "", bnk->imin, &bnk->imin,NULL);
1160:   PetscOptionsReal("-tao_bnk_imax", "(developer) maximum initial Hessian perturbation", "", bnk->imax, &bnk->imax,NULL);
1161:   PetscOptionsReal("-tao_bnk_imfac", "(developer) initial merit factor for Hessian perturbation", "", bnk->imfac, &bnk->imfac,NULL);
1162:   PetscOptionsReal("-tao_bnk_pmin", "(developer) minimum Hessian perturbation", "", bnk->pmin, &bnk->pmin,NULL);
1163:   PetscOptionsReal("-tao_bnk_pmax", "(developer) maximum Hessian perturbation", "", bnk->pmax, &bnk->pmax,NULL);
1164:   PetscOptionsReal("-tao_bnk_pgfac", "(developer) Hessian perturbation growth factor", "", bnk->pgfac, &bnk->pgfac,NULL);
1165:   PetscOptionsReal("-tao_bnk_psfac", "(developer) Hessian perturbation shrink factor", "", bnk->psfac, &bnk->psfac,NULL);
1166:   PetscOptionsReal("-tao_bnk_pmgfac", "(developer) merit growth factor for Hessian perturbation", "", bnk->pmgfac, &bnk->pmgfac,NULL);
1167:   PetscOptionsReal("-tao_bnk_pmsfac", "(developer) merit shrink factor for Hessian perturbation", "", bnk->pmsfac, &bnk->pmsfac,NULL);
1168:   PetscOptionsReal("-tao_bnk_eta1", "(developer) threshold for rejecting step (-tao_bnk_update_type reduction)", "", bnk->eta1, &bnk->eta1,NULL);
1169:   PetscOptionsReal("-tao_bnk_eta2", "(developer) threshold for accepting marginal step (-tao_bnk_update_type reduction)", "", bnk->eta2, &bnk->eta2,NULL);
1170:   PetscOptionsReal("-tao_bnk_eta3", "(developer) threshold for accepting reasonable step (-tao_bnk_update_type reduction)", "", bnk->eta3, &bnk->eta3,NULL);
1171:   PetscOptionsReal("-tao_bnk_eta4", "(developer) threshold for accepting good step (-tao_bnk_update_type reduction)", "", bnk->eta4, &bnk->eta4,NULL);
1172:   PetscOptionsReal("-tao_bnk_alpha1", "(developer) radius reduction factor for rejected step (-tao_bnk_update_type reduction)", "", bnk->alpha1, &bnk->alpha1,NULL);
1173:   PetscOptionsReal("-tao_bnk_alpha2", "(developer) radius reduction factor for marginally accepted bad step (-tao_bnk_update_type reduction)", "", bnk->alpha2, &bnk->alpha2,NULL);
1174:   PetscOptionsReal("-tao_bnk_alpha3", "(developer) radius increase factor for reasonable accepted step (-tao_bnk_update_type reduction)", "", bnk->alpha3, &bnk->alpha3,NULL);
1175:   PetscOptionsReal("-tao_bnk_alpha4", "(developer) radius increase factor for good accepted step (-tao_bnk_update_type reduction)", "", bnk->alpha4, &bnk->alpha4,NULL);
1176:   PetscOptionsReal("-tao_bnk_alpha5", "(developer) radius increase factor for very good accepted step (-tao_bnk_update_type reduction)", "", bnk->alpha5, &bnk->alpha5,NULL);
1177:   PetscOptionsReal("-tao_bnk_nu1", "(developer) threshold for small line-search step length (-tao_bnk_update_type step)", "", bnk->nu1, &bnk->nu1,NULL);
1178:   PetscOptionsReal("-tao_bnk_nu2", "(developer) threshold for reasonable line-search step length (-tao_bnk_update_type step)", "", bnk->nu2, &bnk->nu2,NULL);
1179:   PetscOptionsReal("-tao_bnk_nu3", "(developer) threshold for large line-search step length (-tao_bnk_update_type step)", "", bnk->nu3, &bnk->nu3,NULL);
1180:   PetscOptionsReal("-tao_bnk_nu4", "(developer) threshold for very large line-search step length (-tao_bnk_update_type step)", "", bnk->nu4, &bnk->nu4,NULL);
1181:   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);
1182:   PetscOptionsReal("-tao_bnk_omega2", "(developer) radius reduction factor for small line-search step length (-tao_bnk_update_type step)", "", bnk->omega2, &bnk->omega2,NULL);
1183:   PetscOptionsReal("-tao_bnk_omega3", "(developer) radius factor for decent line-search step length (-tao_bnk_update_type step)", "", bnk->omega3, &bnk->omega3,NULL);
1184:   PetscOptionsReal("-tao_bnk_omega4", "(developer) radius increase factor for large line-search step length (-tao_bnk_update_type step)", "", bnk->omega4, &bnk->omega4,NULL);
1185:   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);
1186:   PetscOptionsReal("-tao_bnk_mu1_i", "(developer) threshold for accepting very good step (-tao_bnk_init_type interpolation)", "", bnk->mu1_i, &bnk->mu1_i,NULL);
1187:   PetscOptionsReal("-tao_bnk_mu2_i", "(developer) threshold for accepting good step (-tao_bnk_init_type interpolation)", "", bnk->mu2_i, &bnk->mu2_i,NULL);
1188:   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);
1189:   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);
1190:   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);
1191:   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);
1192:   PetscOptionsReal("-tao_bnk_theta_i", "(developer) trust region interpolation factor (-tao_bnk_init_type interpolation)", "", bnk->theta_i, &bnk->theta_i,NULL);
1193:   PetscOptionsReal("-tao_bnk_mu1", "(developer) threshold for accepting very good step (-tao_bnk_update_type interpolation)", "", bnk->mu1, &bnk->mu1,NULL);
1194:   PetscOptionsReal("-tao_bnk_mu2", "(developer) threshold for accepting good step (-tao_bnk_update_type interpolation)", "", bnk->mu2, &bnk->mu2,NULL);
1195:   PetscOptionsReal("-tao_bnk_gamma1", "(developer) radius reduction factor for rejected very bad step (-tao_bnk_update_type interpolation)", "", bnk->gamma1, &bnk->gamma1,NULL);
1196:   PetscOptionsReal("-tao_bnk_gamma2", "(developer) radius reduction factor for rejected bad step (-tao_bnk_update_type interpolation)", "", bnk->gamma2, &bnk->gamma2,NULL);
1197:   PetscOptionsReal("-tao_bnk_gamma3", "(developer) radius increase factor for accepted good step (-tao_bnk_update_type interpolation)", "", bnk->gamma3, &bnk->gamma3,NULL);
1198:   PetscOptionsReal("-tao_bnk_gamma4", "(developer) radius increase factor for accepted very good step (-tao_bnk_update_type interpolation)", "", bnk->gamma4, &bnk->gamma4,NULL);
1199:   PetscOptionsReal("-tao_bnk_theta", "(developer) trust region interpolation factor (-tao_bnk_update_type interpolation)", "", bnk->theta, &bnk->theta,NULL);
1200:   PetscOptionsReal("-tao_bnk_min_radius", "(developer) lower bound on initial radius", "", bnk->min_radius, &bnk->min_radius,NULL);
1201:   PetscOptionsReal("-tao_bnk_max_radius", "(developer) upper bound on radius", "", bnk->max_radius, &bnk->max_radius,NULL);
1202:   PetscOptionsReal("-tao_bnk_epsilon", "(developer) tolerance used when computing actual and predicted reduction", "", bnk->epsilon, &bnk->epsilon,NULL);
1203:   PetscOptionsReal("-tao_bnk_as_tol", "(developer) initial tolerance used when estimating actively bounded variables", "", bnk->as_tol, &bnk->as_tol,NULL);
1204:   PetscOptionsReal("-tao_bnk_as_step", "(developer) step length used when estimating actively bounded variables", "", bnk->as_step, &bnk->as_step,NULL);
1205:   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);
1206:   PetscOptionsTail();
1207:   TaoSetFromOptions(bnk->bncg);
1208:   TaoLineSearchSetFromOptions(tao->linesearch);
1209:   KSPSetFromOptions(tao->ksp);
1210:   KSPGetType(tao->ksp,&ksp_type);
1211:   PetscStrcmp(ksp_type,KSPNASH,&bnk->is_nash);
1212:   PetscStrcmp(ksp_type,KSPSTCG,&bnk->is_stcg);
1213:   PetscStrcmp(ksp_type,KSPGLTR,&bnk->is_gltr);
1214:   return(0);
1215: }

1217: /*------------------------------------------------------------*/

1219: PetscErrorCode TaoView_BNK(Tao tao, PetscViewer viewer)
1220: {
1221:   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
1222:   PetscInt       nrejects;
1223:   PetscBool      isascii;

1227:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
1228:   if (isascii) {
1229:     PetscViewerASCIIPushTab(viewer);
1230:     if (bnk->M) {
1231:       MatLMVMGetRejectCount(bnk->M,&nrejects);
1232:       PetscViewerASCIIPrintf(viewer, "Rejected BFGS updates: %D\n",nrejects);
1233:     }
1234:     PetscViewerASCIIPrintf(viewer, "CG steps: %D\n", bnk->tot_cg_its);
1235:     PetscViewerASCIIPrintf(viewer, "Newton steps: %D\n", bnk->newt);
1236:     if (bnk->M) {
1237:       PetscViewerASCIIPrintf(viewer, "BFGS steps: %D\n", bnk->bfgs);
1238:     }
1239:     PetscViewerASCIIPrintf(viewer, "Scaled gradient steps: %D\n", bnk->sgrad);
1240:     PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", bnk->grad);
1241:     PetscViewerASCIIPrintf(viewer, "KSP termination reasons:\n");
1242:     PetscViewerASCIIPrintf(viewer, "  atol: %D\n", bnk->ksp_atol);
1243:     PetscViewerASCIIPrintf(viewer, "  rtol: %D\n", bnk->ksp_rtol);
1244:     PetscViewerASCIIPrintf(viewer, "  ctol: %D\n", bnk->ksp_ctol);
1245:     PetscViewerASCIIPrintf(viewer, "  negc: %D\n", bnk->ksp_negc);
1246:     PetscViewerASCIIPrintf(viewer, "  dtol: %D\n", bnk->ksp_dtol);
1247:     PetscViewerASCIIPrintf(viewer, "  iter: %D\n", bnk->ksp_iter);
1248:     PetscViewerASCIIPrintf(viewer, "  othr: %D\n", bnk->ksp_othr);
1249:     PetscViewerASCIIPopTab(viewer);
1250:   }
1251:   return(0);
1252: }

1254: /* ---------------------------------------------------------- */

1256: /*MC
1257:   TAOBNK - Shared base-type for Bounded Newton-Krylov type algorithms.
1258:   At each iteration, the BNK methods solve the symmetric
1259:   system of equations to obtain the step diretion dk:
1260:               Hk dk = -gk
1261:   for free variables only. The step can be globalized either through 
1262:   trust-region methods, or a line search, or a heuristic mixture of both.

1264:     Options Database Keys:
1265: + -max_cg_its - maximum number of bounded conjugate-gradient iterations taken in each Newton loop
1266: . -init_type - trust radius initialization method ("constant", "direction", "interpolation")
1267: . -update_type - trust radius update method ("step", "direction", "interpolation")
1268: . -as_type - active-set estimation method ("none", "bertsekas")
1269: . -as_tol - (developer) initial tolerance used in estimating bounded active variables (-as_type bertsekas)
1270: . -as_step - (developer) trial step length used in estimating bounded active variables (-as_type bertsekas)
1271: . -sval - (developer) Hessian perturbation starting value
1272: . -imin - (developer) minimum initial Hessian perturbation
1273: . -imax - (developer) maximum initial Hessian perturbation
1274: . -pmin - (developer) minimum Hessian perturbation
1275: . -pmax - (developer) aximum Hessian perturbation
1276: . -pgfac - (developer) Hessian perturbation growth factor
1277: . -psfac - (developer) Hessian perturbation shrink factor
1278: . -imfac - (developer) initial merit factor for Hessian perturbation
1279: . -pmgfac - (developer) merit growth factor for Hessian perturbation
1280: . -pmsfac - (developer) merit shrink factor for Hessian perturbation
1281: . -eta1 - (developer) threshold for rejecting step (-update_type reduction)
1282: . -eta2 - (developer) threshold for accepting marginal step (-update_type reduction)
1283: . -eta3 - (developer) threshold for accepting reasonable step (-update_type reduction)
1284: . -eta4 - (developer) threshold for accepting good step (-update_type reduction)
1285: . -alpha1 - (developer) radius reduction factor for rejected step (-update_type reduction)
1286: . -alpha2 - (developer) radius reduction factor for marginally accepted bad step (-update_type reduction)
1287: . -alpha3 - (developer) radius increase factor for reasonable accepted step (-update_type reduction)
1288: . -alpha4 - (developer) radius increase factor for good accepted step (-update_type reduction)
1289: . -alpha5 - (developer) radius increase factor for very good accepted step (-update_type reduction)
1290: . -epsilon - (developer) tolerance for small pred/actual ratios that trigger automatic step acceptance (-update_type reduction)
1291: . -mu1 - (developer) threshold for accepting very good step (-update_type interpolation)
1292: . -mu2 - (developer) threshold for accepting good step (-update_type interpolation)
1293: . -gamma1 - (developer) radius reduction factor for rejected very bad step (-update_type interpolation)
1294: . -gamma2 - (developer) radius reduction factor for rejected bad step (-update_type interpolation)
1295: . -gamma3 - (developer) radius increase factor for accepted good step (-update_type interpolation)
1296: . -gamma4 - (developer) radius increase factor for accepted very good step (-update_type interpolation)
1297: . -theta - (developer) trust region interpolation factor (-update_type interpolation)
1298: . -nu1 - (developer) threshold for small line-search step length (-update_type step)
1299: . -nu2 - (developer) threshold for reasonable line-search step length (-update_type step)
1300: . -nu3 - (developer) threshold for large line-search step length (-update_type step)
1301: . -nu4 - (developer) threshold for very large line-search step length (-update_type step)
1302: . -omega1 - (developer) radius reduction factor for very small line-search step length (-update_type step)
1303: . -omega2 - (developer) radius reduction factor for small line-search step length (-update_type step)
1304: . -omega3 - (developer) radius factor for decent line-search step length (-update_type step)
1305: . -omega4 - (developer) radius increase factor for large line-search step length (-update_type step)
1306: . -omega5 - (developer) radius increase factor for very large line-search step length (-update_type step)
1307: . -mu1_i -  (developer) threshold for accepting very good step (-init_type interpolation)
1308: . -mu2_i -  (developer) threshold for accepting good step (-init_type interpolation)
1309: . -gamma1_i - (developer) radius reduction factor for rejected very bad step (-init_type interpolation)
1310: . -gamma2_i - (developer) radius reduction factor for rejected bad step (-init_type interpolation)
1311: . -gamma3_i - (developer) radius increase factor for accepted good step (-init_type interpolation)
1312: . -gamma4_i - (developer) radius increase factor for accepted very good step (-init_type interpolation)
1313: - -theta_i - (developer) trust region interpolation factor (-init_type interpolation)

1315:   Level: beginner
1316: M*/

1318: PetscErrorCode TaoCreate_BNK(Tao tao)
1319: {
1320:   TAO_BNK        *bnk;
1321:   const char     *morethuente_type = TAOLINESEARCHMT;
1323:   PC             pc;

1326:   PetscNewLog(tao,&bnk);

1328:   tao->ops->setup = TaoSetUp_BNK;
1329:   tao->ops->view = TaoView_BNK;
1330:   tao->ops->setfromoptions = TaoSetFromOptions_BNK;
1331:   tao->ops->destroy = TaoDestroy_BNK;

1333:   /*  Override default settings (unless already changed) */
1334:   if (!tao->max_it_changed) tao->max_it = 50;
1335:   if (!tao->trust0_changed) tao->trust0 = 100.0;

1337:   tao->data = (void*)bnk;
1338: 
1339:   /*  Hessian shifting parameters */
1340:   bnk->computehessian = TaoBNKComputeHessian;
1341:   bnk->computestep = TaoBNKComputeStep;
1342: 
1343:   bnk->sval   = 0.0;
1344:   bnk->imin   = 1.0e-4;
1345:   bnk->imax   = 1.0e+2;
1346:   bnk->imfac  = 1.0e-1;

1348:   bnk->pmin   = 1.0e-12;
1349:   bnk->pmax   = 1.0e+2;
1350:   bnk->pgfac  = 1.0e+1;
1351:   bnk->psfac  = 4.0e-1;
1352:   bnk->pmgfac = 1.0e-1;
1353:   bnk->pmsfac = 1.0e-1;

1355:   /*  Default values for trust-region radius update based on steplength */
1356:   bnk->nu1 = 0.25;
1357:   bnk->nu2 = 0.50;
1358:   bnk->nu3 = 1.00;
1359:   bnk->nu4 = 1.25;

1361:   bnk->omega1 = 0.25;
1362:   bnk->omega2 = 0.50;
1363:   bnk->omega3 = 1.00;
1364:   bnk->omega4 = 2.00;
1365:   bnk->omega5 = 4.00;

1367:   /*  Default values for trust-region radius update based on reduction */
1368:   bnk->eta1 = 1.0e-4;
1369:   bnk->eta2 = 0.25;
1370:   bnk->eta3 = 0.50;
1371:   bnk->eta4 = 0.90;

1373:   bnk->alpha1 = 0.25;
1374:   bnk->alpha2 = 0.50;
1375:   bnk->alpha3 = 1.00;
1376:   bnk->alpha4 = 2.00;
1377:   bnk->alpha5 = 4.00;

1379:   /*  Default values for trust-region radius update based on interpolation */
1380:   bnk->mu1 = 0.10;
1381:   bnk->mu2 = 0.50;

1383:   bnk->gamma1 = 0.25;
1384:   bnk->gamma2 = 0.50;
1385:   bnk->gamma3 = 2.00;
1386:   bnk->gamma4 = 4.00;

1388:   bnk->theta = 0.05;

1390:   /*  Default values for trust region initialization based on interpolation */
1391:   bnk->mu1_i = 0.35;
1392:   bnk->mu2_i = 0.50;

1394:   bnk->gamma1_i = 0.0625;
1395:   bnk->gamma2_i = 0.5;
1396:   bnk->gamma3_i = 2.0;
1397:   bnk->gamma4_i = 5.0;

1399:   bnk->theta_i = 0.25;

1401:   /*  Remaining parameters */
1402:   bnk->max_cg_its = 0;
1403:   bnk->min_radius = 1.0e-10;
1404:   bnk->max_radius = 1.0e10;
1405:   bnk->epsilon = PetscPowReal(PETSC_MACHINE_EPSILON, 2.0/3.0);
1406:   bnk->as_tol = 1.0e-3;
1407:   bnk->as_step = 1.0e-3;
1408:   bnk->dmin = 1.0e-6;
1409:   bnk->dmax = 1.0e6;
1410: 
1411:   bnk->M               = 0;
1412:   bnk->bfgs_pre        = 0;
1413:   bnk->init_type       = BNK_INIT_INTERPOLATION;
1414:   bnk->update_type     = BNK_UPDATE_REDUCTION;
1415:   bnk->as_type         = BNK_AS_BERTSEKAS;
1416: 
1417:   /* Create the embedded BNCG solver */
1418:   TaoCreate(PetscObjectComm((PetscObject)tao), &bnk->bncg);
1419:   PetscObjectIncrementTabLevel((PetscObject)bnk->bncg, (PetscObject)tao, 1);
1420:   TaoSetOptionsPrefix(bnk->bncg, "tao_bnk_");
1421:   TaoSetType(bnk->bncg, TAOBNCG);

1423:   /* Create the line search */
1424:   TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch);
1425:   PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);
1426:   TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);
1427:   TaoLineSearchSetType(tao->linesearch,morethuente_type);
1428:   TaoLineSearchUseTaoRoutines(tao->linesearch,tao);

1430:   /*  Set linear solver to default for symmetric matrices */
1431:   KSPCreate(((PetscObject)tao)->comm,&tao->ksp);
1432:   PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1);
1433:   KSPSetOptionsPrefix(tao->ksp,"tao_bnk_");
1434:   KSPSetType(tao->ksp,KSPSTCG);
1435:   KSPGetPC(tao->ksp, &pc);
1436:   PCSetType(pc, PCLMVM);
1437:   return(0);
1438: }