xref: /petsc/src/tao/bound/impls/bnk/bnk.c (revision 0875260307e97403c513ee2ad2ca062f01ff5a8a)
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 
44c0f10754SAlp Dener   PetscInt                     n,N;
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 */
5728017e9fSAlp Dener   ierr = VecMedian(tao->XL,tao->solution,tao->XU,tao->solution);CHKERRQ(ierr);
5828017e9fSAlp Dener 
5928017e9fSAlp Dener   /* Check convergence criteria */
6028017e9fSAlp Dener   ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &bnk->f, bnk->unprojected_gradient);CHKERRQ(ierr);
6128017e9fSAlp Dener   ierr = VecBoundGradientProjection(bnk->unprojected_gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr);
62*08752603SAlp Dener   ierr = VecNorm(tao->gradient,NORM_2,&bnk->gnorm);CHKERRQ(ierr);
6328017e9fSAlp Dener   if (PetscIsInfOrNanReal(bnk->f) || PetscIsInfOrNanReal(bnk->gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
6428017e9fSAlp Dener 
65c0f10754SAlp Dener   /* Test the initial point for convergence */
66e031d6f5SAlp Dener   ierr = VecFischer(tao->solution, bnk->unprojected_gradient, tao->XL, tao->XU, bnk->Gwork);CHKERRQ(ierr);
67e031d6f5SAlp Dener   ierr = VecNorm(bnk->Gwork, NORM_2, &resnorm);CHKERRQ(ierr);
68e031d6f5SAlp Dener   ierr = TaoLogConvergenceHistory(tao,bnk->f,resnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
69e031d6f5SAlp Dener   ierr = TaoMonitor(tao,tao->niter,bnk->f,resnorm,0.0,1.0);CHKERRQ(ierr);
70c0f10754SAlp Dener   ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
71c0f10754SAlp Dener   if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
72c0f10754SAlp Dener 
73e031d6f5SAlp Dener   /* Reset KSP stopping reason counters */
74eb910715SAlp Dener   bnk->ksp_atol = 0;
75eb910715SAlp Dener   bnk->ksp_rtol = 0;
76eb910715SAlp Dener   bnk->ksp_dtol = 0;
77eb910715SAlp Dener   bnk->ksp_ctol = 0;
78eb910715SAlp Dener   bnk->ksp_negc = 0;
79eb910715SAlp Dener   bnk->ksp_iter = 0;
80eb910715SAlp Dener   bnk->ksp_othr = 0;
81eb910715SAlp Dener 
82e031d6f5SAlp Dener   /* Reset accepted step type counters */
83e031d6f5SAlp Dener   bnk->tot_cg_its = 0;
84e031d6f5SAlp Dener   bnk->newt = 0;
85e031d6f5SAlp Dener   bnk->bfgs = 0;
86e031d6f5SAlp Dener   bnk->sgrad = 0;
87e031d6f5SAlp Dener   bnk->grad = 0;
88e031d6f5SAlp Dener 
89fed79b8eSAlp Dener   /* Initialize the Hessian perturbation */
90fed79b8eSAlp Dener   bnk->pert = bnk->sval;
91fed79b8eSAlp Dener 
92e031d6f5SAlp Dener   /* Allocate the vectors needed for the BFGS approximation */
93e031d6f5SAlp Dener   if (BNK_PC_BFGS == bnk->pc_type) {
94e031d6f5SAlp Dener     if (!bnk->M) {
95eb910715SAlp Dener       ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr);
96eb910715SAlp Dener       ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr);
97eb910715SAlp Dener       ierr = MatCreateLMVM(((PetscObject)tao)->comm,n,N,&bnk->M);CHKERRQ(ierr);
98eb910715SAlp Dener       ierr = MatLMVMAllocateVectors(bnk->M,tao->solution);CHKERRQ(ierr);
99eb910715SAlp Dener     }
100e031d6f5SAlp Dener     if (bnk->bfgs_scale_type != BFGS_SCALE_BFGS && !bnk->Diag) {
101e031d6f5SAlp Dener       ierr = VecDuplicate(tao->solution,&bnk->Diag);CHKERRQ(ierr);
1025e9b73cbSAlp Dener     }
103e031d6f5SAlp Dener   }
104e031d6f5SAlp Dener 
105e031d6f5SAlp Dener   /* Prepare the min/max vectors for safeguarding diagonal scales */
10662675beeSAlp Dener   ierr = VecSet(bnk->Diag_min, bnk->dmin);CHKERRQ(ierr);
10762675beeSAlp Dener   ierr = VecSet(bnk->Diag_max, bnk->dmax);CHKERRQ(ierr);
108eb910715SAlp Dener 
109eb910715SAlp Dener   /* Modify the preconditioner to use the bfgs approximation */
110eb910715SAlp Dener   ierr = KSPGetPC(tao->ksp, &pc);CHKERRQ(ierr);
111eb910715SAlp Dener   switch(bnk->pc_type) {
112eb910715SAlp Dener   case BNK_PC_NONE:
113eb910715SAlp Dener     ierr = PCSetType(pc, PCNONE);CHKERRQ(ierr);
114eb910715SAlp Dener     ierr = PCSetFromOptions(pc);CHKERRQ(ierr);
115eb910715SAlp Dener     break;
116eb910715SAlp Dener 
117eb910715SAlp Dener   case BNK_PC_AHESS:
118eb910715SAlp Dener     ierr = PCSetType(pc, PCJACOBI);CHKERRQ(ierr);
119eb910715SAlp Dener     ierr = PCSetFromOptions(pc);CHKERRQ(ierr);
120eb910715SAlp Dener     ierr = PCJacobiSetUseAbs(pc,PETSC_TRUE);CHKERRQ(ierr);
121eb910715SAlp Dener     break;
122eb910715SAlp Dener 
123eb910715SAlp Dener   case BNK_PC_BFGS:
124eb910715SAlp Dener     ierr = PCSetType(pc, PCSHELL);CHKERRQ(ierr);
125eb910715SAlp Dener     ierr = PCSetFromOptions(pc);CHKERRQ(ierr);
126eb910715SAlp Dener     ierr = PCShellSetName(pc, "bfgs");CHKERRQ(ierr);
127eb910715SAlp Dener     ierr = PCShellSetContext(pc, bnk->M);CHKERRQ(ierr);
128eb910715SAlp Dener     ierr = PCShellSetApply(pc, MatLMVMSolveShell);CHKERRQ(ierr);
129eb910715SAlp Dener     break;
130eb910715SAlp Dener 
131eb910715SAlp Dener   default:
132eb910715SAlp Dener     /* Use the pc method set by pc_type */
133eb910715SAlp Dener     break;
134eb910715SAlp Dener   }
135eb910715SAlp Dener 
136eb910715SAlp Dener   /* Initialize trust-region radius.  The initialization is only performed
137eb910715SAlp Dener      when we are using Nash, Steihaug-Toint or the Generalized Lanczos method. */
138c0f10754SAlp Dener   *needH = PETSC_TRUE;
139eb910715SAlp Dener   if (bnk->is_nash || bnk->is_stcg || bnk->is_gltr) {
14062675beeSAlp Dener     switch(initType) {
141eb910715SAlp Dener     case BNK_INIT_CONSTANT:
142eb910715SAlp Dener       /* Use the initial radius specified */
143c0f10754SAlp Dener       tao->trust = tao->trust0;
144eb910715SAlp Dener       break;
145eb910715SAlp Dener 
146eb910715SAlp Dener     case BNK_INIT_INTERPOLATION:
147c0f10754SAlp Dener       /* Use interpolation based on the initial Hessian */
148eb910715SAlp Dener       max_radius = 0.0;
149*08752603SAlp Dener       tao->trust = tao->trust0;
150eb910715SAlp Dener       for (j = 0; j < j_max; ++j) {
1510a4511e9SAlp Dener         f_min = bnk->f;
152eb910715SAlp Dener         sigma = 0.0;
153eb910715SAlp Dener 
154c0f10754SAlp Dener         if (*needH) {
15562602cfbSAlp Dener           /* Compute the Hessian at the new step, and extract the inactive subsystem */
156e031d6f5SAlp Dener           ierr = TaoBNKComputeHessian(tao);CHKERRQ(ierr);
157*08752603SAlp Dener           ierr = TaoBNKEstimateActiveSet(tao, BNK_AS_NONE);CHKERRQ(ierr);
15828017e9fSAlp Dener           if (bnk->inactive_idx) {
1592ab2a32cSAlp Dener             ierr = MatCreateSubMatrix(tao->hessian, bnk->inactive_idx, bnk->inactive_idx, MAT_INITIAL_MATRIX, &bnk->H_inactive);CHKERRQ(ierr);
16028017e9fSAlp Dener           } else {
161*08752603SAlp Dener             ierr = MatDuplicate(tao->hessian, MAT_COPY_VALUES, &bnk->H_inactive);CHKERRQ(ierr);
16228017e9fSAlp Dener           }
163c0f10754SAlp Dener           *needH = PETSC_FALSE;
164eb910715SAlp Dener         }
165eb910715SAlp Dener 
166eb910715SAlp Dener         for (i = 0; i < i_max; ++i) {
16762602cfbSAlp Dener           /* Take a steepest descent step and snap it to bounds */
16862602cfbSAlp Dener           ierr = VecCopy(tao->solution, bnk->Xold);CHKERRQ(ierr);
16962602cfbSAlp Dener           ierr = VecAXPY(tao->solution, -tao->trust/bnk->gnorm, tao->gradient);CHKERRQ(ierr);
17062602cfbSAlp Dener           ierr = VecMedian(tao->XL, tao->solution, tao->XU, tao->solution);CHKERRQ(ierr);
171770b7498SAlp Dener           /* Recompute the step after bound snapping so that it can be used in predicted decrease calculation later */
172eb910715SAlp Dener           ierr = VecCopy(tao->solution, bnk->W);CHKERRQ(ierr);
17362602cfbSAlp Dener           ierr = VecAXPY(bnk->W, -1.0, bnk->Xold);CHKERRQ(ierr);
17462602cfbSAlp Dener           /* Compute the objective at the trial */
17562602cfbSAlp Dener           ierr = TaoComputeObjective(tao, tao->solution, &ftrial);CHKERRQ(ierr);
17662602cfbSAlp Dener           ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr);
177eb910715SAlp Dener           if (PetscIsInfOrNanReal(ftrial)) {
178eb910715SAlp Dener             tau = bnk->gamma1_i;
179eb910715SAlp Dener           } else {
1800a4511e9SAlp Dener             if (ftrial < f_min) {
1810a4511e9SAlp Dener               f_min = ftrial;
182eb910715SAlp Dener               sigma = -tao->trust / bnk->gnorm;
183eb910715SAlp Dener             }
184*08752603SAlp Dener 
185770b7498SAlp Dener             /* Compute the predicted and actual reduction */
1862ab2a32cSAlp Dener             if (bnk->inactive_idx) {
187*08752603SAlp Dener               ierr = VecGetSubVector(bnk->W, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr);
188*08752603SAlp Dener               ierr = VecGetSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work);CHKERRQ(ierr);
1892ab2a32cSAlp Dener             } else {
190*08752603SAlp Dener               bnk->X_inactive = bnk->W;
191*08752603SAlp Dener               bnk->inactive_work = bnk->Xwork;
1922ab2a32cSAlp Dener             }
193*08752603SAlp Dener             ierr = MatMult(bnk->H_inactive, bnk->X_inactive, bnk->inactive_work);CHKERRQ(ierr);
194*08752603SAlp Dener             ierr = VecDot(bnk->X_inactive, bnk->inactive_work, &prered);CHKERRQ(ierr);
1952ab2a32cSAlp Dener             if (bnk->inactive_idx) {
196*08752603SAlp Dener               ierr = VecRestoreSubVector(bnk->W, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr);
197*08752603SAlp Dener               ierr = VecRestoreSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work);CHKERRQ(ierr);
1982ab2a32cSAlp Dener             }
199eb910715SAlp Dener             prered = tao->trust * (bnk->gnorm - 0.5 * tao->trust * prered / (bnk->gnorm * bnk->gnorm));
200eb910715SAlp Dener             actred = bnk->f - ftrial;
201*08752603SAlp Dener             if ((PetscAbsScalar(actred) <= bnk->epsilon) &&
202*08752603SAlp Dener                 (PetscAbsScalar(prered) <= bnk->epsilon)) {
203eb910715SAlp Dener               kappa = 1.0;
204*08752603SAlp Dener             }
205*08752603SAlp Dener             else {
206eb910715SAlp Dener               kappa = actred / prered;
207eb910715SAlp Dener             }
208eb910715SAlp Dener 
209eb910715SAlp Dener             tau_1 = bnk->theta_i * bnk->gnorm * tao->trust / (bnk->theta_i * bnk->gnorm * tao->trust + (1.0 - bnk->theta_i) * prered - actred);
210eb910715SAlp Dener             tau_2 = bnk->theta_i * bnk->gnorm * tao->trust / (bnk->theta_i * bnk->gnorm * tao->trust - (1.0 + bnk->theta_i) * prered + actred);
211eb910715SAlp Dener             tau_min = PetscMin(tau_1, tau_2);
212eb910715SAlp Dener             tau_max = PetscMax(tau_1, tau_2);
213eb910715SAlp Dener 
214eb910715SAlp Dener             if (PetscAbsScalar(kappa - 1.0) <= bnk->mu1_i) {
215eb910715SAlp Dener               /*  Great agreement */
216eb910715SAlp Dener               max_radius = PetscMax(max_radius, tao->trust);
217eb910715SAlp Dener 
218eb910715SAlp Dener               if (tau_max < 1.0) {
219eb910715SAlp Dener                 tau = bnk->gamma3_i;
220*08752603SAlp Dener               }
221*08752603SAlp Dener               else if (tau_max > bnk->gamma4_i) {
222eb910715SAlp Dener                 tau = bnk->gamma4_i;
223*08752603SAlp Dener               }
224*08752603SAlp Dener               else {
225eb910715SAlp Dener                 tau = tau_max;
226eb910715SAlp Dener               }
227*08752603SAlp Dener             }
228*08752603SAlp Dener             else if (PetscAbsScalar(kappa - 1.0) <= bnk->mu2_i) {
229eb910715SAlp Dener               /*  Good agreement */
230eb910715SAlp Dener               max_radius = PetscMax(max_radius, tao->trust);
231eb910715SAlp Dener 
232eb910715SAlp Dener               if (tau_max < bnk->gamma2_i) {
233eb910715SAlp Dener                 tau = bnk->gamma2_i;
234eb910715SAlp Dener               } else if (tau_max > bnk->gamma3_i) {
235eb910715SAlp Dener                 tau = bnk->gamma3_i;
236eb910715SAlp Dener               } else {
237eb910715SAlp Dener                 tau = tau_max;
238eb910715SAlp Dener               }
239*08752603SAlp Dener             }
240*08752603SAlp Dener             else {
241eb910715SAlp Dener               /*  Not good agreement */
242eb910715SAlp Dener               if (tau_min > 1.0) {
243eb910715SAlp Dener                 tau = bnk->gamma2_i;
244eb910715SAlp Dener               } else if (tau_max < bnk->gamma1_i) {
245eb910715SAlp Dener                 tau = bnk->gamma1_i;
246eb910715SAlp Dener               } else if ((tau_min < bnk->gamma1_i) && (tau_max >= 1.0)) {
247eb910715SAlp Dener                 tau = bnk->gamma1_i;
248*08752603SAlp Dener               } else if ((tau_1 >= bnk->gamma1_i) && (tau_1 < 1.0) &&
249*08752603SAlp Dener                         ((tau_2 < bnk->gamma1_i) || (tau_2 >= 1.0))) {
250eb910715SAlp Dener                 tau = tau_1;
251*08752603SAlp Dener               } else if ((tau_2 >= bnk->gamma1_i) && (tau_2 < 1.0) &&
252*08752603SAlp Dener                         ((tau_1 < bnk->gamma1_i) || (tau_2 >= 1.0))) {
253eb910715SAlp Dener                 tau = tau_2;
254eb910715SAlp Dener               } else {
255eb910715SAlp Dener                 tau = tau_max;
256eb910715SAlp Dener               }
257eb910715SAlp Dener             }
258eb910715SAlp Dener           }
259eb910715SAlp Dener           tao->trust = tau * tao->trust;
260eb910715SAlp Dener         }
261eb910715SAlp Dener 
2620a4511e9SAlp Dener         if (f_min < bnk->f) {
263e031d6f5SAlp Dener           /* We found a solution better than the initial, so let's test and accept it */
2640a4511e9SAlp Dener           bnk->f = f_min;
265eb910715SAlp Dener           ierr = VecAXPY(tao->solution,sigma,tao->gradient);CHKERRQ(ierr);
26687602d1eSAlp Dener           ierr = VecMedian(tao->XL, tao->solution, tao->XU, tao->solution);CHKERRQ(ierr);
26709164190SAlp Dener           ierr = TaoComputeGradient(tao,tao->solution,bnk->unprojected_gradient);CHKERRQ(ierr);
26809164190SAlp Dener           ierr = VecBoundGradientProjection(bnk->unprojected_gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr);
269eb910715SAlp Dener 
270eb910715SAlp Dener           ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&bnk->gnorm);CHKERRQ(ierr);
271eb910715SAlp Dener           if (PetscIsInfOrNanReal(bnk->gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute gradient generated Inf or NaN");
272c0f10754SAlp Dener           *needH = PETSC_TRUE;
273eb910715SAlp Dener 
274e031d6f5SAlp Dener           ierr = VecFischer(tao->solution, bnk->unprojected_gradient, tao->XL, tao->XU, bnk->Gwork);CHKERRQ(ierr);
275e031d6f5SAlp Dener           ierr = VecNorm(bnk->Gwork, NORM_2, &resnorm);CHKERRQ(ierr);
276e031d6f5SAlp Dener           ierr = TaoLogConvergenceHistory(tao,bnk->f,resnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
277e031d6f5SAlp Dener           ierr = TaoMonitor(tao,tao->niter,bnk->f,resnorm,0.0,1.0);CHKERRQ(ierr);
278eb910715SAlp Dener           ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
279eb910715SAlp Dener           if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
280eb910715SAlp Dener         }
281eb910715SAlp Dener       }
282eb910715SAlp Dener       tao->trust = PetscMax(tao->trust, max_radius);
283e031d6f5SAlp Dener 
284e031d6f5SAlp Dener       /* Ensure that the trust radius is within the limits */
285e031d6f5SAlp Dener       tao->trust = PetscMax(tao->trust, bnk->min_radius);
286e031d6f5SAlp Dener       tao->trust = PetscMin(tao->trust, bnk->max_radius);
287eb910715SAlp Dener       break;
288eb910715SAlp Dener 
289eb910715SAlp Dener     default:
290eb910715SAlp Dener       /* Norm of the first direction will initialize radius */
291eb910715SAlp Dener       tao->trust = 0.0;
292eb910715SAlp Dener       break;
293eb910715SAlp Dener     }
294eb910715SAlp Dener   }
295e031d6f5SAlp Dener 
296eb910715SAlp Dener   /* Set initial scaling for the BFGS preconditioner
297eb910715SAlp Dener      This step is done after computing the initial trust-region radius
298eb910715SAlp Dener      since the function value may have decreased */
299eb910715SAlp Dener   if (BNK_PC_BFGS == bnk->pc_type) {
3009b6ef848SAlp Dener     delta = 2.0 * PetscMax(1.0, PetscAbsScalar(bnk->f)) / (bnk->gnorm*bnk->gnorm);
301eb910715SAlp Dener     ierr = MatLMVMSetDelta(bnk->M,delta);CHKERRQ(ierr);
302eb910715SAlp Dener   }
303eb910715SAlp Dener   PetscFunctionReturn(0);
304eb910715SAlp Dener }
305eb910715SAlp Dener 
306df278d8fSAlp Dener /*------------------------------------------------------------*/
307df278d8fSAlp Dener 
30862675beeSAlp Dener /* Routine for computing the Hessian and preparing the preconditioner at the new iterate */
30962675beeSAlp Dener 
31062675beeSAlp Dener PetscErrorCode TaoBNKComputeHessian(Tao tao)
31162675beeSAlp Dener {
31262675beeSAlp Dener   PetscErrorCode               ierr;
31362675beeSAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
31462675beeSAlp Dener 
31562675beeSAlp Dener   PetscFunctionBegin;
31662675beeSAlp Dener   /* Compute the Hessian */
31762675beeSAlp Dener   ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr);
31862675beeSAlp Dener   /* Add a correction to the BFGS preconditioner */
31962675beeSAlp Dener   if (BNK_PC_BFGS == bnk->pc_type) {
32062675beeSAlp Dener     ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
32162675beeSAlp Dener     /* Update the BFGS diagonal scaling */
32262675beeSAlp Dener     if (BFGS_SCALE_AHESS == bnk->bfgs_scale_type) {
32362675beeSAlp Dener       ierr = MatGetDiagonal(tao->hessian, bnk->Diag);CHKERRQ(ierr);
32462675beeSAlp Dener       ierr = VecAbs(bnk->Diag);CHKERRQ(ierr);
32562675beeSAlp Dener       ierr = VecMedian(bnk->Diag_min, bnk->Diag, bnk->Diag_max, bnk->Diag);CHKERRQ(ierr);
32662675beeSAlp Dener       ierr = VecReciprocal(bnk->Diag);CHKERRQ(ierr);
32762675beeSAlp Dener       ierr = MatLMVMSetScale(bnk->M,bnk->Diag);CHKERRQ(ierr);
32862675beeSAlp Dener     }
32962675beeSAlp Dener   }
33062675beeSAlp Dener   PetscFunctionReturn(0);
33162675beeSAlp Dener }
33262675beeSAlp Dener 
33362675beeSAlp Dener /*------------------------------------------------------------*/
33462675beeSAlp Dener 
3352f75a4aaSAlp Dener /* Routine for estimating the active set */
3362f75a4aaSAlp Dener 
337*08752603SAlp Dener PetscErrorCode TaoBNKEstimateActiveSet(Tao tao, PetscInt asType)
3382f75a4aaSAlp Dener {
3392f75a4aaSAlp Dener   PetscErrorCode               ierr;
3402f75a4aaSAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
3412f75a4aaSAlp Dener 
3422f75a4aaSAlp Dener   PetscFunctionBegin;
343*08752603SAlp Dener   switch (asType) {
3442f75a4aaSAlp Dener   case BNK_AS_NONE:
3452f75a4aaSAlp Dener     ierr = ISDestroy(&bnk->inactive_idx);CHKERRQ(ierr);
3462f75a4aaSAlp Dener     ierr = VecWhichInactive(tao->XL, tao->solution, bnk->unprojected_gradient, tao->XU, PETSC_TRUE, &bnk->inactive_idx);CHKERRQ(ierr);
3472f75a4aaSAlp Dener     ierr = ISDestroy(&bnk->active_idx);CHKERRQ(ierr);
3482f75a4aaSAlp Dener     ierr = ISComplementVec(bnk->inactive_idx, tao->solution, &bnk->active_idx);CHKERRQ(ierr);
3492f75a4aaSAlp Dener     break;
3502f75a4aaSAlp Dener 
3512f75a4aaSAlp Dener   case BNK_AS_BERTSEKAS:
3522f75a4aaSAlp Dener     /* Compute the trial step vector with which we will estimate the active set at the next iteration */
3532f75a4aaSAlp Dener     if (BNK_PC_BFGS == bnk->pc_type) {
3542f75a4aaSAlp Dener       /* If the BFGS preconditioner matrix is available, we will construct a trial step with it */
3555e9b73cbSAlp Dener       ierr = MatLMVMSetInactive(bnk->M, NULL);CHKERRQ(ierr);
3562f75a4aaSAlp Dener       ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, bnk->W);CHKERRQ(ierr);
3572f75a4aaSAlp Dener     } else {
3589b6ef848SAlp Dener       /* BFGS preconditioner doesn't exist so let's invert the absolute diagonal of the Hessian instead onto the gradient */
3592f75a4aaSAlp Dener       ierr = MatGetDiagonal(tao->hessian, bnk->Xwork);CHKERRQ(ierr);
3609b6ef848SAlp Dener       ierr = VecAbs(bnk->Xwork);CHKERRQ(ierr);
3619b6ef848SAlp Dener       ierr = VecMedian(bnk->Diag_min, bnk->Xwork, bnk->Diag_max, bnk->Xwork);CHKERRQ(ierr);
3622f75a4aaSAlp Dener       ierr = VecReciprocal(bnk->Xwork);CHKERRQ(ierr);CHKERRQ(ierr);
3632f75a4aaSAlp Dener       ierr = VecPointwiseMult(bnk->W, bnk->Xwork, bnk->unprojected_gradient);CHKERRQ(ierr);
3642f75a4aaSAlp Dener     }
3652f75a4aaSAlp Dener     ierr = VecScale(bnk->W, -1.0);CHKERRQ(ierr);
3660a4511e9SAlp Dener     ierr = TaoEstimateActiveBounds(tao->solution, tao->XL, tao->XU, bnk->unprojected_gradient, bnk->W, bnk->as_step, &bnk->as_tol, &bnk->active_lower, &bnk->active_upper, &bnk->active_fixed, &bnk->active_idx, &bnk->inactive_idx);CHKERRQ(ierr);
3672f75a4aaSAlp Dener 
3682f75a4aaSAlp Dener   default:
3692f75a4aaSAlp Dener     break;
3702f75a4aaSAlp Dener   }
3712f75a4aaSAlp Dener   PetscFunctionReturn(0);
3722f75a4aaSAlp Dener }
3732f75a4aaSAlp Dener 
3742f75a4aaSAlp Dener /*------------------------------------------------------------*/
3752f75a4aaSAlp Dener 
3762f75a4aaSAlp Dener /* Routine for bounding the step direction */
3772f75a4aaSAlp Dener 
3782f75a4aaSAlp Dener PetscErrorCode TaoBNKBoundStep(Tao tao, Vec step)
3792f75a4aaSAlp Dener {
3802f75a4aaSAlp Dener   PetscErrorCode               ierr;
3812f75a4aaSAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
3822f75a4aaSAlp Dener 
3832f75a4aaSAlp Dener   PetscFunctionBegin;
38428017e9fSAlp Dener   if (bnk->active_idx) {
3852f75a4aaSAlp Dener     switch (bnk->as_type) {
3862f75a4aaSAlp Dener     case BNK_AS_NONE:
3872f75a4aaSAlp Dener       if (bnk->active_idx) {
3882f75a4aaSAlp Dener         ierr = VecGetSubVector(step, bnk->active_idx, &bnk->active_work);CHKERRQ(ierr);
3892f75a4aaSAlp Dener         ierr = VecSet(bnk->active_work, 0.0);CHKERRQ(ierr);
3902f75a4aaSAlp Dener         ierr = VecRestoreSubVector(step, bnk->active_idx, &bnk->active_work);CHKERRQ(ierr);
3912f75a4aaSAlp Dener       }
3922f75a4aaSAlp Dener       break;
3932f75a4aaSAlp Dener 
3942f75a4aaSAlp Dener     case BNK_AS_BERTSEKAS:
3952f75a4aaSAlp Dener       ierr = TaoBoundStep(tao->solution, tao->XL, tao->XU, bnk->active_lower, bnk->active_upper, bnk->active_fixed, step);CHKERRQ(ierr);
3962f75a4aaSAlp Dener       break;
3972f75a4aaSAlp Dener 
3982f75a4aaSAlp Dener     default:
3992f75a4aaSAlp Dener       break;
4002f75a4aaSAlp Dener     }
401e465cd6fSAlp Dener   }
4022f75a4aaSAlp Dener   PetscFunctionReturn(0);
4032f75a4aaSAlp Dener }
4042f75a4aaSAlp Dener 
405e031d6f5SAlp Dener /*------------------------------------------------------------*/
406e031d6f5SAlp Dener 
407e031d6f5SAlp Dener /* Routine for taking a finite number of BNCG iterations to
408e031d6f5SAlp Dener    accelerate Newton convergence.
409e031d6f5SAlp Dener 
410e031d6f5SAlp Dener    In practice, this approach simply trades off Hessian evaluations
411e031d6f5SAlp Dener    for more gradient evaluations.
412e031d6f5SAlp Dener */
413e031d6f5SAlp Dener 
414c0f10754SAlp Dener PetscErrorCode TaoBNKTakeCGSteps(Tao tao, PetscBool *terminate)
415c0f10754SAlp Dener {
416c0f10754SAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
417c0f10754SAlp Dener   PetscErrorCode               ierr;
418c0f10754SAlp Dener 
419c0f10754SAlp Dener   PetscFunctionBegin;
420c0f10754SAlp Dener   *terminate = PETSC_FALSE;
421c0f10754SAlp Dener   if (bnk->max_cg_its > 0) {
422c0f10754SAlp Dener     /* Copy the current solution, unprojected gradient and step info into BNCG */
423c0f10754SAlp Dener     ierr = VecCopy(tao->solution, bnk->bncg->solution);CHKERRQ(ierr);
424c0f10754SAlp Dener     if (tao->niter > 1) {
425c0f10754SAlp Dener       bnk->bncg_ctx->f = bnk->f;
426c0f10754SAlp Dener       ierr = VecCopy(bnk->unprojected_gradient, bnk->bncg_ctx->unprojected_gradient);CHKERRQ(ierr);
427c0f10754SAlp Dener       ierr = VecCopy(tao->stepdirection, bnk->bncg->stepdirection);CHKERRQ(ierr);
428c0f10754SAlp Dener       ierr = TaoBNCGSetRecycleFlag(bnk->bncg, PETSC_TRUE);CHKERRQ(ierr);
429c0f10754SAlp Dener     }
430c0f10754SAlp Dener     /* Take some small finite number of BNCG iterations */
431c0f10754SAlp Dener     ierr = TaoSolve(bnk->bncg);CHKERRQ(ierr);
432c0f10754SAlp Dener     /* Add the number of gradient and function evaluations to the total */
433c0f10754SAlp Dener     tao->nfuncs += bnk->bncg->nfuncs;
434c0f10754SAlp Dener     tao->nfuncgrads += bnk->bncg->nfuncgrads;
435c0f10754SAlp Dener     tao->ngrads += bnk->bncg->ngrads;
436c0f10754SAlp Dener     tao->nhess += bnk->bncg->nhess;
437e031d6f5SAlp Dener     bnk->tot_cg_its += bnk->bncg->niter;
438c0f10754SAlp Dener     /* Extract the BNCG solution out and save it into BNK */
439c0f10754SAlp Dener     bnk->f = bnk->bncg_ctx->f;
440c0f10754SAlp Dener     ierr = VecCopy(bnk->bncg->solution, tao->solution);
441c0f10754SAlp Dener     ierr = VecCopy(bnk->bncg_ctx->unprojected_gradient, bnk->unprojected_gradient);
442c0f10754SAlp Dener     ierr = VecCopy(bnk->bncg->gradient, tao->gradient);
443c0f10754SAlp Dener     /* Check to see if BNCG converged the problem */
444c0f10754SAlp 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) {
445c0f10754SAlp Dener       *terminate = PETSC_TRUE;
446c0f10754SAlp Dener     }
447c0f10754SAlp Dener   }
448c0f10754SAlp Dener   PetscFunctionReturn(0);
449c0f10754SAlp Dener }
450c0f10754SAlp Dener 
4512f75a4aaSAlp Dener /*------------------------------------------------------------*/
4522f75a4aaSAlp Dener 
453c0f10754SAlp Dener /* Routine for computing the Newton step. */
454df278d8fSAlp Dener 
45562675beeSAlp Dener PetscErrorCode TaoBNKComputeStep(Tao tao, PetscBool shift, KSPConvergedReason *ksp_reason)
456eb910715SAlp Dener {
457eb910715SAlp Dener   PetscErrorCode               ierr;
458eb910715SAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
459eb910715SAlp Dener 
460e465cd6fSAlp Dener   PetscReal                    delta;
461eb910715SAlp Dener   PetscInt                     bfgsUpdates = 0;
462eb910715SAlp Dener   PetscInt                     kspits;
463eb910715SAlp Dener 
464eb910715SAlp Dener   PetscFunctionBegin;
4655e9b73cbSAlp Dener   /* Prepare the reduced sub-matrices for the inactive set */
466*08752603SAlp Dener   ierr = TaoBNKEstimateActiveSet(tao, bnk->as_type);CHKERRQ(ierr);
4672f75a4aaSAlp Dener   if (BNK_PC_BFGS == bnk->pc_type) { ierr = MatLMVMSetInactive(bnk->M, bnk->inactive_idx);CHKERRQ(ierr); }
4682f75a4aaSAlp Dener   if (bnk->inactive_idx) {
4695e9b73cbSAlp Dener     ierr = MatDestroy(&bnk->H_inactive);
4705e9b73cbSAlp Dener     ierr = MatCreateSubMatrix(tao->hessian, bnk->inactive_idx, bnk->inactive_idx, MAT_INITIAL_MATRIX, &bnk->H_inactive);CHKERRQ(ierr);
471eb3ba6a7SAlp Dener     if (tao->hessian == tao->hessian_pre) {
472eb3ba6a7SAlp Dener       bnk->Hpre_inactive = bnk->H_inactive;
473eb3ba6a7SAlp Dener     } else {
4745e9b73cbSAlp Dener       ierr = MatDestroy(&bnk->Hpre_inactive);
4755e9b73cbSAlp Dener       ierr = MatCreateSubMatrix(tao->hessian_pre, bnk->inactive_idx, bnk->inactive_idx, MAT_INITIAL_MATRIX, &bnk->Hpre_inactive);CHKERRQ(ierr);
476eb3ba6a7SAlp Dener     }
4772f75a4aaSAlp Dener   } else {
47862675beeSAlp Dener     ierr = MatDestroy(&bnk->H_inactive);
47962675beeSAlp Dener     ierr = MatDuplicate(tao->hessian, MAT_COPY_VALUES, &bnk->H_inactive);
48062675beeSAlp Dener     if (tao->hessian == tao->hessian_pre) {
48162675beeSAlp Dener       bnk->Hpre_inactive = bnk->H_inactive;
48262675beeSAlp Dener     } else {
48362675beeSAlp Dener       ierr = MatDestroy(&bnk->Hpre_inactive);
48462675beeSAlp Dener       ierr = MatDuplicate(tao->hessian_pre, MAT_COPY_VALUES, &bnk->Hpre_inactive);
48562675beeSAlp Dener     }
48662675beeSAlp Dener   }
48762675beeSAlp Dener 
48862675beeSAlp Dener   /* Shift the reduced Hessian matrix */
48962675beeSAlp Dener   if ((shift) && (bnk->pert > 0)) {
49062675beeSAlp Dener     ierr = MatShift(bnk->H_inactive, bnk->pert);CHKERRQ(ierr);
49162675beeSAlp Dener     if (bnk->H_inactive != bnk->Hpre_inactive) {
49262675beeSAlp Dener       ierr = MatShift(bnk->Hpre_inactive, bnk->pert);CHKERRQ(ierr);
49362675beeSAlp Dener     }
49462675beeSAlp Dener   }
49562675beeSAlp Dener 
49662675beeSAlp Dener   /* Update the diagonal scaling for the BFGS preconditioner, this time with the Hessian perturbation */
49762675beeSAlp Dener   if ((BNK_PC_BFGS == bnk->pc_type) && (BFGS_SCALE_PHESS == bnk->bfgs_scale_type)) {
49862675beeSAlp Dener     /* Obtain diagonal for the bfgs preconditioner  */
4995e9b73cbSAlp Dener     ierr = VecSet(bnk->Diag, 1.0);CHKERRQ(ierr);
5005e9b73cbSAlp Dener     if (bnk->inactive_idx) {
5015e9b73cbSAlp Dener       ierr = VecGetSubVector(bnk->Diag, bnk->inactive_idx, &bnk->Diag_red);CHKERRQ(ierr);
5025e9b73cbSAlp Dener     } else {
5035e9b73cbSAlp Dener       bnk->Diag_red = bnk->Diag;
5045e9b73cbSAlp Dener     }
5055e9b73cbSAlp Dener     ierr = MatGetDiagonal(bnk->H_inactive, bnk->Diag_red);CHKERRQ(ierr);
5065e9b73cbSAlp Dener     if (bnk->inactive_idx) {
5075e9b73cbSAlp Dener       ierr = VecRestoreSubVector(bnk->Diag, bnk->inactive_idx, &bnk->Diag_red);CHKERRQ(ierr);
5085e9b73cbSAlp Dener     }
50962675beeSAlp Dener     ierr = VecAbs(bnk->Diag);CHKERRQ(ierr);
51062675beeSAlp Dener     ierr = VecMedian(bnk->Diag_min, bnk->Diag, bnk->Diag_max, bnk->Diag);CHKERRQ(ierr);
51162675beeSAlp Dener     ierr = VecReciprocal(bnk->Diag);CHKERRQ(ierr);
51262675beeSAlp Dener     ierr = MatLMVMSetScale(bnk->M,bnk->Diag);CHKERRQ(ierr);
5132f75a4aaSAlp Dener   }
51409164190SAlp Dener 
515eb910715SAlp Dener   /* Solve the Newton system of equations */
5162f75a4aaSAlp Dener   ierr = VecSet(tao->stepdirection, 0.0);CHKERRQ(ierr);
5175e9b73cbSAlp Dener   ierr = KSPReset(tao->ksp);CHKERRQ(ierr);
51809164190SAlp Dener   ierr = KSPSetOperators(tao->ksp,bnk->H_inactive,bnk->Hpre_inactive);CHKERRQ(ierr);
5195e9b73cbSAlp Dener   ierr = VecCopy(bnk->unprojected_gradient, bnk->Gwork);CHKERRQ(ierr);
5205e9b73cbSAlp Dener   if (bnk->inactive_idx) {
5215e9b73cbSAlp Dener     ierr = VecGetSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive);CHKERRQ(ierr);
5225e9b73cbSAlp Dener     ierr = VecGetSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr);
5235e9b73cbSAlp Dener   } else {
5245e9b73cbSAlp Dener     bnk->G_inactive = bnk->unprojected_gradient;
5255e9b73cbSAlp Dener     bnk->X_inactive = tao->stepdirection;
52628017e9fSAlp Dener   }
527eb910715SAlp Dener   if (bnk->is_nash || bnk->is_stcg || bnk->is_gltr) {
528fed79b8eSAlp Dener     ierr = KSPCGSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr);
5295e9b73cbSAlp Dener     ierr = KSPSolve(tao->ksp, bnk->G_inactive, bnk->X_inactive);CHKERRQ(ierr);
530eb910715SAlp Dener     ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr);
531eb910715SAlp Dener     tao->ksp_its+=kspits;
532eb910715SAlp Dener     tao->ksp_tot_its+=kspits;
533080d2917SAlp Dener     ierr = KSPCGGetNormD(tao->ksp,&bnk->dnorm);CHKERRQ(ierr);
534eb910715SAlp Dener 
535eb910715SAlp Dener     if (0.0 == tao->trust) {
536eb910715SAlp Dener       /* Radius was uninitialized; use the norm of the direction */
537080d2917SAlp Dener       if (bnk->dnorm > 0.0) {
538080d2917SAlp Dener         tao->trust = bnk->dnorm;
539eb910715SAlp Dener 
540eb910715SAlp Dener         /* Modify the radius if it is too large or small */
541eb910715SAlp Dener         tao->trust = PetscMax(tao->trust, bnk->min_radius);
542eb910715SAlp Dener         tao->trust = PetscMin(tao->trust, bnk->max_radius);
543eb910715SAlp Dener       } else {
544eb910715SAlp Dener         /* The direction was bad; set radius to default value and re-solve
545eb910715SAlp Dener            the trust-region subproblem to get a direction */
546eb910715SAlp Dener         tao->trust = tao->trust0;
547eb910715SAlp Dener 
548eb910715SAlp Dener         /* Modify the radius if it is too large or small */
549eb910715SAlp Dener         tao->trust = PetscMax(tao->trust, bnk->min_radius);
550eb910715SAlp Dener         tao->trust = PetscMin(tao->trust, bnk->max_radius);
551eb910715SAlp Dener 
552fed79b8eSAlp Dener         ierr = KSPCGSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr);
5535e9b73cbSAlp Dener         ierr = KSPSolve(tao->ksp, bnk->G_inactive, bnk->X_inactive);CHKERRQ(ierr);
554eb910715SAlp Dener         ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr);
555eb910715SAlp Dener         tao->ksp_its+=kspits;
556eb910715SAlp Dener         tao->ksp_tot_its+=kspits;
557080d2917SAlp Dener         ierr = KSPCGGetNormD(tao->ksp,&bnk->dnorm);CHKERRQ(ierr);
558eb910715SAlp Dener 
559080d2917SAlp Dener         if (bnk->dnorm == 0.0) SETERRQ(PETSC_COMM_SELF,1, "Initial direction zero");
560eb910715SAlp Dener       }
561eb910715SAlp Dener     }
562eb910715SAlp Dener   } else {
5635e9b73cbSAlp Dener     ierr = KSPSolve(tao->ksp, bnk->G_inactive, bnk->X_inactive);CHKERRQ(ierr);
564eb910715SAlp Dener     ierr = KSPGetIterationNumber(tao->ksp, &kspits);CHKERRQ(ierr);
565eb910715SAlp Dener     tao->ksp_its += kspits;
566eb910715SAlp Dener     tao->ksp_tot_its+=kspits;
567eb910715SAlp Dener   }
5685e9b73cbSAlp Dener   /* Restore sub vectors back */
5695e9b73cbSAlp Dener   if (bnk->inactive_idx) {
5705e9b73cbSAlp Dener     ierr = VecRestoreSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive);CHKERRQ(ierr);
5715e9b73cbSAlp Dener     ierr = VecRestoreSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr);
5725e9b73cbSAlp Dener   }
573770b7498SAlp Dener   /* Make sure the safeguarded fall-back step is zero for actively bounded variables */
574fed79b8eSAlp Dener   ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
5752f75a4aaSAlp Dener   ierr = TaoBNKBoundStep(tao, tao->stepdirection);CHKERRQ(ierr);
576770b7498SAlp Dener 
577770b7498SAlp Dener   /* Record convergence reasons */
578e465cd6fSAlp Dener   ierr = KSPGetConvergedReason(tao->ksp, ksp_reason);CHKERRQ(ierr);
579e465cd6fSAlp Dener   if (KSP_CONVERGED_ATOL == *ksp_reason) {
580770b7498SAlp Dener     ++bnk->ksp_atol;
581e465cd6fSAlp Dener   } else if (KSP_CONVERGED_RTOL == *ksp_reason) {
582770b7498SAlp Dener     ++bnk->ksp_rtol;
583e465cd6fSAlp Dener   } else if (KSP_CONVERGED_CG_CONSTRAINED == *ksp_reason) {
584770b7498SAlp Dener     ++bnk->ksp_ctol;
585e465cd6fSAlp Dener   } else if (KSP_CONVERGED_CG_NEG_CURVE == *ksp_reason) {
586770b7498SAlp Dener     ++bnk->ksp_negc;
587e465cd6fSAlp Dener   } else if (KSP_DIVERGED_DTOL == *ksp_reason) {
588770b7498SAlp Dener     ++bnk->ksp_dtol;
589e465cd6fSAlp Dener   } else if (KSP_DIVERGED_ITS == *ksp_reason) {
590770b7498SAlp Dener     ++bnk->ksp_iter;
591770b7498SAlp Dener   } else {
592770b7498SAlp Dener     ++bnk->ksp_othr;
593770b7498SAlp Dener   }
594fed79b8eSAlp Dener 
595fed79b8eSAlp Dener   /* Make sure the BFGS preconditioner is healthy */
596fed79b8eSAlp Dener   if (bnk->pc_type == BNK_PC_BFGS) {
597fed79b8eSAlp Dener     ierr = MatLMVMGetUpdates(bnk->M, &bfgsUpdates);CHKERRQ(ierr);
598e465cd6fSAlp Dener     if ((KSP_DIVERGED_INDEFINITE_PC == *ksp_reason) && (bfgsUpdates > 1)) {
599fed79b8eSAlp Dener       /* Preconditioner is numerically indefinite; reset the approximation. */
6009b6ef848SAlp Dener       delta = 2.0 * PetscMax(1.0, PetscAbsScalar(bnk->f)) / (bnk->gnorm*bnk->gnorm);
601eb910715SAlp Dener       ierr = MatLMVMSetDelta(bnk->M,delta);CHKERRQ(ierr);
602eb910715SAlp Dener       ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr);
60309164190SAlp Dener       ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
604eb910715SAlp Dener     }
605fed79b8eSAlp Dener   }
606e465cd6fSAlp Dener   PetscFunctionReturn(0);
607e465cd6fSAlp Dener }
608eb910715SAlp Dener 
60962675beeSAlp Dener /*------------------------------------------------------------*/
61062675beeSAlp Dener 
6115e9b73cbSAlp Dener /* Routine for recomputing the predicted reduction for a given step vector */
6125e9b73cbSAlp Dener 
6135e9b73cbSAlp Dener PetscErrorCode TaoBNKRecomputePred(Tao tao, Vec S, PetscReal *prered)
6145e9b73cbSAlp Dener {
6155e9b73cbSAlp Dener   PetscErrorCode               ierr;
6165e9b73cbSAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
6175e9b73cbSAlp Dener 
6185e9b73cbSAlp Dener   PetscFunctionBegin;
6195e9b73cbSAlp Dener   /* Extract subvectors associated with the inactive set */
6205e9b73cbSAlp Dener   if (bnk->inactive_idx){
6215e9b73cbSAlp Dener     ierr = VecGetSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr);
6225e9b73cbSAlp Dener     ierr = VecGetSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work);CHKERRQ(ierr);
6235e9b73cbSAlp Dener     ierr = VecGetSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive);CHKERRQ(ierr);
6245e9b73cbSAlp Dener   } else {
6255e9b73cbSAlp Dener     bnk->X_inactive = tao->stepdirection;
6265e9b73cbSAlp Dener     bnk->inactive_work = bnk->Xwork;
6275e9b73cbSAlp Dener     bnk->G_inactive = bnk->Gwork;
6285e9b73cbSAlp Dener   }
6295e9b73cbSAlp Dener   /* Recompute the predicted decrease based on the quadratic model */
6305e9b73cbSAlp Dener   ierr = MatMult(bnk->H_inactive, bnk->X_inactive, bnk->inactive_work);CHKERRQ(ierr);
6315e9b73cbSAlp Dener   ierr = VecAYPX(bnk->inactive_work, -0.5, bnk->G_inactive);CHKERRQ(ierr);
6325e9b73cbSAlp Dener   ierr = VecDot(bnk->inactive_work, bnk->X_inactive, prered);
6335e9b73cbSAlp Dener   /* Restore the sub vectors */
6345e9b73cbSAlp Dener   if (bnk->inactive_idx){
6355e9b73cbSAlp Dener     ierr = VecRestoreSubVector(tao->stepdirection, bnk->inactive_idx, &bnk->X_inactive);CHKERRQ(ierr);
6365e9b73cbSAlp Dener     ierr = VecRestoreSubVector(bnk->Xwork, bnk->inactive_idx, &bnk->inactive_work);CHKERRQ(ierr);
6375e9b73cbSAlp Dener     ierr = VecRestoreSubVector(bnk->Gwork, bnk->inactive_idx, &bnk->G_inactive);CHKERRQ(ierr);
6385e9b73cbSAlp Dener   }
6395e9b73cbSAlp Dener   PetscFunctionReturn(0);
6405e9b73cbSAlp Dener }
6415e9b73cbSAlp Dener 
6425e9b73cbSAlp Dener /*------------------------------------------------------------*/
6435e9b73cbSAlp Dener 
64462675beeSAlp Dener /* Routine for ensuring that the Newton step is a descent direction.
64562675beeSAlp Dener 
64662675beeSAlp Dener    The step direction falls back onto BFGS, scaled gradient and gradient steps
64762675beeSAlp Dener    in the event that the Newton step fails the test.
64862675beeSAlp Dener */
64962675beeSAlp Dener 
650e465cd6fSAlp Dener PetscErrorCode TaoBNKSafeguardStep(Tao tao, KSPConvergedReason ksp_reason, PetscInt *stepType)
651e465cd6fSAlp Dener {
652e465cd6fSAlp Dener   PetscErrorCode               ierr;
653e465cd6fSAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
654e465cd6fSAlp Dener 
655e465cd6fSAlp Dener   PetscReal                    gdx, delta, e_min;
656e465cd6fSAlp Dener   PetscInt                     bfgsUpdates;
657e465cd6fSAlp Dener 
658e465cd6fSAlp Dener   PetscFunctionBegin;
659080d2917SAlp Dener   ierr = VecDot(tao->stepdirection, tao->gradient, &gdx);CHKERRQ(ierr);
660eb910715SAlp Dener   if ((gdx >= 0.0) || PetscIsInfOrNanReal(gdx)) {
661eb910715SAlp Dener     /* Newton step is not descent or direction produced Inf or NaN
662eb910715SAlp Dener        Update the perturbation for next time */
663eb910715SAlp Dener     if (bnk->pert <= 0.0) {
664eb910715SAlp Dener       /* Initialize the perturbation */
665eb910715SAlp Dener       bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
666eb910715SAlp Dener       if (bnk->is_gltr) {
667eb910715SAlp Dener         ierr = KSPCGGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr);
668eb910715SAlp Dener         bnk->pert = PetscMax(bnk->pert, -e_min);
669eb910715SAlp Dener       }
670eb910715SAlp Dener     } else {
671eb910715SAlp Dener       /* Increase the perturbation */
672eb910715SAlp Dener       bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
673eb910715SAlp Dener     }
674eb910715SAlp Dener 
675eb910715SAlp Dener     if (BNK_PC_BFGS != bnk->pc_type) {
676eb910715SAlp Dener       /* We don't have the bfgs matrix around and updated
677eb910715SAlp Dener          Must use gradient direction in this case */
678080d2917SAlp Dener       ierr = VecCopy(tao->gradient, tao->stepdirection);CHKERRQ(ierr);
679eb910715SAlp Dener       *stepType = BNK_GRADIENT;
680eb910715SAlp Dener     } else {
681eb910715SAlp Dener       /* Attempt to use the BFGS direction */
6822ab2a32cSAlp Dener       ierr = MatLMVMSetInactive(bnk->M, NULL);CHKERRQ(ierr);
68309164190SAlp Dener       ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
684eb910715SAlp Dener 
6858d5ead36SAlp Dener       /* Check for success (descent direction)
6868d5ead36SAlp Dener          NOTE: Negative gdx here means not a descent direction because
6878d5ead36SAlp Dener          the fall-back step is missing a negative sign. */
688080d2917SAlp Dener       ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr);
6898d5ead36SAlp Dener       if ((gdx <= 0) || PetscIsInfOrNanReal(gdx)) {
690eb910715SAlp Dener         /* BFGS direction is not descent or direction produced not a number
691eb910715SAlp Dener            We can assert bfgsUpdates > 1 in this case because
692eb910715SAlp Dener            the first solve produces the scaled gradient direction,
693eb910715SAlp Dener            which is guaranteed to be descent */
694eb910715SAlp Dener 
695eb910715SAlp Dener         /* Use steepest descent direction (scaled) */
6969b6ef848SAlp Dener         delta = 2.0 * PetscMax(1.0, PetscAbsScalar(bnk->f)) / (bnk->gnorm*bnk->gnorm);
697eb910715SAlp Dener         ierr = MatLMVMSetDelta(bnk->M, delta);CHKERRQ(ierr);
698eb910715SAlp Dener         ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr);
69909164190SAlp Dener         ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
70009164190SAlp Dener         ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
701eb910715SAlp Dener 
702eb910715SAlp Dener         *stepType = BNK_SCALED_GRADIENT;
703eb910715SAlp Dener       } else {
704770b7498SAlp Dener         ierr = MatLMVMGetUpdates(bnk->M, &bfgsUpdates);CHKERRQ(ierr);
705eb910715SAlp Dener         if (1 == bfgsUpdates) {
706eb910715SAlp Dener           /* The first BFGS direction is always the scaled gradient */
707eb910715SAlp Dener           *stepType = BNK_SCALED_GRADIENT;
708eb910715SAlp Dener         } else {
709eb910715SAlp Dener           *stepType = BNK_BFGS;
710eb910715SAlp Dener         }
711eb910715SAlp Dener       }
712eb910715SAlp Dener     }
7138d5ead36SAlp Dener     /* Make sure the safeguarded fall-back step is zero for actively bounded variables */
7148d5ead36SAlp Dener     ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
7152f75a4aaSAlp Dener     ierr = TaoBNKBoundStep(tao, tao->stepdirection);CHKERRQ(ierr);
716eb910715SAlp Dener   } else {
717eb910715SAlp Dener     /* Computed Newton step is descent */
718eb910715SAlp Dener     switch (ksp_reason) {
719eb910715SAlp Dener     case KSP_DIVERGED_NANORINF:
720eb910715SAlp Dener     case KSP_DIVERGED_BREAKDOWN:
721eb910715SAlp Dener     case KSP_DIVERGED_INDEFINITE_MAT:
722eb910715SAlp Dener     case KSP_DIVERGED_INDEFINITE_PC:
723eb910715SAlp Dener     case KSP_CONVERGED_CG_NEG_CURVE:
724eb910715SAlp Dener       /* Matrix or preconditioner is indefinite; increase perturbation */
725eb910715SAlp Dener       if (bnk->pert <= 0.0) {
726eb910715SAlp Dener         /* Initialize the perturbation */
727eb910715SAlp Dener         bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
728eb910715SAlp Dener         if (bnk->is_gltr) {
729eb910715SAlp Dener           ierr = KSPCGGLTRGetMinEig(tao->ksp, &e_min);CHKERRQ(ierr);
730eb910715SAlp Dener           bnk->pert = PetscMax(bnk->pert, -e_min);
731eb910715SAlp Dener         }
732eb910715SAlp Dener       } else {
733eb910715SAlp Dener         /* Increase the perturbation */
734eb910715SAlp Dener         bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
735eb910715SAlp Dener       }
736eb910715SAlp Dener       break;
737eb910715SAlp Dener 
738eb910715SAlp Dener     default:
739eb910715SAlp Dener       /* Newton step computation is good; decrease perturbation */
740eb910715SAlp Dener       bnk->pert = PetscMin(bnk->psfac * bnk->pert, bnk->pmsfac * bnk->gnorm);
741eb910715SAlp Dener       if (bnk->pert < bnk->pmin) {
742eb910715SAlp Dener         bnk->pert = 0.0;
743eb910715SAlp Dener       }
744eb910715SAlp Dener       break;
745eb910715SAlp Dener     }
746fed79b8eSAlp Dener     *stepType = BNK_NEWTON;
747eb910715SAlp Dener   }
748eb910715SAlp Dener   PetscFunctionReturn(0);
749eb910715SAlp Dener }
750eb910715SAlp Dener 
751df278d8fSAlp Dener /*------------------------------------------------------------*/
752df278d8fSAlp Dener 
753df278d8fSAlp Dener /* Routine for performing a bound-projected More-Thuente line search.
754df278d8fSAlp Dener 
755df278d8fSAlp Dener   Includes fallbacks to BFGS, scaled gradient, and unscaled gradient steps if the
756df278d8fSAlp Dener   Newton step does not produce a valid step length.
757df278d8fSAlp Dener */
758df278d8fSAlp Dener 
759c14b763aSAlp Dener PetscErrorCode TaoBNKPerformLineSearch(Tao tao, PetscInt stepType, PetscReal *steplen, TaoLineSearchConvergedReason *reason)
760c14b763aSAlp Dener {
761c14b763aSAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
762c14b763aSAlp Dener   PetscErrorCode ierr;
763c14b763aSAlp Dener   TaoLineSearchConvergedReason ls_reason;
764c14b763aSAlp Dener 
765c14b763aSAlp Dener   PetscReal      e_min, gdx, delta;
766c14b763aSAlp Dener   PetscInt       bfgsUpdates;
767c14b763aSAlp Dener 
768c14b763aSAlp Dener   PetscFunctionBegin;
769c14b763aSAlp Dener   /* Perform the linesearch */
770c14b763aSAlp Dener   ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, steplen, &ls_reason);CHKERRQ(ierr);
771c14b763aSAlp Dener   ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
772c14b763aSAlp Dener 
7739b6ef848SAlp Dener   while (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER && (stepType != BNK_GRADIENT || stepType !=BNK_SCALED_GRADIENT)) {
774c14b763aSAlp Dener     /* Linesearch failed, revert solution */
775c14b763aSAlp Dener     bnk->f = bnk->fold;
776c14b763aSAlp Dener     ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr);
777c14b763aSAlp Dener     ierr = VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient);CHKERRQ(ierr);
778c14b763aSAlp Dener 
779c14b763aSAlp Dener     switch(stepType) {
780c14b763aSAlp Dener     case BNK_NEWTON:
7818d5ead36SAlp Dener       /* Failed to obtain acceptable iterate with Newton step
782c14b763aSAlp Dener          Update the perturbation for next time */
783c14b763aSAlp Dener       if (bnk->pert <= 0.0) {
784c14b763aSAlp Dener         /* Initialize the perturbation */
785c14b763aSAlp Dener         bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
786c14b763aSAlp Dener         if (bnk->is_gltr) {
787c14b763aSAlp Dener           ierr = KSPCGGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr);
788c14b763aSAlp Dener           bnk->pert = PetscMax(bnk->pert, -e_min);
789c14b763aSAlp Dener         }
790c14b763aSAlp Dener       } else {
791c14b763aSAlp Dener         /* Increase the perturbation */
792c14b763aSAlp Dener         bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
793c14b763aSAlp Dener       }
794c14b763aSAlp Dener 
795c14b763aSAlp Dener       if (BNK_PC_BFGS != bnk->pc_type) {
796c14b763aSAlp Dener         /* We don't have the bfgs matrix around and being updated
797c14b763aSAlp Dener            Must use gradient direction in this case */
798c14b763aSAlp Dener         ierr = VecCopy(tao->gradient, tao->stepdirection);CHKERRQ(ierr);
799c14b763aSAlp Dener         stepType = BNK_GRADIENT;
800c14b763aSAlp Dener       } else {
801c14b763aSAlp Dener         /* Attempt to use the BFGS direction */
802c14b763aSAlp Dener         ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
8038d5ead36SAlp Dener         /* Check for success (descent direction)
8048d5ead36SAlp Dener            NOTE: Negative gdx means not a descent direction because the step here is missing a negative sign. */
805c14b763aSAlp Dener         ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr);
806c14b763aSAlp Dener         if ((gdx <= 0) || PetscIsInfOrNanReal(gdx)) {
807c14b763aSAlp Dener           /* BFGS direction is not descent or direction produced not a number
808c14b763aSAlp Dener              We can assert bfgsUpdates > 1 in this case
809c14b763aSAlp Dener              Use steepest descent direction (scaled) */
8109b6ef848SAlp Dener           delta = 2.0 * PetscMax(1.0, PetscAbsScalar(bnk->f)) / (bnk->gnorm*bnk->gnorm);
811c14b763aSAlp Dener           ierr = MatLMVMSetDelta(bnk->M, delta);CHKERRQ(ierr);
812c14b763aSAlp Dener           ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr);
813c14b763aSAlp Dener           ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
814c14b763aSAlp Dener           ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
815c14b763aSAlp Dener 
816c14b763aSAlp Dener           bfgsUpdates = 1;
817c14b763aSAlp Dener           stepType = BNK_SCALED_GRADIENT;
818c14b763aSAlp Dener         } else {
819c14b763aSAlp Dener           ierr = MatLMVMGetUpdates(bnk->M, &bfgsUpdates);CHKERRQ(ierr);
820c14b763aSAlp Dener           if (1 == bfgsUpdates) {
821c14b763aSAlp Dener             /* The first BFGS direction is always the scaled gradient */
822c14b763aSAlp Dener             stepType = BNK_SCALED_GRADIENT;
823c14b763aSAlp Dener           } else {
824c14b763aSAlp Dener             stepType = BNK_BFGS;
825c14b763aSAlp Dener           }
826c14b763aSAlp Dener         }
827c14b763aSAlp Dener       }
828c14b763aSAlp Dener       break;
829c14b763aSAlp Dener 
830c14b763aSAlp Dener     case BNK_BFGS:
831c14b763aSAlp Dener       /* Can only enter if pc_type == BNK_PC_BFGS
832c14b763aSAlp Dener          Failed to obtain acceptable iterate with BFGS step
833c14b763aSAlp Dener          Attempt to use the scaled gradient direction */
8349b6ef848SAlp Dener       delta = 2.0 * PetscMax(1.0, PetscAbsScalar(bnk->f)) / (bnk->gnorm*bnk->gnorm);
835c14b763aSAlp Dener       ierr = MatLMVMSetDelta(bnk->M, delta);CHKERRQ(ierr);
836c14b763aSAlp Dener       ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr);
837c14b763aSAlp Dener       ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
838c14b763aSAlp Dener       ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
839c14b763aSAlp Dener 
840c14b763aSAlp Dener       bfgsUpdates = 1;
841c14b763aSAlp Dener       stepType = BNK_SCALED_GRADIENT;
842c14b763aSAlp Dener       break;
843c14b763aSAlp Dener 
844c14b763aSAlp Dener     case BNK_SCALED_GRADIENT:
845c14b763aSAlp Dener       /* Can only enter if pc_type == BNK_PC_BFGS
846c14b763aSAlp Dener          The scaled gradient step did not produce a new iterate;
847c14b763aSAlp Dener          attemp to use the gradient direction.
848c14b763aSAlp Dener          Need to make sure we are not using a different diagonal scaling */
849c14b763aSAlp Dener       ierr = MatLMVMSetScale(bnk->M,0);CHKERRQ(ierr);
850c14b763aSAlp Dener       ierr = MatLMVMSetDelta(bnk->M,1.0);CHKERRQ(ierr);
851c14b763aSAlp Dener       ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr);
852c14b763aSAlp Dener       ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
853c14b763aSAlp Dener       ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
854c14b763aSAlp Dener 
855c14b763aSAlp Dener       bfgsUpdates = 1;
856c14b763aSAlp Dener       stepType = BNK_GRADIENT;
857c14b763aSAlp Dener       break;
858c14b763aSAlp Dener     }
8598d5ead36SAlp Dener     /* Make sure the safeguarded fall-back step is zero for actively bounded variables */
8608d5ead36SAlp Dener     ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
8612f75a4aaSAlp Dener     ierr = TaoBNKBoundStep(tao, tao->stepdirection);CHKERRQ(ierr);
862c14b763aSAlp Dener 
8638d5ead36SAlp Dener     /* Perform one last line search with the fall-back step */
864c14b763aSAlp Dener     ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, steplen, &ls_reason);CHKERRQ(ierr);
865c14b763aSAlp Dener     ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
866c14b763aSAlp Dener   }
867c14b763aSAlp Dener   *reason = ls_reason;
868c14b763aSAlp Dener   PetscFunctionReturn(0);
869c14b763aSAlp Dener }
870c14b763aSAlp Dener 
871df278d8fSAlp Dener /*------------------------------------------------------------*/
872df278d8fSAlp Dener 
873df278d8fSAlp Dener /* Routine for updating the trust radius.
874df278d8fSAlp Dener 
875df278d8fSAlp Dener   Function features three different update methods:
876df278d8fSAlp Dener   1) Line-search step length based
877df278d8fSAlp Dener   2) Predicted decrease on the CG quadratic model
878df278d8fSAlp Dener   3) Interpolation
879df278d8fSAlp Dener */
880df278d8fSAlp Dener 
88128017e9fSAlp Dener PetscErrorCode TaoBNKUpdateTrustRadius(Tao tao, PetscReal prered, PetscReal actred, PetscInt updateType, PetscInt stepType, PetscBool *accept)
882080d2917SAlp Dener {
883080d2917SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
884080d2917SAlp Dener   PetscErrorCode ierr;
885080d2917SAlp Dener 
886b1c2d0e3SAlp Dener   PetscReal      step, kappa;
887080d2917SAlp Dener   PetscReal      gdx, tau_1, tau_2, tau_min, tau_max;
888080d2917SAlp Dener 
889080d2917SAlp Dener   PetscFunctionBegin;
890080d2917SAlp Dener   /* Update trust region radius */
891080d2917SAlp Dener   *accept = PETSC_FALSE;
89228017e9fSAlp Dener   switch(updateType) {
893080d2917SAlp Dener   case BNK_UPDATE_STEP:
894c14b763aSAlp Dener     *accept = PETSC_TRUE; /* always accept here because line search succeeded */
895080d2917SAlp Dener     if (stepType == BNK_NEWTON) {
896080d2917SAlp Dener       ierr = TaoLineSearchGetStepLength(tao->linesearch, &step);CHKERRQ(ierr);
897080d2917SAlp Dener       if (step < bnk->nu1) {
898080d2917SAlp Dener         /* Very bad step taken; reduce radius */
899080d2917SAlp Dener         tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust);
900080d2917SAlp Dener       } else if (step < bnk->nu2) {
901080d2917SAlp Dener         /* Reasonably bad step taken; reduce radius */
902080d2917SAlp Dener         tao->trust = bnk->omega2 * PetscMin(bnk->dnorm, tao->trust);
903080d2917SAlp Dener       } else if (step < bnk->nu3) {
904080d2917SAlp Dener         /*  Reasonable step was taken; leave radius alone */
905080d2917SAlp Dener         if (bnk->omega3 < 1.0) {
906080d2917SAlp Dener           tao->trust = bnk->omega3 * PetscMin(bnk->dnorm, tao->trust);
907080d2917SAlp Dener         } else if (bnk->omega3 > 1.0) {
908080d2917SAlp Dener           tao->trust = PetscMax(bnk->omega3 * bnk->dnorm, tao->trust);
909080d2917SAlp Dener         }
910080d2917SAlp Dener       } else if (step < bnk->nu4) {
911080d2917SAlp Dener         /*  Full step taken; increase the radius */
912080d2917SAlp Dener         tao->trust = PetscMax(bnk->omega4 * bnk->dnorm, tao->trust);
913080d2917SAlp Dener       } else {
914080d2917SAlp Dener         /*  More than full step taken; increase the radius */
915080d2917SAlp Dener         tao->trust = PetscMax(bnk->omega5 * bnk->dnorm, tao->trust);
916080d2917SAlp Dener       }
917080d2917SAlp Dener     } else {
918080d2917SAlp Dener       /*  Newton step was not good; reduce the radius */
919080d2917SAlp Dener       tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust);
920080d2917SAlp Dener     }
921080d2917SAlp Dener     break;
922080d2917SAlp Dener 
923080d2917SAlp Dener   case BNK_UPDATE_REDUCTION:
924080d2917SAlp Dener     if (stepType == BNK_NEWTON) {
925b1c2d0e3SAlp Dener       if (prered < 0.0) {
926fed79b8eSAlp Dener         /* The predicted reduction has the wrong sign.  This cannot
927fed79b8eSAlp Dener            happen in infinite precision arithmetic.  Step should
928fed79b8eSAlp Dener            be rejected! */
929080d2917SAlp Dener         tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
930fed79b8eSAlp Dener       }
931fed79b8eSAlp Dener       else {
932b1c2d0e3SAlp Dener         if (PetscIsInfOrNanReal(actred)) {
933080d2917SAlp Dener           tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
934080d2917SAlp Dener         } else {
9350a4511e9SAlp Dener           if ((PetscAbsScalar(actred) <= PetscMax(1.0, PetscAbsScalar(bnk->f))*bnk->epsilon) &&
9360a4511e9SAlp Dener               (PetscAbsScalar(prered) <= PetscMax(1.0, PetscAbsScalar(bnk->f))*bnk->epsilon)) {
937080d2917SAlp Dener             kappa = 1.0;
938fed79b8eSAlp Dener           }
939fed79b8eSAlp Dener           else {
940080d2917SAlp Dener             kappa = actred / prered;
941080d2917SAlp Dener           }
942fed79b8eSAlp Dener 
943fed79b8eSAlp Dener           /* Accept or reject the step and update radius */
944080d2917SAlp Dener           if (kappa < bnk->eta1) {
945fed79b8eSAlp Dener             /* Reject the step */
946080d2917SAlp Dener             tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
947fed79b8eSAlp Dener           }
948fed79b8eSAlp Dener           else {
949fed79b8eSAlp Dener             /* Accept the step */
950c133c014SAlp Dener             *accept = PETSC_TRUE;
951c133c014SAlp Dener             /* Update the trust region radius only if the computed step is at the trust radius boundary */
9528d5ead36SAlp Dener             if (bnk->dnorm == tao->trust) {
953080d2917SAlp Dener               if (kappa < bnk->eta2) {
954080d2917SAlp Dener                 /* Marginal bad step */
955c133c014SAlp Dener                 tao->trust = bnk->alpha2 * tao->trust;
956080d2917SAlp Dener               }
957fed79b8eSAlp Dener               else if (kappa < bnk->eta3) {
958fed79b8eSAlp Dener                 /* Reasonable step */
959fed79b8eSAlp Dener                 tao->trust = bnk->alpha3 * tao->trust;
960fed79b8eSAlp Dener               }
961fed79b8eSAlp Dener               else if (kappa < bnk->eta4) {
962080d2917SAlp Dener                 /* Good step */
963c133c014SAlp Dener                 tao->trust = bnk->alpha4 * tao->trust;
964fed79b8eSAlp Dener               }
965fed79b8eSAlp Dener               else {
966080d2917SAlp Dener                 /* Very good step */
967c133c014SAlp Dener                 tao->trust = bnk->alpha5 * tao->trust;
968080d2917SAlp Dener               }
969c133c014SAlp Dener             }
970080d2917SAlp Dener           }
971080d2917SAlp Dener         }
972080d2917SAlp Dener       }
973080d2917SAlp Dener     } else {
974080d2917SAlp Dener       /*  Newton step was not good; reduce the radius */
975080d2917SAlp Dener       tao->trust = bnk->alpha1 * PetscMin(bnk->dnorm, tao->trust);
976080d2917SAlp Dener     }
977080d2917SAlp Dener     break;
978080d2917SAlp Dener 
979080d2917SAlp Dener   default:
980080d2917SAlp Dener     if (stepType == BNK_NEWTON) {
981b1c2d0e3SAlp Dener       if (prered < 0.0) {
982080d2917SAlp Dener         /*  The predicted reduction has the wrong sign.  This cannot */
983080d2917SAlp Dener         /*  happen in infinite precision arithmetic.  Step should */
984080d2917SAlp Dener         /*  be rejected! */
985080d2917SAlp Dener         tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
986080d2917SAlp Dener       } else {
987b1c2d0e3SAlp Dener         if (PetscIsInfOrNanReal(actred)) {
988080d2917SAlp Dener           tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
989080d2917SAlp Dener         } else {
990080d2917SAlp Dener           if ((PetscAbsScalar(actred) <= bnk->epsilon) && (PetscAbsScalar(prered) <= bnk->epsilon)) {
991080d2917SAlp Dener             kappa = 1.0;
992080d2917SAlp Dener           } else {
993080d2917SAlp Dener             kappa = actred / prered;
994080d2917SAlp Dener           }
995080d2917SAlp Dener 
996080d2917SAlp Dener           ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr);
997080d2917SAlp Dener           tau_1 = bnk->theta * gdx / (bnk->theta * gdx - (1.0 - bnk->theta) * prered + actred);
998080d2917SAlp Dener           tau_2 = bnk->theta * gdx / (bnk->theta * gdx + (1.0 + bnk->theta) * prered - actred);
999080d2917SAlp Dener           tau_min = PetscMin(tau_1, tau_2);
1000080d2917SAlp Dener           tau_max = PetscMax(tau_1, tau_2);
1001080d2917SAlp Dener 
1002080d2917SAlp Dener           if (kappa >= 1.0 - bnk->mu1) {
1003080d2917SAlp Dener             /*  Great agreement */
1004080d2917SAlp Dener             *accept = PETSC_TRUE;
1005080d2917SAlp Dener             if (tau_max < 1.0) {
1006080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm);
1007080d2917SAlp Dener             } else if (tau_max > bnk->gamma4) {
1008080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, bnk->gamma4 * bnk->dnorm);
1009080d2917SAlp Dener             } else {
1010080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm);
1011080d2917SAlp Dener             }
1012080d2917SAlp Dener           } else if (kappa >= 1.0 - bnk->mu2) {
1013080d2917SAlp Dener             /*  Good agreement */
1014080d2917SAlp Dener             *accept = PETSC_TRUE;
1015080d2917SAlp Dener             if (tau_max < bnk->gamma2) {
1016080d2917SAlp Dener               tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm);
1017080d2917SAlp Dener             } else if (tau_max > bnk->gamma3) {
1018080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm);
1019080d2917SAlp Dener             } else if (tau_max < 1.0) {
1020080d2917SAlp Dener               tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm);
1021080d2917SAlp Dener             } else {
1022080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm);
1023080d2917SAlp Dener             }
1024080d2917SAlp Dener           } else {
1025080d2917SAlp Dener             /*  Not good agreement */
1026080d2917SAlp Dener             if (tau_min > 1.0) {
1027080d2917SAlp Dener               tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm);
1028080d2917SAlp Dener             } else if (tau_max < bnk->gamma1) {
1029080d2917SAlp Dener               tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
1030080d2917SAlp Dener             } else if ((tau_min < bnk->gamma1) && (tau_max >= 1.0)) {
1031080d2917SAlp Dener               tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
1032080d2917SAlp Dener             } else if ((tau_1 >= bnk->gamma1) && (tau_1 < 1.0) && ((tau_2 < bnk->gamma1) || (tau_2 >= 1.0))) {
1033080d2917SAlp Dener               tao->trust = tau_1 * PetscMin(tao->trust, bnk->dnorm);
1034080d2917SAlp Dener             } else if ((tau_2 >= bnk->gamma1) && (tau_2 < 1.0) && ((tau_1 < bnk->gamma1) || (tau_2 >= 1.0))) {
1035080d2917SAlp Dener               tao->trust = tau_2 * PetscMin(tao->trust, bnk->dnorm);
1036080d2917SAlp Dener             } else {
1037080d2917SAlp Dener               tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm);
1038080d2917SAlp Dener             }
1039080d2917SAlp Dener           }
1040080d2917SAlp Dener         }
1041080d2917SAlp Dener       }
1042080d2917SAlp Dener     } else {
1043080d2917SAlp Dener       /*  Newton step was not good; reduce the radius */
1044080d2917SAlp Dener       tao->trust = bnk->gamma1 * PetscMin(bnk->dnorm, tao->trust);
1045080d2917SAlp Dener     }
104628017e9fSAlp Dener     break;
1047080d2917SAlp Dener   }
1048c133c014SAlp Dener   /* Make sure the radius does not violate min and max settings */
1049c133c014SAlp Dener   tao->trust = PetscMin(tao->trust, bnk->max_radius);
1050fed79b8eSAlp Dener   tao->trust = PetscMax(tao->trust, bnk->min_radius);
1051080d2917SAlp Dener   PetscFunctionReturn(0);
1052080d2917SAlp Dener }
1053080d2917SAlp Dener 
1054eb910715SAlp Dener /* ---------------------------------------------------------- */
1055df278d8fSAlp Dener 
105662675beeSAlp Dener PetscErrorCode TaoBNKAddStepCounts(Tao tao, PetscInt stepType)
105762675beeSAlp Dener {
105862675beeSAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
105962675beeSAlp Dener 
106062675beeSAlp Dener   PetscFunctionBegin;
106162675beeSAlp Dener   switch (stepType) {
106262675beeSAlp Dener   case BNK_NEWTON:
106362675beeSAlp Dener     ++bnk->newt;
106462675beeSAlp Dener     break;
106562675beeSAlp Dener   case BNK_BFGS:
106662675beeSAlp Dener     ++bnk->bfgs;
106762675beeSAlp Dener     break;
106862675beeSAlp Dener   case BNK_SCALED_GRADIENT:
106962675beeSAlp Dener     ++bnk->sgrad;
107062675beeSAlp Dener     break;
107162675beeSAlp Dener   case BNK_GRADIENT:
107262675beeSAlp Dener     ++bnk->grad;
107362675beeSAlp Dener     break;
107462675beeSAlp Dener   default:
107562675beeSAlp Dener     break;
107662675beeSAlp Dener   }
107762675beeSAlp Dener   PetscFunctionReturn(0);
107862675beeSAlp Dener }
107962675beeSAlp Dener 
108062675beeSAlp Dener /* ---------------------------------------------------------- */
108162675beeSAlp Dener 
10829b6ef848SAlp Dener PetscErrorCode TaoSetUp_BNK(Tao tao)
1083eb910715SAlp Dener {
1084eb910715SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
1085eb910715SAlp Dener   PetscErrorCode ierr;
10869b6ef848SAlp Dener   KSPType        ksp_type;
1087e031d6f5SAlp Dener   PetscInt       i;
1088eb910715SAlp Dener 
1089eb910715SAlp Dener   PetscFunctionBegin;
1090eb910715SAlp Dener   if (!tao->gradient) {ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);}
1091eb910715SAlp Dener   if (!tao->stepdirection) {ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr);}
1092eb910715SAlp Dener   if (!bnk->W) {ierr = VecDuplicate(tao->solution,&bnk->W);CHKERRQ(ierr);}
1093eb910715SAlp Dener   if (!bnk->Xold) {ierr = VecDuplicate(tao->solution,&bnk->Xold);CHKERRQ(ierr);}
1094eb910715SAlp Dener   if (!bnk->Gold) {ierr = VecDuplicate(tao->solution,&bnk->Gold);CHKERRQ(ierr);}
109509164190SAlp Dener   if (!bnk->Xwork) {ierr = VecDuplicate(tao->solution,&bnk->Xwork);CHKERRQ(ierr);}
109609164190SAlp Dener   if (!bnk->Gwork) {ierr = VecDuplicate(tao->solution,&bnk->Gwork);CHKERRQ(ierr);}
109709164190SAlp Dener   if (!bnk->unprojected_gradient) {ierr = VecDuplicate(tao->solution,&bnk->unprojected_gradient);CHKERRQ(ierr);}
109809164190SAlp Dener   if (!bnk->unprojected_gradient_old) {ierr = VecDuplicate(tao->solution,&bnk->unprojected_gradient_old);CHKERRQ(ierr);}
10995e9b73cbSAlp Dener   if (!bnk->Diag_min) {ierr = VecDuplicate(tao->solution,&bnk->Diag_min);CHKERRQ(ierr);}
11005e9b73cbSAlp Dener   if (!bnk->Diag_max) {ierr = VecDuplicate(tao->solution,&bnk->Diag_max);CHKERRQ(ierr);}
1101e031d6f5SAlp Dener   if (bnk->max_cg_its > 0) {
1102e031d6f5SAlp Dener     if (!bnk->bncg_sol) {ierr = VecDuplicate(tao->solution,&bnk->bncg_sol);CHKERRQ(ierr);}
1103e031d6f5SAlp Dener     ierr = TaoSetInitialVector(bnk->bncg, bnk->bncg_sol);CHKERRQ(ierr);
1104e031d6f5SAlp Dener     ierr = TaoSetMaximumIterations(bnk->bncg, bnk->max_cg_its);CHKERRQ(ierr);
1105e031d6f5SAlp Dener     ierr = TaoSetTolerances(bnk->bncg, tao->gatol, tao->grtol, tao->gttol);CHKERRQ(ierr);
1106e031d6f5SAlp Dener     ierr = TaoSetFunctionLowerBound(bnk->bncg, tao->fmin);CHKERRQ(ierr);
1107e031d6f5SAlp Dener 
1108e031d6f5SAlp Dener     ierr = TaoSetObjectiveRoutine(bnk->bncg, tao->ops->computeobjective, tao->user_objP);CHKERRQ(ierr);
1109e031d6f5SAlp Dener     ierr = TaoSetGradientRoutine(bnk->bncg, tao->ops->computegradient, tao->user_gradP);CHKERRQ(ierr);
1110e031d6f5SAlp Dener     ierr = TaoSetObjectiveAndGradientRoutine(bnk->bncg, tao->ops->computeobjectiveandgradient, tao->user_objgradP);CHKERRQ(ierr);
1111e031d6f5SAlp Dener     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)tao, (PetscObject)(bnk->bncg));CHKERRQ(ierr);
1112e031d6f5SAlp Dener     for (i=0; i<tao->numbermonitors; i++) {
1113e031d6f5SAlp Dener       ierr = TaoSetMonitor(bnk->bncg, tao->monitor[i], tao->monitorcontext[i], tao->monitordestroy[i]);CHKERRQ(ierr);
1114e031d6f5SAlp Dener       ierr = PetscObjectReference((PetscObject)(tao->monitorcontext[i]));CHKERRQ(ierr);
1115e031d6f5SAlp Dener     }
1116e031d6f5SAlp Dener     bnk->bncg_ctx = (TAO_BNCG *)bnk->bncg->data;
1117e031d6f5SAlp Dener   }
1118eb910715SAlp Dener   bnk->Diag = 0;
1119c0f10754SAlp Dener   bnk->Diag_red = 0;
1120c0f10754SAlp Dener   bnk->X_inactive = 0;
1121c0f10754SAlp Dener   bnk->G_inactive = 0;
112262675beeSAlp Dener   bnk->inactive_work = 0;
112362675beeSAlp Dener   bnk->active_work = 0;
112462675beeSAlp Dener   bnk->inactive_idx = 0;
112562675beeSAlp Dener   bnk->active_idx = 0;
112662675beeSAlp Dener   bnk->active_lower = 0;
112762675beeSAlp Dener   bnk->active_upper = 0;
112862675beeSAlp Dener   bnk->active_fixed = 0;
1129eb910715SAlp Dener   bnk->M = 0;
113009164190SAlp Dener   bnk->H_inactive = 0;
113109164190SAlp Dener   bnk->Hpre_inactive = 0;
11329b6ef848SAlp Dener   ierr = KSPGetType(tao->ksp,&ksp_type);CHKERRQ(ierr);
11339b6ef848SAlp Dener   ierr = PetscStrcmp(ksp_type,KSPCGNASH,&bnk->is_nash);CHKERRQ(ierr);
11349b6ef848SAlp Dener   ierr = PetscStrcmp(ksp_type,KSPCGSTCG,&bnk->is_stcg);CHKERRQ(ierr);
11359b6ef848SAlp Dener   ierr = PetscStrcmp(ksp_type,KSPCGGLTR,&bnk->is_gltr);CHKERRQ(ierr);
1136eb910715SAlp Dener   PetscFunctionReturn(0);
1137eb910715SAlp Dener }
1138eb910715SAlp Dener 
1139eb910715SAlp Dener /*------------------------------------------------------------*/
1140df278d8fSAlp Dener 
1141eb910715SAlp Dener static PetscErrorCode TaoDestroy_BNK(Tao tao)
1142eb910715SAlp Dener {
1143eb910715SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
1144eb910715SAlp Dener   PetscErrorCode ierr;
1145eb910715SAlp Dener 
1146eb910715SAlp Dener   PetscFunctionBegin;
1147eb910715SAlp Dener   if (tao->setupcalled) {
1148eb910715SAlp Dener     ierr = VecDestroy(&bnk->W);CHKERRQ(ierr);
1149eb910715SAlp Dener     ierr = VecDestroy(&bnk->Xold);CHKERRQ(ierr);
1150eb910715SAlp Dener     ierr = VecDestroy(&bnk->Gold);CHKERRQ(ierr);
115109164190SAlp Dener     ierr = VecDestroy(&bnk->Xwork);CHKERRQ(ierr);
115209164190SAlp Dener     ierr = VecDestroy(&bnk->Gwork);CHKERRQ(ierr);
115309164190SAlp Dener     ierr = VecDestroy(&bnk->unprojected_gradient);CHKERRQ(ierr);
115409164190SAlp Dener     ierr = VecDestroy(&bnk->unprojected_gradient_old);CHKERRQ(ierr);
115562675beeSAlp Dener     ierr = VecDestroy(&bnk->Diag_min);CHKERRQ(ierr);
115662675beeSAlp Dener     ierr = VecDestroy(&bnk->Diag_max);CHKERRQ(ierr);
1157c0f10754SAlp Dener     if (bnk->max_cg_its > 0) {
1158c0f10754SAlp Dener       ierr = TaoDestroy(&bnk->bncg);CHKERRQ(ierr);
1159c0f10754SAlp Dener       ierr = VecDestroy(&bnk->bncg_sol);CHKERRQ(ierr);
1160c0f10754SAlp Dener     }
11615e9b73cbSAlp Dener   }
11625e9b73cbSAlp Dener   ierr = VecDestroy(&bnk->Diag);CHKERRQ(ierr);
1163eb910715SAlp Dener   ierr = MatDestroy(&bnk->M);CHKERRQ(ierr);
1164c0f10754SAlp Dener   if (bnk->Hpre_inactive != tao->hessian_pre && bnk->Hpre_inactive != bnk->H_inactive) {ierr = MatDestroy(&bnk->Hpre_inactive);CHKERRQ(ierr);}
1165c0f10754SAlp Dener   if (bnk->H_inactive != tao->hessian) {ierr = MatDestroy(&bnk->H_inactive);CHKERRQ(ierr);}
1166eb910715SAlp Dener   ierr = PetscFree(tao->data);CHKERRQ(ierr);
1167eb910715SAlp Dener   PetscFunctionReturn(0);
1168eb910715SAlp Dener }
1169eb910715SAlp Dener 
1170eb910715SAlp Dener /*------------------------------------------------------------*/
1171df278d8fSAlp Dener 
1172eb910715SAlp Dener static PetscErrorCode TaoSetFromOptions_BNK(PetscOptionItems *PetscOptionsObject,Tao tao)
1173eb910715SAlp Dener {
1174eb910715SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
1175eb910715SAlp Dener   PetscErrorCode ierr;
1176eb910715SAlp Dener 
1177eb910715SAlp Dener   PetscFunctionBegin;
1178eb910715SAlp Dener   ierr = PetscOptionsHead(PetscOptionsObject,"Newton line search method for unconstrained optimization");CHKERRQ(ierr);
11798d5ead36SAlp 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);
11808d5ead36SAlp 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);
11818d5ead36SAlp 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);
11828d5ead36SAlp 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);
11832f75a4aaSAlp 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);
11848d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_sval", "perturbation starting value", "", bnk->sval, &bnk->sval,NULL);CHKERRQ(ierr);
11858d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_imin", "minimum initial perturbation", "", bnk->imin, &bnk->imin,NULL);CHKERRQ(ierr);
11868d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_imax", "maximum initial perturbation", "", bnk->imax, &bnk->imax,NULL);CHKERRQ(ierr);
11878d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_imfac", "initial merit factor", "", bnk->imfac, &bnk->imfac,NULL);CHKERRQ(ierr);
11888d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_pmin", "minimum perturbation", "", bnk->pmin, &bnk->pmin,NULL);CHKERRQ(ierr);
11898d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_pmax", "maximum perturbation", "", bnk->pmax, &bnk->pmax,NULL);CHKERRQ(ierr);
11908d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_pgfac", "growth factor", "", bnk->pgfac, &bnk->pgfac,NULL);CHKERRQ(ierr);
11918d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_psfac", "shrink factor", "", bnk->psfac, &bnk->psfac,NULL);CHKERRQ(ierr);
11928d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_pmgfac", "merit growth factor", "", bnk->pmgfac, &bnk->pmgfac,NULL);CHKERRQ(ierr);
11938d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_pmsfac", "merit shrink factor", "", bnk->pmsfac, &bnk->pmsfac,NULL);CHKERRQ(ierr);
11948d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_eta1", "poor steplength; reduce radius", "", bnk->eta1, &bnk->eta1,NULL);CHKERRQ(ierr);
11958d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_eta2", "reasonable steplength; leave radius alone", "", bnk->eta2, &bnk->eta2,NULL);CHKERRQ(ierr);
11968d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_eta3", "good steplength; increase radius", "", bnk->eta3, &bnk->eta3,NULL);CHKERRQ(ierr);
11978d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_eta4", "excellent steplength; greatly increase radius", "", bnk->eta4, &bnk->eta4,NULL);CHKERRQ(ierr);
11988d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_alpha1", "", "", bnk->alpha1, &bnk->alpha1,NULL);CHKERRQ(ierr);
11998d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_alpha2", "", "", bnk->alpha2, &bnk->alpha2,NULL);CHKERRQ(ierr);
12008d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_alpha3", "", "", bnk->alpha3, &bnk->alpha3,NULL);CHKERRQ(ierr);
12018d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_alpha4", "", "", bnk->alpha4, &bnk->alpha4,NULL);CHKERRQ(ierr);
12028d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_alpha5", "", "", bnk->alpha5, &bnk->alpha5,NULL);CHKERRQ(ierr);
12038d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_nu1", "poor steplength; reduce radius", "", bnk->nu1, &bnk->nu1,NULL);CHKERRQ(ierr);
12048d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_nu2", "reasonable steplength; leave radius alone", "", bnk->nu2, &bnk->nu2,NULL);CHKERRQ(ierr);
12058d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_nu3", "good steplength; increase radius", "", bnk->nu3, &bnk->nu3,NULL);CHKERRQ(ierr);
12068d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_nu4", "excellent steplength; greatly increase radius", "", bnk->nu4, &bnk->nu4,NULL);CHKERRQ(ierr);
12078d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_omega1", "", "", bnk->omega1, &bnk->omega1,NULL);CHKERRQ(ierr);
12088d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_omega2", "", "", bnk->omega2, &bnk->omega2,NULL);CHKERRQ(ierr);
12098d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_omega3", "", "", bnk->omega3, &bnk->omega3,NULL);CHKERRQ(ierr);
12108d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_omega4", "", "", bnk->omega4, &bnk->omega4,NULL);CHKERRQ(ierr);
12118d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_omega5", "", "", bnk->omega5, &bnk->omega5,NULL);CHKERRQ(ierr);
12128d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_mu1_i", "", "", bnk->mu1_i, &bnk->mu1_i,NULL);CHKERRQ(ierr);
12138d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_mu2_i", "", "", bnk->mu2_i, &bnk->mu2_i,NULL);CHKERRQ(ierr);
12148d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_gamma1_i", "", "", bnk->gamma1_i, &bnk->gamma1_i,NULL);CHKERRQ(ierr);
12158d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_gamma2_i", "", "", bnk->gamma2_i, &bnk->gamma2_i,NULL);CHKERRQ(ierr);
12168d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_gamma3_i", "", "", bnk->gamma3_i, &bnk->gamma3_i,NULL);CHKERRQ(ierr);
12178d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_gamma4_i", "", "", bnk->gamma4_i, &bnk->gamma4_i,NULL);CHKERRQ(ierr);
12188d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_theta_i", "", "", bnk->theta_i, &bnk->theta_i,NULL);CHKERRQ(ierr);
12198d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_mu1", "", "", bnk->mu1, &bnk->mu1,NULL);CHKERRQ(ierr);
12208d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_mu2", "", "", bnk->mu2, &bnk->mu2,NULL);CHKERRQ(ierr);
12218d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_gamma1", "", "", bnk->gamma1, &bnk->gamma1,NULL);CHKERRQ(ierr);
12228d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_gamma2", "", "", bnk->gamma2, &bnk->gamma2,NULL);CHKERRQ(ierr);
12238d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_gamma3", "", "", bnk->gamma3, &bnk->gamma3,NULL);CHKERRQ(ierr);
12248d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_gamma4", "", "", bnk->gamma4, &bnk->gamma4,NULL);CHKERRQ(ierr);
12258d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_theta", "", "", bnk->theta, &bnk->theta,NULL);CHKERRQ(ierr);
12268d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_min_radius", "lower bound on initial radius", "", bnk->min_radius, &bnk->min_radius,NULL);CHKERRQ(ierr);
12278d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_max_radius", "upper bound on radius", "", bnk->max_radius, &bnk->max_radius,NULL);CHKERRQ(ierr);
12288d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_epsilon", "tolerance used when computing actual and predicted reduction", "", bnk->epsilon, &bnk->epsilon,NULL);CHKERRQ(ierr);
12290a4511e9SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_as_tol", "initial tolerance used when estimating actively bounded variables", "", bnk->as_tol, &bnk->as_tol,NULL);CHKERRQ(ierr);
12300a4511e9SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_as_step", "step length used when estimating actively bounded variables", "", bnk->as_step, &bnk->as_step,NULL);CHKERRQ(ierr);
1231c0f10754SAlp 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);
1232eb910715SAlp Dener   ierr = PetscOptionsTail();CHKERRQ(ierr);
1233e031d6f5SAlp Dener   ierr = TaoSetFromOptions(bnk->bncg);
1234eb910715SAlp Dener   ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr);
1235eb910715SAlp Dener   ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
1236eb910715SAlp Dener   PetscFunctionReturn(0);
1237eb910715SAlp Dener }
1238eb910715SAlp Dener 
1239eb910715SAlp Dener /*------------------------------------------------------------*/
1240df278d8fSAlp Dener 
1241eb910715SAlp Dener static PetscErrorCode TaoView_BNK(Tao tao, PetscViewer viewer)
1242eb910715SAlp Dener {
1243eb910715SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
1244eb910715SAlp Dener   PetscInt       nrejects;
1245eb910715SAlp Dener   PetscBool      isascii;
1246eb910715SAlp Dener   PetscErrorCode ierr;
1247eb910715SAlp Dener 
1248eb910715SAlp Dener   PetscFunctionBegin;
1249eb910715SAlp Dener   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
1250eb910715SAlp Dener   if (isascii) {
1251eb910715SAlp Dener     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1252eb910715SAlp Dener     if (BNK_PC_BFGS == bnk->pc_type && bnk->M) {
1253eb910715SAlp Dener       ierr = MatLMVMGetRejects(bnk->M,&nrejects);CHKERRQ(ierr);
1254eb910715SAlp Dener       ierr = PetscViewerASCIIPrintf(viewer, "Rejected matrix updates: %D\n",nrejects);CHKERRQ(ierr);
1255eb910715SAlp Dener     }
1256e031d6f5SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "CG steps: %D\n", bnk->tot_cg_its);CHKERRQ(ierr);
1257eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "Newton steps: %D\n", bnk->newt);CHKERRQ(ierr);
1258eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "BFGS steps: %D\n", bnk->bfgs);CHKERRQ(ierr);
1259eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "Scaled gradient steps: %D\n", bnk->sgrad);CHKERRQ(ierr);
1260eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", bnk->grad);CHKERRQ(ierr);
1261eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "KSP termination reasons:\n");CHKERRQ(ierr);
1262eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  atol: %D\n", bnk->ksp_atol);CHKERRQ(ierr);
1263eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  rtol: %D\n", bnk->ksp_rtol);CHKERRQ(ierr);
1264eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  ctol: %D\n", bnk->ksp_ctol);CHKERRQ(ierr);
1265eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  negc: %D\n", bnk->ksp_negc);CHKERRQ(ierr);
1266eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  dtol: %D\n", bnk->ksp_dtol);CHKERRQ(ierr);
1267eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  iter: %D\n", bnk->ksp_iter);CHKERRQ(ierr);
1268eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  othr: %D\n", bnk->ksp_othr);CHKERRQ(ierr);
1269eb910715SAlp Dener     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1270eb910715SAlp Dener   }
1271eb910715SAlp Dener   PetscFunctionReturn(0);
1272eb910715SAlp Dener }
1273eb910715SAlp Dener 
1274eb910715SAlp Dener /* ---------------------------------------------------------- */
1275df278d8fSAlp Dener 
1276eb910715SAlp Dener /*MC
1277eb910715SAlp Dener   TAOBNK - Shared base-type for Bounded Newton-Krylov type algorithms.
127866ed3702SAlp Dener   At each iteration, the BNK methods solve the symmetric
1279eb910715SAlp Dener   system of equations to obtain the step diretion dk:
1280eb910715SAlp Dener               Hk dk = -gk
12812b97c8d8SAlp Dener   for free variables only. The step can be globalized either through
12822b97c8d8SAlp Dener   trust-region methods, or a line search, or a heuristic mixture of both.
1283eb910715SAlp Dener 
1284eb910715SAlp Dener     Options Database Keys:
12858d5ead36SAlp Dener + -tao_bnk_pc_type - "none","ahess","bfgs","petsc"
12868d5ead36SAlp Dener . -tao_bnk_bfgs_scale_type - "ahess","phess","bfgs"
12878d5ead36SAlp Dener . -tao_bnk_init_type - "constant","direction","interpolation"
12888d5ead36SAlp Dener . -tao_bnk_update_type - "step","direction","interpolation"
12892f75a4aaSAlp Dener . -tao_bnk_as_type - "none","bertsekas"
12908d5ead36SAlp Dener . -tao_bnk_sval - perturbation starting value
12918d5ead36SAlp Dener . -tao_bnk_imin - minimum initial perturbation
12928d5ead36SAlp Dener . -tao_bnk_imax - maximum initial perturbation
12938d5ead36SAlp Dener . -tao_bnk_pmin - minimum perturbation
12948d5ead36SAlp Dener . -tao_bnk_pmax - maximum perturbation
12958d5ead36SAlp Dener . -tao_bnk_pgfac - growth factor
12968d5ead36SAlp Dener . -tao_bnk_psfac - shrink factor
12978d5ead36SAlp Dener . -tao_bnk_imfac - initial merit factor
12988d5ead36SAlp Dener . -tao_bnk_pmgfac - merit growth factor
12998d5ead36SAlp Dener . -tao_bnk_pmsfac - merit shrink factor
13008d5ead36SAlp Dener . -tao_bnk_eta1 - poor steplength; reduce radius
13018d5ead36SAlp Dener . -tao_bnk_eta2 - reasonable steplength; leave radius
13028d5ead36SAlp Dener . -tao_bnk_eta3 - good steplength; increase readius
13038d5ead36SAlp Dener . -tao_bnk_eta4 - excellent steplength; greatly increase radius
13048d5ead36SAlp Dener . -tao_bnk_alpha1 - alpha1 reduction
13058d5ead36SAlp Dener . -tao_bnk_alpha2 - alpha2 reduction
13068d5ead36SAlp Dener . -tao_bnk_alpha3 - alpha3 reduction
13078d5ead36SAlp Dener . -tao_bnk_alpha4 - alpha4 reduction
13088d5ead36SAlp Dener . -tao_bnk_alpha - alpha5 reduction
13098d5ead36SAlp Dener . -tao_bnk_mu1 - mu1 interpolation update
13108d5ead36SAlp Dener . -tao_bnk_mu2 - mu2 interpolation update
13118d5ead36SAlp Dener . -tao_bnk_gamma1 - gamma1 interpolation update
13128d5ead36SAlp Dener . -tao_bnk_gamma2 - gamma2 interpolation update
13138d5ead36SAlp Dener . -tao_bnk_gamma3 - gamma3 interpolation update
13148d5ead36SAlp Dener . -tao_bnk_gamma4 - gamma4 interpolation update
13158d5ead36SAlp Dener . -tao_bnk_theta - theta interpolation update
13168d5ead36SAlp Dener . -tao_bnk_omega1 - omega1 step update
13178d5ead36SAlp Dener . -tao_bnk_omega2 - omega2 step update
13188d5ead36SAlp Dener . -tao_bnk_omega3 - omega3 step update
13198d5ead36SAlp Dener . -tao_bnk_omega4 - omega4 step update
13208d5ead36SAlp Dener . -tao_bnk_omega5 - omega5 step update
13218d5ead36SAlp Dener . -tao_bnk_mu1_i -  mu1 interpolation init factor
13228d5ead36SAlp Dener . -tao_bnk_mu2_i -  mu2 interpolation init factor
13238d5ead36SAlp Dener . -tao_bnk_gamma1_i -  gamma1 interpolation init factor
13248d5ead36SAlp Dener . -tao_bnk_gamma2_i -  gamma2 interpolation init factor
13258d5ead36SAlp Dener . -tao_bnk_gamma3_i -  gamma3 interpolation init factor
13268d5ead36SAlp Dener . -tao_bnk_gamma4_i -  gamma4 interpolation init factor
13272f75a4aaSAlp Dener . -tao_bnk_theta_i -  theta interpolation init factor
13282f75a4aaSAlp Dener - -tao_bnk_bound_tol -  initial tolerance used in estimating bounded active variables
1329eb910715SAlp Dener 
1330eb910715SAlp Dener   Level: beginner
1331eb910715SAlp Dener M*/
1332eb910715SAlp Dener 
1333eb910715SAlp Dener PetscErrorCode TaoCreate_BNK(Tao tao)
1334eb910715SAlp Dener {
1335eb910715SAlp Dener   TAO_BNK        *bnk;
1336eb910715SAlp Dener   const char     *morethuente_type = TAOLINESEARCHMT;
1337eb910715SAlp Dener   PetscErrorCode ierr;
1338eb910715SAlp Dener 
1339eb910715SAlp Dener   PetscFunctionBegin;
1340eb910715SAlp Dener   ierr = PetscNewLog(tao,&bnk);CHKERRQ(ierr);
1341eb910715SAlp Dener 
1342eb910715SAlp Dener   tao->ops->setup = TaoSetUp_BNK;
1343eb910715SAlp Dener   tao->ops->view = TaoView_BNK;
1344eb910715SAlp Dener   tao->ops->setfromoptions = TaoSetFromOptions_BNK;
1345eb910715SAlp Dener   tao->ops->destroy = TaoDestroy_BNK;
1346eb910715SAlp Dener 
1347eb910715SAlp Dener   /*  Override default settings (unless already changed) */
1348eb910715SAlp Dener   if (!tao->max_it_changed) tao->max_it = 50;
1349eb910715SAlp Dener   if (!tao->trust0_changed) tao->trust0 = 100.0;
1350eb910715SAlp Dener 
1351eb910715SAlp Dener   tao->data = (void*)bnk;
1352eb910715SAlp Dener 
135366ed3702SAlp Dener   /*  Hessian shifting parameters */
1354eb910715SAlp Dener   bnk->sval   = 0.0;
1355eb910715SAlp Dener   bnk->imin   = 1.0e-4;
1356eb910715SAlp Dener   bnk->imax   = 1.0e+2;
1357eb910715SAlp Dener   bnk->imfac  = 1.0e-1;
1358eb910715SAlp Dener 
1359eb910715SAlp Dener   bnk->pmin   = 1.0e-12;
1360eb910715SAlp Dener   bnk->pmax   = 1.0e+2;
1361eb910715SAlp Dener   bnk->pgfac  = 1.0e+1;
1362eb910715SAlp Dener   bnk->psfac  = 4.0e-1;
1363eb910715SAlp Dener   bnk->pmgfac = 1.0e-1;
1364eb910715SAlp Dener   bnk->pmsfac = 1.0e-1;
1365eb910715SAlp Dener 
1366eb910715SAlp Dener   /*  Default values for trust-region radius update based on steplength */
1367eb910715SAlp Dener   bnk->nu1 = 0.25;
1368eb910715SAlp Dener   bnk->nu2 = 0.50;
1369eb910715SAlp Dener   bnk->nu3 = 1.00;
1370eb910715SAlp Dener   bnk->nu4 = 1.25;
1371eb910715SAlp Dener 
1372eb910715SAlp Dener   bnk->omega1 = 0.25;
1373eb910715SAlp Dener   bnk->omega2 = 0.50;
1374eb910715SAlp Dener   bnk->omega3 = 1.00;
1375eb910715SAlp Dener   bnk->omega4 = 2.00;
1376eb910715SAlp Dener   bnk->omega5 = 4.00;
1377eb910715SAlp Dener 
1378eb910715SAlp Dener   /*  Default values for trust-region radius update based on reduction */
1379eb910715SAlp Dener   bnk->eta1 = 1.0e-4;
1380eb910715SAlp Dener   bnk->eta2 = 0.25;
1381eb910715SAlp Dener   bnk->eta3 = 0.50;
1382eb910715SAlp Dener   bnk->eta4 = 0.90;
1383eb910715SAlp Dener 
1384eb910715SAlp Dener   bnk->alpha1 = 0.25;
1385eb910715SAlp Dener   bnk->alpha2 = 0.50;
1386eb910715SAlp Dener   bnk->alpha3 = 1.00;
1387eb910715SAlp Dener   bnk->alpha4 = 2.00;
1388eb910715SAlp Dener   bnk->alpha5 = 4.00;
1389eb910715SAlp Dener 
1390eb910715SAlp Dener   /*  Default values for trust-region radius update based on interpolation */
1391eb910715SAlp Dener   bnk->mu1 = 0.10;
1392eb910715SAlp Dener   bnk->mu2 = 0.50;
1393eb910715SAlp Dener 
1394eb910715SAlp Dener   bnk->gamma1 = 0.25;
1395eb910715SAlp Dener   bnk->gamma2 = 0.50;
1396eb910715SAlp Dener   bnk->gamma3 = 2.00;
1397eb910715SAlp Dener   bnk->gamma4 = 4.00;
1398eb910715SAlp Dener 
1399eb910715SAlp Dener   bnk->theta = 0.05;
1400eb910715SAlp Dener 
1401eb910715SAlp Dener   /*  Default values for trust region initialization based on interpolation */
1402eb910715SAlp Dener   bnk->mu1_i = 0.35;
1403eb910715SAlp Dener   bnk->mu2_i = 0.50;
1404eb910715SAlp Dener 
1405eb910715SAlp Dener   bnk->gamma1_i = 0.0625;
1406eb910715SAlp Dener   bnk->gamma2_i = 0.5;
1407eb910715SAlp Dener   bnk->gamma3_i = 2.0;
1408eb910715SAlp Dener   bnk->gamma4_i = 5.0;
1409eb910715SAlp Dener 
1410eb910715SAlp Dener   bnk->theta_i = 0.25;
1411eb910715SAlp Dener 
1412eb910715SAlp Dener   /*  Remaining parameters */
1413c0f10754SAlp Dener   bnk->max_cg_its = 0;
1414eb910715SAlp Dener   bnk->min_radius = 1.0e-10;
1415eb910715SAlp Dener   bnk->max_radius = 1.0e10;
1416770b7498SAlp Dener   bnk->epsilon = PetscPowReal(PETSC_MACHINE_EPSILON, 2.0/3.0);
14170a4511e9SAlp Dener   bnk->as_tol = 1.0e-3;
14180a4511e9SAlp Dener   bnk->as_step = 1.0e-3;
141962675beeSAlp Dener   bnk->dmin = 1.0e-6;
142062675beeSAlp Dener   bnk->dmax = 1.0e6;
1421eb910715SAlp Dener 
1422e031d6f5SAlp Dener   bnk->pc_type         = BNK_PC_AHESS;
1423eb910715SAlp Dener   bnk->bfgs_scale_type = BFGS_SCALE_PHESS;
1424eb910715SAlp Dener   bnk->init_type       = BNK_INIT_INTERPOLATION;
1425fed79b8eSAlp Dener   bnk->update_type     = BNK_UPDATE_INTERPOLATION;
14262f75a4aaSAlp Dener   bnk->as_type         = BNK_AS_BERTSEKAS;
1427eb910715SAlp Dener 
1428e031d6f5SAlp Dener   /* Create the embedded BNCG solver */
1429e031d6f5SAlp Dener   ierr = TaoCreate(PetscObjectComm((PetscObject)tao), &bnk->bncg);CHKERRQ(ierr);
1430e031d6f5SAlp Dener   ierr = PetscObjectIncrementTabLevel((PetscObject)bnk->bncg, (PetscObject)tao, 1);CHKERRQ(ierr);
1431e031d6f5SAlp Dener   ierr = TaoSetOptionsPrefix(bnk->bncg, "tao_bnk_");CHKERRQ(ierr);
1432e031d6f5SAlp Dener   ierr = TaoSetType(bnk->bncg, TAOBNCG);CHKERRQ(ierr);
1433e031d6f5SAlp Dener 
1434c0f10754SAlp Dener   /* Create the line search */
1435eb910715SAlp Dener   ierr = TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch);CHKERRQ(ierr);
1436eb910715SAlp Dener   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr);
1437e031d6f5SAlp Dener   ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr);
1438eb910715SAlp Dener   ierr = TaoLineSearchSetType(tao->linesearch,morethuente_type);CHKERRQ(ierr);
1439eb910715SAlp Dener   ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr);
1440eb910715SAlp Dener 
1441eb910715SAlp Dener   /*  Set linear solver to default for symmetric matrices */
1442eb910715SAlp Dener   ierr = KSPCreate(((PetscObject)tao)->comm,&tao->ksp);CHKERRQ(ierr);
1443eb910715SAlp Dener   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1);CHKERRQ(ierr);
1444eb910715SAlp Dener   ierr = KSPSetOptionsPrefix(tao->ksp,tao->hdr.prefix);CHKERRQ(ierr);
1445eb910715SAlp Dener   ierr = KSPSetType(tao->ksp,KSPCGSTCG);CHKERRQ(ierr);
1446eb910715SAlp Dener   PetscFunctionReturn(0);
1447eb910715SAlp Dener }
1448