xref: /petsc/src/tao/bound/impls/bnk/bnk.c (revision 8f8a4e060205d21e1b9befe0ac7d9fb21177bcd2)
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 */
573b063059SAlp Dener   ierr = TaoBoundSolution(tao->solution, tao->XL,tao->XU, 0.0, &nDiff, tao->solution);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);
1763b063059SAlp Dener           ierr = TaoBoundSolution(tao->solution, tao->XL,tao->XU, 0.0, &nDiff, tao->solution);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               }
230*8f8a4e06SAlp Dener             } else if (PetscAbsScalar(kappa - 1.0) <= bnk->mu2_i) {
231eb910715SAlp Dener               /*  Good agreement */
232eb910715SAlp Dener               max_radius = PetscMax(max_radius, tao->trust);
233eb910715SAlp Dener 
234eb910715SAlp Dener               if (tau_max < bnk->gamma2_i) {
235eb910715SAlp Dener                 tau = bnk->gamma2_i;
236eb910715SAlp Dener               } else if (tau_max > bnk->gamma3_i) {
237eb910715SAlp Dener                 tau = bnk->gamma3_i;
238eb910715SAlp Dener               } else {
239eb910715SAlp Dener                 tau = tau_max;
240eb910715SAlp Dener               }
241*8f8a4e06SAlp Dener             } else {
242eb910715SAlp Dener               /*  Not good agreement */
243eb910715SAlp Dener               if (tau_min > 1.0) {
244eb910715SAlp Dener                 tau = bnk->gamma2_i;
245eb910715SAlp Dener               } else if (tau_max < bnk->gamma1_i) {
246eb910715SAlp Dener                 tau = bnk->gamma1_i;
247eb910715SAlp Dener               } else if ((tau_min < bnk->gamma1_i) && (tau_max >= 1.0)) {
248eb910715SAlp Dener                 tau = bnk->gamma1_i;
2493105154fSTodd Munson               } else if ((tau_1 >= bnk->gamma1_i) && (tau_1 < 1.0) && ((tau_2 < bnk->gamma1_i) || (tau_2 >= 1.0))) {
250eb910715SAlp Dener                 tau = tau_1;
2513105154fSTodd Munson               } else if ((tau_2 >= bnk->gamma1_i) && (tau_2 < 1.0) && ((tau_1 < bnk->gamma1_i) || (tau_2 >= 1.0))) {
252eb910715SAlp Dener                 tau = tau_2;
253eb910715SAlp Dener               } else {
254eb910715SAlp Dener                 tau = tau_max;
255eb910715SAlp Dener               }
256eb910715SAlp Dener             }
257eb910715SAlp Dener           }
258eb910715SAlp Dener           tao->trust = tau * tao->trust;
259eb910715SAlp Dener         }
260eb910715SAlp Dener 
2610a4511e9SAlp Dener         if (f_min < bnk->f) {
262937a31a1SAlp Dener           /* We accidentally found a solution better than the initial, so accept it */
2630a4511e9SAlp Dener           bnk->f = f_min;
264937a31a1SAlp Dener           ierr = VecCopy(tao->solution, bnk->Xold);CHKERRQ(ierr);
265eb910715SAlp Dener           ierr = VecAXPY(tao->solution,sigma,tao->gradient);CHKERRQ(ierr);
2663b063059SAlp Dener           ierr = TaoBoundSolution(tao->solution, tao->XL,tao->XU, 0.0, &nDiff, tao->solution);CHKERRQ(ierr);
267937a31a1SAlp Dener           ierr = VecCopy(tao->solution, tao->stepdirection);CHKERRQ(ierr);
268937a31a1SAlp Dener           ierr = VecAXPY(tao->stepdirection, -1.0, bnk->Xold);CHKERRQ(ierr);
26909164190SAlp Dener           ierr = TaoComputeGradient(tao,tao->solution,bnk->unprojected_gradient);CHKERRQ(ierr);
27061be54a6SAlp Dener           ierr = TaoBNKEstimateActiveSet(tao, bnk->as_type);CHKERRQ(ierr);
27161be54a6SAlp Dener           ierr = VecCopy(bnk->unprojected_gradient, tao->gradient);CHKERRQ(ierr);
27261be54a6SAlp Dener           ierr = VecISSet(tao->gradient, bnk->active_idx, 0.0);CHKERRQ(ierr);
273937a31a1SAlp Dener           /* Compute gradient at the new iterate and flip switch to compute the Hessian later */
274c4b75bccSAlp Dener           ierr = VecNorm(tao->gradient, NORM_2, &bnk->gnorm);CHKERRQ(ierr);
275c0f10754SAlp Dener           *needH = PETSC_TRUE;
276937a31a1SAlp Dener           /* Test the new step for convergence */
27789da521bSAlp Dener           ierr = VecFischer(tao->solution, bnk->unprojected_gradient, tao->XL, tao->XU, bnk->W);CHKERRQ(ierr);
27889da521bSAlp Dener           ierr = VecNorm(bnk->W, NORM_2, &resnorm);CHKERRQ(ierr);
279b4a30f08SAlp Dener           if (PetscIsInfOrNanReal(resnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
280e031d6f5SAlp Dener           ierr = TaoLogConvergenceHistory(tao,bnk->f,resnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
281e031d6f5SAlp Dener           ierr = TaoMonitor(tao,tao->niter,bnk->f,resnorm,0.0,1.0);CHKERRQ(ierr);
282eb910715SAlp Dener           ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
283eb910715SAlp Dener           if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
284937a31a1SAlp Dener           /* active BNCG recycling early because we have a stepdirection computed */
285937a31a1SAlp Dener           ierr = TaoBNCGSetRecycleFlag(bnk->bncg, PETSC_TRUE);CHKERRQ(ierr);
286eb910715SAlp Dener         }
287eb910715SAlp Dener       }
288eb910715SAlp Dener       tao->trust = PetscMax(tao->trust, max_radius);
289e031d6f5SAlp Dener 
290e031d6f5SAlp Dener       /* Ensure that the trust radius is within the limits */
291e031d6f5SAlp Dener       tao->trust = PetscMax(tao->trust, bnk->min_radius);
292e031d6f5SAlp Dener       tao->trust = PetscMin(tao->trust, bnk->max_radius);
293eb910715SAlp Dener       break;
294eb910715SAlp Dener 
295eb910715SAlp Dener     default:
296eb910715SAlp Dener       /* Norm of the first direction will initialize radius */
297eb910715SAlp Dener       tao->trust = 0.0;
298eb910715SAlp Dener       break;
299eb910715SAlp Dener     }
300eb910715SAlp Dener   }
301e031d6f5SAlp Dener 
302eb910715SAlp Dener   /* Set initial scaling for the BFGS preconditioner
303eb910715SAlp Dener      This step is done after computing the initial trust-region radius
304eb910715SAlp Dener      since the function value may have decreased */
305eb910715SAlp Dener   if (BNK_PC_BFGS == bnk->pc_type) {
3069b6ef848SAlp Dener     delta = 2.0 * PetscMax(1.0, PetscAbsScalar(bnk->f)) / (bnk->gnorm*bnk->gnorm);
307eb910715SAlp Dener     ierr = MatLMVMSetDelta(bnk->M,delta);CHKERRQ(ierr);
308eb910715SAlp Dener   }
309eb910715SAlp Dener   PetscFunctionReturn(0);
310eb910715SAlp Dener }
311eb910715SAlp Dener 
312df278d8fSAlp Dener /*------------------------------------------------------------*/
313df278d8fSAlp Dener 
31462675beeSAlp Dener /* Routine for computing the Hessian and preparing the preconditioner at the new iterate */
31562675beeSAlp Dener 
31662675beeSAlp Dener PetscErrorCode TaoBNKComputeHessian(Tao tao)
31762675beeSAlp Dener {
31862675beeSAlp Dener   PetscErrorCode               ierr;
31962675beeSAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
320c4b75bccSAlp Dener   PetscBool                    diagExists;
321c4b75bccSAlp Dener   PetscReal                    delta;
32262675beeSAlp Dener 
32362675beeSAlp Dener   PetscFunctionBegin;
32462675beeSAlp Dener   /* Compute the Hessian */
32562675beeSAlp Dener   ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr);
32662675beeSAlp Dener   /* Add a correction to the BFGS preconditioner */
32762675beeSAlp Dener   if (BNK_PC_BFGS == bnk->pc_type) {
32862675beeSAlp Dener     ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
32962675beeSAlp Dener     /* Update the BFGS diagonal scaling */
330c4b75bccSAlp Dener     ierr = MatHasOperation(tao->hessian, MATOP_GET_DIAGONAL, &diagExists);CHKERRQ(ierr);
331c4b75bccSAlp Dener     if (BFGS_SCALE_AHESS == bnk->bfgs_scale_type && diagExists) {
33262675beeSAlp Dener       ierr = MatGetDiagonal(tao->hessian, bnk->Diag);CHKERRQ(ierr);
33362675beeSAlp Dener       ierr = VecAbs(bnk->Diag);CHKERRQ(ierr);
33462675beeSAlp Dener       ierr = VecMedian(bnk->Diag_min, bnk->Diag, bnk->Diag_max, bnk->Diag);CHKERRQ(ierr);
33562675beeSAlp Dener       ierr = VecReciprocal(bnk->Diag);CHKERRQ(ierr);
33662675beeSAlp Dener       ierr = MatLMVMSetScale(bnk->M,bnk->Diag);CHKERRQ(ierr);
337c4b75bccSAlp Dener     } else {
338c4b75bccSAlp Dener       /* Fall back onto stand-alone BFGS scaling */
339c4b75bccSAlp Dener       delta = 2.0 * PetscMax(1.0, PetscAbsScalar(bnk->f)) / (bnk->gnorm*bnk->gnorm);
340c4b75bccSAlp Dener       ierr = MatLMVMSetDelta(bnk->M,delta);CHKERRQ(ierr);
34162675beeSAlp Dener     }
34262675beeSAlp Dener   }
34362675beeSAlp Dener   PetscFunctionReturn(0);
34462675beeSAlp Dener }
34562675beeSAlp Dener 
34662675beeSAlp Dener /*------------------------------------------------------------*/
34762675beeSAlp Dener 
3482f75a4aaSAlp Dener /* Routine for estimating the active set */
3492f75a4aaSAlp Dener 
35008752603SAlp Dener PetscErrorCode TaoBNKEstimateActiveSet(Tao tao, PetscInt asType)
3512f75a4aaSAlp Dener {
3522f75a4aaSAlp Dener   PetscErrorCode               ierr;
3532f75a4aaSAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
354c4b75bccSAlp Dener   PetscBool                    hessComputed, diagExists;
3552f75a4aaSAlp Dener 
3562f75a4aaSAlp Dener   PetscFunctionBegin;
35708752603SAlp Dener   switch (asType) {
3582f75a4aaSAlp Dener   case BNK_AS_NONE:
3592f75a4aaSAlp Dener     ierr = ISDestroy(&bnk->inactive_idx);CHKERRQ(ierr);
3602f75a4aaSAlp Dener     ierr = VecWhichInactive(tao->XL, tao->solution, bnk->unprojected_gradient, tao->XU, PETSC_TRUE, &bnk->inactive_idx);CHKERRQ(ierr);
3612f75a4aaSAlp Dener     ierr = ISDestroy(&bnk->active_idx);CHKERRQ(ierr);
3622f75a4aaSAlp Dener     ierr = ISComplementVec(bnk->inactive_idx, tao->solution, &bnk->active_idx);CHKERRQ(ierr);
3632f75a4aaSAlp Dener     break;
3642f75a4aaSAlp Dener 
3652f75a4aaSAlp Dener   case BNK_AS_BERTSEKAS:
3662f75a4aaSAlp Dener     /* Compute the trial step vector with which we will estimate the active set at the next iteration */
3672f75a4aaSAlp Dener     if (BNK_PC_BFGS == bnk->pc_type) {
3682f75a4aaSAlp Dener       /* If the BFGS preconditioner matrix is available, we will construct a trial step with it */
3695e9b73cbSAlp Dener       ierr = MatLMVMSetInactive(bnk->M, NULL);CHKERRQ(ierr);
3702f75a4aaSAlp Dener       ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, bnk->W);CHKERRQ(ierr);
3712f75a4aaSAlp Dener     } else {
37261be54a6SAlp Dener       ierr = MatAssembled(tao->hessian, &hessComputed);CHKERRQ(ierr);
373c4b75bccSAlp Dener       ierr = MatHasOperation(tao->hessian, MATOP_GET_DIAGONAL, &diagExists);CHKERRQ(ierr);
374c4b75bccSAlp Dener       if (hessComputed && diagExists) {
3759b6ef848SAlp Dener         /* BFGS preconditioner doesn't exist so let's invert the absolute diagonal of the Hessian instead onto the gradient */
3762f75a4aaSAlp Dener         ierr = MatGetDiagonal(tao->hessian, bnk->Xwork);CHKERRQ(ierr);
3779b6ef848SAlp Dener         ierr = VecAbs(bnk->Xwork);CHKERRQ(ierr);
3789b6ef848SAlp Dener         ierr = VecMedian(bnk->Diag_min, bnk->Xwork, bnk->Diag_max, bnk->Xwork);CHKERRQ(ierr);
3792f75a4aaSAlp Dener         ierr = VecReciprocal(bnk->Xwork);CHKERRQ(ierr);CHKERRQ(ierr);
3802f75a4aaSAlp Dener         ierr = VecPointwiseMult(bnk->W, bnk->Xwork, bnk->unprojected_gradient);CHKERRQ(ierr);
38161be54a6SAlp Dener       } else {
382c4b75bccSAlp Dener         /* If the Hessian or its diagonal does not exist, we will simply use gradient step */
38361be54a6SAlp Dener         ierr = VecCopy(bnk->unprojected_gradient, bnk->W);CHKERRQ(ierr);
38461be54a6SAlp Dener       }
3852f75a4aaSAlp Dener     }
3862f75a4aaSAlp Dener     ierr = VecScale(bnk->W, -1.0);CHKERRQ(ierr);
38789da521bSAlp Dener     ierr = TaoEstimateActiveBounds(tao->solution, tao->XL, tao->XU, bnk->unprojected_gradient, bnk->W, bnk->Xwork, bnk->as_step, &bnk->as_tol,
38889da521bSAlp Dener                                    &bnk->active_lower, &bnk->active_upper, &bnk->active_fixed, &bnk->active_idx, &bnk->inactive_idx);CHKERRQ(ierr);
389c4b75bccSAlp Dener     break;
3902f75a4aaSAlp Dener 
3912f75a4aaSAlp Dener   default:
3922f75a4aaSAlp Dener     break;
3932f75a4aaSAlp Dener   }
3942f75a4aaSAlp Dener   PetscFunctionReturn(0);
3952f75a4aaSAlp Dener }
3962f75a4aaSAlp Dener 
3972f75a4aaSAlp Dener /*------------------------------------------------------------*/
3982f75a4aaSAlp Dener 
3992f75a4aaSAlp Dener /* Routine for bounding the step direction */
4002f75a4aaSAlp Dener 
401a1318120SAlp Dener PetscErrorCode TaoBNKBoundStep(Tao tao, PetscInt asType, Vec step)
4022f75a4aaSAlp Dener {
4032f75a4aaSAlp Dener   PetscErrorCode               ierr;
4042f75a4aaSAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
4052f75a4aaSAlp Dener 
4062f75a4aaSAlp Dener   PetscFunctionBegin;
407a1318120SAlp Dener   switch (asType) {
4082f75a4aaSAlp Dener   case BNK_AS_NONE:
409c4b75bccSAlp Dener     ierr = VecISSet(step, bnk->active_idx, 0.0);CHKERRQ(ierr);
4102f75a4aaSAlp Dener     break;
4112f75a4aaSAlp Dener 
4122f75a4aaSAlp Dener   case BNK_AS_BERTSEKAS:
413c4b75bccSAlp Dener     ierr = TaoBoundStep(tao->solution, tao->XL, tao->XU, bnk->active_lower, bnk->active_upper, bnk->active_fixed, 1.0, step);CHKERRQ(ierr);
4142f75a4aaSAlp Dener     break;
4152f75a4aaSAlp Dener 
4162f75a4aaSAlp Dener   default:
4172f75a4aaSAlp Dener     break;
4182f75a4aaSAlp Dener   }
4192f75a4aaSAlp Dener   PetscFunctionReturn(0);
4202f75a4aaSAlp Dener }
4212f75a4aaSAlp Dener 
422e031d6f5SAlp Dener /*------------------------------------------------------------*/
423e031d6f5SAlp Dener 
424e031d6f5SAlp Dener /* Routine for taking a finite number of BNCG iterations to
425e031d6f5SAlp Dener    accelerate Newton convergence.
426e031d6f5SAlp Dener 
427e031d6f5SAlp Dener    In practice, this approach simply trades off Hessian evaluations
428e031d6f5SAlp Dener    for more gradient evaluations.
429e031d6f5SAlp Dener */
430e031d6f5SAlp Dener 
431c0f10754SAlp Dener PetscErrorCode TaoBNKTakeCGSteps(Tao tao, PetscBool *terminate)
432c0f10754SAlp Dener {
433c0f10754SAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
434c0f10754SAlp Dener   PetscErrorCode               ierr;
435c0f10754SAlp Dener 
436c0f10754SAlp Dener   PetscFunctionBegin;
437c0f10754SAlp Dener   *terminate = PETSC_FALSE;
438c0f10754SAlp Dener   if (bnk->max_cg_its > 0) {
439c4b75bccSAlp Dener     /* Copy the current function value (important vectors are already shared) */
440c0f10754SAlp Dener     bnk->bncg_ctx->f = bnk->f;
441c0f10754SAlp Dener     /* Take some small finite number of BNCG iterations */
442c0f10754SAlp Dener     ierr = TaoSolve(bnk->bncg);CHKERRQ(ierr);
443c0f10754SAlp Dener     /* Add the number of gradient and function evaluations to the total */
444c0f10754SAlp Dener     tao->nfuncs += bnk->bncg->nfuncs;
445c0f10754SAlp Dener     tao->nfuncgrads += bnk->bncg->nfuncgrads;
446c0f10754SAlp Dener     tao->ngrads += bnk->bncg->ngrads;
447c0f10754SAlp Dener     tao->nhess += bnk->bncg->nhess;
448e031d6f5SAlp Dener     bnk->tot_cg_its += bnk->bncg->niter;
449c4b75bccSAlp Dener     /* Extract the BNCG function value out and save it into BNK */
450c0f10754SAlp Dener     bnk->f = bnk->bncg_ctx->f;
451c0f10754SAlp 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) {
452c0f10754SAlp Dener       *terminate = PETSC_TRUE;
45361be54a6SAlp Dener     } else {
45461be54a6SAlp Dener       ierr = TaoBNKEstimateActiveSet(tao, bnk->as_type);
455c0f10754SAlp Dener     }
456c0f10754SAlp Dener   }
457c0f10754SAlp Dener   PetscFunctionReturn(0);
458c0f10754SAlp Dener }
459c0f10754SAlp Dener 
4602f75a4aaSAlp Dener /*------------------------------------------------------------*/
4612f75a4aaSAlp Dener 
462c0f10754SAlp Dener /* Routine for computing the Newton step. */
463df278d8fSAlp Dener 
46462675beeSAlp Dener PetscErrorCode TaoBNKComputeStep(Tao tao, PetscBool shift, KSPConvergedReason *ksp_reason)
465eb910715SAlp Dener {
466eb910715SAlp Dener   PetscErrorCode               ierr;
467eb910715SAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
468eb910715SAlp Dener 
469e465cd6fSAlp Dener   PetscReal                    delta;
470eb910715SAlp Dener   PetscInt                     bfgsUpdates = 0;
471eb910715SAlp Dener   PetscInt                     kspits;
472c4b75bccSAlp Dener   PetscBool                    diagExists;
473eb910715SAlp Dener 
474eb910715SAlp Dener   PetscFunctionBegin;
47589da521bSAlp Dener   /* If there are no inactive variables left, save some computation and return an adjusted zero step
47689da521bSAlp Dener      that has (l-x) and (u-x) for lower and upper bounded variables. */
47789da521bSAlp Dener   if (!bnk->inactive_idx) {
47889da521bSAlp Dener     ierr = VecSet(tao->stepdirection, 0.0);CHKERRQ(ierr);
479a1318120SAlp Dener     ierr = TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection);CHKERRQ(ierr);
48089da521bSAlp Dener     PetscFunctionReturn(0);
48189da521bSAlp Dener   }
48289da521bSAlp Dener 
4835e9b73cbSAlp Dener   /* Prepare the reduced sub-matrices for the inactive set */
4843105154fSTodd Munson   if (BNK_PC_BFGS == bnk->pc_type) {
4853105154fSTodd Munson     ierr = MatLMVMSetInactive(bnk->M, bnk->inactive_idx);CHKERRQ(ierr);
4863105154fSTodd Munson   }
48789da521bSAlp Dener   if (bnk->active_idx) {
4885e9b73cbSAlp Dener     ierr = MatDestroy(&bnk->H_inactive);
4895e9b73cbSAlp Dener     ierr = MatCreateSubMatrix(tao->hessian, bnk->inactive_idx, bnk->inactive_idx, MAT_INITIAL_MATRIX, &bnk->H_inactive);CHKERRQ(ierr);
490eb3ba6a7SAlp Dener     if (tao->hessian == tao->hessian_pre) {
491eb3ba6a7SAlp Dener       bnk->Hpre_inactive = bnk->H_inactive;
492eb3ba6a7SAlp Dener     } else {
4935e9b73cbSAlp Dener       ierr = MatDestroy(&bnk->Hpre_inactive);
4945e9b73cbSAlp Dener       ierr = MatCreateSubMatrix(tao->hessian_pre, bnk->inactive_idx, bnk->inactive_idx, MAT_INITIAL_MATRIX, &bnk->Hpre_inactive);CHKERRQ(ierr);
495eb3ba6a7SAlp Dener     }
4962f75a4aaSAlp Dener   } else {
49762675beeSAlp Dener     ierr = MatDestroy(&bnk->H_inactive);
49862675beeSAlp Dener     ierr = MatDuplicate(tao->hessian, MAT_COPY_VALUES, &bnk->H_inactive);
49962675beeSAlp Dener     if (tao->hessian == tao->hessian_pre) {
50062675beeSAlp Dener       bnk->Hpre_inactive = bnk->H_inactive;
50162675beeSAlp Dener     } else {
50262675beeSAlp Dener       ierr = MatDestroy(&bnk->Hpre_inactive);
50362675beeSAlp Dener       ierr = MatDuplicate(tao->hessian_pre, MAT_COPY_VALUES, &bnk->Hpre_inactive);
50462675beeSAlp Dener     }
50562675beeSAlp Dener   }
50662675beeSAlp Dener 
50762675beeSAlp Dener   /* Shift the reduced Hessian matrix */
50862675beeSAlp Dener   if ((shift) && (bnk->pert > 0)) {
50962675beeSAlp Dener     ierr = MatShift(bnk->H_inactive, bnk->pert);CHKERRQ(ierr);
51062675beeSAlp Dener     if (bnk->H_inactive != bnk->Hpre_inactive) {
51162675beeSAlp Dener       ierr = MatShift(bnk->Hpre_inactive, bnk->pert);CHKERRQ(ierr);
51262675beeSAlp Dener     }
51362675beeSAlp Dener   }
51462675beeSAlp Dener 
51562675beeSAlp Dener   /* Update the diagonal scaling for the BFGS preconditioner, this time with the Hessian perturbation */
51662675beeSAlp Dener   if ((BNK_PC_BFGS == bnk->pc_type) && (BFGS_SCALE_PHESS == bnk->bfgs_scale_type)) {
517c4b75bccSAlp Dener     ierr = MatHasOperation(bnk->H_inactive, MATOP_GET_DIAGONAL, &diagExists);CHKERRQ(ierr);
518c4b75bccSAlp Dener     if (diagExists) {
51962675beeSAlp Dener       /* Obtain diagonal for the bfgs preconditioner  */
5205e9b73cbSAlp Dener       ierr = VecSet(bnk->Diag, 1.0);CHKERRQ(ierr);
52189da521bSAlp Dener       if (bnk->active_idx) {
5225e9b73cbSAlp Dener         ierr = VecGetSubVector(bnk->Diag, bnk->inactive_idx, &bnk->Diag_red);CHKERRQ(ierr);
5235e9b73cbSAlp Dener       } else {
5245e9b73cbSAlp Dener         bnk->Diag_red = bnk->Diag;
5255e9b73cbSAlp Dener       }
5265e9b73cbSAlp Dener       ierr = MatGetDiagonal(bnk->H_inactive, bnk->Diag_red);CHKERRQ(ierr);
52789da521bSAlp Dener       if (bnk->active_idx) {
5285e9b73cbSAlp Dener         ierr = VecRestoreSubVector(bnk->Diag, bnk->inactive_idx, &bnk->Diag_red);CHKERRQ(ierr);
5295e9b73cbSAlp Dener       }
53062675beeSAlp Dener       ierr = VecAbs(bnk->Diag);CHKERRQ(ierr);
53162675beeSAlp Dener       ierr = VecMedian(bnk->Diag_min, bnk->Diag, bnk->Diag_max, bnk->Diag);CHKERRQ(ierr);
53262675beeSAlp Dener       ierr = VecReciprocal(bnk->Diag);CHKERRQ(ierr);
53362675beeSAlp Dener       ierr = MatLMVMSetScale(bnk->M,bnk->Diag);CHKERRQ(ierr);
534c4b75bccSAlp Dener     } else {
535c4b75bccSAlp Dener       /* Fall back onto stand-alone BFGS scaling */
536c4b75bccSAlp Dener       delta = 2.0 * PetscMax(1.0, PetscAbsScalar(bnk->f)) / (bnk->gnorm*bnk->gnorm);
537c4b75bccSAlp Dener       ierr = MatLMVMSetDelta(bnk->M,delta);CHKERRQ(ierr);
538c4b75bccSAlp Dener     }
5392f75a4aaSAlp Dener   }
54009164190SAlp Dener 
541eb910715SAlp Dener   /* Solve the Newton system of equations */
542937a31a1SAlp Dener   tao->ksp_its = 0;
5432f75a4aaSAlp Dener   ierr = VecSet(tao->stepdirection, 0.0);CHKERRQ(ierr);
5445e9b73cbSAlp Dener   ierr = KSPReset(tao->ksp);CHKERRQ(ierr);
54509164190SAlp Dener   ierr = KSPSetOperators(tao->ksp,bnk->H_inactive,bnk->Hpre_inactive);CHKERRQ(ierr);
5465e9b73cbSAlp Dener   ierr = VecCopy(bnk->unprojected_gradient, bnk->Gwork);CHKERRQ(ierr);
54789da521bSAlp Dener   if (bnk->active_idx) {
5485e9b73cbSAlp Dener     ierr = VecGetSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive);CHKERRQ(ierr);
5495e9b73cbSAlp Dener     ierr = VecGetSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr);
5505e9b73cbSAlp Dener   } else {
5515e9b73cbSAlp Dener     bnk->G_inactive = bnk->unprojected_gradient;
5525e9b73cbSAlp Dener     bnk->X_inactive = tao->stepdirection;
55328017e9fSAlp Dener   }
554eb910715SAlp Dener   if (bnk->is_nash || bnk->is_stcg || bnk->is_gltr) {
555fed79b8eSAlp Dener     ierr = KSPCGSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr);
5565e9b73cbSAlp Dener     ierr = KSPSolve(tao->ksp, bnk->G_inactive, bnk->X_inactive);CHKERRQ(ierr);
557eb910715SAlp Dener     ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr);
558eb910715SAlp Dener     tao->ksp_its+=kspits;
559eb910715SAlp Dener     tao->ksp_tot_its+=kspits;
560080d2917SAlp Dener     ierr = KSPCGGetNormD(tao->ksp,&bnk->dnorm);CHKERRQ(ierr);
561eb910715SAlp Dener 
562eb910715SAlp Dener     if (0.0 == tao->trust) {
563eb910715SAlp Dener       /* Radius was uninitialized; use the norm of the direction */
564080d2917SAlp Dener       if (bnk->dnorm > 0.0) {
565080d2917SAlp Dener         tao->trust = bnk->dnorm;
566eb910715SAlp Dener 
567eb910715SAlp Dener         /* Modify the radius if it is too large or small */
568eb910715SAlp Dener         tao->trust = PetscMax(tao->trust, bnk->min_radius);
569eb910715SAlp Dener         tao->trust = PetscMin(tao->trust, bnk->max_radius);
570eb910715SAlp Dener       } else {
571eb910715SAlp Dener         /* The direction was bad; set radius to default value and re-solve
572eb910715SAlp Dener            the trust-region subproblem to get a direction */
573eb910715SAlp Dener         tao->trust = tao->trust0;
574eb910715SAlp Dener 
575eb910715SAlp Dener         /* Modify the radius if it is too large or small */
576eb910715SAlp Dener         tao->trust = PetscMax(tao->trust, bnk->min_radius);
577eb910715SAlp Dener         tao->trust = PetscMin(tao->trust, bnk->max_radius);
578eb910715SAlp Dener 
579fed79b8eSAlp Dener         ierr = KSPCGSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr);
5805e9b73cbSAlp Dener         ierr = KSPSolve(tao->ksp, bnk->G_inactive, bnk->X_inactive);CHKERRQ(ierr);
581eb910715SAlp Dener         ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr);
582eb910715SAlp Dener         tao->ksp_its+=kspits;
583eb910715SAlp Dener         tao->ksp_tot_its+=kspits;
584080d2917SAlp Dener         ierr = KSPCGGetNormD(tao->ksp,&bnk->dnorm);CHKERRQ(ierr);
585eb910715SAlp Dener 
586080d2917SAlp Dener         if (bnk->dnorm == 0.0) SETERRQ(PETSC_COMM_SELF,1, "Initial direction zero");
587eb910715SAlp Dener       }
588eb910715SAlp Dener     }
589eb910715SAlp Dener   } else {
5905e9b73cbSAlp Dener     ierr = KSPSolve(tao->ksp, bnk->G_inactive, bnk->X_inactive);CHKERRQ(ierr);
591eb910715SAlp Dener     ierr = KSPGetIterationNumber(tao->ksp, &kspits);CHKERRQ(ierr);
592eb910715SAlp Dener     tao->ksp_its += kspits;
593eb910715SAlp Dener     tao->ksp_tot_its+=kspits;
594eb910715SAlp Dener   }
5955e9b73cbSAlp Dener   /* Restore sub vectors back */
59689da521bSAlp Dener   if (bnk->active_idx) {
5975e9b73cbSAlp Dener     ierr = VecRestoreSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive);CHKERRQ(ierr);
5985e9b73cbSAlp Dener     ierr = VecRestoreSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr);
5995e9b73cbSAlp Dener   }
600770b7498SAlp Dener   /* Make sure the safeguarded fall-back step is zero for actively bounded variables */
601fed79b8eSAlp Dener   ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
602a1318120SAlp Dener   ierr = TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection);CHKERRQ(ierr);
603770b7498SAlp Dener 
604770b7498SAlp Dener   /* Record convergence reasons */
605e465cd6fSAlp Dener   ierr = KSPGetConvergedReason(tao->ksp, ksp_reason);CHKERRQ(ierr);
606e465cd6fSAlp Dener   if (KSP_CONVERGED_ATOL == *ksp_reason) {
607770b7498SAlp Dener     ++bnk->ksp_atol;
608e465cd6fSAlp Dener   } else if (KSP_CONVERGED_RTOL == *ksp_reason) {
609770b7498SAlp Dener     ++bnk->ksp_rtol;
610e465cd6fSAlp Dener   } else if (KSP_CONVERGED_CG_CONSTRAINED == *ksp_reason) {
611770b7498SAlp Dener     ++bnk->ksp_ctol;
612e465cd6fSAlp Dener   } else if (KSP_CONVERGED_CG_NEG_CURVE == *ksp_reason) {
613770b7498SAlp Dener     ++bnk->ksp_negc;
614e465cd6fSAlp Dener   } else if (KSP_DIVERGED_DTOL == *ksp_reason) {
615770b7498SAlp Dener     ++bnk->ksp_dtol;
616e465cd6fSAlp Dener   } else if (KSP_DIVERGED_ITS == *ksp_reason) {
617770b7498SAlp Dener     ++bnk->ksp_iter;
618770b7498SAlp Dener   } else {
619770b7498SAlp Dener     ++bnk->ksp_othr;
620770b7498SAlp Dener   }
621fed79b8eSAlp Dener 
622fed79b8eSAlp Dener   /* Make sure the BFGS preconditioner is healthy */
623fed79b8eSAlp Dener   if (bnk->pc_type == BNK_PC_BFGS) {
624fed79b8eSAlp Dener     ierr = MatLMVMGetUpdates(bnk->M, &bfgsUpdates);CHKERRQ(ierr);
625e465cd6fSAlp Dener     if ((KSP_DIVERGED_INDEFINITE_PC == *ksp_reason) && (bfgsUpdates > 1)) {
626fed79b8eSAlp Dener       /* Preconditioner is numerically indefinite; reset the approximation. */
6279b6ef848SAlp Dener       delta = 2.0 * PetscMax(1.0, PetscAbsScalar(bnk->f)) / (bnk->gnorm*bnk->gnorm);
628eb910715SAlp Dener       ierr = MatLMVMSetDelta(bnk->M,delta);CHKERRQ(ierr);
629eb910715SAlp Dener       ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr);
63009164190SAlp Dener       ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
631eb910715SAlp Dener     }
632fed79b8eSAlp Dener   }
633e465cd6fSAlp Dener   PetscFunctionReturn(0);
634e465cd6fSAlp Dener }
635eb910715SAlp Dener 
63662675beeSAlp Dener /*------------------------------------------------------------*/
63762675beeSAlp Dener 
6385e9b73cbSAlp Dener /* Routine for recomputing the predicted reduction for a given step vector */
6395e9b73cbSAlp Dener 
6405e9b73cbSAlp Dener PetscErrorCode TaoBNKRecomputePred(Tao tao, Vec S, PetscReal *prered)
6415e9b73cbSAlp Dener {
6425e9b73cbSAlp Dener   PetscErrorCode               ierr;
6435e9b73cbSAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
6445e9b73cbSAlp Dener 
6455e9b73cbSAlp Dener   PetscFunctionBegin;
6465e9b73cbSAlp Dener   /* Extract subvectors associated with the inactive set */
64789da521bSAlp Dener   if (bnk->active_idx){
6485e9b73cbSAlp Dener     ierr = VecGetSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr);
6495e9b73cbSAlp Dener     ierr = VecGetSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work);CHKERRQ(ierr);
6505e9b73cbSAlp Dener     ierr = VecGetSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive);CHKERRQ(ierr);
6515e9b73cbSAlp Dener   } else {
6525e9b73cbSAlp Dener     bnk->X_inactive = tao->stepdirection;
6535e9b73cbSAlp Dener     bnk->inactive_work = bnk->Xwork;
6545e9b73cbSAlp Dener     bnk->G_inactive = bnk->Gwork;
6555e9b73cbSAlp Dener   }
6565e9b73cbSAlp Dener   /* Recompute the predicted decrease based on the quadratic model */
6575e9b73cbSAlp Dener   ierr = MatMult(bnk->H_inactive, bnk->X_inactive, bnk->inactive_work);CHKERRQ(ierr);
6585e9b73cbSAlp Dener   ierr = VecAYPX(bnk->inactive_work, -0.5, bnk->G_inactive);CHKERRQ(ierr);
6595e9b73cbSAlp Dener   ierr = VecDot(bnk->inactive_work, bnk->X_inactive, prered);
6605e9b73cbSAlp Dener   /* Restore the sub vectors */
66189da521bSAlp Dener   if (bnk->active_idx){
6625e9b73cbSAlp Dener     ierr = VecRestoreSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr);
6635e9b73cbSAlp Dener     ierr = VecRestoreSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work);CHKERRQ(ierr);
6645e9b73cbSAlp Dener     ierr = VecRestoreSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive);CHKERRQ(ierr);
6655e9b73cbSAlp Dener   }
6665e9b73cbSAlp Dener   PetscFunctionReturn(0);
6675e9b73cbSAlp Dener }
6685e9b73cbSAlp Dener 
6695e9b73cbSAlp Dener /*------------------------------------------------------------*/
6705e9b73cbSAlp Dener 
67162675beeSAlp Dener /* Routine for ensuring that the Newton step is a descent direction.
67262675beeSAlp Dener 
67362675beeSAlp Dener    The step direction falls back onto BFGS, scaled gradient and gradient steps
67462675beeSAlp Dener    in the event that the Newton step fails the test.
67562675beeSAlp Dener */
67662675beeSAlp Dener 
677e465cd6fSAlp Dener PetscErrorCode TaoBNKSafeguardStep(Tao tao, KSPConvergedReason ksp_reason, PetscInt *stepType)
678e465cd6fSAlp Dener {
679e465cd6fSAlp Dener   PetscErrorCode               ierr;
680e465cd6fSAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
681e465cd6fSAlp Dener 
682e465cd6fSAlp Dener   PetscReal                    gdx, delta, e_min;
683e465cd6fSAlp Dener   PetscInt                     bfgsUpdates;
684e465cd6fSAlp Dener 
685e465cd6fSAlp Dener   PetscFunctionBegin;
686080d2917SAlp Dener   ierr = VecDot(tao->stepdirection, tao->gradient, &gdx);CHKERRQ(ierr);
687eb910715SAlp Dener   if ((gdx >= 0.0) || PetscIsInfOrNanReal(gdx)) {
688eb910715SAlp Dener     /* Newton step is not descent or direction produced Inf or NaN
689eb910715SAlp Dener        Update the perturbation for next time */
690eb910715SAlp Dener     if (bnk->pert <= 0.0) {
691eb910715SAlp Dener       /* Initialize the perturbation */
692eb910715SAlp Dener       bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
693eb910715SAlp Dener       if (bnk->is_gltr) {
694eb910715SAlp Dener         ierr = KSPCGGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr);
695eb910715SAlp Dener         bnk->pert = PetscMax(bnk->pert, -e_min);
696eb910715SAlp Dener       }
697eb910715SAlp Dener     } else {
698eb910715SAlp Dener       /* Increase the perturbation */
699eb910715SAlp Dener       bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
700eb910715SAlp Dener     }
701eb910715SAlp Dener 
702eb910715SAlp Dener     if (BNK_PC_BFGS != bnk->pc_type) {
703eb910715SAlp Dener       /* We don't have the bfgs matrix around and updated
704eb910715SAlp Dener          Must use gradient direction in this case */
705080d2917SAlp Dener       ierr = VecCopy(tao->gradient, tao->stepdirection);CHKERRQ(ierr);
706eb910715SAlp Dener       *stepType = BNK_GRADIENT;
707eb910715SAlp Dener     } else {
708eb910715SAlp Dener       /* Attempt to use the BFGS direction */
7092ab2a32cSAlp Dener       ierr = MatLMVMSetInactive(bnk->M, NULL);CHKERRQ(ierr);
71009164190SAlp Dener       ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
711eb910715SAlp Dener 
7128d5ead36SAlp Dener       /* Check for success (descent direction)
7138d5ead36SAlp Dener          NOTE: Negative gdx here means not a descent direction because
7148d5ead36SAlp Dener          the fall-back step is missing a negative sign. */
715080d2917SAlp Dener       ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr);
7163105154fSTodd Munson       if ((gdx <= 0.0) || PetscIsInfOrNanReal(gdx)) {
717eb910715SAlp Dener         /* BFGS direction is not descent or direction produced not a number
718eb910715SAlp Dener            We can assert bfgsUpdates > 1 in this case because
719eb910715SAlp Dener            the first solve produces the scaled gradient direction,
720eb910715SAlp Dener            which is guaranteed to be descent */
721eb910715SAlp Dener 
722eb910715SAlp Dener         /* Use steepest descent direction (scaled) */
7239b6ef848SAlp Dener         delta = 2.0 * PetscMax(1.0, PetscAbsScalar(bnk->f)) / (bnk->gnorm*bnk->gnorm);
724eb910715SAlp Dener         ierr = MatLMVMSetDelta(bnk->M, delta);CHKERRQ(ierr);
725eb910715SAlp Dener         ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr);
72609164190SAlp Dener         ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
72709164190SAlp Dener         ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
728eb910715SAlp Dener 
729eb910715SAlp Dener         *stepType = BNK_SCALED_GRADIENT;
730eb910715SAlp Dener       } else {
731770b7498SAlp Dener         ierr = MatLMVMGetUpdates(bnk->M, &bfgsUpdates);CHKERRQ(ierr);
732eb910715SAlp Dener         if (1 == bfgsUpdates) {
733eb910715SAlp Dener           /* The first BFGS direction is always the scaled gradient */
734eb910715SAlp Dener           *stepType = BNK_SCALED_GRADIENT;
735eb910715SAlp Dener         } else {
736eb910715SAlp Dener           *stepType = BNK_BFGS;
737eb910715SAlp Dener         }
738eb910715SAlp Dener       }
739eb910715SAlp Dener     }
7408d5ead36SAlp Dener     /* Make sure the safeguarded fall-back step is zero for actively bounded variables */
7418d5ead36SAlp Dener     ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
742a1318120SAlp Dener     ierr = TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection);CHKERRQ(ierr);
743eb910715SAlp Dener   } else {
744eb910715SAlp Dener     /* Computed Newton step is descent */
745eb910715SAlp Dener     switch (ksp_reason) {
746eb910715SAlp Dener     case KSP_DIVERGED_NANORINF:
747eb910715SAlp Dener     case KSP_DIVERGED_BREAKDOWN:
748eb910715SAlp Dener     case KSP_DIVERGED_INDEFINITE_MAT:
749eb910715SAlp Dener     case KSP_DIVERGED_INDEFINITE_PC:
750eb910715SAlp Dener     case KSP_CONVERGED_CG_NEG_CURVE:
751eb910715SAlp Dener       /* Matrix or preconditioner is indefinite; increase perturbation */
752eb910715SAlp Dener       if (bnk->pert <= 0.0) {
753eb910715SAlp Dener         /* Initialize the perturbation */
754eb910715SAlp Dener         bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
755eb910715SAlp Dener         if (bnk->is_gltr) {
756eb910715SAlp Dener           ierr = KSPCGGLTRGetMinEig(tao->ksp, &e_min);CHKERRQ(ierr);
757eb910715SAlp Dener           bnk->pert = PetscMax(bnk->pert, -e_min);
758eb910715SAlp Dener         }
759eb910715SAlp Dener       } else {
760eb910715SAlp Dener         /* Increase the perturbation */
761eb910715SAlp Dener         bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
762eb910715SAlp Dener       }
763eb910715SAlp Dener       break;
764eb910715SAlp Dener 
765eb910715SAlp Dener     default:
766eb910715SAlp Dener       /* Newton step computation is good; decrease perturbation */
767eb910715SAlp Dener       bnk->pert = PetscMin(bnk->psfac * bnk->pert, bnk->pmsfac * bnk->gnorm);
768eb910715SAlp Dener       if (bnk->pert < bnk->pmin) {
769eb910715SAlp Dener         bnk->pert = 0.0;
770eb910715SAlp Dener       }
771eb910715SAlp Dener       break;
772eb910715SAlp Dener     }
773fed79b8eSAlp Dener     *stepType = BNK_NEWTON;
774eb910715SAlp Dener   }
775eb910715SAlp Dener   PetscFunctionReturn(0);
776eb910715SAlp Dener }
777eb910715SAlp Dener 
778df278d8fSAlp Dener /*------------------------------------------------------------*/
779df278d8fSAlp Dener 
780df278d8fSAlp Dener /* Routine for performing a bound-projected More-Thuente line search.
781df278d8fSAlp Dener 
782df278d8fSAlp Dener   Includes fallbacks to BFGS, scaled gradient, and unscaled gradient steps if the
783df278d8fSAlp Dener   Newton step does not produce a valid step length.
784df278d8fSAlp Dener */
785df278d8fSAlp Dener 
786937a31a1SAlp Dener PetscErrorCode TaoBNKPerformLineSearch(Tao tao, PetscInt *stepType, PetscReal *steplen, TaoLineSearchConvergedReason *reason)
787c14b763aSAlp Dener {
788c14b763aSAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
789c14b763aSAlp Dener   PetscErrorCode ierr;
790c14b763aSAlp Dener   TaoLineSearchConvergedReason ls_reason;
791c14b763aSAlp Dener 
792c14b763aSAlp Dener   PetscReal      e_min, gdx, delta;
793c14b763aSAlp Dener   PetscInt       bfgsUpdates;
794c14b763aSAlp Dener 
795c14b763aSAlp Dener   PetscFunctionBegin;
796c14b763aSAlp Dener   /* Perform the linesearch */
797c14b763aSAlp Dener   ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, steplen, &ls_reason);CHKERRQ(ierr);
798c14b763aSAlp Dener   ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
799c14b763aSAlp Dener 
800937a31a1SAlp Dener   while (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER && *stepType != BNK_GRADIENT) {
801c14b763aSAlp Dener     /* Linesearch failed, revert solution */
802c14b763aSAlp Dener     bnk->f = bnk->fold;
803c14b763aSAlp Dener     ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr);
804c14b763aSAlp Dener     ierr = VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient);CHKERRQ(ierr);
805c14b763aSAlp Dener 
806937a31a1SAlp Dener     switch(*stepType) {
807c14b763aSAlp Dener     case BNK_NEWTON:
8088d5ead36SAlp Dener       /* Failed to obtain acceptable iterate with Newton step
809c14b763aSAlp Dener          Update the perturbation for next time */
810c14b763aSAlp Dener       if (bnk->pert <= 0.0) {
811c14b763aSAlp Dener         /* Initialize the perturbation */
812c14b763aSAlp Dener         bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
813c14b763aSAlp Dener         if (bnk->is_gltr) {
814c14b763aSAlp Dener           ierr = KSPCGGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr);
815c14b763aSAlp Dener           bnk->pert = PetscMax(bnk->pert, -e_min);
816c14b763aSAlp Dener         }
817c14b763aSAlp Dener       } else {
818c14b763aSAlp Dener         /* Increase the perturbation */
819c14b763aSAlp Dener         bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
820c14b763aSAlp Dener       }
821c14b763aSAlp Dener 
822c14b763aSAlp Dener       if (BNK_PC_BFGS != bnk->pc_type) {
823c14b763aSAlp Dener         /* We don't have the bfgs matrix around and being updated
824c14b763aSAlp Dener            Must use gradient direction in this case */
825937a31a1SAlp Dener         ierr = VecCopy(bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
826937a31a1SAlp Dener         *stepType = BNK_GRADIENT;
827c14b763aSAlp Dener       } else {
828c14b763aSAlp Dener         /* Attempt to use the BFGS direction */
829937a31a1SAlp Dener         ierr = MatLMVMSetInactive(bnk->M, NULL);CHKERRQ(ierr);
830c14b763aSAlp Dener         ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
8318d5ead36SAlp Dener         /* Check for success (descent direction)
8328d5ead36SAlp Dener            NOTE: Negative gdx means not a descent direction because the step here is missing a negative sign. */
833c14b763aSAlp Dener         ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr);
8343105154fSTodd Munson         if ((gdx <= 0.0) || PetscIsInfOrNanReal(gdx)) {
835c14b763aSAlp Dener           /* BFGS direction is not descent or direction produced not a number
836c14b763aSAlp Dener              We can assert bfgsUpdates > 1 in this case
837c14b763aSAlp Dener              Use steepest descent direction (scaled) */
8389b6ef848SAlp Dener           delta = 2.0 * PetscMax(1.0, PetscAbsScalar(bnk->f)) / (bnk->gnorm*bnk->gnorm);
839c14b763aSAlp Dener           ierr = MatLMVMSetDelta(bnk->M, delta);CHKERRQ(ierr);
840c14b763aSAlp Dener           ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr);
841c14b763aSAlp Dener           ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
842c14b763aSAlp Dener           ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
843c14b763aSAlp Dener 
844c14b763aSAlp Dener           bfgsUpdates = 1;
845937a31a1SAlp Dener           *stepType = BNK_SCALED_GRADIENT;
846c14b763aSAlp Dener         } else {
847c14b763aSAlp Dener           ierr = MatLMVMGetUpdates(bnk->M, &bfgsUpdates);CHKERRQ(ierr);
848c14b763aSAlp Dener           if (1 == bfgsUpdates) {
849c14b763aSAlp Dener             /* The first BFGS direction is always the scaled gradient */
850937a31a1SAlp Dener             *stepType = BNK_SCALED_GRADIENT;
851c14b763aSAlp Dener           } else {
852937a31a1SAlp Dener             *stepType = BNK_BFGS;
853c14b763aSAlp Dener           }
854c14b763aSAlp Dener         }
855c14b763aSAlp Dener       }
856c14b763aSAlp Dener       break;
857c14b763aSAlp Dener 
858c14b763aSAlp Dener     case BNK_BFGS:
859c14b763aSAlp Dener       /* Can only enter if pc_type == BNK_PC_BFGS
860c14b763aSAlp Dener          Failed to obtain acceptable iterate with BFGS step
861c14b763aSAlp Dener          Attempt to use the scaled gradient direction */
8629b6ef848SAlp Dener       delta = 2.0 * PetscMax(1.0, PetscAbsScalar(bnk->f)) / (bnk->gnorm*bnk->gnorm);
863c14b763aSAlp Dener       ierr = MatLMVMSetDelta(bnk->M, delta);CHKERRQ(ierr);
864c14b763aSAlp Dener       ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr);
865c14b763aSAlp Dener       ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
866937a31a1SAlp Dener       ierr = MatLMVMSetInactive(bnk->M, NULL);CHKERRQ(ierr);
867c14b763aSAlp Dener       ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
868c14b763aSAlp Dener 
869c14b763aSAlp Dener       bfgsUpdates = 1;
870937a31a1SAlp Dener       *stepType = BNK_SCALED_GRADIENT;
871c14b763aSAlp Dener       break;
872c14b763aSAlp Dener 
873c14b763aSAlp Dener     case BNK_SCALED_GRADIENT:
874c14b763aSAlp Dener       /* Can only enter if pc_type == BNK_PC_BFGS
875c14b763aSAlp Dener          The scaled gradient step did not produce a new iterate;
876937a31a1SAlp Dener          reset the BFGS matrix and attemp to use the gradient direction. */
877c14b763aSAlp Dener       ierr = MatLMVMSetScale(bnk->M,0);CHKERRQ(ierr);
878c14b763aSAlp Dener       ierr = MatLMVMSetDelta(bnk->M,1.0);CHKERRQ(ierr);
879c14b763aSAlp Dener       ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr);
880c14b763aSAlp Dener       ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
881937a31a1SAlp Dener       ierr = VecCopy(bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
882c14b763aSAlp Dener 
883c14b763aSAlp Dener       bfgsUpdates = 1;
884937a31a1SAlp Dener       *stepType = BNK_GRADIENT;
885c14b763aSAlp Dener       break;
886c14b763aSAlp Dener     }
8878d5ead36SAlp Dener     /* Make sure the safeguarded fall-back step is zero for actively bounded variables */
8888d5ead36SAlp Dener     ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
889a1318120SAlp Dener     ierr = TaoBNKBoundStep(tao, bnk->as_type, tao->stepdirection);CHKERRQ(ierr);
890c14b763aSAlp Dener 
8918d5ead36SAlp Dener     /* Perform one last line search with the fall-back step */
892c14b763aSAlp Dener     ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, steplen, &ls_reason);CHKERRQ(ierr);
893c14b763aSAlp Dener     ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
894c14b763aSAlp Dener   }
895c14b763aSAlp Dener   *reason = ls_reason;
896c14b763aSAlp Dener   PetscFunctionReturn(0);
897c14b763aSAlp Dener }
898c14b763aSAlp Dener 
899df278d8fSAlp Dener /*------------------------------------------------------------*/
900df278d8fSAlp Dener 
901df278d8fSAlp Dener /* Routine for updating the trust radius.
902df278d8fSAlp Dener 
903df278d8fSAlp Dener   Function features three different update methods:
904df278d8fSAlp Dener   1) Line-search step length based
905df278d8fSAlp Dener   2) Predicted decrease on the CG quadratic model
906df278d8fSAlp Dener   3) Interpolation
907df278d8fSAlp Dener */
908df278d8fSAlp Dener 
90928017e9fSAlp Dener PetscErrorCode TaoBNKUpdateTrustRadius(Tao tao, PetscReal prered, PetscReal actred, PetscInt updateType, PetscInt stepType, PetscBool *accept)
910080d2917SAlp Dener {
911080d2917SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
912080d2917SAlp Dener   PetscErrorCode ierr;
913080d2917SAlp Dener 
914b1c2d0e3SAlp Dener   PetscReal      step, kappa;
915080d2917SAlp Dener   PetscReal      gdx, tau_1, tau_2, tau_min, tau_max;
916080d2917SAlp Dener 
917080d2917SAlp Dener   PetscFunctionBegin;
918080d2917SAlp Dener   /* Update trust region radius */
919080d2917SAlp Dener   *accept = PETSC_FALSE;
92028017e9fSAlp Dener   switch(updateType) {
921080d2917SAlp Dener   case BNK_UPDATE_STEP:
922c14b763aSAlp Dener     *accept = PETSC_TRUE; /* always accept here because line search succeeded */
923080d2917SAlp Dener     if (stepType == BNK_NEWTON) {
924080d2917SAlp Dener       ierr = TaoLineSearchGetStepLength(tao->linesearch, &step);CHKERRQ(ierr);
925080d2917SAlp Dener       if (step < bnk->nu1) {
926080d2917SAlp Dener         /* Very bad step taken; reduce radius */
927080d2917SAlp Dener         tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust);
928080d2917SAlp Dener       } else if (step < bnk->nu2) {
929080d2917SAlp Dener         /* Reasonably bad step taken; reduce radius */
930080d2917SAlp Dener         tao->trust = bnk->omega2 * PetscMin(bnk->dnorm, tao->trust);
931080d2917SAlp Dener       } else if (step < bnk->nu3) {
932080d2917SAlp Dener         /*  Reasonable step was taken; leave radius alone */
933080d2917SAlp Dener         if (bnk->omega3 < 1.0) {
934080d2917SAlp Dener           tao->trust = bnk->omega3 * PetscMin(bnk->dnorm, tao->trust);
935080d2917SAlp Dener         } else if (bnk->omega3 > 1.0) {
936080d2917SAlp Dener           tao->trust = PetscMax(bnk->omega3 * bnk->dnorm, tao->trust);
937080d2917SAlp Dener         }
938080d2917SAlp Dener       } else if (step < bnk->nu4) {
939080d2917SAlp Dener         /*  Full step taken; increase the radius */
940080d2917SAlp Dener         tao->trust = PetscMax(bnk->omega4 * bnk->dnorm, tao->trust);
941080d2917SAlp Dener       } else {
942080d2917SAlp Dener         /*  More than full step taken; increase the radius */
943080d2917SAlp Dener         tao->trust = PetscMax(bnk->omega5 * bnk->dnorm, tao->trust);
944080d2917SAlp Dener       }
945080d2917SAlp Dener     } else {
946080d2917SAlp Dener       /*  Newton step was not good; reduce the radius */
947080d2917SAlp Dener       tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust);
948080d2917SAlp Dener     }
949080d2917SAlp Dener     break;
950080d2917SAlp Dener 
951080d2917SAlp Dener   case BNK_UPDATE_REDUCTION:
952080d2917SAlp Dener     if (stepType == BNK_NEWTON) {
953b1c2d0e3SAlp Dener       if (prered < 0.0) {
954fed79b8eSAlp Dener         /* The predicted reduction has the wrong sign.  This cannot
955fed79b8eSAlp Dener            happen in infinite precision arithmetic.  Step should
956fed79b8eSAlp Dener            be rejected! */
957080d2917SAlp Dener         tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
9583105154fSTodd Munson       } else {
959b1c2d0e3SAlp Dener         if (PetscIsInfOrNanReal(actred)) {
960080d2917SAlp Dener           tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
961080d2917SAlp Dener         } else {
9623105154fSTodd Munson           if ((PetscAbsScalar(actred) <= PetscMax(1.0, PetscAbsScalar(bnk->f))*bnk->epsilon) && (PetscAbsScalar(prered) <= PetscMax(1.0, PetscAbsScalar(bnk->f))*bnk->epsilon)) {
963080d2917SAlp Dener             kappa = 1.0;
9643105154fSTodd Munson           } else {
965080d2917SAlp Dener             kappa = actred / prered;
966080d2917SAlp Dener           }
967fed79b8eSAlp Dener 
968fed79b8eSAlp Dener           /* Accept or reject the step and update radius */
969080d2917SAlp Dener           if (kappa < bnk->eta1) {
970fed79b8eSAlp Dener             /* Reject the step */
971080d2917SAlp Dener             tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
9723105154fSTodd Munson           } else {
973fed79b8eSAlp Dener             /* Accept the step */
974c133c014SAlp Dener             *accept = PETSC_TRUE;
975c133c014SAlp Dener             /* Update the trust region radius only if the computed step is at the trust radius boundary */
9768d5ead36SAlp Dener             if (bnk->dnorm == tao->trust) {
977080d2917SAlp Dener               if (kappa < bnk->eta2) {
978080d2917SAlp Dener                 /* Marginal bad step */
979c133c014SAlp Dener                 tao->trust = bnk->alpha2 * tao->trust;
9803105154fSTodd Munson               } else if (kappa < bnk->eta3) {
981fed79b8eSAlp Dener                 /* Reasonable step */
982fed79b8eSAlp Dener                 tao->trust = bnk->alpha3 * tao->trust;
9833105154fSTodd Munson               } else if (kappa < bnk->eta4) {
984080d2917SAlp Dener                 /* Good step */
985c133c014SAlp Dener                 tao->trust = bnk->alpha4 * tao->trust;
9863105154fSTodd Munson               } else {
987080d2917SAlp Dener                 /* Very good step */
988c133c014SAlp Dener                 tao->trust = bnk->alpha5 * tao->trust;
989080d2917SAlp Dener               }
990c133c014SAlp Dener             }
991080d2917SAlp Dener           }
992080d2917SAlp Dener         }
993080d2917SAlp Dener       }
994080d2917SAlp Dener     } else {
995080d2917SAlp Dener       /*  Newton step was not good; reduce the radius */
996080d2917SAlp Dener       tao->trust = bnk->alpha1 * PetscMin(bnk->dnorm, tao->trust);
997080d2917SAlp Dener     }
998080d2917SAlp Dener     break;
999080d2917SAlp Dener 
1000080d2917SAlp Dener   default:
1001080d2917SAlp Dener     if (stepType == BNK_NEWTON) {
1002b1c2d0e3SAlp Dener       if (prered < 0.0) {
1003080d2917SAlp Dener         /*  The predicted reduction has the wrong sign.  This cannot */
1004080d2917SAlp Dener         /*  happen in infinite precision arithmetic.  Step should */
1005080d2917SAlp Dener         /*  be rejected! */
1006080d2917SAlp Dener         tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
1007080d2917SAlp Dener       } else {
1008b1c2d0e3SAlp Dener         if (PetscIsInfOrNanReal(actred)) {
1009080d2917SAlp Dener           tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
1010080d2917SAlp Dener         } else {
1011080d2917SAlp Dener           if ((PetscAbsScalar(actred) <= bnk->epsilon) && (PetscAbsScalar(prered) <= bnk->epsilon)) {
1012080d2917SAlp Dener             kappa = 1.0;
1013080d2917SAlp Dener           } else {
1014080d2917SAlp Dener             kappa = actred / prered;
1015080d2917SAlp Dener           }
1016080d2917SAlp Dener 
1017080d2917SAlp Dener           ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr);
1018080d2917SAlp Dener           tau_1 = bnk->theta * gdx / (bnk->theta * gdx - (1.0 - bnk->theta) * prered + actred);
1019080d2917SAlp Dener           tau_2 = bnk->theta * gdx / (bnk->theta * gdx + (1.0 + bnk->theta) * prered - actred);
1020080d2917SAlp Dener           tau_min = PetscMin(tau_1, tau_2);
1021080d2917SAlp Dener           tau_max = PetscMax(tau_1, tau_2);
1022080d2917SAlp Dener 
1023080d2917SAlp Dener           if (kappa >= 1.0 - bnk->mu1) {
1024080d2917SAlp Dener             /*  Great agreement */
1025080d2917SAlp Dener             *accept = PETSC_TRUE;
1026080d2917SAlp Dener             if (tau_max < 1.0) {
1027080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm);
1028080d2917SAlp Dener             } else if (tau_max > bnk->gamma4) {
1029080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, bnk->gamma4 * bnk->dnorm);
1030080d2917SAlp Dener             } else {
1031080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm);
1032080d2917SAlp Dener             }
1033080d2917SAlp Dener           } else if (kappa >= 1.0 - bnk->mu2) {
1034080d2917SAlp Dener             /*  Good agreement */
1035080d2917SAlp Dener             *accept = PETSC_TRUE;
1036080d2917SAlp Dener             if (tau_max < bnk->gamma2) {
1037080d2917SAlp Dener               tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm);
1038080d2917SAlp Dener             } else if (tau_max > bnk->gamma3) {
1039080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm);
1040080d2917SAlp Dener             } else if (tau_max < 1.0) {
1041080d2917SAlp Dener               tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm);
1042080d2917SAlp Dener             } else {
1043080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm);
1044080d2917SAlp Dener             }
1045080d2917SAlp Dener           } else {
1046080d2917SAlp Dener             /*  Not good agreement */
1047080d2917SAlp Dener             if (tau_min > 1.0) {
1048080d2917SAlp Dener               tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm);
1049080d2917SAlp Dener             } else if (tau_max < bnk->gamma1) {
1050080d2917SAlp Dener               tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
1051080d2917SAlp Dener             } else if ((tau_min < bnk->gamma1) && (tau_max >= 1.0)) {
1052080d2917SAlp Dener               tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
1053080d2917SAlp Dener             } else if ((tau_1 >= bnk->gamma1) && (tau_1 < 1.0) && ((tau_2 < bnk->gamma1) || (tau_2 >= 1.0))) {
1054080d2917SAlp Dener               tao->trust = tau_1 * PetscMin(tao->trust, bnk->dnorm);
1055080d2917SAlp Dener             } else if ((tau_2 >= bnk->gamma1) && (tau_2 < 1.0) && ((tau_1 < bnk->gamma1) || (tau_2 >= 1.0))) {
1056080d2917SAlp Dener               tao->trust = tau_2 * PetscMin(tao->trust, bnk->dnorm);
1057080d2917SAlp Dener             } else {
1058080d2917SAlp Dener               tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm);
1059080d2917SAlp Dener             }
1060080d2917SAlp Dener           }
1061080d2917SAlp Dener         }
1062080d2917SAlp Dener       }
1063080d2917SAlp Dener     } else {
1064080d2917SAlp Dener       /*  Newton step was not good; reduce the radius */
1065080d2917SAlp Dener       tao->trust = bnk->gamma1 * PetscMin(bnk->dnorm, tao->trust);
1066080d2917SAlp Dener     }
106728017e9fSAlp Dener     break;
1068080d2917SAlp Dener   }
1069c133c014SAlp Dener   /* Make sure the radius does not violate min and max settings */
1070c133c014SAlp Dener   tao->trust = PetscMin(tao->trust, bnk->max_radius);
1071fed79b8eSAlp Dener   tao->trust = PetscMax(tao->trust, bnk->min_radius);
1072080d2917SAlp Dener   PetscFunctionReturn(0);
1073080d2917SAlp Dener }
1074080d2917SAlp Dener 
1075eb910715SAlp Dener /* ---------------------------------------------------------- */
1076df278d8fSAlp Dener 
107762675beeSAlp Dener PetscErrorCode TaoBNKAddStepCounts(Tao tao, PetscInt stepType)
107862675beeSAlp Dener {
107962675beeSAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
108062675beeSAlp Dener 
108162675beeSAlp Dener   PetscFunctionBegin;
108262675beeSAlp Dener   switch (stepType) {
108362675beeSAlp Dener   case BNK_NEWTON:
108462675beeSAlp Dener     ++bnk->newt;
108562675beeSAlp Dener     break;
108662675beeSAlp Dener   case BNK_BFGS:
108762675beeSAlp Dener     ++bnk->bfgs;
108862675beeSAlp Dener     break;
108962675beeSAlp Dener   case BNK_SCALED_GRADIENT:
109062675beeSAlp Dener     ++bnk->sgrad;
109162675beeSAlp Dener     break;
109262675beeSAlp Dener   case BNK_GRADIENT:
109362675beeSAlp Dener     ++bnk->grad;
109462675beeSAlp Dener     break;
109562675beeSAlp Dener   default:
109662675beeSAlp Dener     break;
109762675beeSAlp Dener   }
109862675beeSAlp Dener   PetscFunctionReturn(0);
109962675beeSAlp Dener }
110062675beeSAlp Dener 
110162675beeSAlp Dener /* ---------------------------------------------------------- */
110262675beeSAlp Dener 
11039b6ef848SAlp Dener PetscErrorCode TaoSetUp_BNK(Tao tao)
1104eb910715SAlp Dener {
1105eb910715SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
1106eb910715SAlp Dener   PetscErrorCode ierr;
11079b6ef848SAlp Dener   KSPType        ksp_type;
1108e031d6f5SAlp Dener   PetscInt       i;
1109eb910715SAlp Dener 
1110eb910715SAlp Dener   PetscFunctionBegin;
1111c4b75bccSAlp Dener   if (!tao->gradient) {
1112c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);
1113c4b75bccSAlp Dener   }
1114c4b75bccSAlp Dener   if (!tao->stepdirection) {
1115c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr);
1116c4b75bccSAlp Dener   }
1117c4b75bccSAlp Dener   if (!bnk->W) {
1118c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->W);CHKERRQ(ierr);
1119c4b75bccSAlp Dener   }
1120c4b75bccSAlp Dener   if (!bnk->Xold) {
1121c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->Xold);CHKERRQ(ierr);
1122c4b75bccSAlp Dener   }
1123c4b75bccSAlp Dener   if (!bnk->Gold) {
1124c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->Gold);CHKERRQ(ierr);
1125c4b75bccSAlp Dener   }
1126c4b75bccSAlp Dener   if (!bnk->Xwork) {
1127c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->Xwork);CHKERRQ(ierr);
1128c4b75bccSAlp Dener   }
1129c4b75bccSAlp Dener   if (!bnk->Gwork) {
1130c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->Gwork);CHKERRQ(ierr);
1131c4b75bccSAlp Dener   }
1132c4b75bccSAlp Dener   if (!bnk->unprojected_gradient) {
1133c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->unprojected_gradient);CHKERRQ(ierr);
1134c4b75bccSAlp Dener   }
1135c4b75bccSAlp Dener   if (!bnk->unprojected_gradient_old) {
1136c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->unprojected_gradient_old);CHKERRQ(ierr);
1137c4b75bccSAlp Dener   }
1138c4b75bccSAlp Dener   if (!bnk->Diag_min) {
1139c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->Diag_min);CHKERRQ(ierr);
1140c4b75bccSAlp Dener   }
1141c4b75bccSAlp Dener   if (!bnk->Diag_max) {
1142c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->Diag_max);CHKERRQ(ierr);
1143c4b75bccSAlp Dener   }
1144e031d6f5SAlp Dener   if (bnk->max_cg_its > 0) {
1145c4b75bccSAlp Dener     /* Ensure that the important common vectors are shared between BNK and embedded BNCG */
1146c4b75bccSAlp Dener     bnk->bncg_ctx = (TAO_BNCG *)bnk->bncg->data;
114789da521bSAlp Dener     ierr = PetscObjectReference((PetscObject)(bnk->unprojected_gradient_old));CHKERRQ(ierr);
114889da521bSAlp Dener     ierr = VecDestroy(&bnk->bncg_ctx->unprojected_gradient_old);CHKERRQ(ierr);
114989da521bSAlp Dener     bnk->bncg_ctx->unprojected_gradient_old = bnk->unprojected_gradient_old;
1150c4b75bccSAlp Dener     ierr = PetscObjectReference((PetscObject)(bnk->unprojected_gradient));CHKERRQ(ierr);
1151c4b75bccSAlp Dener     ierr = VecDestroy(&bnk->bncg_ctx->unprojected_gradient);CHKERRQ(ierr);
1152c4b75bccSAlp Dener     bnk->bncg_ctx->unprojected_gradient = bnk->unprojected_gradient;
1153c4b75bccSAlp Dener     ierr = PetscObjectReference((PetscObject)(bnk->Gold));CHKERRQ(ierr);
1154c4b75bccSAlp Dener     ierr = VecDestroy(&bnk->bncg_ctx->G_old);CHKERRQ(ierr);
1155c4b75bccSAlp Dener     bnk->bncg_ctx->G_old = bnk->Gold;
1156c4b75bccSAlp Dener     ierr = PetscObjectReference((PetscObject)(tao->gradient));CHKERRQ(ierr);
1157c4b75bccSAlp Dener     ierr = VecDestroy(&bnk->bncg->gradient);CHKERRQ(ierr);
1158c4b75bccSAlp Dener     bnk->bncg->gradient = tao->gradient;
1159c4b75bccSAlp Dener     ierr = PetscObjectReference((PetscObject)(tao->stepdirection));CHKERRQ(ierr);
1160c4b75bccSAlp Dener     ierr = VecDestroy(&bnk->bncg->stepdirection);CHKERRQ(ierr);
1161c4b75bccSAlp Dener     bnk->bncg->stepdirection = tao->stepdirection;
1162c4b75bccSAlp Dener     ierr = TaoSetInitialVector(bnk->bncg, tao->solution);CHKERRQ(ierr);
1163c4b75bccSAlp Dener     /* Copy over some settings from BNK into BNCG */
1164e031d6f5SAlp Dener     ierr = TaoSetMaximumIterations(bnk->bncg, bnk->max_cg_its);CHKERRQ(ierr);
1165e031d6f5SAlp Dener     ierr = TaoSetTolerances(bnk->bncg, tao->gatol, tao->grtol, tao->gttol);CHKERRQ(ierr);
1166e031d6f5SAlp Dener     ierr = TaoSetFunctionLowerBound(bnk->bncg, tao->fmin);CHKERRQ(ierr);
1167937a31a1SAlp Dener     ierr = TaoSetConvergenceTest(bnk->bncg, tao->ops->convergencetest, tao->cnvP);CHKERRQ(ierr);
1168e031d6f5SAlp Dener     ierr = TaoSetObjectiveRoutine(bnk->bncg, tao->ops->computeobjective, tao->user_objP);CHKERRQ(ierr);
1169e031d6f5SAlp Dener     ierr = TaoSetGradientRoutine(bnk->bncg, tao->ops->computegradient, tao->user_gradP);CHKERRQ(ierr);
1170e031d6f5SAlp Dener     ierr = TaoSetObjectiveAndGradientRoutine(bnk->bncg, tao->ops->computeobjectiveandgradient, tao->user_objgradP);CHKERRQ(ierr);
1171e031d6f5SAlp Dener     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)tao, (PetscObject)(bnk->bncg));CHKERRQ(ierr);
1172c4b75bccSAlp Dener     for (i=0; i<tao->numbermonitors; ++i) {
1173e031d6f5SAlp Dener       ierr = TaoSetMonitor(bnk->bncg, tao->monitor[i], tao->monitorcontext[i], tao->monitordestroy[i]);CHKERRQ(ierr);
1174e031d6f5SAlp Dener       ierr = PetscObjectReference((PetscObject)(tao->monitorcontext[i]));CHKERRQ(ierr);
1175e031d6f5SAlp Dener     }
1176e031d6f5SAlp Dener   }
1177eb910715SAlp Dener   bnk->Diag = 0;
1178c0f10754SAlp Dener   bnk->Diag_red = 0;
1179c0f10754SAlp Dener   bnk->X_inactive = 0;
1180c0f10754SAlp Dener   bnk->G_inactive = 0;
118162675beeSAlp Dener   bnk->inactive_work = 0;
118262675beeSAlp Dener   bnk->active_work = 0;
118362675beeSAlp Dener   bnk->inactive_idx = 0;
118462675beeSAlp Dener   bnk->active_idx = 0;
118562675beeSAlp Dener   bnk->active_lower = 0;
118662675beeSAlp Dener   bnk->active_upper = 0;
118762675beeSAlp Dener   bnk->active_fixed = 0;
1188eb910715SAlp Dener   bnk->M = 0;
118909164190SAlp Dener   bnk->H_inactive = 0;
119009164190SAlp Dener   bnk->Hpre_inactive = 0;
11919b6ef848SAlp Dener   ierr = KSPGetType(tao->ksp,&ksp_type);CHKERRQ(ierr);
11929b6ef848SAlp Dener   ierr = PetscStrcmp(ksp_type,KSPCGNASH,&bnk->is_nash);CHKERRQ(ierr);
11939b6ef848SAlp Dener   ierr = PetscStrcmp(ksp_type,KSPCGSTCG,&bnk->is_stcg);CHKERRQ(ierr);
11949b6ef848SAlp Dener   ierr = PetscStrcmp(ksp_type,KSPCGGLTR,&bnk->is_gltr);CHKERRQ(ierr);
1195eb910715SAlp Dener   PetscFunctionReturn(0);
1196eb910715SAlp Dener }
1197eb910715SAlp Dener 
1198eb910715SAlp Dener /*------------------------------------------------------------*/
1199df278d8fSAlp Dener 
1200eb910715SAlp Dener static PetscErrorCode TaoDestroy_BNK(Tao tao)
1201eb910715SAlp Dener {
1202eb910715SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
1203eb910715SAlp Dener   PetscErrorCode ierr;
1204eb910715SAlp Dener 
1205eb910715SAlp Dener   PetscFunctionBegin;
1206eb910715SAlp Dener   if (tao->setupcalled) {
1207eb910715SAlp Dener     ierr = VecDestroy(&bnk->W);CHKERRQ(ierr);
1208eb910715SAlp Dener     ierr = VecDestroy(&bnk->Xold);CHKERRQ(ierr);
1209eb910715SAlp Dener     ierr = VecDestroy(&bnk->Gold);CHKERRQ(ierr);
121009164190SAlp Dener     ierr = VecDestroy(&bnk->Xwork);CHKERRQ(ierr);
121109164190SAlp Dener     ierr = VecDestroy(&bnk->Gwork);CHKERRQ(ierr);
121209164190SAlp Dener     ierr = VecDestroy(&bnk->unprojected_gradient);CHKERRQ(ierr);
121309164190SAlp Dener     ierr = VecDestroy(&bnk->unprojected_gradient_old);CHKERRQ(ierr);
121462675beeSAlp Dener     ierr = VecDestroy(&bnk->Diag_min);CHKERRQ(ierr);
121562675beeSAlp Dener     ierr = VecDestroy(&bnk->Diag_max);CHKERRQ(ierr);
1216c4b75bccSAlp Dener   }
1217ca964c49SAlp Dener   ierr = ISDestroy(&bnk->active_lower);CHKERRQ(ierr);
1218ca964c49SAlp Dener   ierr = ISDestroy(&bnk->active_upper);CHKERRQ(ierr);
1219ca964c49SAlp Dener   ierr = ISDestroy(&bnk->active_fixed);CHKERRQ(ierr);
1220ca964c49SAlp Dener   ierr = ISDestroy(&bnk->active_idx);CHKERRQ(ierr);
1221ca964c49SAlp Dener   ierr = ISDestroy(&bnk->inactive_idx);CHKERRQ(ierr);
12225e9b73cbSAlp Dener   ierr = VecDestroy(&bnk->Diag);CHKERRQ(ierr);
1223eb910715SAlp Dener   ierr = MatDestroy(&bnk->M);CHKERRQ(ierr);
1224c4b75bccSAlp Dener   if (bnk->Hpre_inactive != tao->hessian_pre && bnk->Hpre_inactive != bnk->H_inactive) {
1225c4b75bccSAlp Dener     ierr = MatDestroy(&bnk->Hpre_inactive);CHKERRQ(ierr);
1226c4b75bccSAlp Dener   }
1227c4b75bccSAlp Dener   if (bnk->H_inactive != tao->hessian) {
1228c4b75bccSAlp Dener     ierr = MatDestroy(&bnk->H_inactive);CHKERRQ(ierr);
1229c4b75bccSAlp Dener   }
1230ca964c49SAlp Dener   if (bnk->max_cg_its > 0) {
1231ca964c49SAlp Dener     ierr = TaoDestroy(&bnk->bncg);CHKERRQ(ierr);
1232ca964c49SAlp Dener   }
1233eb910715SAlp Dener   ierr = PetscFree(tao->data);CHKERRQ(ierr);
1234eb910715SAlp Dener   PetscFunctionReturn(0);
1235eb910715SAlp Dener }
1236eb910715SAlp Dener 
1237eb910715SAlp Dener /*------------------------------------------------------------*/
1238df278d8fSAlp Dener 
1239eb910715SAlp Dener static PetscErrorCode TaoSetFromOptions_BNK(PetscOptionItems *PetscOptionsObject,Tao tao)
1240eb910715SAlp Dener {
1241eb910715SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
1242eb910715SAlp Dener   PetscErrorCode ierr;
1243eb910715SAlp Dener 
1244eb910715SAlp Dener   PetscFunctionBegin;
1245eb910715SAlp Dener   ierr = PetscOptionsHead(PetscOptionsObject,"Newton line search method for unconstrained optimization");CHKERRQ(ierr);
12468d5ead36SAlp 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);
12478d5ead36SAlp 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);
12488d5ead36SAlp 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);
12498d5ead36SAlp 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);
12502f75a4aaSAlp 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);
12518d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_sval", "perturbation starting value", "", bnk->sval, &bnk->sval,NULL);CHKERRQ(ierr);
12528d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_imin", "minimum initial perturbation", "", bnk->imin, &bnk->imin,NULL);CHKERRQ(ierr);
12538d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_imax", "maximum initial perturbation", "", bnk->imax, &bnk->imax,NULL);CHKERRQ(ierr);
12548d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_imfac", "initial merit factor", "", bnk->imfac, &bnk->imfac,NULL);CHKERRQ(ierr);
12558d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_pmin", "minimum perturbation", "", bnk->pmin, &bnk->pmin,NULL);CHKERRQ(ierr);
12568d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_pmax", "maximum perturbation", "", bnk->pmax, &bnk->pmax,NULL);CHKERRQ(ierr);
12578d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_pgfac", "growth factor", "", bnk->pgfac, &bnk->pgfac,NULL);CHKERRQ(ierr);
12588d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_psfac", "shrink factor", "", bnk->psfac, &bnk->psfac,NULL);CHKERRQ(ierr);
12598d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_pmgfac", "merit growth factor", "", bnk->pmgfac, &bnk->pmgfac,NULL);CHKERRQ(ierr);
12608d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_pmsfac", "merit shrink factor", "", bnk->pmsfac, &bnk->pmsfac,NULL);CHKERRQ(ierr);
12618d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_eta1", "poor steplength; reduce radius", "", bnk->eta1, &bnk->eta1,NULL);CHKERRQ(ierr);
12628d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_eta2", "reasonable steplength; leave radius alone", "", bnk->eta2, &bnk->eta2,NULL);CHKERRQ(ierr);
12638d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_eta3", "good steplength; increase radius", "", bnk->eta3, &bnk->eta3,NULL);CHKERRQ(ierr);
12648d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_eta4", "excellent steplength; greatly increase radius", "", bnk->eta4, &bnk->eta4,NULL);CHKERRQ(ierr);
12658d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_alpha1", "", "", bnk->alpha1, &bnk->alpha1,NULL);CHKERRQ(ierr);
12668d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_alpha2", "", "", bnk->alpha2, &bnk->alpha2,NULL);CHKERRQ(ierr);
12678d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_alpha3", "", "", bnk->alpha3, &bnk->alpha3,NULL);CHKERRQ(ierr);
12688d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_alpha4", "", "", bnk->alpha4, &bnk->alpha4,NULL);CHKERRQ(ierr);
12698d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_alpha5", "", "", bnk->alpha5, &bnk->alpha5,NULL);CHKERRQ(ierr);
12708d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_nu1", "poor steplength; reduce radius", "", bnk->nu1, &bnk->nu1,NULL);CHKERRQ(ierr);
12718d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_nu2", "reasonable steplength; leave radius alone", "", bnk->nu2, &bnk->nu2,NULL);CHKERRQ(ierr);
12728d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_nu3", "good steplength; increase radius", "", bnk->nu3, &bnk->nu3,NULL);CHKERRQ(ierr);
12738d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_nu4", "excellent steplength; greatly increase radius", "", bnk->nu4, &bnk->nu4,NULL);CHKERRQ(ierr);
12748d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_omega1", "", "", bnk->omega1, &bnk->omega1,NULL);CHKERRQ(ierr);
12758d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_omega2", "", "", bnk->omega2, &bnk->omega2,NULL);CHKERRQ(ierr);
12768d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_omega3", "", "", bnk->omega3, &bnk->omega3,NULL);CHKERRQ(ierr);
12778d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_omega4", "", "", bnk->omega4, &bnk->omega4,NULL);CHKERRQ(ierr);
12788d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_omega5", "", "", bnk->omega5, &bnk->omega5,NULL);CHKERRQ(ierr);
12798d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_mu1_i", "", "", bnk->mu1_i, &bnk->mu1_i,NULL);CHKERRQ(ierr);
12808d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_mu2_i", "", "", bnk->mu2_i, &bnk->mu2_i,NULL);CHKERRQ(ierr);
12818d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_gamma1_i", "", "", bnk->gamma1_i, &bnk->gamma1_i,NULL);CHKERRQ(ierr);
12828d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_gamma2_i", "", "", bnk->gamma2_i, &bnk->gamma2_i,NULL);CHKERRQ(ierr);
12838d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_gamma3_i", "", "", bnk->gamma3_i, &bnk->gamma3_i,NULL);CHKERRQ(ierr);
12848d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_gamma4_i", "", "", bnk->gamma4_i, &bnk->gamma4_i,NULL);CHKERRQ(ierr);
12858d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_theta_i", "", "", bnk->theta_i, &bnk->theta_i,NULL);CHKERRQ(ierr);
12868d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_mu1", "", "", bnk->mu1, &bnk->mu1,NULL);CHKERRQ(ierr);
12878d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_mu2", "", "", bnk->mu2, &bnk->mu2,NULL);CHKERRQ(ierr);
12888d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_gamma1", "", "", bnk->gamma1, &bnk->gamma1,NULL);CHKERRQ(ierr);
12898d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_gamma2", "", "", bnk->gamma2, &bnk->gamma2,NULL);CHKERRQ(ierr);
12908d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_gamma3", "", "", bnk->gamma3, &bnk->gamma3,NULL);CHKERRQ(ierr);
12918d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_gamma4", "", "", bnk->gamma4, &bnk->gamma4,NULL);CHKERRQ(ierr);
12928d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_theta", "", "", bnk->theta, &bnk->theta,NULL);CHKERRQ(ierr);
12938d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_min_radius", "lower bound on initial radius", "", bnk->min_radius, &bnk->min_radius,NULL);CHKERRQ(ierr);
12948d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_max_radius", "upper bound on radius", "", bnk->max_radius, &bnk->max_radius,NULL);CHKERRQ(ierr);
12958d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_epsilon", "tolerance used when computing actual and predicted reduction", "", bnk->epsilon, &bnk->epsilon,NULL);CHKERRQ(ierr);
12960a4511e9SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_as_tol", "initial tolerance used when estimating actively bounded variables", "", bnk->as_tol, &bnk->as_tol,NULL);CHKERRQ(ierr);
12970a4511e9SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_as_step", "step length used when estimating actively bounded variables", "", bnk->as_step, &bnk->as_step,NULL);CHKERRQ(ierr);
1298c0f10754SAlp 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);
1299eb910715SAlp Dener   ierr = PetscOptionsTail();CHKERRQ(ierr);
1300e031d6f5SAlp Dener   ierr = TaoSetFromOptions(bnk->bncg);
1301eb910715SAlp Dener   ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr);
1302eb910715SAlp Dener   ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
1303eb910715SAlp Dener   PetscFunctionReturn(0);
1304eb910715SAlp Dener }
1305eb910715SAlp Dener 
1306eb910715SAlp Dener /*------------------------------------------------------------*/
1307df278d8fSAlp Dener 
1308eb910715SAlp Dener static PetscErrorCode TaoView_BNK(Tao tao, PetscViewer viewer)
1309eb910715SAlp Dener {
1310eb910715SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
1311eb910715SAlp Dener   PetscInt       nrejects;
1312eb910715SAlp Dener   PetscBool      isascii;
1313eb910715SAlp Dener   PetscErrorCode ierr;
1314eb910715SAlp Dener 
1315eb910715SAlp Dener   PetscFunctionBegin;
1316eb910715SAlp Dener   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
1317eb910715SAlp Dener   if (isascii) {
1318eb910715SAlp Dener     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1319eb910715SAlp Dener     if (BNK_PC_BFGS == bnk->pc_type && bnk->M) {
1320eb910715SAlp Dener       ierr = MatLMVMGetRejects(bnk->M,&nrejects);CHKERRQ(ierr);
1321eb910715SAlp Dener       ierr = PetscViewerASCIIPrintf(viewer, "Rejected matrix updates: %D\n",nrejects);CHKERRQ(ierr);
1322eb910715SAlp Dener     }
1323e031d6f5SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "CG steps: %D\n", bnk->tot_cg_its);CHKERRQ(ierr);
1324eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "Newton steps: %D\n", bnk->newt);CHKERRQ(ierr);
1325eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "BFGS steps: %D\n", bnk->bfgs);CHKERRQ(ierr);
1326eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "Scaled gradient steps: %D\n", bnk->sgrad);CHKERRQ(ierr);
1327eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", bnk->grad);CHKERRQ(ierr);
1328eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "KSP termination reasons:\n");CHKERRQ(ierr);
1329eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  atol: %D\n", bnk->ksp_atol);CHKERRQ(ierr);
1330eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  rtol: %D\n", bnk->ksp_rtol);CHKERRQ(ierr);
1331eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  ctol: %D\n", bnk->ksp_ctol);CHKERRQ(ierr);
1332eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  negc: %D\n", bnk->ksp_negc);CHKERRQ(ierr);
1333eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  dtol: %D\n", bnk->ksp_dtol);CHKERRQ(ierr);
1334eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  iter: %D\n", bnk->ksp_iter);CHKERRQ(ierr);
1335eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  othr: %D\n", bnk->ksp_othr);CHKERRQ(ierr);
1336eb910715SAlp Dener     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1337eb910715SAlp Dener   }
1338eb910715SAlp Dener   PetscFunctionReturn(0);
1339eb910715SAlp Dener }
1340eb910715SAlp Dener 
1341eb910715SAlp Dener /* ---------------------------------------------------------- */
1342df278d8fSAlp Dener 
1343eb910715SAlp Dener /*MC
1344eb910715SAlp Dener   TAOBNK - Shared base-type for Bounded Newton-Krylov type algorithms.
134566ed3702SAlp Dener   At each iteration, the BNK methods solve the symmetric
1346eb910715SAlp Dener   system of equations to obtain the step diretion dk:
1347eb910715SAlp Dener               Hk dk = -gk
13482b97c8d8SAlp Dener   for free variables only. The step can be globalized either through
13492b97c8d8SAlp Dener   trust-region methods, or a line search, or a heuristic mixture of both.
1350eb910715SAlp Dener 
1351eb910715SAlp Dener     Options Database Keys:
13528d5ead36SAlp Dener + -tao_bnk_pc_type - "none","ahess","bfgs","petsc"
13538d5ead36SAlp Dener . -tao_bnk_bfgs_scale_type - "ahess","phess","bfgs"
13548d5ead36SAlp Dener . -tao_bnk_init_type - "constant","direction","interpolation"
13558d5ead36SAlp Dener . -tao_bnk_update_type - "step","direction","interpolation"
13562f75a4aaSAlp Dener . -tao_bnk_as_type - "none","bertsekas"
13578d5ead36SAlp Dener . -tao_bnk_sval - perturbation starting value
13588d5ead36SAlp Dener . -tao_bnk_imin - minimum initial perturbation
13598d5ead36SAlp Dener . -tao_bnk_imax - maximum initial perturbation
13608d5ead36SAlp Dener . -tao_bnk_pmin - minimum perturbation
13618d5ead36SAlp Dener . -tao_bnk_pmax - maximum perturbation
13628d5ead36SAlp Dener . -tao_bnk_pgfac - growth factor
13638d5ead36SAlp Dener . -tao_bnk_psfac - shrink factor
13648d5ead36SAlp Dener . -tao_bnk_imfac - initial merit factor
13658d5ead36SAlp Dener . -tao_bnk_pmgfac - merit growth factor
13668d5ead36SAlp Dener . -tao_bnk_pmsfac - merit shrink factor
13678d5ead36SAlp Dener . -tao_bnk_eta1 - poor steplength; reduce radius
13688d5ead36SAlp Dener . -tao_bnk_eta2 - reasonable steplength; leave radius
13698d5ead36SAlp Dener . -tao_bnk_eta3 - good steplength; increase readius
13708d5ead36SAlp Dener . -tao_bnk_eta4 - excellent steplength; greatly increase radius
13718d5ead36SAlp Dener . -tao_bnk_alpha1 - alpha1 reduction
13728d5ead36SAlp Dener . -tao_bnk_alpha2 - alpha2 reduction
13738d5ead36SAlp Dener . -tao_bnk_alpha3 - alpha3 reduction
13748d5ead36SAlp Dener . -tao_bnk_alpha4 - alpha4 reduction
13758d5ead36SAlp Dener . -tao_bnk_alpha - alpha5 reduction
13768d5ead36SAlp Dener . -tao_bnk_mu1 - mu1 interpolation update
13778d5ead36SAlp Dener . -tao_bnk_mu2 - mu2 interpolation update
13788d5ead36SAlp Dener . -tao_bnk_gamma1 - gamma1 interpolation update
13798d5ead36SAlp Dener . -tao_bnk_gamma2 - gamma2 interpolation update
13808d5ead36SAlp Dener . -tao_bnk_gamma3 - gamma3 interpolation update
13818d5ead36SAlp Dener . -tao_bnk_gamma4 - gamma4 interpolation update
13828d5ead36SAlp Dener . -tao_bnk_theta - theta interpolation update
13838d5ead36SAlp Dener . -tao_bnk_omega1 - omega1 step update
13848d5ead36SAlp Dener . -tao_bnk_omega2 - omega2 step update
13858d5ead36SAlp Dener . -tao_bnk_omega3 - omega3 step update
13868d5ead36SAlp Dener . -tao_bnk_omega4 - omega4 step update
13878d5ead36SAlp Dener . -tao_bnk_omega5 - omega5 step update
13888d5ead36SAlp Dener . -tao_bnk_mu1_i -  mu1 interpolation init factor
13898d5ead36SAlp Dener . -tao_bnk_mu2_i -  mu2 interpolation init factor
13908d5ead36SAlp Dener . -tao_bnk_gamma1_i -  gamma1 interpolation init factor
13918d5ead36SAlp Dener . -tao_bnk_gamma2_i -  gamma2 interpolation init factor
13928d5ead36SAlp Dener . -tao_bnk_gamma3_i -  gamma3 interpolation init factor
13938d5ead36SAlp Dener . -tao_bnk_gamma4_i -  gamma4 interpolation init factor
13942f75a4aaSAlp Dener . -tao_bnk_theta_i -  theta interpolation init factor
13952f75a4aaSAlp Dener - -tao_bnk_bound_tol -  initial tolerance used in estimating bounded active variables
1396eb910715SAlp Dener 
1397eb910715SAlp Dener   Level: beginner
1398eb910715SAlp Dener M*/
1399eb910715SAlp Dener 
1400eb910715SAlp Dener PetscErrorCode TaoCreate_BNK(Tao tao)
1401eb910715SAlp Dener {
1402eb910715SAlp Dener   TAO_BNK        *bnk;
1403eb910715SAlp Dener   const char     *morethuente_type = TAOLINESEARCHMT;
1404eb910715SAlp Dener   PetscErrorCode ierr;
1405eb910715SAlp Dener 
1406eb910715SAlp Dener   PetscFunctionBegin;
1407eb910715SAlp Dener   ierr = PetscNewLog(tao,&bnk);CHKERRQ(ierr);
1408eb910715SAlp Dener 
1409eb910715SAlp Dener   tao->ops->setup = TaoSetUp_BNK;
1410eb910715SAlp Dener   tao->ops->view = TaoView_BNK;
1411eb910715SAlp Dener   tao->ops->setfromoptions = TaoSetFromOptions_BNK;
1412eb910715SAlp Dener   tao->ops->destroy = TaoDestroy_BNK;
1413eb910715SAlp Dener 
1414eb910715SAlp Dener   /*  Override default settings (unless already changed) */
1415eb910715SAlp Dener   if (!tao->max_it_changed) tao->max_it = 50;
1416eb910715SAlp Dener   if (!tao->trust0_changed) tao->trust0 = 100.0;
1417eb910715SAlp Dener 
1418eb910715SAlp Dener   tao->data = (void*)bnk;
1419eb910715SAlp Dener 
142066ed3702SAlp Dener   /*  Hessian shifting parameters */
1421eb910715SAlp Dener   bnk->sval   = 0.0;
1422eb910715SAlp Dener   bnk->imin   = 1.0e-4;
1423eb910715SAlp Dener   bnk->imax   = 1.0e+2;
1424eb910715SAlp Dener   bnk->imfac  = 1.0e-1;
1425eb910715SAlp Dener 
1426eb910715SAlp Dener   bnk->pmin   = 1.0e-12;
1427eb910715SAlp Dener   bnk->pmax   = 1.0e+2;
1428eb910715SAlp Dener   bnk->pgfac  = 1.0e+1;
1429eb910715SAlp Dener   bnk->psfac  = 4.0e-1;
1430eb910715SAlp Dener   bnk->pmgfac = 1.0e-1;
1431eb910715SAlp Dener   bnk->pmsfac = 1.0e-1;
1432eb910715SAlp Dener 
1433eb910715SAlp Dener   /*  Default values for trust-region radius update based on steplength */
1434eb910715SAlp Dener   bnk->nu1 = 0.25;
1435eb910715SAlp Dener   bnk->nu2 = 0.50;
1436eb910715SAlp Dener   bnk->nu3 = 1.00;
1437eb910715SAlp Dener   bnk->nu4 = 1.25;
1438eb910715SAlp Dener 
1439eb910715SAlp Dener   bnk->omega1 = 0.25;
1440eb910715SAlp Dener   bnk->omega2 = 0.50;
1441eb910715SAlp Dener   bnk->omega3 = 1.00;
1442eb910715SAlp Dener   bnk->omega4 = 2.00;
1443eb910715SAlp Dener   bnk->omega5 = 4.00;
1444eb910715SAlp Dener 
1445eb910715SAlp Dener   /*  Default values for trust-region radius update based on reduction */
1446eb910715SAlp Dener   bnk->eta1 = 1.0e-4;
1447eb910715SAlp Dener   bnk->eta2 = 0.25;
1448eb910715SAlp Dener   bnk->eta3 = 0.50;
1449eb910715SAlp Dener   bnk->eta4 = 0.90;
1450eb910715SAlp Dener 
1451eb910715SAlp Dener   bnk->alpha1 = 0.25;
1452eb910715SAlp Dener   bnk->alpha2 = 0.50;
1453eb910715SAlp Dener   bnk->alpha3 = 1.00;
1454eb910715SAlp Dener   bnk->alpha4 = 2.00;
1455eb910715SAlp Dener   bnk->alpha5 = 4.00;
1456eb910715SAlp Dener 
1457eb910715SAlp Dener   /*  Default values for trust-region radius update based on interpolation */
1458eb910715SAlp Dener   bnk->mu1 = 0.10;
1459eb910715SAlp Dener   bnk->mu2 = 0.50;
1460eb910715SAlp Dener 
1461eb910715SAlp Dener   bnk->gamma1 = 0.25;
1462eb910715SAlp Dener   bnk->gamma2 = 0.50;
1463eb910715SAlp Dener   bnk->gamma3 = 2.00;
1464eb910715SAlp Dener   bnk->gamma4 = 4.00;
1465eb910715SAlp Dener 
1466eb910715SAlp Dener   bnk->theta = 0.05;
1467eb910715SAlp Dener 
1468eb910715SAlp Dener   /*  Default values for trust region initialization based on interpolation */
1469eb910715SAlp Dener   bnk->mu1_i = 0.35;
1470eb910715SAlp Dener   bnk->mu2_i = 0.50;
1471eb910715SAlp Dener 
1472eb910715SAlp Dener   bnk->gamma1_i = 0.0625;
1473eb910715SAlp Dener   bnk->gamma2_i = 0.5;
1474eb910715SAlp Dener   bnk->gamma3_i = 2.0;
1475eb910715SAlp Dener   bnk->gamma4_i = 5.0;
1476eb910715SAlp Dener 
1477eb910715SAlp Dener   bnk->theta_i = 0.25;
1478eb910715SAlp Dener 
1479eb910715SAlp Dener   /*  Remaining parameters */
1480c0f10754SAlp Dener   bnk->max_cg_its = 0;
1481eb910715SAlp Dener   bnk->min_radius = 1.0e-10;
1482eb910715SAlp Dener   bnk->max_radius = 1.0e10;
1483770b7498SAlp Dener   bnk->epsilon = PetscPowReal(PETSC_MACHINE_EPSILON, 2.0/3.0);
14840a4511e9SAlp Dener   bnk->as_tol = 1.0e-3;
14850a4511e9SAlp Dener   bnk->as_step = 1.0e-3;
148662675beeSAlp Dener   bnk->dmin = 1.0e-6;
148762675beeSAlp Dener   bnk->dmax = 1.0e6;
1488eb910715SAlp Dener 
1489e031d6f5SAlp Dener   bnk->pc_type         = BNK_PC_AHESS;
1490eb910715SAlp Dener   bnk->bfgs_scale_type = BFGS_SCALE_PHESS;
1491eb910715SAlp Dener   bnk->init_type       = BNK_INIT_INTERPOLATION;
1492fed79b8eSAlp Dener   bnk->update_type     = BNK_UPDATE_INTERPOLATION;
14932f75a4aaSAlp Dener   bnk->as_type         = BNK_AS_BERTSEKAS;
1494eb910715SAlp Dener 
1495e031d6f5SAlp Dener   /* Create the embedded BNCG solver */
1496e031d6f5SAlp Dener   ierr = TaoCreate(PetscObjectComm((PetscObject)tao), &bnk->bncg);CHKERRQ(ierr);
1497e031d6f5SAlp Dener   ierr = PetscObjectIncrementTabLevel((PetscObject)bnk->bncg, (PetscObject)tao, 1);CHKERRQ(ierr);
1498e031d6f5SAlp Dener   ierr = TaoSetOptionsPrefix(bnk->bncg, "tao_bnk_");CHKERRQ(ierr);
1499e031d6f5SAlp Dener   ierr = TaoSetType(bnk->bncg, TAOBNCG);CHKERRQ(ierr);
1500e031d6f5SAlp Dener 
1501c0f10754SAlp Dener   /* Create the line search */
1502eb910715SAlp Dener   ierr = TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch);CHKERRQ(ierr);
1503eb910715SAlp Dener   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr);
1504e031d6f5SAlp Dener   ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr);
1505eb910715SAlp Dener   ierr = TaoLineSearchSetType(tao->linesearch,morethuente_type);CHKERRQ(ierr);
1506eb910715SAlp Dener   ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr);
1507eb910715SAlp Dener 
1508eb910715SAlp Dener   /*  Set linear solver to default for symmetric matrices */
1509eb910715SAlp Dener   ierr = KSPCreate(((PetscObject)tao)->comm,&tao->ksp);CHKERRQ(ierr);
1510eb910715SAlp Dener   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1);CHKERRQ(ierr);
1511eb910715SAlp Dener   ierr = KSPSetOptionsPrefix(tao->ksp,tao->hdr.prefix);CHKERRQ(ierr);
1512eb910715SAlp Dener   ierr = KSPSetType(tao->ksp,KSPCGSTCG);CHKERRQ(ierr);
1513eb910715SAlp Dener   PetscFunctionReturn(0);
1514eb910715SAlp Dener }
1515