xref: /petsc/src/tao/bound/impls/bnk/bnk.c (revision ca964c4924045cf1d9552da2e7805dac7209e511)
1eb910715SAlp Dener #include <petsctaolinesearch.h>
2eb910715SAlp Dener #include <../src/tao/bound/impls/bnk/bnk.h>
3eb910715SAlp Dener 
4eb910715SAlp Dener #include <petscksp.h>
5eb910715SAlp Dener 
6e031d6f5SAlp Dener static const char *BNK_PC[64] = {"none", "ahess", "bfgs", "petsc"};
7e031d6f5SAlp Dener static const char *BFGS_SCALE[64] = {"ahess", "phess", "bfgs"};
8e031d6f5SAlp Dener static const char *BNK_INIT[64] = {"constant", "direction", "interpolation"};
9e031d6f5SAlp Dener static const char *BNK_UPDATE[64] = {"step", "reduction", "interpolation"};
10e031d6f5SAlp Dener static const char *BNK_AS[64] = {"none", "bertsekas"};
11e031d6f5SAlp Dener 
12e031d6f5SAlp Dener /*------------------------------------------------------------*/
13e031d6f5SAlp Dener 
14eb910715SAlp Dener /* Routine for BFGS preconditioner */
15eb910715SAlp Dener 
16eb910715SAlp Dener PetscErrorCode MatLMVMSolveShell(PC pc, Vec b, Vec x)
17eb910715SAlp Dener {
18eb910715SAlp Dener   PetscErrorCode ierr;
19eb910715SAlp Dener   Mat            M;
20eb910715SAlp Dener 
21eb910715SAlp Dener   PetscFunctionBegin;
22eb910715SAlp Dener   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
23eb910715SAlp Dener   PetscValidHeaderSpecific(b,VEC_CLASSID,2);
24eb910715SAlp Dener   PetscValidHeaderSpecific(x,VEC_CLASSID,3);
25eb910715SAlp Dener   ierr = PCShellGetContext(pc,(void**)&M);CHKERRQ(ierr);
265e9b73cbSAlp Dener   ierr = MatLMVMSolve(M, b, x);CHKERRQ(ierr);
27eb910715SAlp Dener   PetscFunctionReturn(0);
28eb910715SAlp Dener }
29eb910715SAlp Dener 
30df278d8fSAlp Dener /*------------------------------------------------------------*/
31df278d8fSAlp Dener 
32df278d8fSAlp Dener /* Routine for initializing the KSP solver, the BFGS preconditioner, and the initial trust radius estimation */
33df278d8fSAlp Dener 
34c0f10754SAlp Dener PetscErrorCode TaoBNKInitialize(Tao tao, PetscInt initType, PetscBool *needH)
35eb910715SAlp Dener {
36eb910715SAlp Dener   PetscErrorCode               ierr;
37eb910715SAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
38eb910715SAlp Dener   PC                           pc;
39eb910715SAlp Dener 
4089da521bSAlp Dener   PetscReal                    f_min, ftrial, prered, actred, kappa, sigma, resnorm;
41eb910715SAlp Dener   PetscReal                    tau, tau_1, tau_2, tau_max, tau_min, max_radius;
4289da521bSAlp Dener   PetscReal                    delta;
43eb910715SAlp Dener 
44c4b75bccSAlp Dener   PetscInt                     n,N,nDiff;
45eb910715SAlp Dener 
46eb910715SAlp Dener   PetscInt                     i_max = 5;
47eb910715SAlp Dener   PetscInt                     j_max = 1;
48eb910715SAlp Dener   PetscInt                     i, j;
49eb910715SAlp Dener 
50eb910715SAlp Dener   PetscFunctionBegin;
5128017e9fSAlp Dener   /* Project the current point onto the feasible set */
5228017e9fSAlp Dener   ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr);
53e031d6f5SAlp Dener   ierr = TaoSetVariableBounds(bnk->bncg, tao->XL, tao->XU);CHKERRQ(ierr);
5428017e9fSAlp Dener   ierr = TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU);CHKERRQ(ierr);
5528017e9fSAlp Dener 
5628017e9fSAlp Dener   /* Project the initial point onto the feasible region */
5789da521bSAlp Dener   ierr = TaoBoundSolution(tao->XL, tao->XU, tao->solution, 0.0, &nDiff);CHKERRQ(ierr);
5828017e9fSAlp Dener 
5928017e9fSAlp Dener   /* Check convergence criteria */
6028017e9fSAlp Dener   ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &bnk->f, bnk->unprojected_gradient);CHKERRQ(ierr);
6161be54a6SAlp Dener   ierr = TaoBNKEstimateActiveSet(tao, bnk->as_type);CHKERRQ(ierr);
6261be54a6SAlp Dener   ierr = VecCopy(bnk->unprojected_gradient, tao->gradient);CHKERRQ(ierr);
6361be54a6SAlp Dener   ierr = VecISSet(tao->gradient, bnk->active_idx, 0.0);CHKERRQ(ierr);
6408752603SAlp Dener   ierr = VecNorm(tao->gradient,NORM_2,&bnk->gnorm);CHKERRQ(ierr);
6528017e9fSAlp Dener 
66c0f10754SAlp Dener   /* Test the initial point for convergence */
6789da521bSAlp Dener   ierr = VecFischer(tao->solution, bnk->unprojected_gradient, tao->XL, tao->XU, bnk->W);CHKERRQ(ierr);
6889da521bSAlp Dener   ierr = VecNorm(bnk->W, NORM_2, &resnorm);CHKERRQ(ierr);
69b4a30f08SAlp Dener   if (PetscIsInfOrNanReal(bnk->f) || PetscIsInfOrNanReal(resnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
70e031d6f5SAlp Dener   ierr = TaoLogConvergenceHistory(tao,bnk->f,resnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
71e031d6f5SAlp Dener   ierr = TaoMonitor(tao,tao->niter,bnk->f,resnorm,0.0,1.0);CHKERRQ(ierr);
72c0f10754SAlp Dener   ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
73c0f10754SAlp Dener   if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
74c0f10754SAlp Dener 
75e031d6f5SAlp Dener   /* Reset KSP stopping reason counters */
76eb910715SAlp Dener   bnk->ksp_atol = 0;
77eb910715SAlp Dener   bnk->ksp_rtol = 0;
78eb910715SAlp Dener   bnk->ksp_dtol = 0;
79eb910715SAlp Dener   bnk->ksp_ctol = 0;
80eb910715SAlp Dener   bnk->ksp_negc = 0;
81eb910715SAlp Dener   bnk->ksp_iter = 0;
82eb910715SAlp Dener   bnk->ksp_othr = 0;
83eb910715SAlp Dener 
84e031d6f5SAlp Dener   /* Reset accepted step type counters */
85e031d6f5SAlp Dener   bnk->tot_cg_its = 0;
86e031d6f5SAlp Dener   bnk->newt = 0;
87e031d6f5SAlp Dener   bnk->bfgs = 0;
88e031d6f5SAlp Dener   bnk->sgrad = 0;
89e031d6f5SAlp Dener   bnk->grad = 0;
90e031d6f5SAlp Dener 
91fed79b8eSAlp Dener   /* Initialize the Hessian perturbation */
92fed79b8eSAlp Dener   bnk->pert = bnk->sval;
93fed79b8eSAlp Dener 
94937a31a1SAlp Dener   /* Reset initial steplength to zero (this helps BNCG reset its direction internally) */
95937a31a1SAlp Dener   ierr = VecSet(tao->stepdirection, 0.0);CHKERRQ(ierr);
96937a31a1SAlp Dener 
97e031d6f5SAlp Dener   /* Allocate the vectors needed for the BFGS approximation */
98e031d6f5SAlp Dener   if (BNK_PC_BFGS == bnk->pc_type) {
99e031d6f5SAlp Dener     if (!bnk->M) {
100eb910715SAlp Dener       ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr);
101eb910715SAlp Dener       ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr);
102eb910715SAlp Dener       ierr = MatCreateLMVM(((PetscObject)tao)->comm,n,N,&bnk->M);CHKERRQ(ierr);
103eb910715SAlp Dener       ierr = MatLMVMAllocateVectors(bnk->M,tao->solution);CHKERRQ(ierr);
104eb910715SAlp Dener     }
105e031d6f5SAlp Dener     if (bnk->bfgs_scale_type != BFGS_SCALE_BFGS && !bnk->Diag) {
106e031d6f5SAlp Dener       ierr = VecDuplicate(tao->solution,&bnk->Diag);CHKERRQ(ierr);
1075e9b73cbSAlp Dener     }
108e031d6f5SAlp Dener   }
109e031d6f5SAlp Dener 
110e031d6f5SAlp Dener   /* Prepare the min/max vectors for safeguarding diagonal scales */
11162675beeSAlp Dener   ierr = VecSet(bnk->Diag_min, bnk->dmin);CHKERRQ(ierr);
11262675beeSAlp Dener   ierr = VecSet(bnk->Diag_max, bnk->dmax);CHKERRQ(ierr);
113eb910715SAlp Dener 
114eb910715SAlp Dener   /* Modify the preconditioner to use the bfgs approximation */
115eb910715SAlp Dener   ierr = KSPGetPC(tao->ksp, &pc);CHKERRQ(ierr);
116eb910715SAlp Dener   switch(bnk->pc_type) {
117eb910715SAlp Dener   case BNK_PC_NONE:
118eb910715SAlp Dener     ierr = PCSetType(pc, PCNONE);CHKERRQ(ierr);
119eb910715SAlp Dener     ierr = PCSetFromOptions(pc);CHKERRQ(ierr);
120eb910715SAlp Dener     break;
121eb910715SAlp Dener 
122eb910715SAlp Dener   case BNK_PC_AHESS:
123eb910715SAlp Dener     ierr = PCSetType(pc, PCJACOBI);CHKERRQ(ierr);
124eb910715SAlp Dener     ierr = PCSetFromOptions(pc);CHKERRQ(ierr);
125eb910715SAlp Dener     ierr = PCJacobiSetUseAbs(pc,PETSC_TRUE);CHKERRQ(ierr);
126eb910715SAlp Dener     break;
127eb910715SAlp Dener 
128eb910715SAlp Dener   case BNK_PC_BFGS:
129eb910715SAlp Dener     ierr = PCSetType(pc, PCSHELL);CHKERRQ(ierr);
130eb910715SAlp Dener     ierr = PCSetFromOptions(pc);CHKERRQ(ierr);
131eb910715SAlp Dener     ierr = PCShellSetName(pc, "bfgs");CHKERRQ(ierr);
132eb910715SAlp Dener     ierr = PCShellSetContext(pc, bnk->M);CHKERRQ(ierr);
133eb910715SAlp Dener     ierr = PCShellSetApply(pc, MatLMVMSolveShell);CHKERRQ(ierr);
134eb910715SAlp Dener     break;
135eb910715SAlp Dener 
136eb910715SAlp Dener   default:
137eb910715SAlp Dener     /* Use the pc method set by pc_type */
138eb910715SAlp Dener     break;
139eb910715SAlp Dener   }
140eb910715SAlp Dener 
141eb910715SAlp Dener   /* Initialize trust-region radius.  The initialization is only performed
142eb910715SAlp Dener      when we are using Nash, Steihaug-Toint or the Generalized Lanczos method. */
143c0f10754SAlp Dener   *needH = PETSC_TRUE;
144eb910715SAlp Dener   if (bnk->is_nash || bnk->is_stcg || bnk->is_gltr) {
14562675beeSAlp Dener     switch(initType) {
146eb910715SAlp Dener     case BNK_INIT_CONSTANT:
147eb910715SAlp Dener       /* Use the initial radius specified */
148c0f10754SAlp Dener       tao->trust = tao->trust0;
149eb910715SAlp Dener       break;
150eb910715SAlp Dener 
151eb910715SAlp Dener     case BNK_INIT_INTERPOLATION:
152c0f10754SAlp Dener       /* Use interpolation based on the initial Hessian */
153eb910715SAlp Dener       max_radius = 0.0;
15408752603SAlp Dener       tao->trust = tao->trust0;
155eb910715SAlp Dener       for (j = 0; j < j_max; ++j) {
1560a4511e9SAlp Dener         f_min = bnk->f;
157eb910715SAlp Dener         sigma = 0.0;
158eb910715SAlp Dener 
159c0f10754SAlp Dener         if (*needH) {
16062602cfbSAlp Dener           /* Compute the Hessian at the new step, and extract the inactive subsystem */
161e031d6f5SAlp Dener           ierr = TaoBNKComputeHessian(tao);CHKERRQ(ierr);
16208752603SAlp Dener           ierr = TaoBNKEstimateActiveSet(tao, BNK_AS_NONE);CHKERRQ(ierr);
16389da521bSAlp Dener           ierr = MatDestroy(&bnk->H_inactive);CHKERRQ(ierr);
16489da521bSAlp Dener           if (bnk->active_idx) {
1652ab2a32cSAlp Dener             ierr = MatCreateSubMatrix(tao->hessian, bnk->inactive_idx, bnk->inactive_idx, MAT_INITIAL_MATRIX, &bnk->H_inactive);CHKERRQ(ierr);
16628017e9fSAlp Dener           } else {
16708752603SAlp Dener             ierr = MatDuplicate(tao->hessian, MAT_COPY_VALUES, &bnk->H_inactive);CHKERRQ(ierr);
16828017e9fSAlp Dener           }
169c0f10754SAlp Dener           *needH = PETSC_FALSE;
170eb910715SAlp Dener         }
171eb910715SAlp Dener 
172eb910715SAlp Dener         for (i = 0; i < i_max; ++i) {
17362602cfbSAlp Dener           /* Take a steepest descent step and snap it to bounds */
17462602cfbSAlp Dener           ierr = VecCopy(tao->solution, bnk->Xold);CHKERRQ(ierr);
17562602cfbSAlp Dener           ierr = VecAXPY(tao->solution, -tao->trust/bnk->gnorm, tao->gradient);CHKERRQ(ierr);
17689da521bSAlp Dener           ierr = TaoBoundSolution(tao->XL, tao->XU, tao->solution, 0.0, &nDiff);CHKERRQ(ierr);
17789da521bSAlp Dener           /* Compute the step we actually accepted */
178eb910715SAlp Dener           ierr = VecCopy(tao->solution, bnk->W);CHKERRQ(ierr);
17962602cfbSAlp Dener           ierr = VecAXPY(bnk->W, -1.0, bnk->Xold);CHKERRQ(ierr);
18062602cfbSAlp Dener           /* Compute the objective at the trial */
18162602cfbSAlp Dener           ierr = TaoComputeObjective(tao, tao->solution, &ftrial);CHKERRQ(ierr);
182b4a30f08SAlp Dener           if (PetscIsInfOrNanReal(bnk->f)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
18362602cfbSAlp Dener           ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr);
184eb910715SAlp Dener           if (PetscIsInfOrNanReal(ftrial)) {
185eb910715SAlp Dener             tau = bnk->gamma1_i;
186eb910715SAlp Dener           } else {
1870a4511e9SAlp Dener             if (ftrial < f_min) {
1880a4511e9SAlp Dener               f_min = ftrial;
189eb910715SAlp Dener               sigma = -tao->trust / bnk->gnorm;
190eb910715SAlp Dener             }
19108752603SAlp Dener 
192770b7498SAlp Dener             /* Compute the predicted and actual reduction */
19389da521bSAlp Dener             if (bnk->active_idx) {
19408752603SAlp Dener               ierr = VecGetSubVector(bnk->W, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr);
19508752603SAlp Dener               ierr = VecGetSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work);CHKERRQ(ierr);
1962ab2a32cSAlp Dener             } else {
19708752603SAlp Dener               bnk->X_inactive = bnk->W;
19808752603SAlp Dener               bnk->inactive_work = bnk->Xwork;
1992ab2a32cSAlp Dener             }
20008752603SAlp Dener             ierr = MatMult(bnk->H_inactive, bnk->X_inactive, bnk->inactive_work);CHKERRQ(ierr);
20108752603SAlp Dener             ierr = VecDot(bnk->X_inactive, bnk->inactive_work, &prered);CHKERRQ(ierr);
20289da521bSAlp Dener             if (bnk->active_idx) {
20308752603SAlp Dener               ierr = VecRestoreSubVector(bnk->W, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr);
20408752603SAlp Dener               ierr = VecRestoreSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work);CHKERRQ(ierr);
2052ab2a32cSAlp Dener             }
206eb910715SAlp Dener             prered = tao->trust * (bnk->gnorm - 0.5 * tao->trust * prered / (bnk->gnorm * bnk->gnorm));
207eb910715SAlp Dener             actred = bnk->f - ftrial;
2083105154fSTodd Munson             if ((PetscAbsScalar(actred) <= bnk->epsilon) && (PetscAbsScalar(prered) <= bnk->epsilon)) {
209eb910715SAlp Dener               kappa = 1.0;
2103105154fSTodd Munson             } else {
211eb910715SAlp Dener               kappa = actred / prered;
212eb910715SAlp Dener             }
213eb910715SAlp Dener 
214eb910715SAlp Dener             tau_1 = bnk->theta_i * bnk->gnorm * tao->trust / (bnk->theta_i * bnk->gnorm * tao->trust + (1.0 - bnk->theta_i) * prered - actred);
215eb910715SAlp Dener             tau_2 = bnk->theta_i * bnk->gnorm * tao->trust / (bnk->theta_i * bnk->gnorm * tao->trust - (1.0 + bnk->theta_i) * prered + actred);
216eb910715SAlp Dener             tau_min = PetscMin(tau_1, tau_2);
217eb910715SAlp Dener             tau_max = PetscMax(tau_1, tau_2);
218eb910715SAlp Dener 
219eb910715SAlp Dener             if (PetscAbsScalar(kappa - 1.0) <= bnk->mu1_i) {
220eb910715SAlp Dener               /*  Great agreement */
221eb910715SAlp Dener               max_radius = PetscMax(max_radius, tao->trust);
222eb910715SAlp Dener 
223eb910715SAlp Dener               if (tau_max < 1.0) {
224eb910715SAlp Dener                 tau = bnk->gamma3_i;
2253105154fSTodd Munson               } else if (tau_max > bnk->gamma4_i) {
226eb910715SAlp Dener                 tau = bnk->gamma4_i;
2273105154fSTodd Munson               } else {
228eb910715SAlp Dener                 tau = tau_max;
229eb910715SAlp Dener               }
23008752603SAlp Dener             }
23108752603SAlp Dener             else if (PetscAbsScalar(kappa - 1.0) <= bnk->mu2_i) {
232eb910715SAlp Dener               /*  Good agreement */
233eb910715SAlp Dener               max_radius = PetscMax(max_radius, tao->trust);
234eb910715SAlp Dener 
235eb910715SAlp Dener               if (tau_max < bnk->gamma2_i) {
236eb910715SAlp Dener                 tau = bnk->gamma2_i;
237eb910715SAlp Dener               } else if (tau_max > bnk->gamma3_i) {
238eb910715SAlp Dener                 tau = bnk->gamma3_i;
239eb910715SAlp Dener               } else {
240eb910715SAlp Dener                 tau = tau_max;
241eb910715SAlp Dener               }
24208752603SAlp Dener             }
24308752603SAlp Dener             else {
244eb910715SAlp Dener               /*  Not good agreement */
245eb910715SAlp Dener               if (tau_min > 1.0) {
246eb910715SAlp Dener                 tau = bnk->gamma2_i;
247eb910715SAlp Dener               } else if (tau_max < bnk->gamma1_i) {
248eb910715SAlp Dener                 tau = bnk->gamma1_i;
249eb910715SAlp Dener               } else if ((tau_min < bnk->gamma1_i) && (tau_max >= 1.0)) {
250eb910715SAlp Dener                 tau = bnk->gamma1_i;
2513105154fSTodd Munson               } else if ((tau_1 >= bnk->gamma1_i) && (tau_1 < 1.0) && ((tau_2 < bnk->gamma1_i) || (tau_2 >= 1.0))) {
252eb910715SAlp Dener                 tau = tau_1;
2533105154fSTodd Munson               } else if ((tau_2 >= bnk->gamma1_i) && (tau_2 < 1.0) && ((tau_1 < bnk->gamma1_i) || (tau_2 >= 1.0))) {
254eb910715SAlp Dener                 tau = tau_2;
255eb910715SAlp Dener               } else {
256eb910715SAlp Dener                 tau = tau_max;
257eb910715SAlp Dener               }
258eb910715SAlp Dener             }
259eb910715SAlp Dener           }
260eb910715SAlp Dener           tao->trust = tau * tao->trust;
261eb910715SAlp Dener         }
262eb910715SAlp Dener 
2630a4511e9SAlp Dener         if (f_min < bnk->f) {
264937a31a1SAlp Dener           /* We accidentally found a solution better than the initial, so accept it */
2650a4511e9SAlp Dener           bnk->f = f_min;
266937a31a1SAlp Dener           ierr = VecCopy(tao->solution, bnk->Xold);CHKERRQ(ierr);
267eb910715SAlp Dener           ierr = VecAXPY(tao->solution,sigma,tao->gradient);CHKERRQ(ierr);
26889da521bSAlp Dener           ierr = TaoBoundSolution(tao->XL, tao->XU, tao->solution, 0.0, &nDiff);CHKERRQ(ierr);
269937a31a1SAlp Dener           ierr = VecCopy(tao->solution, tao->stepdirection);CHKERRQ(ierr);
270937a31a1SAlp Dener           ierr = VecAXPY(tao->stepdirection, -1.0, bnk->Xold);CHKERRQ(ierr);
27109164190SAlp Dener           ierr = TaoComputeGradient(tao,tao->solution,bnk->unprojected_gradient);CHKERRQ(ierr);
27261be54a6SAlp Dener           ierr = TaoBNKEstimateActiveSet(tao, bnk->as_type);CHKERRQ(ierr);
27361be54a6SAlp Dener           ierr = VecCopy(bnk->unprojected_gradient, tao->gradient);CHKERRQ(ierr);
27461be54a6SAlp Dener           ierr = VecISSet(tao->gradient, bnk->active_idx, 0.0);CHKERRQ(ierr);
275937a31a1SAlp Dener           /* Compute gradient at the new iterate and flip switch to compute the Hessian later */
276c4b75bccSAlp Dener           ierr = VecNorm(tao->gradient, NORM_2, &bnk->gnorm);CHKERRQ(ierr);
277c0f10754SAlp Dener           *needH = PETSC_TRUE;
278937a31a1SAlp Dener           /* Test the new step for convergence */
27989da521bSAlp Dener           ierr = VecFischer(tao->solution, bnk->unprojected_gradient, tao->XL, tao->XU, bnk->W);CHKERRQ(ierr);
28089da521bSAlp Dener           ierr = VecNorm(bnk->W, NORM_2, &resnorm);CHKERRQ(ierr);
281b4a30f08SAlp Dener           if (PetscIsInfOrNanReal(resnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
282e031d6f5SAlp Dener           ierr = TaoLogConvergenceHistory(tao,bnk->f,resnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
283e031d6f5SAlp Dener           ierr = TaoMonitor(tao,tao->niter,bnk->f,resnorm,0.0,1.0);CHKERRQ(ierr);
284eb910715SAlp Dener           ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
285eb910715SAlp Dener           if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
286937a31a1SAlp Dener           /* active BNCG recycling early because we have a stepdirection computed */
287937a31a1SAlp Dener           ierr = TaoBNCGSetRecycleFlag(bnk->bncg, PETSC_TRUE);CHKERRQ(ierr);
288eb910715SAlp Dener         }
289eb910715SAlp Dener       }
290eb910715SAlp Dener       tao->trust = PetscMax(tao->trust, max_radius);
291e031d6f5SAlp Dener 
292e031d6f5SAlp Dener       /* Ensure that the trust radius is within the limits */
293e031d6f5SAlp Dener       tao->trust = PetscMax(tao->trust, bnk->min_radius);
294e031d6f5SAlp Dener       tao->trust = PetscMin(tao->trust, bnk->max_radius);
295eb910715SAlp Dener       break;
296eb910715SAlp Dener 
297eb910715SAlp Dener     default:
298eb910715SAlp Dener       /* Norm of the first direction will initialize radius */
299eb910715SAlp Dener       tao->trust = 0.0;
300eb910715SAlp Dener       break;
301eb910715SAlp Dener     }
302eb910715SAlp Dener   }
303e031d6f5SAlp Dener 
304eb910715SAlp Dener   /* Set initial scaling for the BFGS preconditioner
305eb910715SAlp Dener      This step is done after computing the initial trust-region radius
306eb910715SAlp Dener      since the function value may have decreased */
307eb910715SAlp Dener   if (BNK_PC_BFGS == bnk->pc_type) {
3089b6ef848SAlp Dener     delta = 2.0 * PetscMax(1.0, PetscAbsScalar(bnk->f)) / (bnk->gnorm*bnk->gnorm);
309eb910715SAlp Dener     ierr = MatLMVMSetDelta(bnk->M,delta);CHKERRQ(ierr);
310eb910715SAlp Dener   }
311eb910715SAlp Dener   PetscFunctionReturn(0);
312eb910715SAlp Dener }
313eb910715SAlp Dener 
314df278d8fSAlp Dener /*------------------------------------------------------------*/
315df278d8fSAlp Dener 
31662675beeSAlp Dener /* Routine for computing the Hessian and preparing the preconditioner at the new iterate */
31762675beeSAlp Dener 
31862675beeSAlp Dener PetscErrorCode TaoBNKComputeHessian(Tao tao)
31962675beeSAlp Dener {
32062675beeSAlp Dener   PetscErrorCode               ierr;
32162675beeSAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
322c4b75bccSAlp Dener   PetscBool                    diagExists;
323c4b75bccSAlp Dener   PetscReal                    delta;
32462675beeSAlp Dener 
32562675beeSAlp Dener   PetscFunctionBegin;
32662675beeSAlp Dener   /* Compute the Hessian */
32762675beeSAlp Dener   ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr);
32862675beeSAlp Dener   /* Add a correction to the BFGS preconditioner */
32962675beeSAlp Dener   if (BNK_PC_BFGS == bnk->pc_type) {
33062675beeSAlp Dener     ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
33162675beeSAlp Dener     /* Update the BFGS diagonal scaling */
332c4b75bccSAlp Dener     ierr = MatHasOperation(tao->hessian, MATOP_GET_DIAGONAL, &diagExists);CHKERRQ(ierr);
333c4b75bccSAlp Dener     if (BFGS_SCALE_AHESS == bnk->bfgs_scale_type && diagExists) {
33462675beeSAlp Dener       ierr = MatGetDiagonal(tao->hessian, bnk->Diag);CHKERRQ(ierr);
33562675beeSAlp Dener       ierr = VecAbs(bnk->Diag);CHKERRQ(ierr);
33662675beeSAlp Dener       ierr = VecMedian(bnk->Diag_min, bnk->Diag, bnk->Diag_max, bnk->Diag);CHKERRQ(ierr);
33762675beeSAlp Dener       ierr = VecReciprocal(bnk->Diag);CHKERRQ(ierr);
33862675beeSAlp Dener       ierr = MatLMVMSetScale(bnk->M,bnk->Diag);CHKERRQ(ierr);
339c4b75bccSAlp Dener     } else {
340c4b75bccSAlp Dener       /* Fall back onto stand-alone BFGS scaling */
341c4b75bccSAlp Dener       delta = 2.0 * PetscMax(1.0, PetscAbsScalar(bnk->f)) / (bnk->gnorm*bnk->gnorm);
342c4b75bccSAlp Dener       ierr = MatLMVMSetDelta(bnk->M,delta);CHKERRQ(ierr);
34362675beeSAlp Dener     }
34462675beeSAlp Dener   }
34562675beeSAlp Dener   PetscFunctionReturn(0);
34662675beeSAlp Dener }
34762675beeSAlp Dener 
34862675beeSAlp Dener /*------------------------------------------------------------*/
34962675beeSAlp Dener 
3502f75a4aaSAlp Dener /* Routine for estimating the active set */
3512f75a4aaSAlp Dener 
35208752603SAlp Dener PetscErrorCode TaoBNKEstimateActiveSet(Tao tao, PetscInt asType)
3532f75a4aaSAlp Dener {
3542f75a4aaSAlp Dener   PetscErrorCode               ierr;
3552f75a4aaSAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
356c4b75bccSAlp Dener   PetscBool                    hessComputed, diagExists;
3572f75a4aaSAlp Dener 
3582f75a4aaSAlp Dener   PetscFunctionBegin;
35908752603SAlp Dener   switch (asType) {
3602f75a4aaSAlp Dener   case BNK_AS_NONE:
3612f75a4aaSAlp Dener     ierr = ISDestroy(&bnk->inactive_idx);CHKERRQ(ierr);
3622f75a4aaSAlp Dener     ierr = VecWhichInactive(tao->XL, tao->solution, bnk->unprojected_gradient, tao->XU, PETSC_TRUE, &bnk->inactive_idx);CHKERRQ(ierr);
3632f75a4aaSAlp Dener     ierr = ISDestroy(&bnk->active_idx);CHKERRQ(ierr);
3642f75a4aaSAlp Dener     ierr = ISComplementVec(bnk->inactive_idx, tao->solution, &bnk->active_idx);CHKERRQ(ierr);
3652f75a4aaSAlp Dener     break;
3662f75a4aaSAlp Dener 
3672f75a4aaSAlp Dener   case BNK_AS_BERTSEKAS:
3682f75a4aaSAlp Dener     /* Compute the trial step vector with which we will estimate the active set at the next iteration */
3692f75a4aaSAlp Dener     if (BNK_PC_BFGS == bnk->pc_type) {
3702f75a4aaSAlp Dener       /* If the BFGS preconditioner matrix is available, we will construct a trial step with it */
3715e9b73cbSAlp Dener       ierr = MatLMVMSetInactive(bnk->M, NULL);CHKERRQ(ierr);
3722f75a4aaSAlp Dener       ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, bnk->W);CHKERRQ(ierr);
3732f75a4aaSAlp Dener     } else {
37461be54a6SAlp Dener       ierr = MatAssembled(tao->hessian, &hessComputed);CHKERRQ(ierr);
375c4b75bccSAlp Dener       ierr = MatHasOperation(tao->hessian, MATOP_GET_DIAGONAL, &diagExists);CHKERRQ(ierr);
376c4b75bccSAlp Dener       if (hessComputed && diagExists) {
3779b6ef848SAlp Dener         /* BFGS preconditioner doesn't exist so let's invert the absolute diagonal of the Hessian instead onto the gradient */
3782f75a4aaSAlp Dener         ierr = MatGetDiagonal(tao->hessian, bnk->Xwork);CHKERRQ(ierr);
3799b6ef848SAlp Dener         ierr = VecAbs(bnk->Xwork);CHKERRQ(ierr);
3809b6ef848SAlp Dener         ierr = VecMedian(bnk->Diag_min, bnk->Xwork, bnk->Diag_max, bnk->Xwork);CHKERRQ(ierr);
3812f75a4aaSAlp Dener         ierr = VecReciprocal(bnk->Xwork);CHKERRQ(ierr);CHKERRQ(ierr);
3822f75a4aaSAlp Dener         ierr = VecPointwiseMult(bnk->W, bnk->Xwork, bnk->unprojected_gradient);CHKERRQ(ierr);
38361be54a6SAlp Dener       } else {
384c4b75bccSAlp Dener         /* If the Hessian or its diagonal does not exist, we will simply use gradient step */
38561be54a6SAlp Dener         ierr = VecCopy(bnk->unprojected_gradient, bnk->W);CHKERRQ(ierr);
38661be54a6SAlp Dener       }
3872f75a4aaSAlp Dener     }
3882f75a4aaSAlp Dener     ierr = VecScale(bnk->W, -1.0);CHKERRQ(ierr);
38989da521bSAlp Dener     ierr = TaoEstimateActiveBounds(tao->solution, tao->XL, tao->XU, bnk->unprojected_gradient, bnk->W, bnk->Xwork, bnk->as_step, &bnk->as_tol,
39089da521bSAlp Dener                                    &bnk->active_lower, &bnk->active_upper, &bnk->active_fixed, &bnk->active_idx, &bnk->inactive_idx);CHKERRQ(ierr);
391c4b75bccSAlp Dener     break;
3922f75a4aaSAlp Dener 
3932f75a4aaSAlp Dener   default:
3942f75a4aaSAlp Dener     break;
3952f75a4aaSAlp Dener   }
3962f75a4aaSAlp Dener   PetscFunctionReturn(0);
3972f75a4aaSAlp Dener }
3982f75a4aaSAlp Dener 
3992f75a4aaSAlp Dener /*------------------------------------------------------------*/
4002f75a4aaSAlp Dener 
4012f75a4aaSAlp Dener /* Routine for bounding the step direction */
4022f75a4aaSAlp Dener 
403a1318120SAlp Dener PetscErrorCode TaoBNKBoundStep(Tao tao, PetscInt asType, Vec step)
4042f75a4aaSAlp Dener {
4052f75a4aaSAlp Dener   PetscErrorCode               ierr;
4062f75a4aaSAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
4072f75a4aaSAlp Dener 
4082f75a4aaSAlp Dener   PetscFunctionBegin;
409a1318120SAlp Dener   switch (asType) {
4102f75a4aaSAlp Dener   case BNK_AS_NONE:
411c4b75bccSAlp Dener     ierr = VecISSet(step, bnk->active_idx, 0.0);CHKERRQ(ierr);
4122f75a4aaSAlp Dener     break;
4132f75a4aaSAlp Dener 
4142f75a4aaSAlp Dener   case BNK_AS_BERTSEKAS:
415c4b75bccSAlp Dener     ierr = TaoBoundStep(tao->solution, tao->XL, tao->XU, bnk->active_lower, bnk->active_upper, bnk->active_fixed, 1.0, step);CHKERRQ(ierr);
4162f75a4aaSAlp Dener     break;
4172f75a4aaSAlp Dener 
4182f75a4aaSAlp Dener   default:
4192f75a4aaSAlp Dener     break;
4202f75a4aaSAlp Dener   }
4212f75a4aaSAlp Dener   PetscFunctionReturn(0);
4222f75a4aaSAlp Dener }
4232f75a4aaSAlp Dener 
424e031d6f5SAlp Dener /*------------------------------------------------------------*/
425e031d6f5SAlp Dener 
426e031d6f5SAlp Dener /* Routine for taking a finite number of BNCG iterations to
427e031d6f5SAlp Dener    accelerate Newton convergence.
428e031d6f5SAlp Dener 
429e031d6f5SAlp Dener    In practice, this approach simply trades off Hessian evaluations
430e031d6f5SAlp Dener    for more gradient evaluations.
431e031d6f5SAlp Dener */
432e031d6f5SAlp Dener 
433c0f10754SAlp Dener PetscErrorCode TaoBNKTakeCGSteps(Tao tao, PetscBool *terminate)
434c0f10754SAlp Dener {
435c0f10754SAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
436c0f10754SAlp Dener   PetscErrorCode               ierr;
437c0f10754SAlp Dener 
438c0f10754SAlp Dener   PetscFunctionBegin;
439c0f10754SAlp Dener   *terminate = PETSC_FALSE;
440c0f10754SAlp Dener   if (bnk->max_cg_its > 0) {
441c4b75bccSAlp Dener     /* Copy the current function value (important vectors are already shared) */
442c0f10754SAlp Dener     bnk->bncg_ctx->f = bnk->f;
443c0f10754SAlp Dener     /* Take some small finite number of BNCG iterations */
444c0f10754SAlp Dener     ierr = TaoSolve(bnk->bncg);CHKERRQ(ierr);
445c0f10754SAlp Dener     /* Add the number of gradient and function evaluations to the total */
446c0f10754SAlp Dener     tao->nfuncs += bnk->bncg->nfuncs;
447c0f10754SAlp Dener     tao->nfuncgrads += bnk->bncg->nfuncgrads;
448c0f10754SAlp Dener     tao->ngrads += bnk->bncg->ngrads;
449c0f10754SAlp Dener     tao->nhess += bnk->bncg->nhess;
450e031d6f5SAlp Dener     bnk->tot_cg_its += bnk->bncg->niter;
451c4b75bccSAlp Dener     /* Extract the BNCG function value out and save it into BNK */
452c0f10754SAlp Dener     bnk->f = bnk->bncg_ctx->f;
453c0f10754SAlp Dener     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) {
454c0f10754SAlp Dener       *terminate = PETSC_TRUE;
45561be54a6SAlp Dener     } else {
45661be54a6SAlp Dener       ierr = TaoBNKEstimateActiveSet(tao, bnk->as_type);
457c0f10754SAlp Dener     }
458c0f10754SAlp Dener   }
459c0f10754SAlp Dener   PetscFunctionReturn(0);
460c0f10754SAlp Dener }
461c0f10754SAlp Dener 
4622f75a4aaSAlp Dener /*------------------------------------------------------------*/
4632f75a4aaSAlp Dener 
464c0f10754SAlp Dener /* Routine for computing the Newton step. */
465df278d8fSAlp Dener 
46662675beeSAlp Dener PetscErrorCode TaoBNKComputeStep(Tao tao, PetscBool shift, KSPConvergedReason *ksp_reason)
467eb910715SAlp Dener {
468eb910715SAlp Dener   PetscErrorCode               ierr;
469eb910715SAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
470eb910715SAlp Dener 
471e465cd6fSAlp Dener   PetscReal                    delta;
472eb910715SAlp Dener   PetscInt                     bfgsUpdates = 0;
473eb910715SAlp Dener   PetscInt                     kspits;
474c4b75bccSAlp Dener   PetscBool                    diagExists;
475eb910715SAlp Dener 
476eb910715SAlp Dener   PetscFunctionBegin;
47789da521bSAlp Dener   /* If there are no inactive variables left, save some computation and return an adjusted zero step
47889da521bSAlp Dener      that has (l-x) and (u-x) for lower and upper bounded variables. */
47989da521bSAlp Dener   if (!bnk->inactive_idx) {
48089da521bSAlp Dener     ierr = VecSet(tao->stepdirection, 0.0);CHKERRQ(ierr);
481a1318120SAlp Dener     ierr = TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection);CHKERRQ(ierr);
48289da521bSAlp Dener     PetscFunctionReturn(0);
48389da521bSAlp Dener   }
48489da521bSAlp Dener 
4855e9b73cbSAlp Dener   /* Prepare the reduced sub-matrices for the inactive set */
4863105154fSTodd Munson   if (BNK_PC_BFGS == bnk->pc_type) {
4873105154fSTodd Munson     ierr = MatLMVMSetInactive(bnk->M, bnk->inactive_idx);CHKERRQ(ierr);
4883105154fSTodd Munson   }
48989da521bSAlp Dener   if (bnk->active_idx) {
4905e9b73cbSAlp Dener     ierr = MatDestroy(&bnk->H_inactive);
4915e9b73cbSAlp Dener     ierr = MatCreateSubMatrix(tao->hessian, bnk->inactive_idx, bnk->inactive_idx, MAT_INITIAL_MATRIX, &bnk->H_inactive);CHKERRQ(ierr);
492eb3ba6a7SAlp Dener     if (tao->hessian == tao->hessian_pre) {
493eb3ba6a7SAlp Dener       bnk->Hpre_inactive = bnk->H_inactive;
494eb3ba6a7SAlp Dener     } else {
4955e9b73cbSAlp Dener       ierr = MatDestroy(&bnk->Hpre_inactive);
4965e9b73cbSAlp Dener       ierr = MatCreateSubMatrix(tao->hessian_pre, bnk->inactive_idx, bnk->inactive_idx, MAT_INITIAL_MATRIX, &bnk->Hpre_inactive);CHKERRQ(ierr);
497eb3ba6a7SAlp Dener     }
4982f75a4aaSAlp Dener   } else {
49962675beeSAlp Dener     ierr = MatDestroy(&bnk->H_inactive);
50062675beeSAlp Dener     ierr = MatDuplicate(tao->hessian, MAT_COPY_VALUES, &bnk->H_inactive);
50162675beeSAlp Dener     if (tao->hessian == tao->hessian_pre) {
50262675beeSAlp Dener       bnk->Hpre_inactive = bnk->H_inactive;
50362675beeSAlp Dener     } else {
50462675beeSAlp Dener       ierr = MatDestroy(&bnk->Hpre_inactive);
50562675beeSAlp Dener       ierr = MatDuplicate(tao->hessian_pre, MAT_COPY_VALUES, &bnk->Hpre_inactive);
50662675beeSAlp Dener     }
50762675beeSAlp Dener   }
50862675beeSAlp Dener 
50962675beeSAlp Dener   /* Shift the reduced Hessian matrix */
51062675beeSAlp Dener   if ((shift) && (bnk->pert > 0)) {
51162675beeSAlp Dener     ierr = MatShift(bnk->H_inactive, bnk->pert);CHKERRQ(ierr);
51262675beeSAlp Dener     if (bnk->H_inactive != bnk->Hpre_inactive) {
51362675beeSAlp Dener       ierr = MatShift(bnk->Hpre_inactive, bnk->pert);CHKERRQ(ierr);
51462675beeSAlp Dener     }
51562675beeSAlp Dener   }
51662675beeSAlp Dener 
51762675beeSAlp Dener   /* Update the diagonal scaling for the BFGS preconditioner, this time with the Hessian perturbation */
51862675beeSAlp Dener   if ((BNK_PC_BFGS == bnk->pc_type) && (BFGS_SCALE_PHESS == bnk->bfgs_scale_type)) {
519c4b75bccSAlp Dener     ierr = MatHasOperation(bnk->H_inactive, MATOP_GET_DIAGONAL, &diagExists);CHKERRQ(ierr);
520c4b75bccSAlp Dener     if (diagExists) {
52162675beeSAlp Dener       /* Obtain diagonal for the bfgs preconditioner  */
5225e9b73cbSAlp Dener       ierr = VecSet(bnk->Diag, 1.0);CHKERRQ(ierr);
52389da521bSAlp Dener       if (bnk->active_idx) {
5245e9b73cbSAlp Dener         ierr = VecGetSubVector(bnk->Diag, bnk->inactive_idx, &bnk->Diag_red);CHKERRQ(ierr);
5255e9b73cbSAlp Dener       } else {
5265e9b73cbSAlp Dener         bnk->Diag_red = bnk->Diag;
5275e9b73cbSAlp Dener       }
5285e9b73cbSAlp Dener       ierr = MatGetDiagonal(bnk->H_inactive, bnk->Diag_red);CHKERRQ(ierr);
52989da521bSAlp Dener       if (bnk->active_idx) {
5305e9b73cbSAlp Dener         ierr = VecRestoreSubVector(bnk->Diag, bnk->inactive_idx, &bnk->Diag_red);CHKERRQ(ierr);
5315e9b73cbSAlp Dener       }
53262675beeSAlp Dener       ierr = VecAbs(bnk->Diag);CHKERRQ(ierr);
53362675beeSAlp Dener       ierr = VecMedian(bnk->Diag_min, bnk->Diag, bnk->Diag_max, bnk->Diag);CHKERRQ(ierr);
53462675beeSAlp Dener       ierr = VecReciprocal(bnk->Diag);CHKERRQ(ierr);
53562675beeSAlp Dener       ierr = MatLMVMSetScale(bnk->M,bnk->Diag);CHKERRQ(ierr);
536c4b75bccSAlp Dener     } else {
537c4b75bccSAlp Dener       /* Fall back onto stand-alone BFGS scaling */
538c4b75bccSAlp Dener       delta = 2.0 * PetscMax(1.0, PetscAbsScalar(bnk->f)) / (bnk->gnorm*bnk->gnorm);
539c4b75bccSAlp Dener       ierr = MatLMVMSetDelta(bnk->M,delta);CHKERRQ(ierr);
540c4b75bccSAlp Dener     }
5412f75a4aaSAlp Dener   }
54209164190SAlp Dener 
543eb910715SAlp Dener   /* Solve the Newton system of equations */
544937a31a1SAlp Dener   tao->ksp_its = 0;
5452f75a4aaSAlp Dener   ierr = VecSet(tao->stepdirection, 0.0);CHKERRQ(ierr);
5465e9b73cbSAlp Dener   ierr = KSPReset(tao->ksp);CHKERRQ(ierr);
54709164190SAlp Dener   ierr = KSPSetOperators(tao->ksp,bnk->H_inactive,bnk->Hpre_inactive);CHKERRQ(ierr);
5485e9b73cbSAlp Dener   ierr = VecCopy(bnk->unprojected_gradient, bnk->Gwork);CHKERRQ(ierr);
54989da521bSAlp Dener   if (bnk->active_idx) {
5505e9b73cbSAlp Dener     ierr = VecGetSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive);CHKERRQ(ierr);
5515e9b73cbSAlp Dener     ierr = VecGetSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr);
5525e9b73cbSAlp Dener   } else {
5535e9b73cbSAlp Dener     bnk->G_inactive = bnk->unprojected_gradient;
5545e9b73cbSAlp Dener     bnk->X_inactive = tao->stepdirection;
55528017e9fSAlp Dener   }
556eb910715SAlp Dener   if (bnk->is_nash || bnk->is_stcg || bnk->is_gltr) {
557fed79b8eSAlp Dener     ierr = KSPCGSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr);
5585e9b73cbSAlp Dener     ierr = KSPSolve(tao->ksp, bnk->G_inactive, bnk->X_inactive);CHKERRQ(ierr);
559eb910715SAlp Dener     ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr);
560eb910715SAlp Dener     tao->ksp_its+=kspits;
561eb910715SAlp Dener     tao->ksp_tot_its+=kspits;
562080d2917SAlp Dener     ierr = KSPCGGetNormD(tao->ksp,&bnk->dnorm);CHKERRQ(ierr);
563eb910715SAlp Dener 
564eb910715SAlp Dener     if (0.0 == tao->trust) {
565eb910715SAlp Dener       /* Radius was uninitialized; use the norm of the direction */
566080d2917SAlp Dener       if (bnk->dnorm > 0.0) {
567080d2917SAlp Dener         tao->trust = bnk->dnorm;
568eb910715SAlp Dener 
569eb910715SAlp Dener         /* Modify the radius if it is too large or small */
570eb910715SAlp Dener         tao->trust = PetscMax(tao->trust, bnk->min_radius);
571eb910715SAlp Dener         tao->trust = PetscMin(tao->trust, bnk->max_radius);
572eb910715SAlp Dener       } else {
573eb910715SAlp Dener         /* The direction was bad; set radius to default value and re-solve
574eb910715SAlp Dener            the trust-region subproblem to get a direction */
575eb910715SAlp Dener         tao->trust = tao->trust0;
576eb910715SAlp Dener 
577eb910715SAlp Dener         /* Modify the radius if it is too large or small */
578eb910715SAlp Dener         tao->trust = PetscMax(tao->trust, bnk->min_radius);
579eb910715SAlp Dener         tao->trust = PetscMin(tao->trust, bnk->max_radius);
580eb910715SAlp Dener 
581fed79b8eSAlp Dener         ierr = KSPCGSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr);
5825e9b73cbSAlp Dener         ierr = KSPSolve(tao->ksp, bnk->G_inactive, bnk->X_inactive);CHKERRQ(ierr);
583eb910715SAlp Dener         ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr);
584eb910715SAlp Dener         tao->ksp_its+=kspits;
585eb910715SAlp Dener         tao->ksp_tot_its+=kspits;
586080d2917SAlp Dener         ierr = KSPCGGetNormD(tao->ksp,&bnk->dnorm);CHKERRQ(ierr);
587eb910715SAlp Dener 
588080d2917SAlp Dener         if (bnk->dnorm == 0.0) SETERRQ(PETSC_COMM_SELF,1, "Initial direction zero");
589eb910715SAlp Dener       }
590eb910715SAlp Dener     }
591eb910715SAlp Dener   } else {
5925e9b73cbSAlp Dener     ierr = KSPSolve(tao->ksp, bnk->G_inactive, bnk->X_inactive);CHKERRQ(ierr);
593eb910715SAlp Dener     ierr = KSPGetIterationNumber(tao->ksp, &kspits);CHKERRQ(ierr);
594eb910715SAlp Dener     tao->ksp_its += kspits;
595eb910715SAlp Dener     tao->ksp_tot_its+=kspits;
596eb910715SAlp Dener   }
5975e9b73cbSAlp Dener   /* Restore sub vectors back */
59889da521bSAlp Dener   if (bnk->active_idx) {
5995e9b73cbSAlp Dener     ierr = VecRestoreSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive);CHKERRQ(ierr);
6005e9b73cbSAlp Dener     ierr = VecRestoreSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr);
6015e9b73cbSAlp Dener   }
602770b7498SAlp Dener   /* Make sure the safeguarded fall-back step is zero for actively bounded variables */
603fed79b8eSAlp Dener   ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
604a1318120SAlp Dener   ierr = TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection);CHKERRQ(ierr);
605770b7498SAlp Dener 
606770b7498SAlp Dener   /* Record convergence reasons */
607e465cd6fSAlp Dener   ierr = KSPGetConvergedReason(tao->ksp, ksp_reason);CHKERRQ(ierr);
608e465cd6fSAlp Dener   if (KSP_CONVERGED_ATOL == *ksp_reason) {
609770b7498SAlp Dener     ++bnk->ksp_atol;
610e465cd6fSAlp Dener   } else if (KSP_CONVERGED_RTOL == *ksp_reason) {
611770b7498SAlp Dener     ++bnk->ksp_rtol;
612e465cd6fSAlp Dener   } else if (KSP_CONVERGED_CG_CONSTRAINED == *ksp_reason) {
613770b7498SAlp Dener     ++bnk->ksp_ctol;
614e465cd6fSAlp Dener   } else if (KSP_CONVERGED_CG_NEG_CURVE == *ksp_reason) {
615770b7498SAlp Dener     ++bnk->ksp_negc;
616e465cd6fSAlp Dener   } else if (KSP_DIVERGED_DTOL == *ksp_reason) {
617770b7498SAlp Dener     ++bnk->ksp_dtol;
618e465cd6fSAlp Dener   } else if (KSP_DIVERGED_ITS == *ksp_reason) {
619770b7498SAlp Dener     ++bnk->ksp_iter;
620770b7498SAlp Dener   } else {
621770b7498SAlp Dener     ++bnk->ksp_othr;
622770b7498SAlp Dener   }
623fed79b8eSAlp Dener 
624fed79b8eSAlp Dener   /* Make sure the BFGS preconditioner is healthy */
625fed79b8eSAlp Dener   if (bnk->pc_type == BNK_PC_BFGS) {
626fed79b8eSAlp Dener     ierr = MatLMVMGetUpdates(bnk->M, &bfgsUpdates);CHKERRQ(ierr);
627e465cd6fSAlp Dener     if ((KSP_DIVERGED_INDEFINITE_PC == *ksp_reason) && (bfgsUpdates > 1)) {
628fed79b8eSAlp Dener       /* Preconditioner is numerically indefinite; reset the approximation. */
6299b6ef848SAlp Dener       delta = 2.0 * PetscMax(1.0, PetscAbsScalar(bnk->f)) / (bnk->gnorm*bnk->gnorm);
630eb910715SAlp Dener       ierr = MatLMVMSetDelta(bnk->M,delta);CHKERRQ(ierr);
631eb910715SAlp Dener       ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr);
63209164190SAlp Dener       ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
633eb910715SAlp Dener     }
634fed79b8eSAlp Dener   }
635e465cd6fSAlp Dener   PetscFunctionReturn(0);
636e465cd6fSAlp Dener }
637eb910715SAlp Dener 
63862675beeSAlp Dener /*------------------------------------------------------------*/
63962675beeSAlp Dener 
6405e9b73cbSAlp Dener /* Routine for recomputing the predicted reduction for a given step vector */
6415e9b73cbSAlp Dener 
6425e9b73cbSAlp Dener PetscErrorCode TaoBNKRecomputePred(Tao tao, Vec S, PetscReal *prered)
6435e9b73cbSAlp Dener {
6445e9b73cbSAlp Dener   PetscErrorCode               ierr;
6455e9b73cbSAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
6465e9b73cbSAlp Dener 
6475e9b73cbSAlp Dener   PetscFunctionBegin;
6485e9b73cbSAlp Dener   /* Extract subvectors associated with the inactive set */
64989da521bSAlp Dener   if (bnk->active_idx){
6505e9b73cbSAlp Dener     ierr = VecGetSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr);
6515e9b73cbSAlp Dener     ierr = VecGetSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work);CHKERRQ(ierr);
6525e9b73cbSAlp Dener     ierr = VecGetSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive);CHKERRQ(ierr);
6535e9b73cbSAlp Dener   } else {
6545e9b73cbSAlp Dener     bnk->X_inactive = tao->stepdirection;
6555e9b73cbSAlp Dener     bnk->inactive_work = bnk->Xwork;
6565e9b73cbSAlp Dener     bnk->G_inactive = bnk->Gwork;
6575e9b73cbSAlp Dener   }
6585e9b73cbSAlp Dener   /* Recompute the predicted decrease based on the quadratic model */
6595e9b73cbSAlp Dener   ierr = MatMult(bnk->H_inactive, bnk->X_inactive, bnk->inactive_work);CHKERRQ(ierr);
6605e9b73cbSAlp Dener   ierr = VecAYPX(bnk->inactive_work, -0.5, bnk->G_inactive);CHKERRQ(ierr);
6615e9b73cbSAlp Dener   ierr = VecDot(bnk->inactive_work, bnk->X_inactive, prered);
6625e9b73cbSAlp Dener   /* Restore the sub vectors */
66389da521bSAlp Dener   if (bnk->active_idx){
6645e9b73cbSAlp Dener     ierr = VecRestoreSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr);
6655e9b73cbSAlp Dener     ierr = VecRestoreSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work);CHKERRQ(ierr);
6665e9b73cbSAlp Dener     ierr = VecRestoreSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive);CHKERRQ(ierr);
6675e9b73cbSAlp Dener   }
6685e9b73cbSAlp Dener   PetscFunctionReturn(0);
6695e9b73cbSAlp Dener }
6705e9b73cbSAlp Dener 
6715e9b73cbSAlp Dener /*------------------------------------------------------------*/
6725e9b73cbSAlp Dener 
67362675beeSAlp Dener /* Routine for ensuring that the Newton step is a descent direction.
67462675beeSAlp Dener 
67562675beeSAlp Dener    The step direction falls back onto BFGS, scaled gradient and gradient steps
67662675beeSAlp Dener    in the event that the Newton step fails the test.
67762675beeSAlp Dener */
67862675beeSAlp Dener 
679e465cd6fSAlp Dener PetscErrorCode TaoBNKSafeguardStep(Tao tao, KSPConvergedReason ksp_reason, PetscInt *stepType)
680e465cd6fSAlp Dener {
681e465cd6fSAlp Dener   PetscErrorCode               ierr;
682e465cd6fSAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
683e465cd6fSAlp Dener 
684e465cd6fSAlp Dener   PetscReal                    gdx, delta, e_min;
685e465cd6fSAlp Dener   PetscInt                     bfgsUpdates;
686e465cd6fSAlp Dener 
687e465cd6fSAlp Dener   PetscFunctionBegin;
688080d2917SAlp Dener   ierr = VecDot(tao->stepdirection, tao->gradient, &gdx);CHKERRQ(ierr);
689eb910715SAlp Dener   if ((gdx >= 0.0) || PetscIsInfOrNanReal(gdx)) {
690eb910715SAlp Dener     /* Newton step is not descent or direction produced Inf or NaN
691eb910715SAlp Dener        Update the perturbation for next time */
692eb910715SAlp Dener     if (bnk->pert <= 0.0) {
693eb910715SAlp Dener       /* Initialize the perturbation */
694eb910715SAlp Dener       bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
695eb910715SAlp Dener       if (bnk->is_gltr) {
696eb910715SAlp Dener         ierr = KSPCGGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr);
697eb910715SAlp Dener         bnk->pert = PetscMax(bnk->pert, -e_min);
698eb910715SAlp Dener       }
699eb910715SAlp Dener     } else {
700eb910715SAlp Dener       /* Increase the perturbation */
701eb910715SAlp Dener       bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
702eb910715SAlp Dener     }
703eb910715SAlp Dener 
704eb910715SAlp Dener     if (BNK_PC_BFGS != bnk->pc_type) {
705eb910715SAlp Dener       /* We don't have the bfgs matrix around and updated
706eb910715SAlp Dener          Must use gradient direction in this case */
707080d2917SAlp Dener       ierr = VecCopy(tao->gradient, tao->stepdirection);CHKERRQ(ierr);
708eb910715SAlp Dener       *stepType = BNK_GRADIENT;
709eb910715SAlp Dener     } else {
710eb910715SAlp Dener       /* Attempt to use the BFGS direction */
7112ab2a32cSAlp Dener       ierr = MatLMVMSetInactive(bnk->M, NULL);CHKERRQ(ierr);
71209164190SAlp Dener       ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
713eb910715SAlp Dener 
7148d5ead36SAlp Dener       /* Check for success (descent direction)
7158d5ead36SAlp Dener          NOTE: Negative gdx here means not a descent direction because
7168d5ead36SAlp Dener          the fall-back step is missing a negative sign. */
717080d2917SAlp Dener       ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr);
7183105154fSTodd Munson       if ((gdx <= 0.0) || PetscIsInfOrNanReal(gdx)) {
719eb910715SAlp Dener         /* BFGS direction is not descent or direction produced not a number
720eb910715SAlp Dener            We can assert bfgsUpdates > 1 in this case because
721eb910715SAlp Dener            the first solve produces the scaled gradient direction,
722eb910715SAlp Dener            which is guaranteed to be descent */
723eb910715SAlp Dener 
724eb910715SAlp Dener         /* Use steepest descent direction (scaled) */
7259b6ef848SAlp Dener         delta = 2.0 * PetscMax(1.0, PetscAbsScalar(bnk->f)) / (bnk->gnorm*bnk->gnorm);
726eb910715SAlp Dener         ierr = MatLMVMSetDelta(bnk->M, delta);CHKERRQ(ierr);
727eb910715SAlp Dener         ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr);
72809164190SAlp Dener         ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
72909164190SAlp Dener         ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
730eb910715SAlp Dener 
731eb910715SAlp Dener         *stepType = BNK_SCALED_GRADIENT;
732eb910715SAlp Dener       } else {
733770b7498SAlp Dener         ierr = MatLMVMGetUpdates(bnk->M, &bfgsUpdates);CHKERRQ(ierr);
734eb910715SAlp Dener         if (1 == bfgsUpdates) {
735eb910715SAlp Dener           /* The first BFGS direction is always the scaled gradient */
736eb910715SAlp Dener           *stepType = BNK_SCALED_GRADIENT;
737eb910715SAlp Dener         } else {
738eb910715SAlp Dener           *stepType = BNK_BFGS;
739eb910715SAlp Dener         }
740eb910715SAlp Dener       }
741eb910715SAlp Dener     }
7428d5ead36SAlp Dener     /* Make sure the safeguarded fall-back step is zero for actively bounded variables */
7438d5ead36SAlp Dener     ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
744a1318120SAlp Dener     ierr = TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection);CHKERRQ(ierr);
745eb910715SAlp Dener   } else {
746eb910715SAlp Dener     /* Computed Newton step is descent */
747eb910715SAlp Dener     switch (ksp_reason) {
748eb910715SAlp Dener     case KSP_DIVERGED_NANORINF:
749eb910715SAlp Dener     case KSP_DIVERGED_BREAKDOWN:
750eb910715SAlp Dener     case KSP_DIVERGED_INDEFINITE_MAT:
751eb910715SAlp Dener     case KSP_DIVERGED_INDEFINITE_PC:
752eb910715SAlp Dener     case KSP_CONVERGED_CG_NEG_CURVE:
753eb910715SAlp Dener       /* Matrix or preconditioner is indefinite; increase perturbation */
754eb910715SAlp Dener       if (bnk->pert <= 0.0) {
755eb910715SAlp Dener         /* Initialize the perturbation */
756eb910715SAlp Dener         bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
757eb910715SAlp Dener         if (bnk->is_gltr) {
758eb910715SAlp Dener           ierr = KSPCGGLTRGetMinEig(tao->ksp, &e_min);CHKERRQ(ierr);
759eb910715SAlp Dener           bnk->pert = PetscMax(bnk->pert, -e_min);
760eb910715SAlp Dener         }
761eb910715SAlp Dener       } else {
762eb910715SAlp Dener         /* Increase the perturbation */
763eb910715SAlp Dener         bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
764eb910715SAlp Dener       }
765eb910715SAlp Dener       break;
766eb910715SAlp Dener 
767eb910715SAlp Dener     default:
768eb910715SAlp Dener       /* Newton step computation is good; decrease perturbation */
769eb910715SAlp Dener       bnk->pert = PetscMin(bnk->psfac * bnk->pert, bnk->pmsfac * bnk->gnorm);
770eb910715SAlp Dener       if (bnk->pert < bnk->pmin) {
771eb910715SAlp Dener         bnk->pert = 0.0;
772eb910715SAlp Dener       }
773eb910715SAlp Dener       break;
774eb910715SAlp Dener     }
775fed79b8eSAlp Dener     *stepType = BNK_NEWTON;
776eb910715SAlp Dener   }
777eb910715SAlp Dener   PetscFunctionReturn(0);
778eb910715SAlp Dener }
779eb910715SAlp Dener 
780df278d8fSAlp Dener /*------------------------------------------------------------*/
781df278d8fSAlp Dener 
782df278d8fSAlp Dener /* Routine for performing a bound-projected More-Thuente line search.
783df278d8fSAlp Dener 
784df278d8fSAlp Dener   Includes fallbacks to BFGS, scaled gradient, and unscaled gradient steps if the
785df278d8fSAlp Dener   Newton step does not produce a valid step length.
786df278d8fSAlp Dener */
787df278d8fSAlp Dener 
788937a31a1SAlp Dener PetscErrorCode TaoBNKPerformLineSearch(Tao tao, PetscInt *stepType, PetscReal *steplen, TaoLineSearchConvergedReason *reason)
789c14b763aSAlp Dener {
790c14b763aSAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
791c14b763aSAlp Dener   PetscErrorCode ierr;
792c14b763aSAlp Dener   TaoLineSearchConvergedReason ls_reason;
793c14b763aSAlp Dener 
794c14b763aSAlp Dener   PetscReal      e_min, gdx, delta;
795c14b763aSAlp Dener   PetscInt       bfgsUpdates;
796c14b763aSAlp Dener 
797c14b763aSAlp Dener   PetscFunctionBegin;
798c14b763aSAlp Dener   /* Perform the linesearch */
799c14b763aSAlp Dener   ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, steplen, &ls_reason);CHKERRQ(ierr);
800c14b763aSAlp Dener   ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
801c14b763aSAlp Dener 
802937a31a1SAlp Dener   while (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER && *stepType != BNK_GRADIENT) {
803c14b763aSAlp Dener     /* Linesearch failed, revert solution */
804c14b763aSAlp Dener     bnk->f = bnk->fold;
805c14b763aSAlp Dener     ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr);
806c14b763aSAlp Dener     ierr = VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient);CHKERRQ(ierr);
807c14b763aSAlp Dener 
808937a31a1SAlp Dener     switch(*stepType) {
809c14b763aSAlp Dener     case BNK_NEWTON:
8108d5ead36SAlp Dener       /* Failed to obtain acceptable iterate with Newton step
811c14b763aSAlp Dener          Update the perturbation for next time */
812c14b763aSAlp Dener       if (bnk->pert <= 0.0) {
813c14b763aSAlp Dener         /* Initialize the perturbation */
814c14b763aSAlp Dener         bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
815c14b763aSAlp Dener         if (bnk->is_gltr) {
816c14b763aSAlp Dener           ierr = KSPCGGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr);
817c14b763aSAlp Dener           bnk->pert = PetscMax(bnk->pert, -e_min);
818c14b763aSAlp Dener         }
819c14b763aSAlp Dener       } else {
820c14b763aSAlp Dener         /* Increase the perturbation */
821c14b763aSAlp Dener         bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
822c14b763aSAlp Dener       }
823c14b763aSAlp Dener 
824c14b763aSAlp Dener       if (BNK_PC_BFGS != bnk->pc_type) {
825c14b763aSAlp Dener         /* We don't have the bfgs matrix around and being updated
826c14b763aSAlp Dener            Must use gradient direction in this case */
827937a31a1SAlp Dener         ierr = VecCopy(bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
828937a31a1SAlp Dener         *stepType = BNK_GRADIENT;
829c14b763aSAlp Dener       } else {
830c14b763aSAlp Dener         /* Attempt to use the BFGS direction */
831937a31a1SAlp Dener         ierr = MatLMVMSetInactive(bnk->M, NULL);CHKERRQ(ierr);
832c14b763aSAlp Dener         ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
8338d5ead36SAlp Dener         /* Check for success (descent direction)
8348d5ead36SAlp Dener            NOTE: Negative gdx means not a descent direction because the step here is missing a negative sign. */
835c14b763aSAlp Dener         ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr);
8363105154fSTodd Munson         if ((gdx <= 0.0) || PetscIsInfOrNanReal(gdx)) {
837c14b763aSAlp Dener           /* BFGS direction is not descent or direction produced not a number
838c14b763aSAlp Dener              We can assert bfgsUpdates > 1 in this case
839c14b763aSAlp Dener              Use steepest descent direction (scaled) */
8409b6ef848SAlp Dener           delta = 2.0 * PetscMax(1.0, PetscAbsScalar(bnk->f)) / (bnk->gnorm*bnk->gnorm);
841c14b763aSAlp Dener           ierr = MatLMVMSetDelta(bnk->M, delta);CHKERRQ(ierr);
842c14b763aSAlp Dener           ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr);
843c14b763aSAlp Dener           ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
844c14b763aSAlp Dener           ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
845c14b763aSAlp Dener 
846c14b763aSAlp Dener           bfgsUpdates = 1;
847937a31a1SAlp Dener           *stepType = BNK_SCALED_GRADIENT;
848c14b763aSAlp Dener         } else {
849c14b763aSAlp Dener           ierr = MatLMVMGetUpdates(bnk->M, &bfgsUpdates);CHKERRQ(ierr);
850c14b763aSAlp Dener           if (1 == bfgsUpdates) {
851c14b763aSAlp Dener             /* The first BFGS direction is always the scaled gradient */
852937a31a1SAlp Dener             *stepType = BNK_SCALED_GRADIENT;
853c14b763aSAlp Dener           } else {
854937a31a1SAlp Dener             *stepType = BNK_BFGS;
855c14b763aSAlp Dener           }
856c14b763aSAlp Dener         }
857c14b763aSAlp Dener       }
858c14b763aSAlp Dener       break;
859c14b763aSAlp Dener 
860c14b763aSAlp Dener     case BNK_BFGS:
861c14b763aSAlp Dener       /* Can only enter if pc_type == BNK_PC_BFGS
862c14b763aSAlp Dener          Failed to obtain acceptable iterate with BFGS step
863c14b763aSAlp Dener          Attempt to use the scaled gradient direction */
8649b6ef848SAlp Dener       delta = 2.0 * PetscMax(1.0, PetscAbsScalar(bnk->f)) / (bnk->gnorm*bnk->gnorm);
865c14b763aSAlp Dener       ierr = MatLMVMSetDelta(bnk->M, delta);CHKERRQ(ierr);
866c14b763aSAlp Dener       ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr);
867c14b763aSAlp Dener       ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
868937a31a1SAlp Dener       ierr = MatLMVMSetInactive(bnk->M, NULL);CHKERRQ(ierr);
869c14b763aSAlp Dener       ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
870c14b763aSAlp Dener 
871c14b763aSAlp Dener       bfgsUpdates = 1;
872937a31a1SAlp Dener       *stepType = BNK_SCALED_GRADIENT;
873c14b763aSAlp Dener       break;
874c14b763aSAlp Dener 
875c14b763aSAlp Dener     case BNK_SCALED_GRADIENT:
876c14b763aSAlp Dener       /* Can only enter if pc_type == BNK_PC_BFGS
877c14b763aSAlp Dener          The scaled gradient step did not produce a new iterate;
878937a31a1SAlp Dener          reset the BFGS matrix and attemp to use the gradient direction. */
879c14b763aSAlp Dener       ierr = MatLMVMSetScale(bnk->M,0);CHKERRQ(ierr);
880c14b763aSAlp Dener       ierr = MatLMVMSetDelta(bnk->M,1.0);CHKERRQ(ierr);
881c14b763aSAlp Dener       ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr);
882c14b763aSAlp Dener       ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
883937a31a1SAlp Dener       ierr = VecCopy(bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
884c14b763aSAlp Dener 
885c14b763aSAlp Dener       bfgsUpdates = 1;
886937a31a1SAlp Dener       *stepType = BNK_GRADIENT;
887c14b763aSAlp Dener       break;
888c14b763aSAlp Dener     }
8898d5ead36SAlp Dener     /* Make sure the safeguarded fall-back step is zero for actively bounded variables */
8908d5ead36SAlp Dener     ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
891a1318120SAlp Dener     ierr = TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection);CHKERRQ(ierr);
892c14b763aSAlp Dener 
8938d5ead36SAlp Dener     /* Perform one last line search with the fall-back step */
894c14b763aSAlp Dener     ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, steplen, &ls_reason);CHKERRQ(ierr);
895c14b763aSAlp Dener     ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
896c14b763aSAlp Dener   }
897c14b763aSAlp Dener   *reason = ls_reason;
898c14b763aSAlp Dener   PetscFunctionReturn(0);
899c14b763aSAlp Dener }
900c14b763aSAlp Dener 
901df278d8fSAlp Dener /*------------------------------------------------------------*/
902df278d8fSAlp Dener 
903df278d8fSAlp Dener /* Routine for updating the trust radius.
904df278d8fSAlp Dener 
905df278d8fSAlp Dener   Function features three different update methods:
906df278d8fSAlp Dener   1) Line-search step length based
907df278d8fSAlp Dener   2) Predicted decrease on the CG quadratic model
908df278d8fSAlp Dener   3) Interpolation
909df278d8fSAlp Dener */
910df278d8fSAlp Dener 
91128017e9fSAlp Dener PetscErrorCode TaoBNKUpdateTrustRadius(Tao tao, PetscReal prered, PetscReal actred, PetscInt updateType, PetscInt stepType, PetscBool *accept)
912080d2917SAlp Dener {
913080d2917SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
914080d2917SAlp Dener   PetscErrorCode ierr;
915080d2917SAlp Dener 
916b1c2d0e3SAlp Dener   PetscReal      step, kappa;
917080d2917SAlp Dener   PetscReal      gdx, tau_1, tau_2, tau_min, tau_max;
918080d2917SAlp Dener 
919080d2917SAlp Dener   PetscFunctionBegin;
920080d2917SAlp Dener   /* Update trust region radius */
921080d2917SAlp Dener   *accept = PETSC_FALSE;
92228017e9fSAlp Dener   switch(updateType) {
923080d2917SAlp Dener   case BNK_UPDATE_STEP:
924c14b763aSAlp Dener     *accept = PETSC_TRUE; /* always accept here because line search succeeded */
925080d2917SAlp Dener     if (stepType == BNK_NEWTON) {
926080d2917SAlp Dener       ierr = TaoLineSearchGetStepLength(tao->linesearch, &step);CHKERRQ(ierr);
927080d2917SAlp Dener       if (step < bnk->nu1) {
928080d2917SAlp Dener         /* Very bad step taken; reduce radius */
929080d2917SAlp Dener         tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust);
930080d2917SAlp Dener       } else if (step < bnk->nu2) {
931080d2917SAlp Dener         /* Reasonably bad step taken; reduce radius */
932080d2917SAlp Dener         tao->trust = bnk->omega2 * PetscMin(bnk->dnorm, tao->trust);
933080d2917SAlp Dener       } else if (step < bnk->nu3) {
934080d2917SAlp Dener         /*  Reasonable step was taken; leave radius alone */
935080d2917SAlp Dener         if (bnk->omega3 < 1.0) {
936080d2917SAlp Dener           tao->trust = bnk->omega3 * PetscMin(bnk->dnorm, tao->trust);
937080d2917SAlp Dener         } else if (bnk->omega3 > 1.0) {
938080d2917SAlp Dener           tao->trust = PetscMax(bnk->omega3 * bnk->dnorm, tao->trust);
939080d2917SAlp Dener         }
940080d2917SAlp Dener       } else if (step < bnk->nu4) {
941080d2917SAlp Dener         /*  Full step taken; increase the radius */
942080d2917SAlp Dener         tao->trust = PetscMax(bnk->omega4 * bnk->dnorm, tao->trust);
943080d2917SAlp Dener       } else {
944080d2917SAlp Dener         /*  More than full step taken; increase the radius */
945080d2917SAlp Dener         tao->trust = PetscMax(bnk->omega5 * bnk->dnorm, tao->trust);
946080d2917SAlp Dener       }
947080d2917SAlp Dener     } else {
948080d2917SAlp Dener       /*  Newton step was not good; reduce the radius */
949080d2917SAlp Dener       tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust);
950080d2917SAlp Dener     }
951080d2917SAlp Dener     break;
952080d2917SAlp Dener 
953080d2917SAlp Dener   case BNK_UPDATE_REDUCTION:
954080d2917SAlp Dener     if (stepType == BNK_NEWTON) {
955b1c2d0e3SAlp Dener       if (prered < 0.0) {
956fed79b8eSAlp Dener         /* The predicted reduction has the wrong sign.  This cannot
957fed79b8eSAlp Dener            happen in infinite precision arithmetic.  Step should
958fed79b8eSAlp Dener            be rejected! */
959080d2917SAlp Dener         tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
9603105154fSTodd Munson       } else {
961b1c2d0e3SAlp Dener         if (PetscIsInfOrNanReal(actred)) {
962080d2917SAlp Dener           tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
963080d2917SAlp Dener         } else {
9643105154fSTodd Munson           if ((PetscAbsScalar(actred) <= PetscMax(1.0, PetscAbsScalar(bnk->f))*bnk->epsilon) && (PetscAbsScalar(prered) <= PetscMax(1.0, PetscAbsScalar(bnk->f))*bnk->epsilon)) {
965080d2917SAlp Dener             kappa = 1.0;
9663105154fSTodd Munson           } else {
967080d2917SAlp Dener             kappa = actred / prered;
968080d2917SAlp Dener           }
969fed79b8eSAlp Dener 
970fed79b8eSAlp Dener           /* Accept or reject the step and update radius */
971080d2917SAlp Dener           if (kappa < bnk->eta1) {
972fed79b8eSAlp Dener             /* Reject the step */
973080d2917SAlp Dener             tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
9743105154fSTodd Munson           } else {
975fed79b8eSAlp Dener             /* Accept the step */
976c133c014SAlp Dener             *accept = PETSC_TRUE;
977c133c014SAlp Dener             /* Update the trust region radius only if the computed step is at the trust radius boundary */
9788d5ead36SAlp Dener             if (bnk->dnorm == tao->trust) {
979080d2917SAlp Dener               if (kappa < bnk->eta2) {
980080d2917SAlp Dener                 /* Marginal bad step */
981c133c014SAlp Dener                 tao->trust = bnk->alpha2 * tao->trust;
9823105154fSTodd Munson               } else if (kappa < bnk->eta3) {
983fed79b8eSAlp Dener                 /* Reasonable step */
984fed79b8eSAlp Dener                 tao->trust = bnk->alpha3 * tao->trust;
9853105154fSTodd Munson               } else if (kappa < bnk->eta4) {
986080d2917SAlp Dener                 /* Good step */
987c133c014SAlp Dener                 tao->trust = bnk->alpha4 * tao->trust;
9883105154fSTodd Munson               } else {
989080d2917SAlp Dener                 /* Very good step */
990c133c014SAlp Dener                 tao->trust = bnk->alpha5 * tao->trust;
991080d2917SAlp Dener               }
992c133c014SAlp Dener             }
993080d2917SAlp Dener           }
994080d2917SAlp Dener         }
995080d2917SAlp Dener       }
996080d2917SAlp Dener     } else {
997080d2917SAlp Dener       /*  Newton step was not good; reduce the radius */
998080d2917SAlp Dener       tao->trust = bnk->alpha1 * PetscMin(bnk->dnorm, tao->trust);
999080d2917SAlp Dener     }
1000080d2917SAlp Dener     break;
1001080d2917SAlp Dener 
1002080d2917SAlp Dener   default:
1003080d2917SAlp Dener     if (stepType == BNK_NEWTON) {
1004b1c2d0e3SAlp Dener       if (prered < 0.0) {
1005080d2917SAlp Dener         /*  The predicted reduction has the wrong sign.  This cannot */
1006080d2917SAlp Dener         /*  happen in infinite precision arithmetic.  Step should */
1007080d2917SAlp Dener         /*  be rejected! */
1008080d2917SAlp Dener         tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
1009080d2917SAlp Dener       } else {
1010b1c2d0e3SAlp Dener         if (PetscIsInfOrNanReal(actred)) {
1011080d2917SAlp Dener           tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
1012080d2917SAlp Dener         } else {
1013080d2917SAlp Dener           if ((PetscAbsScalar(actred) <= bnk->epsilon) && (PetscAbsScalar(prered) <= bnk->epsilon)) {
1014080d2917SAlp Dener             kappa = 1.0;
1015080d2917SAlp Dener           } else {
1016080d2917SAlp Dener             kappa = actred / prered;
1017080d2917SAlp Dener           }
1018080d2917SAlp Dener 
1019080d2917SAlp Dener           ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr);
1020080d2917SAlp Dener           tau_1 = bnk->theta * gdx / (bnk->theta * gdx - (1.0 - bnk->theta) * prered + actred);
1021080d2917SAlp Dener           tau_2 = bnk->theta * gdx / (bnk->theta * gdx + (1.0 + bnk->theta) * prered - actred);
1022080d2917SAlp Dener           tau_min = PetscMin(tau_1, tau_2);
1023080d2917SAlp Dener           tau_max = PetscMax(tau_1, tau_2);
1024080d2917SAlp Dener 
1025080d2917SAlp Dener           if (kappa >= 1.0 - bnk->mu1) {
1026080d2917SAlp Dener             /*  Great agreement */
1027080d2917SAlp Dener             *accept = PETSC_TRUE;
1028080d2917SAlp Dener             if (tau_max < 1.0) {
1029080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm);
1030080d2917SAlp Dener             } else if (tau_max > bnk->gamma4) {
1031080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, bnk->gamma4 * bnk->dnorm);
1032080d2917SAlp Dener             } else {
1033080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm);
1034080d2917SAlp Dener             }
1035080d2917SAlp Dener           } else if (kappa >= 1.0 - bnk->mu2) {
1036080d2917SAlp Dener             /*  Good agreement */
1037080d2917SAlp Dener             *accept = PETSC_TRUE;
1038080d2917SAlp Dener             if (tau_max < bnk->gamma2) {
1039080d2917SAlp Dener               tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm);
1040080d2917SAlp Dener             } else if (tau_max > bnk->gamma3) {
1041080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm);
1042080d2917SAlp Dener             } else if (tau_max < 1.0) {
1043080d2917SAlp Dener               tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm);
1044080d2917SAlp Dener             } else {
1045080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm);
1046080d2917SAlp Dener             }
1047080d2917SAlp Dener           } else {
1048080d2917SAlp Dener             /*  Not good agreement */
1049080d2917SAlp Dener             if (tau_min > 1.0) {
1050080d2917SAlp Dener               tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm);
1051080d2917SAlp Dener             } else if (tau_max < bnk->gamma1) {
1052080d2917SAlp Dener               tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
1053080d2917SAlp Dener             } else if ((tau_min < bnk->gamma1) && (tau_max >= 1.0)) {
1054080d2917SAlp Dener               tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
1055080d2917SAlp Dener             } else if ((tau_1 >= bnk->gamma1) && (tau_1 < 1.0) && ((tau_2 < bnk->gamma1) || (tau_2 >= 1.0))) {
1056080d2917SAlp Dener               tao->trust = tau_1 * PetscMin(tao->trust, bnk->dnorm);
1057080d2917SAlp Dener             } else if ((tau_2 >= bnk->gamma1) && (tau_2 < 1.0) && ((tau_1 < bnk->gamma1) || (tau_2 >= 1.0))) {
1058080d2917SAlp Dener               tao->trust = tau_2 * PetscMin(tao->trust, bnk->dnorm);
1059080d2917SAlp Dener             } else {
1060080d2917SAlp Dener               tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm);
1061080d2917SAlp Dener             }
1062080d2917SAlp Dener           }
1063080d2917SAlp Dener         }
1064080d2917SAlp Dener       }
1065080d2917SAlp Dener     } else {
1066080d2917SAlp Dener       /*  Newton step was not good; reduce the radius */
1067080d2917SAlp Dener       tao->trust = bnk->gamma1 * PetscMin(bnk->dnorm, tao->trust);
1068080d2917SAlp Dener     }
106928017e9fSAlp Dener     break;
1070080d2917SAlp Dener   }
1071c133c014SAlp Dener   /* Make sure the radius does not violate min and max settings */
1072c133c014SAlp Dener   tao->trust = PetscMin(tao->trust, bnk->max_radius);
1073fed79b8eSAlp Dener   tao->trust = PetscMax(tao->trust, bnk->min_radius);
1074080d2917SAlp Dener   PetscFunctionReturn(0);
1075080d2917SAlp Dener }
1076080d2917SAlp Dener 
1077eb910715SAlp Dener /* ---------------------------------------------------------- */
1078df278d8fSAlp Dener 
107962675beeSAlp Dener PetscErrorCode TaoBNKAddStepCounts(Tao tao, PetscInt stepType)
108062675beeSAlp Dener {
108162675beeSAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
108262675beeSAlp Dener 
108362675beeSAlp Dener   PetscFunctionBegin;
108462675beeSAlp Dener   switch (stepType) {
108562675beeSAlp Dener   case BNK_NEWTON:
108662675beeSAlp Dener     ++bnk->newt;
108762675beeSAlp Dener     break;
108862675beeSAlp Dener   case BNK_BFGS:
108962675beeSAlp Dener     ++bnk->bfgs;
109062675beeSAlp Dener     break;
109162675beeSAlp Dener   case BNK_SCALED_GRADIENT:
109262675beeSAlp Dener     ++bnk->sgrad;
109362675beeSAlp Dener     break;
109462675beeSAlp Dener   case BNK_GRADIENT:
109562675beeSAlp Dener     ++bnk->grad;
109662675beeSAlp Dener     break;
109762675beeSAlp Dener   default:
109862675beeSAlp Dener     break;
109962675beeSAlp Dener   }
110062675beeSAlp Dener   PetscFunctionReturn(0);
110162675beeSAlp Dener }
110262675beeSAlp Dener 
110362675beeSAlp Dener /* ---------------------------------------------------------- */
110462675beeSAlp Dener 
11059b6ef848SAlp Dener PetscErrorCode TaoSetUp_BNK(Tao tao)
1106eb910715SAlp Dener {
1107eb910715SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
1108eb910715SAlp Dener   PetscErrorCode ierr;
11099b6ef848SAlp Dener   KSPType        ksp_type;
1110e031d6f5SAlp Dener   PetscInt       i;
1111eb910715SAlp Dener 
1112eb910715SAlp Dener   PetscFunctionBegin;
1113c4b75bccSAlp Dener   if (!tao->gradient) {
1114c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);
1115c4b75bccSAlp Dener   }
1116c4b75bccSAlp Dener   if (!tao->stepdirection) {
1117c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr);
1118c4b75bccSAlp Dener   }
1119c4b75bccSAlp Dener   if (!bnk->W) {
1120c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->W);CHKERRQ(ierr);
1121c4b75bccSAlp Dener   }
1122c4b75bccSAlp Dener   if (!bnk->Xold) {
1123c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->Xold);CHKERRQ(ierr);
1124c4b75bccSAlp Dener   }
1125c4b75bccSAlp Dener   if (!bnk->Gold) {
1126c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->Gold);CHKERRQ(ierr);
1127c4b75bccSAlp Dener   }
1128c4b75bccSAlp Dener   if (!bnk->Xwork) {
1129c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->Xwork);CHKERRQ(ierr);
1130c4b75bccSAlp Dener   }
1131c4b75bccSAlp Dener   if (!bnk->Gwork) {
1132c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->Gwork);CHKERRQ(ierr);
1133c4b75bccSAlp Dener   }
1134c4b75bccSAlp Dener   if (!bnk->unprojected_gradient) {
1135c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->unprojected_gradient);CHKERRQ(ierr);
1136c4b75bccSAlp Dener   }
1137c4b75bccSAlp Dener   if (!bnk->unprojected_gradient_old) {
1138c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->unprojected_gradient_old);CHKERRQ(ierr);
1139c4b75bccSAlp Dener   }
1140c4b75bccSAlp Dener   if (!bnk->Diag_min) {
1141c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->Diag_min);CHKERRQ(ierr);
1142c4b75bccSAlp Dener   }
1143c4b75bccSAlp Dener   if (!bnk->Diag_max) {
1144c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->Diag_max);CHKERRQ(ierr);
1145c4b75bccSAlp Dener   }
1146e031d6f5SAlp Dener   if (bnk->max_cg_its > 0) {
1147c4b75bccSAlp Dener     /* Ensure that the important common vectors are shared between BNK and embedded BNCG */
1148c4b75bccSAlp Dener     bnk->bncg_ctx = (TAO_BNCG *)bnk->bncg->data;
114989da521bSAlp Dener     ierr = PetscObjectReference((PetscObject)(bnk->unprojected_gradient_old));CHKERRQ(ierr);
115089da521bSAlp Dener     ierr = VecDestroy(&bnk->bncg_ctx->unprojected_gradient_old);CHKERRQ(ierr);
115189da521bSAlp Dener     bnk->bncg_ctx->unprojected_gradient_old = bnk->unprojected_gradient_old;
1152c4b75bccSAlp Dener     ierr = PetscObjectReference((PetscObject)(bnk->unprojected_gradient));CHKERRQ(ierr);
1153c4b75bccSAlp Dener     ierr = VecDestroy(&bnk->bncg_ctx->unprojected_gradient);CHKERRQ(ierr);
1154c4b75bccSAlp Dener     bnk->bncg_ctx->unprojected_gradient = bnk->unprojected_gradient;
1155c4b75bccSAlp Dener     ierr = PetscObjectReference((PetscObject)(bnk->Gold));CHKERRQ(ierr);
1156c4b75bccSAlp Dener     ierr = VecDestroy(&bnk->bncg_ctx->G_old);CHKERRQ(ierr);
1157c4b75bccSAlp Dener     bnk->bncg_ctx->G_old = bnk->Gold;
1158c4b75bccSAlp Dener     ierr = PetscObjectReference((PetscObject)(tao->gradient));CHKERRQ(ierr);
1159c4b75bccSAlp Dener     ierr = VecDestroy(&bnk->bncg->gradient);CHKERRQ(ierr);
1160c4b75bccSAlp Dener     bnk->bncg->gradient = tao->gradient;
1161c4b75bccSAlp Dener     ierr = PetscObjectReference((PetscObject)(tao->stepdirection));CHKERRQ(ierr);
1162c4b75bccSAlp Dener     ierr = VecDestroy(&bnk->bncg->stepdirection);CHKERRQ(ierr);
1163c4b75bccSAlp Dener     bnk->bncg->stepdirection = tao->stepdirection;
1164c4b75bccSAlp Dener     ierr = TaoSetInitialVector(bnk->bncg, tao->solution);CHKERRQ(ierr);
1165c4b75bccSAlp Dener     /* Copy over some settings from BNK into BNCG */
1166e031d6f5SAlp Dener     ierr = TaoSetMaximumIterations(bnk->bncg, bnk->max_cg_its);CHKERRQ(ierr);
1167e031d6f5SAlp Dener     ierr = TaoSetTolerances(bnk->bncg, tao->gatol, tao->grtol, tao->gttol);CHKERRQ(ierr);
1168e031d6f5SAlp Dener     ierr = TaoSetFunctionLowerBound(bnk->bncg, tao->fmin);CHKERRQ(ierr);
1169937a31a1SAlp Dener     ierr = TaoSetConvergenceTest(bnk->bncg, tao->ops->convergencetest, tao->cnvP);CHKERRQ(ierr);
1170e031d6f5SAlp Dener     ierr = TaoSetObjectiveRoutine(bnk->bncg, tao->ops->computeobjective, tao->user_objP);CHKERRQ(ierr);
1171e031d6f5SAlp Dener     ierr = TaoSetGradientRoutine(bnk->bncg, tao->ops->computegradient, tao->user_gradP);CHKERRQ(ierr);
1172e031d6f5SAlp Dener     ierr = TaoSetObjectiveAndGradientRoutine(bnk->bncg, tao->ops->computeobjectiveandgradient, tao->user_objgradP);CHKERRQ(ierr);
1173e031d6f5SAlp Dener     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)tao, (PetscObject)(bnk->bncg));CHKERRQ(ierr);
1174c4b75bccSAlp Dener     for (i=0; i<tao->numbermonitors; ++i) {
1175e031d6f5SAlp Dener       ierr = TaoSetMonitor(bnk->bncg, tao->monitor[i], tao->monitorcontext[i], tao->monitordestroy[i]);CHKERRQ(ierr);
1176e031d6f5SAlp Dener       ierr = PetscObjectReference((PetscObject)(tao->monitorcontext[i]));CHKERRQ(ierr);
1177e031d6f5SAlp Dener     }
1178e031d6f5SAlp Dener   }
1179eb910715SAlp Dener   bnk->Diag = 0;
1180c0f10754SAlp Dener   bnk->Diag_red = 0;
1181c0f10754SAlp Dener   bnk->X_inactive = 0;
1182c0f10754SAlp Dener   bnk->G_inactive = 0;
118362675beeSAlp Dener   bnk->inactive_work = 0;
118462675beeSAlp Dener   bnk->active_work = 0;
118562675beeSAlp Dener   bnk->inactive_idx = 0;
118662675beeSAlp Dener   bnk->active_idx = 0;
118762675beeSAlp Dener   bnk->active_lower = 0;
118862675beeSAlp Dener   bnk->active_upper = 0;
118962675beeSAlp Dener   bnk->active_fixed = 0;
1190eb910715SAlp Dener   bnk->M = 0;
119109164190SAlp Dener   bnk->H_inactive = 0;
119209164190SAlp Dener   bnk->Hpre_inactive = 0;
11939b6ef848SAlp Dener   ierr = KSPGetType(tao->ksp,&ksp_type);CHKERRQ(ierr);
11949b6ef848SAlp Dener   ierr = PetscStrcmp(ksp_type,KSPCGNASH,&bnk->is_nash);CHKERRQ(ierr);
11959b6ef848SAlp Dener   ierr = PetscStrcmp(ksp_type,KSPCGSTCG,&bnk->is_stcg);CHKERRQ(ierr);
11969b6ef848SAlp Dener   ierr = PetscStrcmp(ksp_type,KSPCGGLTR,&bnk->is_gltr);CHKERRQ(ierr);
1197eb910715SAlp Dener   PetscFunctionReturn(0);
1198eb910715SAlp Dener }
1199eb910715SAlp Dener 
1200eb910715SAlp Dener /*------------------------------------------------------------*/
1201df278d8fSAlp Dener 
1202eb910715SAlp Dener static PetscErrorCode TaoDestroy_BNK(Tao tao)
1203eb910715SAlp Dener {
1204eb910715SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
1205eb910715SAlp Dener   PetscErrorCode ierr;
1206eb910715SAlp Dener 
1207eb910715SAlp Dener   PetscFunctionBegin;
1208eb910715SAlp Dener   if (tao->setupcalled) {
1209eb910715SAlp Dener     ierr = VecDestroy(&bnk->W);CHKERRQ(ierr);
1210eb910715SAlp Dener     ierr = VecDestroy(&bnk->Xold);CHKERRQ(ierr);
1211eb910715SAlp Dener     ierr = VecDestroy(&bnk->Gold);CHKERRQ(ierr);
121209164190SAlp Dener     ierr = VecDestroy(&bnk->Xwork);CHKERRQ(ierr);
121309164190SAlp Dener     ierr = VecDestroy(&bnk->Gwork);CHKERRQ(ierr);
121409164190SAlp Dener     ierr = VecDestroy(&bnk->unprojected_gradient);CHKERRQ(ierr);
121509164190SAlp Dener     ierr = VecDestroy(&bnk->unprojected_gradient_old);CHKERRQ(ierr);
121662675beeSAlp Dener     ierr = VecDestroy(&bnk->Diag_min);CHKERRQ(ierr);
121762675beeSAlp Dener     ierr = VecDestroy(&bnk->Diag_max);CHKERRQ(ierr);
1218c4b75bccSAlp Dener   }
1219*ca964c49SAlp Dener   ierr = ISDestroy(&bnk->active_lower);CHKERRQ(ierr);
1220*ca964c49SAlp Dener   ierr = ISDestroy(&bnk->active_upper);CHKERRQ(ierr);
1221*ca964c49SAlp Dener   ierr = ISDestroy(&bnk->active_fixed);CHKERRQ(ierr);
1222*ca964c49SAlp Dener   ierr = ISDestroy(&bnk->active_idx);CHKERRQ(ierr);
1223*ca964c49SAlp Dener   ierr = ISDestroy(&bnk->inactive_idx);CHKERRQ(ierr);
12245e9b73cbSAlp Dener   ierr = VecDestroy(&bnk->Diag);CHKERRQ(ierr);
1225eb910715SAlp Dener   ierr = MatDestroy(&bnk->M);CHKERRQ(ierr);
1226c4b75bccSAlp Dener   if (bnk->Hpre_inactive != tao->hessian_pre && bnk->Hpre_inactive != bnk->H_inactive) {
1227c4b75bccSAlp Dener     ierr = MatDestroy(&bnk->Hpre_inactive);CHKERRQ(ierr);
1228c4b75bccSAlp Dener   }
1229c4b75bccSAlp Dener   if (bnk->H_inactive != tao->hessian) {
1230c4b75bccSAlp Dener     ierr = MatDestroy(&bnk->H_inactive);CHKERRQ(ierr);
1231c4b75bccSAlp Dener   }
1232*ca964c49SAlp Dener   if (bnk->max_cg_its > 0) {
1233*ca964c49SAlp Dener     ierr = TaoDestroy(&bnk->bncg);CHKERRQ(ierr);
1234*ca964c49SAlp Dener   }
1235eb910715SAlp Dener   ierr = PetscFree(tao->data);CHKERRQ(ierr);
1236eb910715SAlp Dener   PetscFunctionReturn(0);
1237eb910715SAlp Dener }
1238eb910715SAlp Dener 
1239eb910715SAlp Dener /*------------------------------------------------------------*/
1240df278d8fSAlp Dener 
1241eb910715SAlp Dener static PetscErrorCode TaoSetFromOptions_BNK(PetscOptionItems *PetscOptionsObject,Tao tao)
1242eb910715SAlp Dener {
1243eb910715SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
1244eb910715SAlp Dener   PetscErrorCode ierr;
1245eb910715SAlp Dener 
1246eb910715SAlp Dener   PetscFunctionBegin;
1247eb910715SAlp Dener   ierr = PetscOptionsHead(PetscOptionsObject,"Newton line search method for unconstrained optimization");CHKERRQ(ierr);
12488d5ead36SAlp Dener   ierr = PetscOptionsEList("-tao_bnk_pc_type", "pc type", "", BNK_PC, BNK_PC_TYPES, BNK_PC[bnk->pc_type], &bnk->pc_type, 0);CHKERRQ(ierr);
12498d5ead36SAlp Dener   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);
12508d5ead36SAlp Dener   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);
12518d5ead36SAlp Dener   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);
12522f75a4aaSAlp Dener   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);
12538d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_sval", "perturbation starting value", "", bnk->sval, &bnk->sval,NULL);CHKERRQ(ierr);
12548d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_imin", "minimum initial perturbation", "", bnk->imin, &bnk->imin,NULL);CHKERRQ(ierr);
12558d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_imax", "maximum initial perturbation", "", bnk->imax, &bnk->imax,NULL);CHKERRQ(ierr);
12568d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_imfac", "initial merit factor", "", bnk->imfac, &bnk->imfac,NULL);CHKERRQ(ierr);
12578d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_pmin", "minimum perturbation", "", bnk->pmin, &bnk->pmin,NULL);CHKERRQ(ierr);
12588d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_pmax", "maximum perturbation", "", bnk->pmax, &bnk->pmax,NULL);CHKERRQ(ierr);
12598d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_pgfac", "growth factor", "", bnk->pgfac, &bnk->pgfac,NULL);CHKERRQ(ierr);
12608d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_psfac", "shrink factor", "", bnk->psfac, &bnk->psfac,NULL);CHKERRQ(ierr);
12618d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_pmgfac", "merit growth factor", "", bnk->pmgfac, &bnk->pmgfac,NULL);CHKERRQ(ierr);
12628d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_pmsfac", "merit shrink factor", "", bnk->pmsfac, &bnk->pmsfac,NULL);CHKERRQ(ierr);
12638d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_eta1", "poor steplength; reduce radius", "", bnk->eta1, &bnk->eta1,NULL);CHKERRQ(ierr);
12648d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_eta2", "reasonable steplength; leave radius alone", "", bnk->eta2, &bnk->eta2,NULL);CHKERRQ(ierr);
12658d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_eta3", "good steplength; increase radius", "", bnk->eta3, &bnk->eta3,NULL);CHKERRQ(ierr);
12668d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_eta4", "excellent steplength; greatly increase radius", "", bnk->eta4, &bnk->eta4,NULL);CHKERRQ(ierr);
12678d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_alpha1", "", "", bnk->alpha1, &bnk->alpha1,NULL);CHKERRQ(ierr);
12688d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_alpha2", "", "", bnk->alpha2, &bnk->alpha2,NULL);CHKERRQ(ierr);
12698d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_alpha3", "", "", bnk->alpha3, &bnk->alpha3,NULL);CHKERRQ(ierr);
12708d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_alpha4", "", "", bnk->alpha4, &bnk->alpha4,NULL);CHKERRQ(ierr);
12718d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_alpha5", "", "", bnk->alpha5, &bnk->alpha5,NULL);CHKERRQ(ierr);
12728d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_nu1", "poor steplength; reduce radius", "", bnk->nu1, &bnk->nu1,NULL);CHKERRQ(ierr);
12738d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_nu2", "reasonable steplength; leave radius alone", "", bnk->nu2, &bnk->nu2,NULL);CHKERRQ(ierr);
12748d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_nu3", "good steplength; increase radius", "", bnk->nu3, &bnk->nu3,NULL);CHKERRQ(ierr);
12758d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_nu4", "excellent steplength; greatly increase radius", "", bnk->nu4, &bnk->nu4,NULL);CHKERRQ(ierr);
12768d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_omega1", "", "", bnk->omega1, &bnk->omega1,NULL);CHKERRQ(ierr);
12778d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_omega2", "", "", bnk->omega2, &bnk->omega2,NULL);CHKERRQ(ierr);
12788d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_omega3", "", "", bnk->omega3, &bnk->omega3,NULL);CHKERRQ(ierr);
12798d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_omega4", "", "", bnk->omega4, &bnk->omega4,NULL);CHKERRQ(ierr);
12808d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_omega5", "", "", bnk->omega5, &bnk->omega5,NULL);CHKERRQ(ierr);
12818d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_mu1_i", "", "", bnk->mu1_i, &bnk->mu1_i,NULL);CHKERRQ(ierr);
12828d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_mu2_i", "", "", bnk->mu2_i, &bnk->mu2_i,NULL);CHKERRQ(ierr);
12838d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_gamma1_i", "", "", bnk->gamma1_i, &bnk->gamma1_i,NULL);CHKERRQ(ierr);
12848d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_gamma2_i", "", "", bnk->gamma2_i, &bnk->gamma2_i,NULL);CHKERRQ(ierr);
12858d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_gamma3_i", "", "", bnk->gamma3_i, &bnk->gamma3_i,NULL);CHKERRQ(ierr);
12868d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_gamma4_i", "", "", bnk->gamma4_i, &bnk->gamma4_i,NULL);CHKERRQ(ierr);
12878d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_theta_i", "", "", bnk->theta_i, &bnk->theta_i,NULL);CHKERRQ(ierr);
12888d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_mu1", "", "", bnk->mu1, &bnk->mu1,NULL);CHKERRQ(ierr);
12898d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_mu2", "", "", bnk->mu2, &bnk->mu2,NULL);CHKERRQ(ierr);
12908d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_gamma1", "", "", bnk->gamma1, &bnk->gamma1,NULL);CHKERRQ(ierr);
12918d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_gamma2", "", "", bnk->gamma2, &bnk->gamma2,NULL);CHKERRQ(ierr);
12928d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_gamma3", "", "", bnk->gamma3, &bnk->gamma3,NULL);CHKERRQ(ierr);
12938d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_gamma4", "", "", bnk->gamma4, &bnk->gamma4,NULL);CHKERRQ(ierr);
12948d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_theta", "", "", bnk->theta, &bnk->theta,NULL);CHKERRQ(ierr);
12958d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_min_radius", "lower bound on initial radius", "", bnk->min_radius, &bnk->min_radius,NULL);CHKERRQ(ierr);
12968d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_max_radius", "upper bound on radius", "", bnk->max_radius, &bnk->max_radius,NULL);CHKERRQ(ierr);
12978d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_epsilon", "tolerance used when computing actual and predicted reduction", "", bnk->epsilon, &bnk->epsilon,NULL);CHKERRQ(ierr);
12980a4511e9SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_as_tol", "initial tolerance used when estimating actively bounded variables", "", bnk->as_tol, &bnk->as_tol,NULL);CHKERRQ(ierr);
12990a4511e9SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_as_step", "step length used when estimating actively bounded variables", "", bnk->as_step, &bnk->as_step,NULL);CHKERRQ(ierr);
1300c0f10754SAlp Dener   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);
1301eb910715SAlp Dener   ierr = PetscOptionsTail();CHKERRQ(ierr);
1302e031d6f5SAlp Dener   ierr = TaoSetFromOptions(bnk->bncg);
1303eb910715SAlp Dener   ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr);
1304eb910715SAlp Dener   ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
1305eb910715SAlp Dener   PetscFunctionReturn(0);
1306eb910715SAlp Dener }
1307eb910715SAlp Dener 
1308eb910715SAlp Dener /*------------------------------------------------------------*/
1309df278d8fSAlp Dener 
1310eb910715SAlp Dener static PetscErrorCode TaoView_BNK(Tao tao, PetscViewer viewer)
1311eb910715SAlp Dener {
1312eb910715SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
1313eb910715SAlp Dener   PetscInt       nrejects;
1314eb910715SAlp Dener   PetscBool      isascii;
1315eb910715SAlp Dener   PetscErrorCode ierr;
1316eb910715SAlp Dener 
1317eb910715SAlp Dener   PetscFunctionBegin;
1318eb910715SAlp Dener   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
1319eb910715SAlp Dener   if (isascii) {
1320eb910715SAlp Dener     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1321eb910715SAlp Dener     if (BNK_PC_BFGS == bnk->pc_type && bnk->M) {
1322eb910715SAlp Dener       ierr = MatLMVMGetRejects(bnk->M,&nrejects);CHKERRQ(ierr);
1323eb910715SAlp Dener       ierr = PetscViewerASCIIPrintf(viewer, "Rejected matrix updates: %D\n",nrejects);CHKERRQ(ierr);
1324eb910715SAlp Dener     }
1325e031d6f5SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "CG steps: %D\n", bnk->tot_cg_its);CHKERRQ(ierr);
1326eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "Newton steps: %D\n", bnk->newt);CHKERRQ(ierr);
1327eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "BFGS steps: %D\n", bnk->bfgs);CHKERRQ(ierr);
1328eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "Scaled gradient steps: %D\n", bnk->sgrad);CHKERRQ(ierr);
1329eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", bnk->grad);CHKERRQ(ierr);
1330eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "KSP termination reasons:\n");CHKERRQ(ierr);
1331eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  atol: %D\n", bnk->ksp_atol);CHKERRQ(ierr);
1332eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  rtol: %D\n", bnk->ksp_rtol);CHKERRQ(ierr);
1333eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  ctol: %D\n", bnk->ksp_ctol);CHKERRQ(ierr);
1334eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  negc: %D\n", bnk->ksp_negc);CHKERRQ(ierr);
1335eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  dtol: %D\n", bnk->ksp_dtol);CHKERRQ(ierr);
1336eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  iter: %D\n", bnk->ksp_iter);CHKERRQ(ierr);
1337eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  othr: %D\n", bnk->ksp_othr);CHKERRQ(ierr);
1338eb910715SAlp Dener     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1339eb910715SAlp Dener   }
1340eb910715SAlp Dener   PetscFunctionReturn(0);
1341eb910715SAlp Dener }
1342eb910715SAlp Dener 
1343eb910715SAlp Dener /* ---------------------------------------------------------- */
1344df278d8fSAlp Dener 
1345eb910715SAlp Dener /*MC
1346eb910715SAlp Dener   TAOBNK - Shared base-type for Bounded Newton-Krylov type algorithms.
134766ed3702SAlp Dener   At each iteration, the BNK methods solve the symmetric
1348eb910715SAlp Dener   system of equations to obtain the step diretion dk:
1349eb910715SAlp Dener               Hk dk = -gk
13502b97c8d8SAlp Dener   for free variables only. The step can be globalized either through
13512b97c8d8SAlp Dener   trust-region methods, or a line search, or a heuristic mixture of both.
1352eb910715SAlp Dener 
1353eb910715SAlp Dener     Options Database Keys:
13548d5ead36SAlp Dener + -tao_bnk_pc_type - "none","ahess","bfgs","petsc"
13558d5ead36SAlp Dener . -tao_bnk_bfgs_scale_type - "ahess","phess","bfgs"
13568d5ead36SAlp Dener . -tao_bnk_init_type - "constant","direction","interpolation"
13578d5ead36SAlp Dener . -tao_bnk_update_type - "step","direction","interpolation"
13582f75a4aaSAlp Dener . -tao_bnk_as_type - "none","bertsekas"
13598d5ead36SAlp Dener . -tao_bnk_sval - perturbation starting value
13608d5ead36SAlp Dener . -tao_bnk_imin - minimum initial perturbation
13618d5ead36SAlp Dener . -tao_bnk_imax - maximum initial perturbation
13628d5ead36SAlp Dener . -tao_bnk_pmin - minimum perturbation
13638d5ead36SAlp Dener . -tao_bnk_pmax - maximum perturbation
13648d5ead36SAlp Dener . -tao_bnk_pgfac - growth factor
13658d5ead36SAlp Dener . -tao_bnk_psfac - shrink factor
13668d5ead36SAlp Dener . -tao_bnk_imfac - initial merit factor
13678d5ead36SAlp Dener . -tao_bnk_pmgfac - merit growth factor
13688d5ead36SAlp Dener . -tao_bnk_pmsfac - merit shrink factor
13698d5ead36SAlp Dener . -tao_bnk_eta1 - poor steplength; reduce radius
13708d5ead36SAlp Dener . -tao_bnk_eta2 - reasonable steplength; leave radius
13718d5ead36SAlp Dener . -tao_bnk_eta3 - good steplength; increase readius
13728d5ead36SAlp Dener . -tao_bnk_eta4 - excellent steplength; greatly increase radius
13738d5ead36SAlp Dener . -tao_bnk_alpha1 - alpha1 reduction
13748d5ead36SAlp Dener . -tao_bnk_alpha2 - alpha2 reduction
13758d5ead36SAlp Dener . -tao_bnk_alpha3 - alpha3 reduction
13768d5ead36SAlp Dener . -tao_bnk_alpha4 - alpha4 reduction
13778d5ead36SAlp Dener . -tao_bnk_alpha - alpha5 reduction
13788d5ead36SAlp Dener . -tao_bnk_mu1 - mu1 interpolation update
13798d5ead36SAlp Dener . -tao_bnk_mu2 - mu2 interpolation update
13808d5ead36SAlp Dener . -tao_bnk_gamma1 - gamma1 interpolation update
13818d5ead36SAlp Dener . -tao_bnk_gamma2 - gamma2 interpolation update
13828d5ead36SAlp Dener . -tao_bnk_gamma3 - gamma3 interpolation update
13838d5ead36SAlp Dener . -tao_bnk_gamma4 - gamma4 interpolation update
13848d5ead36SAlp Dener . -tao_bnk_theta - theta interpolation update
13858d5ead36SAlp Dener . -tao_bnk_omega1 - omega1 step update
13868d5ead36SAlp Dener . -tao_bnk_omega2 - omega2 step update
13878d5ead36SAlp Dener . -tao_bnk_omega3 - omega3 step update
13888d5ead36SAlp Dener . -tao_bnk_omega4 - omega4 step update
13898d5ead36SAlp Dener . -tao_bnk_omega5 - omega5 step update
13908d5ead36SAlp Dener . -tao_bnk_mu1_i -  mu1 interpolation init factor
13918d5ead36SAlp Dener . -tao_bnk_mu2_i -  mu2 interpolation init factor
13928d5ead36SAlp Dener . -tao_bnk_gamma1_i -  gamma1 interpolation init factor
13938d5ead36SAlp Dener . -tao_bnk_gamma2_i -  gamma2 interpolation init factor
13948d5ead36SAlp Dener . -tao_bnk_gamma3_i -  gamma3 interpolation init factor
13958d5ead36SAlp Dener . -tao_bnk_gamma4_i -  gamma4 interpolation init factor
13962f75a4aaSAlp Dener . -tao_bnk_theta_i -  theta interpolation init factor
13972f75a4aaSAlp Dener - -tao_bnk_bound_tol -  initial tolerance used in estimating bounded active variables
1398eb910715SAlp Dener 
1399eb910715SAlp Dener   Level: beginner
1400eb910715SAlp Dener M*/
1401eb910715SAlp Dener 
1402eb910715SAlp Dener PetscErrorCode TaoCreate_BNK(Tao tao)
1403eb910715SAlp Dener {
1404eb910715SAlp Dener   TAO_BNK        *bnk;
1405eb910715SAlp Dener   const char     *morethuente_type = TAOLINESEARCHMT;
1406eb910715SAlp Dener   PetscErrorCode ierr;
1407eb910715SAlp Dener 
1408eb910715SAlp Dener   PetscFunctionBegin;
1409eb910715SAlp Dener   ierr = PetscNewLog(tao,&bnk);CHKERRQ(ierr);
1410eb910715SAlp Dener 
1411eb910715SAlp Dener   tao->ops->setup = TaoSetUp_BNK;
1412eb910715SAlp Dener   tao->ops->view = TaoView_BNK;
1413eb910715SAlp Dener   tao->ops->setfromoptions = TaoSetFromOptions_BNK;
1414eb910715SAlp Dener   tao->ops->destroy = TaoDestroy_BNK;
1415eb910715SAlp Dener 
1416eb910715SAlp Dener   /*  Override default settings (unless already changed) */
1417eb910715SAlp Dener   if (!tao->max_it_changed) tao->max_it = 50;
1418eb910715SAlp Dener   if (!tao->trust0_changed) tao->trust0 = 100.0;
1419eb910715SAlp Dener 
1420eb910715SAlp Dener   tao->data = (void*)bnk;
1421eb910715SAlp Dener 
142266ed3702SAlp Dener   /*  Hessian shifting parameters */
1423eb910715SAlp Dener   bnk->sval   = 0.0;
1424eb910715SAlp Dener   bnk->imin   = 1.0e-4;
1425eb910715SAlp Dener   bnk->imax   = 1.0e+2;
1426eb910715SAlp Dener   bnk->imfac  = 1.0e-1;
1427eb910715SAlp Dener 
1428eb910715SAlp Dener   bnk->pmin   = 1.0e-12;
1429eb910715SAlp Dener   bnk->pmax   = 1.0e+2;
1430eb910715SAlp Dener   bnk->pgfac  = 1.0e+1;
1431eb910715SAlp Dener   bnk->psfac  = 4.0e-1;
1432eb910715SAlp Dener   bnk->pmgfac = 1.0e-1;
1433eb910715SAlp Dener   bnk->pmsfac = 1.0e-1;
1434eb910715SAlp Dener 
1435eb910715SAlp Dener   /*  Default values for trust-region radius update based on steplength */
1436eb910715SAlp Dener   bnk->nu1 = 0.25;
1437eb910715SAlp Dener   bnk->nu2 = 0.50;
1438eb910715SAlp Dener   bnk->nu3 = 1.00;
1439eb910715SAlp Dener   bnk->nu4 = 1.25;
1440eb910715SAlp Dener 
1441eb910715SAlp Dener   bnk->omega1 = 0.25;
1442eb910715SAlp Dener   bnk->omega2 = 0.50;
1443eb910715SAlp Dener   bnk->omega3 = 1.00;
1444eb910715SAlp Dener   bnk->omega4 = 2.00;
1445eb910715SAlp Dener   bnk->omega5 = 4.00;
1446eb910715SAlp Dener 
1447eb910715SAlp Dener   /*  Default values for trust-region radius update based on reduction */
1448eb910715SAlp Dener   bnk->eta1 = 1.0e-4;
1449eb910715SAlp Dener   bnk->eta2 = 0.25;
1450eb910715SAlp Dener   bnk->eta3 = 0.50;
1451eb910715SAlp Dener   bnk->eta4 = 0.90;
1452eb910715SAlp Dener 
1453eb910715SAlp Dener   bnk->alpha1 = 0.25;
1454eb910715SAlp Dener   bnk->alpha2 = 0.50;
1455eb910715SAlp Dener   bnk->alpha3 = 1.00;
1456eb910715SAlp Dener   bnk->alpha4 = 2.00;
1457eb910715SAlp Dener   bnk->alpha5 = 4.00;
1458eb910715SAlp Dener 
1459eb910715SAlp Dener   /*  Default values for trust-region radius update based on interpolation */
1460eb910715SAlp Dener   bnk->mu1 = 0.10;
1461eb910715SAlp Dener   bnk->mu2 = 0.50;
1462eb910715SAlp Dener 
1463eb910715SAlp Dener   bnk->gamma1 = 0.25;
1464eb910715SAlp Dener   bnk->gamma2 = 0.50;
1465eb910715SAlp Dener   bnk->gamma3 = 2.00;
1466eb910715SAlp Dener   bnk->gamma4 = 4.00;
1467eb910715SAlp Dener 
1468eb910715SAlp Dener   bnk->theta = 0.05;
1469eb910715SAlp Dener 
1470eb910715SAlp Dener   /*  Default values for trust region initialization based on interpolation */
1471eb910715SAlp Dener   bnk->mu1_i = 0.35;
1472eb910715SAlp Dener   bnk->mu2_i = 0.50;
1473eb910715SAlp Dener 
1474eb910715SAlp Dener   bnk->gamma1_i = 0.0625;
1475eb910715SAlp Dener   bnk->gamma2_i = 0.5;
1476eb910715SAlp Dener   bnk->gamma3_i = 2.0;
1477eb910715SAlp Dener   bnk->gamma4_i = 5.0;
1478eb910715SAlp Dener 
1479eb910715SAlp Dener   bnk->theta_i = 0.25;
1480eb910715SAlp Dener 
1481eb910715SAlp Dener   /*  Remaining parameters */
1482c0f10754SAlp Dener   bnk->max_cg_its = 0;
1483eb910715SAlp Dener   bnk->min_radius = 1.0e-10;
1484eb910715SAlp Dener   bnk->max_radius = 1.0e10;
1485770b7498SAlp Dener   bnk->epsilon = PetscPowReal(PETSC_MACHINE_EPSILON, 2.0/3.0);
14860a4511e9SAlp Dener   bnk->as_tol = 1.0e-3;
14870a4511e9SAlp Dener   bnk->as_step = 1.0e-3;
148862675beeSAlp Dener   bnk->dmin = 1.0e-6;
148962675beeSAlp Dener   bnk->dmax = 1.0e6;
1490eb910715SAlp Dener 
1491e031d6f5SAlp Dener   bnk->pc_type         = BNK_PC_AHESS;
1492eb910715SAlp Dener   bnk->bfgs_scale_type = BFGS_SCALE_PHESS;
1493eb910715SAlp Dener   bnk->init_type       = BNK_INIT_INTERPOLATION;
1494fed79b8eSAlp Dener   bnk->update_type     = BNK_UPDATE_INTERPOLATION;
14952f75a4aaSAlp Dener   bnk->as_type         = BNK_AS_BERTSEKAS;
1496eb910715SAlp Dener 
1497e031d6f5SAlp Dener   /* Create the embedded BNCG solver */
1498e031d6f5SAlp Dener   ierr = TaoCreate(PetscObjectComm((PetscObject)tao), &bnk->bncg);CHKERRQ(ierr);
1499e031d6f5SAlp Dener   ierr = PetscObjectIncrementTabLevel((PetscObject)bnk->bncg, (PetscObject)tao, 1);CHKERRQ(ierr);
1500e031d6f5SAlp Dener   ierr = TaoSetOptionsPrefix(bnk->bncg, "tao_bnk_");CHKERRQ(ierr);
1501e031d6f5SAlp Dener   ierr = TaoSetType(bnk->bncg, TAOBNCG);CHKERRQ(ierr);
1502e031d6f5SAlp Dener 
1503c0f10754SAlp Dener   /* Create the line search */
1504eb910715SAlp Dener   ierr = TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch);CHKERRQ(ierr);
1505eb910715SAlp Dener   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr);
1506e031d6f5SAlp Dener   ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr);
1507eb910715SAlp Dener   ierr = TaoLineSearchSetType(tao->linesearch,morethuente_type);CHKERRQ(ierr);
1508eb910715SAlp Dener   ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr);
1509eb910715SAlp Dener 
1510eb910715SAlp Dener   /*  Set linear solver to default for symmetric matrices */
1511eb910715SAlp Dener   ierr = KSPCreate(((PetscObject)tao)->comm,&tao->ksp);CHKERRQ(ierr);
1512eb910715SAlp Dener   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1);CHKERRQ(ierr);
1513eb910715SAlp Dener   ierr = KSPSetOptionsPrefix(tao->ksp,tao->hdr.prefix);CHKERRQ(ierr);
1514eb910715SAlp Dener   ierr = KSPSetType(tao->ksp,KSPCGSTCG);CHKERRQ(ierr);
1515eb910715SAlp Dener   PetscFunctionReturn(0);
1516eb910715SAlp Dener }
1517