xref: /petsc/src/tao/bound/impls/bnk/bnk.c (revision e465cd6f30f93207df959d4c9fc4f47b6c5c3063)
1 #include <petsctaolinesearch.h>
2 #include <../src/tao/bound/impls/bnk/bnk.h>
3 
4 #include <petscksp.h>
5 
6 /* Routine for BFGS preconditioner */
7 
8 PetscErrorCode MatLMVMSolveShell(PC pc, Vec b, Vec x)
9 {
10   PetscErrorCode ierr;
11   Mat            M;
12 
13   PetscFunctionBegin;
14   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
15   PetscValidHeaderSpecific(b,VEC_CLASSID,2);
16   PetscValidHeaderSpecific(x,VEC_CLASSID,3);
17   ierr = PCShellGetContext(pc,(void**)&M);CHKERRQ(ierr);
18   ierr = MatLMVMSolveInactive(M, b, x);CHKERRQ(ierr);
19   PetscFunctionReturn(0);
20 }
21 
22 /*------------------------------------------------------------*/
23 
24 /* Routine for initializing the KSP solver, the BFGS preconditioner, and the initial trust radius estimation */
25 
26 PetscErrorCode TaoBNKInitialize(Tao tao)
27 {
28   PetscErrorCode               ierr;
29   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
30   KSPType                      ksp_type;
31   PC                           pc;
32 
33   PetscReal                    f_min, ftrial, prered, actred, kappa, sigma;
34   PetscReal                    tau, tau_1, tau_2, tau_max, tau_min, max_radius;
35   PetscReal                    delta, step = 1.0;
36 
37   PetscInt                     n,N,needH = 1;
38 
39   PetscInt                     i_max = 5;
40   PetscInt                     j_max = 1;
41   PetscInt                     i, j;
42 
43   PetscFunctionBegin;
44   /* Number of times ksp stopped because of these reasons */
45   bnk->ksp_atol = 0;
46   bnk->ksp_rtol = 0;
47   bnk->ksp_dtol = 0;
48   bnk->ksp_ctol = 0;
49   bnk->ksp_negc = 0;
50   bnk->ksp_iter = 0;
51   bnk->ksp_othr = 0;
52 
53   /* Initialize the Hessian perturbation */
54   bnk->pert = bnk->sval;
55 
56   /* Initialize trust-region radius when using nash, stcg, or gltr
57      Command automatically ignored for other methods
58      Will be reset during the first iteration
59   */
60   ierr = KSPGetType(tao->ksp,&ksp_type);CHKERRQ(ierr);
61   ierr = PetscStrcmp(ksp_type,KSPCGNASH,&bnk->is_nash);CHKERRQ(ierr);
62   ierr = PetscStrcmp(ksp_type,KSPCGSTCG,&bnk->is_stcg);CHKERRQ(ierr);
63   ierr = PetscStrcmp(ksp_type,KSPCGGLTR,&bnk->is_gltr);CHKERRQ(ierr);
64 
65   ierr = KSPCGSetRadius(tao->ksp,bnk->max_radius);CHKERRQ(ierr);
66 
67   if (bnk->is_nash || bnk->is_stcg || bnk->is_gltr) {
68     if (tao->trust0 < 0.0) SETERRQ(PETSC_COMM_SELF,1,"Initial radius negative");
69     tao->trust = tao->trust0;
70     tao->trust = PetscMax(tao->trust, bnk->min_radius);
71     tao->trust = PetscMin(tao->trust, bnk->max_radius);
72   }
73 
74   /* Get vectors we will need */
75   if (BNK_PC_BFGS == bnk->pc_type && !bnk->M) {
76     ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr);
77     ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr);
78     ierr = MatCreateLMVM(((PetscObject)tao)->comm,n,N,&bnk->M);CHKERRQ(ierr);
79     ierr = MatLMVMAllocateVectors(bnk->M,tao->solution);CHKERRQ(ierr);
80   }
81 
82   /* create vectors for the limited memory preconditioner */
83   if ((BNK_PC_BFGS == bnk->pc_type) && (BFGS_SCALE_BFGS != bnk->bfgs_scale_type)) {
84     if (!bnk->Diag) {
85       ierr = VecDuplicate(tao->solution,&bnk->Diag);CHKERRQ(ierr);
86     }
87   }
88 
89   /* Modify the preconditioner to use the bfgs approximation */
90   ierr = KSPGetPC(tao->ksp, &pc);CHKERRQ(ierr);
91   switch(bnk->pc_type) {
92   case BNK_PC_NONE:
93     ierr = PCSetType(pc, PCNONE);CHKERRQ(ierr);
94     ierr = PCSetFromOptions(pc);CHKERRQ(ierr);
95     break;
96 
97   case BNK_PC_AHESS:
98     ierr = PCSetType(pc, PCJACOBI);CHKERRQ(ierr);
99     ierr = PCSetFromOptions(pc);CHKERRQ(ierr);
100     ierr = PCJacobiSetUseAbs(pc,PETSC_TRUE);CHKERRQ(ierr);
101     break;
102 
103   case BNK_PC_BFGS:
104     ierr = PCSetType(pc, PCSHELL);CHKERRQ(ierr);
105     ierr = PCSetFromOptions(pc);CHKERRQ(ierr);
106     ierr = PCShellSetName(pc, "bfgs");CHKERRQ(ierr);
107     ierr = PCShellSetContext(pc, bnk->M);CHKERRQ(ierr);
108     ierr = PCShellSetApply(pc, MatLMVMSolveShell);CHKERRQ(ierr);
109     break;
110 
111   default:
112     /* Use the pc method set by pc_type */
113     break;
114   }
115 
116   /* Initialize trust-region radius.  The initialization is only performed
117      when we are using Nash, Steihaug-Toint or the Generalized Lanczos method. */
118   if (bnk->is_nash || bnk->is_stcg || bnk->is_gltr) {
119     switch(bnk->init_type) {
120     case BNK_INIT_CONSTANT:
121       /* Use the initial radius specified */
122       break;
123 
124     case BNK_INIT_INTERPOLATION:
125       /* Use the initial radius specified */
126       max_radius = 0.0;
127 
128       for (j = 0; j < j_max; ++j) {
129         f_min = bnk->f;
130         sigma = 0.0;
131 
132         if (needH) {
133           /* Compute the Hessian at the new step, and extract the inactive subsystem */
134           ierr = TaoComputeHessian(tao, tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr);
135           ierr = ISDestroy(&bnk->inactive_idx);CHKERRQ(ierr);
136           ierr = VecWhichInactive(tao->XL,tao->solution,bnk->unprojected_gradient,tao->XU,PETSC_TRUE,&bnk->inactive_idx);CHKERRQ(ierr);
137           ierr = TaoMatGetSubMat(tao->hessian, bnk->inactive_idx, bnk->Xwork, TAO_SUBSET_MASK, &bnk->H_inactive);CHKERRQ(ierr);
138           needH = 0;
139         }
140 
141         for (i = 0; i < i_max; ++i) {
142           /* Take a steepest descent step and snap it to bounds */
143           ierr = VecCopy(tao->solution, bnk->Xold);CHKERRQ(ierr);
144           ierr = VecAXPY(tao->solution, -tao->trust/bnk->gnorm, tao->gradient);CHKERRQ(ierr);
145           ierr = VecMedian(tao->XL, tao->solution, tao->XU, tao->solution);CHKERRQ(ierr);
146           /* Recompute the step after bound snapping so that it can be used in predicted decrease calculation later */
147           ierr = VecCopy(tao->solution, bnk->W);CHKERRQ(ierr);
148           ierr = VecAXPY(bnk->W, -1.0, bnk->Xold);CHKERRQ(ierr);
149           /* Compute the objective at the trial */
150           ierr = TaoComputeObjective(tao, tao->solution, &ftrial);CHKERRQ(ierr);
151           ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr);
152           if (PetscIsInfOrNanReal(ftrial)) {
153             tau = bnk->gamma1_i;
154           } else {
155             if (ftrial < f_min) {
156               f_min = ftrial;
157               sigma = -tao->trust / bnk->gnorm;
158             }
159             /* Compute the predicted and actual reduction */
160             ierr = MatMult(bnk->H_inactive, tao->gradient, bnk->W);CHKERRQ(ierr);
161             ierr = VecDot(tao->gradient, bnk->W, &prered);CHKERRQ(ierr);
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             }
169 
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);
174 
175             if (PetscAbsScalar(kappa - 1.0) <= bnk->mu1_i) {
176               /* Great agreement */
177               max_radius = PetscMax(max_radius, tao->trust);
178 
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 if (tau_1 >= 1.0 && tau_1 <= bnk->gamma4_i && tau_2 < 1.0) {
184                 tau = tau_1;
185               } else if (tau_2 >= 1.0 && tau_2 <= bnk->gamma4_i && tau_1 < 1.0) {
186                 tau = tau_2;
187               } else {
188                 tau = tau_max;
189               }
190             } else if (PetscAbsScalar(kappa - 1.0) <= bnk->mu2_i) {
191               /* Good agreement */
192               max_radius = PetscMax(max_radius, tao->trust);
193 
194               if (tau_max < bnk->gamma2_i) {
195                 tau = bnk->gamma2_i;
196               } else if (tau_max > bnk->gamma3_i) {
197                 tau = bnk->gamma3_i;
198               } else {
199                 tau = tau_max;
200               }
201             } else {
202               /* Not good agreement */
203               if (tau_min > 1.0) {
204                 tau = bnk->gamma2_i;
205               } else if (tau_max < bnk->gamma1_i) {
206                 tau = bnk->gamma1_i;
207               } else if ((tau_min < bnk->gamma1_i) && (tau_max >= 1.0)) {
208                 tau = bnk->gamma1_i;
209               } else if ((tau_1 >= bnk->gamma1_i) && (tau_1 < 1.0) && ((tau_2 < bnk->gamma1_i) || (tau_2 >= 1.0))) {
210                 tau = tau_1;
211               } else if ((tau_2 >= bnk->gamma1_i) && (tau_2 < 1.0) && ((tau_1 < bnk->gamma1_i) || (tau_2 >= 1.0))) {
212                 tau = tau_2;
213               } else {
214                 tau = tau_max;
215               }
216             }
217           }
218           tao->trust = tau * tao->trust;
219         }
220 
221         if (f_min < bnk->f) {
222           bnk->f = f_min;
223           ierr = VecAXPY(tao->solution,sigma,tao->gradient);CHKERRQ(ierr);
224           ierr = VecMedian(tao->XL, tao->solution, tao->XU, tao->solution);CHKERRQ(ierr);
225           ierr = TaoComputeGradient(tao,tao->solution,bnk->unprojected_gradient);CHKERRQ(ierr);
226           ierr = VecBoundGradientProjection(bnk->unprojected_gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr);
227 
228           ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&bnk->gnorm);CHKERRQ(ierr);
229           if (PetscIsInfOrNanReal(bnk->gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute gradient generated Inf or NaN");
230           needH = 1;
231 
232           ierr = TaoLogConvergenceHistory(tao,bnk->f,bnk->gnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
233           ierr = TaoMonitor(tao,tao->niter,bnk->f,bnk->gnorm,0.0,step);CHKERRQ(ierr);
234           ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
235           if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
236         }
237       }
238       tao->trust = PetscMax(tao->trust, max_radius);
239 
240       /* Modify the radius if it is too large or small */
241       tao->trust = PetscMax(tao->trust, bnk->min_radius);
242       tao->trust = PetscMin(tao->trust, bnk->max_radius);
243       break;
244 
245     default:
246       /* Norm of the first direction will initialize radius */
247       tao->trust = 0.0;
248       break;
249     }
250   }
251 
252   /* Set initial scaling for the BFGS preconditioner
253      This step is done after computing the initial trust-region radius
254      since the function value may have decreased */
255   if (BNK_PC_BFGS == bnk->pc_type) {
256     if (bnk->f != 0.0) {
257       delta = 2.0 * PetscAbsScalar(bnk->f) / (bnk->gnorm*bnk->gnorm);
258     } else {
259       delta = 2.0 / (bnk->gnorm*bnk->gnorm);
260     }
261     ierr = MatLMVMSetDelta(bnk->M,delta);CHKERRQ(ierr);
262   }
263 
264   /* Set counter for gradient/reset steps*/
265   bnk->newt = 0;
266   bnk->bfgs = 0;
267   bnk->sgrad = 0;
268   bnk->grad = 0;
269   PetscFunctionReturn(0);
270 }
271 
272 /*------------------------------------------------------------*/
273 
274 /* Routine for estimating the active set */
275 
276 PetscErrorCode TaoBNKEstimateActiveSet(Tao tao)
277 {
278   PetscErrorCode               ierr;
279   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
280 
281   PetscFunctionBegin;
282   switch (bnk->as_type) {
283   case BNK_AS_NONE:
284     ierr = ISDestroy(&bnk->inactive_idx);CHKERRQ(ierr);
285     ierr = VecWhichInactive(tao->XL, tao->solution, bnk->unprojected_gradient, tao->XU, PETSC_TRUE, &bnk->inactive_idx);CHKERRQ(ierr);
286     ierr = ISDestroy(&bnk->active_idx);CHKERRQ(ierr);
287     ierr = ISComplementVec(bnk->inactive_idx, tao->solution, &bnk->active_idx);CHKERRQ(ierr);
288     break;
289 
290   case BNK_AS_BERTSEKAS:
291     /* Compute the trial step vector with which we will estimate the active set at the next iteration */
292     if (BNK_PC_BFGS == bnk->pc_type) {
293       /* If the BFGS preconditioner matrix is available, we will construct a trial step with it */
294       ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, bnk->W);CHKERRQ(ierr);
295     } else {
296       /* BFGS preconditioner doesn't exist so let's invert the diagonal of the Hessian instead onto the gradient*/
297       ierr = MatGetDiagonal(tao->hessian, bnk->Xwork);CHKERRQ(ierr);
298       ierr = VecReciprocal(bnk->Xwork);CHKERRQ(ierr);CHKERRQ(ierr);
299       ierr = VecPointwiseMult(bnk->W, bnk->Xwork, bnk->unprojected_gradient);CHKERRQ(ierr);
300     }
301     ierr = VecScale(bnk->W, -1.0);CHKERRQ(ierr);
302     ierr = TaoEstimateActiveBounds(tao->solution, tao->XL, tao->XU, bnk->unprojected_gradient, bnk->W, bnk->as_step, &bnk->as_tol, &bnk->active_lower, &bnk->active_upper, &bnk->active_fixed, &bnk->active_idx, &bnk->inactive_idx);CHKERRQ(ierr);
303 
304   default:
305     break;
306   }
307   if (bnk->active_idx) {
308     ierr = ISGetSize(bnk->active_idx, &bnk->num_active);CHKERRQ(ierr);
309   } else {
310     bnk->num_active = 0;
311   }
312   PetscFunctionReturn(0);
313 }
314 
315 /*------------------------------------------------------------*/
316 
317 /* Routine for bounding the step direction */
318 
319 PetscErrorCode TaoBNKBoundStep(Tao tao, Vec step)
320 {
321   PetscErrorCode               ierr;
322   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
323 
324   PetscFunctionBegin;
325   if (bnk->num_active > 0) {
326     switch (bnk->as_type) {
327     case BNK_AS_NONE:
328       if (bnk->active_idx) {
329         ierr = VecGetSubVector(step, bnk->active_idx, &bnk->active_work);CHKERRQ(ierr);
330         ierr = VecSet(bnk->active_work, 0.0);CHKERRQ(ierr);
331         ierr = VecRestoreSubVector(step, bnk->active_idx, &bnk->active_work);CHKERRQ(ierr);
332       }
333       break;
334 
335     case BNK_AS_BERTSEKAS:
336       ierr = TaoBoundStep(tao->solution, tao->XL, tao->XU, bnk->active_lower, bnk->active_upper, bnk->active_fixed, step);CHKERRQ(ierr);
337       break;
338 
339     default:
340       break;
341     }
342   }
343   PetscFunctionReturn(0);
344 }
345 
346 /*------------------------------------------------------------*/
347 
348 /* Routine for computing the Newton step.
349 
350   If the safeguard is enabled, the Newton step is verified to be a
351   descent direction, with fallbacks onto BFGS, scaled gradient, and unscaled
352   gradient steps if/when necessary.
353 
354   The function reports back on which type of step has ultimately been stored
355   under tao->stepdirection.
356 */
357 
358 PetscErrorCode TaoBNKComputeStep(Tao tao, KSPConvergedReason *ksp_reason)
359 {
360   PetscErrorCode               ierr;
361   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
362 
363   PetscReal                    delta;
364   PetscInt                     bfgsUpdates = 0;
365   PetscInt                     kspits;
366 
367   PetscFunctionBegin;
368   /* Determine the active and inactive sets */
369   ierr = TaoBNKEstimateActiveSet(tao);CHKERRQ(ierr);
370 
371   /* Prepare masked matrices for the inactive set */
372   if (BNK_PC_BFGS == bnk->pc_type) { ierr = MatLMVMSetInactive(bnk->M, bnk->inactive_idx);CHKERRQ(ierr); }
373   if (bnk->inactive_idx) {
374     ierr = TaoMatGetSubMat(tao->hessian, bnk->inactive_idx, bnk->Xwork, TAO_SUBSET_MASK, &bnk->H_inactive);CHKERRQ(ierr);
375     if (tao->hessian == tao->hessian_pre) {
376       bnk->Hpre_inactive = bnk->H_inactive;
377     } else {
378       ierr = TaoMatGetSubMat(tao->hessian_pre, bnk->inactive_idx, bnk->Xwork, TAO_SUBSET_MASK, &bnk->Hpre_inactive);CHKERRQ(ierr);
379     }
380   } else {
381     bnk->H_inactive = tao->hessian;
382     bnk->Hpre_inactive = tao->hessian_pre;
383   }
384 
385   /* Solve the Newton system of equations */
386   ierr = VecSet(tao->stepdirection, 0.0);CHKERRQ(ierr);
387   ierr = KSPSetOperators(tao->ksp,bnk->H_inactive,bnk->Hpre_inactive);CHKERRQ(ierr);
388   if (bnk->is_nash || bnk->is_stcg || bnk->is_gltr) {
389     ierr = KSPCGSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr);
390     ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr);
391     ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr);
392     tao->ksp_its+=kspits;
393     tao->ksp_tot_its+=kspits;
394     ierr = KSPCGGetNormD(tao->ksp,&bnk->dnorm);CHKERRQ(ierr);
395 
396     if (0.0 == tao->trust) {
397       /* Radius was uninitialized; use the norm of the direction */
398       if (bnk->dnorm > 0.0) {
399         tao->trust = bnk->dnorm;
400 
401         /* Modify the radius if it is too large or small */
402         tao->trust = PetscMax(tao->trust, bnk->min_radius);
403         tao->trust = PetscMin(tao->trust, bnk->max_radius);
404       } else {
405         /* The direction was bad; set radius to default value and re-solve
406            the trust-region subproblem to get a direction */
407         tao->trust = tao->trust0;
408 
409         /* Modify the radius if it is too large or small */
410         tao->trust = PetscMax(tao->trust, bnk->min_radius);
411         tao->trust = PetscMin(tao->trust, bnk->max_radius);
412 
413         ierr = KSPCGSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr);
414         ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr);
415         ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr);
416         tao->ksp_its+=kspits;
417         tao->ksp_tot_its+=kspits;
418         ierr = KSPCGGetNormD(tao->ksp,&bnk->dnorm);CHKERRQ(ierr);
419 
420         if (bnk->dnorm == 0.0) SETERRQ(PETSC_COMM_SELF,1, "Initial direction zero");
421       }
422     }
423   } else {
424     ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr);
425     ierr = KSPGetIterationNumber(tao->ksp, &kspits);CHKERRQ(ierr);
426     tao->ksp_its += kspits;
427     tao->ksp_tot_its+=kspits;
428   }
429   /* Make sure the safeguarded fall-back step is zero for actively bounded variables */
430   ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
431   ierr = TaoBNKBoundStep(tao, tao->stepdirection);CHKERRQ(ierr);
432 
433   /* Record convergence reasons */
434   ierr = KSPGetConvergedReason(tao->ksp, ksp_reason);CHKERRQ(ierr);
435   if (KSP_CONVERGED_ATOL == *ksp_reason) {
436     ++bnk->ksp_atol;
437   } else if (KSP_CONVERGED_RTOL == *ksp_reason) {
438     ++bnk->ksp_rtol;
439   } else if (KSP_CONVERGED_CG_CONSTRAINED == *ksp_reason) {
440     ++bnk->ksp_ctol;
441   } else if (KSP_CONVERGED_CG_NEG_CURVE == *ksp_reason) {
442     ++bnk->ksp_negc;
443   } else if (KSP_DIVERGED_DTOL == *ksp_reason) {
444     ++bnk->ksp_dtol;
445   } else if (KSP_DIVERGED_ITS == *ksp_reason) {
446     ++bnk->ksp_iter;
447   } else {
448     ++bnk->ksp_othr;
449   }
450 
451   /* Make sure the BFGS preconditioner is healthy */
452   if (bnk->pc_type == BNK_PC_BFGS) {
453     ierr = MatLMVMGetUpdates(bnk->M, &bfgsUpdates);CHKERRQ(ierr);
454     if ((KSP_DIVERGED_INDEFINITE_PC == *ksp_reason) && (bfgsUpdates > 1)) {
455       /* Preconditioner is numerically indefinite; reset the approximation. */
456       if (bnk->f != 0.0) {
457         delta = 2.0 * PetscAbsScalar(bnk->f) / (bnk->gnorm*bnk->gnorm);
458       } else {
459         delta = 2.0 / (bnk->gnorm*bnk->gnorm);
460       }
461       ierr = MatLMVMSetDelta(bnk->M,delta);CHKERRQ(ierr);
462       ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr);
463       ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
464     }
465   }
466   PetscFunctionReturn(0);
467 }
468 
469 PetscErrorCode TaoBNKSafeguardStep(Tao tao, KSPConvergedReason ksp_reason, PetscInt *stepType)
470 {
471   PetscErrorCode               ierr;
472   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
473 
474   PetscReal                    gdx, delta, e_min;
475   PetscInt                     bfgsUpdates;
476 
477   PetscFunctionBegin;
478   ierr = VecDot(tao->stepdirection, tao->gradient, &gdx);CHKERRQ(ierr);
479   if ((gdx >= 0.0) || PetscIsInfOrNanReal(gdx)) {
480     /* Newton step is not descent or direction produced Inf or NaN
481        Update the perturbation for next time */
482     if (bnk->pert <= 0.0) {
483       /* Initialize the perturbation */
484       bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
485       if (bnk->is_gltr) {
486         ierr = KSPCGGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr);
487         bnk->pert = PetscMax(bnk->pert, -e_min);
488       }
489     } else {
490       /* Increase the perturbation */
491       bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
492     }
493 
494     if (BNK_PC_BFGS != bnk->pc_type) {
495       /* We don't have the bfgs matrix around and updated
496          Must use gradient direction in this case */
497       ierr = VecCopy(tao->gradient, tao->stepdirection);CHKERRQ(ierr);
498       *stepType = BNK_GRADIENT;
499     } else {
500       /* Attempt to use the BFGS direction */
501       ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
502 
503       /* Check for success (descent direction)
504          NOTE: Negative gdx here means not a descent direction because
505          the fall-back step is missing a negative sign. */
506       ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr);
507       if ((gdx <= 0) || PetscIsInfOrNanReal(gdx)) {
508         /* BFGS direction is not descent or direction produced not a number
509            We can assert bfgsUpdates > 1 in this case because
510            the first solve produces the scaled gradient direction,
511            which is guaranteed to be descent */
512 
513         /* Use steepest descent direction (scaled) */
514         if (bnk->f != 0.0) {
515           delta = 2.0 * PetscAbsScalar(bnk->f) / (bnk->gnorm*bnk->gnorm);
516         } else {
517           delta = 2.0 / (bnk->gnorm*bnk->gnorm);
518         }
519         ierr = MatLMVMSetDelta(bnk->M, delta);CHKERRQ(ierr);
520         ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr);
521         ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
522         ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
523 
524         *stepType = BNK_SCALED_GRADIENT;
525       } else {
526         ierr = MatLMVMGetUpdates(bnk->M, &bfgsUpdates);CHKERRQ(ierr);
527         if (1 == bfgsUpdates) {
528           /* The first BFGS direction is always the scaled gradient */
529           *stepType = BNK_SCALED_GRADIENT;
530         } else {
531           *stepType = BNK_BFGS;
532         }
533       }
534     }
535     /* Make sure the safeguarded fall-back step is zero for actively bounded variables */
536     ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
537     ierr = TaoBNKBoundStep(tao, tao->stepdirection);CHKERRQ(ierr);
538   } else {
539     /* Computed Newton step is descent */
540     switch (ksp_reason) {
541     case KSP_DIVERGED_NANORINF:
542     case KSP_DIVERGED_BREAKDOWN:
543     case KSP_DIVERGED_INDEFINITE_MAT:
544     case KSP_DIVERGED_INDEFINITE_PC:
545     case KSP_CONVERGED_CG_NEG_CURVE:
546       /* Matrix or preconditioner is indefinite; increase perturbation */
547       if (bnk->pert <= 0.0) {
548         /* Initialize the perturbation */
549         bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
550         if (bnk->is_gltr) {
551           ierr = KSPCGGLTRGetMinEig(tao->ksp, &e_min);CHKERRQ(ierr);
552           bnk->pert = PetscMax(bnk->pert, -e_min);
553         }
554       } else {
555         /* Increase the perturbation */
556         bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
557       }
558       break;
559 
560     default:
561       /* Newton step computation is good; decrease perturbation */
562       bnk->pert = PetscMin(bnk->psfac * bnk->pert, bnk->pmsfac * bnk->gnorm);
563       if (bnk->pert < bnk->pmin) {
564         bnk->pert = 0.0;
565       }
566       break;
567     }
568     *stepType = BNK_NEWTON;
569   }
570   PetscFunctionReturn(0);
571 }
572 
573 /*------------------------------------------------------------*/
574 
575 /* Routine for performing a bound-projected More-Thuente line search.
576 
577   Includes fallbacks to BFGS, scaled gradient, and unscaled gradient steps if the
578   Newton step does not produce a valid step length.
579 */
580 
581 PetscErrorCode TaoBNKPerformLineSearch(Tao tao, PetscInt stepType, PetscReal *steplen, TaoLineSearchConvergedReason *reason)
582 {
583   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
584   PetscErrorCode ierr;
585   TaoLineSearchConvergedReason ls_reason;
586 
587   PetscReal      e_min, gdx, delta;
588   PetscInt       bfgsUpdates;
589 
590   PetscFunctionBegin;
591   /* Perform the linesearch */
592   *steplen = 1.0;
593   ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, steplen, &ls_reason);CHKERRQ(ierr);
594   ierr = VecBoundGradientProjection(bnk->unprojected_gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr);
595   ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
596 
597   while (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER && stepType != BNK_GRADIENT) {
598     /* Linesearch failed, revert solution */
599     bnk->f = bnk->fold;
600     ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr);
601     ierr = VecCopy(bnk->Gold, tao->gradient);CHKERRQ(ierr);
602     ierr = VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient);CHKERRQ(ierr);
603 
604     switch(stepType) {
605     case BNK_NEWTON:
606       /* Failed to obtain acceptable iterate with Newton step
607          Update the perturbation for next time */
608       if (bnk->pert <= 0.0) {
609         /* Initialize the perturbation */
610         bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
611         if (bnk->is_gltr) {
612           ierr = KSPCGGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr);
613           bnk->pert = PetscMax(bnk->pert, -e_min);
614         }
615       } else {
616         /* Increase the perturbation */
617         bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
618       }
619 
620       if (BNK_PC_BFGS != bnk->pc_type) {
621         /* We don't have the bfgs matrix around and being updated
622            Must use gradient direction in this case */
623         ierr = VecCopy(tao->gradient, tao->stepdirection);CHKERRQ(ierr);
624         stepType = BNK_GRADIENT;
625       } else {
626         /* Attempt to use the BFGS direction */
627         ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
628         /* Check for success (descent direction)
629            NOTE: Negative gdx means not a descent direction because the step here is missing a negative sign. */
630         ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr);
631         if ((gdx <= 0) || PetscIsInfOrNanReal(gdx)) {
632           /* BFGS direction is not descent or direction produced not a number
633              We can assert bfgsUpdates > 1 in this case
634              Use steepest descent direction (scaled) */
635           if (bnk->f != 0.0) {
636             delta = 2.0 * PetscAbsScalar(bnk->f) / (bnk->gnorm*bnk->gnorm);
637           } else {
638             delta = 2.0 / (bnk->gnorm*bnk->gnorm);
639           }
640           ierr = MatLMVMSetDelta(bnk->M, delta);CHKERRQ(ierr);
641           ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr);
642           ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
643           ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
644 
645           bfgsUpdates = 1;
646           stepType = BNK_SCALED_GRADIENT;
647         } else {
648           ierr = MatLMVMGetUpdates(bnk->M, &bfgsUpdates);CHKERRQ(ierr);
649           if (1 == bfgsUpdates) {
650             /* The first BFGS direction is always the scaled gradient */
651             stepType = BNK_SCALED_GRADIENT;
652           } else {
653             stepType = BNK_BFGS;
654           }
655         }
656       }
657       break;
658 
659     case BNK_BFGS:
660       /* Can only enter if pc_type == BNK_PC_BFGS
661          Failed to obtain acceptable iterate with BFGS step
662          Attempt to use the scaled gradient direction */
663       if (bnk->f != 0.0) {
664         delta = 2.0 * PetscAbsScalar(bnk->f) / (bnk->gnorm*bnk->gnorm);
665       } else {
666         delta = 2.0 / (bnk->gnorm*bnk->gnorm);
667       }
668       ierr = MatLMVMSetDelta(bnk->M, delta);CHKERRQ(ierr);
669       ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr);
670       ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
671       ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
672 
673       bfgsUpdates = 1;
674       stepType = BNK_SCALED_GRADIENT;
675       break;
676 
677     case BNK_SCALED_GRADIENT:
678       /* Can only enter if pc_type == BNK_PC_BFGS
679          The scaled gradient step did not produce a new iterate;
680          attemp to use the gradient direction.
681          Need to make sure we are not using a different diagonal scaling */
682       ierr = MatLMVMSetScale(bnk->M,0);CHKERRQ(ierr);
683       ierr = MatLMVMSetDelta(bnk->M,1.0);CHKERRQ(ierr);
684       ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr);
685       ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
686       ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
687 
688       bfgsUpdates = 1;
689       stepType = BNK_GRADIENT;
690       break;
691     }
692     /* Make sure the safeguarded fall-back step is zero for actively bounded variables */
693     ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
694     ierr = TaoBNKBoundStep(tao, tao->stepdirection);CHKERRQ(ierr);
695 
696     /* Perform one last line search with the fall-back step */
697     *steplen = 1.0;
698     ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, steplen, &ls_reason);CHKERRQ(ierr);
699     ierr = VecBoundGradientProjection(bnk->unprojected_gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr);
700     ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
701   }
702   *reason = ls_reason;
703   PetscFunctionReturn(0);
704 }
705 
706 /*------------------------------------------------------------*/
707 
708 /* Routine for updating the trust radius.
709 
710   Function features three different update methods:
711   1) Line-search step length based
712   2) Predicted decrease on the CG quadratic model
713   3) Interpolation
714 */
715 
716 PetscErrorCode TaoBNKUpdateTrustRadius(Tao tao, PetscReal prered, PetscReal actred, PetscInt stepType, PetscBool *accept)
717 {
718   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
719   PetscErrorCode ierr;
720 
721   PetscReal      step, kappa;
722   PetscReal      gdx, tau_1, tau_2, tau_min, tau_max;
723 
724   PetscFunctionBegin;
725   /* Update trust region radius */
726   *accept = PETSC_FALSE;
727   switch(bnk->update_type) {
728   case BNK_UPDATE_STEP:
729     *accept = PETSC_TRUE; /* always accept here because line search succeeded */
730     if (stepType == BNK_NEWTON) {
731       ierr = TaoLineSearchGetStepLength(tao->linesearch, &step);CHKERRQ(ierr);
732       if (step < bnk->nu1) {
733         /* Very bad step taken; reduce radius */
734         tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust);
735       } else if (step < bnk->nu2) {
736         /* Reasonably bad step taken; reduce radius */
737         tao->trust = bnk->omega2 * PetscMin(bnk->dnorm, tao->trust);
738       } else if (step < bnk->nu3) {
739         /*  Reasonable step was taken; leave radius alone */
740         if (bnk->omega3 < 1.0) {
741           tao->trust = bnk->omega3 * PetscMin(bnk->dnorm, tao->trust);
742         } else if (bnk->omega3 > 1.0) {
743           tao->trust = PetscMax(bnk->omega3 * bnk->dnorm, tao->trust);
744         }
745       } else if (step < bnk->nu4) {
746         /*  Full step taken; increase the radius */
747         tao->trust = PetscMax(bnk->omega4 * bnk->dnorm, tao->trust);
748       } else {
749         /*  More than full step taken; increase the radius */
750         tao->trust = PetscMax(bnk->omega5 * bnk->dnorm, tao->trust);
751       }
752     } else {
753       /*  Newton step was not good; reduce the radius */
754       tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust);
755     }
756     break;
757 
758   case BNK_UPDATE_REDUCTION:
759     if (stepType == BNK_NEWTON) {
760       if (prered < 0.0) {
761         /* The predicted reduction has the wrong sign.  This cannot
762            happen in infinite precision arithmetic.  Step should
763            be rejected! */
764         tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
765       }
766       else {
767         if (PetscIsInfOrNanReal(actred)) {
768           tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
769         } else {
770           if ((PetscAbsScalar(actred) <= PetscMax(1.0, PetscAbsScalar(bnk->f))*bnk->epsilon) &&
771               (PetscAbsScalar(prered) <= PetscMax(1.0, PetscAbsScalar(bnk->f))*bnk->epsilon)) {
772             kappa = 1.0;
773           }
774           else {
775             kappa = actred / prered;
776           }
777 
778           /* Accept or reject the step and update radius */
779           if (kappa < bnk->eta1) {
780             /* Reject the step */
781             tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
782           }
783           else {
784             /* Accept the step */
785             *accept = PETSC_TRUE;
786             /* Update the trust region radius only if the computed step is at the trust radius boundary */
787             if (bnk->dnorm == tao->trust) {
788               if (kappa < bnk->eta2) {
789                 /* Marginal bad step */
790                 tao->trust = bnk->alpha2 * tao->trust;
791               }
792               else if (kappa < bnk->eta3) {
793                 /* Reasonable step */
794                 tao->trust = bnk->alpha3 * tao->trust;
795               }
796               else if (kappa < bnk->eta4) {
797                 /* Good step */
798                 tao->trust = bnk->alpha4 * tao->trust;
799               }
800               else {
801                 /* Very good step */
802                 tao->trust = bnk->alpha5 * tao->trust;
803               }
804             }
805           }
806         }
807       }
808     } else {
809       /*  Newton step was not good; reduce the radius */
810       tao->trust = bnk->alpha1 * PetscMin(bnk->dnorm, tao->trust);
811     }
812     break;
813 
814   default:
815     if (stepType == BNK_NEWTON) {
816       if (prered < 0.0) {
817         /*  The predicted reduction has the wrong sign.  This cannot */
818         /*  happen in infinite precision arithmetic.  Step should */
819         /*  be rejected! */
820         tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
821       } else {
822         if (PetscIsInfOrNanReal(actred)) {
823           tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
824         } else {
825           if ((PetscAbsScalar(actred) <= bnk->epsilon) && (PetscAbsScalar(prered) <= bnk->epsilon)) {
826             kappa = 1.0;
827           } else {
828             kappa = actred / prered;
829           }
830 
831           ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr);
832           tau_1 = bnk->theta * gdx / (bnk->theta * gdx - (1.0 - bnk->theta) * prered + actred);
833           tau_2 = bnk->theta * gdx / (bnk->theta * gdx + (1.0 + bnk->theta) * prered - actred);
834           tau_min = PetscMin(tau_1, tau_2);
835           tau_max = PetscMax(tau_1, tau_2);
836 
837           if (kappa >= 1.0 - bnk->mu1) {
838             /*  Great agreement */
839             *accept = PETSC_TRUE;
840             if (tau_max < 1.0) {
841               tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm);
842             } else if (tau_max > bnk->gamma4) {
843               tao->trust = PetscMax(tao->trust, bnk->gamma4 * bnk->dnorm);
844             } else {
845               tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm);
846             }
847           } else if (kappa >= 1.0 - bnk->mu2) {
848             /*  Good agreement */
849             *accept = PETSC_TRUE;
850             if (tau_max < bnk->gamma2) {
851               tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm);
852             } else if (tau_max > bnk->gamma3) {
853               tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm);
854             } else if (tau_max < 1.0) {
855               tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm);
856             } else {
857               tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm);
858             }
859           } else {
860             /*  Not good agreement */
861             if (tau_min > 1.0) {
862               tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm);
863             } else if (tau_max < bnk->gamma1) {
864               tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
865             } else if ((tau_min < bnk->gamma1) && (tau_max >= 1.0)) {
866               tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
867             } else if ((tau_1 >= bnk->gamma1) && (tau_1 < 1.0) && ((tau_2 < bnk->gamma1) || (tau_2 >= 1.0))) {
868               tao->trust = tau_1 * PetscMin(tao->trust, bnk->dnorm);
869             } else if ((tau_2 >= bnk->gamma1) && (tau_2 < 1.0) && ((tau_1 < bnk->gamma1) || (tau_2 >= 1.0))) {
870               tao->trust = tau_2 * PetscMin(tao->trust, bnk->dnorm);
871             } else {
872               tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm);
873             }
874           }
875         }
876       }
877     } else {
878       /*  Newton step was not good; reduce the radius */
879       tao->trust = bnk->gamma1 * PetscMin(bnk->dnorm, tao->trust);
880     }
881   }
882   /* Make sure the radius does not violate min and max settings */
883   tao->trust = PetscMin(tao->trust, bnk->max_radius);
884   tao->trust = PetscMax(tao->trust, bnk->min_radius);
885   PetscFunctionReturn(0);
886 }
887 
888 /* ---------------------------------------------------------- */
889 
890 static PetscErrorCode TaoSetUp_BNK(Tao tao)
891 {
892   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
893   PetscErrorCode ierr;
894 
895   PetscFunctionBegin;
896   if (!tao->gradient) {ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);}
897   if (!tao->stepdirection) {ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr);}
898   if (!bnk->W) {ierr = VecDuplicate(tao->solution,&bnk->W);CHKERRQ(ierr);}
899   if (!bnk->Xold) {ierr = VecDuplicate(tao->solution,&bnk->Xold);CHKERRQ(ierr);}
900   if (!bnk->Gold) {ierr = VecDuplicate(tao->solution,&bnk->Gold);CHKERRQ(ierr);}
901   if (!bnk->Xwork) {ierr = VecDuplicate(tao->solution,&bnk->Xwork);CHKERRQ(ierr);}
902   if (!bnk->Gwork) {ierr = VecDuplicate(tao->solution,&bnk->Gwork);CHKERRQ(ierr);}
903   if (!bnk->unprojected_gradient) {ierr = VecDuplicate(tao->solution,&bnk->unprojected_gradient);CHKERRQ(ierr);}
904   if (!bnk->unprojected_gradient_old) {ierr = VecDuplicate(tao->solution,&bnk->unprojected_gradient_old);CHKERRQ(ierr);}
905   bnk->Diag = 0;
906   bnk->M = 0;
907   bnk->H_inactive = 0;
908   bnk->Hpre_inactive = 0;
909   PetscFunctionReturn(0);
910 }
911 
912 /*------------------------------------------------------------*/
913 
914 static PetscErrorCode TaoDestroy_BNK(Tao tao)
915 {
916   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
917   PetscErrorCode ierr;
918 
919   PetscFunctionBegin;
920   if (tao->setupcalled) {
921     ierr = VecDestroy(&bnk->W);CHKERRQ(ierr);
922     ierr = VecDestroy(&bnk->Xold);CHKERRQ(ierr);
923     ierr = VecDestroy(&bnk->Gold);CHKERRQ(ierr);
924     ierr = VecDestroy(&bnk->Xwork);CHKERRQ(ierr);
925     ierr = VecDestroy(&bnk->Gwork);CHKERRQ(ierr);
926     ierr = VecDestroy(&bnk->unprojected_gradient);CHKERRQ(ierr);
927     ierr = VecDestroy(&bnk->unprojected_gradient_old);CHKERRQ(ierr);
928   }
929   ierr = VecDestroy(&bnk->Diag);CHKERRQ(ierr);
930   ierr = MatDestroy(&bnk->M);CHKERRQ(ierr);
931   ierr = PetscFree(tao->data);CHKERRQ(ierr);
932   PetscFunctionReturn(0);
933 }
934 
935 /*------------------------------------------------------------*/
936 
937 static PetscErrorCode TaoSetFromOptions_BNK(PetscOptionItems *PetscOptionsObject,Tao tao)
938 {
939   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
940   PetscErrorCode ierr;
941 
942   PetscFunctionBegin;
943   ierr = PetscOptionsHead(PetscOptionsObject,"Newton line search method for unconstrained optimization");CHKERRQ(ierr);
944   ierr = PetscOptionsEList("-tao_bnk_pc_type", "pc type", "", BNK_PC, BNK_PC_TYPES, BNK_PC[bnk->pc_type], &bnk->pc_type, 0);CHKERRQ(ierr);
945   ierr = PetscOptionsEList("-tao_bnk_bfgs_scale_type", "bfgs scale type", "", BFGS_SCALE, BFGS_SCALE_TYPES, BFGS_SCALE[bnk->bfgs_scale_type], &bnk->bfgs_scale_type, 0);CHKERRQ(ierr);
946   ierr = PetscOptionsEList("-tao_bnk_init_type", "radius initialization type", "", BNK_INIT, BNK_INIT_TYPES, BNK_INIT[bnk->init_type], &bnk->init_type, 0);CHKERRQ(ierr);
947   ierr = PetscOptionsEList("-tao_bnk_update_type", "radius update type", "", BNK_UPDATE, BNK_UPDATE_TYPES, BNK_UPDATE[bnk->update_type], &bnk->update_type, 0);CHKERRQ(ierr);
948   ierr = PetscOptionsEList("-tao_bnk_as_type", "active set estimation method", "", BNK_AS, BNK_AS_TYPES, BNK_AS[bnk->as_type], &bnk->as_type, 0);CHKERRQ(ierr);
949   ierr = PetscOptionsReal("-tao_bnk_sval", "perturbation starting value", "", bnk->sval, &bnk->sval,NULL);CHKERRQ(ierr);
950   ierr = PetscOptionsReal("-tao_bnk_imin", "minimum initial perturbation", "", bnk->imin, &bnk->imin,NULL);CHKERRQ(ierr);
951   ierr = PetscOptionsReal("-tao_bnk_imax", "maximum initial perturbation", "", bnk->imax, &bnk->imax,NULL);CHKERRQ(ierr);
952   ierr = PetscOptionsReal("-tao_bnk_imfac", "initial merit factor", "", bnk->imfac, &bnk->imfac,NULL);CHKERRQ(ierr);
953   ierr = PetscOptionsReal("-tao_bnk_pmin", "minimum perturbation", "", bnk->pmin, &bnk->pmin,NULL);CHKERRQ(ierr);
954   ierr = PetscOptionsReal("-tao_bnk_pmax", "maximum perturbation", "", bnk->pmax, &bnk->pmax,NULL);CHKERRQ(ierr);
955   ierr = PetscOptionsReal("-tao_bnk_pgfac", "growth factor", "", bnk->pgfac, &bnk->pgfac,NULL);CHKERRQ(ierr);
956   ierr = PetscOptionsReal("-tao_bnk_psfac", "shrink factor", "", bnk->psfac, &bnk->psfac,NULL);CHKERRQ(ierr);
957   ierr = PetscOptionsReal("-tao_bnk_pmgfac", "merit growth factor", "", bnk->pmgfac, &bnk->pmgfac,NULL);CHKERRQ(ierr);
958   ierr = PetscOptionsReal("-tao_bnk_pmsfac", "merit shrink factor", "", bnk->pmsfac, &bnk->pmsfac,NULL);CHKERRQ(ierr);
959   ierr = PetscOptionsReal("-tao_bnk_eta1", "poor steplength; reduce radius", "", bnk->eta1, &bnk->eta1,NULL);CHKERRQ(ierr);
960   ierr = PetscOptionsReal("-tao_bnk_eta2", "reasonable steplength; leave radius alone", "", bnk->eta2, &bnk->eta2,NULL);CHKERRQ(ierr);
961   ierr = PetscOptionsReal("-tao_bnk_eta3", "good steplength; increase radius", "", bnk->eta3, &bnk->eta3,NULL);CHKERRQ(ierr);
962   ierr = PetscOptionsReal("-tao_bnk_eta4", "excellent steplength; greatly increase radius", "", bnk->eta4, &bnk->eta4,NULL);CHKERRQ(ierr);
963   ierr = PetscOptionsReal("-tao_bnk_alpha1", "", "", bnk->alpha1, &bnk->alpha1,NULL);CHKERRQ(ierr);
964   ierr = PetscOptionsReal("-tao_bnk_alpha2", "", "", bnk->alpha2, &bnk->alpha2,NULL);CHKERRQ(ierr);
965   ierr = PetscOptionsReal("-tao_bnk_alpha3", "", "", bnk->alpha3, &bnk->alpha3,NULL);CHKERRQ(ierr);
966   ierr = PetscOptionsReal("-tao_bnk_alpha4", "", "", bnk->alpha4, &bnk->alpha4,NULL);CHKERRQ(ierr);
967   ierr = PetscOptionsReal("-tao_bnk_alpha5", "", "", bnk->alpha5, &bnk->alpha5,NULL);CHKERRQ(ierr);
968   ierr = PetscOptionsReal("-tao_bnk_nu1", "poor steplength; reduce radius", "", bnk->nu1, &bnk->nu1,NULL);CHKERRQ(ierr);
969   ierr = PetscOptionsReal("-tao_bnk_nu2", "reasonable steplength; leave radius alone", "", bnk->nu2, &bnk->nu2,NULL);CHKERRQ(ierr);
970   ierr = PetscOptionsReal("-tao_bnk_nu3", "good steplength; increase radius", "", bnk->nu3, &bnk->nu3,NULL);CHKERRQ(ierr);
971   ierr = PetscOptionsReal("-tao_bnk_nu4", "excellent steplength; greatly increase radius", "", bnk->nu4, &bnk->nu4,NULL);CHKERRQ(ierr);
972   ierr = PetscOptionsReal("-tao_bnk_omega1", "", "", bnk->omega1, &bnk->omega1,NULL);CHKERRQ(ierr);
973   ierr = PetscOptionsReal("-tao_bnk_omega2", "", "", bnk->omega2, &bnk->omega2,NULL);CHKERRQ(ierr);
974   ierr = PetscOptionsReal("-tao_bnk_omega3", "", "", bnk->omega3, &bnk->omega3,NULL);CHKERRQ(ierr);
975   ierr = PetscOptionsReal("-tao_bnk_omega4", "", "", bnk->omega4, &bnk->omega4,NULL);CHKERRQ(ierr);
976   ierr = PetscOptionsReal("-tao_bnk_omega5", "", "", bnk->omega5, &bnk->omega5,NULL);CHKERRQ(ierr);
977   ierr = PetscOptionsReal("-tao_bnk_mu1_i", "", "", bnk->mu1_i, &bnk->mu1_i,NULL);CHKERRQ(ierr);
978   ierr = PetscOptionsReal("-tao_bnk_mu2_i", "", "", bnk->mu2_i, &bnk->mu2_i,NULL);CHKERRQ(ierr);
979   ierr = PetscOptionsReal("-tao_bnk_gamma1_i", "", "", bnk->gamma1_i, &bnk->gamma1_i,NULL);CHKERRQ(ierr);
980   ierr = PetscOptionsReal("-tao_bnk_gamma2_i", "", "", bnk->gamma2_i, &bnk->gamma2_i,NULL);CHKERRQ(ierr);
981   ierr = PetscOptionsReal("-tao_bnk_gamma3_i", "", "", bnk->gamma3_i, &bnk->gamma3_i,NULL);CHKERRQ(ierr);
982   ierr = PetscOptionsReal("-tao_bnk_gamma4_i", "", "", bnk->gamma4_i, &bnk->gamma4_i,NULL);CHKERRQ(ierr);
983   ierr = PetscOptionsReal("-tao_bnk_theta_i", "", "", bnk->theta_i, &bnk->theta_i,NULL);CHKERRQ(ierr);
984   ierr = PetscOptionsReal("-tao_bnk_mu1", "", "", bnk->mu1, &bnk->mu1,NULL);CHKERRQ(ierr);
985   ierr = PetscOptionsReal("-tao_bnk_mu2", "", "", bnk->mu2, &bnk->mu2,NULL);CHKERRQ(ierr);
986   ierr = PetscOptionsReal("-tao_bnk_gamma1", "", "", bnk->gamma1, &bnk->gamma1,NULL);CHKERRQ(ierr);
987   ierr = PetscOptionsReal("-tao_bnk_gamma2", "", "", bnk->gamma2, &bnk->gamma2,NULL);CHKERRQ(ierr);
988   ierr = PetscOptionsReal("-tao_bnk_gamma3", "", "", bnk->gamma3, &bnk->gamma3,NULL);CHKERRQ(ierr);
989   ierr = PetscOptionsReal("-tao_bnk_gamma4", "", "", bnk->gamma4, &bnk->gamma4,NULL);CHKERRQ(ierr);
990   ierr = PetscOptionsReal("-tao_bnk_theta", "", "", bnk->theta, &bnk->theta,NULL);CHKERRQ(ierr);
991   ierr = PetscOptionsReal("-tao_bnk_min_radius", "lower bound on initial radius", "", bnk->min_radius, &bnk->min_radius,NULL);CHKERRQ(ierr);
992   ierr = PetscOptionsReal("-tao_bnk_max_radius", "upper bound on radius", "", bnk->max_radius, &bnk->max_radius,NULL);CHKERRQ(ierr);
993   ierr = PetscOptionsReal("-tao_bnk_epsilon", "tolerance used when computing actual and predicted reduction", "", bnk->epsilon, &bnk->epsilon,NULL);CHKERRQ(ierr);
994   ierr = PetscOptionsReal("-tao_bnk_as_tol", "initial tolerance used when estimating actively bounded variables", "", bnk->as_tol, &bnk->as_tol,NULL);CHKERRQ(ierr);
995   ierr = PetscOptionsReal("-tao_bnk_as_step", "step length used when estimating actively bounded variables", "", bnk->as_step, &bnk->as_step,NULL);CHKERRQ(ierr);
996   ierr = PetscOptionsTail();CHKERRQ(ierr);
997   ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr);
998   ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
999   PetscFunctionReturn(0);
1000 }
1001 
1002 /*------------------------------------------------------------*/
1003 
1004 static PetscErrorCode TaoView_BNK(Tao tao, PetscViewer viewer)
1005 {
1006   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
1007   PetscInt       nrejects;
1008   PetscBool      isascii;
1009   PetscErrorCode ierr;
1010 
1011   PetscFunctionBegin;
1012   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
1013   if (isascii) {
1014     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1015     if (BNK_PC_BFGS == bnk->pc_type && bnk->M) {
1016       ierr = MatLMVMGetRejects(bnk->M,&nrejects);CHKERRQ(ierr);
1017       ierr = PetscViewerASCIIPrintf(viewer, "Rejected matrix updates: %D\n",nrejects);CHKERRQ(ierr);
1018     }
1019     ierr = PetscViewerASCIIPrintf(viewer, "Newton steps: %D\n", bnk->newt);CHKERRQ(ierr);
1020     ierr = PetscViewerASCIIPrintf(viewer, "BFGS steps: %D\n", bnk->bfgs);CHKERRQ(ierr);
1021     ierr = PetscViewerASCIIPrintf(viewer, "Scaled gradient steps: %D\n", bnk->sgrad);CHKERRQ(ierr);
1022     ierr = PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", bnk->grad);CHKERRQ(ierr);
1023     ierr = PetscViewerASCIIPrintf(viewer, "KSP termination reasons:\n");CHKERRQ(ierr);
1024     ierr = PetscViewerASCIIPrintf(viewer, "  atol: %D\n", bnk->ksp_atol);CHKERRQ(ierr);
1025     ierr = PetscViewerASCIIPrintf(viewer, "  rtol: %D\n", bnk->ksp_rtol);CHKERRQ(ierr);
1026     ierr = PetscViewerASCIIPrintf(viewer, "  ctol: %D\n", bnk->ksp_ctol);CHKERRQ(ierr);
1027     ierr = PetscViewerASCIIPrintf(viewer, "  negc: %D\n", bnk->ksp_negc);CHKERRQ(ierr);
1028     ierr = PetscViewerASCIIPrintf(viewer, "  dtol: %D\n", bnk->ksp_dtol);CHKERRQ(ierr);
1029     ierr = PetscViewerASCIIPrintf(viewer, "  iter: %D\n", bnk->ksp_iter);CHKERRQ(ierr);
1030     ierr = PetscViewerASCIIPrintf(viewer, "  othr: %D\n", bnk->ksp_othr);CHKERRQ(ierr);
1031     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1032   }
1033   PetscFunctionReturn(0);
1034 }
1035 
1036 /* ---------------------------------------------------------- */
1037 
1038 /*MC
1039   TAOBNK - Shared base-type for Bounded Newton-Krylov type algorithms.
1040   At each iteration, the BNK methods solve the symmetric
1041   system of equations to obtain the step diretion dk:
1042               Hk dk = -gk
1043   for free variables only. The step can be globalized either through
1044   trust-region methods, or a line search, or a heuristic mixture of both.
1045 
1046     Options Database Keys:
1047 + -tao_bnk_pc_type - "none","ahess","bfgs","petsc"
1048 . -tao_bnk_bfgs_scale_type - "ahess","phess","bfgs"
1049 . -tao_bnk_init_type - "constant","direction","interpolation"
1050 . -tao_bnk_update_type - "step","direction","interpolation"
1051 . -tao_bnk_as_type - "none","bertsekas"
1052 . -tao_bnk_sval - perturbation starting value
1053 . -tao_bnk_imin - minimum initial perturbation
1054 . -tao_bnk_imax - maximum initial perturbation
1055 . -tao_bnk_pmin - minimum perturbation
1056 . -tao_bnk_pmax - maximum perturbation
1057 . -tao_bnk_pgfac - growth factor
1058 . -tao_bnk_psfac - shrink factor
1059 . -tao_bnk_imfac - initial merit factor
1060 . -tao_bnk_pmgfac - merit growth factor
1061 . -tao_bnk_pmsfac - merit shrink factor
1062 . -tao_bnk_eta1 - poor steplength; reduce radius
1063 . -tao_bnk_eta2 - reasonable steplength; leave radius
1064 . -tao_bnk_eta3 - good steplength; increase readius
1065 . -tao_bnk_eta4 - excellent steplength; greatly increase radius
1066 . -tao_bnk_alpha1 - alpha1 reduction
1067 . -tao_bnk_alpha2 - alpha2 reduction
1068 . -tao_bnk_alpha3 - alpha3 reduction
1069 . -tao_bnk_alpha4 - alpha4 reduction
1070 . -tao_bnk_alpha - alpha5 reduction
1071 . -tao_bnk_mu1 - mu1 interpolation update
1072 . -tao_bnk_mu2 - mu2 interpolation update
1073 . -tao_bnk_gamma1 - gamma1 interpolation update
1074 . -tao_bnk_gamma2 - gamma2 interpolation update
1075 . -tao_bnk_gamma3 - gamma3 interpolation update
1076 . -tao_bnk_gamma4 - gamma4 interpolation update
1077 . -tao_bnk_theta - theta interpolation update
1078 . -tao_bnk_omega1 - omega1 step update
1079 . -tao_bnk_omega2 - omega2 step update
1080 . -tao_bnk_omega3 - omega3 step update
1081 . -tao_bnk_omega4 - omega4 step update
1082 . -tao_bnk_omega5 - omega5 step update
1083 . -tao_bnk_mu1_i -  mu1 interpolation init factor
1084 . -tao_bnk_mu2_i -  mu2 interpolation init factor
1085 . -tao_bnk_gamma1_i -  gamma1 interpolation init factor
1086 . -tao_bnk_gamma2_i -  gamma2 interpolation init factor
1087 . -tao_bnk_gamma3_i -  gamma3 interpolation init factor
1088 . -tao_bnk_gamma4_i -  gamma4 interpolation init factor
1089 . -tao_bnk_theta_i -  theta interpolation init factor
1090 - -tao_bnk_bound_tol -  initial tolerance used in estimating bounded active variables
1091 
1092   Level: beginner
1093 M*/
1094 
1095 PetscErrorCode TaoCreate_BNK(Tao tao)
1096 {
1097   TAO_BNK        *bnk;
1098   const char     *morethuente_type = TAOLINESEARCHMT;
1099   PetscErrorCode ierr;
1100 
1101   PetscFunctionBegin;
1102   ierr = PetscNewLog(tao,&bnk);CHKERRQ(ierr);
1103 
1104   tao->ops->setup = TaoSetUp_BNK;
1105   tao->ops->view = TaoView_BNK;
1106   tao->ops->setfromoptions = TaoSetFromOptions_BNK;
1107   tao->ops->destroy = TaoDestroy_BNK;
1108 
1109   /*  Override default settings (unless already changed) */
1110   if (!tao->max_it_changed) tao->max_it = 50;
1111   if (!tao->trust0_changed) tao->trust0 = 100.0;
1112 
1113   tao->data = (void*)bnk;
1114 
1115   /*  Hessian shifting parameters */
1116   bnk->sval   = 0.0;
1117   bnk->imin   = 1.0e-4;
1118   bnk->imax   = 1.0e+2;
1119   bnk->imfac  = 1.0e-1;
1120 
1121   bnk->pmin   = 1.0e-12;
1122   bnk->pmax   = 1.0e+2;
1123   bnk->pgfac  = 1.0e+1;
1124   bnk->psfac  = 4.0e-1;
1125   bnk->pmgfac = 1.0e-1;
1126   bnk->pmsfac = 1.0e-1;
1127 
1128   /*  Default values for trust-region radius update based on steplength */
1129   bnk->nu1 = 0.25;
1130   bnk->nu2 = 0.50;
1131   bnk->nu3 = 1.00;
1132   bnk->nu4 = 1.25;
1133 
1134   bnk->omega1 = 0.25;
1135   bnk->omega2 = 0.50;
1136   bnk->omega3 = 1.00;
1137   bnk->omega4 = 2.00;
1138   bnk->omega5 = 4.00;
1139 
1140   /*  Default values for trust-region radius update based on reduction */
1141   bnk->eta1 = 1.0e-4;
1142   bnk->eta2 = 0.25;
1143   bnk->eta3 = 0.50;
1144   bnk->eta4 = 0.90;
1145 
1146   bnk->alpha1 = 0.25;
1147   bnk->alpha2 = 0.50;
1148   bnk->alpha3 = 1.00;
1149   bnk->alpha4 = 2.00;
1150   bnk->alpha5 = 4.00;
1151 
1152   /*  Default values for trust-region radius update based on interpolation */
1153   bnk->mu1 = 0.10;
1154   bnk->mu2 = 0.50;
1155 
1156   bnk->gamma1 = 0.25;
1157   bnk->gamma2 = 0.50;
1158   bnk->gamma3 = 2.00;
1159   bnk->gamma4 = 4.00;
1160 
1161   bnk->theta = 0.05;
1162 
1163   /*  Default values for trust region initialization based on interpolation */
1164   bnk->mu1_i = 0.35;
1165   bnk->mu2_i = 0.50;
1166 
1167   bnk->gamma1_i = 0.0625;
1168   bnk->gamma2_i = 0.5;
1169   bnk->gamma3_i = 2.0;
1170   bnk->gamma4_i = 5.0;
1171 
1172   bnk->theta_i = 0.25;
1173 
1174   /*  Remaining parameters */
1175   bnk->min_radius = 1.0e-10;
1176   bnk->max_radius = 1.0e10;
1177   bnk->epsilon = PetscPowReal(PETSC_MACHINE_EPSILON, 2.0/3.0);
1178   bnk->as_tol = 1.0e-3;
1179   bnk->as_step = 1.0e-3;
1180 
1181   bnk->pc_type         = BNK_PC_BFGS;
1182   bnk->bfgs_scale_type = BFGS_SCALE_PHESS;
1183   bnk->init_type       = BNK_INIT_INTERPOLATION;
1184   bnk->update_type     = BNK_UPDATE_INTERPOLATION;
1185   bnk->as_type         = BNK_AS_BERTSEKAS;
1186 
1187   ierr = TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch);CHKERRQ(ierr);
1188   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr);
1189   ierr = TaoLineSearchSetType(tao->linesearch,morethuente_type);CHKERRQ(ierr);
1190   ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr);
1191   ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr);
1192 
1193   /*  Set linear solver to default for symmetric matrices */
1194   ierr = KSPCreate(((PetscObject)tao)->comm,&tao->ksp);CHKERRQ(ierr);
1195   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1);CHKERRQ(ierr);
1196   ierr = KSPSetOptionsPrefix(tao->ksp,tao->hdr.prefix);CHKERRQ(ierr);
1197   ierr = KSPSetType(tao->ksp,KSPCGSTCG);CHKERRQ(ierr);
1198   PetscFunctionReturn(0);
1199 }
1200