xref: /petsc/src/tao/bound/impls/bnk/bnk.c (revision 080d2917cae42b0a43dd96fa397634a71dcc282e)
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 = MatLMVMSolve(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,tao->gradient);CHKERRQ(ierr);
208 
209           ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&bnk->gnorm);CHKERRQ(ierr);
210           if (PetscIsInfOrNanReal(bnk->gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute gradient generated Inf or NaN");
211           needH = 1;
212 
213           ierr = TaoLogConvergenceHistory(tao,bnk->f,bnk->gnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
214           ierr = TaoMonitor(tao,tao->niter,bnk->f,bnk->gnorm,0.0,step);CHKERRQ(ierr);
215           ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
216           if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
217         }
218       }
219       tao->trust = PetscMax(tao->trust, max_radius);
220 
221       /* Modify the radius if it is too large or small */
222       tao->trust = PetscMax(tao->trust, bnk->min_radius);
223       tao->trust = PetscMin(tao->trust, bnk->max_radius);
224       break;
225 
226     default:
227       /* Norm of the first direction will initialize radius */
228       tao->trust = 0.0;
229       break;
230     }
231   }
232 
233   /* Set initial scaling for the BFGS preconditioner
234      This step is done after computing the initial trust-region radius
235      since the function value may have decreased */
236   if (BNK_PC_BFGS == bnk->pc_type) {
237     if (bnk->f != 0.0) {
238       delta = 2.0 * PetscAbsScalar(bnk->f) / (bnk->gnorm*bnk->gnorm);
239     } else {
240       delta = 2.0 / (bnk->gnorm*bnk->gnorm);
241     }
242     ierr = MatLMVMSetDelta(bnk->M,delta);CHKERRQ(ierr);
243   }
244 
245   /* Set counter for gradient/reset steps*/
246   bnk->newt = 0;
247   bnk->bfgs = 0;
248   bnk->sgrad = 0;
249   bnk->grad = 0;
250   PetscFunctionReturn(0);
251 }
252 
253 PetscErrorCode TaoBNKComputeStep(Tao tao, PetscInt *stepType)
254 {
255   PetscErrorCode               ierr;
256   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
257   KSPConvergedReason           ksp_reason;
258 
259   PetscReal                    gdx, delta, e_min;
260 
261   PetscInt                     bfgsUpdates = 0;
262   PetscInt                     kspits;
263   PetscInt                     needH = 1;
264 
265   PetscFunctionBegin;
266   /* Compute the Hessian */
267   if (needH) {
268     ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr);
269   }
270 
271   if ((BNK_PC_BFGS == bnk->pc_type) && (BFGS_SCALE_AHESS == bnk->bfgs_scale_type)) {
272     /* Obtain diagonal for the bfgs preconditioner  */
273     ierr = MatGetDiagonal(tao->hessian, bnk->Diag);CHKERRQ(ierr);
274     ierr = VecAbs(bnk->Diag);CHKERRQ(ierr);
275     ierr = VecReciprocal(bnk->Diag);CHKERRQ(ierr);
276     ierr = MatLMVMSetScale(bnk->M,bnk->Diag);CHKERRQ(ierr);
277   }
278 
279   /* Shift the Hessian matrix */
280   bnk->pert = bnk->sval;
281   if (bnk->pert > 0) {
282     ierr = MatShift(tao->hessian, bnk->pert);CHKERRQ(ierr);
283     if (tao->hessian != tao->hessian_pre) {
284       ierr = MatShift(tao->hessian_pre, bnk->pert);CHKERRQ(ierr);
285     }
286   }
287 
288   if (BNK_PC_BFGS == bnk->pc_type) {
289     if (BFGS_SCALE_PHESS == bnk->bfgs_scale_type) {
290       /* Obtain diagonal for the bfgs preconditioner  */
291       ierr = MatGetDiagonal(tao->hessian, bnk->Diag);CHKERRQ(ierr);
292       ierr = VecAbs(bnk->Diag);CHKERRQ(ierr);
293       ierr = VecReciprocal(bnk->Diag);CHKERRQ(ierr);
294       ierr = MatLMVMSetScale(bnk->M,bnk->Diag);CHKERRQ(ierr);
295     }
296     /* Update the limited memory preconditioner and get existing # of updates */
297     ierr = MatLMVMUpdate(bnk->M, tao->solution, tao->gradient);CHKERRQ(ierr);
298     ierr = MatLMVMGetUpdates(bnk->M, &bfgsUpdates);CHKERRQ(ierr);
299   }
300 
301   /* Solve the Newton system of equations */
302   ierr = KSPSetOperators(tao->ksp,tao->hessian,tao->hessian_pre);CHKERRQ(ierr);
303   if (bnk->is_nash || bnk->is_stcg || bnk->is_gltr) {
304     ierr = KSPCGSetRadius(tao->ksp,bnk->max_radius);CHKERRQ(ierr);
305     ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr);
306     ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr);
307     tao->ksp_its+=kspits;
308     tao->ksp_tot_its+=kspits;
309     ierr = KSPCGGetNormD(tao->ksp,&bnk->dnorm);CHKERRQ(ierr);
310 
311     if (0.0 == tao->trust) {
312       /* Radius was uninitialized; use the norm of the direction */
313       if (bnk->dnorm > 0.0) {
314         tao->trust = bnk->dnorm;
315 
316         /* Modify the radius if it is too large or small */
317         tao->trust = PetscMax(tao->trust, bnk->min_radius);
318         tao->trust = PetscMin(tao->trust, bnk->max_radius);
319       } else {
320         /* The direction was bad; set radius to default value and re-solve
321            the trust-region subproblem to get a direction */
322         tao->trust = tao->trust0;
323 
324         /* Modify the radius if it is too large or small */
325         tao->trust = PetscMax(tao->trust, bnk->min_radius);
326         tao->trust = PetscMin(tao->trust, bnk->max_radius);
327 
328         ierr = KSPCGSetRadius(tao->ksp,bnk->max_radius);CHKERRQ(ierr);
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         ierr = KSPCGGetNormD(tao->ksp,&bnk->dnorm);CHKERRQ(ierr);
334 
335         if (bnk->dnorm == 0.0) SETERRQ(PETSC_COMM_SELF,1, "Initial direction zero");
336       }
337     }
338   } else {
339     ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr);
340     ierr = KSPGetIterationNumber(tao->ksp, &kspits);CHKERRQ(ierr);
341     tao->ksp_its += kspits;
342     tao->ksp_tot_its+=kspits;
343   }
344   ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
345   ierr = KSPGetConvergedReason(tao->ksp, &ksp_reason);CHKERRQ(ierr);
346   if ((KSP_DIVERGED_INDEFINITE_PC == ksp_reason) &&  (BNK_PC_BFGS == bnk->pc_type) && (bfgsUpdates > 1)) {
347     /* Preconditioner is numerically indefinite; reset the
348        approximate if using BFGS preconditioning. */
349 
350     if (bnk->f != 0.0) {
351       delta = 2.0 * PetscAbsScalar(bnk->f) / (bnk->gnorm*bnk->gnorm);
352     } else {
353       delta = 2.0 / (bnk->gnorm*bnk->gnorm);
354     }
355     ierr = MatLMVMSetDelta(bnk->M,delta);CHKERRQ(ierr);
356     ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr);
357     ierr = MatLMVMUpdate(bnk->M, tao->solution, tao->gradient);CHKERRQ(ierr);
358     bfgsUpdates = 1;
359   }
360 
361   if (KSP_CONVERGED_ATOL == ksp_reason) {
362     ++bnk->ksp_atol;
363   } else if (KSP_CONVERGED_RTOL == ksp_reason) {
364     ++bnk->ksp_rtol;
365   } else if (KSP_CONVERGED_CG_CONSTRAINED == ksp_reason) {
366     ++bnk->ksp_ctol;
367   } else if (KSP_CONVERGED_CG_NEG_CURVE == ksp_reason) {
368     ++bnk->ksp_negc;
369   } else if (KSP_DIVERGED_DTOL == ksp_reason) {
370     ++bnk->ksp_dtol;
371   } else if (KSP_DIVERGED_ITS == ksp_reason) {
372     ++bnk->ksp_iter;
373   } else {
374     ++bnk->ksp_othr;
375   }
376 
377   /* Check for success (descent direction) */
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, tao->gradient, tao->stepdirection);CHKERRQ(ierr);
404       ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
405 
406       /* Check for success (descent direction) */
407       ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr);
408       if ((gdx >= 0) || PetscIsInfOrNanReal(gdx)) {
409         /* BFGS direction is not descent or direction produced not a number
410            We can assert bfgsUpdates > 1 in this case because
411            the first solve produces the scaled gradient direction,
412            which is guaranteed to be descent */
413 
414         /* Use steepest descent direction (scaled) */
415 
416         if (bnk->f != 0.0) {
417           delta = 2.0 * PetscAbsScalar(bnk->f) / (bnk->gnorm*bnk->gnorm);
418         } else {
419           delta = 2.0 / (bnk->gnorm*bnk->gnorm);
420         }
421         ierr = MatLMVMSetDelta(bnk->M, delta);CHKERRQ(ierr);
422         ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr);
423         ierr = MatLMVMUpdate(bnk->M, tao->solution, tao->gradient);CHKERRQ(ierr);
424         ierr = MatLMVMSolve(bnk->M, tao->gradient, tao->stepdirection);CHKERRQ(ierr);
425         ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
426 
427         bfgsUpdates = 1;
428         ++bnk->sgrad;
429         *stepType = BNK_SCALED_GRADIENT;
430       } else {
431         if (1 == bfgsUpdates) {
432           /* The first BFGS direction is always the scaled gradient */
433           ++bnk->sgrad;
434           *stepType = BNK_SCALED_GRADIENT;
435         } else {
436           ++bnk->bfgs;
437           *stepType = BNK_BFGS;
438         }
439       }
440     }
441   } else {
442     /* Computed Newton step is descent */
443     switch (ksp_reason) {
444     case KSP_DIVERGED_NANORINF:
445     case KSP_DIVERGED_BREAKDOWN:
446     case KSP_DIVERGED_INDEFINITE_MAT:
447     case KSP_DIVERGED_INDEFINITE_PC:
448     case KSP_CONVERGED_CG_NEG_CURVE:
449       /* Matrix or preconditioner is indefinite; increase perturbation */
450       if (bnk->pert <= 0.0) {
451         /* Initialize the perturbation */
452         bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
453         if (bnk->is_gltr) {
454           ierr = KSPCGGLTRGetMinEig(tao->ksp, &e_min);CHKERRQ(ierr);
455           bnk->pert = PetscMax(bnk->pert, -e_min);
456         }
457       } else {
458         /* Increase the perturbation */
459         bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
460       }
461       break;
462 
463     default:
464       /* Newton step computation is good; decrease perturbation */
465       bnk->pert = PetscMin(bnk->psfac * bnk->pert, bnk->pmsfac * bnk->gnorm);
466       if (bnk->pert < bnk->pmin) {
467         bnk->pert = 0.0;
468       }
469       break;
470     }
471 
472     ++bnk->newt;
473     stepType = BNK_NEWTON;
474   }
475   PetscFunctionReturn(0);
476 }
477 
478 PetscErrorCode TaoBNKUpdateTrustRadius(Tao tao, PetscReal fold, PetscReal fnew, PetscInt stepType, PetscBool *accept)
479 {
480   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
481   PetscErrorCode ierr;
482 
483   PetscReal      step, prered, actred, kappa;
484   PetscReal      gdx, tau_1, tau_2, tau_min, tau_max;
485 
486   PetscFunctionBegin;
487   /* Update trust region radius */
488   *accept = PETSC_FALSE;
489   switch(bnk->update_type) {
490   case BNK_UPDATE_STEP:
491     if (stepType == BNK_NEWTON) {
492       ierr = TaoLineSearchGetStepLength(tao->linesearch, &step);CHKERRQ(ierr);
493       if (step < bnk->nu1) {
494         /* Very bad step taken; reduce radius */
495         tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust);
496       } else if (step < bnk->nu2) {
497         /* Reasonably bad step taken; reduce radius */
498         tao->trust = bnk->omega2 * PetscMin(bnk->dnorm, tao->trust);
499       } else if (step < bnk->nu3) {
500         /*  Reasonable step was taken; leave radius alone */
501         if (bnk->omega3 < 1.0) {
502           tao->trust = bnk->omega3 * PetscMin(bnk->dnorm, tao->trust);
503         } else if (bnk->omega3 > 1.0) {
504           tao->trust = PetscMax(bnk->omega3 * bnk->dnorm, tao->trust);
505         }
506       } else if (step < bnk->nu4) {
507         /*  Full step taken; increase the radius */
508         tao->trust = PetscMax(bnk->omega4 * bnk->dnorm, tao->trust);
509       } else {
510         /*  More than full step taken; increase the radius */
511         tao->trust = PetscMax(bnk->omega5 * bnk->dnorm, tao->trust);
512       }
513     } else {
514       /*  Newton step was not good; reduce the radius */
515       tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust);
516     }
517     break;
518 
519   case BNK_UPDATE_REDUCTION:
520     if (stepType == BNK_NEWTON) {
521       /*  Get predicted reduction */
522       ierr = KSPCGGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr);
523       if (prered >= 0.0) {
524         /*  The predicted reduction has the wrong sign.  This cannot */
525         /*  happen in infinite precision arithmetic.  Step should */
526         /*  be rejected! */
527         tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
528       } else {
529         if (PetscIsInfOrNanReal(fnew)) {
530           tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
531         } else {
532           /*  Compute and actual reduction */
533           actred = fold - fnew;
534           prered = -prered;
535           if ((PetscAbsScalar(actred) <= bnk->epsilon) &&
536               (PetscAbsScalar(prered) <= bnk->epsilon)) {
537             kappa = 1.0;
538           } else {
539             kappa = actred / prered;
540           }
541           /*  Accept of reject the step and update radius */
542           if (kappa < bnk->eta1) {
543             /*  Very bad step, rejected */
544             tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
545           } else {
546             /* Accept step here */
547             *accept = PETSC_TRUE;
548             /* Adjust trust radius */
549             if (kappa < bnk->eta2) {
550               /*  Marginal bad step */
551               tao->trust = bnk->alpha2 * PetscMin(tao->trust, bnk->dnorm);
552             } else if (kappa < bnk->eta3) {
553               /*  Reasonable step */
554               if (bnk->alpha3 < 1.0) {
555                 tao->trust = bnk->alpha3 * PetscMin(bnk->dnorm, tao->trust);
556               } else if (bnk->alpha3 > 1.0) {
557                 tao->trust = PetscMax(bnk->alpha3 * bnk->dnorm, tao->trust);
558               }
559             } else if (kappa < bnk->eta4) {
560               /*  Good step */
561               tao->trust = PetscMax(bnk->alpha4 * bnk->dnorm, tao->trust);
562             } else {
563               /*  Very good step */
564               tao->trust = PetscMax(bnk->alpha5 * bnk->dnorm, tao->trust);
565             }
566           }
567         }
568       }
569     } else {
570       /*  Newton step was not good; reduce the radius */
571       tao->trust = bnk->alpha1 * PetscMin(bnk->dnorm, tao->trust);
572     }
573     break;
574 
575   default:
576     if (stepType == BNK_NEWTON) {
577       ierr = KSPCGGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr);
578       if (prered >= 0.0) {
579         /*  The predicted reduction has the wrong sign.  This cannot */
580         /*  happen in infinite precision arithmetic.  Step should */
581         /*  be rejected! */
582         tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
583       } else {
584         if (PetscIsInfOrNanReal(fnew)) {
585           tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
586         } else {
587           actred = fold - fnew;
588           prered = -prered;
589           if ((PetscAbsScalar(actred) <= bnk->epsilon) && (PetscAbsScalar(prered) <= bnk->epsilon)) {
590             kappa = 1.0;
591           } else {
592             kappa = actred / prered;
593           }
594 
595           ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr);
596           tau_1 = bnk->theta * gdx / (bnk->theta * gdx - (1.0 - bnk->theta) * prered + actred);
597           tau_2 = bnk->theta * gdx / (bnk->theta * gdx + (1.0 + bnk->theta) * prered - actred);
598           tau_min = PetscMin(tau_1, tau_2);
599           tau_max = PetscMax(tau_1, tau_2);
600 
601           if (kappa >= 1.0 - bnk->mu1) {
602             /*  Great agreement */
603             *accept = PETSC_TRUE;
604             if (tau_max < 1.0) {
605               tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm);
606             } else if (tau_max > bnk->gamma4) {
607               tao->trust = PetscMax(tao->trust, bnk->gamma4 * bnk->dnorm);
608             } else {
609               tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm);
610             }
611           } else if (kappa >= 1.0 - bnk->mu2) {
612             /*  Good agreement */
613             *accept = PETSC_TRUE;
614             if (tau_max < bnk->gamma2) {
615               tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm);
616             } else if (tau_max > bnk->gamma3) {
617               tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm);
618             } else if (tau_max < 1.0) {
619               tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm);
620             } else {
621               tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm);
622             }
623           } else {
624             /*  Not good agreement */
625             if (tau_min > 1.0) {
626               tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm);
627             } else if (tau_max < bnk->gamma1) {
628               tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
629             } else if ((tau_min < bnk->gamma1) && (tau_max >= 1.0)) {
630               tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
631             } else if ((tau_1 >= bnk->gamma1) && (tau_1 < 1.0) && ((tau_2 < bnk->gamma1) || (tau_2 >= 1.0))) {
632               tao->trust = tau_1 * PetscMin(tao->trust, bnk->dnorm);
633             } else if ((tau_2 >= bnk->gamma1) && (tau_2 < 1.0) && ((tau_1 < bnk->gamma1) || (tau_2 >= 1.0))) {
634               tao->trust = tau_2 * PetscMin(tao->trust, bnk->dnorm);
635             } else {
636               tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm);
637             }
638           }
639         }
640       }
641     } else {
642       /*  Newton step was not good; reduce the radius */
643       tao->trust = bnk->gamma1 * PetscMin(bnk->dnorm, tao->trust);
644     }
645     /*  The radius may have been increased; modify if it is too large */
646     tao->trust = PetscMin(tao->trust, bnk->max_radius);
647   }
648   PetscFunctionReturn(0);
649 }
650 
651 /* ---------------------------------------------------------- */
652 static PetscErrorCode TaoSetUp_BNK(Tao tao)
653 {
654   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
655   PetscErrorCode ierr;
656 
657   PetscFunctionBegin;
658   if (!tao->gradient) {ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);}
659   if (!tao->stepdirection) {ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr);}
660   if (!bnk->W) {ierr = VecDuplicate(tao->solution,&bnk->W);CHKERRQ(ierr);}
661   if (!bnk->Xold) {ierr = VecDuplicate(tao->solution,&bnk->Xold);CHKERRQ(ierr);}
662   if (!bnk->Gold) {ierr = VecDuplicate(tao->solution,&bnk->Gold);CHKERRQ(ierr);}
663   bnk->Diag = 0;
664   bnk->M = 0;
665   PetscFunctionReturn(0);
666 }
667 
668 /*------------------------------------------------------------*/
669 static PetscErrorCode TaoDestroy_BNK(Tao tao)
670 {
671   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
672   PetscErrorCode ierr;
673 
674   PetscFunctionBegin;
675   if (tao->setupcalled) {
676     ierr = VecDestroy(&bnk->W);CHKERRQ(ierr);
677     ierr = VecDestroy(&bnk->Xold);CHKERRQ(ierr);
678     ierr = VecDestroy(&bnk->Gold);CHKERRQ(ierr);
679   }
680   ierr = VecDestroy(&bnk->Diag);CHKERRQ(ierr);
681   ierr = MatDestroy(&bnk->M);CHKERRQ(ierr);
682   ierr = PetscFree(tao->data);CHKERRQ(ierr);
683   PetscFunctionReturn(0);
684 }
685 
686 /*------------------------------------------------------------*/
687 static PetscErrorCode TaoSetFromOptions_BNK(PetscOptionItems *PetscOptionsObject,Tao tao)
688 {
689   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
690   PetscErrorCode ierr;
691 
692   PetscFunctionBegin;
693   ierr = PetscOptionsHead(PetscOptionsObject,"Newton line search method for unconstrained optimization");CHKERRQ(ierr);
694   ierr = PetscOptionsEList("-tao_BNK_pc_type", "pc type", "", BNK_PC, BNK_PC_TYPES, BNK_PC[bnk->pc_type], &bnk->pc_type, 0);CHKERRQ(ierr);
695   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);
696   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);
697   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);
698   ierr = PetscOptionsReal("-tao_BNK_sval", "perturbation starting value", "", bnk->sval, &bnk->sval,NULL);CHKERRQ(ierr);
699   ierr = PetscOptionsReal("-tao_BNK_imin", "minimum initial perturbation", "", bnk->imin, &bnk->imin,NULL);CHKERRQ(ierr);
700   ierr = PetscOptionsReal("-tao_BNK_imax", "maximum initial perturbation", "", bnk->imax, &bnk->imax,NULL);CHKERRQ(ierr);
701   ierr = PetscOptionsReal("-tao_BNK_imfac", "initial merit factor", "", bnk->imfac, &bnk->imfac,NULL);CHKERRQ(ierr);
702   ierr = PetscOptionsReal("-tao_BNK_pmin", "minimum perturbation", "", bnk->pmin, &bnk->pmin,NULL);CHKERRQ(ierr);
703   ierr = PetscOptionsReal("-tao_BNK_pmax", "maximum perturbation", "", bnk->pmax, &bnk->pmax,NULL);CHKERRQ(ierr);
704   ierr = PetscOptionsReal("-tao_BNK_pgfac", "growth factor", "", bnk->pgfac, &bnk->pgfac,NULL);CHKERRQ(ierr);
705   ierr = PetscOptionsReal("-tao_BNK_psfac", "shrink factor", "", bnk->psfac, &bnk->psfac,NULL);CHKERRQ(ierr);
706   ierr = PetscOptionsReal("-tao_BNK_pmgfac", "merit growth factor", "", bnk->pmgfac, &bnk->pmgfac,NULL);CHKERRQ(ierr);
707   ierr = PetscOptionsReal("-tao_BNK_pmsfac", "merit shrink factor", "", bnk->pmsfac, &bnk->pmsfac,NULL);CHKERRQ(ierr);
708   ierr = PetscOptionsReal("-tao_BNK_eta1", "poor steplength; reduce radius", "", bnk->eta1, &bnk->eta1,NULL);CHKERRQ(ierr);
709   ierr = PetscOptionsReal("-tao_BNK_eta2", "reasonable steplength; leave radius alone", "", bnk->eta2, &bnk->eta2,NULL);CHKERRQ(ierr);
710   ierr = PetscOptionsReal("-tao_BNK_eta3", "good steplength; increase radius", "", bnk->eta3, &bnk->eta3,NULL);CHKERRQ(ierr);
711   ierr = PetscOptionsReal("-tao_BNK_eta4", "excellent steplength; greatly increase radius", "", bnk->eta4, &bnk->eta4,NULL);CHKERRQ(ierr);
712   ierr = PetscOptionsReal("-tao_BNK_alpha1", "", "", bnk->alpha1, &bnk->alpha1,NULL);CHKERRQ(ierr);
713   ierr = PetscOptionsReal("-tao_BNK_alpha2", "", "", bnk->alpha2, &bnk->alpha2,NULL);CHKERRQ(ierr);
714   ierr = PetscOptionsReal("-tao_BNK_alpha3", "", "", bnk->alpha3, &bnk->alpha3,NULL);CHKERRQ(ierr);
715   ierr = PetscOptionsReal("-tao_BNK_alpha4", "", "", bnk->alpha4, &bnk->alpha4,NULL);CHKERRQ(ierr);
716   ierr = PetscOptionsReal("-tao_BNK_alpha5", "", "", bnk->alpha5, &bnk->alpha5,NULL);CHKERRQ(ierr);
717   ierr = PetscOptionsReal("-tao_BNK_nu1", "poor steplength; reduce radius", "", bnk->nu1, &bnk->nu1,NULL);CHKERRQ(ierr);
718   ierr = PetscOptionsReal("-tao_BNK_nu2", "reasonable steplength; leave radius alone", "", bnk->nu2, &bnk->nu2,NULL);CHKERRQ(ierr);
719   ierr = PetscOptionsReal("-tao_BNK_nu3", "good steplength; increase radius", "", bnk->nu3, &bnk->nu3,NULL);CHKERRQ(ierr);
720   ierr = PetscOptionsReal("-tao_BNK_nu4", "excellent steplength; greatly increase radius", "", bnk->nu4, &bnk->nu4,NULL);CHKERRQ(ierr);
721   ierr = PetscOptionsReal("-tao_BNK_omega1", "", "", bnk->omega1, &bnk->omega1,NULL);CHKERRQ(ierr);
722   ierr = PetscOptionsReal("-tao_BNK_omega2", "", "", bnk->omega2, &bnk->omega2,NULL);CHKERRQ(ierr);
723   ierr = PetscOptionsReal("-tao_BNK_omega3", "", "", bnk->omega3, &bnk->omega3,NULL);CHKERRQ(ierr);
724   ierr = PetscOptionsReal("-tao_BNK_omega4", "", "", bnk->omega4, &bnk->omega4,NULL);CHKERRQ(ierr);
725   ierr = PetscOptionsReal("-tao_BNK_omega5", "", "", bnk->omega5, &bnk->omega5,NULL);CHKERRQ(ierr);
726   ierr = PetscOptionsReal("-tao_BNK_mu1_i", "", "", bnk->mu1_i, &bnk->mu1_i,NULL);CHKERRQ(ierr);
727   ierr = PetscOptionsReal("-tao_BNK_mu2_i", "", "", bnk->mu2_i, &bnk->mu2_i,NULL);CHKERRQ(ierr);
728   ierr = PetscOptionsReal("-tao_BNK_gamma1_i", "", "", bnk->gamma1_i, &bnk->gamma1_i,NULL);CHKERRQ(ierr);
729   ierr = PetscOptionsReal("-tao_BNK_gamma2_i", "", "", bnk->gamma2_i, &bnk->gamma2_i,NULL);CHKERRQ(ierr);
730   ierr = PetscOptionsReal("-tao_BNK_gamma3_i", "", "", bnk->gamma3_i, &bnk->gamma3_i,NULL);CHKERRQ(ierr);
731   ierr = PetscOptionsReal("-tao_BNK_gamma4_i", "", "", bnk->gamma4_i, &bnk->gamma4_i,NULL);CHKERRQ(ierr);
732   ierr = PetscOptionsReal("-tao_BNK_theta_i", "", "", bnk->theta_i, &bnk->theta_i,NULL);CHKERRQ(ierr);
733   ierr = PetscOptionsReal("-tao_BNK_mu1", "", "", bnk->mu1, &bnk->mu1,NULL);CHKERRQ(ierr);
734   ierr = PetscOptionsReal("-tao_BNK_mu2", "", "", bnk->mu2, &bnk->mu2,NULL);CHKERRQ(ierr);
735   ierr = PetscOptionsReal("-tao_BNK_gamma1", "", "", bnk->gamma1, &bnk->gamma1,NULL);CHKERRQ(ierr);
736   ierr = PetscOptionsReal("-tao_BNK_gamma2", "", "", bnk->gamma2, &bnk->gamma2,NULL);CHKERRQ(ierr);
737   ierr = PetscOptionsReal("-tao_BNK_gamma3", "", "", bnk->gamma3, &bnk->gamma3,NULL);CHKERRQ(ierr);
738   ierr = PetscOptionsReal("-tao_BNK_gamma4", "", "", bnk->gamma4, &bnk->gamma4,NULL);CHKERRQ(ierr);
739   ierr = PetscOptionsReal("-tao_BNK_theta", "", "", bnk->theta, &bnk->theta,NULL);CHKERRQ(ierr);
740   ierr = PetscOptionsReal("-tao_BNK_min_radius", "lower bound on initial radius", "", bnk->min_radius, &bnk->min_radius,NULL);CHKERRQ(ierr);
741   ierr = PetscOptionsReal("-tao_BNK_max_radius", "upper bound on radius", "", bnk->max_radius, &bnk->max_radius,NULL);CHKERRQ(ierr);
742   ierr = PetscOptionsReal("-tao_BNK_epsilon", "tolerance used when computing actual and predicted reduction", "", bnk->epsilon, &bnk->epsilon,NULL);CHKERRQ(ierr);
743   ierr = PetscOptionsTail();CHKERRQ(ierr);
744   ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr);
745   ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
746   PetscFunctionReturn(0);
747 }
748 
749 /*------------------------------------------------------------*/
750 static PetscErrorCode TaoView_BNK(Tao tao, PetscViewer viewer)
751 {
752   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
753   PetscInt       nrejects;
754   PetscBool      isascii;
755   PetscErrorCode ierr;
756 
757   PetscFunctionBegin;
758   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
759   if (isascii) {
760     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
761     if (BNK_PC_BFGS == bnk->pc_type && bnk->M) {
762       ierr = MatLMVMGetRejects(bnk->M,&nrejects);CHKERRQ(ierr);
763       ierr = PetscViewerASCIIPrintf(viewer, "Rejected matrix updates: %D\n",nrejects);CHKERRQ(ierr);
764     }
765     ierr = PetscViewerASCIIPrintf(viewer, "Newton steps: %D\n", bnk->newt);CHKERRQ(ierr);
766     ierr = PetscViewerASCIIPrintf(viewer, "BFGS steps: %D\n", bnk->bfgs);CHKERRQ(ierr);
767     ierr = PetscViewerASCIIPrintf(viewer, "Scaled gradient steps: %D\n", bnk->sgrad);CHKERRQ(ierr);
768     ierr = PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", bnk->grad);CHKERRQ(ierr);
769     ierr = PetscViewerASCIIPrintf(viewer, "KSP termination reasons:\n");CHKERRQ(ierr);
770     ierr = PetscViewerASCIIPrintf(viewer, "  atol: %D\n", bnk->ksp_atol);CHKERRQ(ierr);
771     ierr = PetscViewerASCIIPrintf(viewer, "  rtol: %D\n", bnk->ksp_rtol);CHKERRQ(ierr);
772     ierr = PetscViewerASCIIPrintf(viewer, "  ctol: %D\n", bnk->ksp_ctol);CHKERRQ(ierr);
773     ierr = PetscViewerASCIIPrintf(viewer, "  negc: %D\n", bnk->ksp_negc);CHKERRQ(ierr);
774     ierr = PetscViewerASCIIPrintf(viewer, "  dtol: %D\n", bnk->ksp_dtol);CHKERRQ(ierr);
775     ierr = PetscViewerASCIIPrintf(viewer, "  iter: %D\n", bnk->ksp_iter);CHKERRQ(ierr);
776     ierr = PetscViewerASCIIPrintf(viewer, "  othr: %D\n", bnk->ksp_othr);CHKERRQ(ierr);
777     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
778   }
779   PetscFunctionReturn(0);
780 }
781 
782 /* ---------------------------------------------------------- */
783 /*MC
784   TAOBNK - Shared base-type for Bounded Newton-Krylov type algorithms.
785   At each iteration, the BNK method solves the symmetric
786   system of equations to obtain the step diretion dk:
787               Hk dk = -gk
788   at which point the step can be globalized either through trust-region
789   methods, or a line search, or a heuristic mixture of both.
790 
791     Options Database Keys:
792 + -tao_BNK_pc_type - "none","ahess","bfgs","petsc"
793 . -tao_BNK_bfgs_scale_type - "ahess","phess","bfgs"
794 . -tao_BNK_init_type - "constant","direction","interpolation"
795 . -tao_BNK_update_type - "step","direction","interpolation"
796 . -tao_BNK_sval - perturbation starting value
797 . -tao_BNK_imin - minimum initial perturbation
798 . -tao_BNK_imax - maximum initial perturbation
799 . -tao_BNK_pmin - minimum perturbation
800 . -tao_BNK_pmax - maximum perturbation
801 . -tao_BNK_pgfac - growth factor
802 . -tao_BNK_psfac - shrink factor
803 . -tao_BNK_imfac - initial merit factor
804 . -tao_BNK_pmgfac - merit growth factor
805 . -tao_BNK_pmsfac - merit shrink factor
806 . -tao_BNK_eta1 - poor steplength; reduce radius
807 . -tao_BNK_eta2 - reasonable steplength; leave radius
808 . -tao_BNK_eta3 - good steplength; increase readius
809 . -tao_BNK_eta4 - excellent steplength; greatly increase radius
810 . -tao_BNK_alpha1 - alpha1 reduction
811 . -tao_BNK_alpha2 - alpha2 reduction
812 . -tao_BNK_alpha3 - alpha3 reduction
813 . -tao_BNK_alpha4 - alpha4 reduction
814 . -tao_BNK_alpha - alpha5 reduction
815 . -tao_BNK_mu1 - mu1 interpolation update
816 . -tao_BNK_mu2 - mu2 interpolation update
817 . -tao_BNK_gamma1 - gamma1 interpolation update
818 . -tao_BNK_gamma2 - gamma2 interpolation update
819 . -tao_BNK_gamma3 - gamma3 interpolation update
820 . -tao_BNK_gamma4 - gamma4 interpolation update
821 . -tao_BNK_theta - theta interpolation update
822 . -tao_BNK_omega1 - omega1 step update
823 . -tao_BNK_omega2 - omega2 step update
824 . -tao_BNK_omega3 - omega3 step update
825 . -tao_BNK_omega4 - omega4 step update
826 . -tao_BNK_omega5 - omega5 step update
827 . -tao_BNK_mu1_i -  mu1 interpolation init factor
828 . -tao_BNK_mu2_i -  mu2 interpolation init factor
829 . -tao_BNK_gamma1_i -  gamma1 interpolation init factor
830 . -tao_BNK_gamma2_i -  gamma2 interpolation init factor
831 . -tao_BNK_gamma3_i -  gamma3 interpolation init factor
832 . -tao_BNK_gamma4_i -  gamma4 interpolation init factor
833 - -tao_BNK_theta_i -  theta interpolation init factor
834 
835   Level: beginner
836 M*/
837 
838 PetscErrorCode TaoCreate_BNK(Tao tao)
839 {
840   TAO_BNK        *bnk;
841   const char     *morethuente_type = TAOLINESEARCHMT;
842   PetscErrorCode ierr;
843 
844   PetscFunctionBegin;
845   ierr = PetscNewLog(tao,&bnk);CHKERRQ(ierr);
846 
847   tao->ops->setup = TaoSetUp_BNK;
848   tao->ops->view = TaoView_BNK;
849   tao->ops->setfromoptions = TaoSetFromOptions_BNK;
850   tao->ops->destroy = TaoDestroy_BNK;
851 
852   /* Override default settings (unless already changed) */
853   if (!tao->max_it_changed) tao->max_it = 50;
854   if (!tao->trust0_changed) tao->trust0 = 100.0;
855 
856   tao->data = (void*)bnk;
857 
858   bnk->sval   = 0.0;
859   bnk->imin   = 1.0e-4;
860   bnk->imax   = 1.0e+2;
861   bnk->imfac  = 1.0e-1;
862 
863   bnk->pmin   = 1.0e-12;
864   bnk->pmax   = 1.0e+2;
865   bnk->pgfac  = 1.0e+1;
866   bnk->psfac  = 4.0e-1;
867   bnk->pmgfac = 1.0e-1;
868   bnk->pmsfac = 1.0e-1;
869 
870   /*  Default values for trust-region radius update based on steplength */
871   bnk->nu1 = 0.25;
872   bnk->nu2 = 0.50;
873   bnk->nu3 = 1.00;
874   bnk->nu4 = 1.25;
875 
876   bnk->omega1 = 0.25;
877   bnk->omega2 = 0.50;
878   bnk->omega3 = 1.00;
879   bnk->omega4 = 2.00;
880   bnk->omega5 = 4.00;
881 
882   /*  Default values for trust-region radius update based on reduction */
883   bnk->eta1 = 1.0e-4;
884   bnk->eta2 = 0.25;
885   bnk->eta3 = 0.50;
886   bnk->eta4 = 0.90;
887 
888   bnk->alpha1 = 0.25;
889   bnk->alpha2 = 0.50;
890   bnk->alpha3 = 1.00;
891   bnk->alpha4 = 2.00;
892   bnk->alpha5 = 4.00;
893 
894   /*  Default values for trust-region radius update based on interpolation */
895   bnk->mu1 = 0.10;
896   bnk->mu2 = 0.50;
897 
898   bnk->gamma1 = 0.25;
899   bnk->gamma2 = 0.50;
900   bnk->gamma3 = 2.00;
901   bnk->gamma4 = 4.00;
902 
903   bnk->theta = 0.05;
904 
905   /*  Default values for trust region initialization based on interpolation */
906   bnk->mu1_i = 0.35;
907   bnk->mu2_i = 0.50;
908 
909   bnk->gamma1_i = 0.0625;
910   bnk->gamma2_i = 0.5;
911   bnk->gamma3_i = 2.0;
912   bnk->gamma4_i = 5.0;
913 
914   bnk->theta_i = 0.25;
915 
916   /*  Remaining parameters */
917   bnk->min_radius = 1.0e-10;
918   bnk->max_radius = 1.0e10;
919   bnk->epsilon = 1.0e-6;
920 
921   bnk->pc_type         = BNK_PC_BFGS;
922   bnk->bfgs_scale_type = BFGS_SCALE_PHESS;
923   bnk->init_type       = BNK_INIT_INTERPOLATION;
924   bnk->update_type     = BNK_UPDATE_STEP;
925 
926   ierr = TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch);CHKERRQ(ierr);
927   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr);
928   ierr = TaoLineSearchSetType(tao->linesearch,morethuente_type);CHKERRQ(ierr);
929   ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr);
930   ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr);
931 
932   /*  Set linear solver to default for symmetric matrices */
933   ierr = KSPCreate(((PetscObject)tao)->comm,&tao->ksp);CHKERRQ(ierr);
934   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1);CHKERRQ(ierr);
935   ierr = KSPSetOptionsPrefix(tao->ksp,tao->hdr.prefix);CHKERRQ(ierr);
936   ierr = KSPSetType(tao->ksp,KSPCGSTCG);CHKERRQ(ierr);
937   PetscFunctionReturn(0);
938 }
939