xref: /petsc/src/tao/bound/impls/bnk/bnk.c (revision 62675bee48cd89cd8757a6de6d7b87f2edd3afcc)
1eb910715SAlp Dener #include <petsctaolinesearch.h>
2eb910715SAlp Dener #include <../src/tao/bound/impls/bnk/bnk.h>
3eb910715SAlp Dener 
4eb910715SAlp Dener #include <petscksp.h>
5eb910715SAlp Dener 
6eb910715SAlp Dener /* Routine for BFGS preconditioner */
7eb910715SAlp Dener 
8eb910715SAlp Dener PetscErrorCode MatLMVMSolveShell(PC pc, Vec b, Vec x)
9eb910715SAlp Dener {
10eb910715SAlp Dener   PetscErrorCode ierr;
11eb910715SAlp Dener   Mat            M;
12eb910715SAlp Dener 
13eb910715SAlp Dener   PetscFunctionBegin;
14eb910715SAlp Dener   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
15eb910715SAlp Dener   PetscValidHeaderSpecific(b,VEC_CLASSID,2);
16eb910715SAlp Dener   PetscValidHeaderSpecific(x,VEC_CLASSID,3);
17eb910715SAlp Dener   ierr = PCShellGetContext(pc,(void**)&M);CHKERRQ(ierr);
1809164190SAlp Dener   ierr = MatLMVMSolveInactive(M, b, x);CHKERRQ(ierr);
19eb910715SAlp Dener   PetscFunctionReturn(0);
20eb910715SAlp Dener }
21eb910715SAlp Dener 
22df278d8fSAlp Dener /*------------------------------------------------------------*/
23df278d8fSAlp Dener 
24df278d8fSAlp Dener /* Routine for initializing the KSP solver, the BFGS preconditioner, and the initial trust radius estimation */
25df278d8fSAlp Dener 
26*62675beeSAlp Dener PetscErrorCode TaoBNKInitialize(Tao tao, PetscInt initType)
27eb910715SAlp Dener {
28eb910715SAlp Dener   PetscErrorCode               ierr;
29eb910715SAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
30eb910715SAlp Dener   KSPType                      ksp_type;
31eb910715SAlp Dener   PC                           pc;
32eb910715SAlp Dener 
330a4511e9SAlp Dener   PetscReal                    f_min, ftrial, prered, actred, kappa, sigma;
34eb910715SAlp Dener   PetscReal                    tau, tau_1, tau_2, tau_max, tau_min, max_radius;
3528017e9fSAlp Dener   PetscReal                    delta;
36eb910715SAlp Dener 
37eb910715SAlp Dener   PetscInt                     n,N,needH = 1;
38eb910715SAlp Dener 
39eb910715SAlp Dener   PetscInt                     i_max = 5;
40eb910715SAlp Dener   PetscInt                     j_max = 1;
41eb910715SAlp Dener   PetscInt                     i, j;
42eb910715SAlp Dener 
43eb910715SAlp Dener   PetscFunctionBegin;
4428017e9fSAlp Dener   /*   Project the current point onto the feasible set */
4528017e9fSAlp Dener   ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr);
4628017e9fSAlp Dener   ierr = TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU);CHKERRQ(ierr);
4728017e9fSAlp Dener 
4828017e9fSAlp Dener   /* Project the initial point onto the feasible region */
4928017e9fSAlp Dener   ierr = VecMedian(tao->XL,tao->solution,tao->XU,tao->solution);CHKERRQ(ierr);
5028017e9fSAlp Dener 
5128017e9fSAlp Dener   /* Check convergence criteria */
5228017e9fSAlp Dener   ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &bnk->f, bnk->unprojected_gradient);CHKERRQ(ierr);
5328017e9fSAlp Dener   ierr = VecBoundGradientProjection(bnk->unprojected_gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr);
5428017e9fSAlp Dener   ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&bnk->gnorm);CHKERRQ(ierr);
5528017e9fSAlp Dener   if (PetscIsInfOrNanReal(bnk->f) || PetscIsInfOrNanReal(bnk->gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
5628017e9fSAlp Dener 
57eb910715SAlp Dener   /* Number of times ksp stopped because of these reasons */
58eb910715SAlp Dener   bnk->ksp_atol = 0;
59eb910715SAlp Dener   bnk->ksp_rtol = 0;
60eb910715SAlp Dener   bnk->ksp_dtol = 0;
61eb910715SAlp Dener   bnk->ksp_ctol = 0;
62eb910715SAlp Dener   bnk->ksp_negc = 0;
63eb910715SAlp Dener   bnk->ksp_iter = 0;
64eb910715SAlp Dener   bnk->ksp_othr = 0;
65eb910715SAlp Dener 
66fed79b8eSAlp Dener   /* Initialize the Hessian perturbation */
67fed79b8eSAlp Dener   bnk->pert = bnk->sval;
68fed79b8eSAlp Dener 
69eb910715SAlp Dener   /* Initialize trust-region radius when using nash, stcg, or gltr
70eb910715SAlp Dener      Command automatically ignored for other methods
71eb910715SAlp Dener      Will be reset during the first iteration
72eb910715SAlp Dener   */
73eb910715SAlp Dener   ierr = KSPGetType(tao->ksp,&ksp_type);CHKERRQ(ierr);
74eb910715SAlp Dener   ierr = PetscStrcmp(ksp_type,KSPCGNASH,&bnk->is_nash);CHKERRQ(ierr);
75eb910715SAlp Dener   ierr = PetscStrcmp(ksp_type,KSPCGSTCG,&bnk->is_stcg);CHKERRQ(ierr);
76eb910715SAlp Dener   ierr = PetscStrcmp(ksp_type,KSPCGGLTR,&bnk->is_gltr);CHKERRQ(ierr);
77eb910715SAlp Dener 
78eb910715SAlp Dener   ierr = KSPCGSetRadius(tao->ksp,bnk->max_radius);CHKERRQ(ierr);
79eb910715SAlp Dener 
80eb910715SAlp Dener   if (bnk->is_nash || bnk->is_stcg || bnk->is_gltr) {
81eb910715SAlp Dener     if (tao->trust0 < 0.0) SETERRQ(PETSC_COMM_SELF,1,"Initial radius negative");
82eb910715SAlp Dener     tao->trust = tao->trust0;
83eb910715SAlp Dener     tao->trust = PetscMax(tao->trust, bnk->min_radius);
84eb910715SAlp Dener     tao->trust = PetscMin(tao->trust, bnk->max_radius);
85eb910715SAlp Dener   }
86eb910715SAlp Dener 
8728017e9fSAlp Dener   ierr = TaoLogConvergenceHistory(tao,bnk->f,bnk->gnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
8828017e9fSAlp Dener   ierr = TaoMonitor(tao,tao->niter,bnk->f,bnk->gnorm,0.0,1.0);CHKERRQ(ierr);
8928017e9fSAlp Dener   ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
9028017e9fSAlp Dener   if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
9128017e9fSAlp Dener 
92eb910715SAlp Dener   /* Get vectors we will need */
93eb910715SAlp Dener   if (BNK_PC_BFGS == bnk->pc_type && !bnk->M) {
94eb910715SAlp Dener     ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr);
95eb910715SAlp Dener     ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr);
96eb910715SAlp Dener     ierr = MatCreateLMVM(((PetscObject)tao)->comm,n,N,&bnk->M);CHKERRQ(ierr);
97eb910715SAlp Dener     ierr = MatLMVMAllocateVectors(bnk->M,tao->solution);CHKERRQ(ierr);
98eb910715SAlp Dener   }
99eb910715SAlp Dener 
100eb910715SAlp Dener   /* create vectors for the limited memory preconditioner */
101eb910715SAlp Dener   if ((BNK_PC_BFGS == bnk->pc_type) && (BFGS_SCALE_BFGS != bnk->bfgs_scale_type)) {
102*62675beeSAlp Dener     if (!bnk->Diag) {ierr = VecDuplicate(tao->solution,&bnk->Diag);CHKERRQ(ierr);}
103*62675beeSAlp Dener     if (!bnk->Diag_min) {ierr = VecDuplicate(tao->solution,&bnk->Diag_min);CHKERRQ(ierr);}
104*62675beeSAlp Dener     if (!bnk->Diag_max) {ierr = VecDuplicate(tao->solution,&bnk->Diag_max);CHKERRQ(ierr);}
105*62675beeSAlp Dener     ierr = VecSet(bnk->Diag_min, bnk->dmin);CHKERRQ(ierr);
106*62675beeSAlp Dener     ierr = VecSet(bnk->Diag_max, bnk->dmax);CHKERRQ(ierr);
107eb910715SAlp Dener   }
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. */
138eb910715SAlp Dener   if (bnk->is_nash || bnk->is_stcg || bnk->is_gltr) {
139*62675beeSAlp Dener     switch(initType) {
140eb910715SAlp Dener     case BNK_INIT_CONSTANT:
141eb910715SAlp Dener       /* Use the initial radius specified */
142eb910715SAlp Dener       break;
143eb910715SAlp Dener 
144eb910715SAlp Dener     case BNK_INIT_INTERPOLATION:
145eb910715SAlp Dener       /* Use the initial radius specified */
146eb910715SAlp Dener       max_radius = 0.0;
147eb910715SAlp Dener 
148eb910715SAlp Dener       for (j = 0; j < j_max; ++j) {
1490a4511e9SAlp Dener         f_min = bnk->f;
150eb910715SAlp Dener         sigma = 0.0;
151eb910715SAlp Dener 
152eb910715SAlp Dener         if (needH) {
15362602cfbSAlp Dener           /* Compute the Hessian at the new step, and extract the inactive subsystem */
154eb910715SAlp Dener           ierr = TaoComputeHessian(tao, tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr);
15528017e9fSAlp Dener           if (bnk->inactive_idx) {
15662602cfbSAlp Dener             ierr = ISDestroy(&bnk->inactive_idx);CHKERRQ(ierr);
15762602cfbSAlp Dener             ierr = VecWhichInactive(tao->XL,tao->solution,bnk->unprojected_gradient,tao->XU,PETSC_TRUE,&bnk->inactive_idx);CHKERRQ(ierr);
15862602cfbSAlp Dener             ierr = TaoMatGetSubMat(tao->hessian, bnk->inactive_idx, bnk->Xwork, TAO_SUBSET_MASK, &bnk->H_inactive);CHKERRQ(ierr);
15928017e9fSAlp Dener           } else {
160*62675beeSAlp Dener             ierr = MatDestroy(&bnk->H_inactive);
161*62675beeSAlp Dener             ierr = MatDuplicate(tao->hessian, MAT_COPY_VALUES, &bnk->H_inactive);
16228017e9fSAlp Dener           }
163eb910715SAlp Dener           needH = 0;
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             }
184770b7498SAlp Dener             /* Compute the predicted and actual reduction */
18562602cfbSAlp Dener             ierr = MatMult(bnk->H_inactive, tao->gradient, bnk->W);CHKERRQ(ierr);
18662602cfbSAlp Dener             ierr = VecDot(tao->gradient, bnk->W, &prered);CHKERRQ(ierr);
187eb910715SAlp Dener             prered = tao->trust * (bnk->gnorm - 0.5 * tao->trust * prered / (bnk->gnorm * bnk->gnorm));
188eb910715SAlp Dener             actred = bnk->f - ftrial;
189eb910715SAlp Dener             if ((PetscAbsScalar(actred) <= bnk->epsilon) && (PetscAbsScalar(prered) <= bnk->epsilon)) {
190eb910715SAlp Dener               kappa = 1.0;
191eb910715SAlp Dener             } else {
192eb910715SAlp Dener               kappa = actred / prered;
193eb910715SAlp Dener             }
194eb910715SAlp Dener 
195eb910715SAlp Dener             tau_1 = bnk->theta_i * bnk->gnorm * tao->trust / (bnk->theta_i * bnk->gnorm * tao->trust + (1.0 - bnk->theta_i) * prered - actred);
196eb910715SAlp Dener             tau_2 = bnk->theta_i * bnk->gnorm * tao->trust / (bnk->theta_i * bnk->gnorm * tao->trust - (1.0 + bnk->theta_i) * prered + actred);
197eb910715SAlp Dener             tau_min = PetscMin(tau_1, tau_2);
198eb910715SAlp Dener             tau_max = PetscMax(tau_1, tau_2);
199eb910715SAlp Dener 
200eb910715SAlp Dener             if (PetscAbsScalar(kappa - 1.0) <= bnk->mu1_i) {
201eb910715SAlp Dener               /* Great agreement */
202eb910715SAlp Dener               max_radius = PetscMax(max_radius, tao->trust);
203eb910715SAlp Dener 
204eb910715SAlp Dener               if (tau_max < 1.0) {
205eb910715SAlp Dener                 tau = bnk->gamma3_i;
206eb910715SAlp Dener               } else if (tau_max > bnk->gamma4_i) {
207eb910715SAlp Dener                 tau = bnk->gamma4_i;
208eb910715SAlp Dener               } else if (tau_1 >= 1.0 && tau_1 <= bnk->gamma4_i && tau_2 < 1.0) {
209eb910715SAlp Dener                 tau = tau_1;
210eb910715SAlp Dener               } else if (tau_2 >= 1.0 && tau_2 <= bnk->gamma4_i && tau_1 < 1.0) {
211eb910715SAlp Dener                 tau = tau_2;
212eb910715SAlp Dener               } else {
213eb910715SAlp Dener                 tau = tau_max;
214eb910715SAlp Dener               }
215eb910715SAlp Dener             } else if (PetscAbsScalar(kappa - 1.0) <= bnk->mu2_i) {
216eb910715SAlp Dener               /* Good agreement */
217eb910715SAlp Dener               max_radius = PetscMax(max_radius, tao->trust);
218eb910715SAlp Dener 
219eb910715SAlp Dener               if (tau_max < bnk->gamma2_i) {
220eb910715SAlp Dener                 tau = bnk->gamma2_i;
221eb910715SAlp Dener               } else if (tau_max > bnk->gamma3_i) {
222eb910715SAlp Dener                 tau = bnk->gamma3_i;
223eb910715SAlp Dener               } else {
224eb910715SAlp Dener                 tau = tau_max;
225eb910715SAlp Dener               }
226eb910715SAlp Dener             } else {
227eb910715SAlp Dener               /* Not good agreement */
228eb910715SAlp Dener               if (tau_min > 1.0) {
229eb910715SAlp Dener                 tau = bnk->gamma2_i;
230eb910715SAlp Dener               } else if (tau_max < bnk->gamma1_i) {
231eb910715SAlp Dener                 tau = bnk->gamma1_i;
232eb910715SAlp Dener               } else if ((tau_min < bnk->gamma1_i) && (tau_max >= 1.0)) {
233eb910715SAlp Dener                 tau = bnk->gamma1_i;
234eb910715SAlp Dener               } else if ((tau_1 >= bnk->gamma1_i) && (tau_1 < 1.0) && ((tau_2 < bnk->gamma1_i) || (tau_2 >= 1.0))) {
235eb910715SAlp Dener                 tau = tau_1;
236eb910715SAlp Dener               } else if ((tau_2 >= bnk->gamma1_i) && (tau_2 < 1.0) && ((tau_1 < bnk->gamma1_i) || (tau_2 >= 1.0))) {
237eb910715SAlp Dener                 tau = tau_2;
238eb910715SAlp Dener               } else {
239eb910715SAlp Dener                 tau = tau_max;
240eb910715SAlp Dener               }
241eb910715SAlp Dener             }
242eb910715SAlp Dener           }
243eb910715SAlp Dener           tao->trust = tau * tao->trust;
244eb910715SAlp Dener         }
245eb910715SAlp Dener 
2460a4511e9SAlp Dener         if (f_min < bnk->f) {
2470a4511e9SAlp Dener           bnk->f = f_min;
248eb910715SAlp Dener           ierr = VecAXPY(tao->solution,sigma,tao->gradient);CHKERRQ(ierr);
24987602d1eSAlp Dener           ierr = VecMedian(tao->XL, tao->solution, tao->XU, tao->solution);CHKERRQ(ierr);
25009164190SAlp Dener           ierr = TaoComputeGradient(tao,tao->solution,bnk->unprojected_gradient);CHKERRQ(ierr);
25109164190SAlp Dener           ierr = VecBoundGradientProjection(bnk->unprojected_gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr);
252eb910715SAlp Dener 
253eb910715SAlp Dener           ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&bnk->gnorm);CHKERRQ(ierr);
254eb910715SAlp Dener           if (PetscIsInfOrNanReal(bnk->gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute gradient generated Inf or NaN");
255eb910715SAlp Dener           needH = 1;
256eb910715SAlp Dener 
257eb910715SAlp Dener           ierr = TaoLogConvergenceHistory(tao,bnk->f,bnk->gnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
25828017e9fSAlp Dener           ierr = TaoMonitor(tao,tao->niter,bnk->f,bnk->gnorm,0.0,1.0);CHKERRQ(ierr);
259eb910715SAlp Dener           ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
260eb910715SAlp Dener           if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
261eb910715SAlp Dener         }
262eb910715SAlp Dener       }
263eb910715SAlp Dener       tao->trust = PetscMax(tao->trust, max_radius);
264eb910715SAlp Dener 
265eb910715SAlp Dener       /* Modify the radius if it is too large or small */
266eb910715SAlp Dener       tao->trust = PetscMax(tao->trust, bnk->min_radius);
267eb910715SAlp Dener       tao->trust = PetscMin(tao->trust, bnk->max_radius);
268eb910715SAlp Dener       break;
269eb910715SAlp Dener 
270eb910715SAlp Dener     default:
271eb910715SAlp Dener       /* Norm of the first direction will initialize radius */
272eb910715SAlp Dener       tao->trust = 0.0;
273eb910715SAlp Dener       break;
274eb910715SAlp Dener     }
275eb910715SAlp Dener   }
276eb910715SAlp Dener 
277eb910715SAlp Dener   /* Set initial scaling for the BFGS preconditioner
278eb910715SAlp Dener      This step is done after computing the initial trust-region radius
279eb910715SAlp Dener      since the function value may have decreased */
280eb910715SAlp Dener   if (BNK_PC_BFGS == bnk->pc_type) {
281eb910715SAlp Dener     if (bnk->f != 0.0) {
282eb910715SAlp Dener       delta = 2.0 * PetscAbsScalar(bnk->f) / (bnk->gnorm*bnk->gnorm);
283eb910715SAlp Dener     } else {
284eb910715SAlp Dener       delta = 2.0 / (bnk->gnorm*bnk->gnorm);
285eb910715SAlp Dener     }
286eb910715SAlp Dener     ierr = MatLMVMSetDelta(bnk->M,delta);CHKERRQ(ierr);
287eb910715SAlp Dener   }
288eb910715SAlp Dener 
289eb910715SAlp Dener   /* Set counter for gradient/reset steps*/
290eb910715SAlp Dener   bnk->newt = 0;
291eb910715SAlp Dener   bnk->bfgs = 0;
292eb910715SAlp Dener   bnk->sgrad = 0;
293eb910715SAlp Dener   bnk->grad = 0;
294eb910715SAlp Dener   PetscFunctionReturn(0);
295eb910715SAlp Dener }
296eb910715SAlp Dener 
297df278d8fSAlp Dener /*------------------------------------------------------------*/
298df278d8fSAlp Dener 
299*62675beeSAlp Dener /* Routine for computing the Hessian and preparing the preconditioner at the new iterate */
300*62675beeSAlp Dener 
301*62675beeSAlp Dener PetscErrorCode TaoBNKComputeHessian(Tao tao)
302*62675beeSAlp Dener {
303*62675beeSAlp Dener   PetscErrorCode               ierr;
304*62675beeSAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
305*62675beeSAlp Dener 
306*62675beeSAlp Dener   PetscFunctionBegin;
307*62675beeSAlp Dener   /* Compute the Hessian */
308*62675beeSAlp Dener   ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr);
309*62675beeSAlp Dener   /* Add a correction to the BFGS preconditioner */
310*62675beeSAlp Dener   if (BNK_PC_BFGS == bnk->pc_type) {
311*62675beeSAlp Dener     ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
312*62675beeSAlp Dener     /* Update the BFGS diagonal scaling */
313*62675beeSAlp Dener     if (BFGS_SCALE_AHESS == bnk->bfgs_scale_type) {
314*62675beeSAlp Dener       ierr = MatGetDiagonal(tao->hessian, bnk->Diag);CHKERRQ(ierr);
315*62675beeSAlp Dener       ierr = VecAbs(bnk->Diag);CHKERRQ(ierr);
316*62675beeSAlp Dener       ierr = VecMedian(bnk->Diag_min, bnk->Diag, bnk->Diag_max, bnk->Diag);CHKERRQ(ierr);
317*62675beeSAlp Dener       ierr = VecReciprocal(bnk->Diag);CHKERRQ(ierr);
318*62675beeSAlp Dener       ierr = MatLMVMSetScale(bnk->M,bnk->Diag);CHKERRQ(ierr);
319*62675beeSAlp Dener     }
320*62675beeSAlp Dener   }
321*62675beeSAlp Dener   PetscFunctionReturn(0);
322*62675beeSAlp Dener }
323*62675beeSAlp Dener 
324*62675beeSAlp Dener /*------------------------------------------------------------*/
325*62675beeSAlp Dener 
3262f75a4aaSAlp Dener /* Routine for estimating the active set */
3272f75a4aaSAlp Dener 
3282f75a4aaSAlp Dener PetscErrorCode TaoBNKEstimateActiveSet(Tao tao)
3292f75a4aaSAlp Dener {
3302f75a4aaSAlp Dener   PetscErrorCode               ierr;
3312f75a4aaSAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
3322f75a4aaSAlp Dener 
3332f75a4aaSAlp Dener   PetscFunctionBegin;
3342f75a4aaSAlp Dener   switch (bnk->as_type) {
3352f75a4aaSAlp Dener   case BNK_AS_NONE:
3362f75a4aaSAlp Dener     ierr = ISDestroy(&bnk->inactive_idx);CHKERRQ(ierr);
3372f75a4aaSAlp Dener     ierr = VecWhichInactive(tao->XL, tao->solution, bnk->unprojected_gradient, tao->XU, PETSC_TRUE, &bnk->inactive_idx);CHKERRQ(ierr);
3382f75a4aaSAlp Dener     ierr = ISDestroy(&bnk->active_idx);CHKERRQ(ierr);
3392f75a4aaSAlp Dener     ierr = ISComplementVec(bnk->inactive_idx, tao->solution, &bnk->active_idx);CHKERRQ(ierr);
3402f75a4aaSAlp Dener     break;
3412f75a4aaSAlp Dener 
3422f75a4aaSAlp Dener   case BNK_AS_BERTSEKAS:
3432f75a4aaSAlp Dener     /* Compute the trial step vector with which we will estimate the active set at the next iteration */
3442f75a4aaSAlp Dener     if (BNK_PC_BFGS == bnk->pc_type) {
3452f75a4aaSAlp Dener       /* If the BFGS preconditioner matrix is available, we will construct a trial step with it */
3462f75a4aaSAlp Dener       ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, bnk->W);CHKERRQ(ierr);
3472f75a4aaSAlp Dener     } else {
3482f75a4aaSAlp Dener       /* BFGS preconditioner doesn't exist so let's invert the diagonal of the Hessian instead onto the gradient*/
3492f75a4aaSAlp Dener       ierr = MatGetDiagonal(tao->hessian, bnk->Xwork);CHKERRQ(ierr);
3502f75a4aaSAlp Dener       ierr = VecReciprocal(bnk->Xwork);CHKERRQ(ierr);CHKERRQ(ierr);
3512f75a4aaSAlp Dener       ierr = VecPointwiseMult(bnk->W, bnk->Xwork, bnk->unprojected_gradient);CHKERRQ(ierr);
3522f75a4aaSAlp Dener     }
3532f75a4aaSAlp Dener     ierr = VecScale(bnk->W, -1.0);CHKERRQ(ierr);
3540a4511e9SAlp 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);
3552f75a4aaSAlp Dener 
3562f75a4aaSAlp Dener   default:
3572f75a4aaSAlp Dener     break;
3582f75a4aaSAlp Dener   }
3592f75a4aaSAlp Dener   PetscFunctionReturn(0);
3602f75a4aaSAlp Dener }
3612f75a4aaSAlp Dener 
3622f75a4aaSAlp Dener /*------------------------------------------------------------*/
3632f75a4aaSAlp Dener 
3642f75a4aaSAlp Dener /* Routine for bounding the step direction */
3652f75a4aaSAlp Dener 
3662f75a4aaSAlp Dener PetscErrorCode TaoBNKBoundStep(Tao tao, Vec step)
3672f75a4aaSAlp Dener {
3682f75a4aaSAlp Dener   PetscErrorCode               ierr;
3692f75a4aaSAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
3702f75a4aaSAlp Dener 
3712f75a4aaSAlp Dener   PetscFunctionBegin;
37228017e9fSAlp Dener   if (bnk->active_idx) {
3732f75a4aaSAlp Dener     switch (bnk->as_type) {
3742f75a4aaSAlp Dener     case BNK_AS_NONE:
3752f75a4aaSAlp Dener       if (bnk->active_idx) {
3762f75a4aaSAlp Dener         ierr = VecGetSubVector(step, bnk->active_idx, &bnk->active_work);CHKERRQ(ierr);
3772f75a4aaSAlp Dener         ierr = VecSet(bnk->active_work, 0.0);CHKERRQ(ierr);
3782f75a4aaSAlp Dener         ierr = VecRestoreSubVector(step, bnk->active_idx, &bnk->active_work);CHKERRQ(ierr);
3792f75a4aaSAlp Dener       }
3802f75a4aaSAlp Dener       break;
3812f75a4aaSAlp Dener 
3822f75a4aaSAlp Dener     case BNK_AS_BERTSEKAS:
3832f75a4aaSAlp Dener       ierr = TaoBoundStep(tao->solution, tao->XL, tao->XU, bnk->active_lower, bnk->active_upper, bnk->active_fixed, step);CHKERRQ(ierr);
3842f75a4aaSAlp Dener       break;
3852f75a4aaSAlp Dener 
3862f75a4aaSAlp Dener     default:
3872f75a4aaSAlp Dener       break;
3882f75a4aaSAlp Dener     }
389e465cd6fSAlp Dener   }
3902f75a4aaSAlp Dener   PetscFunctionReturn(0);
3912f75a4aaSAlp Dener }
3922f75a4aaSAlp Dener 
3932f75a4aaSAlp Dener /*------------------------------------------------------------*/
3942f75a4aaSAlp Dener 
395df278d8fSAlp Dener /* Routine for computing the Newton step.
396df278d8fSAlp Dener 
397df278d8fSAlp Dener   If the safeguard is enabled, the Newton step is verified to be a
398df278d8fSAlp Dener   descent direction, with fallbacks onto BFGS, scaled gradient, and unscaled
399df278d8fSAlp Dener   gradient steps if/when necessary.
400df278d8fSAlp Dener 
401df278d8fSAlp Dener   The function reports back on which type of step has ultimately been stored
402df278d8fSAlp Dener   under tao->stepdirection.
403df278d8fSAlp Dener */
404df278d8fSAlp Dener 
405*62675beeSAlp Dener PetscErrorCode TaoBNKComputeStep(Tao tao, PetscBool shift, KSPConvergedReason *ksp_reason)
406eb910715SAlp Dener {
407eb910715SAlp Dener   PetscErrorCode               ierr;
408eb910715SAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
409eb910715SAlp Dener 
410e465cd6fSAlp Dener   PetscReal                    delta;
411eb910715SAlp Dener   PetscInt                     bfgsUpdates = 0;
412eb910715SAlp Dener   PetscInt                     kspits;
413eb910715SAlp Dener 
414eb910715SAlp Dener   PetscFunctionBegin;
4152f75a4aaSAlp Dener   /* Determine the active and inactive sets */
4162f75a4aaSAlp Dener   ierr = TaoBNKEstimateActiveSet(tao);CHKERRQ(ierr);
41709164190SAlp Dener 
418eb3ba6a7SAlp Dener   /* Prepare masked matrices for the inactive set */
4192f75a4aaSAlp Dener   if (BNK_PC_BFGS == bnk->pc_type) { ierr = MatLMVMSetInactive(bnk->M, bnk->inactive_idx);CHKERRQ(ierr); }
4202f75a4aaSAlp Dener   if (bnk->inactive_idx) {
421eb3ba6a7SAlp Dener     ierr = TaoMatGetSubMat(tao->hessian, bnk->inactive_idx, bnk->Xwork, TAO_SUBSET_MASK, &bnk->H_inactive);CHKERRQ(ierr);
422eb3ba6a7SAlp Dener     if (tao->hessian == tao->hessian_pre) {
423eb3ba6a7SAlp Dener       bnk->Hpre_inactive = bnk->H_inactive;
424eb3ba6a7SAlp Dener     } else {
425eb3ba6a7SAlp Dener       ierr = TaoMatGetSubMat(tao->hessian_pre, bnk->inactive_idx, bnk->Xwork, TAO_SUBSET_MASK, &bnk->Hpre_inactive);CHKERRQ(ierr);
426eb3ba6a7SAlp Dener     }
4272f75a4aaSAlp Dener   } else {
428*62675beeSAlp Dener     ierr = MatDestroy(&bnk->H_inactive);
429*62675beeSAlp Dener     ierr = MatDuplicate(tao->hessian, MAT_COPY_VALUES, &bnk->H_inactive);
430*62675beeSAlp Dener     if (tao->hessian == tao->hessian_pre) {
431*62675beeSAlp Dener       bnk->Hpre_inactive = bnk->H_inactive;
432*62675beeSAlp Dener     } else {
433*62675beeSAlp Dener       ierr = MatDestroy(&bnk->Hpre_inactive);
434*62675beeSAlp Dener       ierr = MatDuplicate(tao->hessian_pre, MAT_COPY_VALUES, &bnk->Hpre_inactive);
435*62675beeSAlp Dener     }
436*62675beeSAlp Dener   }
437*62675beeSAlp Dener 
438*62675beeSAlp Dener   /* Shift the reduced Hessian matrix */
439*62675beeSAlp Dener   if ((shift) && (bnk->pert > 0)) {
440*62675beeSAlp Dener     ierr = MatShift(bnk->H_inactive, bnk->pert);CHKERRQ(ierr);
441*62675beeSAlp Dener     if (bnk->H_inactive != bnk->Hpre_inactive) {
442*62675beeSAlp Dener       ierr = MatShift(bnk->Hpre_inactive, bnk->pert);CHKERRQ(ierr);
443*62675beeSAlp Dener     }
444*62675beeSAlp Dener   }
445*62675beeSAlp Dener 
446*62675beeSAlp Dener   /* Update the diagonal scaling for the BFGS preconditioner, this time with the Hessian perturbation */
447*62675beeSAlp Dener   if ((BNK_PC_BFGS == bnk->pc_type) && (BFGS_SCALE_PHESS == bnk->bfgs_scale_type)) {
448*62675beeSAlp Dener     /* Obtain diagonal for the bfgs preconditioner  */
449*62675beeSAlp Dener     ierr = MatGetDiagonal(bnk->H_inactive, bnk->Diag);CHKERRQ(ierr);
450*62675beeSAlp Dener     ierr = VecAbs(bnk->Diag);CHKERRQ(ierr);
451*62675beeSAlp Dener     ierr = VecMedian(bnk->Diag_min, bnk->Diag, bnk->Diag_max, bnk->Diag);CHKERRQ(ierr);
452*62675beeSAlp Dener     ierr = VecReciprocal(bnk->Diag);CHKERRQ(ierr);
453*62675beeSAlp Dener     ierr = MatLMVMSetScale(bnk->M,bnk->Diag);CHKERRQ(ierr);
4542f75a4aaSAlp Dener   }
45509164190SAlp Dener 
456eb910715SAlp Dener   /* Solve the Newton system of equations */
4572f75a4aaSAlp Dener   ierr = VecSet(tao->stepdirection, 0.0);CHKERRQ(ierr);
45809164190SAlp Dener   ierr = KSPSetOperators(tao->ksp,bnk->H_inactive,bnk->Hpre_inactive);CHKERRQ(ierr);
45928017e9fSAlp Dener   ierr = VecCopy(bnk->unprojected_gradient, bnk->Gwork);CHKERRQ(ierr);
46028017e9fSAlp Dener   if (bnk->active_idx) {
46128017e9fSAlp Dener     ierr = VecGetSubVector(bnk->Gwork, bnk->active_idx, &bnk->active_work);CHKERRQ(ierr);
46228017e9fSAlp Dener     ierr = VecSet(bnk->active_work, 0.0);CHKERRQ(ierr);
46328017e9fSAlp Dener     ierr = VecRestoreSubVector(bnk->Gwork, bnk->active_idx, &bnk->active_work);CHKERRQ(ierr);
46428017e9fSAlp Dener   }
465eb910715SAlp Dener   if (bnk->is_nash || bnk->is_stcg || bnk->is_gltr) {
466fed79b8eSAlp Dener     ierr = KSPCGSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr);
46728017e9fSAlp Dener     ierr = KSPSolve(tao->ksp, bnk->Gwork, tao->stepdirection);CHKERRQ(ierr);
468eb910715SAlp Dener     ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr);
469eb910715SAlp Dener     tao->ksp_its+=kspits;
470eb910715SAlp Dener     tao->ksp_tot_its+=kspits;
471080d2917SAlp Dener     ierr = KSPCGGetNormD(tao->ksp,&bnk->dnorm);CHKERRQ(ierr);
472eb910715SAlp Dener 
473eb910715SAlp Dener     if (0.0 == tao->trust) {
474eb910715SAlp Dener       /* Radius was uninitialized; use the norm of the direction */
475080d2917SAlp Dener       if (bnk->dnorm > 0.0) {
476080d2917SAlp Dener         tao->trust = bnk->dnorm;
477eb910715SAlp Dener 
478eb910715SAlp Dener         /* Modify the radius if it is too large or small */
479eb910715SAlp Dener         tao->trust = PetscMax(tao->trust, bnk->min_radius);
480eb910715SAlp Dener         tao->trust = PetscMin(tao->trust, bnk->max_radius);
481eb910715SAlp Dener       } else {
482eb910715SAlp Dener         /* The direction was bad; set radius to default value and re-solve
483eb910715SAlp Dener            the trust-region subproblem to get a direction */
484eb910715SAlp Dener         tao->trust = tao->trust0;
485eb910715SAlp Dener 
486eb910715SAlp Dener         /* Modify the radius if it is too large or small */
487eb910715SAlp Dener         tao->trust = PetscMax(tao->trust, bnk->min_radius);
488eb910715SAlp Dener         tao->trust = PetscMin(tao->trust, bnk->max_radius);
489eb910715SAlp Dener 
490fed79b8eSAlp Dener         ierr = KSPCGSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr);
49128017e9fSAlp Dener         ierr = KSPSolve(tao->ksp, bnk->Gwork, tao->stepdirection);CHKERRQ(ierr);
492eb910715SAlp Dener         ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr);
493eb910715SAlp Dener         tao->ksp_its+=kspits;
494eb910715SAlp Dener         tao->ksp_tot_its+=kspits;
495080d2917SAlp Dener         ierr = KSPCGGetNormD(tao->ksp,&bnk->dnorm);CHKERRQ(ierr);
496eb910715SAlp Dener 
497080d2917SAlp Dener         if (bnk->dnorm == 0.0) SETERRQ(PETSC_COMM_SELF,1, "Initial direction zero");
498eb910715SAlp Dener       }
499eb910715SAlp Dener     }
500eb910715SAlp Dener   } else {
50128017e9fSAlp Dener     ierr = KSPSolve(tao->ksp, bnk->Gwork, tao->stepdirection);CHKERRQ(ierr);
502eb910715SAlp Dener     ierr = KSPGetIterationNumber(tao->ksp, &kspits);CHKERRQ(ierr);
503eb910715SAlp Dener     tao->ksp_its += kspits;
504eb910715SAlp Dener     tao->ksp_tot_its+=kspits;
505eb910715SAlp Dener   }
506770b7498SAlp Dener   /* Make sure the safeguarded fall-back step is zero for actively bounded variables */
507fed79b8eSAlp Dener   ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
5082f75a4aaSAlp Dener   ierr = TaoBNKBoundStep(tao, tao->stepdirection);CHKERRQ(ierr);
509770b7498SAlp Dener 
510770b7498SAlp Dener   /* Record convergence reasons */
511e465cd6fSAlp Dener   ierr = KSPGetConvergedReason(tao->ksp, ksp_reason);CHKERRQ(ierr);
512e465cd6fSAlp Dener   if (KSP_CONVERGED_ATOL == *ksp_reason) {
513770b7498SAlp Dener     ++bnk->ksp_atol;
514e465cd6fSAlp Dener   } else if (KSP_CONVERGED_RTOL == *ksp_reason) {
515770b7498SAlp Dener     ++bnk->ksp_rtol;
516e465cd6fSAlp Dener   } else if (KSP_CONVERGED_CG_CONSTRAINED == *ksp_reason) {
517770b7498SAlp Dener     ++bnk->ksp_ctol;
518e465cd6fSAlp Dener   } else if (KSP_CONVERGED_CG_NEG_CURVE == *ksp_reason) {
519770b7498SAlp Dener     ++bnk->ksp_negc;
520e465cd6fSAlp Dener   } else if (KSP_DIVERGED_DTOL == *ksp_reason) {
521770b7498SAlp Dener     ++bnk->ksp_dtol;
522e465cd6fSAlp Dener   } else if (KSP_DIVERGED_ITS == *ksp_reason) {
523770b7498SAlp Dener     ++bnk->ksp_iter;
524770b7498SAlp Dener   } else {
525770b7498SAlp Dener     ++bnk->ksp_othr;
526770b7498SAlp Dener   }
527fed79b8eSAlp Dener 
528fed79b8eSAlp Dener   /* Make sure the BFGS preconditioner is healthy */
529fed79b8eSAlp Dener   if (bnk->pc_type == BNK_PC_BFGS) {
530fed79b8eSAlp Dener     ierr = MatLMVMGetUpdates(bnk->M, &bfgsUpdates);CHKERRQ(ierr);
531e465cd6fSAlp Dener     if ((KSP_DIVERGED_INDEFINITE_PC == *ksp_reason) && (bfgsUpdates > 1)) {
532fed79b8eSAlp Dener       /* Preconditioner is numerically indefinite; reset the approximation. */
533eb910715SAlp Dener       if (bnk->f != 0.0) {
534eb910715SAlp Dener         delta = 2.0 * PetscAbsScalar(bnk->f) / (bnk->gnorm*bnk->gnorm);
535eb910715SAlp Dener       } else {
536eb910715SAlp Dener         delta = 2.0 / (bnk->gnorm*bnk->gnorm);
537eb910715SAlp Dener       }
538eb910715SAlp Dener       ierr = MatLMVMSetDelta(bnk->M,delta);CHKERRQ(ierr);
539eb910715SAlp Dener       ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr);
54009164190SAlp Dener       ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
541eb910715SAlp Dener     }
542fed79b8eSAlp Dener   }
543e465cd6fSAlp Dener   PetscFunctionReturn(0);
544e465cd6fSAlp Dener }
545eb910715SAlp Dener 
546*62675beeSAlp Dener /*------------------------------------------------------------*/
547*62675beeSAlp Dener 
548*62675beeSAlp Dener /* Routine for ensuring that the Newton step is a descent direction.
549*62675beeSAlp Dener 
550*62675beeSAlp Dener    The step direction falls back onto BFGS, scaled gradient and gradient steps
551*62675beeSAlp Dener    in the event that the Newton step fails the test.
552*62675beeSAlp Dener */
553*62675beeSAlp Dener 
554e465cd6fSAlp Dener PetscErrorCode TaoBNKSafeguardStep(Tao tao, KSPConvergedReason ksp_reason, PetscInt *stepType)
555e465cd6fSAlp Dener {
556e465cd6fSAlp Dener   PetscErrorCode               ierr;
557e465cd6fSAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
558e465cd6fSAlp Dener 
559e465cd6fSAlp Dener   PetscReal                    gdx, delta, e_min;
560e465cd6fSAlp Dener   PetscInt                     bfgsUpdates;
561e465cd6fSAlp Dener 
562e465cd6fSAlp Dener   PetscFunctionBegin;
563080d2917SAlp Dener   ierr = VecDot(tao->stepdirection, tao->gradient, &gdx);CHKERRQ(ierr);
564eb910715SAlp Dener   if ((gdx >= 0.0) || PetscIsInfOrNanReal(gdx)) {
565eb910715SAlp Dener     /* Newton step is not descent or direction produced Inf or NaN
566eb910715SAlp Dener        Update the perturbation for next time */
567eb910715SAlp Dener     if (bnk->pert <= 0.0) {
568eb910715SAlp Dener       /* Initialize the perturbation */
569eb910715SAlp Dener       bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
570eb910715SAlp Dener       if (bnk->is_gltr) {
571eb910715SAlp Dener         ierr = KSPCGGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr);
572eb910715SAlp Dener         bnk->pert = PetscMax(bnk->pert, -e_min);
573eb910715SAlp Dener       }
574eb910715SAlp Dener     } else {
575eb910715SAlp Dener       /* Increase the perturbation */
576eb910715SAlp Dener       bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
577eb910715SAlp Dener     }
578eb910715SAlp Dener 
579eb910715SAlp Dener     if (BNK_PC_BFGS != bnk->pc_type) {
580eb910715SAlp Dener       /* We don't have the bfgs matrix around and updated
581eb910715SAlp Dener          Must use gradient direction in this case */
582080d2917SAlp Dener       ierr = VecCopy(tao->gradient, tao->stepdirection);CHKERRQ(ierr);
583eb910715SAlp Dener       *stepType = BNK_GRADIENT;
584eb910715SAlp Dener     } else {
585eb910715SAlp Dener       /* Attempt to use the BFGS direction */
58609164190SAlp Dener       ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
587eb910715SAlp Dener 
5888d5ead36SAlp Dener       /* Check for success (descent direction)
5898d5ead36SAlp Dener          NOTE: Negative gdx here means not a descent direction because
5908d5ead36SAlp Dener          the fall-back step is missing a negative sign. */
591080d2917SAlp Dener       ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr);
5928d5ead36SAlp Dener       if ((gdx <= 0) || PetscIsInfOrNanReal(gdx)) {
593eb910715SAlp Dener         /* BFGS direction is not descent or direction produced not a number
594eb910715SAlp Dener            We can assert bfgsUpdates > 1 in this case because
595eb910715SAlp Dener            the first solve produces the scaled gradient direction,
596eb910715SAlp Dener            which is guaranteed to be descent */
597eb910715SAlp Dener 
598eb910715SAlp Dener         /* Use steepest descent direction (scaled) */
599eb910715SAlp Dener         if (bnk->f != 0.0) {
600eb910715SAlp Dener           delta = 2.0 * PetscAbsScalar(bnk->f) / (bnk->gnorm*bnk->gnorm);
601eb910715SAlp Dener         } else {
602eb910715SAlp Dener           delta = 2.0 / (bnk->gnorm*bnk->gnorm);
603eb910715SAlp Dener         }
604eb910715SAlp Dener         ierr = MatLMVMSetDelta(bnk->M, delta);CHKERRQ(ierr);
605eb910715SAlp Dener         ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr);
60609164190SAlp Dener         ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
60709164190SAlp Dener         ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
608eb910715SAlp Dener 
609eb910715SAlp Dener         *stepType = BNK_SCALED_GRADIENT;
610eb910715SAlp Dener       } else {
611770b7498SAlp Dener         ierr = MatLMVMGetUpdates(bnk->M, &bfgsUpdates);CHKERRQ(ierr);
612eb910715SAlp Dener         if (1 == bfgsUpdates) {
613eb910715SAlp Dener           /* The first BFGS direction is always the scaled gradient */
614eb910715SAlp Dener           *stepType = BNK_SCALED_GRADIENT;
615eb910715SAlp Dener         } else {
616eb910715SAlp Dener           *stepType = BNK_BFGS;
617eb910715SAlp Dener         }
618eb910715SAlp Dener       }
619eb910715SAlp Dener     }
6208d5ead36SAlp Dener     /* Make sure the safeguarded fall-back step is zero for actively bounded variables */
6218d5ead36SAlp Dener     ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
6222f75a4aaSAlp Dener     ierr = TaoBNKBoundStep(tao, tao->stepdirection);CHKERRQ(ierr);
623eb910715SAlp Dener   } else {
624eb910715SAlp Dener     /* Computed Newton step is descent */
625eb910715SAlp Dener     switch (ksp_reason) {
626eb910715SAlp Dener     case KSP_DIVERGED_NANORINF:
627eb910715SAlp Dener     case KSP_DIVERGED_BREAKDOWN:
628eb910715SAlp Dener     case KSP_DIVERGED_INDEFINITE_MAT:
629eb910715SAlp Dener     case KSP_DIVERGED_INDEFINITE_PC:
630eb910715SAlp Dener     case KSP_CONVERGED_CG_NEG_CURVE:
631eb910715SAlp Dener       /* Matrix or preconditioner is indefinite; increase perturbation */
632eb910715SAlp Dener       if (bnk->pert <= 0.0) {
633eb910715SAlp Dener         /* Initialize the perturbation */
634eb910715SAlp Dener         bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
635eb910715SAlp Dener         if (bnk->is_gltr) {
636eb910715SAlp Dener           ierr = KSPCGGLTRGetMinEig(tao->ksp, &e_min);CHKERRQ(ierr);
637eb910715SAlp Dener           bnk->pert = PetscMax(bnk->pert, -e_min);
638eb910715SAlp Dener         }
639eb910715SAlp Dener       } else {
640eb910715SAlp Dener         /* Increase the perturbation */
641eb910715SAlp Dener         bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
642eb910715SAlp Dener       }
643eb910715SAlp Dener       break;
644eb910715SAlp Dener 
645eb910715SAlp Dener     default:
646eb910715SAlp Dener       /* Newton step computation is good; decrease perturbation */
647eb910715SAlp Dener       bnk->pert = PetscMin(bnk->psfac * bnk->pert, bnk->pmsfac * bnk->gnorm);
648eb910715SAlp Dener       if (bnk->pert < bnk->pmin) {
649eb910715SAlp Dener         bnk->pert = 0.0;
650eb910715SAlp Dener       }
651eb910715SAlp Dener       break;
652eb910715SAlp Dener     }
653fed79b8eSAlp Dener     *stepType = BNK_NEWTON;
654eb910715SAlp Dener   }
655eb910715SAlp Dener   PetscFunctionReturn(0);
656eb910715SAlp Dener }
657eb910715SAlp Dener 
658df278d8fSAlp Dener /*------------------------------------------------------------*/
659df278d8fSAlp Dener 
660df278d8fSAlp Dener /* Routine for performing a bound-projected More-Thuente line search.
661df278d8fSAlp Dener 
662df278d8fSAlp Dener   Includes fallbacks to BFGS, scaled gradient, and unscaled gradient steps if the
663df278d8fSAlp Dener   Newton step does not produce a valid step length.
664df278d8fSAlp Dener */
665df278d8fSAlp Dener 
666c14b763aSAlp Dener PetscErrorCode TaoBNKPerformLineSearch(Tao tao, PetscInt stepType, PetscReal *steplen, TaoLineSearchConvergedReason *reason)
667c14b763aSAlp Dener {
668c14b763aSAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
669c14b763aSAlp Dener   PetscErrorCode ierr;
670c14b763aSAlp Dener   TaoLineSearchConvergedReason ls_reason;
671c14b763aSAlp Dener 
672c14b763aSAlp Dener   PetscReal      e_min, gdx, delta;
673c14b763aSAlp Dener   PetscInt       bfgsUpdates;
674c14b763aSAlp Dener 
675c14b763aSAlp Dener   PetscFunctionBegin;
676c14b763aSAlp Dener   /* Perform the linesearch */
677770b7498SAlp Dener   *steplen = 1.0;
678c14b763aSAlp Dener   ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, steplen, &ls_reason);CHKERRQ(ierr);
679c14b763aSAlp Dener   ierr = VecBoundGradientProjection(bnk->unprojected_gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr);
680c14b763aSAlp Dener   ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
681c14b763aSAlp Dener 
682c14b763aSAlp Dener   while (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER && stepType != BNK_GRADIENT) {
683c14b763aSAlp Dener     /* Linesearch failed, revert solution */
684c14b763aSAlp Dener     bnk->f = bnk->fold;
685c14b763aSAlp Dener     ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr);
686c14b763aSAlp Dener     ierr = VecCopy(bnk->Gold, tao->gradient);CHKERRQ(ierr);
687c14b763aSAlp Dener     ierr = VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient);CHKERRQ(ierr);
688c14b763aSAlp Dener 
689c14b763aSAlp Dener     switch(stepType) {
690c14b763aSAlp Dener     case BNK_NEWTON:
6918d5ead36SAlp Dener       /* Failed to obtain acceptable iterate with Newton step
692c14b763aSAlp Dener          Update the perturbation for next time */
693c14b763aSAlp Dener       if (bnk->pert <= 0.0) {
694c14b763aSAlp Dener         /* Initialize the perturbation */
695c14b763aSAlp Dener         bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
696c14b763aSAlp Dener         if (bnk->is_gltr) {
697c14b763aSAlp Dener           ierr = KSPCGGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr);
698c14b763aSAlp Dener           bnk->pert = PetscMax(bnk->pert, -e_min);
699c14b763aSAlp Dener         }
700c14b763aSAlp Dener       } else {
701c14b763aSAlp Dener         /* Increase the perturbation */
702c14b763aSAlp Dener         bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
703c14b763aSAlp Dener       }
704c14b763aSAlp Dener 
705c14b763aSAlp Dener       if (BNK_PC_BFGS != bnk->pc_type) {
706c14b763aSAlp Dener         /* We don't have the bfgs matrix around and being updated
707c14b763aSAlp Dener            Must use gradient direction in this case */
708c14b763aSAlp Dener         ierr = VecCopy(tao->gradient, tao->stepdirection);CHKERRQ(ierr);
709c14b763aSAlp Dener         stepType = BNK_GRADIENT;
710c14b763aSAlp Dener       } else {
711c14b763aSAlp Dener         /* Attempt to use the BFGS direction */
712c14b763aSAlp Dener         ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
7138d5ead36SAlp Dener         /* Check for success (descent direction)
7148d5ead36SAlp Dener            NOTE: Negative gdx means not a descent direction because the step here is missing a negative sign. */
715c14b763aSAlp Dener         ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr);
716c14b763aSAlp Dener         if ((gdx <= 0) || PetscIsInfOrNanReal(gdx)) {
717c14b763aSAlp Dener           /* BFGS direction is not descent or direction produced not a number
718c14b763aSAlp Dener              We can assert bfgsUpdates > 1 in this case
719c14b763aSAlp Dener              Use steepest descent direction (scaled) */
720c14b763aSAlp Dener           if (bnk->f != 0.0) {
721c14b763aSAlp Dener             delta = 2.0 * PetscAbsScalar(bnk->f) / (bnk->gnorm*bnk->gnorm);
722c14b763aSAlp Dener           } else {
723c14b763aSAlp Dener             delta = 2.0 / (bnk->gnorm*bnk->gnorm);
724c14b763aSAlp Dener           }
725c14b763aSAlp Dener           ierr = MatLMVMSetDelta(bnk->M, delta);CHKERRQ(ierr);
726c14b763aSAlp Dener           ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr);
727c14b763aSAlp Dener           ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
728c14b763aSAlp Dener           ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
729c14b763aSAlp Dener 
730c14b763aSAlp Dener           bfgsUpdates = 1;
731c14b763aSAlp Dener           stepType = BNK_SCALED_GRADIENT;
732c14b763aSAlp Dener         } else {
733c14b763aSAlp Dener           ierr = MatLMVMGetUpdates(bnk->M, &bfgsUpdates);CHKERRQ(ierr);
734c14b763aSAlp Dener           if (1 == bfgsUpdates) {
735c14b763aSAlp Dener             /* The first BFGS direction is always the scaled gradient */
736c14b763aSAlp Dener             stepType = BNK_SCALED_GRADIENT;
737c14b763aSAlp Dener           } else {
738c14b763aSAlp Dener             stepType = BNK_BFGS;
739c14b763aSAlp Dener           }
740c14b763aSAlp Dener         }
741c14b763aSAlp Dener       }
742c14b763aSAlp Dener       break;
743c14b763aSAlp Dener 
744c14b763aSAlp Dener     case BNK_BFGS:
745c14b763aSAlp Dener       /* Can only enter if pc_type == BNK_PC_BFGS
746c14b763aSAlp Dener          Failed to obtain acceptable iterate with BFGS step
747c14b763aSAlp Dener          Attempt to use the scaled gradient direction */
748c14b763aSAlp Dener       if (bnk->f != 0.0) {
749c14b763aSAlp Dener         delta = 2.0 * PetscAbsScalar(bnk->f) / (bnk->gnorm*bnk->gnorm);
750c14b763aSAlp Dener       } else {
751c14b763aSAlp Dener         delta = 2.0 / (bnk->gnorm*bnk->gnorm);
752c14b763aSAlp Dener       }
753c14b763aSAlp Dener       ierr = MatLMVMSetDelta(bnk->M, delta);CHKERRQ(ierr);
754c14b763aSAlp Dener       ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr);
755c14b763aSAlp Dener       ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
756c14b763aSAlp Dener       ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
757c14b763aSAlp Dener 
758c14b763aSAlp Dener       bfgsUpdates = 1;
759c14b763aSAlp Dener       stepType = BNK_SCALED_GRADIENT;
760c14b763aSAlp Dener       break;
761c14b763aSAlp Dener 
762c14b763aSAlp Dener     case BNK_SCALED_GRADIENT:
763c14b763aSAlp Dener       /* Can only enter if pc_type == BNK_PC_BFGS
764c14b763aSAlp Dener          The scaled gradient step did not produce a new iterate;
765c14b763aSAlp Dener          attemp to use the gradient direction.
766c14b763aSAlp Dener          Need to make sure we are not using a different diagonal scaling */
767c14b763aSAlp Dener       ierr = MatLMVMSetScale(bnk->M,0);CHKERRQ(ierr);
768c14b763aSAlp Dener       ierr = MatLMVMSetDelta(bnk->M,1.0);CHKERRQ(ierr);
769c14b763aSAlp Dener       ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr);
770c14b763aSAlp Dener       ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
771c14b763aSAlp Dener       ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
772c14b763aSAlp Dener 
773c14b763aSAlp Dener       bfgsUpdates = 1;
774c14b763aSAlp Dener       stepType = BNK_GRADIENT;
775c14b763aSAlp Dener       break;
776c14b763aSAlp Dener     }
7778d5ead36SAlp Dener     /* Make sure the safeguarded fall-back step is zero for actively bounded variables */
7788d5ead36SAlp Dener     ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
7792f75a4aaSAlp Dener     ierr = TaoBNKBoundStep(tao, tao->stepdirection);CHKERRQ(ierr);
780c14b763aSAlp Dener 
7818d5ead36SAlp Dener     /* Perform one last line search with the fall-back step */
782770b7498SAlp Dener     *steplen = 1.0;
783c14b763aSAlp Dener     ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, steplen, &ls_reason);CHKERRQ(ierr);
784c14b763aSAlp Dener     ierr = VecBoundGradientProjection(bnk->unprojected_gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr);
785c14b763aSAlp Dener     ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
786c14b763aSAlp Dener   }
787c14b763aSAlp Dener   *reason = ls_reason;
788c14b763aSAlp Dener   PetscFunctionReturn(0);
789c14b763aSAlp Dener }
790c14b763aSAlp Dener 
791df278d8fSAlp Dener /*------------------------------------------------------------*/
792df278d8fSAlp Dener 
793df278d8fSAlp Dener /* Routine for updating the trust radius.
794df278d8fSAlp Dener 
795df278d8fSAlp Dener   Function features three different update methods:
796df278d8fSAlp Dener   1) Line-search step length based
797df278d8fSAlp Dener   2) Predicted decrease on the CG quadratic model
798df278d8fSAlp Dener   3) Interpolation
799df278d8fSAlp Dener */
800df278d8fSAlp Dener 
80128017e9fSAlp Dener PetscErrorCode TaoBNKUpdateTrustRadius(Tao tao, PetscReal prered, PetscReal actred, PetscInt updateType, PetscInt stepType, PetscBool *accept)
802080d2917SAlp Dener {
803080d2917SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
804080d2917SAlp Dener   PetscErrorCode ierr;
805080d2917SAlp Dener 
806b1c2d0e3SAlp Dener   PetscReal      step, kappa;
807080d2917SAlp Dener   PetscReal      gdx, tau_1, tau_2, tau_min, tau_max;
808080d2917SAlp Dener 
809080d2917SAlp Dener   PetscFunctionBegin;
810080d2917SAlp Dener   /* Update trust region radius */
811080d2917SAlp Dener   *accept = PETSC_FALSE;
81228017e9fSAlp Dener   switch(updateType) {
813080d2917SAlp Dener   case BNK_UPDATE_STEP:
814c14b763aSAlp Dener     *accept = PETSC_TRUE; /* always accept here because line search succeeded */
815080d2917SAlp Dener     if (stepType == BNK_NEWTON) {
816080d2917SAlp Dener       ierr = TaoLineSearchGetStepLength(tao->linesearch, &step);CHKERRQ(ierr);
817080d2917SAlp Dener       if (step < bnk->nu1) {
818080d2917SAlp Dener         /* Very bad step taken; reduce radius */
819080d2917SAlp Dener         tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust);
820080d2917SAlp Dener       } else if (step < bnk->nu2) {
821080d2917SAlp Dener         /* Reasonably bad step taken; reduce radius */
822080d2917SAlp Dener         tao->trust = bnk->omega2 * PetscMin(bnk->dnorm, tao->trust);
823080d2917SAlp Dener       } else if (step < bnk->nu3) {
824080d2917SAlp Dener         /*  Reasonable step was taken; leave radius alone */
825080d2917SAlp Dener         if (bnk->omega3 < 1.0) {
826080d2917SAlp Dener           tao->trust = bnk->omega3 * PetscMin(bnk->dnorm, tao->trust);
827080d2917SAlp Dener         } else if (bnk->omega3 > 1.0) {
828080d2917SAlp Dener           tao->trust = PetscMax(bnk->omega3 * bnk->dnorm, tao->trust);
829080d2917SAlp Dener         }
830080d2917SAlp Dener       } else if (step < bnk->nu4) {
831080d2917SAlp Dener         /*  Full step taken; increase the radius */
832080d2917SAlp Dener         tao->trust = PetscMax(bnk->omega4 * bnk->dnorm, tao->trust);
833080d2917SAlp Dener       } else {
834080d2917SAlp Dener         /*  More than full step taken; increase the radius */
835080d2917SAlp Dener         tao->trust = PetscMax(bnk->omega5 * bnk->dnorm, tao->trust);
836080d2917SAlp Dener       }
837080d2917SAlp Dener     } else {
838080d2917SAlp Dener       /*  Newton step was not good; reduce the radius */
839080d2917SAlp Dener       tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust);
840080d2917SAlp Dener     }
841080d2917SAlp Dener     break;
842080d2917SAlp Dener 
843080d2917SAlp Dener   case BNK_UPDATE_REDUCTION:
844080d2917SAlp Dener     if (stepType == BNK_NEWTON) {
845b1c2d0e3SAlp Dener       if (prered < 0.0) {
846fed79b8eSAlp Dener         /* The predicted reduction has the wrong sign.  This cannot
847fed79b8eSAlp Dener            happen in infinite precision arithmetic.  Step should
848fed79b8eSAlp Dener            be rejected! */
849080d2917SAlp Dener         tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
850fed79b8eSAlp Dener       }
851fed79b8eSAlp Dener       else {
852b1c2d0e3SAlp Dener         if (PetscIsInfOrNanReal(actred)) {
853080d2917SAlp Dener           tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
854080d2917SAlp Dener         } else {
8550a4511e9SAlp Dener           if ((PetscAbsScalar(actred) <= PetscMax(1.0, PetscAbsScalar(bnk->f))*bnk->epsilon) &&
8560a4511e9SAlp Dener               (PetscAbsScalar(prered) <= PetscMax(1.0, PetscAbsScalar(bnk->f))*bnk->epsilon)) {
857080d2917SAlp Dener             kappa = 1.0;
858fed79b8eSAlp Dener           }
859fed79b8eSAlp Dener           else {
860080d2917SAlp Dener             kappa = actred / prered;
861080d2917SAlp Dener           }
862fed79b8eSAlp Dener 
863fed79b8eSAlp Dener           /* Accept or reject the step and update radius */
864080d2917SAlp Dener           if (kappa < bnk->eta1) {
865fed79b8eSAlp Dener             /* Reject the step */
866080d2917SAlp Dener             tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
867fed79b8eSAlp Dener           }
868fed79b8eSAlp Dener           else {
869fed79b8eSAlp Dener             /* Accept the step */
870c133c014SAlp Dener             *accept = PETSC_TRUE;
871c133c014SAlp Dener             /* Update the trust region radius only if the computed step is at the trust radius boundary */
8728d5ead36SAlp Dener             if (bnk->dnorm == tao->trust) {
873080d2917SAlp Dener               if (kappa < bnk->eta2) {
874080d2917SAlp Dener                 /* Marginal bad step */
875c133c014SAlp Dener                 tao->trust = bnk->alpha2 * tao->trust;
876080d2917SAlp Dener               }
877fed79b8eSAlp Dener               else if (kappa < bnk->eta3) {
878fed79b8eSAlp Dener                 /* Reasonable step */
879fed79b8eSAlp Dener                 tao->trust = bnk->alpha3 * tao->trust;
880fed79b8eSAlp Dener               }
881fed79b8eSAlp Dener               else if (kappa < bnk->eta4) {
882080d2917SAlp Dener                 /* Good step */
883c133c014SAlp Dener                 tao->trust = bnk->alpha4 * tao->trust;
884fed79b8eSAlp Dener               }
885fed79b8eSAlp Dener               else {
886080d2917SAlp Dener                 /* Very good step */
887c133c014SAlp Dener                 tao->trust = bnk->alpha5 * tao->trust;
888080d2917SAlp Dener               }
889c133c014SAlp Dener             }
890080d2917SAlp Dener           }
891080d2917SAlp Dener         }
892080d2917SAlp Dener       }
893080d2917SAlp Dener     } else {
894080d2917SAlp Dener       /*  Newton step was not good; reduce the radius */
895080d2917SAlp Dener       tao->trust = bnk->alpha1 * PetscMin(bnk->dnorm, tao->trust);
896080d2917SAlp Dener     }
897080d2917SAlp Dener     break;
898080d2917SAlp Dener 
899080d2917SAlp Dener   default:
900080d2917SAlp Dener     if (stepType == BNK_NEWTON) {
901b1c2d0e3SAlp Dener       if (prered < 0.0) {
902080d2917SAlp Dener         /*  The predicted reduction has the wrong sign.  This cannot */
903080d2917SAlp Dener         /*  happen in infinite precision arithmetic.  Step should */
904080d2917SAlp Dener         /*  be rejected! */
905080d2917SAlp Dener         tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
906080d2917SAlp Dener       } else {
907b1c2d0e3SAlp Dener         if (PetscIsInfOrNanReal(actred)) {
908080d2917SAlp Dener           tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
909080d2917SAlp Dener         } else {
910080d2917SAlp Dener           if ((PetscAbsScalar(actred) <= bnk->epsilon) && (PetscAbsScalar(prered) <= bnk->epsilon)) {
911080d2917SAlp Dener             kappa = 1.0;
912080d2917SAlp Dener           } else {
913080d2917SAlp Dener             kappa = actred / prered;
914080d2917SAlp Dener           }
915080d2917SAlp Dener 
916080d2917SAlp Dener           ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr);
917080d2917SAlp Dener           tau_1 = bnk->theta * gdx / (bnk->theta * gdx - (1.0 - bnk->theta) * prered + actred);
918080d2917SAlp Dener           tau_2 = bnk->theta * gdx / (bnk->theta * gdx + (1.0 + bnk->theta) * prered - actred);
919080d2917SAlp Dener           tau_min = PetscMin(tau_1, tau_2);
920080d2917SAlp Dener           tau_max = PetscMax(tau_1, tau_2);
921080d2917SAlp Dener 
922080d2917SAlp Dener           if (kappa >= 1.0 - bnk->mu1) {
923080d2917SAlp Dener             /*  Great agreement */
924080d2917SAlp Dener             *accept = PETSC_TRUE;
925080d2917SAlp Dener             if (tau_max < 1.0) {
926080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm);
927080d2917SAlp Dener             } else if (tau_max > bnk->gamma4) {
928080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, bnk->gamma4 * bnk->dnorm);
929080d2917SAlp Dener             } else {
930080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm);
931080d2917SAlp Dener             }
932080d2917SAlp Dener           } else if (kappa >= 1.0 - bnk->mu2) {
933080d2917SAlp Dener             /*  Good agreement */
934080d2917SAlp Dener             *accept = PETSC_TRUE;
935080d2917SAlp Dener             if (tau_max < bnk->gamma2) {
936080d2917SAlp Dener               tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm);
937080d2917SAlp Dener             } else if (tau_max > bnk->gamma3) {
938080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm);
939080d2917SAlp Dener             } else if (tau_max < 1.0) {
940080d2917SAlp Dener               tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm);
941080d2917SAlp Dener             } else {
942080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm);
943080d2917SAlp Dener             }
944080d2917SAlp Dener           } else {
945080d2917SAlp Dener             /*  Not good agreement */
946080d2917SAlp Dener             if (tau_min > 1.0) {
947080d2917SAlp Dener               tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm);
948080d2917SAlp Dener             } else if (tau_max < bnk->gamma1) {
949080d2917SAlp Dener               tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
950080d2917SAlp Dener             } else if ((tau_min < bnk->gamma1) && (tau_max >= 1.0)) {
951080d2917SAlp Dener               tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
952080d2917SAlp Dener             } else if ((tau_1 >= bnk->gamma1) && (tau_1 < 1.0) && ((tau_2 < bnk->gamma1) || (tau_2 >= 1.0))) {
953080d2917SAlp Dener               tao->trust = tau_1 * PetscMin(tao->trust, bnk->dnorm);
954080d2917SAlp Dener             } else if ((tau_2 >= bnk->gamma1) && (tau_2 < 1.0) && ((tau_1 < bnk->gamma1) || (tau_2 >= 1.0))) {
955080d2917SAlp Dener               tao->trust = tau_2 * PetscMin(tao->trust, bnk->dnorm);
956080d2917SAlp Dener             } else {
957080d2917SAlp Dener               tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm);
958080d2917SAlp Dener             }
959080d2917SAlp Dener           }
960080d2917SAlp Dener         }
961080d2917SAlp Dener       }
962080d2917SAlp Dener     } else {
963080d2917SAlp Dener       /*  Newton step was not good; reduce the radius */
964080d2917SAlp Dener       tao->trust = bnk->gamma1 * PetscMin(bnk->dnorm, tao->trust);
965080d2917SAlp Dener     }
96628017e9fSAlp Dener     break;
967080d2917SAlp Dener   }
968c133c014SAlp Dener   /* Make sure the radius does not violate min and max settings */
969c133c014SAlp Dener   tao->trust = PetscMin(tao->trust, bnk->max_radius);
970fed79b8eSAlp Dener   tao->trust = PetscMax(tao->trust, bnk->min_radius);
971080d2917SAlp Dener   PetscFunctionReturn(0);
972080d2917SAlp Dener }
973080d2917SAlp Dener 
974eb910715SAlp Dener /* ---------------------------------------------------------- */
975df278d8fSAlp Dener 
976*62675beeSAlp Dener PetscErrorCode TaoBNKAddStepCounts(Tao tao, PetscInt stepType)
977*62675beeSAlp Dener {
978*62675beeSAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
979*62675beeSAlp Dener 
980*62675beeSAlp Dener   PetscFunctionBegin;
981*62675beeSAlp Dener   switch (stepType) {
982*62675beeSAlp Dener   case BNK_NEWTON:
983*62675beeSAlp Dener     ++bnk->newt;
984*62675beeSAlp Dener     break;
985*62675beeSAlp Dener   case BNK_BFGS:
986*62675beeSAlp Dener     ++bnk->bfgs;
987*62675beeSAlp Dener     break;
988*62675beeSAlp Dener   case BNK_SCALED_GRADIENT:
989*62675beeSAlp Dener     ++bnk->sgrad;
990*62675beeSAlp Dener     break;
991*62675beeSAlp Dener   case BNK_GRADIENT:
992*62675beeSAlp Dener     ++bnk->grad;
993*62675beeSAlp Dener     break;
994*62675beeSAlp Dener   default:
995*62675beeSAlp Dener     break;
996*62675beeSAlp Dener   }
997*62675beeSAlp Dener   PetscFunctionReturn(0);
998*62675beeSAlp Dener }
999*62675beeSAlp Dener 
1000*62675beeSAlp Dener /* ---------------------------------------------------------- */
1001*62675beeSAlp Dener 
1002eb910715SAlp Dener static PetscErrorCode TaoSetUp_BNK(Tao tao)
1003eb910715SAlp Dener {
1004eb910715SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
1005eb910715SAlp Dener   PetscErrorCode ierr;
1006eb910715SAlp Dener 
1007eb910715SAlp Dener   PetscFunctionBegin;
1008eb910715SAlp Dener   if (!tao->gradient) {ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);}
1009eb910715SAlp Dener   if (!tao->stepdirection) {ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr);}
1010eb910715SAlp Dener   if (!bnk->W) {ierr = VecDuplicate(tao->solution,&bnk->W);CHKERRQ(ierr);}
1011eb910715SAlp Dener   if (!bnk->Xold) {ierr = VecDuplicate(tao->solution,&bnk->Xold);CHKERRQ(ierr);}
1012eb910715SAlp Dener   if (!bnk->Gold) {ierr = VecDuplicate(tao->solution,&bnk->Gold);CHKERRQ(ierr);}
101309164190SAlp Dener   if (!bnk->Xwork) {ierr = VecDuplicate(tao->solution,&bnk->Xwork);CHKERRQ(ierr);}
101409164190SAlp Dener   if (!bnk->Gwork) {ierr = VecDuplicate(tao->solution,&bnk->Gwork);CHKERRQ(ierr);}
101509164190SAlp Dener   if (!bnk->unprojected_gradient) {ierr = VecDuplicate(tao->solution,&bnk->unprojected_gradient);CHKERRQ(ierr);}
101609164190SAlp Dener   if (!bnk->unprojected_gradient_old) {ierr = VecDuplicate(tao->solution,&bnk->unprojected_gradient_old);CHKERRQ(ierr);}
1017eb910715SAlp Dener   bnk->Diag = 0;
1018*62675beeSAlp Dener   bnk->Diag_min = 0;
1019*62675beeSAlp Dener   bnk->Diag_max = 0;
1020*62675beeSAlp Dener   bnk->inactive_work = 0;
1021*62675beeSAlp Dener   bnk->active_work = 0;
1022*62675beeSAlp Dener   bnk->inactive_idx = 0;
1023*62675beeSAlp Dener   bnk->active_idx = 0;
1024*62675beeSAlp Dener   bnk->active_lower = 0;
1025*62675beeSAlp Dener   bnk->active_upper = 0;
1026*62675beeSAlp Dener   bnk->active_fixed = 0;
1027eb910715SAlp Dener   bnk->M = 0;
102809164190SAlp Dener   bnk->H_inactive = 0;
102909164190SAlp Dener   bnk->Hpre_inactive = 0;
1030eb910715SAlp Dener   PetscFunctionReturn(0);
1031eb910715SAlp Dener }
1032eb910715SAlp Dener 
1033eb910715SAlp Dener /*------------------------------------------------------------*/
1034df278d8fSAlp Dener 
1035eb910715SAlp Dener static PetscErrorCode TaoDestroy_BNK(Tao tao)
1036eb910715SAlp Dener {
1037eb910715SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
1038eb910715SAlp Dener   PetscErrorCode ierr;
1039eb910715SAlp Dener 
1040eb910715SAlp Dener   PetscFunctionBegin;
1041eb910715SAlp Dener   if (tao->setupcalled) {
1042eb910715SAlp Dener     ierr = VecDestroy(&bnk->W);CHKERRQ(ierr);
1043eb910715SAlp Dener     ierr = VecDestroy(&bnk->Xold);CHKERRQ(ierr);
1044eb910715SAlp Dener     ierr = VecDestroy(&bnk->Gold);CHKERRQ(ierr);
104509164190SAlp Dener     ierr = VecDestroy(&bnk->Xwork);CHKERRQ(ierr);
104609164190SAlp Dener     ierr = VecDestroy(&bnk->Gwork);CHKERRQ(ierr);
104709164190SAlp Dener     ierr = VecDestroy(&bnk->unprojected_gradient);CHKERRQ(ierr);
104809164190SAlp Dener     ierr = VecDestroy(&bnk->unprojected_gradient_old);CHKERRQ(ierr);
1049eb910715SAlp Dener   }
1050eb910715SAlp Dener   ierr = VecDestroy(&bnk->Diag);CHKERRQ(ierr);
1051*62675beeSAlp Dener   ierr = VecDestroy(&bnk->Diag_min);CHKERRQ(ierr);
1052*62675beeSAlp Dener   ierr = VecDestroy(&bnk->Diag_max);CHKERRQ(ierr);
1053eb910715SAlp Dener   ierr = MatDestroy(&bnk->M);CHKERRQ(ierr);
1054*62675beeSAlp Dener   if (bnk->Hpre_inactive != bnk->H_inactive) {ierr = MatDestroy(&bnk->Hpre_inactive);CHKERRQ(ierr);}
1055*62675beeSAlp Dener   ierr = MatDestroy(&bnk->H_inactive);CHKERRQ(ierr);
1056eb910715SAlp Dener   ierr = PetscFree(tao->data);CHKERRQ(ierr);
1057eb910715SAlp Dener   PetscFunctionReturn(0);
1058eb910715SAlp Dener }
1059eb910715SAlp Dener 
1060eb910715SAlp Dener /*------------------------------------------------------------*/
1061df278d8fSAlp Dener 
1062eb910715SAlp Dener static PetscErrorCode TaoSetFromOptions_BNK(PetscOptionItems *PetscOptionsObject,Tao tao)
1063eb910715SAlp Dener {
1064eb910715SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
1065eb910715SAlp Dener   PetscErrorCode ierr;
1066eb910715SAlp Dener 
1067eb910715SAlp Dener   PetscFunctionBegin;
1068eb910715SAlp Dener   ierr = PetscOptionsHead(PetscOptionsObject,"Newton line search method for unconstrained optimization");CHKERRQ(ierr);
10698d5ead36SAlp 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);
10708d5ead36SAlp 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);
10718d5ead36SAlp 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);
10728d5ead36SAlp 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);
10732f75a4aaSAlp 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);
10748d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_sval", "perturbation starting value", "", bnk->sval, &bnk->sval,NULL);CHKERRQ(ierr);
10758d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_imin", "minimum initial perturbation", "", bnk->imin, &bnk->imin,NULL);CHKERRQ(ierr);
10768d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_imax", "maximum initial perturbation", "", bnk->imax, &bnk->imax,NULL);CHKERRQ(ierr);
10778d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_imfac", "initial merit factor", "", bnk->imfac, &bnk->imfac,NULL);CHKERRQ(ierr);
10788d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_pmin", "minimum perturbation", "", bnk->pmin, &bnk->pmin,NULL);CHKERRQ(ierr);
10798d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_pmax", "maximum perturbation", "", bnk->pmax, &bnk->pmax,NULL);CHKERRQ(ierr);
10808d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_pgfac", "growth factor", "", bnk->pgfac, &bnk->pgfac,NULL);CHKERRQ(ierr);
10818d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_psfac", "shrink factor", "", bnk->psfac, &bnk->psfac,NULL);CHKERRQ(ierr);
10828d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_pmgfac", "merit growth factor", "", bnk->pmgfac, &bnk->pmgfac,NULL);CHKERRQ(ierr);
10838d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_pmsfac", "merit shrink factor", "", bnk->pmsfac, &bnk->pmsfac,NULL);CHKERRQ(ierr);
10848d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_eta1", "poor steplength; reduce radius", "", bnk->eta1, &bnk->eta1,NULL);CHKERRQ(ierr);
10858d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_eta2", "reasonable steplength; leave radius alone", "", bnk->eta2, &bnk->eta2,NULL);CHKERRQ(ierr);
10868d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_eta3", "good steplength; increase radius", "", bnk->eta3, &bnk->eta3,NULL);CHKERRQ(ierr);
10878d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_eta4", "excellent steplength; greatly increase radius", "", bnk->eta4, &bnk->eta4,NULL);CHKERRQ(ierr);
10888d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_alpha1", "", "", bnk->alpha1, &bnk->alpha1,NULL);CHKERRQ(ierr);
10898d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_alpha2", "", "", bnk->alpha2, &bnk->alpha2,NULL);CHKERRQ(ierr);
10908d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_alpha3", "", "", bnk->alpha3, &bnk->alpha3,NULL);CHKERRQ(ierr);
10918d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_alpha4", "", "", bnk->alpha4, &bnk->alpha4,NULL);CHKERRQ(ierr);
10928d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_alpha5", "", "", bnk->alpha5, &bnk->alpha5,NULL);CHKERRQ(ierr);
10938d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_nu1", "poor steplength; reduce radius", "", bnk->nu1, &bnk->nu1,NULL);CHKERRQ(ierr);
10948d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_nu2", "reasonable steplength; leave radius alone", "", bnk->nu2, &bnk->nu2,NULL);CHKERRQ(ierr);
10958d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_nu3", "good steplength; increase radius", "", bnk->nu3, &bnk->nu3,NULL);CHKERRQ(ierr);
10968d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_nu4", "excellent steplength; greatly increase radius", "", bnk->nu4, &bnk->nu4,NULL);CHKERRQ(ierr);
10978d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_omega1", "", "", bnk->omega1, &bnk->omega1,NULL);CHKERRQ(ierr);
10988d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_omega2", "", "", bnk->omega2, &bnk->omega2,NULL);CHKERRQ(ierr);
10998d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_omega3", "", "", bnk->omega3, &bnk->omega3,NULL);CHKERRQ(ierr);
11008d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_omega4", "", "", bnk->omega4, &bnk->omega4,NULL);CHKERRQ(ierr);
11018d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_omega5", "", "", bnk->omega5, &bnk->omega5,NULL);CHKERRQ(ierr);
11028d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_mu1_i", "", "", bnk->mu1_i, &bnk->mu1_i,NULL);CHKERRQ(ierr);
11038d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_mu2_i", "", "", bnk->mu2_i, &bnk->mu2_i,NULL);CHKERRQ(ierr);
11048d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_gamma1_i", "", "", bnk->gamma1_i, &bnk->gamma1_i,NULL);CHKERRQ(ierr);
11058d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_gamma2_i", "", "", bnk->gamma2_i, &bnk->gamma2_i,NULL);CHKERRQ(ierr);
11068d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_gamma3_i", "", "", bnk->gamma3_i, &bnk->gamma3_i,NULL);CHKERRQ(ierr);
11078d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_gamma4_i", "", "", bnk->gamma4_i, &bnk->gamma4_i,NULL);CHKERRQ(ierr);
11088d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_theta_i", "", "", bnk->theta_i, &bnk->theta_i,NULL);CHKERRQ(ierr);
11098d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_mu1", "", "", bnk->mu1, &bnk->mu1,NULL);CHKERRQ(ierr);
11108d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_mu2", "", "", bnk->mu2, &bnk->mu2,NULL);CHKERRQ(ierr);
11118d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_gamma1", "", "", bnk->gamma1, &bnk->gamma1,NULL);CHKERRQ(ierr);
11128d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_gamma2", "", "", bnk->gamma2, &bnk->gamma2,NULL);CHKERRQ(ierr);
11138d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_gamma3", "", "", bnk->gamma3, &bnk->gamma3,NULL);CHKERRQ(ierr);
11148d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_gamma4", "", "", bnk->gamma4, &bnk->gamma4,NULL);CHKERRQ(ierr);
11158d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_theta", "", "", bnk->theta, &bnk->theta,NULL);CHKERRQ(ierr);
11168d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_min_radius", "lower bound on initial radius", "", bnk->min_radius, &bnk->min_radius,NULL);CHKERRQ(ierr);
11178d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_max_radius", "upper bound on radius", "", bnk->max_radius, &bnk->max_radius,NULL);CHKERRQ(ierr);
11188d5ead36SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_epsilon", "tolerance used when computing actual and predicted reduction", "", bnk->epsilon, &bnk->epsilon,NULL);CHKERRQ(ierr);
11190a4511e9SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_as_tol", "initial tolerance used when estimating actively bounded variables", "", bnk->as_tol, &bnk->as_tol,NULL);CHKERRQ(ierr);
11200a4511e9SAlp Dener   ierr = PetscOptionsReal("-tao_bnk_as_step", "step length used when estimating actively bounded variables", "", bnk->as_step, &bnk->as_step,NULL);CHKERRQ(ierr);
1121eb910715SAlp Dener   ierr = PetscOptionsTail();CHKERRQ(ierr);
1122eb910715SAlp Dener   ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr);
1123eb910715SAlp Dener   ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
1124eb910715SAlp Dener   PetscFunctionReturn(0);
1125eb910715SAlp Dener }
1126eb910715SAlp Dener 
1127eb910715SAlp Dener /*------------------------------------------------------------*/
1128df278d8fSAlp Dener 
1129eb910715SAlp Dener static PetscErrorCode TaoView_BNK(Tao tao, PetscViewer viewer)
1130eb910715SAlp Dener {
1131eb910715SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
1132eb910715SAlp Dener   PetscInt       nrejects;
1133eb910715SAlp Dener   PetscBool      isascii;
1134eb910715SAlp Dener   PetscErrorCode ierr;
1135eb910715SAlp Dener 
1136eb910715SAlp Dener   PetscFunctionBegin;
1137eb910715SAlp Dener   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
1138eb910715SAlp Dener   if (isascii) {
1139eb910715SAlp Dener     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1140eb910715SAlp Dener     if (BNK_PC_BFGS == bnk->pc_type && bnk->M) {
1141eb910715SAlp Dener       ierr = MatLMVMGetRejects(bnk->M,&nrejects);CHKERRQ(ierr);
1142eb910715SAlp Dener       ierr = PetscViewerASCIIPrintf(viewer, "Rejected matrix updates: %D\n",nrejects);CHKERRQ(ierr);
1143eb910715SAlp Dener     }
1144eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "Newton steps: %D\n", bnk->newt);CHKERRQ(ierr);
1145eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "BFGS steps: %D\n", bnk->bfgs);CHKERRQ(ierr);
1146eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "Scaled gradient steps: %D\n", bnk->sgrad);CHKERRQ(ierr);
1147eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", bnk->grad);CHKERRQ(ierr);
1148eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "KSP termination reasons:\n");CHKERRQ(ierr);
1149eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  atol: %D\n", bnk->ksp_atol);CHKERRQ(ierr);
1150eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  rtol: %D\n", bnk->ksp_rtol);CHKERRQ(ierr);
1151eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  ctol: %D\n", bnk->ksp_ctol);CHKERRQ(ierr);
1152eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  negc: %D\n", bnk->ksp_negc);CHKERRQ(ierr);
1153eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  dtol: %D\n", bnk->ksp_dtol);CHKERRQ(ierr);
1154eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  iter: %D\n", bnk->ksp_iter);CHKERRQ(ierr);
1155eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  othr: %D\n", bnk->ksp_othr);CHKERRQ(ierr);
1156eb910715SAlp Dener     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1157eb910715SAlp Dener   }
1158eb910715SAlp Dener   PetscFunctionReturn(0);
1159eb910715SAlp Dener }
1160eb910715SAlp Dener 
1161eb910715SAlp Dener /* ---------------------------------------------------------- */
1162df278d8fSAlp Dener 
1163eb910715SAlp Dener /*MC
1164eb910715SAlp Dener   TAOBNK - Shared base-type for Bounded Newton-Krylov type algorithms.
116566ed3702SAlp Dener   At each iteration, the BNK methods solve the symmetric
1166eb910715SAlp Dener   system of equations to obtain the step diretion dk:
1167eb910715SAlp Dener               Hk dk = -gk
11682b97c8d8SAlp Dener   for free variables only. The step can be globalized either through
11692b97c8d8SAlp Dener   trust-region methods, or a line search, or a heuristic mixture of both.
1170eb910715SAlp Dener 
1171eb910715SAlp Dener     Options Database Keys:
11728d5ead36SAlp Dener + -tao_bnk_pc_type - "none","ahess","bfgs","petsc"
11738d5ead36SAlp Dener . -tao_bnk_bfgs_scale_type - "ahess","phess","bfgs"
11748d5ead36SAlp Dener . -tao_bnk_init_type - "constant","direction","interpolation"
11758d5ead36SAlp Dener . -tao_bnk_update_type - "step","direction","interpolation"
11762f75a4aaSAlp Dener . -tao_bnk_as_type - "none","bertsekas"
11778d5ead36SAlp Dener . -tao_bnk_sval - perturbation starting value
11788d5ead36SAlp Dener . -tao_bnk_imin - minimum initial perturbation
11798d5ead36SAlp Dener . -tao_bnk_imax - maximum initial perturbation
11808d5ead36SAlp Dener . -tao_bnk_pmin - minimum perturbation
11818d5ead36SAlp Dener . -tao_bnk_pmax - maximum perturbation
11828d5ead36SAlp Dener . -tao_bnk_pgfac - growth factor
11838d5ead36SAlp Dener . -tao_bnk_psfac - shrink factor
11848d5ead36SAlp Dener . -tao_bnk_imfac - initial merit factor
11858d5ead36SAlp Dener . -tao_bnk_pmgfac - merit growth factor
11868d5ead36SAlp Dener . -tao_bnk_pmsfac - merit shrink factor
11878d5ead36SAlp Dener . -tao_bnk_eta1 - poor steplength; reduce radius
11888d5ead36SAlp Dener . -tao_bnk_eta2 - reasonable steplength; leave radius
11898d5ead36SAlp Dener . -tao_bnk_eta3 - good steplength; increase readius
11908d5ead36SAlp Dener . -tao_bnk_eta4 - excellent steplength; greatly increase radius
11918d5ead36SAlp Dener . -tao_bnk_alpha1 - alpha1 reduction
11928d5ead36SAlp Dener . -tao_bnk_alpha2 - alpha2 reduction
11938d5ead36SAlp Dener . -tao_bnk_alpha3 - alpha3 reduction
11948d5ead36SAlp Dener . -tao_bnk_alpha4 - alpha4 reduction
11958d5ead36SAlp Dener . -tao_bnk_alpha - alpha5 reduction
11968d5ead36SAlp Dener . -tao_bnk_mu1 - mu1 interpolation update
11978d5ead36SAlp Dener . -tao_bnk_mu2 - mu2 interpolation update
11988d5ead36SAlp Dener . -tao_bnk_gamma1 - gamma1 interpolation update
11998d5ead36SAlp Dener . -tao_bnk_gamma2 - gamma2 interpolation update
12008d5ead36SAlp Dener . -tao_bnk_gamma3 - gamma3 interpolation update
12018d5ead36SAlp Dener . -tao_bnk_gamma4 - gamma4 interpolation update
12028d5ead36SAlp Dener . -tao_bnk_theta - theta interpolation update
12038d5ead36SAlp Dener . -tao_bnk_omega1 - omega1 step update
12048d5ead36SAlp Dener . -tao_bnk_omega2 - omega2 step update
12058d5ead36SAlp Dener . -tao_bnk_omega3 - omega3 step update
12068d5ead36SAlp Dener . -tao_bnk_omega4 - omega4 step update
12078d5ead36SAlp Dener . -tao_bnk_omega5 - omega5 step update
12088d5ead36SAlp Dener . -tao_bnk_mu1_i -  mu1 interpolation init factor
12098d5ead36SAlp Dener . -tao_bnk_mu2_i -  mu2 interpolation init factor
12108d5ead36SAlp Dener . -tao_bnk_gamma1_i -  gamma1 interpolation init factor
12118d5ead36SAlp Dener . -tao_bnk_gamma2_i -  gamma2 interpolation init factor
12128d5ead36SAlp Dener . -tao_bnk_gamma3_i -  gamma3 interpolation init factor
12138d5ead36SAlp Dener . -tao_bnk_gamma4_i -  gamma4 interpolation init factor
12142f75a4aaSAlp Dener . -tao_bnk_theta_i -  theta interpolation init factor
12152f75a4aaSAlp Dener - -tao_bnk_bound_tol -  initial tolerance used in estimating bounded active variables
1216eb910715SAlp Dener 
1217eb910715SAlp Dener   Level: beginner
1218eb910715SAlp Dener M*/
1219eb910715SAlp Dener 
1220eb910715SAlp Dener PetscErrorCode TaoCreate_BNK(Tao tao)
1221eb910715SAlp Dener {
1222eb910715SAlp Dener   TAO_BNK        *bnk;
1223eb910715SAlp Dener   const char     *morethuente_type = TAOLINESEARCHMT;
1224eb910715SAlp Dener   PetscErrorCode ierr;
1225eb910715SAlp Dener 
1226eb910715SAlp Dener   PetscFunctionBegin;
1227eb910715SAlp Dener   ierr = PetscNewLog(tao,&bnk);CHKERRQ(ierr);
1228eb910715SAlp Dener 
1229eb910715SAlp Dener   tao->ops->setup = TaoSetUp_BNK;
1230eb910715SAlp Dener   tao->ops->view = TaoView_BNK;
1231eb910715SAlp Dener   tao->ops->setfromoptions = TaoSetFromOptions_BNK;
1232eb910715SAlp Dener   tao->ops->destroy = TaoDestroy_BNK;
1233eb910715SAlp Dener 
1234eb910715SAlp Dener   /*  Override default settings (unless already changed) */
1235eb910715SAlp Dener   if (!tao->max_it_changed) tao->max_it = 50;
1236eb910715SAlp Dener   if (!tao->trust0_changed) tao->trust0 = 100.0;
1237eb910715SAlp Dener 
1238eb910715SAlp Dener   tao->data = (void*)bnk;
1239eb910715SAlp Dener 
124066ed3702SAlp Dener   /*  Hessian shifting parameters */
1241eb910715SAlp Dener   bnk->sval   = 0.0;
1242eb910715SAlp Dener   bnk->imin   = 1.0e-4;
1243eb910715SAlp Dener   bnk->imax   = 1.0e+2;
1244eb910715SAlp Dener   bnk->imfac  = 1.0e-1;
1245eb910715SAlp Dener 
1246eb910715SAlp Dener   bnk->pmin   = 1.0e-12;
1247eb910715SAlp Dener   bnk->pmax   = 1.0e+2;
1248eb910715SAlp Dener   bnk->pgfac  = 1.0e+1;
1249eb910715SAlp Dener   bnk->psfac  = 4.0e-1;
1250eb910715SAlp Dener   bnk->pmgfac = 1.0e-1;
1251eb910715SAlp Dener   bnk->pmsfac = 1.0e-1;
1252eb910715SAlp Dener 
1253eb910715SAlp Dener   /*  Default values for trust-region radius update based on steplength */
1254eb910715SAlp Dener   bnk->nu1 = 0.25;
1255eb910715SAlp Dener   bnk->nu2 = 0.50;
1256eb910715SAlp Dener   bnk->nu3 = 1.00;
1257eb910715SAlp Dener   bnk->nu4 = 1.25;
1258eb910715SAlp Dener 
1259eb910715SAlp Dener   bnk->omega1 = 0.25;
1260eb910715SAlp Dener   bnk->omega2 = 0.50;
1261eb910715SAlp Dener   bnk->omega3 = 1.00;
1262eb910715SAlp Dener   bnk->omega4 = 2.00;
1263eb910715SAlp Dener   bnk->omega5 = 4.00;
1264eb910715SAlp Dener 
1265eb910715SAlp Dener   /*  Default values for trust-region radius update based on reduction */
1266eb910715SAlp Dener   bnk->eta1 = 1.0e-4;
1267eb910715SAlp Dener   bnk->eta2 = 0.25;
1268eb910715SAlp Dener   bnk->eta3 = 0.50;
1269eb910715SAlp Dener   bnk->eta4 = 0.90;
1270eb910715SAlp Dener 
1271eb910715SAlp Dener   bnk->alpha1 = 0.25;
1272eb910715SAlp Dener   bnk->alpha2 = 0.50;
1273eb910715SAlp Dener   bnk->alpha3 = 1.00;
1274eb910715SAlp Dener   bnk->alpha4 = 2.00;
1275eb910715SAlp Dener   bnk->alpha5 = 4.00;
1276eb910715SAlp Dener 
1277eb910715SAlp Dener   /*  Default values for trust-region radius update based on interpolation */
1278eb910715SAlp Dener   bnk->mu1 = 0.10;
1279eb910715SAlp Dener   bnk->mu2 = 0.50;
1280eb910715SAlp Dener 
1281eb910715SAlp Dener   bnk->gamma1 = 0.25;
1282eb910715SAlp Dener   bnk->gamma2 = 0.50;
1283eb910715SAlp Dener   bnk->gamma3 = 2.00;
1284eb910715SAlp Dener   bnk->gamma4 = 4.00;
1285eb910715SAlp Dener 
1286eb910715SAlp Dener   bnk->theta = 0.05;
1287eb910715SAlp Dener 
1288eb910715SAlp Dener   /*  Default values for trust region initialization based on interpolation */
1289eb910715SAlp Dener   bnk->mu1_i = 0.35;
1290eb910715SAlp Dener   bnk->mu2_i = 0.50;
1291eb910715SAlp Dener 
1292eb910715SAlp Dener   bnk->gamma1_i = 0.0625;
1293eb910715SAlp Dener   bnk->gamma2_i = 0.5;
1294eb910715SAlp Dener   bnk->gamma3_i = 2.0;
1295eb910715SAlp Dener   bnk->gamma4_i = 5.0;
1296eb910715SAlp Dener 
1297eb910715SAlp Dener   bnk->theta_i = 0.25;
1298eb910715SAlp Dener 
1299eb910715SAlp Dener   /*  Remaining parameters */
1300eb910715SAlp Dener   bnk->min_radius = 1.0e-10;
1301eb910715SAlp Dener   bnk->max_radius = 1.0e10;
1302770b7498SAlp Dener   bnk->epsilon = PetscPowReal(PETSC_MACHINE_EPSILON, 2.0/3.0);
13030a4511e9SAlp Dener   bnk->as_tol = 1.0e-3;
13040a4511e9SAlp Dener   bnk->as_step = 1.0e-3;
1305*62675beeSAlp Dener   bnk->dmin = 1.0e-6;
1306*62675beeSAlp Dener   bnk->dmax = 1.0e6;
1307eb910715SAlp Dener 
1308eb910715SAlp Dener   bnk->pc_type         = BNK_PC_BFGS;
1309eb910715SAlp Dener   bnk->bfgs_scale_type = BFGS_SCALE_PHESS;
1310eb910715SAlp Dener   bnk->init_type       = BNK_INIT_INTERPOLATION;
1311fed79b8eSAlp Dener   bnk->update_type     = BNK_UPDATE_INTERPOLATION;
13122f75a4aaSAlp Dener   bnk->as_type         = BNK_AS_BERTSEKAS;
1313eb910715SAlp Dener 
1314eb910715SAlp Dener   ierr = TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch);CHKERRQ(ierr);
1315eb910715SAlp Dener   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr);
1316eb910715SAlp Dener   ierr = TaoLineSearchSetType(tao->linesearch,morethuente_type);CHKERRQ(ierr);
1317eb910715SAlp Dener   ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr);
1318eb910715SAlp Dener   ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr);
1319eb910715SAlp Dener 
1320eb910715SAlp Dener   /*  Set linear solver to default for symmetric matrices */
1321eb910715SAlp Dener   ierr = KSPCreate(((PetscObject)tao)->comm,&tao->ksp);CHKERRQ(ierr);
1322eb910715SAlp Dener   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1);CHKERRQ(ierr);
1323eb910715SAlp Dener   ierr = KSPSetOptionsPrefix(tao->ksp,tao->hdr.prefix);CHKERRQ(ierr);
1324eb910715SAlp Dener   ierr = KSPSetType(tao->ksp,KSPCGSTCG);CHKERRQ(ierr);
1325eb910715SAlp Dener   PetscFunctionReturn(0);
1326eb910715SAlp Dener }
1327