xref: /petsc/src/tao/bound/impls/bnk/bnk.c (revision fa2f3a98544cffa8147e57031dd6780dd55215d2)
1 #include <petsctaolinesearch.h>
2 #include <../src/tao/bound/impls/bnk/bnk.h>
3 
4 #include <petscksp.h>
5 
6 static const char *BNK_PC[64] = {"none", "ahess", "bfgs", "petsc"};
7 static const char *BFGS_SCALE[64] = {"ahess", "phess", "bfgs"};
8 static const char *BNK_INIT[64] = {"constant", "direction", "interpolation"};
9 static const char *BNK_UPDATE[64] = {"step", "reduction", "interpolation"};
10 static const char *BNK_AS[64] = {"none", "bertsekas"};
11 
12 /*------------------------------------------------------------*/
13 
14 /* Routine for BFGS preconditioner */
15 
16 PetscErrorCode MatLMVMSolveShell(PC pc, Vec b, Vec x)
17 {
18   PetscErrorCode ierr;
19   Mat            M;
20 
21   PetscFunctionBegin;
22   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
23   PetscValidHeaderSpecific(b,VEC_CLASSID,2);
24   PetscValidHeaderSpecific(x,VEC_CLASSID,3);
25   ierr = PCShellGetContext(pc,(void**)&M);CHKERRQ(ierr);
26   ierr = MatLMVMSolve(M, b, x);CHKERRQ(ierr);
27   PetscFunctionReturn(0);
28 }
29 
30 /*------------------------------------------------------------*/
31 
32 /* Routine for initializing the KSP solver, the BFGS preconditioner, and the initial trust radius estimation */
33 
34 PetscErrorCode TaoBNKInitialize(Tao tao, PetscInt initType, PetscBool *needH)
35 {
36   PetscErrorCode               ierr;
37   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
38   PC                           pc;
39 
40   PetscReal                    f_min, ftrial, prered, actred, kappa, sigma, resnorm;
41   PetscReal                    tau, tau_1, tau_2, tau_max, tau_min, max_radius;
42   PetscReal                    delta;
43 
44   PetscInt                     n,N,nDiff;
45 
46   PetscInt                     i_max = 5;
47   PetscInt                     j_max = 1;
48   PetscInt                     i, j;
49 
50   PetscFunctionBegin;
51   /* Project the current point onto the feasible set */
52   ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr);
53   ierr = TaoSetVariableBounds(bnk->bncg, tao->XL, tao->XU);CHKERRQ(ierr);
54   ierr = TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU);CHKERRQ(ierr);
55 
56   /* Project the initial point onto the feasible region */
57   ierr = TaoBoundSolution(tao->solution, tao->XL,tao->XU, 0.0, &nDiff, tao->solution);CHKERRQ(ierr);
58 
59   /* Check convergence criteria */
60   ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &bnk->f, bnk->unprojected_gradient);CHKERRQ(ierr);
61   ierr = TaoBNKEstimateActiveSet(tao, bnk->as_type);CHKERRQ(ierr);
62   ierr = VecCopy(bnk->unprojected_gradient, tao->gradient);CHKERRQ(ierr);
63   ierr = VecISSet(tao->gradient, bnk->active_idx, 0.0);CHKERRQ(ierr);
64   ierr = VecNorm(tao->gradient,NORM_2,&bnk->gnorm);CHKERRQ(ierr);
65 
66   /* Test the initial point for convergence */
67   ierr = VecFischer(tao->solution, bnk->unprojected_gradient, tao->XL, tao->XU, bnk->W);CHKERRQ(ierr);
68   ierr = VecNorm(bnk->W, NORM_2, &resnorm);CHKERRQ(ierr);
69   if (PetscIsInfOrNanReal(bnk->f) || PetscIsInfOrNanReal(resnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
70   ierr = TaoLogConvergenceHistory(tao,bnk->f,resnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
71   ierr = TaoMonitor(tao,tao->niter,bnk->f,resnorm,0.0,1.0);CHKERRQ(ierr);
72   ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
73   if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
74 
75   /* Reset KSP stopping reason counters */
76   bnk->ksp_atol = 0;
77   bnk->ksp_rtol = 0;
78   bnk->ksp_dtol = 0;
79   bnk->ksp_ctol = 0;
80   bnk->ksp_negc = 0;
81   bnk->ksp_iter = 0;
82   bnk->ksp_othr = 0;
83 
84   /* Reset accepted step type counters */
85   bnk->tot_cg_its = 0;
86   bnk->newt = 0;
87   bnk->bfgs = 0;
88   bnk->sgrad = 0;
89   bnk->grad = 0;
90 
91   /* Initialize the Hessian perturbation */
92   bnk->pert = bnk->sval;
93 
94   /* Reset initial steplength to zero (this helps BNCG reset its direction internally) */
95   ierr = VecSet(tao->stepdirection, 0.0);CHKERRQ(ierr);
96 
97   /* Allocate the vectors needed for the BFGS approximation */
98   if (BNK_PC_BFGS == bnk->pc_type) {
99     if (!bnk->M) {
100       ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr);
101       ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr);
102       ierr = MatCreateLMVM(((PetscObject)tao)->comm,n,N,&bnk->M);CHKERRQ(ierr);
103       ierr = MatLMVMAllocateVectors(bnk->M,tao->solution);CHKERRQ(ierr);
104     }
105     if (bnk->bfgs_scale_type != BFGS_SCALE_BFGS && !bnk->Diag) {
106       ierr = VecDuplicate(tao->solution,&bnk->Diag);CHKERRQ(ierr);
107     }
108   }
109 
110   /* Prepare the min/max vectors for safeguarding diagonal scales */
111   ierr = VecSet(bnk->Diag_min, bnk->dmin);CHKERRQ(ierr);
112   ierr = VecSet(bnk->Diag_max, bnk->dmax);CHKERRQ(ierr);
113 
114   /* Modify the preconditioner to use the bfgs approximation */
115   ierr = KSPGetPC(tao->ksp, &pc);CHKERRQ(ierr);
116   switch(bnk->pc_type) {
117   case BNK_PC_NONE:
118     ierr = PCSetType(pc, PCNONE);CHKERRQ(ierr);
119     ierr = PCSetFromOptions(pc);CHKERRQ(ierr);
120     break;
121 
122   case BNK_PC_AHESS:
123     ierr = PCSetType(pc, PCJACOBI);CHKERRQ(ierr);
124     ierr = PCSetFromOptions(pc);CHKERRQ(ierr);
125     ierr = PCJacobiSetUseAbs(pc,PETSC_TRUE);CHKERRQ(ierr);
126     break;
127 
128   case BNK_PC_BFGS:
129     ierr = PCSetType(pc, PCSHELL);CHKERRQ(ierr);
130     ierr = PCSetFromOptions(pc);CHKERRQ(ierr);
131     ierr = PCShellSetName(pc, "bfgs");CHKERRQ(ierr);
132     ierr = PCShellSetContext(pc, bnk->M);CHKERRQ(ierr);
133     ierr = PCShellSetApply(pc, MatLMVMSolveShell);CHKERRQ(ierr);
134     break;
135 
136   default:
137     /* Use the pc method set by pc_type */
138     break;
139   }
140 
141   /* Initialize trust-region radius.  The initialization is only performed
142      when we are using Nash, Steihaug-Toint or the Generalized Lanczos method. */
143   *needH = PETSC_TRUE;
144   if (bnk->is_nash || bnk->is_stcg || bnk->is_gltr) {
145     switch(initType) {
146     case BNK_INIT_CONSTANT:
147       /* Use the initial radius specified */
148       tao->trust = tao->trust0;
149       break;
150 
151     case BNK_INIT_INTERPOLATION:
152       /* Use interpolation based on the initial Hessian */
153       max_radius = 0.0;
154       tao->trust = tao->trust0;
155       for (j = 0; j < j_max; ++j) {
156         f_min = bnk->f;
157         sigma = 0.0;
158 
159         if (*needH) {
160           /* Compute the Hessian at the new step, and extract the inactive subsystem */
161           ierr = TaoBNKComputeHessian(tao);CHKERRQ(ierr);
162           ierr = TaoBNKEstimateActiveSet(tao, BNK_AS_NONE);CHKERRQ(ierr);
163           ierr = MatDestroy(&bnk->H_inactive);CHKERRQ(ierr);
164           if (bnk->active_idx) {
165             ierr = MatCreateSubMatrix(tao->hessian, bnk->inactive_idx, bnk->inactive_idx, MAT_INITIAL_MATRIX, &bnk->H_inactive);CHKERRQ(ierr);
166           } else {
167             ierr = MatDuplicate(tao->hessian, MAT_COPY_VALUES, &bnk->H_inactive);CHKERRQ(ierr);
168           }
169           *needH = PETSC_FALSE;
170         }
171 
172         for (i = 0; i < i_max; ++i) {
173           /* Take a steepest descent step and snap it to bounds */
174           ierr = VecCopy(tao->solution, bnk->Xold);CHKERRQ(ierr);
175           ierr = VecAXPY(tao->solution, -tao->trust/bnk->gnorm, tao->gradient);CHKERRQ(ierr);
176           ierr = TaoBoundSolution(tao->solution, tao->XL,tao->XU, 0.0, &nDiff, tao->solution);CHKERRQ(ierr);
177           /* Compute the step we actually accepted */
178           ierr = VecCopy(tao->solution, bnk->W);CHKERRQ(ierr);
179           ierr = VecAXPY(bnk->W, -1.0, bnk->Xold);CHKERRQ(ierr);
180           /* Compute the objective at the trial */
181           ierr = TaoComputeObjective(tao, tao->solution, &ftrial);CHKERRQ(ierr);
182           if (PetscIsInfOrNanReal(bnk->f)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
183           ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr);
184           if (PetscIsInfOrNanReal(ftrial)) {
185             tau = bnk->gamma1_i;
186           } else {
187             if (ftrial < f_min) {
188               f_min = ftrial;
189               sigma = -tao->trust / bnk->gnorm;
190             }
191 
192             /* Compute the predicted and actual reduction */
193             if (bnk->active_idx) {
194               ierr = VecGetSubVector(bnk->W, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr);
195               ierr = VecGetSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work);CHKERRQ(ierr);
196             } else {
197               bnk->X_inactive = bnk->W;
198               bnk->inactive_work = bnk->Xwork;
199             }
200             ierr = MatMult(bnk->H_inactive, bnk->X_inactive, bnk->inactive_work);CHKERRQ(ierr);
201             ierr = VecDot(bnk->X_inactive, bnk->inactive_work, &prered);CHKERRQ(ierr);
202             if (bnk->active_idx) {
203               ierr = VecRestoreSubVector(bnk->W, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr);
204               ierr = VecRestoreSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work);CHKERRQ(ierr);
205             }
206             prered = tao->trust * (bnk->gnorm - 0.5 * tao->trust * prered / (bnk->gnorm * bnk->gnorm));
207             actred = bnk->f - ftrial;
208             if ((PetscAbsScalar(actred) <= bnk->epsilon) && (PetscAbsScalar(prered) <= bnk->epsilon)) {
209               kappa = 1.0;
210             } else {
211               kappa = actred / prered;
212             }
213 
214             tau_1 = bnk->theta_i * bnk->gnorm * tao->trust / (bnk->theta_i * bnk->gnorm * tao->trust + (1.0 - bnk->theta_i) * prered - actred);
215             tau_2 = bnk->theta_i * bnk->gnorm * tao->trust / (bnk->theta_i * bnk->gnorm * tao->trust - (1.0 + bnk->theta_i) * prered + actred);
216             tau_min = PetscMin(tau_1, tau_2);
217             tau_max = PetscMax(tau_1, tau_2);
218 
219             if (PetscAbsScalar(kappa - 1.0) <= bnk->mu1_i) {
220               /*  Great agreement */
221               max_radius = PetscMax(max_radius, tao->trust);
222 
223               if (tau_max < 1.0) {
224                 tau = bnk->gamma3_i;
225               } else if (tau_max > bnk->gamma4_i) {
226                 tau = bnk->gamma4_i;
227               } else {
228                 tau = tau_max;
229               }
230             } else if (PetscAbsScalar(kappa - 1.0) <= bnk->mu2_i) {
231               /*  Good agreement */
232               max_radius = PetscMax(max_radius, tao->trust);
233 
234               if (tau_max < bnk->gamma2_i) {
235                 tau = bnk->gamma2_i;
236               } else if (tau_max > bnk->gamma3_i) {
237                 tau = bnk->gamma3_i;
238               } else {
239                 tau = tau_max;
240               }
241             } else {
242               /*  Not good agreement */
243               if (tau_min > 1.0) {
244                 tau = bnk->gamma2_i;
245               } else if (tau_max < bnk->gamma1_i) {
246                 tau = bnk->gamma1_i;
247               } else if ((tau_min < bnk->gamma1_i) && (tau_max >= 1.0)) {
248                 tau = bnk->gamma1_i;
249               } else if ((tau_1 >= bnk->gamma1_i) && (tau_1 < 1.0) && ((tau_2 < bnk->gamma1_i) || (tau_2 >= 1.0))) {
250                 tau = tau_1;
251               } else if ((tau_2 >= bnk->gamma1_i) && (tau_2 < 1.0) && ((tau_1 < bnk->gamma1_i) || (tau_2 >= 1.0))) {
252                 tau = tau_2;
253               } else {
254                 tau = tau_max;
255               }
256             }
257           }
258           tao->trust = tau * tao->trust;
259         }
260 
261         if (f_min < bnk->f) {
262           /* We accidentally found a solution better than the initial, so accept it */
263           bnk->f = f_min;
264           ierr = VecCopy(tao->solution, bnk->Xold);CHKERRQ(ierr);
265           ierr = VecAXPY(tao->solution,sigma,tao->gradient);CHKERRQ(ierr);
266           ierr = TaoBoundSolution(tao->solution, tao->XL,tao->XU, 0.0, &nDiff, tao->solution);CHKERRQ(ierr);
267           ierr = VecCopy(tao->solution, tao->stepdirection);CHKERRQ(ierr);
268           ierr = VecAXPY(tao->stepdirection, -1.0, bnk->Xold);CHKERRQ(ierr);
269           ierr = TaoComputeGradient(tao,tao->solution,bnk->unprojected_gradient);CHKERRQ(ierr);
270           ierr = TaoBNKEstimateActiveSet(tao, bnk->as_type);CHKERRQ(ierr);
271           ierr = VecCopy(bnk->unprojected_gradient, tao->gradient);CHKERRQ(ierr);
272           ierr = VecISSet(tao->gradient, bnk->active_idx, 0.0);CHKERRQ(ierr);
273           /* Compute gradient at the new iterate and flip switch to compute the Hessian later */
274           ierr = VecNorm(tao->gradient, NORM_2, &bnk->gnorm);CHKERRQ(ierr);
275           *needH = PETSC_TRUE;
276           /* Test the new step for convergence */
277           ierr = VecFischer(tao->solution, bnk->unprojected_gradient, tao->XL, tao->XU, bnk->W);CHKERRQ(ierr);
278           ierr = VecNorm(bnk->W, NORM_2, &resnorm);CHKERRQ(ierr);
279           if (PetscIsInfOrNanReal(resnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
280           ierr = TaoLogConvergenceHistory(tao,bnk->f,resnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
281           ierr = TaoMonitor(tao,tao->niter,bnk->f,resnorm,0.0,1.0);CHKERRQ(ierr);
282           ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
283           if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
284           /* active BNCG recycling early because we have a stepdirection computed */
285           ierr = TaoBNCGSetRecycleFlag(bnk->bncg, PETSC_TRUE);CHKERRQ(ierr);
286         }
287       }
288       tao->trust = PetscMax(tao->trust, max_radius);
289 
290       /* Ensure that the trust radius is within the limits */
291       tao->trust = PetscMax(tao->trust, bnk->min_radius);
292       tao->trust = PetscMin(tao->trust, bnk->max_radius);
293       break;
294 
295     default:
296       /* Norm of the first direction will initialize radius */
297       tao->trust = 0.0;
298       break;
299     }
300   }
301 
302   /* Set initial scaling for the BFGS preconditioner
303      This step is done after computing the initial trust-region radius
304      since the function value may have decreased */
305   if (BNK_PC_BFGS == bnk->pc_type) {
306     delta = 2.0 * PetscMax(1.0, PetscAbsScalar(bnk->f)) / (bnk->gnorm*bnk->gnorm);
307     ierr = MatLMVMSetDelta(bnk->M,delta);CHKERRQ(ierr);
308   }
309   PetscFunctionReturn(0);
310 }
311 
312 /*------------------------------------------------------------*/
313 
314 /* Routine for computing the Hessian and preparing the preconditioner at the new iterate */
315 
316 PetscErrorCode TaoBNKComputeHessian(Tao tao)
317 {
318   PetscErrorCode               ierr;
319   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
320   PetscBool                    diagExists;
321   PetscReal                    delta;
322 
323   PetscFunctionBegin;
324   /* Compute the Hessian */
325   ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr);
326   /* Add a correction to the BFGS preconditioner */
327   if (BNK_PC_BFGS == bnk->pc_type && bnk->M) {
328     ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
329     /* Update the BFGS diagonal scaling */
330     ierr = MatHasOperation(tao->hessian, MATOP_GET_DIAGONAL, &diagExists);CHKERRQ(ierr);
331     if (BFGS_SCALE_AHESS == bnk->bfgs_scale_type && diagExists) {
332       ierr = MatGetDiagonal(tao->hessian, bnk->Diag);CHKERRQ(ierr);
333       ierr = VecAbs(bnk->Diag);CHKERRQ(ierr);
334       ierr = VecMedian(bnk->Diag_min, bnk->Diag, bnk->Diag_max, bnk->Diag);CHKERRQ(ierr);
335       ierr = VecReciprocal(bnk->Diag);CHKERRQ(ierr);
336       ierr = MatLMVMSetScale(bnk->M,bnk->Diag);CHKERRQ(ierr);
337     } else {
338       /* Fall back onto stand-alone BFGS scaling */
339       delta = 2.0 * PetscMax(1.0, PetscAbsScalar(bnk->f)) / (bnk->gnorm*bnk->gnorm);
340       ierr = MatLMVMSetDelta(bnk->M,delta);CHKERRQ(ierr);
341     }
342   }
343   PetscFunctionReturn(0);
344 }
345 
346 /*------------------------------------------------------------*/
347 
348 /* Routine for estimating the active set */
349 
350 PetscErrorCode TaoBNKEstimateActiveSet(Tao tao, PetscInt asType)
351 {
352   PetscErrorCode               ierr;
353   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
354   PetscBool                    hessComputed, diagExists;
355 
356   PetscFunctionBegin;
357   switch (asType) {
358   case BNK_AS_NONE:
359     ierr = ISDestroy(&bnk->inactive_idx);CHKERRQ(ierr);
360     ierr = VecWhichInactive(tao->XL, tao->solution, bnk->unprojected_gradient, tao->XU, PETSC_TRUE, &bnk->inactive_idx);CHKERRQ(ierr);
361     ierr = ISDestroy(&bnk->active_idx);CHKERRQ(ierr);
362     ierr = ISComplementVec(bnk->inactive_idx, tao->solution, &bnk->active_idx);CHKERRQ(ierr);
363     break;
364 
365   case BNK_AS_BERTSEKAS:
366     /* Compute the trial step vector with which we will estimate the active set at the next iteration */
367     if (BNK_PC_BFGS == bnk->pc_type && bnk->M) {
368       /* If the BFGS preconditioner matrix is available, we will construct a trial step with it */
369       ierr = MatLMVMSetInactive(bnk->M, NULL);CHKERRQ(ierr);
370       ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, bnk->W);CHKERRQ(ierr);
371     } else {
372       ierr = MatAssembled(tao->hessian, &hessComputed);CHKERRQ(ierr);
373       ierr = MatHasOperation(tao->hessian, MATOP_GET_DIAGONAL, &diagExists);CHKERRQ(ierr);
374       if (hessComputed && diagExists) {
375         /* BFGS preconditioner doesn't exist so let's invert the absolute diagonal of the Hessian instead onto the gradient */
376         ierr = MatGetDiagonal(tao->hessian, bnk->Xwork);CHKERRQ(ierr);
377         ierr = VecAbs(bnk->Xwork);CHKERRQ(ierr);
378         ierr = VecMedian(bnk->Diag_min, bnk->Xwork, bnk->Diag_max, bnk->Xwork);CHKERRQ(ierr);
379         ierr = VecReciprocal(bnk->Xwork);CHKERRQ(ierr);CHKERRQ(ierr);
380         ierr = VecPointwiseMult(bnk->W, bnk->Xwork, bnk->unprojected_gradient);CHKERRQ(ierr);
381       } else {
382         /* If the Hessian or its diagonal does not exist, we will simply use gradient step */
383         ierr = VecCopy(bnk->unprojected_gradient, bnk->W);CHKERRQ(ierr);
384       }
385     }
386     ierr = VecScale(bnk->W, -1.0);CHKERRQ(ierr);
387     ierr = TaoEstimateActiveBounds(tao->solution, tao->XL, tao->XU, bnk->unprojected_gradient, bnk->W, bnk->Xwork, bnk->as_step, &bnk->as_tol,
388                                    &bnk->active_lower, &bnk->active_upper, &bnk->active_fixed, &bnk->active_idx, &bnk->inactive_idx);CHKERRQ(ierr);
389     break;
390 
391   default:
392     break;
393   }
394   PetscFunctionReturn(0);
395 }
396 
397 /*------------------------------------------------------------*/
398 
399 /* Routine for bounding the step direction */
400 
401 PetscErrorCode TaoBNKBoundStep(Tao tao, PetscInt asType, Vec step)
402 {
403   PetscErrorCode               ierr;
404   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
405 
406   PetscFunctionBegin;
407   switch (asType) {
408   case BNK_AS_NONE:
409     ierr = VecISSet(step, bnk->active_idx, 0.0);CHKERRQ(ierr);
410     break;
411 
412   case BNK_AS_BERTSEKAS:
413     ierr = TaoBoundStep(tao->solution, tao->XL, tao->XU, bnk->active_lower, bnk->active_upper, bnk->active_fixed, 1.0, step);CHKERRQ(ierr);
414     break;
415 
416   default:
417     break;
418   }
419   PetscFunctionReturn(0);
420 }
421 
422 /*------------------------------------------------------------*/
423 
424 /* Routine for taking a finite number of BNCG iterations to
425    accelerate Newton convergence.
426 
427    In practice, this approach simply trades off Hessian evaluations
428    for more gradient evaluations.
429 */
430 
431 PetscErrorCode TaoBNKTakeCGSteps(Tao tao, PetscBool *terminate)
432 {
433   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
434   PetscErrorCode               ierr;
435 
436   PetscFunctionBegin;
437   *terminate = PETSC_FALSE;
438   if (bnk->max_cg_its > 0) {
439     /* Copy the current function value (important vectors are already shared) */
440     bnk->bncg_ctx->f = bnk->f;
441     /* Take some small finite number of BNCG iterations */
442     ierr = TaoSolve(bnk->bncg);CHKERRQ(ierr);
443     /* Add the number of gradient and function evaluations to the total */
444     tao->nfuncs += bnk->bncg->nfuncs;
445     tao->nfuncgrads += bnk->bncg->nfuncgrads;
446     tao->ngrads += bnk->bncg->ngrads;
447     tao->nhess += bnk->bncg->nhess;
448     bnk->tot_cg_its += bnk->bncg->niter;
449     /* Extract the BNCG function value out and save it into BNK */
450     bnk->f = bnk->bncg_ctx->f;
451     if (bnk->bncg->reason == TAO_CONVERGED_GATOL || bnk->bncg->reason == TAO_CONVERGED_GRTOL || bnk->bncg->reason == TAO_CONVERGED_GTTOL || bnk->bncg->reason == TAO_CONVERGED_MINF) {
452       *terminate = PETSC_TRUE;
453     } else {
454       ierr = TaoBNKEstimateActiveSet(tao, bnk->as_type);CHKERRQ(ierr);
455     }
456   }
457   PetscFunctionReturn(0);
458 }
459 
460 /*------------------------------------------------------------*/
461 
462 /* Routine for computing the Newton step. */
463 
464 PetscErrorCode TaoBNKComputeStep(Tao tao, PetscBool shift, KSPConvergedReason *ksp_reason)
465 {
466   PetscErrorCode               ierr;
467   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
468 
469   PetscReal                    delta;
470   PetscInt                     bfgsUpdates = 0;
471   PetscInt                     kspits;
472   PetscBool                    diagExists;
473 
474   PetscFunctionBegin;
475   /* If there are no inactive variables left, save some computation and return an adjusted zero step
476      that has (l-x) and (u-x) for lower and upper bounded variables. */
477   if (!bnk->inactive_idx) {
478     ierr = VecSet(tao->stepdirection, 0.0);CHKERRQ(ierr);
479     ierr = TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection);CHKERRQ(ierr);
480     PetscFunctionReturn(0);
481   }
482 
483   /* Prepare the reduced sub-matrices for the inactive set */
484   if (BNK_PC_BFGS == bnk->pc_type) {
485     ierr = MatLMVMSetInactive(bnk->M, bnk->inactive_idx);CHKERRQ(ierr);
486   }
487   if (bnk->active_idx) {
488     ierr = MatDestroy(&bnk->H_inactive);CHKERRQ(ierr);
489     ierr = MatCreateSubMatrix(tao->hessian, bnk->inactive_idx, bnk->inactive_idx, MAT_INITIAL_MATRIX, &bnk->H_inactive);CHKERRQ(ierr);
490     if (tao->hessian == tao->hessian_pre) {
491       bnk->Hpre_inactive = bnk->H_inactive;
492     } else {
493       ierr = MatDestroy(&bnk->Hpre_inactive);CHKERRQ(ierr);
494       ierr = MatCreateSubMatrix(tao->hessian_pre, bnk->inactive_idx, bnk->inactive_idx, MAT_INITIAL_MATRIX, &bnk->Hpre_inactive);CHKERRQ(ierr);
495     }
496   } else {
497     ierr = MatDestroy(&bnk->H_inactive);CHKERRQ(ierr);
498     ierr = MatDuplicate(tao->hessian, MAT_COPY_VALUES, &bnk->H_inactive);CHKERRQ(ierr);
499     if (tao->hessian == tao->hessian_pre) {
500       bnk->Hpre_inactive = bnk->H_inactive;
501     } else {
502       ierr = MatDestroy(&bnk->Hpre_inactive);CHKERRQ(ierr);
503       ierr = MatDuplicate(tao->hessian_pre, MAT_COPY_VALUES, &bnk->Hpre_inactive);CHKERRQ(ierr);
504     }
505   }
506 
507   /* Shift the reduced Hessian matrix */
508   if ((shift) && (bnk->pert > 0)) {
509     ierr = MatShift(bnk->H_inactive, bnk->pert);CHKERRQ(ierr);
510     if (bnk->H_inactive != bnk->Hpre_inactive) {
511       ierr = MatShift(bnk->Hpre_inactive, bnk->pert);CHKERRQ(ierr);
512     }
513   }
514 
515   /* Update the diagonal scaling for the BFGS preconditioner, this time with the Hessian perturbation */
516   if ((BNK_PC_BFGS == bnk->pc_type) && (BFGS_SCALE_PHESS == bnk->bfgs_scale_type)) {
517     ierr = MatHasOperation(bnk->H_inactive, MATOP_GET_DIAGONAL, &diagExists);CHKERRQ(ierr);
518     if (diagExists) {
519       /* Obtain diagonal for the bfgs preconditioner  */
520       ierr = VecSet(bnk->Diag, 1.0);CHKERRQ(ierr);
521       if (bnk->active_idx) {
522         ierr = VecGetSubVector(bnk->Diag, bnk->inactive_idx, &bnk->Diag_red);CHKERRQ(ierr);
523       } else {
524         bnk->Diag_red = bnk->Diag;
525       }
526       ierr = MatGetDiagonal(bnk->H_inactive, bnk->Diag_red);CHKERRQ(ierr);
527       if (bnk->active_idx) {
528         ierr = VecRestoreSubVector(bnk->Diag, bnk->inactive_idx, &bnk->Diag_red);CHKERRQ(ierr);
529       }
530       ierr = VecAbs(bnk->Diag);CHKERRQ(ierr);
531       ierr = VecMedian(bnk->Diag_min, bnk->Diag, bnk->Diag_max, bnk->Diag);CHKERRQ(ierr);
532       ierr = VecReciprocal(bnk->Diag);CHKERRQ(ierr);
533       ierr = MatLMVMSetScale(bnk->M,bnk->Diag);CHKERRQ(ierr);
534     } else {
535       /* Fall back onto stand-alone BFGS scaling */
536       delta = 2.0 * PetscMax(1.0, PetscAbsScalar(bnk->f)) / (bnk->gnorm*bnk->gnorm);
537       ierr = MatLMVMSetDelta(bnk->M,delta);CHKERRQ(ierr);
538     }
539   }
540 
541   /* Solve the Newton system of equations */
542   tao->ksp_its = 0;
543   ierr = VecSet(tao->stepdirection, 0.0);CHKERRQ(ierr);
544   ierr = KSPReset(tao->ksp);CHKERRQ(ierr);
545   ierr = KSPSetOperators(tao->ksp,bnk->H_inactive,bnk->Hpre_inactive);CHKERRQ(ierr);
546   ierr = VecCopy(bnk->unprojected_gradient, bnk->Gwork);CHKERRQ(ierr);
547   if (bnk->active_idx) {
548     ierr = VecGetSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive);CHKERRQ(ierr);
549     ierr = VecGetSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr);
550   } else {
551     bnk->G_inactive = bnk->unprojected_gradient;
552     bnk->X_inactive = tao->stepdirection;
553   }
554   if (bnk->is_nash || bnk->is_stcg || bnk->is_gltr) {
555     ierr = KSPCGSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr);
556     ierr = KSPSolve(tao->ksp, bnk->G_inactive, bnk->X_inactive);CHKERRQ(ierr);
557     ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr);
558     tao->ksp_its+=kspits;
559     tao->ksp_tot_its+=kspits;
560     ierr = KSPCGGetNormD(tao->ksp,&bnk->dnorm);CHKERRQ(ierr);
561 
562     if (0.0 == tao->trust) {
563       /* Radius was uninitialized; use the norm of the direction */
564       if (bnk->dnorm > 0.0) {
565         tao->trust = bnk->dnorm;
566 
567         /* Modify the radius if it is too large or small */
568         tao->trust = PetscMax(tao->trust, bnk->min_radius);
569         tao->trust = PetscMin(tao->trust, bnk->max_radius);
570       } else {
571         /* The direction was bad; set radius to default value and re-solve
572            the trust-region subproblem to get a direction */
573         tao->trust = tao->trust0;
574 
575         /* Modify the radius if it is too large or small */
576         tao->trust = PetscMax(tao->trust, bnk->min_radius);
577         tao->trust = PetscMin(tao->trust, bnk->max_radius);
578 
579         ierr = KSPCGSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr);
580         ierr = KSPSolve(tao->ksp, bnk->G_inactive, bnk->X_inactive);CHKERRQ(ierr);
581         ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr);
582         tao->ksp_its+=kspits;
583         tao->ksp_tot_its+=kspits;
584         ierr = KSPCGGetNormD(tao->ksp,&bnk->dnorm);CHKERRQ(ierr);
585 
586         if (bnk->dnorm == 0.0) SETERRQ(PETSC_COMM_SELF,1, "Initial direction zero");
587       }
588     }
589   } else {
590     ierr = KSPSolve(tao->ksp, bnk->G_inactive, bnk->X_inactive);CHKERRQ(ierr);
591     ierr = KSPGetIterationNumber(tao->ksp, &kspits);CHKERRQ(ierr);
592     tao->ksp_its += kspits;
593     tao->ksp_tot_its+=kspits;
594   }
595   /* Restore sub vectors back */
596   if (bnk->active_idx) {
597     ierr = VecRestoreSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive);CHKERRQ(ierr);
598     ierr = VecRestoreSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr);
599   }
600   /* Make sure the safeguarded fall-back step is zero for actively bounded variables */
601   ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
602   ierr = TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection);CHKERRQ(ierr);
603 
604   /* Record convergence reasons */
605   ierr = KSPGetConvergedReason(tao->ksp, ksp_reason);CHKERRQ(ierr);
606   if (KSP_CONVERGED_ATOL == *ksp_reason) {
607     ++bnk->ksp_atol;
608   } else if (KSP_CONVERGED_RTOL == *ksp_reason) {
609     ++bnk->ksp_rtol;
610   } else if (KSP_CONVERGED_CG_CONSTRAINED == *ksp_reason) {
611     ++bnk->ksp_ctol;
612   } else if (KSP_CONVERGED_CG_NEG_CURVE == *ksp_reason) {
613     ++bnk->ksp_negc;
614   } else if (KSP_DIVERGED_DTOL == *ksp_reason) {
615     ++bnk->ksp_dtol;
616   } else if (KSP_DIVERGED_ITS == *ksp_reason) {
617     ++bnk->ksp_iter;
618   } else {
619     ++bnk->ksp_othr;
620   }
621 
622   /* Make sure the BFGS preconditioner is healthy */
623   if (bnk->pc_type == BNK_PC_BFGS) {
624     ierr = MatLMVMGetUpdates(bnk->M, &bfgsUpdates);CHKERRQ(ierr);
625     if ((KSP_DIVERGED_INDEFINITE_PC == *ksp_reason) && (bfgsUpdates > 1)) {
626       /* Preconditioner is numerically indefinite; reset the approximation. */
627       delta = 2.0 * PetscMax(1.0, PetscAbsScalar(bnk->f)) / (bnk->gnorm*bnk->gnorm);
628       ierr = MatLMVMSetDelta(bnk->M,delta);CHKERRQ(ierr);
629       ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr);
630       ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
631     }
632   }
633   PetscFunctionReturn(0);
634 }
635 
636 /*------------------------------------------------------------*/
637 
638 /* Routine for recomputing the predicted reduction for a given step vector */
639 
640 PetscErrorCode TaoBNKRecomputePred(Tao tao, Vec S, PetscReal *prered)
641 {
642   PetscErrorCode               ierr;
643   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
644 
645   PetscFunctionBegin;
646   /* Extract subvectors associated with the inactive set */
647   if (bnk->active_idx){
648     ierr = VecGetSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr);
649     ierr = VecGetSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work);CHKERRQ(ierr);
650     ierr = VecGetSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive);CHKERRQ(ierr);
651   } else {
652     bnk->X_inactive = tao->stepdirection;
653     bnk->inactive_work = bnk->Xwork;
654     bnk->G_inactive = bnk->Gwork;
655   }
656   /* Recompute the predicted decrease based on the quadratic model */
657   ierr = MatMult(bnk->H_inactive, bnk->X_inactive, bnk->inactive_work);CHKERRQ(ierr);
658   ierr = VecAYPX(bnk->inactive_work, -0.5, bnk->G_inactive);CHKERRQ(ierr);
659   ierr = VecDot(bnk->inactive_work, bnk->X_inactive, prered);CHKERRQ(ierr);
660   /* Restore the sub vectors */
661   if (bnk->active_idx){
662     ierr = VecRestoreSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr);
663     ierr = VecRestoreSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work);CHKERRQ(ierr);
664     ierr = VecRestoreSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive);CHKERRQ(ierr);
665   }
666   PetscFunctionReturn(0);
667 }
668 
669 /*------------------------------------------------------------*/
670 
671 /* Routine for ensuring that the Newton step is a descent direction.
672 
673    The step direction falls back onto BFGS, scaled gradient and gradient steps
674    in the event that the Newton step fails the test.
675 */
676 
677 PetscErrorCode TaoBNKSafeguardStep(Tao tao, KSPConvergedReason ksp_reason, PetscInt *stepType)
678 {
679   PetscErrorCode               ierr;
680   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
681 
682   PetscReal                    gdx, delta, e_min;
683   PetscInt                     bfgsUpdates;
684 
685   PetscFunctionBegin;
686   ierr = VecDot(tao->stepdirection, tao->gradient, &gdx);CHKERRQ(ierr);
687   if ((gdx >= 0.0) || PetscIsInfOrNanReal(gdx)) {
688     /* Newton step is not descent or direction produced Inf or NaN
689        Update the perturbation for next time */
690     if (bnk->pert <= 0.0) {
691       /* Initialize the perturbation */
692       bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
693       if (bnk->is_gltr) {
694         ierr = KSPCGGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr);
695         bnk->pert = PetscMax(bnk->pert, -e_min);
696       }
697     } else {
698       /* Increase the perturbation */
699       bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
700     }
701 
702     if (BNK_PC_BFGS != bnk->pc_type) {
703       /* We don't have the bfgs matrix around and updated
704          Must use gradient direction in this case */
705       ierr = VecCopy(tao->gradient, tao->stepdirection);CHKERRQ(ierr);
706       *stepType = BNK_GRADIENT;
707     } else {
708       /* Attempt to use the BFGS direction */
709       ierr = MatLMVMSetInactive(bnk->M, NULL);CHKERRQ(ierr);
710       ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
711 
712       /* Check for success (descent direction)
713          NOTE: Negative gdx here means not a descent direction because
714          the fall-back step is missing a negative sign. */
715       ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr);
716       if ((gdx <= 0.0) || PetscIsInfOrNanReal(gdx)) {
717         /* BFGS direction is not descent or direction produced not a number
718            We can assert bfgsUpdates > 1 in this case because
719            the first solve produces the scaled gradient direction,
720            which is guaranteed to be descent */
721 
722         /* Use steepest descent direction (scaled) */
723         delta = 2.0 * PetscMax(1.0, PetscAbsScalar(bnk->f)) / (bnk->gnorm*bnk->gnorm);
724         ierr = MatLMVMSetDelta(bnk->M, delta);CHKERRQ(ierr);
725         ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr);
726         ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
727         ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
728 
729         *stepType = BNK_SCALED_GRADIENT;
730       } else {
731         ierr = MatLMVMGetUpdates(bnk->M, &bfgsUpdates);CHKERRQ(ierr);
732         if (1 == bfgsUpdates) {
733           /* The first BFGS direction is always the scaled gradient */
734           *stepType = BNK_SCALED_GRADIENT;
735         } else {
736           *stepType = BNK_BFGS;
737         }
738       }
739     }
740     /* Make sure the safeguarded fall-back step is zero for actively bounded variables */
741     ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
742     ierr = TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection);CHKERRQ(ierr);
743   } else {
744     /* Computed Newton step is descent */
745     switch (ksp_reason) {
746     case KSP_DIVERGED_NANORINF:
747     case KSP_DIVERGED_BREAKDOWN:
748     case KSP_DIVERGED_INDEFINITE_MAT:
749     case KSP_DIVERGED_INDEFINITE_PC:
750     case KSP_CONVERGED_CG_NEG_CURVE:
751       /* Matrix or preconditioner is indefinite; increase perturbation */
752       if (bnk->pert <= 0.0) {
753         /* Initialize the perturbation */
754         bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
755         if (bnk->is_gltr) {
756           ierr = KSPCGGLTRGetMinEig(tao->ksp, &e_min);CHKERRQ(ierr);
757           bnk->pert = PetscMax(bnk->pert, -e_min);
758         }
759       } else {
760         /* Increase the perturbation */
761         bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
762       }
763       break;
764 
765     default:
766       /* Newton step computation is good; decrease perturbation */
767       bnk->pert = PetscMin(bnk->psfac * bnk->pert, bnk->pmsfac * bnk->gnorm);
768       if (bnk->pert < bnk->pmin) {
769         bnk->pert = 0.0;
770       }
771       break;
772     }
773     *stepType = BNK_NEWTON;
774   }
775   PetscFunctionReturn(0);
776 }
777 
778 /*------------------------------------------------------------*/
779 
780 /* Routine for performing a bound-projected More-Thuente line search.
781 
782   Includes fallbacks to BFGS, scaled gradient, and unscaled gradient steps if the
783   Newton step does not produce a valid step length.
784 */
785 
786 PetscErrorCode TaoBNKPerformLineSearch(Tao tao, PetscInt *stepType, PetscReal *steplen, TaoLineSearchConvergedReason *reason)
787 {
788   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
789   PetscErrorCode ierr;
790   TaoLineSearchConvergedReason ls_reason;
791 
792   PetscReal      e_min, gdx, delta;
793   PetscInt       bfgsUpdates;
794 
795   PetscFunctionBegin;
796   /* Perform the linesearch */
797   ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, steplen, &ls_reason);CHKERRQ(ierr);
798   ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
799 
800   while (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER && *stepType != BNK_GRADIENT) {
801     /* Linesearch failed, revert solution */
802     bnk->f = bnk->fold;
803     ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr);
804     ierr = VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient);CHKERRQ(ierr);
805 
806     switch(*stepType) {
807     case BNK_NEWTON:
808       /* Failed to obtain acceptable iterate with Newton step
809          Update the perturbation for next time */
810       if (bnk->pert <= 0.0) {
811         /* Initialize the perturbation */
812         bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
813         if (bnk->is_gltr) {
814           ierr = KSPCGGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr);
815           bnk->pert = PetscMax(bnk->pert, -e_min);
816         }
817       } else {
818         /* Increase the perturbation */
819         bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
820       }
821 
822       if (BNK_PC_BFGS != bnk->pc_type) {
823         /* We don't have the bfgs matrix around and being updated
824            Must use gradient direction in this case */
825         ierr = VecCopy(bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
826         *stepType = BNK_GRADIENT;
827       } else {
828         /* Attempt to use the BFGS direction */
829         ierr = MatLMVMSetInactive(bnk->M, NULL);CHKERRQ(ierr);
830         ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
831         /* Check for success (descent direction)
832            NOTE: Negative gdx means not a descent direction because the step here is missing a negative sign. */
833         ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr);
834         if ((gdx <= 0.0) || PetscIsInfOrNanReal(gdx)) {
835           /* BFGS direction is not descent or direction produced not a number
836              We can assert bfgsUpdates > 1 in this case
837              Use steepest descent direction (scaled) */
838           delta = 2.0 * PetscMax(1.0, PetscAbsScalar(bnk->f)) / (bnk->gnorm*bnk->gnorm);
839           ierr = MatLMVMSetDelta(bnk->M, delta);CHKERRQ(ierr);
840           ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr);
841           ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
842           ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
843 
844           bfgsUpdates = 1;
845           *stepType = BNK_SCALED_GRADIENT;
846         } else {
847           ierr = MatLMVMGetUpdates(bnk->M, &bfgsUpdates);CHKERRQ(ierr);
848           if (1 == bfgsUpdates) {
849             /* The first BFGS direction is always the scaled gradient */
850             *stepType = BNK_SCALED_GRADIENT;
851           } else {
852             *stepType = BNK_BFGS;
853           }
854         }
855       }
856       break;
857 
858     case BNK_BFGS:
859       /* Can only enter if pc_type == BNK_PC_BFGS
860          Failed to obtain acceptable iterate with BFGS step
861          Attempt to use the scaled gradient direction */
862       delta = 2.0 * PetscMax(1.0, PetscAbsScalar(bnk->f)) / (bnk->gnorm*bnk->gnorm);
863       ierr = MatLMVMSetDelta(bnk->M, delta);CHKERRQ(ierr);
864       ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr);
865       ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
866       ierr = MatLMVMSetInactive(bnk->M, NULL);CHKERRQ(ierr);
867       ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
868 
869       bfgsUpdates = 1;
870       *stepType = BNK_SCALED_GRADIENT;
871       break;
872 
873     case BNK_SCALED_GRADIENT:
874       /* Can only enter if pc_type == BNK_PC_BFGS
875          The scaled gradient step did not produce a new iterate;
876          reset the BFGS matrix and attemp to use the gradient direction. */
877       ierr = MatLMVMSetScale(bnk->M,0);CHKERRQ(ierr);
878       ierr = MatLMVMSetDelta(bnk->M,1.0);CHKERRQ(ierr);
879       ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr);
880       ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
881       ierr = VecCopy(bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
882 
883       bfgsUpdates = 1;
884       *stepType = BNK_GRADIENT;
885       break;
886     }
887     /* Make sure the safeguarded fall-back step is zero for actively bounded variables */
888     ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
889     ierr = TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection);CHKERRQ(ierr);
890 
891     /* Perform one last line search with the fall-back step */
892     ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, steplen, &ls_reason);CHKERRQ(ierr);
893     ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
894   }
895   *reason = ls_reason;
896   PetscFunctionReturn(0);
897 }
898 
899 /*------------------------------------------------------------*/
900 
901 /* Routine for updating the trust radius.
902 
903   Function features three different update methods:
904   1) Line-search step length based
905   2) Predicted decrease on the CG quadratic model
906   3) Interpolation
907 */
908 
909 PetscErrorCode TaoBNKUpdateTrustRadius(Tao tao, PetscReal prered, PetscReal actred, PetscInt updateType, PetscInt stepType, PetscBool *accept)
910 {
911   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
912   PetscErrorCode ierr;
913 
914   PetscReal      step, kappa;
915   PetscReal      gdx, tau_1, tau_2, tau_min, tau_max;
916 
917   PetscFunctionBegin;
918   /* Update trust region radius */
919   *accept = PETSC_FALSE;
920   switch(updateType) {
921   case BNK_UPDATE_STEP:
922     *accept = PETSC_TRUE; /* always accept here because line search succeeded */
923     if (stepType == BNK_NEWTON) {
924       ierr = TaoLineSearchGetStepLength(tao->linesearch, &step);CHKERRQ(ierr);
925       if (step < bnk->nu1) {
926         /* Very bad step taken; reduce radius */
927         tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust);
928       } else if (step < bnk->nu2) {
929         /* Reasonably bad step taken; reduce radius */
930         tao->trust = bnk->omega2 * PetscMin(bnk->dnorm, tao->trust);
931       } else if (step < bnk->nu3) {
932         /*  Reasonable step was taken; leave radius alone */
933         if (bnk->omega3 < 1.0) {
934           tao->trust = bnk->omega3 * PetscMin(bnk->dnorm, tao->trust);
935         } else if (bnk->omega3 > 1.0) {
936           tao->trust = PetscMax(bnk->omega3 * bnk->dnorm, tao->trust);
937         }
938       } else if (step < bnk->nu4) {
939         /*  Full step taken; increase the radius */
940         tao->trust = PetscMax(bnk->omega4 * bnk->dnorm, tao->trust);
941       } else {
942         /*  More than full step taken; increase the radius */
943         tao->trust = PetscMax(bnk->omega5 * bnk->dnorm, tao->trust);
944       }
945     } else {
946       /*  Newton step was not good; reduce the radius */
947       tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust);
948     }
949     break;
950 
951   case BNK_UPDATE_REDUCTION:
952     if (stepType == BNK_NEWTON) {
953       if (prered < 0.0) {
954         /* The predicted reduction has the wrong sign.  This cannot
955            happen in infinite precision arithmetic.  Step should
956            be rejected! */
957         tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
958       } else {
959         if (PetscIsInfOrNanReal(actred)) {
960           tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
961         } else {
962           if ((PetscAbsScalar(actred) <= PetscMax(1.0, PetscAbsScalar(bnk->f))*bnk->epsilon) && (PetscAbsScalar(prered) <= PetscMax(1.0, PetscAbsScalar(bnk->f))*bnk->epsilon)) {
963             kappa = 1.0;
964           } else {
965             kappa = actred / prered;
966           }
967 
968           /* Accept or reject the step and update radius */
969           if (kappa < bnk->eta1) {
970             /* Reject the step */
971             tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
972           } else {
973             /* Accept the step */
974             *accept = PETSC_TRUE;
975             /* Update the trust region radius only if the computed step is at the trust radius boundary */
976             if (bnk->dnorm == tao->trust) {
977               if (kappa < bnk->eta2) {
978                 /* Marginal bad step */
979                 tao->trust = bnk->alpha2 * tao->trust;
980               } else if (kappa < bnk->eta3) {
981                 /* Reasonable step */
982                 tao->trust = bnk->alpha3 * tao->trust;
983               } else if (kappa < bnk->eta4) {
984                 /* Good step */
985                 tao->trust = bnk->alpha4 * tao->trust;
986               } else {
987                 /* Very good step */
988                 tao->trust = bnk->alpha5 * tao->trust;
989               }
990             }
991           }
992         }
993       }
994     } else {
995       /*  Newton step was not good; reduce the radius */
996       tao->trust = bnk->alpha1 * PetscMin(bnk->dnorm, tao->trust);
997     }
998     break;
999 
1000   default:
1001     if (stepType == BNK_NEWTON) {
1002       if (prered < 0.0) {
1003         /*  The predicted reduction has the wrong sign.  This cannot */
1004         /*  happen in infinite precision arithmetic.  Step should */
1005         /*  be rejected! */
1006         tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
1007       } else {
1008         if (PetscIsInfOrNanReal(actred)) {
1009           tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
1010         } else {
1011           if ((PetscAbsScalar(actred) <= bnk->epsilon) && (PetscAbsScalar(prered) <= bnk->epsilon)) {
1012             kappa = 1.0;
1013           } else {
1014             kappa = actred / prered;
1015           }
1016 
1017           ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr);
1018           tau_1 = bnk->theta * gdx / (bnk->theta * gdx - (1.0 - bnk->theta) * prered + actred);
1019           tau_2 = bnk->theta * gdx / (bnk->theta * gdx + (1.0 + bnk->theta) * prered - actred);
1020           tau_min = PetscMin(tau_1, tau_2);
1021           tau_max = PetscMax(tau_1, tau_2);
1022 
1023           if (kappa >= 1.0 - bnk->mu1) {
1024             /*  Great agreement */
1025             *accept = PETSC_TRUE;
1026             if (tau_max < 1.0) {
1027               tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm);
1028             } else if (tau_max > bnk->gamma4) {
1029               tao->trust = PetscMax(tao->trust, bnk->gamma4 * bnk->dnorm);
1030             } else {
1031               tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm);
1032             }
1033           } else if (kappa >= 1.0 - bnk->mu2) {
1034             /*  Good agreement */
1035             *accept = PETSC_TRUE;
1036             if (tau_max < bnk->gamma2) {
1037               tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm);
1038             } else if (tau_max > bnk->gamma3) {
1039               tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm);
1040             } else if (tau_max < 1.0) {
1041               tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm);
1042             } else {
1043               tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm);
1044             }
1045           } else {
1046             /*  Not good agreement */
1047             if (tau_min > 1.0) {
1048               tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm);
1049             } else if (tau_max < bnk->gamma1) {
1050               tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
1051             } else if ((tau_min < bnk->gamma1) && (tau_max >= 1.0)) {
1052               tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
1053             } else if ((tau_1 >= bnk->gamma1) && (tau_1 < 1.0) && ((tau_2 < bnk->gamma1) || (tau_2 >= 1.0))) {
1054               tao->trust = tau_1 * PetscMin(tao->trust, bnk->dnorm);
1055             } else if ((tau_2 >= bnk->gamma1) && (tau_2 < 1.0) && ((tau_1 < bnk->gamma1) || (tau_2 >= 1.0))) {
1056               tao->trust = tau_2 * PetscMin(tao->trust, bnk->dnorm);
1057             } else {
1058               tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm);
1059             }
1060           }
1061         }
1062       }
1063     } else {
1064       /*  Newton step was not good; reduce the radius */
1065       tao->trust = bnk->gamma1 * PetscMin(bnk->dnorm, tao->trust);
1066     }
1067     break;
1068   }
1069   /* Make sure the radius does not violate min and max settings */
1070   tao->trust = PetscMin(tao->trust, bnk->max_radius);
1071   tao->trust = PetscMax(tao->trust, bnk->min_radius);
1072   PetscFunctionReturn(0);
1073 }
1074 
1075 /* ---------------------------------------------------------- */
1076 
1077 PetscErrorCode TaoBNKAddStepCounts(Tao tao, PetscInt stepType)
1078 {
1079   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
1080 
1081   PetscFunctionBegin;
1082   switch (stepType) {
1083   case BNK_NEWTON:
1084     ++bnk->newt;
1085     break;
1086   case BNK_BFGS:
1087     ++bnk->bfgs;
1088     break;
1089   case BNK_SCALED_GRADIENT:
1090     ++bnk->sgrad;
1091     break;
1092   case BNK_GRADIENT:
1093     ++bnk->grad;
1094     break;
1095   default:
1096     break;
1097   }
1098   PetscFunctionReturn(0);
1099 }
1100 
1101 /* ---------------------------------------------------------- */
1102 
1103 PetscErrorCode TaoSetUp_BNK(Tao tao)
1104 {
1105   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
1106   PetscErrorCode ierr;
1107   KSPType        ksp_type;
1108   PetscInt       i;
1109 
1110   PetscFunctionBegin;
1111   if (!tao->gradient) {
1112     ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);
1113   }
1114   if (!tao->stepdirection) {
1115     ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr);
1116   }
1117   if (!bnk->W) {
1118     ierr = VecDuplicate(tao->solution,&bnk->W);CHKERRQ(ierr);
1119   }
1120   if (!bnk->Xold) {
1121     ierr = VecDuplicate(tao->solution,&bnk->Xold);CHKERRQ(ierr);
1122   }
1123   if (!bnk->Gold) {
1124     ierr = VecDuplicate(tao->solution,&bnk->Gold);CHKERRQ(ierr);
1125   }
1126   if (!bnk->Xwork) {
1127     ierr = VecDuplicate(tao->solution,&bnk->Xwork);CHKERRQ(ierr);
1128   }
1129   if (!bnk->Gwork) {
1130     ierr = VecDuplicate(tao->solution,&bnk->Gwork);CHKERRQ(ierr);
1131   }
1132   if (!bnk->unprojected_gradient) {
1133     ierr = VecDuplicate(tao->solution,&bnk->unprojected_gradient);CHKERRQ(ierr);
1134   }
1135   if (!bnk->unprojected_gradient_old) {
1136     ierr = VecDuplicate(tao->solution,&bnk->unprojected_gradient_old);CHKERRQ(ierr);
1137   }
1138   if (!bnk->Diag_min) {
1139     ierr = VecDuplicate(tao->solution,&bnk->Diag_min);CHKERRQ(ierr);
1140   }
1141   if (!bnk->Diag_max) {
1142     ierr = VecDuplicate(tao->solution,&bnk->Diag_max);CHKERRQ(ierr);
1143   }
1144   if (bnk->max_cg_its > 0) {
1145     /* Ensure that the important common vectors are shared between BNK and embedded BNCG */
1146     bnk->bncg_ctx = (TAO_BNCG *)bnk->bncg->data;
1147     ierr = PetscObjectReference((PetscObject)(bnk->unprojected_gradient_old));CHKERRQ(ierr);
1148     ierr = VecDestroy(&bnk->bncg_ctx->unprojected_gradient_old);CHKERRQ(ierr);
1149     bnk->bncg_ctx->unprojected_gradient_old = bnk->unprojected_gradient_old;
1150     ierr = PetscObjectReference((PetscObject)(bnk->unprojected_gradient));CHKERRQ(ierr);
1151     ierr = VecDestroy(&bnk->bncg_ctx->unprojected_gradient);CHKERRQ(ierr);
1152     bnk->bncg_ctx->unprojected_gradient = bnk->unprojected_gradient;
1153     ierr = PetscObjectReference((PetscObject)(bnk->Gold));CHKERRQ(ierr);
1154     ierr = VecDestroy(&bnk->bncg_ctx->G_old);CHKERRQ(ierr);
1155     bnk->bncg_ctx->G_old = bnk->Gold;
1156     ierr = PetscObjectReference((PetscObject)(tao->gradient));CHKERRQ(ierr);
1157     ierr = VecDestroy(&bnk->bncg->gradient);CHKERRQ(ierr);
1158     bnk->bncg->gradient = tao->gradient;
1159     ierr = PetscObjectReference((PetscObject)(tao->stepdirection));CHKERRQ(ierr);
1160     ierr = VecDestroy(&bnk->bncg->stepdirection);CHKERRQ(ierr);
1161     bnk->bncg->stepdirection = tao->stepdirection;
1162     ierr = TaoSetInitialVector(bnk->bncg, tao->solution);CHKERRQ(ierr);
1163     /* Copy over some settings from BNK into BNCG */
1164     ierr = TaoSetMaximumIterations(bnk->bncg, bnk->max_cg_its);CHKERRQ(ierr);
1165     ierr = TaoSetTolerances(bnk->bncg, tao->gatol, tao->grtol, tao->gttol);CHKERRQ(ierr);
1166     ierr = TaoSetFunctionLowerBound(bnk->bncg, tao->fmin);CHKERRQ(ierr);
1167     ierr = TaoSetConvergenceTest(bnk->bncg, tao->ops->convergencetest, tao->cnvP);CHKERRQ(ierr);
1168     ierr = TaoSetObjectiveRoutine(bnk->bncg, tao->ops->computeobjective, tao->user_objP);CHKERRQ(ierr);
1169     ierr = TaoSetGradientRoutine(bnk->bncg, tao->ops->computegradient, tao->user_gradP);CHKERRQ(ierr);
1170     ierr = TaoSetObjectiveAndGradientRoutine(bnk->bncg, tao->ops->computeobjectiveandgradient, tao->user_objgradP);CHKERRQ(ierr);
1171     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)tao, (PetscObject)(bnk->bncg));CHKERRQ(ierr);
1172     for (i=0; i<tao->numbermonitors; ++i) {
1173       ierr = TaoSetMonitor(bnk->bncg, tao->monitor[i], tao->monitorcontext[i], tao->monitordestroy[i]);CHKERRQ(ierr);
1174       ierr = PetscObjectReference((PetscObject)(tao->monitorcontext[i]));CHKERRQ(ierr);
1175     }
1176   }
1177   bnk->Diag = 0;
1178   bnk->Diag_red = 0;
1179   bnk->X_inactive = 0;
1180   bnk->G_inactive = 0;
1181   bnk->inactive_work = 0;
1182   bnk->active_work = 0;
1183   bnk->inactive_idx = 0;
1184   bnk->active_idx = 0;
1185   bnk->active_lower = 0;
1186   bnk->active_upper = 0;
1187   bnk->active_fixed = 0;
1188   bnk->M = 0;
1189   bnk->H_inactive = 0;
1190   bnk->Hpre_inactive = 0;
1191   ierr = KSPGetType(tao->ksp,&ksp_type);CHKERRQ(ierr);
1192   ierr = PetscStrcmp(ksp_type,KSPCGNASH,&bnk->is_nash);CHKERRQ(ierr);
1193   ierr = PetscStrcmp(ksp_type,KSPCGSTCG,&bnk->is_stcg);CHKERRQ(ierr);
1194   ierr = PetscStrcmp(ksp_type,KSPCGGLTR,&bnk->is_gltr);CHKERRQ(ierr);
1195   PetscFunctionReturn(0);
1196 }
1197 
1198 /*------------------------------------------------------------*/
1199 
1200 static PetscErrorCode TaoDestroy_BNK(Tao tao)
1201 {
1202   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
1203   PetscErrorCode ierr;
1204 
1205   PetscFunctionBegin;
1206   if (tao->setupcalled) {
1207     ierr = VecDestroy(&bnk->W);CHKERRQ(ierr);
1208     ierr = VecDestroy(&bnk->Xold);CHKERRQ(ierr);
1209     ierr = VecDestroy(&bnk->Gold);CHKERRQ(ierr);
1210     ierr = VecDestroy(&bnk->Xwork);CHKERRQ(ierr);
1211     ierr = VecDestroy(&bnk->Gwork);CHKERRQ(ierr);
1212     ierr = VecDestroy(&bnk->unprojected_gradient);CHKERRQ(ierr);
1213     ierr = VecDestroy(&bnk->unprojected_gradient_old);CHKERRQ(ierr);
1214     ierr = VecDestroy(&bnk->Diag_min);CHKERRQ(ierr);
1215     ierr = VecDestroy(&bnk->Diag_max);CHKERRQ(ierr);
1216   }
1217   ierr = ISDestroy(&bnk->active_lower);CHKERRQ(ierr);
1218   ierr = ISDestroy(&bnk->active_upper);CHKERRQ(ierr);
1219   ierr = ISDestroy(&bnk->active_fixed);CHKERRQ(ierr);
1220   ierr = ISDestroy(&bnk->active_idx);CHKERRQ(ierr);
1221   ierr = ISDestroy(&bnk->inactive_idx);CHKERRQ(ierr);
1222   ierr = VecDestroy(&bnk->Diag);CHKERRQ(ierr);
1223   ierr = MatDestroy(&bnk->M);CHKERRQ(ierr);
1224   if (bnk->Hpre_inactive != tao->hessian_pre && bnk->Hpre_inactive != bnk->H_inactive) {
1225     ierr = MatDestroy(&bnk->Hpre_inactive);CHKERRQ(ierr);
1226   }
1227   if (bnk->H_inactive != tao->hessian) {
1228     ierr = MatDestroy(&bnk->H_inactive);CHKERRQ(ierr);
1229   }
1230   ierr = TaoDestroy(&bnk->bncg);CHKERRQ(ierr);
1231   ierr = PetscFree(tao->data);CHKERRQ(ierr);
1232   PetscFunctionReturn(0);
1233 }
1234 
1235 /*------------------------------------------------------------*/
1236 
1237 static PetscErrorCode TaoSetFromOptions_BNK(PetscOptionItems *PetscOptionsObject,Tao tao)
1238 {
1239   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
1240   PetscErrorCode ierr;
1241 
1242   PetscFunctionBegin;
1243   ierr = PetscOptionsHead(PetscOptionsObject,"Newton line search method for unconstrained optimization");CHKERRQ(ierr);
1244   ierr = PetscOptionsEList("-tao_bnk_pc_type", "pc type", "", BNK_PC, BNK_PC_TYPES, BNK_PC[bnk->pc_type], &bnk->pc_type, 0);CHKERRQ(ierr);
1245   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);
1246   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);
1247   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);
1248   ierr = PetscOptionsEList("-tao_bnk_as_type", "active set estimation method", "", BNK_AS, BNK_AS_TYPES, BNK_AS[bnk->as_type], &bnk->as_type, 0);CHKERRQ(ierr);
1249   ierr = PetscOptionsReal("-tao_bnk_sval", "(developer) Hessian perturbation starting value", "", bnk->sval, &bnk->sval,NULL);CHKERRQ(ierr);
1250   ierr = PetscOptionsReal("-tao_bnk_imin", "(developer) minimum initial Hessian perturbation", "", bnk->imin, &bnk->imin,NULL);CHKERRQ(ierr);
1251   ierr = PetscOptionsReal("-tao_bnk_imax", "(developer) maximum initial Hessian perturbation", "", bnk->imax, &bnk->imax,NULL);CHKERRQ(ierr);
1252   ierr = PetscOptionsReal("-tao_bnk_imfac", "(developer) initial merit factor for Hessian perturbation", "", bnk->imfac, &bnk->imfac,NULL);CHKERRQ(ierr);
1253   ierr = PetscOptionsReal("-tao_bnk_pmin", "(developer) minimum Hessian perturbation", "", bnk->pmin, &bnk->pmin,NULL);CHKERRQ(ierr);
1254   ierr = PetscOptionsReal("-tao_bnk_pmax", "(developer) maximum Hessian perturbation", "", bnk->pmax, &bnk->pmax,NULL);CHKERRQ(ierr);
1255   ierr = PetscOptionsReal("-tao_bnk_pgfac", "(developer) Hessian perturbation growth factor", "", bnk->pgfac, &bnk->pgfac,NULL);CHKERRQ(ierr);
1256   ierr = PetscOptionsReal("-tao_bnk_psfac", "(developer) Hessian perturbation shrink factor", "", bnk->psfac, &bnk->psfac,NULL);CHKERRQ(ierr);
1257   ierr = PetscOptionsReal("-tao_bnk_pmgfac", "(developer) merit growth factor for Hessian perturbation", "", bnk->pmgfac, &bnk->pmgfac,NULL);CHKERRQ(ierr);
1258   ierr = PetscOptionsReal("-tao_bnk_pmsfac", "(developer) merit shrink factor for Hessian perturbation", "", bnk->pmsfac, &bnk->pmsfac,NULL);CHKERRQ(ierr);
1259   ierr = PetscOptionsReal("-tao_bnk_eta1", "(developer) threshold for rejecting step (-tao_bnk_update_type reduction)", "", bnk->eta1, &bnk->eta1,NULL);CHKERRQ(ierr);
1260   ierr = PetscOptionsReal("-tao_bnk_eta2", "(developer) threshold for accepting marginal step (-tao_bnk_update_type reduction)", "", bnk->eta2, &bnk->eta2,NULL);CHKERRQ(ierr);
1261   ierr = PetscOptionsReal("-tao_bnk_eta3", "(developer) threshold for accepting reasonable step (-tao_bnk_update_type reduction)", "", bnk->eta3, &bnk->eta3,NULL);CHKERRQ(ierr);
1262   ierr = PetscOptionsReal("-tao_bnk_eta4", "(developer) threshold for accepting good step (-tao_bnk_update_type reduction)", "", bnk->eta4, &bnk->eta4,NULL);CHKERRQ(ierr);
1263   ierr = PetscOptionsReal("-tao_bnk_alpha1", "(developer) radius reduction factor for rejected step (-tao_bnk_update_type reduction)", "", bnk->alpha1, &bnk->alpha1,NULL);CHKERRQ(ierr);
1264   ierr = PetscOptionsReal("-tao_bnk_alpha2", "(developer) radius reduction factor for marginally accepted bad step (-tao_bnk_update_type reduction)", "", bnk->alpha2, &bnk->alpha2,NULL);CHKERRQ(ierr);
1265   ierr = PetscOptionsReal("-tao_bnk_alpha3", "(developer) radius increase factor for reasonable accepted step (-tao_bnk_update_type reduction)", "", bnk->alpha3, &bnk->alpha3,NULL);CHKERRQ(ierr);
1266   ierr = PetscOptionsReal("-tao_bnk_alpha4", "(developer) radius increase factor for good accepted step (-tao_bnk_update_type reduction)", "", bnk->alpha4, &bnk->alpha4,NULL);CHKERRQ(ierr);
1267   ierr = PetscOptionsReal("-tao_bnk_alpha5", "(developer) radius increase factor for very good accepted step (-tao_bnk_update_type reduction)", "", bnk->alpha5, &bnk->alpha5,NULL);CHKERRQ(ierr);
1268   ierr = PetscOptionsReal("-tao_bnk_nu1", "(developer) threshold for small line-search step length (-tao_bnk_update_type step)", "", bnk->nu1, &bnk->nu1,NULL);CHKERRQ(ierr);
1269   ierr = PetscOptionsReal("-tao_bnk_nu2", "(developer) threshold for reasonable line-search step length (-tao_bnk_update_type step)", "", bnk->nu2, &bnk->nu2,NULL);CHKERRQ(ierr);
1270   ierr = PetscOptionsReal("-tao_bnk_nu3", "(developer) threshold for large line-search step length (-tao_bnk_update_type step)", "", bnk->nu3, &bnk->nu3,NULL);CHKERRQ(ierr);
1271   ierr = PetscOptionsReal("-tao_bnk_nu4", "(developer) threshold for very large line-search step length (-tao_bnk_update_type step)", "", bnk->nu4, &bnk->nu4,NULL);CHKERRQ(ierr);
1272   ierr = PetscOptionsReal("-tao_bnk_omega1", "(developer) radius reduction factor for very small line-search step length (-tao_bnk_update_type step)", "", bnk->omega1, &bnk->omega1,NULL);CHKERRQ(ierr);
1273   ierr = PetscOptionsReal("-tao_bnk_omega2", "(developer) radius reduction factor for small line-search step length (-tao_bnk_update_type step)", "", bnk->omega2, &bnk->omega2,NULL);CHKERRQ(ierr);
1274   ierr = PetscOptionsReal("-tao_bnk_omega3", "(developer) radius factor for decent line-search step length (-tao_bnk_update_type step)", "", bnk->omega3, &bnk->omega3,NULL);CHKERRQ(ierr);
1275   ierr = PetscOptionsReal("-tao_bnk_omega4", "(developer) radius increase factor for large line-search step length (-tao_bnk_update_type step)", "", bnk->omega4, &bnk->omega4,NULL);CHKERRQ(ierr);
1276   ierr = PetscOptionsReal("-tao_bnk_omega5", "(developer) radius increase factor for very large line-search step length (-tao_bnk_update_type step)", "", bnk->omega5, &bnk->omega5,NULL);CHKERRQ(ierr);
1277   ierr = PetscOptionsReal("-tao_bnk_mu1_i", "(developer) threshold for accepting very good step (-tao_bnk_init_type interpolation)", "", bnk->mu1_i, &bnk->mu1_i,NULL);CHKERRQ(ierr);
1278   ierr = PetscOptionsReal("-tao_bnk_mu2_i", "(developer) threshold for accepting good step (-tao_bnk_init_type interpolation)", "", bnk->mu2_i, &bnk->mu2_i,NULL);CHKERRQ(ierr);
1279   ierr = PetscOptionsReal("-tao_bnk_gamma1_i", "(developer) radius reduction factor for rejected very bad step (-tao_bnk_init_type interpolation)", "", bnk->gamma1_i, &bnk->gamma1_i,NULL);CHKERRQ(ierr);
1280   ierr = PetscOptionsReal("-tao_bnk_gamma2_i", "(developer) radius reduction factor for rejected bad step (-tao_bnk_init_type interpolation)", "", bnk->gamma2_i, &bnk->gamma2_i,NULL);CHKERRQ(ierr);
1281   ierr = PetscOptionsReal("-tao_bnk_gamma3_i", "(developer) radius increase factor for accepted good step (-tao_bnk_init_type interpolation)", "", bnk->gamma3_i, &bnk->gamma3_i,NULL);CHKERRQ(ierr);
1282   ierr = PetscOptionsReal("-tao_bnk_gamma4_i", "(developer) radius increase factor for accepted very good step (-tao_bnk_init_type interpolation)", "", bnk->gamma4_i, &bnk->gamma4_i,NULL);CHKERRQ(ierr);
1283   ierr = PetscOptionsReal("-tao_bnk_theta_i", "(developer) trust region interpolation factor (-tao_bnk_init_type interpolation)", "", bnk->theta_i, &bnk->theta_i,NULL);CHKERRQ(ierr);
1284   ierr = PetscOptionsReal("-tao_bnk_mu1", "(developer) threshold for accepting very good step (-tao_bnk_update_type interpolation)", "", bnk->mu1, &bnk->mu1,NULL);CHKERRQ(ierr);
1285   ierr = PetscOptionsReal("-tao_bnk_mu2", "(developer) threshold for accepting good step (-tao_bnk_update_type interpolation)", "", bnk->mu2, &bnk->mu2,NULL);CHKERRQ(ierr);
1286   ierr = PetscOptionsReal("-tao_bnk_gamma1", "(developer) radius reduction factor for rejected very bad step (-tao_bnk_update_type interpolation)", "", bnk->gamma1, &bnk->gamma1,NULL);CHKERRQ(ierr);
1287   ierr = PetscOptionsReal("-tao_bnk_gamma2", "(developer) radius reduction factor for rejected bad step (-tao_bnk_update_type interpolation)", "", bnk->gamma2, &bnk->gamma2,NULL);CHKERRQ(ierr);
1288   ierr = PetscOptionsReal("-tao_bnk_gamma3", "(developer) radius increase factor for accepted good step (-tao_bnk_update_type interpolation)", "", bnk->gamma3, &bnk->gamma3,NULL);CHKERRQ(ierr);
1289   ierr = PetscOptionsReal("-tao_bnk_gamma4", "(developer) radius increase factor for accepted very good step (-tao_bnk_update_type interpolation)", "", bnk->gamma4, &bnk->gamma4,NULL);CHKERRQ(ierr);
1290   ierr = PetscOptionsReal("-tao_bnk_theta", "(developer) trust region interpolation factor (-tao_bnk_update_type interpolation)", "", bnk->theta, &bnk->theta,NULL);CHKERRQ(ierr);
1291   ierr = PetscOptionsReal("-tao_bnk_min_radius", "(developer) lower bound on initial radius", "", bnk->min_radius, &bnk->min_radius,NULL);CHKERRQ(ierr);
1292   ierr = PetscOptionsReal("-tao_bnk_max_radius", "(developer) upper bound on radius", "", bnk->max_radius, &bnk->max_radius,NULL);CHKERRQ(ierr);
1293   ierr = PetscOptionsReal("-tao_bnk_epsilon", "(developer) tolerance used when computing actual and predicted reduction", "", bnk->epsilon, &bnk->epsilon,NULL);CHKERRQ(ierr);
1294   ierr = PetscOptionsReal("-tao_bnk_as_tol", "(developer) initial tolerance used when estimating actively bounded variables", "", bnk->as_tol, &bnk->as_tol,NULL);CHKERRQ(ierr);
1295   ierr = PetscOptionsReal("-tao_bnk_as_step", "(developer) step length used when estimating actively bounded variables", "", bnk->as_step, &bnk->as_step,NULL);CHKERRQ(ierr);
1296   ierr = PetscOptionsInt("-tao_bnk_max_cg_its", "number of BNCG iterations to take for each Newton step", "", bnk->max_cg_its, &bnk->max_cg_its,NULL);CHKERRQ(ierr);
1297   ierr = PetscOptionsTail();CHKERRQ(ierr);
1298   ierr = TaoSetFromOptions(bnk->bncg);CHKERRQ(ierr);
1299   ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr);
1300   ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
1301   PetscFunctionReturn(0);
1302 }
1303 
1304 /*------------------------------------------------------------*/
1305 
1306 static PetscErrorCode TaoView_BNK(Tao tao, PetscViewer viewer)
1307 {
1308   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
1309   PetscInt       nrejects;
1310   PetscBool      isascii;
1311   PetscErrorCode ierr;
1312 
1313   PetscFunctionBegin;
1314   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
1315   if (isascii) {
1316     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1317     if (BNK_PC_BFGS == bnk->pc_type && bnk->M) {
1318       ierr = MatLMVMGetRejects(bnk->M,&nrejects);CHKERRQ(ierr);
1319       ierr = PetscViewerASCIIPrintf(viewer, "Rejected matrix updates: %D\n",nrejects);CHKERRQ(ierr);
1320     }
1321     ierr = PetscViewerASCIIPrintf(viewer, "CG steps: %D\n", bnk->tot_cg_its);CHKERRQ(ierr);
1322     ierr = PetscViewerASCIIPrintf(viewer, "Newton steps: %D\n", bnk->newt);CHKERRQ(ierr);
1323     ierr = PetscViewerASCIIPrintf(viewer, "BFGS steps: %D\n", bnk->bfgs);CHKERRQ(ierr);
1324     ierr = PetscViewerASCIIPrintf(viewer, "Scaled gradient steps: %D\n", bnk->sgrad);CHKERRQ(ierr);
1325     ierr = PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", bnk->grad);CHKERRQ(ierr);
1326     ierr = PetscViewerASCIIPrintf(viewer, "KSP termination reasons:\n");CHKERRQ(ierr);
1327     ierr = PetscViewerASCIIPrintf(viewer, "  atol: %D\n", bnk->ksp_atol);CHKERRQ(ierr);
1328     ierr = PetscViewerASCIIPrintf(viewer, "  rtol: %D\n", bnk->ksp_rtol);CHKERRQ(ierr);
1329     ierr = PetscViewerASCIIPrintf(viewer, "  ctol: %D\n", bnk->ksp_ctol);CHKERRQ(ierr);
1330     ierr = PetscViewerASCIIPrintf(viewer, "  negc: %D\n", bnk->ksp_negc);CHKERRQ(ierr);
1331     ierr = PetscViewerASCIIPrintf(viewer, "  dtol: %D\n", bnk->ksp_dtol);CHKERRQ(ierr);
1332     ierr = PetscViewerASCIIPrintf(viewer, "  iter: %D\n", bnk->ksp_iter);CHKERRQ(ierr);
1333     ierr = PetscViewerASCIIPrintf(viewer, "  othr: %D\n", bnk->ksp_othr);CHKERRQ(ierr);
1334     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1335   }
1336   PetscFunctionReturn(0);
1337 }
1338 
1339 /* ---------------------------------------------------------- */
1340 
1341 /*MC
1342   TAOBNK - Shared base-type for Bounded Newton-Krylov type algorithms.
1343   At each iteration, the BNK methods solve the symmetric
1344   system of equations to obtain the step diretion dk:
1345               Hk dk = -gk
1346   for free variables only. The step can be globalized either through
1347   trust-region methods, or a line search, or a heuristic mixture of both.
1348 
1349     Options Database Keys:
1350 + -tao_max_cg_its - maximum number of bounded conjugate-gradient iterations taken in each Newton loop
1351 . -tao_bnk_pc_type - preconditioner type ("none", "ahess", "bfgs", "petsc")
1352 . -tao_bnk_bfgs_scale_type - BFGS preconditioner diagonal scaling type ("ahess", "phess", "bfgs")
1353 . -tao_bnk_init_type - trust radius initialization method ("constant", "direction", "interpolation")
1354 . -tao_bnk_update_type - trust radius update method ("step", "direction", "interpolation")
1355 . -tao_bnk_as_type - active-set estimation method ("none", "bertsekas")
1356 . -tao_bnk_as_tol - (developer) initial tolerance used in estimating bounded active variables (-tao_bnk_as_type bertsekas)
1357 . -tao_bnk_as_step - (developer) trial step length used in estimating bounded active variables (-tao_bnk_as_type bertsekas)
1358 . -tao_bnk_sval - (developer) Hessian perturbation starting value
1359 . -tao_bnk_imin - (developer) minimum initial Hessian perturbation
1360 . -tao_bnk_imax - (developer) maximum initial Hessian perturbation
1361 . -tao_bnk_pmin - (developer) minimum Hessian perturbation
1362 . -tao_bnk_pmax - (developer) aximum Hessian perturbation
1363 . -tao_bnk_pgfac - (developer) Hessian perturbation growth factor
1364 . -tao_bnk_psfac - (developer) Hessian perturbation shrink factor
1365 . -tao_bnk_imfac - (developer) initial merit factor for Hessian perturbation
1366 . -tao_bnk_pmgfac - (developer) merit growth factor for Hessian perturbation
1367 . -tao_bnk_pmsfac - (developer) merit shrink factor for Hessian perturbation
1368 . -tao_bnk_eta1 - (developer) threshold for rejecting step (-tao_bnk_update_type reduction)
1369 . -tao_bnk_eta2 - (developer) threshold for accepting marginal step (-tao_bnk_update_type reduction)
1370 . -tao_bnk_eta3 - (developer) threshold for accepting reasonable step (-tao_bnk_update_type reduction)
1371 . -tao_bnk_eta4 - (developer) threshold for accepting good step (-tao_bnk_update_type reduction)
1372 . -tao_bnk_alpha1 - (developer) radius reduction factor for rejected step (-tao_bnk_update_type reduction)
1373 . -tao_bnk_alpha2 - (developer) radius reduction factor for marginally accepted bad step (-tao_bnk_update_type reduction)
1374 . -tao_bnk_alpha3 - (developer) radius increase factor for reasonable accepted step (-tao_bnk_update_type reduction)
1375 . -tao_bnk_alpha4 - (developer) radius increase factor for good accepted step (-tao_bnk_update_type reduction)
1376 . -tao_bnk_alpha5 - (developer) radius increase factor for very good accepted step (-tao_bnk_update_type reduction)
1377 . -tao_bnk_epsilon - (developer) tolerance for small pred/actual ratios that trigger automatic step acceptance (-tao_bnk_update_type reduction)
1378 . -tao_bnk_mu1 - (developer) threshold for accepting very good step (-tao_bnk_update_type interpolation)
1379 . -tao_bnk_mu2 - (developer) threshold for accepting good step (-tao_bnk_update_type interpolation)
1380 . -tao_bnk_gamma1 - (developer) radius reduction factor for rejected very bad step (-tao_bnk_update_type interpolation)
1381 . -tao_bnk_gamma2 - (developer) radius reduction factor for rejected bad step (-tao_bnk_update_type interpolation)
1382 . -tao_bnk_gamma3 - (developer) radius increase factor for accepted good step (-tao_bnk_update_type interpolation)
1383 . -tao_bnk_gamma4 - (developer) radius increase factor for accepted very good step (-tao_bnk_update_type interpolation)
1384 . -tao_bnk_theta - (developer) trust region interpolation factor (-tao_bnk_update_type interpolation)
1385 . -tao_bnk_nu1 - (developer) threshold for small line-search step length (-tao_bnk_update_type step)
1386 . -tao_bnk_nu2 - (developer) threshold for reasonable line-search step length (-tao_bnk_update_type step)
1387 . -tao_bnk_nu3 - (developer) threshold for large line-search step length (-tao_bnk_update_type step)
1388 . -tao_bnk_nu4 - (developer) threshold for very large line-search step length (-tao_bnk_update_type step)
1389 . -tao_bnk_omega1 - (developer) radius reduction factor for very small line-search step length (-tao_bnk_update_type step)
1390 . -tao_bnk_omega2 - (developer) radius reduction factor for small line-search step length (-tao_bnk_update_type step)
1391 . -tao_bnk_omega3 - (developer) radius factor for decent line-search step length (-tao_bnk_update_type step)
1392 . -tao_bnk_omega4 - (developer) radius increase factor for large line-search step length (-tao_bnk_update_type step)
1393 . -tao_bnk_omega5 - (developer) radius increase factor for very large line-search step length (-tao_bnk_update_type step)
1394 . -tao_bnk_mu1_i -  (developer) threshold for accepting very good step (-tao_bnk_init_type interpolation)
1395 . -tao_bnk_mu2_i -  (developer) threshold for accepting good step (-tao_bnk_init_type interpolation)
1396 . -tao_bnk_gamma1_i - (developer) radius reduction factor for rejected very bad step (-tao_bnk_init_type interpolation)
1397 . -tao_bnk_gamma2_i - (developer) radius reduction factor for rejected bad step (-tao_bnk_init_type interpolation)
1398 . -tao_bnk_gamma3_i - (developer) radius increase factor for accepted good step (-tao_bnk_init_type interpolation)
1399 . -tao_bnk_gamma4_i - (developer) radius increase factor for accepted very good step (-tao_bnk_init_type interpolation)
1400 - -tao_bnk_theta_i - (developer) trust region interpolation factor (-tao_bnk_init_type interpolation)
1401 
1402   Level: beginner
1403 M*/
1404 
1405 PetscErrorCode TaoCreate_BNK(Tao tao)
1406 {
1407   TAO_BNK        *bnk;
1408   const char     *morethuente_type = TAOLINESEARCHMT;
1409   PetscErrorCode ierr;
1410 
1411   PetscFunctionBegin;
1412   ierr = PetscNewLog(tao,&bnk);CHKERRQ(ierr);
1413 
1414   tao->ops->setup = TaoSetUp_BNK;
1415   tao->ops->view = TaoView_BNK;
1416   tao->ops->setfromoptions = TaoSetFromOptions_BNK;
1417   tao->ops->destroy = TaoDestroy_BNK;
1418 
1419   /*  Override default settings (unless already changed) */
1420   if (!tao->max_it_changed) tao->max_it = 50;
1421   if (!tao->trust0_changed) tao->trust0 = 100.0;
1422 
1423   tao->data = (void*)bnk;
1424 
1425   /*  Hessian shifting parameters */
1426   bnk->sval   = 0.0;
1427   bnk->imin   = 1.0e-4;
1428   bnk->imax   = 1.0e+2;
1429   bnk->imfac  = 1.0e-1;
1430 
1431   bnk->pmin   = 1.0e-12;
1432   bnk->pmax   = 1.0e+2;
1433   bnk->pgfac  = 1.0e+1;
1434   bnk->psfac  = 4.0e-1;
1435   bnk->pmgfac = 1.0e-1;
1436   bnk->pmsfac = 1.0e-1;
1437 
1438   /*  Default values for trust-region radius update based on steplength */
1439   bnk->nu1 = 0.25;
1440   bnk->nu2 = 0.50;
1441   bnk->nu3 = 1.00;
1442   bnk->nu4 = 1.25;
1443 
1444   bnk->omega1 = 0.25;
1445   bnk->omega2 = 0.50;
1446   bnk->omega3 = 1.00;
1447   bnk->omega4 = 2.00;
1448   bnk->omega5 = 4.00;
1449 
1450   /*  Default values for trust-region radius update based on reduction */
1451   bnk->eta1 = 1.0e-4;
1452   bnk->eta2 = 0.25;
1453   bnk->eta3 = 0.50;
1454   bnk->eta4 = 0.90;
1455 
1456   bnk->alpha1 = 0.25;
1457   bnk->alpha2 = 0.50;
1458   bnk->alpha3 = 1.00;
1459   bnk->alpha4 = 2.00;
1460   bnk->alpha5 = 4.00;
1461 
1462   /*  Default values for trust-region radius update based on interpolation */
1463   bnk->mu1 = 0.10;
1464   bnk->mu2 = 0.50;
1465 
1466   bnk->gamma1 = 0.25;
1467   bnk->gamma2 = 0.50;
1468   bnk->gamma3 = 2.00;
1469   bnk->gamma4 = 4.00;
1470 
1471   bnk->theta = 0.05;
1472 
1473   /*  Default values for trust region initialization based on interpolation */
1474   bnk->mu1_i = 0.35;
1475   bnk->mu2_i = 0.50;
1476 
1477   bnk->gamma1_i = 0.0625;
1478   bnk->gamma2_i = 0.5;
1479   bnk->gamma3_i = 2.0;
1480   bnk->gamma4_i = 5.0;
1481 
1482   bnk->theta_i = 0.25;
1483 
1484   /*  Remaining parameters */
1485   bnk->max_cg_its = 0;
1486   bnk->min_radius = 1.0e-10;
1487   bnk->max_radius = 1.0e10;
1488   bnk->epsilon = PetscPowReal(PETSC_MACHINE_EPSILON, 2.0/3.0);
1489   bnk->as_tol = 1.0e-3;
1490   bnk->as_step = 1.0e-3;
1491   bnk->dmin = 1.0e-6;
1492   bnk->dmax = 1.0e6;
1493 
1494   bnk->pc_type         = BNK_PC_AHESS;
1495   bnk->bfgs_scale_type = BFGS_SCALE_PHESS;
1496   bnk->init_type       = BNK_INIT_INTERPOLATION;
1497   bnk->update_type     = BNK_UPDATE_INTERPOLATION;
1498   bnk->as_type         = BNK_AS_BERTSEKAS;
1499 
1500   /* Create the embedded BNCG solver */
1501   ierr = TaoCreate(PetscObjectComm((PetscObject)tao), &bnk->bncg);CHKERRQ(ierr);
1502   ierr = PetscObjectIncrementTabLevel((PetscObject)bnk->bncg, (PetscObject)tao, 1);CHKERRQ(ierr);
1503   ierr = TaoSetOptionsPrefix(bnk->bncg, "tao_bnk_");CHKERRQ(ierr);
1504   ierr = TaoSetType(bnk->bncg, TAOBNCG);CHKERRQ(ierr);
1505 
1506   /* Create the line search */
1507   ierr = TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch);CHKERRQ(ierr);
1508   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr);
1509   ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr);
1510   ierr = TaoLineSearchSetType(tao->linesearch,morethuente_type);CHKERRQ(ierr);
1511   ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr);
1512 
1513   /*  Set linear solver to default for symmetric matrices */
1514   ierr = KSPCreate(((PetscObject)tao)->comm,&tao->ksp);CHKERRQ(ierr);
1515   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1);CHKERRQ(ierr);
1516   ierr = KSPSetOptionsPrefix(tao->ksp,tao->hdr.prefix);CHKERRQ(ierr);
1517   ierr = KSPSetType(tao->ksp,KSPCGSTCG);CHKERRQ(ierr);
1518   PetscFunctionReturn(0);
1519 }
1520