xref: /petsc/src/tao/bound/impls/bnk/bnk.c (revision 09164190726e6f6172ca3d0df2bdee7b8aa3e300)
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 PetscErrorCode TaoBNKInitialize(Tao tao)
23 {
24   PetscErrorCode               ierr;
25   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
26   KSPType                      ksp_type;
27   PC                           pc;
28 
29   PetscReal                    fmin, ftrial, prered, actred, kappa, sigma;
30   PetscReal                    tau, tau_1, tau_2, tau_max, tau_min, max_radius;
31   PetscReal                    delta, step = 1.0;
32 
33   PetscInt                     n,N,needH = 1;
34 
35   PetscInt                     i_max = 5;
36   PetscInt                     j_max = 1;
37   PetscInt                     i, j;
38 
39   PetscFunctionBegin;
40   /* Number of times ksp stopped because of these reasons */
41   bnk->ksp_atol = 0;
42   bnk->ksp_rtol = 0;
43   bnk->ksp_dtol = 0;
44   bnk->ksp_ctol = 0;
45   bnk->ksp_negc = 0;
46   bnk->ksp_iter = 0;
47   bnk->ksp_othr = 0;
48 
49   /* Initialize trust-region radius when using nash, stcg, or gltr
50      Command automatically ignored for other methods
51      Will be reset during the first iteration
52   */
53   ierr = KSPGetType(tao->ksp,&ksp_type);CHKERRQ(ierr);
54   ierr = PetscStrcmp(ksp_type,KSPCGNASH,&bnk->is_nash);CHKERRQ(ierr);
55   ierr = PetscStrcmp(ksp_type,KSPCGSTCG,&bnk->is_stcg);CHKERRQ(ierr);
56   ierr = PetscStrcmp(ksp_type,KSPCGGLTR,&bnk->is_gltr);CHKERRQ(ierr);
57 
58   ierr = KSPCGSetRadius(tao->ksp,bnk->max_radius);CHKERRQ(ierr);
59 
60   if (bnk->is_nash || bnk->is_stcg || bnk->is_gltr) {
61     if (tao->trust0 < 0.0) SETERRQ(PETSC_COMM_SELF,1,"Initial radius negative");
62     tao->trust = tao->trust0;
63     tao->trust = PetscMax(tao->trust, bnk->min_radius);
64     tao->trust = PetscMin(tao->trust, bnk->max_radius);
65   }
66 
67   /* Get vectors we will need */
68   if (BNK_PC_BFGS == bnk->pc_type && !bnk->M) {
69     ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr);
70     ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr);
71     ierr = MatCreateLMVM(((PetscObject)tao)->comm,n,N,&bnk->M);CHKERRQ(ierr);
72     ierr = MatLMVMAllocateVectors(bnk->M,tao->solution);CHKERRQ(ierr);
73   }
74 
75   /* create vectors for the limited memory preconditioner */
76   if ((BNK_PC_BFGS == bnk->pc_type) && (BFGS_SCALE_BFGS != bnk->bfgs_scale_type)) {
77     if (!bnk->Diag) {
78       ierr = VecDuplicate(tao->solution,&bnk->Diag);CHKERRQ(ierr);
79     }
80   }
81 
82   /* Modify the preconditioner to use the bfgs approximation */
83   ierr = KSPGetPC(tao->ksp, &pc);CHKERRQ(ierr);
84   switch(bnk->pc_type) {
85   case BNK_PC_NONE:
86     ierr = PCSetType(pc, PCNONE);CHKERRQ(ierr);
87     ierr = PCSetFromOptions(pc);CHKERRQ(ierr);
88     break;
89 
90   case BNK_PC_AHESS:
91     ierr = PCSetType(pc, PCJACOBI);CHKERRQ(ierr);
92     ierr = PCSetFromOptions(pc);CHKERRQ(ierr);
93     ierr = PCJacobiSetUseAbs(pc,PETSC_TRUE);CHKERRQ(ierr);
94     break;
95 
96   case BNK_PC_BFGS:
97     ierr = PCSetType(pc, PCSHELL);CHKERRQ(ierr);
98     ierr = PCSetFromOptions(pc);CHKERRQ(ierr);
99     ierr = PCShellSetName(pc, "bfgs");CHKERRQ(ierr);
100     ierr = PCShellSetContext(pc, bnk->M);CHKERRQ(ierr);
101     ierr = PCShellSetApply(pc, MatLMVMSolveShell);CHKERRQ(ierr);
102     break;
103 
104   default:
105     /* Use the pc method set by pc_type */
106     break;
107   }
108 
109   /* Initialize trust-region radius.  The initialization is only performed
110      when we are using Nash, Steihaug-Toint or the Generalized Lanczos method. */
111   if (bnk->is_nash || bnk->is_stcg || bnk->is_gltr) {
112     switch(bnk->init_type) {
113     case BNK_INIT_CONSTANT:
114       /* Use the initial radius specified */
115       break;
116 
117     case BNK_INIT_INTERPOLATION:
118       /* Use the initial radius specified */
119       max_radius = 0.0;
120 
121       for (j = 0; j < j_max; ++j) {
122         fmin = bnk->f;
123         sigma = 0.0;
124 
125         if (needH) {
126           ierr  = TaoComputeHessian(tao, tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr);
127           needH = 0;
128         }
129 
130         for (i = 0; i < i_max; ++i) {
131           ierr = VecCopy(tao->solution,bnk->W);CHKERRQ(ierr);
132           ierr = VecAXPY(bnk->W,-tao->trust/bnk->gnorm,tao->gradient);CHKERRQ(ierr);
133           ierr = TaoComputeObjective(tao, bnk->W, &ftrial);CHKERRQ(ierr);
134           if (PetscIsInfOrNanReal(ftrial)) {
135             tau = bnk->gamma1_i;
136           } else {
137             if (ftrial < fmin) {
138               fmin = ftrial;
139               sigma = -tao->trust / bnk->gnorm;
140             }
141 
142             ierr = MatMult(tao->hessian, tao->gradient, tao->stepdirection);CHKERRQ(ierr);
143             ierr = VecDot(tao->gradient, tao->stepdirection, &prered);CHKERRQ(ierr);
144 
145             prered = tao->trust * (bnk->gnorm - 0.5 * tao->trust * prered / (bnk->gnorm * bnk->gnorm));
146             actred = bnk->f - ftrial;
147             if ((PetscAbsScalar(actred) <= bnk->epsilon) && (PetscAbsScalar(prered) <= bnk->epsilon)) {
148               kappa = 1.0;
149             } else {
150               kappa = actred / prered;
151             }
152 
153             tau_1 = bnk->theta_i * bnk->gnorm * tao->trust / (bnk->theta_i * bnk->gnorm * tao->trust + (1.0 - bnk->theta_i) * prered - actred);
154             tau_2 = bnk->theta_i * bnk->gnorm * tao->trust / (bnk->theta_i * bnk->gnorm * tao->trust - (1.0 + bnk->theta_i) * prered + actred);
155             tau_min = PetscMin(tau_1, tau_2);
156             tau_max = PetscMax(tau_1, tau_2);
157 
158             if (PetscAbsScalar(kappa - 1.0) <= bnk->mu1_i) {
159               /* Great agreement */
160               max_radius = PetscMax(max_radius, tao->trust);
161 
162               if (tau_max < 1.0) {
163                 tau = bnk->gamma3_i;
164               } else if (tau_max > bnk->gamma4_i) {
165                 tau = bnk->gamma4_i;
166               } else if (tau_1 >= 1.0 && tau_1 <= bnk->gamma4_i && tau_2 < 1.0) {
167                 tau = tau_1;
168               } else if (tau_2 >= 1.0 && tau_2 <= bnk->gamma4_i && tau_1 < 1.0) {
169                 tau = tau_2;
170               } else {
171                 tau = tau_max;
172               }
173             } else if (PetscAbsScalar(kappa - 1.0) <= bnk->mu2_i) {
174               /* Good agreement */
175               max_radius = PetscMax(max_radius, tao->trust);
176 
177               if (tau_max < bnk->gamma2_i) {
178                 tau = bnk->gamma2_i;
179               } else if (tau_max > bnk->gamma3_i) {
180                 tau = bnk->gamma3_i;
181               } else {
182                 tau = tau_max;
183               }
184             } else {
185               /* Not good agreement */
186               if (tau_min > 1.0) {
187                 tau = bnk->gamma2_i;
188               } else if (tau_max < bnk->gamma1_i) {
189                 tau = bnk->gamma1_i;
190               } else if ((tau_min < bnk->gamma1_i) && (tau_max >= 1.0)) {
191                 tau = bnk->gamma1_i;
192               } else if ((tau_1 >= bnk->gamma1_i) && (tau_1 < 1.0) && ((tau_2 < bnk->gamma1_i) || (tau_2 >= 1.0))) {
193                 tau = tau_1;
194               } else if ((tau_2 >= bnk->gamma1_i) && (tau_2 < 1.0) && ((tau_1 < bnk->gamma1_i) || (tau_2 >= 1.0))) {
195                 tau = tau_2;
196               } else {
197                 tau = tau_max;
198               }
199             }
200           }
201           tao->trust = tau * tao->trust;
202         }
203 
204         if (fmin < bnk->f) {
205           bnk->f = fmin;
206           ierr = VecAXPY(tao->solution,sigma,tao->gradient);CHKERRQ(ierr);
207           ierr = TaoComputeGradient(tao,tao->solution,bnk->unprojected_gradient);CHKERRQ(ierr);
208           ierr = VecBoundGradientProjection(bnk->unprojected_gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr);
209 
210           ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&bnk->gnorm);CHKERRQ(ierr);
211           if (PetscIsInfOrNanReal(bnk->gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute gradient generated Inf or NaN");
212           needH = 1;
213 
214           ierr = TaoLogConvergenceHistory(tao,bnk->f,bnk->gnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
215           ierr = TaoMonitor(tao,tao->niter,bnk->f,bnk->gnorm,0.0,step);CHKERRQ(ierr);
216           ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
217           if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
218         }
219       }
220       tao->trust = PetscMax(tao->trust, max_radius);
221 
222       /* Modify the radius if it is too large or small */
223       tao->trust = PetscMax(tao->trust, bnk->min_radius);
224       tao->trust = PetscMin(tao->trust, bnk->max_radius);
225       break;
226 
227     default:
228       /* Norm of the first direction will initialize radius */
229       tao->trust = 0.0;
230       break;
231     }
232   }
233 
234   /* Set initial scaling for the BFGS preconditioner
235      This step is done after computing the initial trust-region radius
236      since the function value may have decreased */
237   if (BNK_PC_BFGS == bnk->pc_type) {
238     if (bnk->f != 0.0) {
239       delta = 2.0 * PetscAbsScalar(bnk->f) / (bnk->gnorm*bnk->gnorm);
240     } else {
241       delta = 2.0 / (bnk->gnorm*bnk->gnorm);
242     }
243     ierr = MatLMVMSetDelta(bnk->M,delta);CHKERRQ(ierr);
244   }
245 
246   /* Set counter for gradient/reset steps*/
247   bnk->newt = 0;
248   bnk->bfgs = 0;
249   bnk->sgrad = 0;
250   bnk->grad = 0;
251   PetscFunctionReturn(0);
252 }
253 
254 PetscErrorCode TaoBNKComputeStep(Tao tao, PetscInt *stepType)
255 {
256   PetscErrorCode               ierr;
257   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
258   KSPConvergedReason           ksp_reason;
259 
260   PetscReal                    gdx, delta, e_min;
261 
262   PetscInt                     bfgsUpdates = 0;
263   PetscInt                     kspits;
264 
265   PetscFunctionBegin;
266   if ((BNK_PC_BFGS == bnk->pc_type) && (BFGS_SCALE_AHESS == bnk->bfgs_scale_type)) {
267     /* Obtain diagonal for the bfgs preconditioner  */
268     ierr = MatGetDiagonal(tao->hessian, bnk->Diag);CHKERRQ(ierr);
269     ierr = VecAbs(bnk->Diag);CHKERRQ(ierr);
270     ierr = VecReciprocal(bnk->Diag);CHKERRQ(ierr);
271     ierr = MatLMVMSetScale(bnk->M,bnk->Diag);CHKERRQ(ierr);
272   }
273 
274   /* Shift the Hessian matrix */
275   bnk->pert = bnk->sval;
276   if (bnk->pert > 0) {
277     ierr = MatShift(tao->hessian, bnk->pert);CHKERRQ(ierr);
278     if (tao->hessian != tao->hessian_pre) {
279       ierr = MatShift(tao->hessian_pre, bnk->pert);CHKERRQ(ierr);
280     }
281   }
282 
283   if (BNK_PC_BFGS == bnk->pc_type) {
284     if (BFGS_SCALE_PHESS == bnk->bfgs_scale_type) {
285       /* Obtain diagonal for the bfgs preconditioner  */
286       ierr = MatGetDiagonal(tao->hessian, bnk->Diag);CHKERRQ(ierr);
287       ierr = VecAbs(bnk->Diag);CHKERRQ(ierr);
288       ierr = VecReciprocal(bnk->Diag);CHKERRQ(ierr);
289       ierr = MatLMVMSetScale(bnk->M,bnk->Diag);CHKERRQ(ierr);
290     }
291     /* Update the limited memory preconditioner and get existing # of updates */
292     ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
293     ierr = MatLMVMGetUpdates(bnk->M, &bfgsUpdates);CHKERRQ(ierr);
294   }
295 
296   /* Determine the inactive set */
297   ierr = ISDestroy(&bnk->inactive_idx);CHKERRQ(ierr);
298   ierr = VecWhichInactive(tao->XL,tao->solution,bnk->unprojected_gradient,tao->XU,PETSC_TRUE,&bnk->inactive_idx);CHKERRQ(ierr);
299 
300   /* Extract the inactive sub-system to solve */
301   ierr = MatLMVMSetInactive(bnk->M, bnk->inactive_idx);CHKERRQ(ierr);
302   ierr = VecSet(tao->stepdirection, 0.0);CHKERRQ(ierr);
303   ierr = VecGetSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->inactive_step);CHKERRQ(ierr);
304   ierr = VecGetSubVector(tao->gradient, bnk->inactive_idx, &bnk->inactive_gradient);CHKERRQ(ierr);
305   ierr = MatCreateSubMatrix(tao->hessian, bnk->inactive_idx, bnk->inactive_idx, MAT_INITIAL_MATRIX, &bnk->H_inactive);CHKERRQ(ierr);
306   ierr = MatCreateSubMatrix(tao->hessian_pre, bnk->inactive_idx, bnk->inactive_idx, MAT_INITIAL_MATRIX, &bnk->Hpre_inactive);CHKERRQ(ierr);
307 
308   /* Solve the Newton system of equations */
309   ierr = KSPSetOperators(tao->ksp,bnk->H_inactive,bnk->Hpre_inactive);CHKERRQ(ierr);
310   if (bnk->is_nash || bnk->is_stcg || bnk->is_gltr) {
311     ierr = KSPCGSetRadius(tao->ksp,bnk->max_radius);CHKERRQ(ierr);
312     ierr = KSPSolve(tao->ksp, bnk->inactive_gradient, bnk->inactive_step);CHKERRQ(ierr);
313     ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr);
314     tao->ksp_its+=kspits;
315     tao->ksp_tot_its+=kspits;
316     ierr = KSPCGGetNormD(tao->ksp,&bnk->dnorm);CHKERRQ(ierr);
317 
318     if (0.0 == tao->trust) {
319       /* Radius was uninitialized; use the norm of the direction */
320       if (bnk->dnorm > 0.0) {
321         tao->trust = bnk->dnorm;
322 
323         /* Modify the radius if it is too large or small */
324         tao->trust = PetscMax(tao->trust, bnk->min_radius);
325         tao->trust = PetscMin(tao->trust, bnk->max_radius);
326       } else {
327         /* The direction was bad; set radius to default value and re-solve
328            the trust-region subproblem to get a direction */
329         tao->trust = tao->trust0;
330 
331         /* Modify the radius if it is too large or small */
332         tao->trust = PetscMax(tao->trust, bnk->min_radius);
333         tao->trust = PetscMin(tao->trust, bnk->max_radius);
334 
335         ierr = KSPCGSetRadius(tao->ksp,bnk->max_radius);CHKERRQ(ierr);
336         ierr = KSPSolve(tao->ksp, bnk->inactive_gradient, bnk->inactive_step);CHKERRQ(ierr);
337         ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr);
338         tao->ksp_its+=kspits;
339         tao->ksp_tot_its+=kspits;
340         ierr = KSPCGGetNormD(tao->ksp,&bnk->dnorm);CHKERRQ(ierr);
341 
342         if (bnk->dnorm == 0.0) SETERRQ(PETSC_COMM_SELF,1, "Initial direction zero");
343       }
344     }
345   } else {
346     ierr = KSPSolve(tao->ksp, bnk->inactive_gradient, bnk->inactive_step);CHKERRQ(ierr);
347     ierr = KSPGetIterationNumber(tao->ksp, &kspits);CHKERRQ(ierr);
348     tao->ksp_its += kspits;
349     tao->ksp_tot_its+=kspits;
350   }
351   /* Restore submatrices and subvecs */
352   ierr = MatDestroy(&bnk->H_inactive);CHKERRQ(ierr);
353   ierr = MatDestroy(&bnk->Hpre_inactive);CHKERRQ(ierr);
354   ierr = VecRestoreSubVector(tao->gradient, bnk->inactive_idx, &bnk->inactive_gradient);CHKERRQ(ierr);
355   ierr = VecRestoreSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->inactive_step);CHKERRQ(ierr);
356 
357   ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
358   ierr = KSPGetConvergedReason(tao->ksp, &ksp_reason);CHKERRQ(ierr);
359   if ((KSP_DIVERGED_INDEFINITE_PC == ksp_reason) &&  (BNK_PC_BFGS == bnk->pc_type) && (bfgsUpdates > 1)) {
360     /* Preconditioner is numerically indefinite; reset the
361        approximate if using BFGS preconditioning. */
362 
363     if (bnk->f != 0.0) {
364       delta = 2.0 * PetscAbsScalar(bnk->f) / (bnk->gnorm*bnk->gnorm);
365     } else {
366       delta = 2.0 / (bnk->gnorm*bnk->gnorm);
367     }
368     ierr = MatLMVMSetDelta(bnk->M,delta);CHKERRQ(ierr);
369     ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr);
370     ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
371     bfgsUpdates = 1;
372   }
373 
374   if (KSP_CONVERGED_ATOL == ksp_reason) {
375     ++bnk->ksp_atol;
376   } else if (KSP_CONVERGED_RTOL == ksp_reason) {
377     ++bnk->ksp_rtol;
378   } else if (KSP_CONVERGED_CG_CONSTRAINED == ksp_reason) {
379     ++bnk->ksp_ctol;
380   } else if (KSP_CONVERGED_CG_NEG_CURVE == ksp_reason) {
381     ++bnk->ksp_negc;
382   } else if (KSP_DIVERGED_DTOL == ksp_reason) {
383     ++bnk->ksp_dtol;
384   } else if (KSP_DIVERGED_ITS == ksp_reason) {
385     ++bnk->ksp_iter;
386   } else {
387     ++bnk->ksp_othr;
388   }
389 
390   /* Check for success (descent direction) */
391   ierr = VecDot(tao->stepdirection, tao->gradient, &gdx);CHKERRQ(ierr);
392   if ((gdx >= 0.0) || PetscIsInfOrNanReal(gdx)) {
393     /* Newton step is not descent or direction produced Inf or NaN
394        Update the perturbation for next time */
395     if (bnk->pert <= 0.0) {
396       /* Initialize the perturbation */
397       bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
398       if (bnk->is_gltr) {
399         ierr = KSPCGGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr);
400         bnk->pert = PetscMax(bnk->pert, -e_min);
401       }
402     } else {
403       /* Increase the perturbation */
404       bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
405     }
406 
407     if (BNK_PC_BFGS != bnk->pc_type) {
408       /* We don't have the bfgs matrix around and updated
409          Must use gradient direction in this case */
410       ierr = VecCopy(tao->gradient, tao->stepdirection);CHKERRQ(ierr);
411       ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
412       ++bnk->grad;
413       *stepType = BNK_GRADIENT;
414     } else {
415       /* Attempt to use the BFGS direction */
416       ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
417       ierr = VecBoundGradientProjection(tao->stepdirection,tao->solution,tao->XL,tao->XU,tao->stepdirection);CHKERRQ(ierr);
418       ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
419 
420       /* Check for success (descent direction) */
421       ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr);
422       if ((gdx >= 0) || PetscIsInfOrNanReal(gdx)) {
423         /* BFGS direction is not descent or direction produced not a number
424            We can assert bfgsUpdates > 1 in this case because
425            the first solve produces the scaled gradient direction,
426            which is guaranteed to be descent */
427 
428         /* Use steepest descent direction (scaled) */
429 
430         if (bnk->f != 0.0) {
431           delta = 2.0 * PetscAbsScalar(bnk->f) / (bnk->gnorm*bnk->gnorm);
432         } else {
433           delta = 2.0 / (bnk->gnorm*bnk->gnorm);
434         }
435         ierr = MatLMVMSetDelta(bnk->M, delta);CHKERRQ(ierr);
436         ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr);
437         ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
438         ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
439         ierr = VecBoundGradientProjection(tao->stepdirection,tao->solution,tao->XL,tao->XU,tao->stepdirection);CHKERRQ(ierr);
440         ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
441 
442         bfgsUpdates = 1;
443         ++bnk->sgrad;
444         *stepType = BNK_SCALED_GRADIENT;
445       } else {
446         if (1 == bfgsUpdates) {
447           /* The first BFGS direction is always the scaled gradient */
448           ++bnk->sgrad;
449           *stepType = BNK_SCALED_GRADIENT;
450         } else {
451           ++bnk->bfgs;
452           *stepType = BNK_BFGS;
453         }
454       }
455     }
456   } else {
457     /* Computed Newton step is descent */
458     switch (ksp_reason) {
459     case KSP_DIVERGED_NANORINF:
460     case KSP_DIVERGED_BREAKDOWN:
461     case KSP_DIVERGED_INDEFINITE_MAT:
462     case KSP_DIVERGED_INDEFINITE_PC:
463     case KSP_CONVERGED_CG_NEG_CURVE:
464       /* Matrix or preconditioner is indefinite; increase perturbation */
465       if (bnk->pert <= 0.0) {
466         /* Initialize the perturbation */
467         bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
468         if (bnk->is_gltr) {
469           ierr = KSPCGGLTRGetMinEig(tao->ksp, &e_min);CHKERRQ(ierr);
470           bnk->pert = PetscMax(bnk->pert, -e_min);
471         }
472       } else {
473         /* Increase the perturbation */
474         bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
475       }
476       break;
477 
478     default:
479       /* Newton step computation is good; decrease perturbation */
480       bnk->pert = PetscMin(bnk->psfac * bnk->pert, bnk->pmsfac * bnk->gnorm);
481       if (bnk->pert < bnk->pmin) {
482         bnk->pert = 0.0;
483       }
484       break;
485     }
486 
487     ++bnk->newt;
488     stepType = BNK_NEWTON;
489   }
490   PetscFunctionReturn(0);
491 }
492 
493 PetscErrorCode TaoBNKUpdateTrustRadius(Tao tao, PetscReal fold, PetscReal fnew, PetscInt stepType, PetscBool *accept)
494 {
495   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
496   PetscErrorCode ierr;
497 
498   PetscReal      step, prered, actred, kappa;
499   PetscReal      gdx, tau_1, tau_2, tau_min, tau_max;
500 
501   PetscFunctionBegin;
502   /* Update trust region radius */
503   *accept = PETSC_FALSE;
504   switch(bnk->update_type) {
505   case BNK_UPDATE_STEP:
506     if (stepType == BNK_NEWTON) {
507       ierr = TaoLineSearchGetStepLength(tao->linesearch, &step);CHKERRQ(ierr);
508       if (step < bnk->nu1) {
509         /* Very bad step taken; reduce radius */
510         tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust);
511       } else if (step < bnk->nu2) {
512         /* Reasonably bad step taken; reduce radius */
513         tao->trust = bnk->omega2 * PetscMin(bnk->dnorm, tao->trust);
514       } else if (step < bnk->nu3) {
515         /*  Reasonable step was taken; leave radius alone */
516         if (bnk->omega3 < 1.0) {
517           tao->trust = bnk->omega3 * PetscMin(bnk->dnorm, tao->trust);
518         } else if (bnk->omega3 > 1.0) {
519           tao->trust = PetscMax(bnk->omega3 * bnk->dnorm, tao->trust);
520         }
521       } else if (step < bnk->nu4) {
522         /*  Full step taken; increase the radius */
523         tao->trust = PetscMax(bnk->omega4 * bnk->dnorm, tao->trust);
524       } else {
525         /*  More than full step taken; increase the radius */
526         tao->trust = PetscMax(bnk->omega5 * bnk->dnorm, tao->trust);
527       }
528     } else {
529       /*  Newton step was not good; reduce the radius */
530       tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust);
531     }
532     break;
533 
534   case BNK_UPDATE_REDUCTION:
535     if (stepType == BNK_NEWTON) {
536       /*  Get predicted reduction */
537       ierr = KSPCGGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr);
538       if (prered >= 0.0) {
539         /*  The predicted reduction has the wrong sign.  This cannot */
540         /*  happen in infinite precision arithmetic.  Step should */
541         /*  be rejected! */
542         tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
543       } else {
544         if (PetscIsInfOrNanReal(fnew)) {
545           tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
546         } else {
547           /*  Compute and actual reduction */
548           actred = fold - fnew;
549           prered = -prered;
550           if ((PetscAbsScalar(actred) <= bnk->epsilon) &&
551               (PetscAbsScalar(prered) <= bnk->epsilon)) {
552             kappa = 1.0;
553           } else {
554             kappa = actred / prered;
555           }
556           /*  Accept of reject the step and update radius */
557           if (kappa < bnk->eta1) {
558             /*  Very bad step, rejected */
559             tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
560           } else {
561             /* Accept step here */
562             *accept = PETSC_TRUE;
563             /* Adjust trust radius */
564             if (kappa < bnk->eta2) {
565               /*  Marginal bad step */
566               tao->trust = bnk->alpha2 * PetscMin(tao->trust, bnk->dnorm);
567             } else if (kappa < bnk->eta3) {
568               /*  Reasonable step */
569               if (bnk->alpha3 < 1.0) {
570                 tao->trust = bnk->alpha3 * PetscMin(bnk->dnorm, tao->trust);
571               } else if (bnk->alpha3 > 1.0) {
572                 tao->trust = PetscMax(bnk->alpha3 * bnk->dnorm, tao->trust);
573               }
574             } else if (kappa < bnk->eta4) {
575               /*  Good step */
576               tao->trust = PetscMax(bnk->alpha4 * bnk->dnorm, tao->trust);
577             } else {
578               /*  Very good step */
579               tao->trust = PetscMax(bnk->alpha5 * bnk->dnorm, tao->trust);
580             }
581           }
582         }
583       }
584     } else {
585       /*  Newton step was not good; reduce the radius */
586       tao->trust = bnk->alpha1 * PetscMin(bnk->dnorm, tao->trust);
587     }
588     break;
589 
590   default:
591     if (stepType == BNK_NEWTON) {
592       ierr = KSPCGGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr);
593       if (prered >= 0.0) {
594         /*  The predicted reduction has the wrong sign.  This cannot */
595         /*  happen in infinite precision arithmetic.  Step should */
596         /*  be rejected! */
597         tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
598       } else {
599         if (PetscIsInfOrNanReal(fnew)) {
600           tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
601         } else {
602           actred = fold - fnew;
603           prered = -prered;
604           if ((PetscAbsScalar(actred) <= bnk->epsilon) && (PetscAbsScalar(prered) <= bnk->epsilon)) {
605             kappa = 1.0;
606           } else {
607             kappa = actred / prered;
608           }
609 
610           ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr);
611           tau_1 = bnk->theta * gdx / (bnk->theta * gdx - (1.0 - bnk->theta) * prered + actred);
612           tau_2 = bnk->theta * gdx / (bnk->theta * gdx + (1.0 + bnk->theta) * prered - actred);
613           tau_min = PetscMin(tau_1, tau_2);
614           tau_max = PetscMax(tau_1, tau_2);
615 
616           if (kappa >= 1.0 - bnk->mu1) {
617             /*  Great agreement */
618             *accept = PETSC_TRUE;
619             if (tau_max < 1.0) {
620               tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm);
621             } else if (tau_max > bnk->gamma4) {
622               tao->trust = PetscMax(tao->trust, bnk->gamma4 * bnk->dnorm);
623             } else {
624               tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm);
625             }
626           } else if (kappa >= 1.0 - bnk->mu2) {
627             /*  Good agreement */
628             *accept = PETSC_TRUE;
629             if (tau_max < bnk->gamma2) {
630               tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm);
631             } else if (tau_max > bnk->gamma3) {
632               tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm);
633             } else if (tau_max < 1.0) {
634               tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm);
635             } else {
636               tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm);
637             }
638           } else {
639             /*  Not good agreement */
640             if (tau_min > 1.0) {
641               tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm);
642             } else if (tau_max < bnk->gamma1) {
643               tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
644             } else if ((tau_min < bnk->gamma1) && (tau_max >= 1.0)) {
645               tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
646             } else if ((tau_1 >= bnk->gamma1) && (tau_1 < 1.0) && ((tau_2 < bnk->gamma1) || (tau_2 >= 1.0))) {
647               tao->trust = tau_1 * PetscMin(tao->trust, bnk->dnorm);
648             } else if ((tau_2 >= bnk->gamma1) && (tau_2 < 1.0) && ((tau_1 < bnk->gamma1) || (tau_2 >= 1.0))) {
649               tao->trust = tau_2 * PetscMin(tao->trust, bnk->dnorm);
650             } else {
651               tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm);
652             }
653           }
654         }
655       }
656     } else {
657       /*  Newton step was not good; reduce the radius */
658       tao->trust = bnk->gamma1 * PetscMin(bnk->dnorm, tao->trust);
659     }
660     /*  The radius may have been increased; modify if it is too large */
661     tao->trust = PetscMin(tao->trust, bnk->max_radius);
662   }
663   PetscFunctionReturn(0);
664 }
665 
666 /* ---------------------------------------------------------- */
667 static PetscErrorCode TaoSetUp_BNK(Tao tao)
668 {
669   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
670   PetscErrorCode ierr;
671 
672   PetscFunctionBegin;
673   if (!tao->gradient) {ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);}
674   if (!tao->stepdirection) {ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr);}
675   if (!bnk->W) {ierr = VecDuplicate(tao->solution,&bnk->W);CHKERRQ(ierr);}
676   if (!bnk->Xold) {ierr = VecDuplicate(tao->solution,&bnk->Xold);CHKERRQ(ierr);}
677   if (!bnk->Gold) {ierr = VecDuplicate(tao->solution,&bnk->Gold);CHKERRQ(ierr);}
678   if (!bnk->Xwork) {ierr = VecDuplicate(tao->solution,&bnk->Xwork);CHKERRQ(ierr);}
679   if (!bnk->Gwork) {ierr = VecDuplicate(tao->solution,&bnk->Gwork);CHKERRQ(ierr);}
680   if (!bnk->unprojected_gradient) {ierr = VecDuplicate(tao->solution,&bnk->unprojected_gradient);CHKERRQ(ierr);}
681   if (!bnk->unprojected_gradient_old) {ierr = VecDuplicate(tao->solution,&bnk->unprojected_gradient_old);CHKERRQ(ierr);}
682   bnk->Diag = 0;
683   bnk->M = 0;
684   bnk->H_inactive = 0;
685   bnk->Hpre_inactive = 0;
686   PetscFunctionReturn(0);
687 }
688 
689 /*------------------------------------------------------------*/
690 static PetscErrorCode TaoDestroy_BNK(Tao tao)
691 {
692   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
693   PetscErrorCode ierr;
694 
695   PetscFunctionBegin;
696   if (tao->setupcalled) {
697     ierr = VecDestroy(&bnk->W);CHKERRQ(ierr);
698     ierr = VecDestroy(&bnk->Xold);CHKERRQ(ierr);
699     ierr = VecDestroy(&bnk->Gold);CHKERRQ(ierr);
700     ierr = VecDestroy(&bnk->Xwork);CHKERRQ(ierr);
701     ierr = VecDestroy(&bnk->Gwork);CHKERRQ(ierr);
702     ierr = VecDestroy(&bnk->unprojected_gradient);CHKERRQ(ierr);
703     ierr = VecDestroy(&bnk->unprojected_gradient_old);CHKERRQ(ierr);
704   }
705   ierr = VecDestroy(&bnk->Diag);CHKERRQ(ierr);
706   ierr = MatDestroy(&bnk->M);CHKERRQ(ierr);
707   ierr = MatDestroy(&bnk->H_inactive);CHKERRQ(ierr);
708   ierr = MatDestroy(&bnk->Hpre_inactive);CHKERRQ(ierr);
709   ierr = PetscFree(tao->data);CHKERRQ(ierr);
710   PetscFunctionReturn(0);
711 }
712 
713 /*------------------------------------------------------------*/
714 static PetscErrorCode TaoSetFromOptions_BNK(PetscOptionItems *PetscOptionsObject,Tao tao)
715 {
716   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
717   PetscErrorCode ierr;
718 
719   PetscFunctionBegin;
720   ierr = PetscOptionsHead(PetscOptionsObject,"Newton line search method for unconstrained optimization");CHKERRQ(ierr);
721   ierr = PetscOptionsEList("-tao_BNK_pc_type", "pc type", "", BNK_PC, BNK_PC_TYPES, BNK_PC[bnk->pc_type], &bnk->pc_type, 0);CHKERRQ(ierr);
722   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);
723   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);
724   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);
725   ierr = PetscOptionsReal("-tao_BNK_sval", "perturbation starting value", "", bnk->sval, &bnk->sval,NULL);CHKERRQ(ierr);
726   ierr = PetscOptionsReal("-tao_BNK_imin", "minimum initial perturbation", "", bnk->imin, &bnk->imin,NULL);CHKERRQ(ierr);
727   ierr = PetscOptionsReal("-tao_BNK_imax", "maximum initial perturbation", "", bnk->imax, &bnk->imax,NULL);CHKERRQ(ierr);
728   ierr = PetscOptionsReal("-tao_BNK_imfac", "initial merit factor", "", bnk->imfac, &bnk->imfac,NULL);CHKERRQ(ierr);
729   ierr = PetscOptionsReal("-tao_BNK_pmin", "minimum perturbation", "", bnk->pmin, &bnk->pmin,NULL);CHKERRQ(ierr);
730   ierr = PetscOptionsReal("-tao_BNK_pmax", "maximum perturbation", "", bnk->pmax, &bnk->pmax,NULL);CHKERRQ(ierr);
731   ierr = PetscOptionsReal("-tao_BNK_pgfac", "growth factor", "", bnk->pgfac, &bnk->pgfac,NULL);CHKERRQ(ierr);
732   ierr = PetscOptionsReal("-tao_BNK_psfac", "shrink factor", "", bnk->psfac, &bnk->psfac,NULL);CHKERRQ(ierr);
733   ierr = PetscOptionsReal("-tao_BNK_pmgfac", "merit growth factor", "", bnk->pmgfac, &bnk->pmgfac,NULL);CHKERRQ(ierr);
734   ierr = PetscOptionsReal("-tao_BNK_pmsfac", "merit shrink factor", "", bnk->pmsfac, &bnk->pmsfac,NULL);CHKERRQ(ierr);
735   ierr = PetscOptionsReal("-tao_BNK_eta1", "poor steplength; reduce radius", "", bnk->eta1, &bnk->eta1,NULL);CHKERRQ(ierr);
736   ierr = PetscOptionsReal("-tao_BNK_eta2", "reasonable steplength; leave radius alone", "", bnk->eta2, &bnk->eta2,NULL);CHKERRQ(ierr);
737   ierr = PetscOptionsReal("-tao_BNK_eta3", "good steplength; increase radius", "", bnk->eta3, &bnk->eta3,NULL);CHKERRQ(ierr);
738   ierr = PetscOptionsReal("-tao_BNK_eta4", "excellent steplength; greatly increase radius", "", bnk->eta4, &bnk->eta4,NULL);CHKERRQ(ierr);
739   ierr = PetscOptionsReal("-tao_BNK_alpha1", "", "", bnk->alpha1, &bnk->alpha1,NULL);CHKERRQ(ierr);
740   ierr = PetscOptionsReal("-tao_BNK_alpha2", "", "", bnk->alpha2, &bnk->alpha2,NULL);CHKERRQ(ierr);
741   ierr = PetscOptionsReal("-tao_BNK_alpha3", "", "", bnk->alpha3, &bnk->alpha3,NULL);CHKERRQ(ierr);
742   ierr = PetscOptionsReal("-tao_BNK_alpha4", "", "", bnk->alpha4, &bnk->alpha4,NULL);CHKERRQ(ierr);
743   ierr = PetscOptionsReal("-tao_BNK_alpha5", "", "", bnk->alpha5, &bnk->alpha5,NULL);CHKERRQ(ierr);
744   ierr = PetscOptionsReal("-tao_BNK_nu1", "poor steplength; reduce radius", "", bnk->nu1, &bnk->nu1,NULL);CHKERRQ(ierr);
745   ierr = PetscOptionsReal("-tao_BNK_nu2", "reasonable steplength; leave radius alone", "", bnk->nu2, &bnk->nu2,NULL);CHKERRQ(ierr);
746   ierr = PetscOptionsReal("-tao_BNK_nu3", "good steplength; increase radius", "", bnk->nu3, &bnk->nu3,NULL);CHKERRQ(ierr);
747   ierr = PetscOptionsReal("-tao_BNK_nu4", "excellent steplength; greatly increase radius", "", bnk->nu4, &bnk->nu4,NULL);CHKERRQ(ierr);
748   ierr = PetscOptionsReal("-tao_BNK_omega1", "", "", bnk->omega1, &bnk->omega1,NULL);CHKERRQ(ierr);
749   ierr = PetscOptionsReal("-tao_BNK_omega2", "", "", bnk->omega2, &bnk->omega2,NULL);CHKERRQ(ierr);
750   ierr = PetscOptionsReal("-tao_BNK_omega3", "", "", bnk->omega3, &bnk->omega3,NULL);CHKERRQ(ierr);
751   ierr = PetscOptionsReal("-tao_BNK_omega4", "", "", bnk->omega4, &bnk->omega4,NULL);CHKERRQ(ierr);
752   ierr = PetscOptionsReal("-tao_BNK_omega5", "", "", bnk->omega5, &bnk->omega5,NULL);CHKERRQ(ierr);
753   ierr = PetscOptionsReal("-tao_BNK_mu1_i", "", "", bnk->mu1_i, &bnk->mu1_i,NULL);CHKERRQ(ierr);
754   ierr = PetscOptionsReal("-tao_BNK_mu2_i", "", "", bnk->mu2_i, &bnk->mu2_i,NULL);CHKERRQ(ierr);
755   ierr = PetscOptionsReal("-tao_BNK_gamma1_i", "", "", bnk->gamma1_i, &bnk->gamma1_i,NULL);CHKERRQ(ierr);
756   ierr = PetscOptionsReal("-tao_BNK_gamma2_i", "", "", bnk->gamma2_i, &bnk->gamma2_i,NULL);CHKERRQ(ierr);
757   ierr = PetscOptionsReal("-tao_BNK_gamma3_i", "", "", bnk->gamma3_i, &bnk->gamma3_i,NULL);CHKERRQ(ierr);
758   ierr = PetscOptionsReal("-tao_BNK_gamma4_i", "", "", bnk->gamma4_i, &bnk->gamma4_i,NULL);CHKERRQ(ierr);
759   ierr = PetscOptionsReal("-tao_BNK_theta_i", "", "", bnk->theta_i, &bnk->theta_i,NULL);CHKERRQ(ierr);
760   ierr = PetscOptionsReal("-tao_BNK_mu1", "", "", bnk->mu1, &bnk->mu1,NULL);CHKERRQ(ierr);
761   ierr = PetscOptionsReal("-tao_BNK_mu2", "", "", bnk->mu2, &bnk->mu2,NULL);CHKERRQ(ierr);
762   ierr = PetscOptionsReal("-tao_BNK_gamma1", "", "", bnk->gamma1, &bnk->gamma1,NULL);CHKERRQ(ierr);
763   ierr = PetscOptionsReal("-tao_BNK_gamma2", "", "", bnk->gamma2, &bnk->gamma2,NULL);CHKERRQ(ierr);
764   ierr = PetscOptionsReal("-tao_BNK_gamma3", "", "", bnk->gamma3, &bnk->gamma3,NULL);CHKERRQ(ierr);
765   ierr = PetscOptionsReal("-tao_BNK_gamma4", "", "", bnk->gamma4, &bnk->gamma4,NULL);CHKERRQ(ierr);
766   ierr = PetscOptionsReal("-tao_BNK_theta", "", "", bnk->theta, &bnk->theta,NULL);CHKERRQ(ierr);
767   ierr = PetscOptionsReal("-tao_BNK_min_radius", "lower bound on initial radius", "", bnk->min_radius, &bnk->min_radius,NULL);CHKERRQ(ierr);
768   ierr = PetscOptionsReal("-tao_BNK_max_radius", "upper bound on radius", "", bnk->max_radius, &bnk->max_radius,NULL);CHKERRQ(ierr);
769   ierr = PetscOptionsReal("-tao_BNK_epsilon", "tolerance used when computing actual and predicted reduction", "", bnk->epsilon, &bnk->epsilon,NULL);CHKERRQ(ierr);
770   ierr = PetscOptionsTail();CHKERRQ(ierr);
771   ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr);
772   ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
773   PetscFunctionReturn(0);
774 }
775 
776 /*------------------------------------------------------------*/
777 static PetscErrorCode TaoView_BNK(Tao tao, PetscViewer viewer)
778 {
779   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
780   PetscInt       nrejects;
781   PetscBool      isascii;
782   PetscErrorCode ierr;
783 
784   PetscFunctionBegin;
785   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
786   if (isascii) {
787     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
788     if (BNK_PC_BFGS == bnk->pc_type && bnk->M) {
789       ierr = MatLMVMGetRejects(bnk->M,&nrejects);CHKERRQ(ierr);
790       ierr = PetscViewerASCIIPrintf(viewer, "Rejected matrix updates: %D\n",nrejects);CHKERRQ(ierr);
791     }
792     ierr = PetscViewerASCIIPrintf(viewer, "Newton steps: %D\n", bnk->newt);CHKERRQ(ierr);
793     ierr = PetscViewerASCIIPrintf(viewer, "BFGS steps: %D\n", bnk->bfgs);CHKERRQ(ierr);
794     ierr = PetscViewerASCIIPrintf(viewer, "Scaled gradient steps: %D\n", bnk->sgrad);CHKERRQ(ierr);
795     ierr = PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", bnk->grad);CHKERRQ(ierr);
796     ierr = PetscViewerASCIIPrintf(viewer, "KSP termination reasons:\n");CHKERRQ(ierr);
797     ierr = PetscViewerASCIIPrintf(viewer, "  atol: %D\n", bnk->ksp_atol);CHKERRQ(ierr);
798     ierr = PetscViewerASCIIPrintf(viewer, "  rtol: %D\n", bnk->ksp_rtol);CHKERRQ(ierr);
799     ierr = PetscViewerASCIIPrintf(viewer, "  ctol: %D\n", bnk->ksp_ctol);CHKERRQ(ierr);
800     ierr = PetscViewerASCIIPrintf(viewer, "  negc: %D\n", bnk->ksp_negc);CHKERRQ(ierr);
801     ierr = PetscViewerASCIIPrintf(viewer, "  dtol: %D\n", bnk->ksp_dtol);CHKERRQ(ierr);
802     ierr = PetscViewerASCIIPrintf(viewer, "  iter: %D\n", bnk->ksp_iter);CHKERRQ(ierr);
803     ierr = PetscViewerASCIIPrintf(viewer, "  othr: %D\n", bnk->ksp_othr);CHKERRQ(ierr);
804     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
805   }
806   PetscFunctionReturn(0);
807 }
808 
809 /* ---------------------------------------------------------- */
810 /*MC
811   TAOBNK - Shared base-type for Bounded Newton-Krylov type algorithms.
812   At each iteration, the BNK method solves the symmetric
813   system of equations to obtain the step diretion dk:
814               Hk dk = -gk
815   at which point the step can be globalized either through trust-region
816   methods, or a line search, or a heuristic mixture of both.
817 
818     Options Database Keys:
819 + -tao_BNK_pc_type - "none","ahess","bfgs","petsc"
820 . -tao_BNK_bfgs_scale_type - "ahess","phess","bfgs"
821 . -tao_BNK_init_type - "constant","direction","interpolation"
822 . -tao_BNK_update_type - "step","direction","interpolation"
823 . -tao_BNK_sval - perturbation starting value
824 . -tao_BNK_imin - minimum initial perturbation
825 . -tao_BNK_imax - maximum initial perturbation
826 . -tao_BNK_pmin - minimum perturbation
827 . -tao_BNK_pmax - maximum perturbation
828 . -tao_BNK_pgfac - growth factor
829 . -tao_BNK_psfac - shrink factor
830 . -tao_BNK_imfac - initial merit factor
831 . -tao_BNK_pmgfac - merit growth factor
832 . -tao_BNK_pmsfac - merit shrink factor
833 . -tao_BNK_eta1 - poor steplength; reduce radius
834 . -tao_BNK_eta2 - reasonable steplength; leave radius
835 . -tao_BNK_eta3 - good steplength; increase readius
836 . -tao_BNK_eta4 - excellent steplength; greatly increase radius
837 . -tao_BNK_alpha1 - alpha1 reduction
838 . -tao_BNK_alpha2 - alpha2 reduction
839 . -tao_BNK_alpha3 - alpha3 reduction
840 . -tao_BNK_alpha4 - alpha4 reduction
841 . -tao_BNK_alpha - alpha5 reduction
842 . -tao_BNK_mu1 - mu1 interpolation update
843 . -tao_BNK_mu2 - mu2 interpolation update
844 . -tao_BNK_gamma1 - gamma1 interpolation update
845 . -tao_BNK_gamma2 - gamma2 interpolation update
846 . -tao_BNK_gamma3 - gamma3 interpolation update
847 . -tao_BNK_gamma4 - gamma4 interpolation update
848 . -tao_BNK_theta - theta interpolation update
849 . -tao_BNK_omega1 - omega1 step update
850 . -tao_BNK_omega2 - omega2 step update
851 . -tao_BNK_omega3 - omega3 step update
852 . -tao_BNK_omega4 - omega4 step update
853 . -tao_BNK_omega5 - omega5 step update
854 . -tao_BNK_mu1_i -  mu1 interpolation init factor
855 . -tao_BNK_mu2_i -  mu2 interpolation init factor
856 . -tao_BNK_gamma1_i -  gamma1 interpolation init factor
857 . -tao_BNK_gamma2_i -  gamma2 interpolation init factor
858 . -tao_BNK_gamma3_i -  gamma3 interpolation init factor
859 . -tao_BNK_gamma4_i -  gamma4 interpolation init factor
860 - -tao_BNK_theta_i -  theta interpolation init factor
861 
862   Level: beginner
863 M*/
864 
865 PetscErrorCode TaoComputeBounds_BNK(Tao tao, Vec XL, Vec XU, void *ctx)
866 {
867   PetscFunctionBegin;
868   PetscFunctionReturn(0);
869 }
870 
871 PetscErrorCode TaoCreate_BNK(Tao tao)
872 {
873   TAO_BNK        *bnk;
874   const char     *morethuente_type = TAOLINESEARCHMT;
875   PetscErrorCode ierr;
876 
877   PetscFunctionBegin;
878   ierr = PetscNewLog(tao,&bnk);CHKERRQ(ierr);
879 
880   tao->ops->setup = TaoSetUp_BNK;
881   tao->ops->view = TaoView_BNK;
882   tao->ops->setfromoptions = TaoSetFromOptions_BNK;
883   tao->ops->destroy = TaoDestroy_BNK;
884   tao->ops->computebounds = TaoComputeBounds_BNK;
885 
886   /* Override default settings (unless already changed) */
887   if (!tao->max_it_changed) tao->max_it = 50;
888   if (!tao->trust0_changed) tao->trust0 = 100.0;
889 
890   tao->data = (void*)bnk;
891 
892   bnk->sval   = 0.0;
893   bnk->imin   = 1.0e-4;
894   bnk->imax   = 1.0e+2;
895   bnk->imfac  = 1.0e-1;
896 
897   bnk->pmin   = 1.0e-12;
898   bnk->pmax   = 1.0e+2;
899   bnk->pgfac  = 1.0e+1;
900   bnk->psfac  = 4.0e-1;
901   bnk->pmgfac = 1.0e-1;
902   bnk->pmsfac = 1.0e-1;
903 
904   /*  Default values for trust-region radius update based on steplength */
905   bnk->nu1 = 0.25;
906   bnk->nu2 = 0.50;
907   bnk->nu3 = 1.00;
908   bnk->nu4 = 1.25;
909 
910   bnk->omega1 = 0.25;
911   bnk->omega2 = 0.50;
912   bnk->omega3 = 1.00;
913   bnk->omega4 = 2.00;
914   bnk->omega5 = 4.00;
915 
916   /*  Default values for trust-region radius update based on reduction */
917   bnk->eta1 = 1.0e-4;
918   bnk->eta2 = 0.25;
919   bnk->eta3 = 0.50;
920   bnk->eta4 = 0.90;
921 
922   bnk->alpha1 = 0.25;
923   bnk->alpha2 = 0.50;
924   bnk->alpha3 = 1.00;
925   bnk->alpha4 = 2.00;
926   bnk->alpha5 = 4.00;
927 
928   /*  Default values for trust-region radius update based on interpolation */
929   bnk->mu1 = 0.10;
930   bnk->mu2 = 0.50;
931 
932   bnk->gamma1 = 0.25;
933   bnk->gamma2 = 0.50;
934   bnk->gamma3 = 2.00;
935   bnk->gamma4 = 4.00;
936 
937   bnk->theta = 0.05;
938 
939   /*  Default values for trust region initialization based on interpolation */
940   bnk->mu1_i = 0.35;
941   bnk->mu2_i = 0.50;
942 
943   bnk->gamma1_i = 0.0625;
944   bnk->gamma2_i = 0.5;
945   bnk->gamma3_i = 2.0;
946   bnk->gamma4_i = 5.0;
947 
948   bnk->theta_i = 0.25;
949 
950   /*  Remaining parameters */
951   bnk->min_radius = 1.0e-10;
952   bnk->max_radius = 1.0e10;
953   bnk->epsilon = 1.0e-6;
954 
955   bnk->pc_type         = BNK_PC_BFGS;
956   bnk->bfgs_scale_type = BFGS_SCALE_PHESS;
957   bnk->init_type       = BNK_INIT_INTERPOLATION;
958   bnk->update_type     = BNK_UPDATE_STEP;
959 
960   ierr = TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch);CHKERRQ(ierr);
961   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr);
962   ierr = TaoLineSearchSetType(tao->linesearch,morethuente_type);CHKERRQ(ierr);
963   ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr);
964   ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr);
965 
966   /*  Set linear solver to default for symmetric matrices */
967   ierr = KSPCreate(((PetscObject)tao)->comm,&tao->ksp);CHKERRQ(ierr);
968   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1);CHKERRQ(ierr);
969   ierr = KSPSetOptionsPrefix(tao->ksp,tao->hdr.prefix);CHKERRQ(ierr);
970   ierr = KSPSetType(tao->ksp,KSPCGSTCG);CHKERRQ(ierr);
971   PetscFunctionReturn(0);
972 }
973