xref: /petsc/src/tao/bound/impls/bnk/bnk.c (revision 3105154f8b5773d8bc272595da50e6ebbb4dfd4f)
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 
400a4511e9SAlp Dener   PetscReal                    f_min, ftrial, prered, actred, kappa, sigma;
41eb910715SAlp Dener   PetscReal                    tau, tau_1, tau_2, tau_max, tau_min, max_radius;
42e031d6f5SAlp Dener   PetscReal                    resnorm, 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 */
57c4b75bccSAlp Dener   ierr = TaoBoundSolution(tao->XL, tao->XU, tao->solution, &nDiff);CHKERRQ(ierr);
5828017e9fSAlp Dener 
5928017e9fSAlp Dener   /* Check convergence criteria */
6028017e9fSAlp Dener   ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &bnk->f, bnk->unprojected_gradient);CHKERRQ(ierr);
6161be54a6SAlp Dener   ierr = TaoBNKEstimateActiveSet(tao, bnk->as_type);CHKERRQ(ierr);
6261be54a6SAlp Dener   ierr = VecCopy(bnk->unprojected_gradient, tao->gradient);CHKERRQ(ierr);
637529f6b4SAlp Dener   if (bnk->active_idx) {
6461be54a6SAlp Dener     ierr = VecISSet(tao->gradient, bnk->active_idx, 0.0);CHKERRQ(ierr);
657529f6b4SAlp Dener   }
6608752603SAlp Dener   ierr = VecNorm(tao->gradient,NORM_2,&bnk->gnorm);CHKERRQ(ierr);
6728017e9fSAlp Dener   if (PetscIsInfOrNanReal(bnk->f) || PetscIsInfOrNanReal(bnk->gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
6828017e9fSAlp Dener 
69c0f10754SAlp Dener   /* Test the initial point for convergence */
70e031d6f5SAlp Dener   ierr = VecFischer(tao->solution, bnk->unprojected_gradient, tao->XL, tao->XU, bnk->Gwork);CHKERRQ(ierr);
71e031d6f5SAlp Dener   ierr = VecNorm(bnk->Gwork, NORM_2, &resnorm);CHKERRQ(ierr);
72e031d6f5SAlp Dener   ierr = TaoLogConvergenceHistory(tao,bnk->f,resnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
73e031d6f5SAlp Dener   ierr = TaoMonitor(tao,tao->niter,bnk->f,resnorm,0.0,1.0);CHKERRQ(ierr);
74c0f10754SAlp Dener   ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
75c0f10754SAlp Dener   if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
76c0f10754SAlp Dener 
77e031d6f5SAlp Dener   /* Reset KSP stopping reason counters */
78eb910715SAlp Dener   bnk->ksp_atol = 0;
79eb910715SAlp Dener   bnk->ksp_rtol = 0;
80eb910715SAlp Dener   bnk->ksp_dtol = 0;
81eb910715SAlp Dener   bnk->ksp_ctol = 0;
82eb910715SAlp Dener   bnk->ksp_negc = 0;
83eb910715SAlp Dener   bnk->ksp_iter = 0;
84eb910715SAlp Dener   bnk->ksp_othr = 0;
85eb910715SAlp Dener 
86e031d6f5SAlp Dener   /* Reset accepted step type counters */
87e031d6f5SAlp Dener   bnk->tot_cg_its = 0;
88e031d6f5SAlp Dener   bnk->newt = 0;
89e031d6f5SAlp Dener   bnk->bfgs = 0;
90e031d6f5SAlp Dener   bnk->sgrad = 0;
91e031d6f5SAlp Dener   bnk->grad = 0;
92e031d6f5SAlp Dener 
93fed79b8eSAlp Dener   /* Initialize the Hessian perturbation */
94fed79b8eSAlp Dener   bnk->pert = bnk->sval;
95fed79b8eSAlp Dener 
96937a31a1SAlp Dener   /* Reset initial steplength to zero (this helps BNCG reset its direction internally) */
97937a31a1SAlp Dener   ierr = VecSet(tao->stepdirection, 0.0);CHKERRQ(ierr);
98937a31a1SAlp Dener 
99e031d6f5SAlp Dener   /* Allocate the vectors needed for the BFGS approximation */
100e031d6f5SAlp Dener   if (BNK_PC_BFGS == bnk->pc_type) {
101e031d6f5SAlp Dener     if (!bnk->M) {
102eb910715SAlp Dener       ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr);
103eb910715SAlp Dener       ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr);
104eb910715SAlp Dener       ierr = MatCreateLMVM(((PetscObject)tao)->comm,n,N,&bnk->M);CHKERRQ(ierr);
105eb910715SAlp Dener       ierr = MatLMVMAllocateVectors(bnk->M,tao->solution);CHKERRQ(ierr);
106eb910715SAlp Dener     }
107e031d6f5SAlp Dener     if (bnk->bfgs_scale_type != BFGS_SCALE_BFGS && !bnk->Diag) {
108e031d6f5SAlp Dener       ierr = VecDuplicate(tao->solution,&bnk->Diag);CHKERRQ(ierr);
1095e9b73cbSAlp Dener     }
110e031d6f5SAlp Dener   }
111e031d6f5SAlp Dener 
112e031d6f5SAlp Dener   /* Prepare the min/max vectors for safeguarding diagonal scales */
11362675beeSAlp Dener   ierr = VecSet(bnk->Diag_min, bnk->dmin);CHKERRQ(ierr);
11462675beeSAlp Dener   ierr = VecSet(bnk->Diag_max, bnk->dmax);CHKERRQ(ierr);
115eb910715SAlp Dener 
116eb910715SAlp Dener   /* Modify the preconditioner to use the bfgs approximation */
117eb910715SAlp Dener   ierr = KSPGetPC(tao->ksp, &pc);CHKERRQ(ierr);
118eb910715SAlp Dener   switch(bnk->pc_type) {
119eb910715SAlp Dener   case BNK_PC_NONE:
120eb910715SAlp Dener     ierr = PCSetType(pc, PCNONE);CHKERRQ(ierr);
121eb910715SAlp Dener     ierr = PCSetFromOptions(pc);CHKERRQ(ierr);
122eb910715SAlp Dener     break;
123eb910715SAlp Dener 
124eb910715SAlp Dener   case BNK_PC_AHESS:
125eb910715SAlp Dener     ierr = PCSetType(pc, PCJACOBI);CHKERRQ(ierr);
126eb910715SAlp Dener     ierr = PCSetFromOptions(pc);CHKERRQ(ierr);
127eb910715SAlp Dener     ierr = PCJacobiSetUseAbs(pc,PETSC_TRUE);CHKERRQ(ierr);
128eb910715SAlp Dener     break;
129eb910715SAlp Dener 
130eb910715SAlp Dener   case BNK_PC_BFGS:
131eb910715SAlp Dener     ierr = PCSetType(pc, PCSHELL);CHKERRQ(ierr);
132eb910715SAlp Dener     ierr = PCSetFromOptions(pc);CHKERRQ(ierr);
133eb910715SAlp Dener     ierr = PCShellSetName(pc, "bfgs");CHKERRQ(ierr);
134eb910715SAlp Dener     ierr = PCShellSetContext(pc, bnk->M);CHKERRQ(ierr);
135eb910715SAlp Dener     ierr = PCShellSetApply(pc, MatLMVMSolveShell);CHKERRQ(ierr);
136eb910715SAlp Dener     break;
137eb910715SAlp Dener 
138eb910715SAlp Dener   default:
139eb910715SAlp Dener     /* Use the pc method set by pc_type */
140eb910715SAlp Dener     break;
141eb910715SAlp Dener   }
142eb910715SAlp Dener 
143eb910715SAlp Dener   /* Initialize trust-region radius.  The initialization is only performed
144eb910715SAlp Dener      when we are using Nash, Steihaug-Toint or the Generalized Lanczos method. */
145c0f10754SAlp Dener   *needH = PETSC_TRUE;
146eb910715SAlp Dener   if (bnk->is_nash || bnk->is_stcg || bnk->is_gltr) {
14762675beeSAlp Dener     switch(initType) {
148eb910715SAlp Dener     case BNK_INIT_CONSTANT:
149eb910715SAlp Dener       /* Use the initial radius specified */
150c0f10754SAlp Dener       tao->trust = tao->trust0;
151eb910715SAlp Dener       break;
152eb910715SAlp Dener 
153eb910715SAlp Dener     case BNK_INIT_INTERPOLATION:
154c0f10754SAlp Dener       /* Use interpolation based on the initial Hessian */
155eb910715SAlp Dener       max_radius = 0.0;
15608752603SAlp Dener       tao->trust = tao->trust0;
157eb910715SAlp Dener       for (j = 0; j < j_max; ++j) {
1580a4511e9SAlp Dener         f_min = bnk->f;
159eb910715SAlp Dener         sigma = 0.0;
160eb910715SAlp Dener 
161c0f10754SAlp Dener         if (*needH) {
16262602cfbSAlp Dener           /* Compute the Hessian at the new step, and extract the inactive subsystem */
163e031d6f5SAlp Dener           ierr = TaoBNKComputeHessian(tao);CHKERRQ(ierr);
16408752603SAlp Dener           ierr = TaoBNKEstimateActiveSet(tao, BNK_AS_NONE);CHKERRQ(ierr);
16528017e9fSAlp Dener           if (bnk->inactive_idx) {
1662ab2a32cSAlp Dener             ierr = MatCreateSubMatrix(tao->hessian, bnk->inactive_idx, bnk->inactive_idx, MAT_INITIAL_MATRIX, &bnk->H_inactive);CHKERRQ(ierr);
16728017e9fSAlp Dener           } else {
16808752603SAlp Dener             ierr = MatDuplicate(tao->hessian, MAT_COPY_VALUES, &bnk->H_inactive);CHKERRQ(ierr);
16928017e9fSAlp Dener           }
170c0f10754SAlp Dener           *needH = PETSC_FALSE;
171eb910715SAlp Dener         }
172eb910715SAlp Dener 
173eb910715SAlp Dener         for (i = 0; i < i_max; ++i) {
17462602cfbSAlp Dener           /* Take a steepest descent step and snap it to bounds */
17562602cfbSAlp Dener           ierr = VecCopy(tao->solution, bnk->Xold);CHKERRQ(ierr);
17662602cfbSAlp Dener           ierr = VecAXPY(tao->solution, -tao->trust/bnk->gnorm, tao->gradient);CHKERRQ(ierr);
177c4b75bccSAlp Dener           ierr = TaoBoundSolution(tao->XL, tao->XU, tao->solution,&nDiff);CHKERRQ(ierr);
178c4b75bccSAlp Dener           if (nDiff > 0) {
179c4b75bccSAlp Dener             /* Solution bounding changed variables so let's recompute the step we actually accepted */
180eb910715SAlp Dener             ierr = VecCopy(tao->solution, bnk->W);CHKERRQ(ierr);
18162602cfbSAlp Dener             ierr = VecAXPY(bnk->W, -1.0, bnk->Xold);CHKERRQ(ierr);
182c4b75bccSAlp Dener           }
18362602cfbSAlp Dener           /* Compute the objective at the trial */
18462602cfbSAlp Dener           ierr = TaoComputeObjective(tao, tao->solution, &ftrial);CHKERRQ(ierr);
18562602cfbSAlp Dener           ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr);
186eb910715SAlp Dener           if (PetscIsInfOrNanReal(ftrial)) {
187eb910715SAlp Dener             tau = bnk->gamma1_i;
188eb910715SAlp Dener           } else {
1890a4511e9SAlp Dener             if (ftrial < f_min) {
1900a4511e9SAlp Dener               f_min = ftrial;
191eb910715SAlp Dener               sigma = -tao->trust / bnk->gnorm;
192eb910715SAlp Dener             }
19308752603SAlp Dener 
194770b7498SAlp Dener             /* Compute the predicted and actual reduction */
1952ab2a32cSAlp Dener             if (bnk->inactive_idx) {
19608752603SAlp Dener               ierr = VecGetSubVector(bnk->W, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr);
19708752603SAlp Dener               ierr = VecGetSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work);CHKERRQ(ierr);
1982ab2a32cSAlp Dener             } else {
19908752603SAlp Dener               bnk->X_inactive = bnk->W;
20008752603SAlp Dener               bnk->inactive_work = bnk->Xwork;
2012ab2a32cSAlp Dener             }
20208752603SAlp Dener             ierr = MatMult(bnk->H_inactive, bnk->X_inactive, bnk->inactive_work);CHKERRQ(ierr);
20308752603SAlp Dener             ierr = VecDot(bnk->X_inactive, bnk->inactive_work, &prered);CHKERRQ(ierr);
2042ab2a32cSAlp Dener             if (bnk->inactive_idx) {
20508752603SAlp Dener               ierr = VecRestoreSubVector(bnk->W, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr);
20608752603SAlp Dener               ierr = VecRestoreSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work);CHKERRQ(ierr);
2072ab2a32cSAlp Dener             }
208eb910715SAlp Dener             prered = tao->trust * (bnk->gnorm - 0.5 * tao->trust * prered / (bnk->gnorm * bnk->gnorm));
209eb910715SAlp Dener             actred = bnk->f - ftrial;
210*3105154fSTodd Munson             if ((PetscAbsScalar(actred) <= bnk->epsilon) && (PetscAbsScalar(prered) <= bnk->epsilon)) {
211eb910715SAlp Dener               kappa = 1.0;
212*3105154fSTodd Munson             } else {
213eb910715SAlp Dener               kappa = actred / prered;
214eb910715SAlp Dener             }
215eb910715SAlp Dener 
216eb910715SAlp Dener             tau_1 = bnk->theta_i * bnk->gnorm * tao->trust / (bnk->theta_i * bnk->gnorm * tao->trust + (1.0 - bnk->theta_i) * prered - actred);
217eb910715SAlp Dener             tau_2 = bnk->theta_i * bnk->gnorm * tao->trust / (bnk->theta_i * bnk->gnorm * tao->trust - (1.0 + bnk->theta_i) * prered + actred);
218eb910715SAlp Dener             tau_min = PetscMin(tau_1, tau_2);
219eb910715SAlp Dener             tau_max = PetscMax(tau_1, tau_2);
220eb910715SAlp Dener 
221eb910715SAlp Dener             if (PetscAbsScalar(kappa - 1.0) <= bnk->mu1_i) {
222eb910715SAlp Dener               /*  Great agreement */
223eb910715SAlp Dener               max_radius = PetscMax(max_radius, tao->trust);
224eb910715SAlp Dener 
225eb910715SAlp Dener               if (tau_max < 1.0) {
226eb910715SAlp Dener                 tau = bnk->gamma3_i;
227*3105154fSTodd Munson               } else if (tau_max > bnk->gamma4_i) {
228eb910715SAlp Dener                 tau = bnk->gamma4_i;
229*3105154fSTodd Munson               } else {
230eb910715SAlp Dener                 tau = tau_max;
231eb910715SAlp Dener               }
23208752603SAlp Dener             }
23308752603SAlp Dener             else if (PetscAbsScalar(kappa - 1.0) <= bnk->mu2_i) {
234eb910715SAlp Dener               /*  Good agreement */
235eb910715SAlp Dener               max_radius = PetscMax(max_radius, tao->trust);
236eb910715SAlp Dener 
237eb910715SAlp Dener               if (tau_max < bnk->gamma2_i) {
238eb910715SAlp Dener                 tau = bnk->gamma2_i;
239eb910715SAlp Dener               } else if (tau_max > bnk->gamma3_i) {
240eb910715SAlp Dener                 tau = bnk->gamma3_i;
241eb910715SAlp Dener               } else {
242eb910715SAlp Dener                 tau = tau_max;
243eb910715SAlp Dener               }
24408752603SAlp Dener             }
24508752603SAlp Dener             else {
246eb910715SAlp Dener               /*  Not good agreement */
247eb910715SAlp Dener               if (tau_min > 1.0) {
248eb910715SAlp Dener                 tau = bnk->gamma2_i;
249eb910715SAlp Dener               } else if (tau_max < bnk->gamma1_i) {
250eb910715SAlp Dener                 tau = bnk->gamma1_i;
251eb910715SAlp Dener               } else if ((tau_min < bnk->gamma1_i) && (tau_max >= 1.0)) {
252eb910715SAlp Dener                 tau = bnk->gamma1_i;
253*3105154fSTodd Munson               } else if ((tau_1 >= bnk->gamma1_i) && (tau_1 < 1.0) && ((tau_2 < bnk->gamma1_i) || (tau_2 >= 1.0))) {
254eb910715SAlp Dener                 tau = tau_1;
255*3105154fSTodd Munson               } else if ((tau_2 >= bnk->gamma1_i) && (tau_2 < 1.0) && ((tau_1 < bnk->gamma1_i) || (tau_2 >= 1.0))) {
256eb910715SAlp Dener                 tau = tau_2;
257eb910715SAlp Dener               } else {
258eb910715SAlp Dener                 tau = tau_max;
259eb910715SAlp Dener               }
260eb910715SAlp Dener             }
261eb910715SAlp Dener           }
262eb910715SAlp Dener           tao->trust = tau * tao->trust;
263eb910715SAlp Dener         }
264eb910715SAlp Dener 
2650a4511e9SAlp Dener         if (f_min < bnk->f) {
266937a31a1SAlp Dener           /* We accidentally found a solution better than the initial, so accept it */
2670a4511e9SAlp Dener           bnk->f = f_min;
268937a31a1SAlp Dener           ierr = VecCopy(tao->solution, bnk->Xold);CHKERRQ(ierr);
269eb910715SAlp Dener           ierr = VecAXPY(tao->solution,sigma,tao->gradient);CHKERRQ(ierr);
270c4b75bccSAlp Dener           ierr = TaoBoundSolution(tao->XL, tao->XU, tao->solution, &nDiff);CHKERRQ(ierr);
271937a31a1SAlp Dener           ierr = VecCopy(tao->solution, tao->stepdirection);CHKERRQ(ierr);
272937a31a1SAlp Dener           ierr = VecAXPY(tao->stepdirection, -1.0, bnk->Xold);CHKERRQ(ierr);
27309164190SAlp Dener           ierr = TaoComputeGradient(tao,tao->solution,bnk->unprojected_gradient);CHKERRQ(ierr);
27461be54a6SAlp Dener           ierr = TaoBNKEstimateActiveSet(tao, bnk->as_type);CHKERRQ(ierr);
27561be54a6SAlp Dener           ierr = VecCopy(bnk->unprojected_gradient, tao->gradient);CHKERRQ(ierr);
2767529f6b4SAlp Dener           if (bnk->active_idx) {
27761be54a6SAlp Dener             ierr = VecISSet(tao->gradient, bnk->active_idx, 0.0);CHKERRQ(ierr);
2787529f6b4SAlp Dener           }
279937a31a1SAlp Dener           /* Compute gradient at the new iterate and flip switch to compute the Hessian later */
280c4b75bccSAlp Dener           ierr = VecNorm(tao->gradient, NORM_2, &bnk->gnorm);CHKERRQ(ierr);
281eb910715SAlp Dener           if (PetscIsInfOrNanReal(bnk->gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute gradient generated Inf or NaN");
282c0f10754SAlp Dener           *needH = PETSC_TRUE;
283937a31a1SAlp Dener           /* Test the new step for convergence */
284e031d6f5SAlp Dener           ierr = VecFischer(tao->solution, bnk->unprojected_gradient, tao->XL, tao->XU, bnk->Gwork);CHKERRQ(ierr);
285e031d6f5SAlp Dener           ierr = VecNorm(bnk->Gwork, NORM_2, &resnorm);CHKERRQ(ierr);
286e031d6f5SAlp Dener           ierr = TaoLogConvergenceHistory(tao,bnk->f,resnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
287e031d6f5SAlp Dener           ierr = TaoMonitor(tao,tao->niter,bnk->f,resnorm,0.0,1.0);CHKERRQ(ierr);
288eb910715SAlp Dener           ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
289eb910715SAlp Dener           if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
290937a31a1SAlp Dener           /* active BNCG recycling early because we have a stepdirection computed */
291937a31a1SAlp Dener           ierr = TaoBNCGSetRecycleFlag(bnk->bncg, PETSC_TRUE);CHKERRQ(ierr);
292eb910715SAlp Dener         }
293eb910715SAlp Dener       }
294eb910715SAlp Dener       tao->trust = PetscMax(tao->trust, max_radius);
295e031d6f5SAlp Dener 
296e031d6f5SAlp Dener       /* Ensure that the trust radius is within the limits */
297e031d6f5SAlp Dener       tao->trust = PetscMax(tao->trust, bnk->min_radius);
298e031d6f5SAlp Dener       tao->trust = PetscMin(tao->trust, bnk->max_radius);
299eb910715SAlp Dener       break;
300eb910715SAlp Dener 
301eb910715SAlp Dener     default:
302eb910715SAlp Dener       /* Norm of the first direction will initialize radius */
303eb910715SAlp Dener       tao->trust = 0.0;
304eb910715SAlp Dener       break;
305eb910715SAlp Dener     }
306eb910715SAlp Dener   }
307e031d6f5SAlp Dener 
308eb910715SAlp Dener   /* Set initial scaling for the BFGS preconditioner
309eb910715SAlp Dener      This step is done after computing the initial trust-region radius
310eb910715SAlp Dener      since the function value may have decreased */
311eb910715SAlp Dener   if (BNK_PC_BFGS == bnk->pc_type) {
3129b6ef848SAlp Dener     delta = 2.0 * PetscMax(1.0, PetscAbsScalar(bnk->f)) / (bnk->gnorm*bnk->gnorm);
313eb910715SAlp Dener     ierr = MatLMVMSetDelta(bnk->M,delta);CHKERRQ(ierr);
314eb910715SAlp Dener   }
315eb910715SAlp Dener   PetscFunctionReturn(0);
316eb910715SAlp Dener }
317eb910715SAlp Dener 
318df278d8fSAlp Dener /*------------------------------------------------------------*/
319df278d8fSAlp Dener 
32062675beeSAlp Dener /* Routine for computing the Hessian and preparing the preconditioner at the new iterate */
32162675beeSAlp Dener 
32262675beeSAlp Dener PetscErrorCode TaoBNKComputeHessian(Tao tao)
32362675beeSAlp Dener {
32462675beeSAlp Dener   PetscErrorCode               ierr;
32562675beeSAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
326c4b75bccSAlp Dener   PetscBool                    diagExists;
327c4b75bccSAlp Dener   PetscReal                    delta;
32862675beeSAlp Dener 
32962675beeSAlp Dener   PetscFunctionBegin;
33062675beeSAlp Dener   /* Compute the Hessian */
33162675beeSAlp Dener   ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr);
33262675beeSAlp Dener   /* Add a correction to the BFGS preconditioner */
33362675beeSAlp Dener   if (BNK_PC_BFGS == bnk->pc_type) {
33462675beeSAlp Dener     ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
33562675beeSAlp Dener     /* Update the BFGS diagonal scaling */
336c4b75bccSAlp Dener     ierr = MatHasOperation(tao->hessian, MATOP_GET_DIAGONAL, &diagExists);CHKERRQ(ierr);
337c4b75bccSAlp Dener     if (BFGS_SCALE_AHESS == bnk->bfgs_scale_type && diagExists) {
33862675beeSAlp Dener       ierr = MatGetDiagonal(tao->hessian, bnk->Diag);CHKERRQ(ierr);
33962675beeSAlp Dener       ierr = VecAbs(bnk->Diag);CHKERRQ(ierr);
34062675beeSAlp Dener       ierr = VecMedian(bnk->Diag_min, bnk->Diag, bnk->Diag_max, bnk->Diag);CHKERRQ(ierr);
34162675beeSAlp Dener       ierr = VecReciprocal(bnk->Diag);CHKERRQ(ierr);
34262675beeSAlp Dener       ierr = MatLMVMSetScale(bnk->M,bnk->Diag);CHKERRQ(ierr);
343c4b75bccSAlp Dener     } else {
344c4b75bccSAlp Dener       /* Fall back onto stand-alone BFGS scaling */
345c4b75bccSAlp Dener       delta = 2.0 * PetscMax(1.0, PetscAbsScalar(bnk->f)) / (bnk->gnorm*bnk->gnorm);
346c4b75bccSAlp Dener       ierr = MatLMVMSetDelta(bnk->M,delta);CHKERRQ(ierr);
34762675beeSAlp Dener     }
34862675beeSAlp Dener   }
34962675beeSAlp Dener   PetscFunctionReturn(0);
35062675beeSAlp Dener }
35162675beeSAlp Dener 
35262675beeSAlp Dener /*------------------------------------------------------------*/
35362675beeSAlp Dener 
3542f75a4aaSAlp Dener /* Routine for estimating the active set */
3552f75a4aaSAlp Dener 
35608752603SAlp Dener PetscErrorCode TaoBNKEstimateActiveSet(Tao tao, PetscInt asType)
3572f75a4aaSAlp Dener {
3582f75a4aaSAlp Dener   PetscErrorCode               ierr;
3592f75a4aaSAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
360c4b75bccSAlp Dener   PetscBool                    hessComputed, diagExists;
3612f75a4aaSAlp Dener 
3622f75a4aaSAlp Dener   PetscFunctionBegin;
36308752603SAlp Dener   switch (asType) {
3642f75a4aaSAlp Dener   case BNK_AS_NONE:
3652f75a4aaSAlp Dener     ierr = ISDestroy(&bnk->inactive_idx);CHKERRQ(ierr);
3662f75a4aaSAlp Dener     ierr = VecWhichInactive(tao->XL, tao->solution, bnk->unprojected_gradient, tao->XU, PETSC_TRUE, &bnk->inactive_idx);CHKERRQ(ierr);
3672f75a4aaSAlp Dener     ierr = ISDestroy(&bnk->active_idx);CHKERRQ(ierr);
3682f75a4aaSAlp Dener     ierr = ISComplementVec(bnk->inactive_idx, tao->solution, &bnk->active_idx);CHKERRQ(ierr);
3692f75a4aaSAlp Dener     break;
3702f75a4aaSAlp Dener 
3712f75a4aaSAlp Dener   case BNK_AS_BERTSEKAS:
3722f75a4aaSAlp Dener     /* Compute the trial step vector with which we will estimate the active set at the next iteration */
3732f75a4aaSAlp Dener     if (BNK_PC_BFGS == bnk->pc_type) {
3742f75a4aaSAlp Dener       /* If the BFGS preconditioner matrix is available, we will construct a trial step with it */
3755e9b73cbSAlp Dener       ierr = MatLMVMSetInactive(bnk->M, NULL);CHKERRQ(ierr);
3762f75a4aaSAlp Dener       ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, bnk->W);CHKERRQ(ierr);
3772f75a4aaSAlp Dener     } else {
37861be54a6SAlp Dener       ierr = MatAssembled(tao->hessian, &hessComputed);CHKERRQ(ierr);
379c4b75bccSAlp Dener       ierr = MatHasOperation(tao->hessian, MATOP_GET_DIAGONAL, &diagExists);CHKERRQ(ierr);
380c4b75bccSAlp Dener       if (hessComputed && diagExists) {
3819b6ef848SAlp Dener         /* BFGS preconditioner doesn't exist so let's invert the absolute diagonal of the Hessian instead onto the gradient */
3822f75a4aaSAlp Dener         ierr = MatGetDiagonal(tao->hessian, bnk->Xwork);CHKERRQ(ierr);
3839b6ef848SAlp Dener         ierr = VecAbs(bnk->Xwork);CHKERRQ(ierr);
3849b6ef848SAlp Dener         ierr = VecMedian(bnk->Diag_min, bnk->Xwork, bnk->Diag_max, bnk->Xwork);CHKERRQ(ierr);
3852f75a4aaSAlp Dener         ierr = VecReciprocal(bnk->Xwork);CHKERRQ(ierr);CHKERRQ(ierr);
3862f75a4aaSAlp Dener         ierr = VecPointwiseMult(bnk->W, bnk->Xwork, bnk->unprojected_gradient);CHKERRQ(ierr);
38761be54a6SAlp Dener       } else {
388c4b75bccSAlp Dener         /* If the Hessian or its diagonal does not exist, we will simply use gradient step */
38961be54a6SAlp Dener         ierr = VecCopy(bnk->unprojected_gradient, bnk->W);CHKERRQ(ierr);
39061be54a6SAlp Dener       }
3912f75a4aaSAlp Dener     }
3922f75a4aaSAlp Dener     ierr = VecScale(bnk->W, -1.0);CHKERRQ(ierr);
393c4b75bccSAlp Dener     ierr = TaoEstimateActiveBounds(tao->solution, tao->XL, tao->XU, bnk->unprojected_gradient, bnk->W, bnk->Xwork, bnk->as_step, &bnk->as_tol, &bnk->active_lower, &bnk->active_upper, &bnk->active_fixed, &bnk->active_idx, &bnk->inactive_idx);CHKERRQ(ierr);
394c4b75bccSAlp Dener     break;
3952f75a4aaSAlp Dener 
3962f75a4aaSAlp Dener   default:
3972f75a4aaSAlp Dener     break;
3982f75a4aaSAlp Dener   }
3992f75a4aaSAlp Dener   PetscFunctionReturn(0);
4002f75a4aaSAlp Dener }
4012f75a4aaSAlp Dener 
4022f75a4aaSAlp Dener /*------------------------------------------------------------*/
4032f75a4aaSAlp Dener 
4042f75a4aaSAlp Dener /* Routine for bounding the step direction */
4052f75a4aaSAlp Dener 
4062f75a4aaSAlp Dener PetscErrorCode TaoBNKBoundStep(Tao tao, Vec step)
4072f75a4aaSAlp Dener {
4082f75a4aaSAlp Dener   PetscErrorCode               ierr;
4092f75a4aaSAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
4102f75a4aaSAlp Dener 
4112f75a4aaSAlp Dener   PetscFunctionBegin;
4122f75a4aaSAlp Dener   switch (bnk->as_type) {
4132f75a4aaSAlp Dener   case BNK_AS_NONE:
414c4b75bccSAlp Dener     if (bnk->active_idx) {
415c4b75bccSAlp Dener       ierr = VecISSet(step, bnk->active_idx, 0.0);CHKERRQ(ierr);
416c4b75bccSAlp Dener     }
4172f75a4aaSAlp Dener     break;
4182f75a4aaSAlp Dener 
4192f75a4aaSAlp Dener   case BNK_AS_BERTSEKAS:
420c4b75bccSAlp Dener     ierr = TaoBoundStep(tao->solution, tao->XL, tao->XU, bnk->active_lower, bnk->active_upper, bnk->active_fixed, 1.0, step);CHKERRQ(ierr);
4212f75a4aaSAlp Dener     break;
4222f75a4aaSAlp Dener 
4232f75a4aaSAlp Dener   default:
4242f75a4aaSAlp Dener     break;
4252f75a4aaSAlp Dener   }
4262f75a4aaSAlp Dener   PetscFunctionReturn(0);
4272f75a4aaSAlp Dener }
4282f75a4aaSAlp Dener 
429e031d6f5SAlp Dener /*------------------------------------------------------------*/
430e031d6f5SAlp Dener 
431e031d6f5SAlp Dener /* Routine for taking a finite number of BNCG iterations to
432e031d6f5SAlp Dener    accelerate Newton convergence.
433e031d6f5SAlp Dener 
434e031d6f5SAlp Dener    In practice, this approach simply trades off Hessian evaluations
435e031d6f5SAlp Dener    for more gradient evaluations.
436e031d6f5SAlp Dener */
437e031d6f5SAlp Dener 
438c0f10754SAlp Dener PetscErrorCode TaoBNKTakeCGSteps(Tao tao, PetscBool *terminate)
439c0f10754SAlp Dener {
440c0f10754SAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
441c0f10754SAlp Dener   PetscErrorCode               ierr;
442c0f10754SAlp Dener 
443c0f10754SAlp Dener   PetscFunctionBegin;
444c0f10754SAlp Dener   *terminate = PETSC_FALSE;
445c0f10754SAlp Dener   if (bnk->max_cg_its > 0) {
446c4b75bccSAlp Dener     /* Copy the current function value (important vectors are already shared) */
447c0f10754SAlp Dener     bnk->bncg_ctx->f = bnk->f;
448c0f10754SAlp Dener     /* Take some small finite number of BNCG iterations */
449c0f10754SAlp Dener     ierr = TaoSolve(bnk->bncg);CHKERRQ(ierr);
450c0f10754SAlp Dener     /* Add the number of gradient and function evaluations to the total */
451c0f10754SAlp Dener     tao->nfuncs += bnk->bncg->nfuncs;
452c0f10754SAlp Dener     tao->nfuncgrads += bnk->bncg->nfuncgrads;
453c0f10754SAlp Dener     tao->ngrads += bnk->bncg->ngrads;
454c0f10754SAlp Dener     tao->nhess += bnk->bncg->nhess;
455e031d6f5SAlp Dener     bnk->tot_cg_its += bnk->bncg->niter;
456c4b75bccSAlp Dener     /* Extract the BNCG function value out and save it into BNK */
457c0f10754SAlp Dener     bnk->f = bnk->bncg_ctx->f;
458c0f10754SAlp 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) {
459c0f10754SAlp Dener       *terminate = PETSC_TRUE;
46061be54a6SAlp Dener     } else {
46161be54a6SAlp Dener       ierr = TaoBNKEstimateActiveSet(tao, bnk->as_type);
462c0f10754SAlp Dener     }
463c0f10754SAlp Dener   }
464c0f10754SAlp Dener   PetscFunctionReturn(0);
465c0f10754SAlp Dener }
466c0f10754SAlp Dener 
4672f75a4aaSAlp Dener /*------------------------------------------------------------*/
4682f75a4aaSAlp Dener 
469c0f10754SAlp Dener /* Routine for computing the Newton step. */
470df278d8fSAlp Dener 
47162675beeSAlp Dener PetscErrorCode TaoBNKComputeStep(Tao tao, PetscBool shift, KSPConvergedReason *ksp_reason)
472eb910715SAlp Dener {
473eb910715SAlp Dener   PetscErrorCode               ierr;
474eb910715SAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
475eb910715SAlp Dener 
476e465cd6fSAlp Dener   PetscReal                    delta;
477eb910715SAlp Dener   PetscInt                     bfgsUpdates = 0;
478eb910715SAlp Dener   PetscInt                     kspits;
479c4b75bccSAlp Dener   PetscBool                    diagExists;
480eb910715SAlp Dener 
481eb910715SAlp Dener   PetscFunctionBegin;
4825e9b73cbSAlp Dener   /* Prepare the reduced sub-matrices for the inactive set */
483*3105154fSTodd Munson   if (BNK_PC_BFGS == bnk->pc_type) {
484*3105154fSTodd Munson     ierr = MatLMVMSetInactive(bnk->M, bnk->inactive_idx);CHKERRQ(ierr);
485*3105154fSTodd Munson   }
4862f75a4aaSAlp Dener   if (bnk->inactive_idx) {
4875e9b73cbSAlp Dener     ierr = MatDestroy(&bnk->H_inactive);
4885e9b73cbSAlp Dener     ierr = MatCreateSubMatrix(tao->hessian, bnk->inactive_idx, bnk->inactive_idx, MAT_INITIAL_MATRIX, &bnk->H_inactive);CHKERRQ(ierr);
489eb3ba6a7SAlp Dener     if (tao->hessian == tao->hessian_pre) {
490eb3ba6a7SAlp Dener       bnk->Hpre_inactive = bnk->H_inactive;
491eb3ba6a7SAlp Dener     } else {
4925e9b73cbSAlp Dener       ierr = MatDestroy(&bnk->Hpre_inactive);
4935e9b73cbSAlp Dener       ierr = MatCreateSubMatrix(tao->hessian_pre, bnk->inactive_idx, bnk->inactive_idx, MAT_INITIAL_MATRIX, &bnk->Hpre_inactive);CHKERRQ(ierr);
494eb3ba6a7SAlp Dener     }
4952f75a4aaSAlp Dener   } else {
49662675beeSAlp Dener     ierr = MatDestroy(&bnk->H_inactive);
49762675beeSAlp Dener     ierr = MatDuplicate(tao->hessian, MAT_COPY_VALUES, &bnk->H_inactive);
49862675beeSAlp Dener     if (tao->hessian == tao->hessian_pre) {
49962675beeSAlp Dener       bnk->Hpre_inactive = bnk->H_inactive;
50062675beeSAlp Dener     } else {
50162675beeSAlp Dener       ierr = MatDestroy(&bnk->Hpre_inactive);
50262675beeSAlp Dener       ierr = MatDuplicate(tao->hessian_pre, MAT_COPY_VALUES, &bnk->Hpre_inactive);
50362675beeSAlp Dener     }
50462675beeSAlp Dener   }
50562675beeSAlp Dener 
50662675beeSAlp Dener   /* Shift the reduced Hessian matrix */
50762675beeSAlp Dener   if ((shift) && (bnk->pert > 0)) {
50862675beeSAlp Dener     ierr = MatShift(bnk->H_inactive, bnk->pert);CHKERRQ(ierr);
50962675beeSAlp Dener     if (bnk->H_inactive != bnk->Hpre_inactive) {
51062675beeSAlp Dener       ierr = MatShift(bnk->Hpre_inactive, bnk->pert);CHKERRQ(ierr);
51162675beeSAlp Dener     }
51262675beeSAlp Dener   }
51362675beeSAlp Dener 
51462675beeSAlp Dener   /* Update the diagonal scaling for the BFGS preconditioner, this time with the Hessian perturbation */
51562675beeSAlp Dener   if ((BNK_PC_BFGS == bnk->pc_type) && (BFGS_SCALE_PHESS == bnk->bfgs_scale_type)) {
516c4b75bccSAlp Dener     ierr = MatHasOperation(bnk->H_inactive, MATOP_GET_DIAGONAL, &diagExists);CHKERRQ(ierr);
517c4b75bccSAlp Dener     if (diagExists) {
51862675beeSAlp Dener       /* Obtain diagonal for the bfgs preconditioner  */
5195e9b73cbSAlp Dener       ierr = VecSet(bnk->Diag, 1.0);CHKERRQ(ierr);
5205e9b73cbSAlp Dener       if (bnk->inactive_idx) {
5215e9b73cbSAlp Dener         ierr = VecGetSubVector(bnk->Diag, bnk->inactive_idx, &bnk->Diag_red);CHKERRQ(ierr);
5225e9b73cbSAlp Dener       } else {
5235e9b73cbSAlp Dener         bnk->Diag_red = bnk->Diag;
5245e9b73cbSAlp Dener       }
5255e9b73cbSAlp Dener       ierr = MatGetDiagonal(bnk->H_inactive, bnk->Diag_red);CHKERRQ(ierr);
5265e9b73cbSAlp Dener       if (bnk->inactive_idx) {
5275e9b73cbSAlp Dener         ierr = VecRestoreSubVector(bnk->Diag, bnk->inactive_idx, &bnk->Diag_red);CHKERRQ(ierr);
5285e9b73cbSAlp Dener       }
52962675beeSAlp Dener       ierr = VecAbs(bnk->Diag);CHKERRQ(ierr);
53062675beeSAlp Dener       ierr = VecMedian(bnk->Diag_min, bnk->Diag, bnk->Diag_max, bnk->Diag);CHKERRQ(ierr);
53162675beeSAlp Dener       ierr = VecReciprocal(bnk->Diag);CHKERRQ(ierr);
53262675beeSAlp Dener       ierr = MatLMVMSetScale(bnk->M,bnk->Diag);CHKERRQ(ierr);
533c4b75bccSAlp Dener     } else {
534c4b75bccSAlp Dener       /* Fall back onto stand-alone BFGS scaling */
535c4b75bccSAlp Dener       delta = 2.0 * PetscMax(1.0, PetscAbsScalar(bnk->f)) / (bnk->gnorm*bnk->gnorm);
536c4b75bccSAlp Dener       ierr = MatLMVMSetDelta(bnk->M,delta);CHKERRQ(ierr);
537c4b75bccSAlp Dener     }
5382f75a4aaSAlp Dener   }
53909164190SAlp Dener 
540eb910715SAlp Dener   /* Solve the Newton system of equations */
541937a31a1SAlp Dener   tao->ksp_its = 0;
5422f75a4aaSAlp Dener   ierr = VecSet(tao->stepdirection, 0.0);CHKERRQ(ierr);
5435e9b73cbSAlp Dener   ierr = KSPReset(tao->ksp);CHKERRQ(ierr);
54409164190SAlp Dener   ierr = KSPSetOperators(tao->ksp,bnk->H_inactive,bnk->Hpre_inactive);CHKERRQ(ierr);
5455e9b73cbSAlp Dener   ierr = VecCopy(bnk->unprojected_gradient, bnk->Gwork);CHKERRQ(ierr);
5465e9b73cbSAlp Dener   if (bnk->inactive_idx) {
5475e9b73cbSAlp Dener     ierr = VecGetSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive);CHKERRQ(ierr);
5485e9b73cbSAlp Dener     ierr = VecGetSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr);
5495e9b73cbSAlp Dener   } else {
5505e9b73cbSAlp Dener     bnk->G_inactive = bnk->unprojected_gradient;
5515e9b73cbSAlp Dener     bnk->X_inactive = tao->stepdirection;
55228017e9fSAlp Dener   }
553eb910715SAlp Dener   if (bnk->is_nash || bnk->is_stcg || bnk->is_gltr) {
554fed79b8eSAlp Dener     ierr = KSPCGSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr);
5555e9b73cbSAlp Dener     ierr = KSPSolve(tao->ksp, bnk->G_inactive, bnk->X_inactive);CHKERRQ(ierr);
556eb910715SAlp Dener     ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr);
557eb910715SAlp Dener     tao->ksp_its+=kspits;
558eb910715SAlp Dener     tao->ksp_tot_its+=kspits;
559080d2917SAlp Dener     ierr = KSPCGGetNormD(tao->ksp,&bnk->dnorm);CHKERRQ(ierr);
560eb910715SAlp Dener 
561eb910715SAlp Dener     if (0.0 == tao->trust) {
562eb910715SAlp Dener       /* Radius was uninitialized; use the norm of the direction */
563080d2917SAlp Dener       if (bnk->dnorm > 0.0) {
564080d2917SAlp Dener         tao->trust = bnk->dnorm;
565eb910715SAlp Dener 
566eb910715SAlp Dener         /* Modify the radius if it is too large or small */
567eb910715SAlp Dener         tao->trust = PetscMax(tao->trust, bnk->min_radius);
568eb910715SAlp Dener         tao->trust = PetscMin(tao->trust, bnk->max_radius);
569eb910715SAlp Dener       } else {
570eb910715SAlp Dener         /* The direction was bad; set radius to default value and re-solve
571eb910715SAlp Dener            the trust-region subproblem to get a direction */
572eb910715SAlp Dener         tao->trust = tao->trust0;
573eb910715SAlp Dener 
574eb910715SAlp Dener         /* Modify the radius if it is too large or small */
575eb910715SAlp Dener         tao->trust = PetscMax(tao->trust, bnk->min_radius);
576eb910715SAlp Dener         tao->trust = PetscMin(tao->trust, bnk->max_radius);
577eb910715SAlp Dener 
578fed79b8eSAlp Dener         ierr = KSPCGSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr);
5795e9b73cbSAlp Dener         ierr = KSPSolve(tao->ksp, bnk->G_inactive, bnk->X_inactive);CHKERRQ(ierr);
580eb910715SAlp Dener         ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr);
581eb910715SAlp Dener         tao->ksp_its+=kspits;
582eb910715SAlp Dener         tao->ksp_tot_its+=kspits;
583080d2917SAlp Dener         ierr = KSPCGGetNormD(tao->ksp,&bnk->dnorm);CHKERRQ(ierr);
584eb910715SAlp Dener 
585080d2917SAlp Dener         if (bnk->dnorm == 0.0) SETERRQ(PETSC_COMM_SELF,1, "Initial direction zero");
586eb910715SAlp Dener       }
587eb910715SAlp Dener     }
588eb910715SAlp Dener   } else {
5895e9b73cbSAlp Dener     ierr = KSPSolve(tao->ksp, bnk->G_inactive, bnk->X_inactive);CHKERRQ(ierr);
590eb910715SAlp Dener     ierr = KSPGetIterationNumber(tao->ksp, &kspits);CHKERRQ(ierr);
591eb910715SAlp Dener     tao->ksp_its += kspits;
592eb910715SAlp Dener     tao->ksp_tot_its+=kspits;
593eb910715SAlp Dener   }
5945e9b73cbSAlp Dener   /* Restore sub vectors back */
5955e9b73cbSAlp Dener   if (bnk->inactive_idx) {
5965e9b73cbSAlp Dener     ierr = VecRestoreSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive);CHKERRQ(ierr);
5975e9b73cbSAlp Dener     ierr = VecRestoreSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr);
5985e9b73cbSAlp Dener   }
599770b7498SAlp Dener   /* Make sure the safeguarded fall-back step is zero for actively bounded variables */
600fed79b8eSAlp Dener   ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
6012f75a4aaSAlp Dener   ierr = TaoBNKBoundStep(tao, tao->stepdirection);CHKERRQ(ierr);
602770b7498SAlp Dener 
603770b7498SAlp Dener   /* Record convergence reasons */
604e465cd6fSAlp Dener   ierr = KSPGetConvergedReason(tao->ksp, ksp_reason);CHKERRQ(ierr);
605e465cd6fSAlp Dener   if (KSP_CONVERGED_ATOL == *ksp_reason) {
606770b7498SAlp Dener     ++bnk->ksp_atol;
607e465cd6fSAlp Dener   } else if (KSP_CONVERGED_RTOL == *ksp_reason) {
608770b7498SAlp Dener     ++bnk->ksp_rtol;
609e465cd6fSAlp Dener   } else if (KSP_CONVERGED_CG_CONSTRAINED == *ksp_reason) {
610770b7498SAlp Dener     ++bnk->ksp_ctol;
611e465cd6fSAlp Dener   } else if (KSP_CONVERGED_CG_NEG_CURVE == *ksp_reason) {
612770b7498SAlp Dener     ++bnk->ksp_negc;
613e465cd6fSAlp Dener   } else if (KSP_DIVERGED_DTOL == *ksp_reason) {
614770b7498SAlp Dener     ++bnk->ksp_dtol;
615e465cd6fSAlp Dener   } else if (KSP_DIVERGED_ITS == *ksp_reason) {
616770b7498SAlp Dener     ++bnk->ksp_iter;
617770b7498SAlp Dener   } else {
618770b7498SAlp Dener     ++bnk->ksp_othr;
619770b7498SAlp Dener   }
620fed79b8eSAlp Dener 
621fed79b8eSAlp Dener   /* Make sure the BFGS preconditioner is healthy */
622fed79b8eSAlp Dener   if (bnk->pc_type == BNK_PC_BFGS) {
623fed79b8eSAlp Dener     ierr = MatLMVMGetUpdates(bnk->M, &bfgsUpdates);CHKERRQ(ierr);
624e465cd6fSAlp Dener     if ((KSP_DIVERGED_INDEFINITE_PC == *ksp_reason) && (bfgsUpdates > 1)) {
625fed79b8eSAlp Dener       /* Preconditioner is numerically indefinite; reset the approximation. */
6269b6ef848SAlp Dener       delta = 2.0 * PetscMax(1.0, PetscAbsScalar(bnk->f)) / (bnk->gnorm*bnk->gnorm);
627eb910715SAlp Dener       ierr = MatLMVMSetDelta(bnk->M,delta);CHKERRQ(ierr);
628eb910715SAlp Dener       ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr);
62909164190SAlp Dener       ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
630eb910715SAlp Dener     }
631fed79b8eSAlp Dener   }
632e465cd6fSAlp Dener   PetscFunctionReturn(0);
633e465cd6fSAlp Dener }
634eb910715SAlp Dener 
63562675beeSAlp Dener /*------------------------------------------------------------*/
63662675beeSAlp Dener 
6375e9b73cbSAlp Dener /* Routine for recomputing the predicted reduction for a given step vector */
6385e9b73cbSAlp Dener 
6395e9b73cbSAlp Dener PetscErrorCode TaoBNKRecomputePred(Tao tao, Vec S, PetscReal *prered)
6405e9b73cbSAlp Dener {
6415e9b73cbSAlp Dener   PetscErrorCode               ierr;
6425e9b73cbSAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
6435e9b73cbSAlp Dener 
6445e9b73cbSAlp Dener   PetscFunctionBegin;
6455e9b73cbSAlp Dener   /* Extract subvectors associated with the inactive set */
6465e9b73cbSAlp Dener   if (bnk->inactive_idx){
6475e9b73cbSAlp Dener     ierr = VecGetSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr);
6485e9b73cbSAlp Dener     ierr = VecGetSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work);CHKERRQ(ierr);
6495e9b73cbSAlp Dener     ierr = VecGetSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive);CHKERRQ(ierr);
6505e9b73cbSAlp Dener   } else {
6515e9b73cbSAlp Dener     bnk->X_inactive = tao->stepdirection;
6525e9b73cbSAlp Dener     bnk->inactive_work = bnk->Xwork;
6535e9b73cbSAlp Dener     bnk->G_inactive = bnk->Gwork;
6545e9b73cbSAlp Dener   }
6555e9b73cbSAlp Dener   /* Recompute the predicted decrease based on the quadratic model */
6565e9b73cbSAlp Dener   ierr = MatMult(bnk->H_inactive, bnk->X_inactive, bnk->inactive_work);CHKERRQ(ierr);
6575e9b73cbSAlp Dener   ierr = VecAYPX(bnk->inactive_work, -0.5, bnk->G_inactive);CHKERRQ(ierr);
6585e9b73cbSAlp Dener   ierr = VecDot(bnk->inactive_work, bnk->X_inactive, prered);
6595e9b73cbSAlp Dener   /* Restore the sub vectors */
6605e9b73cbSAlp Dener   if (bnk->inactive_idx){
6615e9b73cbSAlp Dener     ierr = VecRestoreSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr);
6625e9b73cbSAlp Dener     ierr = VecRestoreSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work);CHKERRQ(ierr);
6635e9b73cbSAlp Dener     ierr = VecRestoreSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive);CHKERRQ(ierr);
6645e9b73cbSAlp Dener   }
6655e9b73cbSAlp Dener   PetscFunctionReturn(0);
6665e9b73cbSAlp Dener }
6675e9b73cbSAlp Dener 
6685e9b73cbSAlp Dener /*------------------------------------------------------------*/
6695e9b73cbSAlp Dener 
67062675beeSAlp Dener /* Routine for ensuring that the Newton step is a descent direction.
67162675beeSAlp Dener 
67262675beeSAlp Dener    The step direction falls back onto BFGS, scaled gradient and gradient steps
67362675beeSAlp Dener    in the event that the Newton step fails the test.
67462675beeSAlp Dener */
67562675beeSAlp Dener 
676e465cd6fSAlp Dener PetscErrorCode TaoBNKSafeguardStep(Tao tao, KSPConvergedReason ksp_reason, PetscInt *stepType)
677e465cd6fSAlp Dener {
678e465cd6fSAlp Dener   PetscErrorCode               ierr;
679e465cd6fSAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
680e465cd6fSAlp Dener 
681e465cd6fSAlp Dener   PetscReal                    gdx, delta, e_min;
682e465cd6fSAlp Dener   PetscInt                     bfgsUpdates;
683e465cd6fSAlp Dener 
684e465cd6fSAlp Dener   PetscFunctionBegin;
685080d2917SAlp Dener   ierr = VecDot(tao->stepdirection, tao->gradient, &gdx);CHKERRQ(ierr);
686eb910715SAlp Dener   if ((gdx >= 0.0) || PetscIsInfOrNanReal(gdx)) {
687eb910715SAlp Dener     /* Newton step is not descent or direction produced Inf or NaN
688eb910715SAlp Dener        Update the perturbation for next time */
689eb910715SAlp Dener     if (bnk->pert <= 0.0) {
690eb910715SAlp Dener       /* Initialize the perturbation */
691eb910715SAlp Dener       bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
692eb910715SAlp Dener       if (bnk->is_gltr) {
693eb910715SAlp Dener         ierr = KSPCGGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr);
694eb910715SAlp Dener         bnk->pert = PetscMax(bnk->pert, -e_min);
695eb910715SAlp Dener       }
696eb910715SAlp Dener     } else {
697eb910715SAlp Dener       /* Increase the perturbation */
698eb910715SAlp Dener       bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
699eb910715SAlp Dener     }
700eb910715SAlp Dener 
701eb910715SAlp Dener     if (BNK_PC_BFGS != bnk->pc_type) {
702eb910715SAlp Dener       /* We don't have the bfgs matrix around and updated
703eb910715SAlp Dener          Must use gradient direction in this case */
704080d2917SAlp Dener       ierr = VecCopy(tao->gradient, tao->stepdirection);CHKERRQ(ierr);
705eb910715SAlp Dener       *stepType = BNK_GRADIENT;
706eb910715SAlp Dener     } else {
707eb910715SAlp Dener       /* Attempt to use the BFGS direction */
7082ab2a32cSAlp Dener       ierr = MatLMVMSetInactive(bnk->M, NULL);CHKERRQ(ierr);
70909164190SAlp Dener       ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
710eb910715SAlp Dener 
7118d5ead36SAlp Dener       /* Check for success (descent direction)
7128d5ead36SAlp Dener          NOTE: Negative gdx here means not a descent direction because
7138d5ead36SAlp Dener          the fall-back step is missing a negative sign. */
714080d2917SAlp Dener       ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr);
715*3105154fSTodd Munson       if ((gdx <= 0.0) || PetscIsInfOrNanReal(gdx)) {
716eb910715SAlp Dener         /* BFGS direction is not descent or direction produced not a number
717eb910715SAlp Dener            We can assert bfgsUpdates > 1 in this case because
718eb910715SAlp Dener            the first solve produces the scaled gradient direction,
719eb910715SAlp Dener            which is guaranteed to be descent */
720eb910715SAlp Dener 
721eb910715SAlp Dener         /* Use steepest descent direction (scaled) */
7229b6ef848SAlp Dener         delta = 2.0 * PetscMax(1.0, PetscAbsScalar(bnk->f)) / (bnk->gnorm*bnk->gnorm);
723eb910715SAlp Dener         ierr = MatLMVMSetDelta(bnk->M, delta);CHKERRQ(ierr);
724eb910715SAlp Dener         ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr);
72509164190SAlp Dener         ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
72609164190SAlp Dener         ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
727eb910715SAlp Dener 
728eb910715SAlp Dener         *stepType = BNK_SCALED_GRADIENT;
729eb910715SAlp Dener       } else {
730770b7498SAlp Dener         ierr = MatLMVMGetUpdates(bnk->M, &bfgsUpdates);CHKERRQ(ierr);
731eb910715SAlp Dener         if (1 == bfgsUpdates) {
732eb910715SAlp Dener           /* The first BFGS direction is always the scaled gradient */
733eb910715SAlp Dener           *stepType = BNK_SCALED_GRADIENT;
734eb910715SAlp Dener         } else {
735eb910715SAlp Dener           *stepType = BNK_BFGS;
736eb910715SAlp Dener         }
737eb910715SAlp Dener       }
738eb910715SAlp Dener     }
7398d5ead36SAlp Dener     /* Make sure the safeguarded fall-back step is zero for actively bounded variables */
7408d5ead36SAlp Dener     ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
7412f75a4aaSAlp Dener     ierr = TaoBNKBoundStep(tao, tao->stepdirection);CHKERRQ(ierr);
742eb910715SAlp Dener   } else {
743eb910715SAlp Dener     /* Computed Newton step is descent */
744eb910715SAlp Dener     switch (ksp_reason) {
745eb910715SAlp Dener     case KSP_DIVERGED_NANORINF:
746eb910715SAlp Dener     case KSP_DIVERGED_BREAKDOWN:
747eb910715SAlp Dener     case KSP_DIVERGED_INDEFINITE_MAT:
748eb910715SAlp Dener     case KSP_DIVERGED_INDEFINITE_PC:
749eb910715SAlp Dener     case KSP_CONVERGED_CG_NEG_CURVE:
750eb910715SAlp Dener       /* Matrix or preconditioner is indefinite; increase perturbation */
751eb910715SAlp Dener       if (bnk->pert <= 0.0) {
752eb910715SAlp Dener         /* Initialize the perturbation */
753eb910715SAlp Dener         bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
754eb910715SAlp Dener         if (bnk->is_gltr) {
755eb910715SAlp Dener           ierr = KSPCGGLTRGetMinEig(tao->ksp, &e_min);CHKERRQ(ierr);
756eb910715SAlp Dener           bnk->pert = PetscMax(bnk->pert, -e_min);
757eb910715SAlp Dener         }
758eb910715SAlp Dener       } else {
759eb910715SAlp Dener         /* Increase the perturbation */
760eb910715SAlp Dener         bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
761eb910715SAlp Dener       }
762eb910715SAlp Dener       break;
763eb910715SAlp Dener 
764eb910715SAlp Dener     default:
765eb910715SAlp Dener       /* Newton step computation is good; decrease perturbation */
766eb910715SAlp Dener       bnk->pert = PetscMin(bnk->psfac * bnk->pert, bnk->pmsfac * bnk->gnorm);
767eb910715SAlp Dener       if (bnk->pert < bnk->pmin) {
768eb910715SAlp Dener         bnk->pert = 0.0;
769eb910715SAlp Dener       }
770eb910715SAlp Dener       break;
771eb910715SAlp Dener     }
772fed79b8eSAlp Dener     *stepType = BNK_NEWTON;
773eb910715SAlp Dener   }
774eb910715SAlp Dener   PetscFunctionReturn(0);
775eb910715SAlp Dener }
776eb910715SAlp Dener 
777df278d8fSAlp Dener /*------------------------------------------------------------*/
778df278d8fSAlp Dener 
779df278d8fSAlp Dener /* Routine for performing a bound-projected More-Thuente line search.
780df278d8fSAlp Dener 
781df278d8fSAlp Dener   Includes fallbacks to BFGS, scaled gradient, and unscaled gradient steps if the
782df278d8fSAlp Dener   Newton step does not produce a valid step length.
783df278d8fSAlp Dener */
784df278d8fSAlp Dener 
785937a31a1SAlp Dener PetscErrorCode TaoBNKPerformLineSearch(Tao tao, PetscInt *stepType, PetscReal *steplen, TaoLineSearchConvergedReason *reason)
786c14b763aSAlp Dener {
787c14b763aSAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
788c14b763aSAlp Dener   PetscErrorCode ierr;
789c14b763aSAlp Dener   TaoLineSearchConvergedReason ls_reason;
790c14b763aSAlp Dener 
791c14b763aSAlp Dener   PetscReal      e_min, gdx, delta;
792c14b763aSAlp Dener   PetscInt       bfgsUpdates;
793c14b763aSAlp Dener 
794c14b763aSAlp Dener   PetscFunctionBegin;
795c14b763aSAlp Dener   /* Perform the linesearch */
796c14b763aSAlp Dener   ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, steplen, &ls_reason);CHKERRQ(ierr);
797c14b763aSAlp Dener   ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
798c14b763aSAlp Dener 
799937a31a1SAlp Dener   while (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER && *stepType != BNK_GRADIENT) {
800c14b763aSAlp Dener     /* Linesearch failed, revert solution */
801c14b763aSAlp Dener     bnk->f = bnk->fold;
802c14b763aSAlp Dener     ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr);
803c14b763aSAlp Dener     ierr = VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient);CHKERRQ(ierr);
804c14b763aSAlp Dener 
805937a31a1SAlp Dener     switch(*stepType) {
806c14b763aSAlp Dener     case BNK_NEWTON:
8078d5ead36SAlp Dener       /* Failed to obtain acceptable iterate with Newton step
808c14b763aSAlp Dener          Update the perturbation for next time */
809c14b763aSAlp Dener       if (bnk->pert <= 0.0) {
810c14b763aSAlp Dener         /* Initialize the perturbation */
811c14b763aSAlp Dener         bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
812c14b763aSAlp Dener         if (bnk->is_gltr) {
813c14b763aSAlp Dener           ierr = KSPCGGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr);
814c14b763aSAlp Dener           bnk->pert = PetscMax(bnk->pert, -e_min);
815c14b763aSAlp Dener         }
816c14b763aSAlp Dener       } else {
817c14b763aSAlp Dener         /* Increase the perturbation */
818c14b763aSAlp Dener         bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
819c14b763aSAlp Dener       }
820c14b763aSAlp Dener 
821c14b763aSAlp Dener       if (BNK_PC_BFGS != bnk->pc_type) {
822c14b763aSAlp Dener         /* We don't have the bfgs matrix around and being updated
823c14b763aSAlp Dener            Must use gradient direction in this case */
824937a31a1SAlp Dener         ierr = VecCopy(bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
825937a31a1SAlp Dener         *stepType = BNK_GRADIENT;
826c14b763aSAlp Dener       } else {
827c14b763aSAlp Dener         /* Attempt to use the BFGS direction */
828937a31a1SAlp Dener         ierr = MatLMVMSetInactive(bnk->M, NULL);CHKERRQ(ierr);
829c14b763aSAlp Dener         ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
8308d5ead36SAlp Dener         /* Check for success (descent direction)
8318d5ead36SAlp Dener            NOTE: Negative gdx means not a descent direction because the step here is missing a negative sign. */
832c14b763aSAlp Dener         ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr);
833*3105154fSTodd Munson         if ((gdx <= 0.0) || PetscIsInfOrNanReal(gdx)) {
834c14b763aSAlp Dener           /* BFGS direction is not descent or direction produced not a number
835c14b763aSAlp Dener              We can assert bfgsUpdates > 1 in this case
836c14b763aSAlp Dener              Use steepest descent direction (scaled) */
8379b6ef848SAlp Dener           delta = 2.0 * PetscMax(1.0, PetscAbsScalar(bnk->f)) / (bnk->gnorm*bnk->gnorm);
838c14b763aSAlp Dener           ierr = MatLMVMSetDelta(bnk->M, delta);CHKERRQ(ierr);
839c14b763aSAlp Dener           ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr);
840c14b763aSAlp Dener           ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
841c14b763aSAlp Dener           ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
842c14b763aSAlp Dener 
843c14b763aSAlp Dener           bfgsUpdates = 1;
844937a31a1SAlp Dener           *stepType = BNK_SCALED_GRADIENT;
845c14b763aSAlp Dener         } else {
846c14b763aSAlp Dener           ierr = MatLMVMGetUpdates(bnk->M, &bfgsUpdates);CHKERRQ(ierr);
847c14b763aSAlp Dener           if (1 == bfgsUpdates) {
848c14b763aSAlp Dener             /* The first BFGS direction is always the scaled gradient */
849937a31a1SAlp Dener             *stepType = BNK_SCALED_GRADIENT;
850c14b763aSAlp Dener           } else {
851937a31a1SAlp Dener             *stepType = BNK_BFGS;
852c14b763aSAlp Dener           }
853c14b763aSAlp Dener         }
854c14b763aSAlp Dener       }
855c14b763aSAlp Dener       break;
856c14b763aSAlp Dener 
857c14b763aSAlp Dener     case BNK_BFGS:
858c14b763aSAlp Dener       /* Can only enter if pc_type == BNK_PC_BFGS
859c14b763aSAlp Dener          Failed to obtain acceptable iterate with BFGS step
860c14b763aSAlp Dener          Attempt to use the scaled gradient direction */
8619b6ef848SAlp Dener       delta = 2.0 * PetscMax(1.0, PetscAbsScalar(bnk->f)) / (bnk->gnorm*bnk->gnorm);
862c14b763aSAlp Dener       ierr = MatLMVMSetDelta(bnk->M, delta);CHKERRQ(ierr);
863c14b763aSAlp Dener       ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr);
864c14b763aSAlp Dener       ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
865937a31a1SAlp Dener       ierr = MatLMVMSetInactive(bnk->M, NULL);CHKERRQ(ierr);
866c14b763aSAlp Dener       ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
867c14b763aSAlp Dener 
868c14b763aSAlp Dener       bfgsUpdates = 1;
869937a31a1SAlp Dener       *stepType = BNK_SCALED_GRADIENT;
870c14b763aSAlp Dener       break;
871c14b763aSAlp Dener 
872c14b763aSAlp Dener     case BNK_SCALED_GRADIENT:
873c14b763aSAlp Dener       /* Can only enter if pc_type == BNK_PC_BFGS
874c14b763aSAlp Dener          The scaled gradient step did not produce a new iterate;
875937a31a1SAlp Dener          reset the BFGS matrix and attemp to use the gradient direction. */
876c14b763aSAlp Dener       ierr = MatLMVMSetScale(bnk->M,0);CHKERRQ(ierr);
877c14b763aSAlp Dener       ierr = MatLMVMSetDelta(bnk->M,1.0);CHKERRQ(ierr);
878c14b763aSAlp Dener       ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr);
879c14b763aSAlp Dener       ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
880937a31a1SAlp Dener       ierr = VecCopy(bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
881c14b763aSAlp Dener 
882c14b763aSAlp Dener       bfgsUpdates = 1;
883937a31a1SAlp Dener       *stepType = BNK_GRADIENT;
884c14b763aSAlp Dener       break;
885c14b763aSAlp Dener     }
8868d5ead36SAlp Dener     /* Make sure the safeguarded fall-back step is zero for actively bounded variables */
8878d5ead36SAlp Dener     ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
8882f75a4aaSAlp Dener     ierr = TaoBNKBoundStep(tao, tao->stepdirection);CHKERRQ(ierr);
889c14b763aSAlp Dener 
8908d5ead36SAlp Dener     /* Perform one last line search with the fall-back step */
891c14b763aSAlp Dener     ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, steplen, &ls_reason);CHKERRQ(ierr);
892c14b763aSAlp Dener     ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
893c14b763aSAlp Dener   }
894c14b763aSAlp Dener   *reason = ls_reason;
895c14b763aSAlp Dener   PetscFunctionReturn(0);
896c14b763aSAlp Dener }
897c14b763aSAlp Dener 
898df278d8fSAlp Dener /*------------------------------------------------------------*/
899df278d8fSAlp Dener 
900df278d8fSAlp Dener /* Routine for updating the trust radius.
901df278d8fSAlp Dener 
902df278d8fSAlp Dener   Function features three different update methods:
903df278d8fSAlp Dener   1) Line-search step length based
904df278d8fSAlp Dener   2) Predicted decrease on the CG quadratic model
905df278d8fSAlp Dener   3) Interpolation
906df278d8fSAlp Dener */
907df278d8fSAlp Dener 
90828017e9fSAlp Dener PetscErrorCode TaoBNKUpdateTrustRadius(Tao tao, PetscReal prered, PetscReal actred, PetscInt updateType, PetscInt stepType, PetscBool *accept)
909080d2917SAlp Dener {
910080d2917SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
911080d2917SAlp Dener   PetscErrorCode ierr;
912080d2917SAlp Dener 
913b1c2d0e3SAlp Dener   PetscReal      step, kappa;
914080d2917SAlp Dener   PetscReal      gdx, tau_1, tau_2, tau_min, tau_max;
915080d2917SAlp Dener 
916080d2917SAlp Dener   PetscFunctionBegin;
917080d2917SAlp Dener   /* Update trust region radius */
918080d2917SAlp Dener   *accept = PETSC_FALSE;
91928017e9fSAlp Dener   switch(updateType) {
920080d2917SAlp Dener   case BNK_UPDATE_STEP:
921c14b763aSAlp Dener     *accept = PETSC_TRUE; /* always accept here because line search succeeded */
922080d2917SAlp Dener     if (stepType == BNK_NEWTON) {
923080d2917SAlp Dener       ierr = TaoLineSearchGetStepLength(tao->linesearch, &step);CHKERRQ(ierr);
924080d2917SAlp Dener       if (step < bnk->nu1) {
925080d2917SAlp Dener         /* Very bad step taken; reduce radius */
926080d2917SAlp Dener         tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust);
927080d2917SAlp Dener       } else if (step < bnk->nu2) {
928080d2917SAlp Dener         /* Reasonably bad step taken; reduce radius */
929080d2917SAlp Dener         tao->trust = bnk->omega2 * PetscMin(bnk->dnorm, tao->trust);
930080d2917SAlp Dener       } else if (step < bnk->nu3) {
931080d2917SAlp Dener         /*  Reasonable step was taken; leave radius alone */
932080d2917SAlp Dener         if (bnk->omega3 < 1.0) {
933080d2917SAlp Dener           tao->trust = bnk->omega3 * PetscMin(bnk->dnorm, tao->trust);
934080d2917SAlp Dener         } else if (bnk->omega3 > 1.0) {
935080d2917SAlp Dener           tao->trust = PetscMax(bnk->omega3 * bnk->dnorm, tao->trust);
936080d2917SAlp Dener         }
937080d2917SAlp Dener       } else if (step < bnk->nu4) {
938080d2917SAlp Dener         /*  Full step taken; increase the radius */
939080d2917SAlp Dener         tao->trust = PetscMax(bnk->omega4 * bnk->dnorm, tao->trust);
940080d2917SAlp Dener       } else {
941080d2917SAlp Dener         /*  More than full step taken; increase the radius */
942080d2917SAlp Dener         tao->trust = PetscMax(bnk->omega5 * bnk->dnorm, tao->trust);
943080d2917SAlp Dener       }
944080d2917SAlp Dener     } else {
945080d2917SAlp Dener       /*  Newton step was not good; reduce the radius */
946080d2917SAlp Dener       tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust);
947080d2917SAlp Dener     }
948080d2917SAlp Dener     break;
949080d2917SAlp Dener 
950080d2917SAlp Dener   case BNK_UPDATE_REDUCTION:
951080d2917SAlp Dener     if (stepType == BNK_NEWTON) {
952b1c2d0e3SAlp Dener       if (prered < 0.0) {
953fed79b8eSAlp Dener         /* The predicted reduction has the wrong sign.  This cannot
954fed79b8eSAlp Dener            happen in infinite precision arithmetic.  Step should
955fed79b8eSAlp Dener            be rejected! */
956080d2917SAlp Dener         tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
957*3105154fSTodd Munson       } else {
958b1c2d0e3SAlp Dener         if (PetscIsInfOrNanReal(actred)) {
959080d2917SAlp Dener           tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
960080d2917SAlp Dener         } else {
961*3105154fSTodd Munson           if ((PetscAbsScalar(actred) <= PetscMax(1.0, PetscAbsScalar(bnk->f))*bnk->epsilon) && (PetscAbsScalar(prered) <= PetscMax(1.0, PetscAbsScalar(bnk->f))*bnk->epsilon)) {
962080d2917SAlp Dener             kappa = 1.0;
963*3105154fSTodd Munson           } else {
964080d2917SAlp Dener             kappa = actred / prered;
965080d2917SAlp Dener           }
966fed79b8eSAlp Dener 
967fed79b8eSAlp Dener           /* Accept or reject the step and update radius */
968080d2917SAlp Dener           if (kappa < bnk->eta1) {
969fed79b8eSAlp Dener             /* Reject the step */
970080d2917SAlp Dener             tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
971*3105154fSTodd Munson           } else {
972fed79b8eSAlp Dener             /* Accept the step */
973c133c014SAlp Dener             *accept = PETSC_TRUE;
974c133c014SAlp Dener             /* Update the trust region radius only if the computed step is at the trust radius boundary */
9758d5ead36SAlp Dener             if (bnk->dnorm == tao->trust) {
976080d2917SAlp Dener               if (kappa < bnk->eta2) {
977080d2917SAlp Dener                 /* Marginal bad step */
978c133c014SAlp Dener                 tao->trust = bnk->alpha2 * tao->trust;
979*3105154fSTodd Munson               } else if (kappa < bnk->eta3) {
980fed79b8eSAlp Dener                 /* Reasonable step */
981fed79b8eSAlp Dener                 tao->trust = bnk->alpha3 * tao->trust;
982*3105154fSTodd Munson               } else if (kappa < bnk->eta4) {
983080d2917SAlp Dener                 /* Good step */
984c133c014SAlp Dener                 tao->trust = bnk->alpha4 * tao->trust;
985*3105154fSTodd Munson               } else {
986080d2917SAlp Dener                 /* Very good step */
987c133c014SAlp Dener                 tao->trust = bnk->alpha5 * tao->trust;
988080d2917SAlp Dener               }
989c133c014SAlp Dener             }
990080d2917SAlp Dener           }
991080d2917SAlp Dener         }
992080d2917SAlp Dener       }
993080d2917SAlp Dener     } else {
994080d2917SAlp Dener       /*  Newton step was not good; reduce the radius */
995080d2917SAlp Dener       tao->trust = bnk->alpha1 * PetscMin(bnk->dnorm, tao->trust);
996080d2917SAlp Dener     }
997080d2917SAlp Dener     break;
998080d2917SAlp Dener 
999080d2917SAlp Dener   default:
1000080d2917SAlp Dener     if (stepType == BNK_NEWTON) {
1001b1c2d0e3SAlp Dener       if (prered < 0.0) {
1002080d2917SAlp Dener         /*  The predicted reduction has the wrong sign.  This cannot */
1003080d2917SAlp Dener         /*  happen in infinite precision arithmetic.  Step should */
1004080d2917SAlp Dener         /*  be rejected! */
1005080d2917SAlp Dener         tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
1006080d2917SAlp Dener       } else {
1007b1c2d0e3SAlp Dener         if (PetscIsInfOrNanReal(actred)) {
1008080d2917SAlp Dener           tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
1009080d2917SAlp Dener         } else {
1010080d2917SAlp Dener           if ((PetscAbsScalar(actred) <= bnk->epsilon) && (PetscAbsScalar(prered) <= bnk->epsilon)) {
1011080d2917SAlp Dener             kappa = 1.0;
1012080d2917SAlp Dener           } else {
1013080d2917SAlp Dener             kappa = actred / prered;
1014080d2917SAlp Dener           }
1015080d2917SAlp Dener 
1016080d2917SAlp Dener           ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr);
1017080d2917SAlp Dener           tau_1 = bnk->theta * gdx / (bnk->theta * gdx - (1.0 - bnk->theta) * prered + actred);
1018080d2917SAlp Dener           tau_2 = bnk->theta * gdx / (bnk->theta * gdx + (1.0 + bnk->theta) * prered - actred);
1019080d2917SAlp Dener           tau_min = PetscMin(tau_1, tau_2);
1020080d2917SAlp Dener           tau_max = PetscMax(tau_1, tau_2);
1021080d2917SAlp Dener 
1022080d2917SAlp Dener           if (kappa >= 1.0 - bnk->mu1) {
1023080d2917SAlp Dener             /*  Great agreement */
1024080d2917SAlp Dener             *accept = PETSC_TRUE;
1025080d2917SAlp Dener             if (tau_max < 1.0) {
1026080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm);
1027080d2917SAlp Dener             } else if (tau_max > bnk->gamma4) {
1028080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, bnk->gamma4 * bnk->dnorm);
1029080d2917SAlp Dener             } else {
1030080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm);
1031080d2917SAlp Dener             }
1032080d2917SAlp Dener           } else if (kappa >= 1.0 - bnk->mu2) {
1033080d2917SAlp Dener             /*  Good agreement */
1034080d2917SAlp Dener             *accept = PETSC_TRUE;
1035080d2917SAlp Dener             if (tau_max < bnk->gamma2) {
1036080d2917SAlp Dener               tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm);
1037080d2917SAlp Dener             } else if (tau_max > bnk->gamma3) {
1038080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm);
1039080d2917SAlp Dener             } else if (tau_max < 1.0) {
1040080d2917SAlp Dener               tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm);
1041080d2917SAlp Dener             } else {
1042080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm);
1043080d2917SAlp Dener             }
1044080d2917SAlp Dener           } else {
1045080d2917SAlp Dener             /*  Not good agreement */
1046080d2917SAlp Dener             if (tau_min > 1.0) {
1047080d2917SAlp Dener               tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm);
1048080d2917SAlp Dener             } else if (tau_max < bnk->gamma1) {
1049080d2917SAlp Dener               tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
1050080d2917SAlp Dener             } else if ((tau_min < bnk->gamma1) && (tau_max >= 1.0)) {
1051080d2917SAlp Dener               tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
1052080d2917SAlp Dener             } else if ((tau_1 >= bnk->gamma1) && (tau_1 < 1.0) && ((tau_2 < bnk->gamma1) || (tau_2 >= 1.0))) {
1053080d2917SAlp Dener               tao->trust = tau_1 * PetscMin(tao->trust, bnk->dnorm);
1054080d2917SAlp Dener             } else if ((tau_2 >= bnk->gamma1) && (tau_2 < 1.0) && ((tau_1 < bnk->gamma1) || (tau_2 >= 1.0))) {
1055080d2917SAlp Dener               tao->trust = tau_2 * PetscMin(tao->trust, bnk->dnorm);
1056080d2917SAlp Dener             } else {
1057080d2917SAlp Dener               tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm);
1058080d2917SAlp Dener             }
1059080d2917SAlp Dener           }
1060080d2917SAlp Dener         }
1061080d2917SAlp Dener       }
1062080d2917SAlp Dener     } else {
1063080d2917SAlp Dener       /*  Newton step was not good; reduce the radius */
1064080d2917SAlp Dener       tao->trust = bnk->gamma1 * PetscMin(bnk->dnorm, tao->trust);
1065080d2917SAlp Dener     }
106628017e9fSAlp Dener     break;
1067080d2917SAlp Dener   }
1068c133c014SAlp Dener   /* Make sure the radius does not violate min and max settings */
1069c133c014SAlp Dener   tao->trust = PetscMin(tao->trust, bnk->max_radius);
1070fed79b8eSAlp Dener   tao->trust = PetscMax(tao->trust, bnk->min_radius);
1071080d2917SAlp Dener   PetscFunctionReturn(0);
1072080d2917SAlp Dener }
1073080d2917SAlp Dener 
1074eb910715SAlp Dener /* ---------------------------------------------------------- */
1075df278d8fSAlp Dener 
107662675beeSAlp Dener PetscErrorCode TaoBNKAddStepCounts(Tao tao, PetscInt stepType)
107762675beeSAlp Dener {
107862675beeSAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
107962675beeSAlp Dener 
108062675beeSAlp Dener   PetscFunctionBegin;
108162675beeSAlp Dener   switch (stepType) {
108262675beeSAlp Dener   case BNK_NEWTON:
108362675beeSAlp Dener     ++bnk->newt;
108462675beeSAlp Dener     break;
108562675beeSAlp Dener   case BNK_BFGS:
108662675beeSAlp Dener     ++bnk->bfgs;
108762675beeSAlp Dener     break;
108862675beeSAlp Dener   case BNK_SCALED_GRADIENT:
108962675beeSAlp Dener     ++bnk->sgrad;
109062675beeSAlp Dener     break;
109162675beeSAlp Dener   case BNK_GRADIENT:
109262675beeSAlp Dener     ++bnk->grad;
109362675beeSAlp Dener     break;
109462675beeSAlp Dener   default:
109562675beeSAlp Dener     break;
109662675beeSAlp Dener   }
109762675beeSAlp Dener   PetscFunctionReturn(0);
109862675beeSAlp Dener }
109962675beeSAlp Dener 
110062675beeSAlp Dener /* ---------------------------------------------------------- */
110162675beeSAlp Dener 
11029b6ef848SAlp Dener PetscErrorCode TaoSetUp_BNK(Tao tao)
1103eb910715SAlp Dener {
1104eb910715SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
1105eb910715SAlp Dener   PetscErrorCode ierr;
11069b6ef848SAlp Dener   KSPType        ksp_type;
1107e031d6f5SAlp Dener   PetscInt       i;
1108eb910715SAlp Dener 
1109eb910715SAlp Dener   PetscFunctionBegin;
1110c4b75bccSAlp Dener   if (!tao->gradient) {
1111c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);
1112c4b75bccSAlp Dener   }
1113c4b75bccSAlp Dener   if (!tao->stepdirection) {
1114c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr);
1115c4b75bccSAlp Dener   }
1116c4b75bccSAlp Dener   if (!bnk->W) {
1117c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->W);CHKERRQ(ierr);
1118c4b75bccSAlp Dener   }
1119c4b75bccSAlp Dener   if (!bnk->Xold) {
1120c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->Xold);CHKERRQ(ierr);
1121c4b75bccSAlp Dener   }
1122c4b75bccSAlp Dener   if (!bnk->Gold) {
1123c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->Gold);CHKERRQ(ierr);
1124c4b75bccSAlp Dener   }
1125c4b75bccSAlp Dener   if (!bnk->Xwork) {
1126c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->Xwork);CHKERRQ(ierr);
1127c4b75bccSAlp Dener   }
1128c4b75bccSAlp Dener   if (!bnk->Gwork) {
1129c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->Gwork);CHKERRQ(ierr);
1130c4b75bccSAlp Dener   }
1131c4b75bccSAlp Dener   if (!bnk->unprojected_gradient) {
1132c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->unprojected_gradient);CHKERRQ(ierr);
1133c4b75bccSAlp Dener   }
1134c4b75bccSAlp Dener   if (!bnk->unprojected_gradient_old) {
1135c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->unprojected_gradient_old);CHKERRQ(ierr);
1136c4b75bccSAlp Dener   }
1137c4b75bccSAlp Dener   if (!bnk->Diag_min) {
1138c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->Diag_min);CHKERRQ(ierr);
1139c4b75bccSAlp Dener   }
1140c4b75bccSAlp Dener   if (!bnk->Diag_max) {
1141c4b75bccSAlp Dener     ierr = VecDuplicate(tao->solution,&bnk->Diag_max);CHKERRQ(ierr);
1142c4b75bccSAlp Dener   }
1143e031d6f5SAlp Dener   if (bnk->max_cg_its > 0) {
1144c4b75bccSAlp Dener     /* Ensure that the important common vectors are shared between BNK and embedded BNCG */
1145c4b75bccSAlp Dener     bnk->bncg_ctx = (TAO_BNCG *)bnk->bncg->data;
1146c4b75bccSAlp Dener     ierr = PetscObjectReference((PetscObject)(bnk->unprojected_gradient));CHKERRQ(ierr);
1147c4b75bccSAlp Dener     ierr = VecDestroy(&bnk->bncg_ctx->unprojected_gradient);CHKERRQ(ierr);
1148c4b75bccSAlp Dener     bnk->bncg_ctx->unprojected_gradient = bnk->unprojected_gradient;
1149c4b75bccSAlp Dener     ierr = PetscObjectReference((PetscObject)(bnk->Gold));CHKERRQ(ierr);
1150c4b75bccSAlp Dener     ierr = VecDestroy(&bnk->bncg_ctx->G_old);CHKERRQ(ierr);
1151c4b75bccSAlp Dener     bnk->bncg_ctx->G_old = bnk->Gold;
1152c4b75bccSAlp Dener     ierr = PetscObjectReference((PetscObject)(tao->gradient));CHKERRQ(ierr);
1153c4b75bccSAlp Dener     ierr = VecDestroy(&bnk->bncg->gradient);CHKERRQ(ierr);
1154c4b75bccSAlp Dener     bnk->bncg->gradient = tao->gradient;
1155c4b75bccSAlp Dener     ierr = PetscObjectReference((PetscObject)(tao->stepdirection));CHKERRQ(ierr);
1156c4b75bccSAlp Dener     ierr = VecDestroy(&bnk->bncg->stepdirection);CHKERRQ(ierr);
1157c4b75bccSAlp Dener     bnk->bncg->stepdirection = tao->stepdirection;
1158c4b75bccSAlp Dener     ierr = TaoSetInitialVector(bnk->bncg, tao->solution);CHKERRQ(ierr);
1159c4b75bccSAlp Dener     /* Copy over some settings from BNK into BNCG */
1160e031d6f5SAlp Dener     ierr = TaoSetMaximumIterations(bnk->bncg, bnk->max_cg_its);CHKERRQ(ierr);
1161e031d6f5SAlp Dener     ierr = TaoSetTolerances(bnk->bncg, tao->gatol, tao->grtol, tao->gttol);CHKERRQ(ierr);
1162e031d6f5SAlp Dener     ierr = TaoSetFunctionLowerBound(bnk->bncg, tao->fmin);CHKERRQ(ierr);
1163937a31a1SAlp Dener     ierr = TaoSetConvergenceTest(bnk->bncg, tao->ops->convergencetest, tao->cnvP);CHKERRQ(ierr);
1164e031d6f5SAlp Dener     ierr = TaoSetObjectiveRoutine(bnk->bncg, tao->ops->computeobjective, tao->user_objP);CHKERRQ(ierr);
1165e031d6f5SAlp Dener     ierr = TaoSetGradientRoutine(bnk->bncg, tao->ops->computegradient, tao->user_gradP);CHKERRQ(ierr);
1166e031d6f5SAlp Dener     ierr = TaoSetObjectiveAndGradientRoutine(bnk->bncg, tao->ops->computeobjectiveandgradient, tao->user_objgradP);CHKERRQ(ierr);
1167e031d6f5SAlp Dener     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)tao, (PetscObject)(bnk->bncg));CHKERRQ(ierr);
1168c4b75bccSAlp Dener     for (i=0; i<tao->numbermonitors; ++i) {
1169e031d6f5SAlp Dener       ierr = TaoSetMonitor(bnk->bncg, tao->monitor[i], tao->monitorcontext[i], tao->monitordestroy[i]);CHKERRQ(ierr);
1170e031d6f5SAlp Dener       ierr = PetscObjectReference((PetscObject)(tao->monitorcontext[i]));CHKERRQ(ierr);
1171e031d6f5SAlp Dener     }
1172e031d6f5SAlp Dener   }
1173eb910715SAlp Dener   bnk->Diag = 0;
1174c0f10754SAlp Dener   bnk->Diag_red = 0;
1175c0f10754SAlp Dener   bnk->X_inactive = 0;
1176c0f10754SAlp Dener   bnk->G_inactive = 0;
117762675beeSAlp Dener   bnk->inactive_work = 0;
117862675beeSAlp Dener   bnk->active_work = 0;
117962675beeSAlp Dener   bnk->inactive_idx = 0;
118062675beeSAlp Dener   bnk->active_idx = 0;
118162675beeSAlp Dener   bnk->active_lower = 0;
118262675beeSAlp Dener   bnk->active_upper = 0;
118362675beeSAlp Dener   bnk->active_fixed = 0;
1184eb910715SAlp Dener   bnk->M = 0;
118509164190SAlp Dener   bnk->H_inactive = 0;
118609164190SAlp Dener   bnk->Hpre_inactive = 0;
11879b6ef848SAlp Dener   ierr = KSPGetType(tao->ksp,&ksp_type);CHKERRQ(ierr);
11889b6ef848SAlp Dener   ierr = PetscStrcmp(ksp_type,KSPCGNASH,&bnk->is_nash);CHKERRQ(ierr);
11899b6ef848SAlp Dener   ierr = PetscStrcmp(ksp_type,KSPCGSTCG,&bnk->is_stcg);CHKERRQ(ierr);
11909b6ef848SAlp Dener   ierr = PetscStrcmp(ksp_type,KSPCGGLTR,&bnk->is_gltr);CHKERRQ(ierr);
1191eb910715SAlp Dener   PetscFunctionReturn(0);
1192eb910715SAlp Dener }
1193eb910715SAlp Dener 
1194eb910715SAlp Dener /*------------------------------------------------------------*/
1195df278d8fSAlp Dener 
1196eb910715SAlp Dener static PetscErrorCode TaoDestroy_BNK(Tao tao)
1197eb910715SAlp Dener {
1198eb910715SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
1199eb910715SAlp Dener   PetscErrorCode ierr;
1200eb910715SAlp Dener 
1201eb910715SAlp Dener   PetscFunctionBegin;
1202eb910715SAlp Dener   if (tao->setupcalled) {
1203eb910715SAlp Dener     ierr = VecDestroy(&bnk->W);CHKERRQ(ierr);
1204eb910715SAlp Dener     ierr = VecDestroy(&bnk->Xold);CHKERRQ(ierr);
1205eb910715SAlp Dener     ierr = VecDestroy(&bnk->Gold);CHKERRQ(ierr);
120609164190SAlp Dener     ierr = VecDestroy(&bnk->Xwork);CHKERRQ(ierr);
120709164190SAlp Dener     ierr = VecDestroy(&bnk->Gwork);CHKERRQ(ierr);
120809164190SAlp Dener     ierr = VecDestroy(&bnk->unprojected_gradient);CHKERRQ(ierr);
120909164190SAlp Dener     ierr = VecDestroy(&bnk->unprojected_gradient_old);CHKERRQ(ierr);
121062675beeSAlp Dener     ierr = VecDestroy(&bnk->Diag_min);CHKERRQ(ierr);
121162675beeSAlp Dener     ierr = VecDestroy(&bnk->Diag_max);CHKERRQ(ierr);
1212c4b75bccSAlp Dener   }
1213c0f10754SAlp Dener   if (bnk->max_cg_its > 0) {
1214c0f10754SAlp Dener     ierr = TaoDestroy(&bnk->bncg);CHKERRQ(ierr);
12155e9b73cbSAlp Dener   }
12165e9b73cbSAlp Dener   ierr = VecDestroy(&bnk->Diag);CHKERRQ(ierr);
1217eb910715SAlp Dener   ierr = MatDestroy(&bnk->M);CHKERRQ(ierr);
1218c4b75bccSAlp Dener   if (bnk->Hpre_inactive != tao->hessian_pre && bnk->Hpre_inactive != bnk->H_inactive) {
1219c4b75bccSAlp Dener     ierr = MatDestroy(&bnk->Hpre_inactive);CHKERRQ(ierr);
1220c4b75bccSAlp Dener   }
1221c4b75bccSAlp Dener   if (bnk->H_inactive != tao->hessian) {
1222c4b75bccSAlp Dener     ierr = MatDestroy(&bnk->H_inactive);CHKERRQ(ierr);
1223c4b75bccSAlp Dener   }
1224eb910715SAlp Dener   ierr = PetscFree(tao->data);CHKERRQ(ierr);
1225eb910715SAlp Dener   PetscFunctionReturn(0);
1226eb910715SAlp Dener }
1227eb910715SAlp Dener 
1228eb910715SAlp Dener /*------------------------------------------------------------*/
1229df278d8fSAlp Dener 
1230eb910715SAlp Dener static PetscErrorCode TaoSetFromOptions_BNK(PetscOptionItems *PetscOptionsObject,Tao tao)
1231eb910715SAlp Dener {
1232eb910715SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
1233eb910715SAlp Dener   PetscErrorCode ierr;
1234eb910715SAlp Dener 
1235eb910715SAlp Dener   PetscFunctionBegin;
1236eb910715SAlp Dener   ierr = PetscOptionsHead(PetscOptionsObject,"Newton line search method for unconstrained optimization");CHKERRQ(ierr);
12378d5ead36SAlp 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);
12388d5ead36SAlp 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);
12398d5ead36SAlp 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);
12408d5ead36SAlp 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);
12412f75a4aaSAlp 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);
12428d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_sval", "perturbation starting value", "", bnk->sval, &bnk->sval,NULL);CHKERRQ(ierr);
12438d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_imin", "minimum initial perturbation", "", bnk->imin, &bnk->imin,NULL);CHKERRQ(ierr);
12448d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_imax", "maximum initial perturbation", "", bnk->imax, &bnk->imax,NULL);CHKERRQ(ierr);
12458d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_imfac", "initial merit factor", "", bnk->imfac, &bnk->imfac,NULL);CHKERRQ(ierr);
12468d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_pmin", "minimum perturbation", "", bnk->pmin, &bnk->pmin,NULL);CHKERRQ(ierr);
12478d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_pmax", "maximum perturbation", "", bnk->pmax, &bnk->pmax,NULL);CHKERRQ(ierr);
12488d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_pgfac", "growth factor", "", bnk->pgfac, &bnk->pgfac,NULL);CHKERRQ(ierr);
12498d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_psfac", "shrink factor", "", bnk->psfac, &bnk->psfac,NULL);CHKERRQ(ierr);
12508d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_pmgfac", "merit growth factor", "", bnk->pmgfac, &bnk->pmgfac,NULL);CHKERRQ(ierr);
12518d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_pmsfac", "merit shrink factor", "", bnk->pmsfac, &bnk->pmsfac,NULL);CHKERRQ(ierr);
12528d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_eta1", "poor steplength; reduce radius", "", bnk->eta1, &bnk->eta1,NULL);CHKERRQ(ierr);
12538d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_eta2", "reasonable steplength; leave radius alone", "", bnk->eta2, &bnk->eta2,NULL);CHKERRQ(ierr);
12548d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_eta3", "good steplength; increase radius", "", bnk->eta3, &bnk->eta3,NULL);CHKERRQ(ierr);
12558d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_eta4", "excellent steplength; greatly increase radius", "", bnk->eta4, &bnk->eta4,NULL);CHKERRQ(ierr);
12568d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_alpha1", "", "", bnk->alpha1, &bnk->alpha1,NULL);CHKERRQ(ierr);
12578d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_alpha2", "", "", bnk->alpha2, &bnk->alpha2,NULL);CHKERRQ(ierr);
12588d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_alpha3", "", "", bnk->alpha3, &bnk->alpha3,NULL);CHKERRQ(ierr);
12598d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_alpha4", "", "", bnk->alpha4, &bnk->alpha4,NULL);CHKERRQ(ierr);
12608d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_alpha5", "", "", bnk->alpha5, &bnk->alpha5,NULL);CHKERRQ(ierr);
12618d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_nu1", "poor steplength; reduce radius", "", bnk->nu1, &bnk->nu1,NULL);CHKERRQ(ierr);
12628d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_nu2", "reasonable steplength; leave radius alone", "", bnk->nu2, &bnk->nu2,NULL);CHKERRQ(ierr);
12638d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_nu3", "good steplength; increase radius", "", bnk->nu3, &bnk->nu3,NULL);CHKERRQ(ierr);
12648d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_nu4", "excellent steplength; greatly increase radius", "", bnk->nu4, &bnk->nu4,NULL);CHKERRQ(ierr);
12658d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_omega1", "", "", bnk->omega1, &bnk->omega1,NULL);CHKERRQ(ierr);
12668d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_omega2", "", "", bnk->omega2, &bnk->omega2,NULL);CHKERRQ(ierr);
12678d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_omega3", "", "", bnk->omega3, &bnk->omega3,NULL);CHKERRQ(ierr);
12688d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_omega4", "", "", bnk->omega4, &bnk->omega4,NULL);CHKERRQ(ierr);
12698d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_omega5", "", "", bnk->omega5, &bnk->omega5,NULL);CHKERRQ(ierr);
12708d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_mu1_i", "", "", bnk->mu1_i, &bnk->mu1_i,NULL);CHKERRQ(ierr);
12718d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_mu2_i", "", "", bnk->mu2_i, &bnk->mu2_i,NULL);CHKERRQ(ierr);
12728d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_gamma1_i", "", "", bnk->gamma1_i, &bnk->gamma1_i,NULL);CHKERRQ(ierr);
12738d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_gamma2_i", "", "", bnk->gamma2_i, &bnk->gamma2_i,NULL);CHKERRQ(ierr);
12748d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_gamma3_i", "", "", bnk->gamma3_i, &bnk->gamma3_i,NULL);CHKERRQ(ierr);
12758d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_gamma4_i", "", "", bnk->gamma4_i, &bnk->gamma4_i,NULL);CHKERRQ(ierr);
12768d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_theta_i", "", "", bnk->theta_i, &bnk->theta_i,NULL);CHKERRQ(ierr);
12778d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_mu1", "", "", bnk->mu1, &bnk->mu1,NULL);CHKERRQ(ierr);
12788d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_mu2", "", "", bnk->mu2, &bnk->mu2,NULL);CHKERRQ(ierr);
12798d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_gamma1", "", "", bnk->gamma1, &bnk->gamma1,NULL);CHKERRQ(ierr);
12808d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_gamma2", "", "", bnk->gamma2, &bnk->gamma2,NULL);CHKERRQ(ierr);
12818d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_gamma3", "", "", bnk->gamma3, &bnk->gamma3,NULL);CHKERRQ(ierr);
12828d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_gamma4", "", "", bnk->gamma4, &bnk->gamma4,NULL);CHKERRQ(ierr);
12838d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_theta", "", "", bnk->theta, &bnk->theta,NULL);CHKERRQ(ierr);
12848d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_min_radius", "lower bound on initial radius", "", bnk->min_radius, &bnk->min_radius,NULL);CHKERRQ(ierr);
12858d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_max_radius", "upper bound on radius", "", bnk->max_radius, &bnk->max_radius,NULL);CHKERRQ(ierr);
12868d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_epsilon", "tolerance used when computing actual and predicted reduction", "", bnk->epsilon, &bnk->epsilon,NULL);CHKERRQ(ierr);
12870a4511e9SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_as_tol", "initial tolerance used when estimating actively bounded variables", "", bnk->as_tol, &bnk->as_tol,NULL);CHKERRQ(ierr);
12880a4511e9SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_as_step", "step length used when estimating actively bounded variables", "", bnk->as_step, &bnk->as_step,NULL);CHKERRQ(ierr);
1289c0f10754SAlp 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);
1290eb910715SAlp Dener   ierr = PetscOptionsTail();CHKERRQ(ierr);
1291e031d6f5SAlp Dener   ierr = TaoSetFromOptions(bnk->bncg);
1292eb910715SAlp Dener   ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr);
1293eb910715SAlp Dener   ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
1294eb910715SAlp Dener   PetscFunctionReturn(0);
1295eb910715SAlp Dener }
1296eb910715SAlp Dener 
1297eb910715SAlp Dener /*------------------------------------------------------------*/
1298df278d8fSAlp Dener 
1299eb910715SAlp Dener static PetscErrorCode TaoView_BNK(Tao tao, PetscViewer viewer)
1300eb910715SAlp Dener {
1301eb910715SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
1302eb910715SAlp Dener   PetscInt       nrejects;
1303eb910715SAlp Dener   PetscBool      isascii;
1304eb910715SAlp Dener   PetscErrorCode ierr;
1305eb910715SAlp Dener 
1306eb910715SAlp Dener   PetscFunctionBegin;
1307eb910715SAlp Dener   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
1308eb910715SAlp Dener   if (isascii) {
1309eb910715SAlp Dener     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1310eb910715SAlp Dener     if (BNK_PC_BFGS == bnk->pc_type && bnk->M) {
1311eb910715SAlp Dener       ierr = MatLMVMGetRejects(bnk->M,&nrejects);CHKERRQ(ierr);
1312eb910715SAlp Dener       ierr = PetscViewerASCIIPrintf(viewer, "Rejected matrix updates: %D\n",nrejects);CHKERRQ(ierr);
1313eb910715SAlp Dener     }
1314e031d6f5SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "CG steps: %D\n", bnk->tot_cg_its);CHKERRQ(ierr);
1315eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "Newton steps: %D\n", bnk->newt);CHKERRQ(ierr);
1316eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "BFGS steps: %D\n", bnk->bfgs);CHKERRQ(ierr);
1317eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "Scaled gradient steps: %D\n", bnk->sgrad);CHKERRQ(ierr);
1318eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", bnk->grad);CHKERRQ(ierr);
1319eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "KSP termination reasons:\n");CHKERRQ(ierr);
1320eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  atol: %D\n", bnk->ksp_atol);CHKERRQ(ierr);
1321eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  rtol: %D\n", bnk->ksp_rtol);CHKERRQ(ierr);
1322eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  ctol: %D\n", bnk->ksp_ctol);CHKERRQ(ierr);
1323eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  negc: %D\n", bnk->ksp_negc);CHKERRQ(ierr);
1324eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  dtol: %D\n", bnk->ksp_dtol);CHKERRQ(ierr);
1325eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  iter: %D\n", bnk->ksp_iter);CHKERRQ(ierr);
1326eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  othr: %D\n", bnk->ksp_othr);CHKERRQ(ierr);
1327eb910715SAlp Dener     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1328eb910715SAlp Dener   }
1329eb910715SAlp Dener   PetscFunctionReturn(0);
1330eb910715SAlp Dener }
1331eb910715SAlp Dener 
1332eb910715SAlp Dener /* ---------------------------------------------------------- */
1333df278d8fSAlp Dener 
1334eb910715SAlp Dener /*MC
1335eb910715SAlp Dener   TAOBNK - Shared base-type for Bounded Newton-Krylov type algorithms.
133666ed3702SAlp Dener   At each iteration, the BNK methods solve the symmetric
1337eb910715SAlp Dener   system of equations to obtain the step diretion dk:
1338eb910715SAlp Dener               Hk dk = -gk
13392b97c8d8SAlp Dener   for free variables only. The step can be globalized either through
13402b97c8d8SAlp Dener   trust-region methods, or a line search, or a heuristic mixture of both.
1341eb910715SAlp Dener 
1342eb910715SAlp Dener     Options Database Keys:
13438d5ead36SAlp Dener + -tao_bnk_pc_type - "none","ahess","bfgs","petsc"
13448d5ead36SAlp Dener . -tao_bnk_bfgs_scale_type - "ahess","phess","bfgs"
13458d5ead36SAlp Dener . -tao_bnk_init_type - "constant","direction","interpolation"
13468d5ead36SAlp Dener . -tao_bnk_update_type - "step","direction","interpolation"
13472f75a4aaSAlp Dener . -tao_bnk_as_type - "none","bertsekas"
13488d5ead36SAlp Dener . -tao_bnk_sval - perturbation starting value
13498d5ead36SAlp Dener . -tao_bnk_imin - minimum initial perturbation
13508d5ead36SAlp Dener . -tao_bnk_imax - maximum initial perturbation
13518d5ead36SAlp Dener . -tao_bnk_pmin - minimum perturbation
13528d5ead36SAlp Dener . -tao_bnk_pmax - maximum perturbation
13538d5ead36SAlp Dener . -tao_bnk_pgfac - growth factor
13548d5ead36SAlp Dener . -tao_bnk_psfac - shrink factor
13558d5ead36SAlp Dener . -tao_bnk_imfac - initial merit factor
13568d5ead36SAlp Dener . -tao_bnk_pmgfac - merit growth factor
13578d5ead36SAlp Dener . -tao_bnk_pmsfac - merit shrink factor
13588d5ead36SAlp Dener . -tao_bnk_eta1 - poor steplength; reduce radius
13598d5ead36SAlp Dener . -tao_bnk_eta2 - reasonable steplength; leave radius
13608d5ead36SAlp Dener . -tao_bnk_eta3 - good steplength; increase readius
13618d5ead36SAlp Dener . -tao_bnk_eta4 - excellent steplength; greatly increase radius
13628d5ead36SAlp Dener . -tao_bnk_alpha1 - alpha1 reduction
13638d5ead36SAlp Dener . -tao_bnk_alpha2 - alpha2 reduction
13648d5ead36SAlp Dener . -tao_bnk_alpha3 - alpha3 reduction
13658d5ead36SAlp Dener . -tao_bnk_alpha4 - alpha4 reduction
13668d5ead36SAlp Dener . -tao_bnk_alpha - alpha5 reduction
13678d5ead36SAlp Dener . -tao_bnk_mu1 - mu1 interpolation update
13688d5ead36SAlp Dener . -tao_bnk_mu2 - mu2 interpolation update
13698d5ead36SAlp Dener . -tao_bnk_gamma1 - gamma1 interpolation update
13708d5ead36SAlp Dener . -tao_bnk_gamma2 - gamma2 interpolation update
13718d5ead36SAlp Dener . -tao_bnk_gamma3 - gamma3 interpolation update
13728d5ead36SAlp Dener . -tao_bnk_gamma4 - gamma4 interpolation update
13738d5ead36SAlp Dener . -tao_bnk_theta - theta interpolation update
13748d5ead36SAlp Dener . -tao_bnk_omega1 - omega1 step update
13758d5ead36SAlp Dener . -tao_bnk_omega2 - omega2 step update
13768d5ead36SAlp Dener . -tao_bnk_omega3 - omega3 step update
13778d5ead36SAlp Dener . -tao_bnk_omega4 - omega4 step update
13788d5ead36SAlp Dener . -tao_bnk_omega5 - omega5 step update
13798d5ead36SAlp Dener . -tao_bnk_mu1_i -  mu1 interpolation init factor
13808d5ead36SAlp Dener . -tao_bnk_mu2_i -  mu2 interpolation init factor
13818d5ead36SAlp Dener . -tao_bnk_gamma1_i -  gamma1 interpolation init factor
13828d5ead36SAlp Dener . -tao_bnk_gamma2_i -  gamma2 interpolation init factor
13838d5ead36SAlp Dener . -tao_bnk_gamma3_i -  gamma3 interpolation init factor
13848d5ead36SAlp Dener . -tao_bnk_gamma4_i -  gamma4 interpolation init factor
13852f75a4aaSAlp Dener . -tao_bnk_theta_i -  theta interpolation init factor
13862f75a4aaSAlp Dener - -tao_bnk_bound_tol -  initial tolerance used in estimating bounded active variables
1387eb910715SAlp Dener 
1388eb910715SAlp Dener   Level: beginner
1389eb910715SAlp Dener M*/
1390eb910715SAlp Dener 
1391eb910715SAlp Dener PetscErrorCode TaoCreate_BNK(Tao tao)
1392eb910715SAlp Dener {
1393eb910715SAlp Dener   TAO_BNK        *bnk;
1394eb910715SAlp Dener   const char     *morethuente_type = TAOLINESEARCHMT;
1395eb910715SAlp Dener   PetscErrorCode ierr;
1396eb910715SAlp Dener 
1397eb910715SAlp Dener   PetscFunctionBegin;
1398eb910715SAlp Dener   ierr = PetscNewLog(tao,&bnk);CHKERRQ(ierr);
1399eb910715SAlp Dener 
1400eb910715SAlp Dener   tao->ops->setup = TaoSetUp_BNK;
1401eb910715SAlp Dener   tao->ops->view = TaoView_BNK;
1402eb910715SAlp Dener   tao->ops->setfromoptions = TaoSetFromOptions_BNK;
1403eb910715SAlp Dener   tao->ops->destroy = TaoDestroy_BNK;
1404eb910715SAlp Dener 
1405eb910715SAlp Dener   /*  Override default settings (unless already changed) */
1406eb910715SAlp Dener   if (!tao->max_it_changed) tao->max_it = 50;
1407eb910715SAlp Dener   if (!tao->trust0_changed) tao->trust0 = 100.0;
1408eb910715SAlp Dener 
1409eb910715SAlp Dener   tao->data = (void*)bnk;
1410eb910715SAlp Dener 
141166ed3702SAlp Dener   /*  Hessian shifting parameters */
1412eb910715SAlp Dener   bnk->sval   = 0.0;
1413eb910715SAlp Dener   bnk->imin   = 1.0e-4;
1414eb910715SAlp Dener   bnk->imax   = 1.0e+2;
1415eb910715SAlp Dener   bnk->imfac  = 1.0e-1;
1416eb910715SAlp Dener 
1417eb910715SAlp Dener   bnk->pmin   = 1.0e-12;
1418eb910715SAlp Dener   bnk->pmax   = 1.0e+2;
1419eb910715SAlp Dener   bnk->pgfac  = 1.0e+1;
1420eb910715SAlp Dener   bnk->psfac  = 4.0e-1;
1421eb910715SAlp Dener   bnk->pmgfac = 1.0e-1;
1422eb910715SAlp Dener   bnk->pmsfac = 1.0e-1;
1423eb910715SAlp Dener 
1424eb910715SAlp Dener   /*  Default values for trust-region radius update based on steplength */
1425eb910715SAlp Dener   bnk->nu1 = 0.25;
1426eb910715SAlp Dener   bnk->nu2 = 0.50;
1427eb910715SAlp Dener   bnk->nu3 = 1.00;
1428eb910715SAlp Dener   bnk->nu4 = 1.25;
1429eb910715SAlp Dener 
1430eb910715SAlp Dener   bnk->omega1 = 0.25;
1431eb910715SAlp Dener   bnk->omega2 = 0.50;
1432eb910715SAlp Dener   bnk->omega3 = 1.00;
1433eb910715SAlp Dener   bnk->omega4 = 2.00;
1434eb910715SAlp Dener   bnk->omega5 = 4.00;
1435eb910715SAlp Dener 
1436eb910715SAlp Dener   /*  Default values for trust-region radius update based on reduction */
1437eb910715SAlp Dener   bnk->eta1 = 1.0e-4;
1438eb910715SAlp Dener   bnk->eta2 = 0.25;
1439eb910715SAlp Dener   bnk->eta3 = 0.50;
1440eb910715SAlp Dener   bnk->eta4 = 0.90;
1441eb910715SAlp Dener 
1442eb910715SAlp Dener   bnk->alpha1 = 0.25;
1443eb910715SAlp Dener   bnk->alpha2 = 0.50;
1444eb910715SAlp Dener   bnk->alpha3 = 1.00;
1445eb910715SAlp Dener   bnk->alpha4 = 2.00;
1446eb910715SAlp Dener   bnk->alpha5 = 4.00;
1447eb910715SAlp Dener 
1448eb910715SAlp Dener   /*  Default values for trust-region radius update based on interpolation */
1449eb910715SAlp Dener   bnk->mu1 = 0.10;
1450eb910715SAlp Dener   bnk->mu2 = 0.50;
1451eb910715SAlp Dener 
1452eb910715SAlp Dener   bnk->gamma1 = 0.25;
1453eb910715SAlp Dener   bnk->gamma2 = 0.50;
1454eb910715SAlp Dener   bnk->gamma3 = 2.00;
1455eb910715SAlp Dener   bnk->gamma4 = 4.00;
1456eb910715SAlp Dener 
1457eb910715SAlp Dener   bnk->theta = 0.05;
1458eb910715SAlp Dener 
1459eb910715SAlp Dener   /*  Default values for trust region initialization based on interpolation */
1460eb910715SAlp Dener   bnk->mu1_i = 0.35;
1461eb910715SAlp Dener   bnk->mu2_i = 0.50;
1462eb910715SAlp Dener 
1463eb910715SAlp Dener   bnk->gamma1_i = 0.0625;
1464eb910715SAlp Dener   bnk->gamma2_i = 0.5;
1465eb910715SAlp Dener   bnk->gamma3_i = 2.0;
1466eb910715SAlp Dener   bnk->gamma4_i = 5.0;
1467eb910715SAlp Dener 
1468eb910715SAlp Dener   bnk->theta_i = 0.25;
1469eb910715SAlp Dener 
1470eb910715SAlp Dener   /*  Remaining parameters */
1471c0f10754SAlp Dener   bnk->max_cg_its = 0;
1472eb910715SAlp Dener   bnk->min_radius = 1.0e-10;
1473eb910715SAlp Dener   bnk->max_radius = 1.0e10;
1474770b7498SAlp Dener   bnk->epsilon = PetscPowReal(PETSC_MACHINE_EPSILON, 2.0/3.0);
14750a4511e9SAlp Dener   bnk->as_tol = 1.0e-3;
14760a4511e9SAlp Dener   bnk->as_step = 1.0e-3;
147762675beeSAlp Dener   bnk->dmin = 1.0e-6;
147862675beeSAlp Dener   bnk->dmax = 1.0e6;
1479eb910715SAlp Dener 
1480e031d6f5SAlp Dener   bnk->pc_type         = BNK_PC_AHESS;
1481eb910715SAlp Dener   bnk->bfgs_scale_type = BFGS_SCALE_PHESS;
1482eb910715SAlp Dener   bnk->init_type       = BNK_INIT_INTERPOLATION;
1483fed79b8eSAlp Dener   bnk->update_type     = BNK_UPDATE_INTERPOLATION;
14842f75a4aaSAlp Dener   bnk->as_type         = BNK_AS_BERTSEKAS;
1485eb910715SAlp Dener 
1486e031d6f5SAlp Dener   /* Create the embedded BNCG solver */
1487e031d6f5SAlp Dener   ierr = TaoCreate(PetscObjectComm((PetscObject)tao), &bnk->bncg);CHKERRQ(ierr);
1488e031d6f5SAlp Dener   ierr = PetscObjectIncrementTabLevel((PetscObject)bnk->bncg, (PetscObject)tao, 1);CHKERRQ(ierr);
1489e031d6f5SAlp Dener   ierr = TaoSetOptionsPrefix(bnk->bncg, "tao_bnk_");CHKERRQ(ierr);
1490e031d6f5SAlp Dener   ierr = TaoSetType(bnk->bncg, TAOBNCG);CHKERRQ(ierr);
1491e031d6f5SAlp Dener 
1492c0f10754SAlp Dener   /* Create the line search */
1493eb910715SAlp Dener   ierr = TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch);CHKERRQ(ierr);
1494eb910715SAlp Dener   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr);
1495e031d6f5SAlp Dener   ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr);
1496eb910715SAlp Dener   ierr = TaoLineSearchSetType(tao->linesearch,morethuente_type);CHKERRQ(ierr);
1497eb910715SAlp Dener   ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr);
1498eb910715SAlp Dener 
1499eb910715SAlp Dener   /*  Set linear solver to default for symmetric matrices */
1500eb910715SAlp Dener   ierr = KSPCreate(((PetscObject)tao)->comm,&tao->ksp);CHKERRQ(ierr);
1501eb910715SAlp Dener   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1);CHKERRQ(ierr);
1502eb910715SAlp Dener   ierr = KSPSetOptionsPrefix(tao->ksp,tao->hdr.prefix);CHKERRQ(ierr);
1503eb910715SAlp Dener   ierr = KSPSetType(tao->ksp,KSPCGSTCG);CHKERRQ(ierr);
1504eb910715SAlp Dener   PetscFunctionReturn(0);
1505eb910715SAlp Dener }
1506