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