xref: /petsc/src/tao/bound/impls/bnk/bnk.c (revision df278d8f2b91157cf6d5980446bdf033f772efe7)
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 
22*df278d8fSAlp Dener /*------------------------------------------------------------*/
23*df278d8fSAlp Dener 
24*df278d8fSAlp Dener /* Routine for initializing the KSP solver, the BFGS preconditioner, and the initial trust radius estimation */
25*df278d8fSAlp Dener 
26eb910715SAlp Dener PetscErrorCode TaoBNKInitialize(Tao tao)
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 
33eb910715SAlp Dener   PetscReal                    fmin, ftrial, prered, actred, kappa, sigma;
34eb910715SAlp Dener   PetscReal                    tau, tau_1, tau_2, tau_max, tau_min, max_radius;
35eb910715SAlp Dener   PetscReal                    delta, step = 1.0;
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;
44eb910715SAlp Dener   /* Number of times ksp stopped because of these reasons */
45eb910715SAlp Dener   bnk->ksp_atol = 0;
46eb910715SAlp Dener   bnk->ksp_rtol = 0;
47eb910715SAlp Dener   bnk->ksp_dtol = 0;
48eb910715SAlp Dener   bnk->ksp_ctol = 0;
49eb910715SAlp Dener   bnk->ksp_negc = 0;
50eb910715SAlp Dener   bnk->ksp_iter = 0;
51eb910715SAlp Dener   bnk->ksp_othr = 0;
52eb910715SAlp Dener 
53fed79b8eSAlp Dener   /* Initialize the Hessian perturbation */
54fed79b8eSAlp Dener   bnk->pert = bnk->sval;
55fed79b8eSAlp Dener 
56eb910715SAlp Dener   /* Initialize trust-region radius when using nash, stcg, or gltr
57eb910715SAlp Dener      Command automatically ignored for other methods
58eb910715SAlp Dener      Will be reset during the first iteration
59eb910715SAlp Dener   */
60eb910715SAlp Dener   ierr = KSPGetType(tao->ksp,&ksp_type);CHKERRQ(ierr);
61eb910715SAlp Dener   ierr = PetscStrcmp(ksp_type,KSPCGNASH,&bnk->is_nash);CHKERRQ(ierr);
62eb910715SAlp Dener   ierr = PetscStrcmp(ksp_type,KSPCGSTCG,&bnk->is_stcg);CHKERRQ(ierr);
63eb910715SAlp Dener   ierr = PetscStrcmp(ksp_type,KSPCGGLTR,&bnk->is_gltr);CHKERRQ(ierr);
64eb910715SAlp Dener 
65eb910715SAlp Dener   ierr = KSPCGSetRadius(tao->ksp,bnk->max_radius);CHKERRQ(ierr);
66eb910715SAlp Dener 
67eb910715SAlp Dener   if (bnk->is_nash || bnk->is_stcg || bnk->is_gltr) {
68eb910715SAlp Dener     if (tao->trust0 < 0.0) SETERRQ(PETSC_COMM_SELF,1,"Initial radius negative");
69eb910715SAlp Dener     tao->trust = tao->trust0;
70eb910715SAlp Dener     tao->trust = PetscMax(tao->trust, bnk->min_radius);
71eb910715SAlp Dener     tao->trust = PetscMin(tao->trust, bnk->max_radius);
72eb910715SAlp Dener   }
73eb910715SAlp Dener 
74eb910715SAlp Dener   /* Get vectors we will need */
75eb910715SAlp Dener   if (BNK_PC_BFGS == bnk->pc_type && !bnk->M) {
76eb910715SAlp Dener     ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr);
77eb910715SAlp Dener     ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr);
78eb910715SAlp Dener     ierr = MatCreateLMVM(((PetscObject)tao)->comm,n,N,&bnk->M);CHKERRQ(ierr);
79eb910715SAlp Dener     ierr = MatLMVMAllocateVectors(bnk->M,tao->solution);CHKERRQ(ierr);
80eb910715SAlp Dener   }
81eb910715SAlp Dener 
82eb910715SAlp Dener   /* create vectors for the limited memory preconditioner */
83eb910715SAlp Dener   if ((BNK_PC_BFGS == bnk->pc_type) && (BFGS_SCALE_BFGS != bnk->bfgs_scale_type)) {
84eb910715SAlp Dener     if (!bnk->Diag) {
85eb910715SAlp Dener       ierr = VecDuplicate(tao->solution,&bnk->Diag);CHKERRQ(ierr);
86eb910715SAlp Dener     }
87eb910715SAlp Dener   }
88eb910715SAlp Dener 
89eb910715SAlp Dener   /* Modify the preconditioner to use the bfgs approximation */
90eb910715SAlp Dener   ierr = KSPGetPC(tao->ksp, &pc);CHKERRQ(ierr);
91eb910715SAlp Dener   switch(bnk->pc_type) {
92eb910715SAlp Dener   case BNK_PC_NONE:
93eb910715SAlp Dener     ierr = PCSetType(pc, PCNONE);CHKERRQ(ierr);
94eb910715SAlp Dener     ierr = PCSetFromOptions(pc);CHKERRQ(ierr);
95eb910715SAlp Dener     break;
96eb910715SAlp Dener 
97eb910715SAlp Dener   case BNK_PC_AHESS:
98eb910715SAlp Dener     ierr = PCSetType(pc, PCJACOBI);CHKERRQ(ierr);
99eb910715SAlp Dener     ierr = PCSetFromOptions(pc);CHKERRQ(ierr);
100eb910715SAlp Dener     ierr = PCJacobiSetUseAbs(pc,PETSC_TRUE);CHKERRQ(ierr);
101eb910715SAlp Dener     break;
102eb910715SAlp Dener 
103eb910715SAlp Dener   case BNK_PC_BFGS:
104eb910715SAlp Dener     ierr = PCSetType(pc, PCSHELL);CHKERRQ(ierr);
105eb910715SAlp Dener     ierr = PCSetFromOptions(pc);CHKERRQ(ierr);
106eb910715SAlp Dener     ierr = PCShellSetName(pc, "bfgs");CHKERRQ(ierr);
107eb910715SAlp Dener     ierr = PCShellSetContext(pc, bnk->M);CHKERRQ(ierr);
108eb910715SAlp Dener     ierr = PCShellSetApply(pc, MatLMVMSolveShell);CHKERRQ(ierr);
109eb910715SAlp Dener     break;
110eb910715SAlp Dener 
111eb910715SAlp Dener   default:
112eb910715SAlp Dener     /* Use the pc method set by pc_type */
113eb910715SAlp Dener     break;
114eb910715SAlp Dener   }
115eb910715SAlp Dener 
116eb910715SAlp Dener   /* Initialize trust-region radius.  The initialization is only performed
117eb910715SAlp Dener      when we are using Nash, Steihaug-Toint or the Generalized Lanczos method. */
118eb910715SAlp Dener   if (bnk->is_nash || bnk->is_stcg || bnk->is_gltr) {
119eb910715SAlp Dener     switch(bnk->init_type) {
120eb910715SAlp Dener     case BNK_INIT_CONSTANT:
121eb910715SAlp Dener       /* Use the initial radius specified */
122eb910715SAlp Dener       break;
123eb910715SAlp Dener 
124eb910715SAlp Dener     case BNK_INIT_INTERPOLATION:
125eb910715SAlp Dener       /* Use the initial radius specified */
126eb910715SAlp Dener       max_radius = 0.0;
127eb910715SAlp Dener 
128eb910715SAlp Dener       for (j = 0; j < j_max; ++j) {
129eb910715SAlp Dener         fmin = bnk->f;
130eb910715SAlp Dener         sigma = 0.0;
131eb910715SAlp Dener 
132eb910715SAlp Dener         if (needH) {
133eb910715SAlp Dener           ierr  = TaoComputeHessian(tao, tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr);
134eb910715SAlp Dener           needH = 0;
135eb910715SAlp Dener         }
136eb910715SAlp Dener 
137eb910715SAlp Dener         for (i = 0; i < i_max; ++i) {
138eb910715SAlp Dener           ierr = VecCopy(tao->solution,bnk->W);CHKERRQ(ierr);
139eb910715SAlp Dener           ierr = VecAXPY(bnk->W,-tao->trust/bnk->gnorm,tao->gradient);CHKERRQ(ierr);
140eb910715SAlp Dener           ierr = TaoComputeObjective(tao, bnk->W, &ftrial);CHKERRQ(ierr);
141eb910715SAlp Dener           if (PetscIsInfOrNanReal(ftrial)) {
142eb910715SAlp Dener             tau = bnk->gamma1_i;
143eb910715SAlp Dener           } else {
144eb910715SAlp Dener             if (ftrial < fmin) {
145eb910715SAlp Dener               fmin = ftrial;
146eb910715SAlp Dener               sigma = -tao->trust / bnk->gnorm;
147eb910715SAlp Dener             }
148eb910715SAlp Dener 
149080d2917SAlp Dener             ierr = MatMult(tao->hessian, tao->gradient, tao->stepdirection);CHKERRQ(ierr);
150080d2917SAlp Dener             ierr = VecDot(tao->gradient, tao->stepdirection, &prered);CHKERRQ(ierr);
151eb910715SAlp Dener 
152eb910715SAlp Dener             prered = tao->trust * (bnk->gnorm - 0.5 * tao->trust * prered / (bnk->gnorm * bnk->gnorm));
153eb910715SAlp Dener             actred = bnk->f - ftrial;
154eb910715SAlp Dener             if ((PetscAbsScalar(actred) <= bnk->epsilon) && (PetscAbsScalar(prered) <= bnk->epsilon)) {
155eb910715SAlp Dener               kappa = 1.0;
156eb910715SAlp Dener             } else {
157eb910715SAlp Dener               kappa = actred / prered;
158eb910715SAlp Dener             }
159eb910715SAlp Dener 
160eb910715SAlp Dener             tau_1 = bnk->theta_i * bnk->gnorm * tao->trust / (bnk->theta_i * bnk->gnorm * tao->trust + (1.0 - bnk->theta_i) * prered - actred);
161eb910715SAlp Dener             tau_2 = bnk->theta_i * bnk->gnorm * tao->trust / (bnk->theta_i * bnk->gnorm * tao->trust - (1.0 + bnk->theta_i) * prered + actred);
162eb910715SAlp Dener             tau_min = PetscMin(tau_1, tau_2);
163eb910715SAlp Dener             tau_max = PetscMax(tau_1, tau_2);
164eb910715SAlp Dener 
165eb910715SAlp Dener             if (PetscAbsScalar(kappa - 1.0) <= bnk->mu1_i) {
166eb910715SAlp Dener               /* Great agreement */
167eb910715SAlp Dener               max_radius = PetscMax(max_radius, tao->trust);
168eb910715SAlp Dener 
169eb910715SAlp Dener               if (tau_max < 1.0) {
170eb910715SAlp Dener                 tau = bnk->gamma3_i;
171eb910715SAlp Dener               } else if (tau_max > bnk->gamma4_i) {
172eb910715SAlp Dener                 tau = bnk->gamma4_i;
173eb910715SAlp Dener               } else if (tau_1 >= 1.0 && tau_1 <= bnk->gamma4_i && tau_2 < 1.0) {
174eb910715SAlp Dener                 tau = tau_1;
175eb910715SAlp Dener               } else if (tau_2 >= 1.0 && tau_2 <= bnk->gamma4_i && tau_1 < 1.0) {
176eb910715SAlp Dener                 tau = tau_2;
177eb910715SAlp Dener               } else {
178eb910715SAlp Dener                 tau = tau_max;
179eb910715SAlp Dener               }
180eb910715SAlp Dener             } else if (PetscAbsScalar(kappa - 1.0) <= bnk->mu2_i) {
181eb910715SAlp Dener               /* Good agreement */
182eb910715SAlp Dener               max_radius = PetscMax(max_radius, tao->trust);
183eb910715SAlp Dener 
184eb910715SAlp Dener               if (tau_max < bnk->gamma2_i) {
185eb910715SAlp Dener                 tau = bnk->gamma2_i;
186eb910715SAlp Dener               } else if (tau_max > bnk->gamma3_i) {
187eb910715SAlp Dener                 tau = bnk->gamma3_i;
188eb910715SAlp Dener               } else {
189eb910715SAlp Dener                 tau = tau_max;
190eb910715SAlp Dener               }
191eb910715SAlp Dener             } else {
192eb910715SAlp Dener               /* Not good agreement */
193eb910715SAlp Dener               if (tau_min > 1.0) {
194eb910715SAlp Dener                 tau = bnk->gamma2_i;
195eb910715SAlp Dener               } else if (tau_max < bnk->gamma1_i) {
196eb910715SAlp Dener                 tau = bnk->gamma1_i;
197eb910715SAlp Dener               } else if ((tau_min < bnk->gamma1_i) && (tau_max >= 1.0)) {
198eb910715SAlp Dener                 tau = bnk->gamma1_i;
199eb910715SAlp Dener               } else if ((tau_1 >= bnk->gamma1_i) && (tau_1 < 1.0) && ((tau_2 < bnk->gamma1_i) || (tau_2 >= 1.0))) {
200eb910715SAlp Dener                 tau = tau_1;
201eb910715SAlp Dener               } else if ((tau_2 >= bnk->gamma1_i) && (tau_2 < 1.0) && ((tau_1 < bnk->gamma1_i) || (tau_2 >= 1.0))) {
202eb910715SAlp Dener                 tau = tau_2;
203eb910715SAlp Dener               } else {
204eb910715SAlp Dener                 tau = tau_max;
205eb910715SAlp Dener               }
206eb910715SAlp Dener             }
207eb910715SAlp Dener           }
208eb910715SAlp Dener           tao->trust = tau * tao->trust;
209eb910715SAlp Dener         }
210eb910715SAlp Dener 
211eb910715SAlp Dener         if (fmin < bnk->f) {
212eb910715SAlp Dener           bnk->f = fmin;
213eb910715SAlp Dener           ierr = VecAXPY(tao->solution,sigma,tao->gradient);CHKERRQ(ierr);
21409164190SAlp Dener           ierr = TaoComputeGradient(tao,tao->solution,bnk->unprojected_gradient);CHKERRQ(ierr);
21509164190SAlp Dener           ierr = VecBoundGradientProjection(bnk->unprojected_gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr);
216eb910715SAlp Dener 
217eb910715SAlp Dener           ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&bnk->gnorm);CHKERRQ(ierr);
218eb910715SAlp Dener           if (PetscIsInfOrNanReal(bnk->gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute gradient generated Inf or NaN");
219eb910715SAlp Dener           needH = 1;
220eb910715SAlp Dener 
221eb910715SAlp Dener           ierr = TaoLogConvergenceHistory(tao,bnk->f,bnk->gnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
222eb910715SAlp Dener           ierr = TaoMonitor(tao,tao->niter,bnk->f,bnk->gnorm,0.0,step);CHKERRQ(ierr);
223eb910715SAlp Dener           ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
224eb910715SAlp Dener           if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
225eb910715SAlp Dener         }
226eb910715SAlp Dener       }
227eb910715SAlp Dener       tao->trust = PetscMax(tao->trust, max_radius);
228eb910715SAlp Dener 
229eb910715SAlp Dener       /* Modify the radius if it is too large or small */
230eb910715SAlp Dener       tao->trust = PetscMax(tao->trust, bnk->min_radius);
231eb910715SAlp Dener       tao->trust = PetscMin(tao->trust, bnk->max_radius);
232eb910715SAlp Dener       break;
233eb910715SAlp Dener 
234eb910715SAlp Dener     default:
235eb910715SAlp Dener       /* Norm of the first direction will initialize radius */
236eb910715SAlp Dener       tao->trust = 0.0;
237eb910715SAlp Dener       break;
238eb910715SAlp Dener     }
239eb910715SAlp Dener   }
240eb910715SAlp Dener 
241eb910715SAlp Dener   /* Set initial scaling for the BFGS preconditioner
242eb910715SAlp Dener      This step is done after computing the initial trust-region radius
243eb910715SAlp Dener      since the function value may have decreased */
244eb910715SAlp Dener   if (BNK_PC_BFGS == bnk->pc_type) {
245eb910715SAlp Dener     if (bnk->f != 0.0) {
246eb910715SAlp Dener       delta = 2.0 * PetscAbsScalar(bnk->f) / (bnk->gnorm*bnk->gnorm);
247eb910715SAlp Dener     } else {
248eb910715SAlp Dener       delta = 2.0 / (bnk->gnorm*bnk->gnorm);
249eb910715SAlp Dener     }
250eb910715SAlp Dener     ierr = MatLMVMSetDelta(bnk->M,delta);CHKERRQ(ierr);
251eb910715SAlp Dener   }
252eb910715SAlp Dener 
253eb910715SAlp Dener   /* Set counter for gradient/reset steps*/
254eb910715SAlp Dener   bnk->newt = 0;
255eb910715SAlp Dener   bnk->bfgs = 0;
256eb910715SAlp Dener   bnk->sgrad = 0;
257eb910715SAlp Dener   bnk->grad = 0;
258eb910715SAlp Dener   PetscFunctionReturn(0);
259eb910715SAlp Dener }
260eb910715SAlp Dener 
261*df278d8fSAlp Dener /*------------------------------------------------------------*/
262*df278d8fSAlp Dener 
263*df278d8fSAlp Dener /* Routine for computing the Newton step.
264*df278d8fSAlp Dener 
265*df278d8fSAlp Dener   If the safeguard is enabled, the Newton step is verified to be a
266*df278d8fSAlp Dener   descent direction, with fallbacks onto BFGS, scaled gradient, and unscaled
267*df278d8fSAlp Dener   gradient steps if/when necessary.
268*df278d8fSAlp Dener 
269*df278d8fSAlp Dener   The function reports back on which type of step has ultimately been stored
270*df278d8fSAlp Dener   under tao->stepdirection.
271*df278d8fSAlp Dener */
272*df278d8fSAlp Dener 
273fed79b8eSAlp Dener PetscErrorCode TaoBNKComputeStep(Tao tao, PetscBool safeguard, PetscInt *stepType)
274eb910715SAlp Dener {
275eb910715SAlp Dener   PetscErrorCode               ierr;
276eb910715SAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
277eb910715SAlp Dener   KSPConvergedReason           ksp_reason;
278eb910715SAlp Dener 
279080d2917SAlp Dener   PetscReal                    gdx, delta, e_min;
280eb910715SAlp Dener 
281eb910715SAlp Dener   PetscInt                     bfgsUpdates = 0;
282eb910715SAlp Dener   PetscInt                     kspits;
283eb910715SAlp Dener 
284eb910715SAlp Dener   PetscFunctionBegin;
285eb910715SAlp Dener   /* Shift the Hessian matrix */
286eb910715SAlp Dener   if (bnk->pert > 0) {
287eb910715SAlp Dener     ierr = MatShift(tao->hessian, bnk->pert);CHKERRQ(ierr);
288eb910715SAlp Dener     if (tao->hessian != tao->hessian_pre) {
289eb910715SAlp Dener       ierr = MatShift(tao->hessian_pre, bnk->pert);CHKERRQ(ierr);
290eb910715SAlp Dener     }
291eb910715SAlp Dener   }
292eb910715SAlp Dener 
29309164190SAlp Dener   /* Determine the inactive set */
29409164190SAlp Dener   ierr = ISDestroy(&bnk->inactive_idx);CHKERRQ(ierr);
29509164190SAlp Dener   ierr = VecWhichInactive(tao->XL,tao->solution,bnk->unprojected_gradient,tao->XU,PETSC_TRUE,&bnk->inactive_idx);CHKERRQ(ierr);
29609164190SAlp Dener 
297eb3ba6a7SAlp Dener   /* Prepare masked matrices for the inactive set */
29809164190SAlp Dener   ierr = MatLMVMSetInactive(bnk->M, bnk->inactive_idx);CHKERRQ(ierr);
29909164190SAlp Dener   ierr = VecSet(tao->stepdirection, 0.0);CHKERRQ(ierr);
300eb3ba6a7SAlp Dener   ierr = TaoMatGetSubMat(tao->hessian, bnk->inactive_idx, bnk->Xwork, TAO_SUBSET_MASK, &bnk->H_inactive);CHKERRQ(ierr);
301eb3ba6a7SAlp Dener   if (tao->hessian == tao->hessian_pre) {
302eb3ba6a7SAlp Dener     bnk->Hpre_inactive = bnk->H_inactive;
303eb3ba6a7SAlp Dener   } else {
304eb3ba6a7SAlp Dener     ierr = TaoMatGetSubMat(tao->hessian_pre, bnk->inactive_idx, bnk->Xwork, TAO_SUBSET_MASK, &bnk->Hpre_inactive);CHKERRQ(ierr);
305eb3ba6a7SAlp Dener   }
30609164190SAlp Dener 
307eb910715SAlp Dener   /* Solve the Newton system of equations */
30809164190SAlp Dener   ierr = KSPSetOperators(tao->ksp,bnk->H_inactive,bnk->Hpre_inactive);CHKERRQ(ierr);
309eb910715SAlp Dener   if (bnk->is_nash || bnk->is_stcg || bnk->is_gltr) {
310fed79b8eSAlp Dener     ierr = KSPCGSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr);
311eb3ba6a7SAlp Dener     ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr);
312eb910715SAlp Dener     ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr);
313eb910715SAlp Dener     tao->ksp_its+=kspits;
314eb910715SAlp Dener     tao->ksp_tot_its+=kspits;
315080d2917SAlp Dener     ierr = KSPCGGetNormD(tao->ksp,&bnk->dnorm);CHKERRQ(ierr);
316eb910715SAlp Dener 
317eb910715SAlp Dener     if (0.0 == tao->trust) {
318eb910715SAlp Dener       /* Radius was uninitialized; use the norm of the direction */
319080d2917SAlp Dener       if (bnk->dnorm > 0.0) {
320080d2917SAlp Dener         tao->trust = bnk->dnorm;
321eb910715SAlp Dener 
322eb910715SAlp Dener         /* Modify the radius if it is too large or small */
323eb910715SAlp Dener         tao->trust = PetscMax(tao->trust, bnk->min_radius);
324eb910715SAlp Dener         tao->trust = PetscMin(tao->trust, bnk->max_radius);
325eb910715SAlp Dener       } else {
326eb910715SAlp Dener         /* The direction was bad; set radius to default value and re-solve
327eb910715SAlp Dener            the trust-region subproblem to get a direction */
328eb910715SAlp Dener         tao->trust = tao->trust0;
329eb910715SAlp Dener 
330eb910715SAlp Dener         /* Modify the radius if it is too large or small */
331eb910715SAlp Dener         tao->trust = PetscMax(tao->trust, bnk->min_radius);
332eb910715SAlp Dener         tao->trust = PetscMin(tao->trust, bnk->max_radius);
333eb910715SAlp Dener 
334fed79b8eSAlp Dener         ierr = KSPCGSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr);
335eb3ba6a7SAlp Dener         ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr);
336eb910715SAlp Dener         ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr);
337eb910715SAlp Dener         tao->ksp_its+=kspits;
338eb910715SAlp Dener         tao->ksp_tot_its+=kspits;
339080d2917SAlp Dener         ierr = KSPCGGetNormD(tao->ksp,&bnk->dnorm);CHKERRQ(ierr);
340eb910715SAlp Dener 
341080d2917SAlp Dener         if (bnk->dnorm == 0.0) SETERRQ(PETSC_COMM_SELF,1, "Initial direction zero");
342eb910715SAlp Dener       }
343eb910715SAlp Dener     }
344eb910715SAlp Dener   } else {
345eb3ba6a7SAlp Dener     ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr);
346eb910715SAlp Dener     ierr = KSPGetIterationNumber(tao->ksp, &kspits);CHKERRQ(ierr);
347eb910715SAlp Dener     tao->ksp_its += kspits;
348eb910715SAlp Dener     tao->ksp_tot_its+=kspits;
349eb910715SAlp Dener   }
350fed79b8eSAlp Dener   ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
351fed79b8eSAlp Dener   ierr = KSPGetConvergedReason(tao->ksp, &ksp_reason);CHKERRQ(ierr);
352fed79b8eSAlp Dener 
353eb3ba6a7SAlp Dener   /* Destroy masked matrices */
354eb3ba6a7SAlp Dener   if (bnk->H_inactive != bnk->Hpre_inactive) {
35509164190SAlp Dener     ierr = MatDestroy(&bnk->Hpre_inactive);CHKERRQ(ierr);
356eb3ba6a7SAlp Dener   }
357eb3ba6a7SAlp Dener   ierr = MatDestroy(&bnk->H_inactive);CHKERRQ(ierr);
35809164190SAlp Dener 
359fed79b8eSAlp Dener   /* Make sure the BFGS preconditioner is healthy */
360fed79b8eSAlp Dener   if (bnk->pc_type == BNK_PC_BFGS) {
361fed79b8eSAlp Dener     ierr = MatLMVMGetUpdates(bnk->M, &bfgsUpdates);CHKERRQ(ierr);
362fed79b8eSAlp Dener     if ((KSP_DIVERGED_INDEFINITE_PC == ksp_reason) && (bfgsUpdates > 1)) {
363fed79b8eSAlp Dener       /* Preconditioner is numerically indefinite; reset the approximation. */
364eb910715SAlp Dener       if (bnk->f != 0.0) {
365eb910715SAlp Dener         delta = 2.0 * PetscAbsScalar(bnk->f) / (bnk->gnorm*bnk->gnorm);
366eb910715SAlp Dener       } else {
367eb910715SAlp Dener         delta = 2.0 / (bnk->gnorm*bnk->gnorm);
368eb910715SAlp Dener       }
369eb910715SAlp Dener       ierr = MatLMVMSetDelta(bnk->M,delta);CHKERRQ(ierr);
370eb910715SAlp Dener       ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr);
37109164190SAlp Dener       ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
372eb910715SAlp Dener       bfgsUpdates = 1;
373eb910715SAlp Dener     }
374fed79b8eSAlp Dener   }
375eb910715SAlp Dener 
376eb910715SAlp Dener   if (KSP_CONVERGED_ATOL == ksp_reason) {
377eb910715SAlp Dener     ++bnk->ksp_atol;
378eb910715SAlp Dener   } else if (KSP_CONVERGED_RTOL == ksp_reason) {
379eb910715SAlp Dener     ++bnk->ksp_rtol;
380eb910715SAlp Dener   } else if (KSP_CONVERGED_CG_CONSTRAINED == ksp_reason) {
381eb910715SAlp Dener     ++bnk->ksp_ctol;
382eb910715SAlp Dener   } else if (KSP_CONVERGED_CG_NEG_CURVE == ksp_reason) {
383eb910715SAlp Dener     ++bnk->ksp_negc;
384eb910715SAlp Dener   } else if (KSP_DIVERGED_DTOL == ksp_reason) {
385eb910715SAlp Dener     ++bnk->ksp_dtol;
386eb910715SAlp Dener   } else if (KSP_DIVERGED_ITS == ksp_reason) {
387eb910715SAlp Dener     ++bnk->ksp_iter;
388eb910715SAlp Dener   } else {
389eb910715SAlp Dener     ++bnk->ksp_othr;
390eb910715SAlp Dener   }
391eb910715SAlp Dener 
392eb910715SAlp Dener   /* Check for success (descent direction) */
393fed79b8eSAlp Dener   if (safeguard) {
394080d2917SAlp Dener     ierr = VecDot(tao->stepdirection, tao->gradient, &gdx);CHKERRQ(ierr);
395eb910715SAlp Dener     if ((gdx >= 0.0) || PetscIsInfOrNanReal(gdx)) {
396eb910715SAlp Dener       /* Newton step is not descent or direction produced Inf or NaN
397eb910715SAlp Dener          Update the perturbation for next time */
398eb910715SAlp Dener       if (bnk->pert <= 0.0) {
399eb910715SAlp Dener         /* Initialize the perturbation */
400eb910715SAlp Dener         bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
401eb910715SAlp Dener         if (bnk->is_gltr) {
402eb910715SAlp Dener           ierr = KSPCGGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr);
403eb910715SAlp Dener           bnk->pert = PetscMax(bnk->pert, -e_min);
404eb910715SAlp Dener         }
405eb910715SAlp Dener       } else {
406eb910715SAlp Dener         /* Increase the perturbation */
407eb910715SAlp Dener         bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
408eb910715SAlp Dener       }
409eb910715SAlp Dener 
410eb910715SAlp Dener       if (BNK_PC_BFGS != bnk->pc_type) {
411eb910715SAlp Dener         /* We don't have the bfgs matrix around and updated
412eb910715SAlp Dener            Must use gradient direction in this case */
413080d2917SAlp Dener         ierr = VecCopy(tao->gradient, tao->stepdirection);CHKERRQ(ierr);
414080d2917SAlp Dener         ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
415eb910715SAlp Dener         ++bnk->grad;
416eb910715SAlp Dener         *stepType = BNK_GRADIENT;
417eb910715SAlp Dener       } else {
418eb910715SAlp Dener         /* Attempt to use the BFGS direction */
41909164190SAlp Dener         ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
42009164190SAlp Dener         ierr = VecBoundGradientProjection(tao->stepdirection,tao->solution,tao->XL,tao->XU,tao->stepdirection);CHKERRQ(ierr);
421080d2917SAlp Dener         ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
422eb910715SAlp Dener 
423eb910715SAlp Dener         /* Check for success (descent direction) */
424080d2917SAlp Dener         ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr);
425eb910715SAlp Dener         if ((gdx >= 0) || PetscIsInfOrNanReal(gdx)) {
426eb910715SAlp Dener           /* BFGS direction is not descent or direction produced not a number
427eb910715SAlp Dener              We can assert bfgsUpdates > 1 in this case because
428eb910715SAlp Dener              the first solve produces the scaled gradient direction,
429eb910715SAlp Dener              which is guaranteed to be descent */
430eb910715SAlp Dener 
431eb910715SAlp Dener           /* Use steepest descent direction (scaled) */
432eb910715SAlp Dener 
433eb910715SAlp Dener           if (bnk->f != 0.0) {
434eb910715SAlp Dener             delta = 2.0 * PetscAbsScalar(bnk->f) / (bnk->gnorm*bnk->gnorm);
435eb910715SAlp Dener           } else {
436eb910715SAlp Dener             delta = 2.0 / (bnk->gnorm*bnk->gnorm);
437eb910715SAlp Dener           }
438eb910715SAlp Dener           ierr = MatLMVMSetDelta(bnk->M, delta);CHKERRQ(ierr);
439eb910715SAlp Dener           ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr);
44009164190SAlp Dener           ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
44109164190SAlp Dener           ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
44209164190SAlp Dener           ierr = VecBoundGradientProjection(tao->stepdirection,tao->solution,tao->XL,tao->XU,tao->stepdirection);CHKERRQ(ierr);
443080d2917SAlp Dener           ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
444eb910715SAlp Dener 
445eb910715SAlp Dener           bfgsUpdates = 1;
446eb910715SAlp Dener           ++bnk->sgrad;
447eb910715SAlp Dener           *stepType = BNK_SCALED_GRADIENT;
448eb910715SAlp Dener         } else {
449eb910715SAlp Dener           if (1 == bfgsUpdates) {
450eb910715SAlp Dener             /* The first BFGS direction is always the scaled gradient */
451eb910715SAlp Dener             ++bnk->sgrad;
452eb910715SAlp Dener             *stepType = BNK_SCALED_GRADIENT;
453eb910715SAlp Dener           } else {
454eb910715SAlp Dener             ++bnk->bfgs;
455eb910715SAlp Dener             *stepType = BNK_BFGS;
456eb910715SAlp Dener           }
457eb910715SAlp Dener         }
458eb910715SAlp Dener       }
459eb910715SAlp Dener     } else {
460eb910715SAlp Dener       /* Computed Newton step is descent */
461eb910715SAlp Dener       switch (ksp_reason) {
462eb910715SAlp Dener       case KSP_DIVERGED_NANORINF:
463eb910715SAlp Dener       case KSP_DIVERGED_BREAKDOWN:
464eb910715SAlp Dener       case KSP_DIVERGED_INDEFINITE_MAT:
465eb910715SAlp Dener       case KSP_DIVERGED_INDEFINITE_PC:
466eb910715SAlp Dener       case KSP_CONVERGED_CG_NEG_CURVE:
467eb910715SAlp Dener         /* Matrix or preconditioner is indefinite; increase perturbation */
468eb910715SAlp Dener         if (bnk->pert <= 0.0) {
469eb910715SAlp Dener           /* Initialize the perturbation */
470eb910715SAlp Dener           bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
471eb910715SAlp Dener           if (bnk->is_gltr) {
472eb910715SAlp Dener             ierr = KSPCGGLTRGetMinEig(tao->ksp, &e_min);CHKERRQ(ierr);
473eb910715SAlp Dener             bnk->pert = PetscMax(bnk->pert, -e_min);
474eb910715SAlp Dener           }
475eb910715SAlp Dener         } else {
476eb910715SAlp Dener           /* Increase the perturbation */
477eb910715SAlp Dener           bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
478eb910715SAlp Dener         }
479eb910715SAlp Dener         break;
480eb910715SAlp Dener 
481eb910715SAlp Dener       default:
482eb910715SAlp Dener         /* Newton step computation is good; decrease perturbation */
483eb910715SAlp Dener         bnk->pert = PetscMin(bnk->psfac * bnk->pert, bnk->pmsfac * bnk->gnorm);
484eb910715SAlp Dener         if (bnk->pert < bnk->pmin) {
485eb910715SAlp Dener           bnk->pert = 0.0;
486eb910715SAlp Dener         }
487eb910715SAlp Dener         break;
488eb910715SAlp Dener       }
489eb910715SAlp Dener 
490eb910715SAlp Dener       ++bnk->newt;
491fed79b8eSAlp Dener       *stepType = BNK_NEWTON;
492fed79b8eSAlp Dener     }
493fed79b8eSAlp Dener   } else {
494fed79b8eSAlp Dener     ++bnk->newt;
495fed79b8eSAlp Dener     bnk->pert = 0.0;
496fed79b8eSAlp Dener     *stepType = BNK_NEWTON;
497eb910715SAlp Dener   }
498eb910715SAlp Dener   PetscFunctionReturn(0);
499eb910715SAlp Dener }
500eb910715SAlp Dener 
501*df278d8fSAlp Dener /*------------------------------------------------------------*/
502*df278d8fSAlp Dener 
503*df278d8fSAlp Dener /* Routine for performing a bound-projected More-Thuente line search.
504*df278d8fSAlp Dener 
505*df278d8fSAlp Dener   Includes fallbacks to BFGS, scaled gradient, and unscaled gradient steps if the
506*df278d8fSAlp Dener   Newton step does not produce a valid step length.
507*df278d8fSAlp Dener */
508*df278d8fSAlp Dener 
509c14b763aSAlp Dener PetscErrorCode TaoBNKPerformLineSearch(Tao tao, PetscInt stepType, PetscReal *steplen, TaoLineSearchConvergedReason *reason)
510c14b763aSAlp Dener {
511c14b763aSAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
512c14b763aSAlp Dener   PetscErrorCode ierr;
513c14b763aSAlp Dener   TaoLineSearchConvergedReason ls_reason;
514c14b763aSAlp Dener 
515c14b763aSAlp Dener   PetscReal      e_min, gdx, delta;
516c14b763aSAlp Dener   PetscInt       bfgsUpdates;
517c14b763aSAlp Dener 
518c14b763aSAlp Dener   PetscFunctionBegin;
519c14b763aSAlp Dener   /* Perform the linesearch */
520c14b763aSAlp Dener   ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, steplen, &ls_reason);CHKERRQ(ierr);
521c14b763aSAlp Dener   ierr = VecBoundGradientProjection(bnk->unprojected_gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr);
522c14b763aSAlp Dener   ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
523c14b763aSAlp Dener 
524c14b763aSAlp Dener   while (ls_reason != TAOLINESEARCH_SUCCESS && ls_reason != TAOLINESEARCH_SUCCESS_USER && stepType != BNK_GRADIENT) {
525c14b763aSAlp Dener     /* Linesearch failed, revert solution */
526c14b763aSAlp Dener     bnk->f = bnk->fold;
527c14b763aSAlp Dener     ierr = VecCopy(bnk->Xold, tao->solution);CHKERRQ(ierr);
528c14b763aSAlp Dener     ierr = VecCopy(bnk->Gold, tao->gradient);CHKERRQ(ierr);
529c14b763aSAlp Dener     ierr = VecCopy(bnk->unprojected_gradient_old, bnk->unprojected_gradient);CHKERRQ(ierr);
530c14b763aSAlp Dener 
531c14b763aSAlp Dener     switch(stepType) {
532c14b763aSAlp Dener     case BNK_NEWTON:
533c14b763aSAlp Dener       /* Failed to obtain acceptable iterate with Newton 1step
534c14b763aSAlp Dener          Update the perturbation for next time */
535c14b763aSAlp Dener       --bnk->newt;
536c14b763aSAlp Dener       if (bnk->pert <= 0.0) {
537c14b763aSAlp Dener         /* Initialize the perturbation */
538c14b763aSAlp Dener         bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
539c14b763aSAlp Dener         if (bnk->is_gltr) {
540c14b763aSAlp Dener           ierr = KSPCGGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr);
541c14b763aSAlp Dener           bnk->pert = PetscMax(bnk->pert, -e_min);
542c14b763aSAlp Dener         }
543c14b763aSAlp Dener       } else {
544c14b763aSAlp Dener         /* Increase the perturbation */
545c14b763aSAlp Dener         bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
546c14b763aSAlp Dener       }
547c14b763aSAlp Dener 
548c14b763aSAlp Dener       if (BNK_PC_BFGS != bnk->pc_type) {
549c14b763aSAlp Dener         /* We don't have the bfgs matrix around and being updated
550c14b763aSAlp Dener            Must use gradient direction in this case */
551c14b763aSAlp Dener         ierr = VecCopy(tao->gradient, tao->stepdirection);CHKERRQ(ierr);
552c14b763aSAlp Dener         ++bnk->grad;
553c14b763aSAlp Dener         stepType = BNK_GRADIENT;
554c14b763aSAlp Dener       } else {
555c14b763aSAlp Dener         /* Attempt to use the BFGS direction */
556c14b763aSAlp Dener         ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
557c14b763aSAlp Dener         ierr = VecBoundGradientProjection(tao->stepdirection,tao->solution,tao->XL,tao->XU,tao->stepdirection);CHKERRQ(ierr);
558c14b763aSAlp Dener         /* Check for success (descent direction) */
559c14b763aSAlp Dener         ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr);
560c14b763aSAlp Dener         if ((gdx <= 0) || PetscIsInfOrNanReal(gdx)) {
561c14b763aSAlp Dener           /* BFGS direction is not descent or direction produced not a number
562c14b763aSAlp Dener              We can assert bfgsUpdates > 1 in this case
563c14b763aSAlp Dener              Use steepest descent direction (scaled) */
564c14b763aSAlp Dener 
565c14b763aSAlp Dener           if (bnk->f != 0.0) {
566c14b763aSAlp Dener             delta = 2.0 * PetscAbsScalar(bnk->f) / (bnk->gnorm*bnk->gnorm);
567c14b763aSAlp Dener           } else {
568c14b763aSAlp Dener             delta = 2.0 / (bnk->gnorm*bnk->gnorm);
569c14b763aSAlp Dener           }
570c14b763aSAlp Dener           ierr = MatLMVMSetDelta(bnk->M, delta);CHKERRQ(ierr);
571c14b763aSAlp Dener           ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr);
572c14b763aSAlp Dener           ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
573c14b763aSAlp Dener           ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
574c14b763aSAlp Dener           ierr = VecBoundGradientProjection(tao->stepdirection,tao->solution,tao->XL,tao->XU,tao->stepdirection);CHKERRQ(ierr);
575c14b763aSAlp Dener 
576c14b763aSAlp Dener           bfgsUpdates = 1;
577c14b763aSAlp Dener           ++bnk->sgrad;
578c14b763aSAlp Dener           stepType = BNK_SCALED_GRADIENT;
579c14b763aSAlp Dener         } else {
580c14b763aSAlp Dener           ierr = MatLMVMGetUpdates(bnk->M, &bfgsUpdates);CHKERRQ(ierr);
581c14b763aSAlp Dener           if (1 == bfgsUpdates) {
582c14b763aSAlp Dener             /* The first BFGS direction is always the scaled gradient */
583c14b763aSAlp Dener             ++bnk->sgrad;
584c14b763aSAlp Dener             stepType = BNK_SCALED_GRADIENT;
585c14b763aSAlp Dener           } else {
586c14b763aSAlp Dener             ++bnk->bfgs;
587c14b763aSAlp Dener             stepType = BNK_BFGS;
588c14b763aSAlp Dener           }
589c14b763aSAlp Dener         }
590c14b763aSAlp Dener       }
591c14b763aSAlp Dener       break;
592c14b763aSAlp Dener 
593c14b763aSAlp Dener     case BNK_BFGS:
594c14b763aSAlp Dener       /* Can only enter if pc_type == BNK_PC_BFGS
595c14b763aSAlp Dener          Failed to obtain acceptable iterate with BFGS step
596c14b763aSAlp Dener          Attempt to use the scaled gradient direction */
597c14b763aSAlp Dener 
598c14b763aSAlp Dener       if (bnk->f != 0.0) {
599c14b763aSAlp Dener         delta = 2.0 * PetscAbsScalar(bnk->f) / (bnk->gnorm*bnk->gnorm);
600c14b763aSAlp Dener       } else {
601c14b763aSAlp Dener         delta = 2.0 / (bnk->gnorm*bnk->gnorm);
602c14b763aSAlp Dener       }
603c14b763aSAlp Dener       ierr = MatLMVMSetDelta(bnk->M, delta);CHKERRQ(ierr);
604c14b763aSAlp Dener       ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr);
605c14b763aSAlp Dener       ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
606c14b763aSAlp Dener       ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
607c14b763aSAlp Dener       ierr = VecBoundGradientProjection(tao->stepdirection,tao->solution,tao->XL,tao->XU,tao->stepdirection);CHKERRQ(ierr);
608c14b763aSAlp Dener 
609c14b763aSAlp Dener       bfgsUpdates = 1;
610c14b763aSAlp Dener       ++bnk->sgrad;
611c14b763aSAlp Dener       stepType = BNK_SCALED_GRADIENT;
612c14b763aSAlp Dener       break;
613c14b763aSAlp Dener 
614c14b763aSAlp Dener     case BNK_SCALED_GRADIENT:
615c14b763aSAlp Dener       /* Can only enter if pc_type == BNK_PC_BFGS
616c14b763aSAlp Dener          The scaled gradient step did not produce a new iterate;
617c14b763aSAlp Dener          attemp to use the gradient direction.
618c14b763aSAlp Dener          Need to make sure we are not using a different diagonal scaling */
619c14b763aSAlp Dener 
620c14b763aSAlp Dener       ierr = MatLMVMSetScale(bnk->M,0);CHKERRQ(ierr);
621c14b763aSAlp Dener       ierr = MatLMVMSetDelta(bnk->M,1.0);CHKERRQ(ierr);
622c14b763aSAlp Dener       ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr);
623c14b763aSAlp Dener       ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
624c14b763aSAlp Dener       ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
625c14b763aSAlp Dener       ierr = VecBoundGradientProjection(tao->stepdirection,tao->solution,tao->XL,tao->XU,tao->stepdirection);CHKERRQ(ierr);
626c14b763aSAlp Dener 
627c14b763aSAlp Dener       bfgsUpdates = 1;
628c14b763aSAlp Dener       ++bnk->grad;
629c14b763aSAlp Dener       stepType = BNK_GRADIENT;
630c14b763aSAlp Dener       break;
631c14b763aSAlp Dener     }
632c14b763aSAlp Dener     ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
633c14b763aSAlp Dener 
634c14b763aSAlp Dener     ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &bnk->f, bnk->unprojected_gradient, tao->stepdirection, steplen, &ls_reason);CHKERRQ(ierr);
635c14b763aSAlp Dener     ierr = VecBoundGradientProjection(bnk->unprojected_gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr);
636c14b763aSAlp Dener     ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
637c14b763aSAlp Dener   }
638c14b763aSAlp Dener   *reason = ls_reason;
639c14b763aSAlp Dener   PetscFunctionReturn(0);
640c14b763aSAlp Dener }
641c14b763aSAlp Dener 
642*df278d8fSAlp Dener /*------------------------------------------------------------*/
643*df278d8fSAlp Dener 
644*df278d8fSAlp Dener /* Routine for updating the trust radius.
645*df278d8fSAlp Dener 
646*df278d8fSAlp Dener   Function features three different update methods:
647*df278d8fSAlp Dener   1) Line-search step length based
648*df278d8fSAlp Dener   2) Predicted decrease on the CG quadratic model
649*df278d8fSAlp Dener   3) Interpolation
650*df278d8fSAlp Dener */
651*df278d8fSAlp Dener 
652b1c2d0e3SAlp Dener PetscErrorCode TaoBNKUpdateTrustRadius(Tao tao, PetscReal prered, PetscReal actred, PetscInt stepType, PetscBool *accept)
653080d2917SAlp Dener {
654080d2917SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
655080d2917SAlp Dener   PetscErrorCode ierr;
656080d2917SAlp Dener 
657b1c2d0e3SAlp Dener   PetscReal      step, kappa;
658080d2917SAlp Dener   PetscReal      gdx, tau_1, tau_2, tau_min, tau_max;
659080d2917SAlp Dener 
660080d2917SAlp Dener   PetscFunctionBegin;
661080d2917SAlp Dener   /* Update trust region radius */
662080d2917SAlp Dener   *accept = PETSC_FALSE;
663080d2917SAlp Dener   switch(bnk->update_type) {
664080d2917SAlp Dener   case BNK_UPDATE_STEP:
665c14b763aSAlp Dener     *accept = PETSC_TRUE; /* always accept here because line search succeeded */
666080d2917SAlp Dener     if (stepType == BNK_NEWTON) {
667080d2917SAlp Dener       ierr = TaoLineSearchGetStepLength(tao->linesearch, &step);CHKERRQ(ierr);
668080d2917SAlp Dener       if (step < bnk->nu1) {
669080d2917SAlp Dener         /* Very bad step taken; reduce radius */
670080d2917SAlp Dener         tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust);
671080d2917SAlp Dener       } else if (step < bnk->nu2) {
672080d2917SAlp Dener         /* Reasonably bad step taken; reduce radius */
673080d2917SAlp Dener         tao->trust = bnk->omega2 * PetscMin(bnk->dnorm, tao->trust);
674080d2917SAlp Dener       } else if (step < bnk->nu3) {
675080d2917SAlp Dener         /*  Reasonable step was taken; leave radius alone */
676080d2917SAlp Dener         if (bnk->omega3 < 1.0) {
677080d2917SAlp Dener           tao->trust = bnk->omega3 * PetscMin(bnk->dnorm, tao->trust);
678080d2917SAlp Dener         } else if (bnk->omega3 > 1.0) {
679080d2917SAlp Dener           tao->trust = PetscMax(bnk->omega3 * bnk->dnorm, tao->trust);
680080d2917SAlp Dener         }
681080d2917SAlp Dener       } else if (step < bnk->nu4) {
682080d2917SAlp Dener         /*  Full step taken; increase the radius */
683080d2917SAlp Dener         tao->trust = PetscMax(bnk->omega4 * bnk->dnorm, tao->trust);
684080d2917SAlp Dener       } else {
685080d2917SAlp Dener         /*  More than full step taken; increase the radius */
686080d2917SAlp Dener         tao->trust = PetscMax(bnk->omega5 * bnk->dnorm, tao->trust);
687080d2917SAlp Dener       }
688080d2917SAlp Dener     } else {
689080d2917SAlp Dener       /*  Newton step was not good; reduce the radius */
690080d2917SAlp Dener       tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust);
691080d2917SAlp Dener     }
692080d2917SAlp Dener     break;
693080d2917SAlp Dener 
694080d2917SAlp Dener   case BNK_UPDATE_REDUCTION:
695080d2917SAlp Dener     if (stepType == BNK_NEWTON) {
696b1c2d0e3SAlp Dener       if (prered < 0.0) {
697fed79b8eSAlp Dener         /* The predicted reduction has the wrong sign.  This cannot
698fed79b8eSAlp Dener            happen in infinite precision arithmetic.  Step should
699fed79b8eSAlp Dener            be rejected! */
700080d2917SAlp Dener         tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
701fed79b8eSAlp Dener       }
702fed79b8eSAlp Dener       else {
703b1c2d0e3SAlp Dener         if (PetscIsInfOrNanReal(actred)) {
704080d2917SAlp Dener           tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
705080d2917SAlp Dener         } else {
706080d2917SAlp Dener           if ((PetscAbsScalar(actred) <= bnk->epsilon) &&
707080d2917SAlp Dener               (PetscAbsScalar(prered) <= bnk->epsilon)) {
708080d2917SAlp Dener             kappa = 1.0;
709fed79b8eSAlp Dener           }
710fed79b8eSAlp Dener           else {
711080d2917SAlp Dener             kappa = actred / prered;
712080d2917SAlp Dener           }
713fed79b8eSAlp Dener 
714fed79b8eSAlp Dener           /* Accept or reject the step and update radius */
715080d2917SAlp Dener           if (kappa < bnk->eta1) {
716fed79b8eSAlp Dener             /* Reject the step */
717080d2917SAlp Dener             tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
718fed79b8eSAlp Dener           }
719fed79b8eSAlp Dener           else {
720fed79b8eSAlp Dener             /* Accept the step */
721080d2917SAlp Dener             if (kappa < bnk->eta2) {
722080d2917SAlp Dener               /* Marginal bad step */
723080d2917SAlp Dener               tao->trust = bnk->alpha2 * PetscMin(tao->trust, bnk->dnorm);
724080d2917SAlp Dener             }
725fed79b8eSAlp Dener             else if (kappa < bnk->eta3) {
726fed79b8eSAlp Dener               /* Reasonable step */
727fed79b8eSAlp Dener               tao->trust = bnk->alpha3 * tao->trust;
728fed79b8eSAlp Dener             }
729fed79b8eSAlp Dener             else if (kappa < bnk->eta4) {
730080d2917SAlp Dener               /* Good step */
731080d2917SAlp Dener               tao->trust = PetscMax(bnk->alpha4 * bnk->dnorm, tao->trust);
732fed79b8eSAlp Dener             }
733fed79b8eSAlp Dener             else {
734080d2917SAlp Dener               /* Very good step */
735080d2917SAlp Dener               tao->trust = PetscMax(bnk->alpha5 * bnk->dnorm, tao->trust);
736080d2917SAlp Dener             }
737fed79b8eSAlp Dener             *accept = PETSC_TRUE;
738080d2917SAlp Dener           }
739080d2917SAlp Dener         }
740080d2917SAlp Dener       }
741080d2917SAlp Dener     } else {
742080d2917SAlp Dener       /*  Newton step was not good; reduce the radius */
743080d2917SAlp Dener       tao->trust = bnk->alpha1 * PetscMin(bnk->dnorm, tao->trust);
744080d2917SAlp Dener     }
745080d2917SAlp Dener     break;
746080d2917SAlp Dener 
747080d2917SAlp Dener   default:
748080d2917SAlp Dener     if (stepType == BNK_NEWTON) {
749b1c2d0e3SAlp Dener       if (prered < 0.0) {
750080d2917SAlp Dener         /*  The predicted reduction has the wrong sign.  This cannot */
751080d2917SAlp Dener         /*  happen in infinite precision arithmetic.  Step should */
752080d2917SAlp Dener         /*  be rejected! */
753080d2917SAlp Dener         tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
754080d2917SAlp Dener       } else {
755b1c2d0e3SAlp Dener         if (PetscIsInfOrNanReal(actred)) {
756080d2917SAlp Dener           tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
757080d2917SAlp Dener         } else {
758080d2917SAlp Dener           if ((PetscAbsScalar(actred) <= bnk->epsilon) && (PetscAbsScalar(prered) <= bnk->epsilon)) {
759080d2917SAlp Dener             kappa = 1.0;
760080d2917SAlp Dener           } else {
761080d2917SAlp Dener             kappa = actred / prered;
762080d2917SAlp Dener           }
763080d2917SAlp Dener 
764080d2917SAlp Dener           ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr);
765080d2917SAlp Dener           tau_1 = bnk->theta * gdx / (bnk->theta * gdx - (1.0 - bnk->theta) * prered + actred);
766080d2917SAlp Dener           tau_2 = bnk->theta * gdx / (bnk->theta * gdx + (1.0 + bnk->theta) * prered - actred);
767080d2917SAlp Dener           tau_min = PetscMin(tau_1, tau_2);
768080d2917SAlp Dener           tau_max = PetscMax(tau_1, tau_2);
769080d2917SAlp Dener 
770080d2917SAlp Dener           if (kappa >= 1.0 - bnk->mu1) {
771080d2917SAlp Dener             /*  Great agreement */
772080d2917SAlp Dener             *accept = PETSC_TRUE;
773080d2917SAlp Dener             if (tau_max < 1.0) {
774080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm);
775080d2917SAlp Dener             } else if (tau_max > bnk->gamma4) {
776080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, bnk->gamma4 * bnk->dnorm);
777080d2917SAlp Dener             } else {
778080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm);
779080d2917SAlp Dener             }
780080d2917SAlp Dener           } else if (kappa >= 1.0 - bnk->mu2) {
781080d2917SAlp Dener             /*  Good agreement */
782080d2917SAlp Dener             *accept = PETSC_TRUE;
783080d2917SAlp Dener             if (tau_max < bnk->gamma2) {
784080d2917SAlp Dener               tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm);
785080d2917SAlp Dener             } else if (tau_max > bnk->gamma3) {
786080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm);
787080d2917SAlp Dener             } else if (tau_max < 1.0) {
788080d2917SAlp Dener               tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm);
789080d2917SAlp Dener             } else {
790080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm);
791080d2917SAlp Dener             }
792080d2917SAlp Dener           } else {
793080d2917SAlp Dener             /*  Not good agreement */
794080d2917SAlp Dener             if (tau_min > 1.0) {
795080d2917SAlp Dener               tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm);
796080d2917SAlp Dener             } else if (tau_max < bnk->gamma1) {
797080d2917SAlp Dener               tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
798080d2917SAlp Dener             } else if ((tau_min < bnk->gamma1) && (tau_max >= 1.0)) {
799080d2917SAlp Dener               tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
800080d2917SAlp Dener             } else if ((tau_1 >= bnk->gamma1) && (tau_1 < 1.0) && ((tau_2 < bnk->gamma1) || (tau_2 >= 1.0))) {
801080d2917SAlp Dener               tao->trust = tau_1 * PetscMin(tao->trust, bnk->dnorm);
802080d2917SAlp Dener             } else if ((tau_2 >= bnk->gamma1) && (tau_2 < 1.0) && ((tau_1 < bnk->gamma1) || (tau_2 >= 1.0))) {
803080d2917SAlp Dener               tao->trust = tau_2 * PetscMin(tao->trust, bnk->dnorm);
804080d2917SAlp Dener             } else {
805080d2917SAlp Dener               tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm);
806080d2917SAlp Dener             }
807080d2917SAlp Dener           }
808080d2917SAlp Dener         }
809080d2917SAlp Dener       }
810080d2917SAlp Dener     } else {
811080d2917SAlp Dener       /*  Newton step was not good; reduce the radius */
812080d2917SAlp Dener       tao->trust = bnk->gamma1 * PetscMin(bnk->dnorm, tao->trust);
813080d2917SAlp Dener     }
814080d2917SAlp Dener     /*  The radius may have been increased; modify if it is too large */
815080d2917SAlp Dener     tao->trust = PetscMin(tao->trust, bnk->max_radius);
816080d2917SAlp Dener   }
817fed79b8eSAlp Dener   tao->trust = PetscMax(tao->trust, bnk->min_radius);
818080d2917SAlp Dener   PetscFunctionReturn(0);
819080d2917SAlp Dener }
820080d2917SAlp Dener 
821eb910715SAlp Dener /* ---------------------------------------------------------- */
822*df278d8fSAlp Dener 
823eb910715SAlp Dener static PetscErrorCode TaoSetUp_BNK(Tao tao)
824eb910715SAlp Dener {
825eb910715SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
826eb910715SAlp Dener   PetscErrorCode ierr;
827eb910715SAlp Dener 
828eb910715SAlp Dener   PetscFunctionBegin;
829eb910715SAlp Dener   if (!tao->gradient) {ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);}
830eb910715SAlp Dener   if (!tao->stepdirection) {ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr);}
831eb910715SAlp Dener   if (!bnk->W) {ierr = VecDuplicate(tao->solution,&bnk->W);CHKERRQ(ierr);}
832eb910715SAlp Dener   if (!bnk->Xold) {ierr = VecDuplicate(tao->solution,&bnk->Xold);CHKERRQ(ierr);}
833eb910715SAlp Dener   if (!bnk->Gold) {ierr = VecDuplicate(tao->solution,&bnk->Gold);CHKERRQ(ierr);}
83409164190SAlp Dener   if (!bnk->Xwork) {ierr = VecDuplicate(tao->solution,&bnk->Xwork);CHKERRQ(ierr);}
83509164190SAlp Dener   if (!bnk->Gwork) {ierr = VecDuplicate(tao->solution,&bnk->Gwork);CHKERRQ(ierr);}
83609164190SAlp Dener   if (!bnk->unprojected_gradient) {ierr = VecDuplicate(tao->solution,&bnk->unprojected_gradient);CHKERRQ(ierr);}
83709164190SAlp Dener   if (!bnk->unprojected_gradient_old) {ierr = VecDuplicate(tao->solution,&bnk->unprojected_gradient_old);CHKERRQ(ierr);}
838eb910715SAlp Dener   bnk->Diag = 0;
839eb910715SAlp Dener   bnk->M = 0;
84009164190SAlp Dener   bnk->H_inactive = 0;
84109164190SAlp Dener   bnk->Hpre_inactive = 0;
842eb910715SAlp Dener   PetscFunctionReturn(0);
843eb910715SAlp Dener }
844eb910715SAlp Dener 
845eb910715SAlp Dener /*------------------------------------------------------------*/
846*df278d8fSAlp Dener 
847eb910715SAlp Dener static PetscErrorCode TaoDestroy_BNK(Tao tao)
848eb910715SAlp Dener {
849eb910715SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
850eb910715SAlp Dener   PetscErrorCode ierr;
851eb910715SAlp Dener 
852eb910715SAlp Dener   PetscFunctionBegin;
853eb910715SAlp Dener   if (tao->setupcalled) {
854eb910715SAlp Dener     ierr = VecDestroy(&bnk->W);CHKERRQ(ierr);
855eb910715SAlp Dener     ierr = VecDestroy(&bnk->Xold);CHKERRQ(ierr);
856eb910715SAlp Dener     ierr = VecDestroy(&bnk->Gold);CHKERRQ(ierr);
85709164190SAlp Dener     ierr = VecDestroy(&bnk->Xwork);CHKERRQ(ierr);
85809164190SAlp Dener     ierr = VecDestroy(&bnk->Gwork);CHKERRQ(ierr);
85909164190SAlp Dener     ierr = VecDestroy(&bnk->unprojected_gradient);CHKERRQ(ierr);
86009164190SAlp Dener     ierr = VecDestroy(&bnk->unprojected_gradient_old);CHKERRQ(ierr);
861eb910715SAlp Dener   }
862eb910715SAlp Dener   ierr = VecDestroy(&bnk->Diag);CHKERRQ(ierr);
863eb910715SAlp Dener   ierr = MatDestroy(&bnk->M);CHKERRQ(ierr);
864eb910715SAlp Dener   ierr = PetscFree(tao->data);CHKERRQ(ierr);
865eb910715SAlp Dener   PetscFunctionReturn(0);
866eb910715SAlp Dener }
867eb910715SAlp Dener 
868eb910715SAlp Dener /*------------------------------------------------------------*/
869*df278d8fSAlp Dener 
870eb910715SAlp Dener static PetscErrorCode TaoSetFromOptions_BNK(PetscOptionItems *PetscOptionsObject,Tao tao)
871eb910715SAlp Dener {
872eb910715SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
873eb910715SAlp Dener   PetscErrorCode ierr;
874eb910715SAlp Dener 
875eb910715SAlp Dener   PetscFunctionBegin;
876eb910715SAlp Dener   ierr = PetscOptionsHead(PetscOptionsObject,"Newton line search method for unconstrained optimization");CHKERRQ(ierr);
877eb910715SAlp 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);
878eb910715SAlp 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);
879eb910715SAlp 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);
880eb910715SAlp 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);
881eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_sval", "perturbation starting value", "", bnk->sval, &bnk->sval,NULL);CHKERRQ(ierr);
882eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_imin", "minimum initial perturbation", "", bnk->imin, &bnk->imin,NULL);CHKERRQ(ierr);
883eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_imax", "maximum initial perturbation", "", bnk->imax, &bnk->imax,NULL);CHKERRQ(ierr);
884eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_imfac", "initial merit factor", "", bnk->imfac, &bnk->imfac,NULL);CHKERRQ(ierr);
885eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_pmin", "minimum perturbation", "", bnk->pmin, &bnk->pmin,NULL);CHKERRQ(ierr);
886eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_pmax", "maximum perturbation", "", bnk->pmax, &bnk->pmax,NULL);CHKERRQ(ierr);
887eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_pgfac", "growth factor", "", bnk->pgfac, &bnk->pgfac,NULL);CHKERRQ(ierr);
888eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_psfac", "shrink factor", "", bnk->psfac, &bnk->psfac,NULL);CHKERRQ(ierr);
889eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_pmgfac", "merit growth factor", "", bnk->pmgfac, &bnk->pmgfac,NULL);CHKERRQ(ierr);
890eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_pmsfac", "merit shrink factor", "", bnk->pmsfac, &bnk->pmsfac,NULL);CHKERRQ(ierr);
891eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_eta1", "poor steplength; reduce radius", "", bnk->eta1, &bnk->eta1,NULL);CHKERRQ(ierr);
892eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_eta2", "reasonable steplength; leave radius alone", "", bnk->eta2, &bnk->eta2,NULL);CHKERRQ(ierr);
893eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_eta3", "good steplength; increase radius", "", bnk->eta3, &bnk->eta3,NULL);CHKERRQ(ierr);
894eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_eta4", "excellent steplength; greatly increase radius", "", bnk->eta4, &bnk->eta4,NULL);CHKERRQ(ierr);
895eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_alpha1", "", "", bnk->alpha1, &bnk->alpha1,NULL);CHKERRQ(ierr);
896eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_alpha2", "", "", bnk->alpha2, &bnk->alpha2,NULL);CHKERRQ(ierr);
897eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_alpha3", "", "", bnk->alpha3, &bnk->alpha3,NULL);CHKERRQ(ierr);
898eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_alpha4", "", "", bnk->alpha4, &bnk->alpha4,NULL);CHKERRQ(ierr);
899eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_alpha5", "", "", bnk->alpha5, &bnk->alpha5,NULL);CHKERRQ(ierr);
900eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_nu1", "poor steplength; reduce radius", "", bnk->nu1, &bnk->nu1,NULL);CHKERRQ(ierr);
901eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_nu2", "reasonable steplength; leave radius alone", "", bnk->nu2, &bnk->nu2,NULL);CHKERRQ(ierr);
902eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_nu3", "good steplength; increase radius", "", bnk->nu3, &bnk->nu3,NULL);CHKERRQ(ierr);
903eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_nu4", "excellent steplength; greatly increase radius", "", bnk->nu4, &bnk->nu4,NULL);CHKERRQ(ierr);
904eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_omega1", "", "", bnk->omega1, &bnk->omega1,NULL);CHKERRQ(ierr);
905eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_omega2", "", "", bnk->omega2, &bnk->omega2,NULL);CHKERRQ(ierr);
906eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_omega3", "", "", bnk->omega3, &bnk->omega3,NULL);CHKERRQ(ierr);
907eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_omega4", "", "", bnk->omega4, &bnk->omega4,NULL);CHKERRQ(ierr);
908eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_omega5", "", "", bnk->omega5, &bnk->omega5,NULL);CHKERRQ(ierr);
909eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_mu1_i", "", "", bnk->mu1_i, &bnk->mu1_i,NULL);CHKERRQ(ierr);
910eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_mu2_i", "", "", bnk->mu2_i, &bnk->mu2_i,NULL);CHKERRQ(ierr);
911eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_gamma1_i", "", "", bnk->gamma1_i, &bnk->gamma1_i,NULL);CHKERRQ(ierr);
912eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_gamma2_i", "", "", bnk->gamma2_i, &bnk->gamma2_i,NULL);CHKERRQ(ierr);
913eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_gamma3_i", "", "", bnk->gamma3_i, &bnk->gamma3_i,NULL);CHKERRQ(ierr);
914eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_gamma4_i", "", "", bnk->gamma4_i, &bnk->gamma4_i,NULL);CHKERRQ(ierr);
915eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_theta_i", "", "", bnk->theta_i, &bnk->theta_i,NULL);CHKERRQ(ierr);
916eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_mu1", "", "", bnk->mu1, &bnk->mu1,NULL);CHKERRQ(ierr);
917eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_mu2", "", "", bnk->mu2, &bnk->mu2,NULL);CHKERRQ(ierr);
918eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_gamma1", "", "", bnk->gamma1, &bnk->gamma1,NULL);CHKERRQ(ierr);
919eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_gamma2", "", "", bnk->gamma2, &bnk->gamma2,NULL);CHKERRQ(ierr);
920eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_gamma3", "", "", bnk->gamma3, &bnk->gamma3,NULL);CHKERRQ(ierr);
921eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_gamma4", "", "", bnk->gamma4, &bnk->gamma4,NULL);CHKERRQ(ierr);
922eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_theta", "", "", bnk->theta, &bnk->theta,NULL);CHKERRQ(ierr);
923eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_min_radius", "lower bound on initial radius", "", bnk->min_radius, &bnk->min_radius,NULL);CHKERRQ(ierr);
924eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_max_radius", "upper bound on radius", "", bnk->max_radius, &bnk->max_radius,NULL);CHKERRQ(ierr);
925eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_epsilon", "tolerance used when computing actual and predicted reduction", "", bnk->epsilon, &bnk->epsilon,NULL);CHKERRQ(ierr);
926eb910715SAlp Dener   ierr = PetscOptionsTail();CHKERRQ(ierr);
927eb910715SAlp Dener   ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr);
928eb910715SAlp Dener   ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
929eb910715SAlp Dener   PetscFunctionReturn(0);
930eb910715SAlp Dener }
931eb910715SAlp Dener 
932eb910715SAlp Dener /*------------------------------------------------------------*/
933*df278d8fSAlp Dener 
934eb910715SAlp Dener static PetscErrorCode TaoView_BNK(Tao tao, PetscViewer viewer)
935eb910715SAlp Dener {
936eb910715SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
937eb910715SAlp Dener   PetscInt       nrejects;
938eb910715SAlp Dener   PetscBool      isascii;
939eb910715SAlp Dener   PetscErrorCode ierr;
940eb910715SAlp Dener 
941eb910715SAlp Dener   PetscFunctionBegin;
942eb910715SAlp Dener   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
943eb910715SAlp Dener   if (isascii) {
944eb910715SAlp Dener     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
945eb910715SAlp Dener     if (BNK_PC_BFGS == bnk->pc_type && bnk->M) {
946eb910715SAlp Dener       ierr = MatLMVMGetRejects(bnk->M,&nrejects);CHKERRQ(ierr);
947eb910715SAlp Dener       ierr = PetscViewerASCIIPrintf(viewer, "Rejected matrix updates: %D\n",nrejects);CHKERRQ(ierr);
948eb910715SAlp Dener     }
949eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "Newton steps: %D\n", bnk->newt);CHKERRQ(ierr);
950eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "BFGS steps: %D\n", bnk->bfgs);CHKERRQ(ierr);
951eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "Scaled gradient steps: %D\n", bnk->sgrad);CHKERRQ(ierr);
952eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", bnk->grad);CHKERRQ(ierr);
953eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "KSP termination reasons:\n");CHKERRQ(ierr);
954eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  atol: %D\n", bnk->ksp_atol);CHKERRQ(ierr);
955eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  rtol: %D\n", bnk->ksp_rtol);CHKERRQ(ierr);
956eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  ctol: %D\n", bnk->ksp_ctol);CHKERRQ(ierr);
957eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  negc: %D\n", bnk->ksp_negc);CHKERRQ(ierr);
958eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  dtol: %D\n", bnk->ksp_dtol);CHKERRQ(ierr);
959eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  iter: %D\n", bnk->ksp_iter);CHKERRQ(ierr);
960eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  othr: %D\n", bnk->ksp_othr);CHKERRQ(ierr);
961eb910715SAlp Dener     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
962eb910715SAlp Dener   }
963eb910715SAlp Dener   PetscFunctionReturn(0);
964eb910715SAlp Dener }
965eb910715SAlp Dener 
966eb910715SAlp Dener /* ---------------------------------------------------------- */
967*df278d8fSAlp Dener 
968eb910715SAlp Dener /*MC
969eb910715SAlp Dener   TAOBNK - Shared base-type for Bounded Newton-Krylov type algorithms.
97066ed3702SAlp Dener   At each iteration, the BNK methods solve the symmetric
971eb910715SAlp Dener   system of equations to obtain the step diretion dk:
972eb910715SAlp Dener               Hk dk = -gk
9732b97c8d8SAlp Dener   for free variables only. The step can be globalized either through
9742b97c8d8SAlp Dener   trust-region methods, or a line search, or a heuristic mixture of both.
975eb910715SAlp Dener 
976eb910715SAlp Dener     Options Database Keys:
977eb910715SAlp Dener + -tao_BNK_pc_type - "none","ahess","bfgs","petsc"
978eb910715SAlp Dener . -tao_BNK_bfgs_scale_type - "ahess","phess","bfgs"
979eb910715SAlp Dener . -tao_BNK_init_type - "constant","direction","interpolation"
980eb910715SAlp Dener . -tao_BNK_update_type - "step","direction","interpolation"
981eb910715SAlp Dener . -tao_BNK_sval - perturbation starting value
982eb910715SAlp Dener . -tao_BNK_imin - minimum initial perturbation
983eb910715SAlp Dener . -tao_BNK_imax - maximum initial perturbation
984eb910715SAlp Dener . -tao_BNK_pmin - minimum perturbation
985eb910715SAlp Dener . -tao_BNK_pmax - maximum perturbation
986eb910715SAlp Dener . -tao_BNK_pgfac - growth factor
987eb910715SAlp Dener . -tao_BNK_psfac - shrink factor
988eb910715SAlp Dener . -tao_BNK_imfac - initial merit factor
989eb910715SAlp Dener . -tao_BNK_pmgfac - merit growth factor
990eb910715SAlp Dener . -tao_BNK_pmsfac - merit shrink factor
991eb910715SAlp Dener . -tao_BNK_eta1 - poor steplength; reduce radius
992eb910715SAlp Dener . -tao_BNK_eta2 - reasonable steplength; leave radius
993eb910715SAlp Dener . -tao_BNK_eta3 - good steplength; increase readius
994eb910715SAlp Dener . -tao_BNK_eta4 - excellent steplength; greatly increase radius
995eb910715SAlp Dener . -tao_BNK_alpha1 - alpha1 reduction
996eb910715SAlp Dener . -tao_BNK_alpha2 - alpha2 reduction
997eb910715SAlp Dener . -tao_BNK_alpha3 - alpha3 reduction
998eb910715SAlp Dener . -tao_BNK_alpha4 - alpha4 reduction
999eb910715SAlp Dener . -tao_BNK_alpha - alpha5 reduction
1000eb910715SAlp Dener . -tao_BNK_mu1 - mu1 interpolation update
1001eb910715SAlp Dener . -tao_BNK_mu2 - mu2 interpolation update
1002eb910715SAlp Dener . -tao_BNK_gamma1 - gamma1 interpolation update
1003eb910715SAlp Dener . -tao_BNK_gamma2 - gamma2 interpolation update
1004eb910715SAlp Dener . -tao_BNK_gamma3 - gamma3 interpolation update
1005eb910715SAlp Dener . -tao_BNK_gamma4 - gamma4 interpolation update
1006eb910715SAlp Dener . -tao_BNK_theta - theta interpolation update
1007eb910715SAlp Dener . -tao_BNK_omega1 - omega1 step update
1008eb910715SAlp Dener . -tao_BNK_omega2 - omega2 step update
1009eb910715SAlp Dener . -tao_BNK_omega3 - omega3 step update
1010eb910715SAlp Dener . -tao_BNK_omega4 - omega4 step update
1011eb910715SAlp Dener . -tao_BNK_omega5 - omega5 step update
1012eb910715SAlp Dener . -tao_BNK_mu1_i -  mu1 interpolation init factor
1013eb910715SAlp Dener . -tao_BNK_mu2_i -  mu2 interpolation init factor
1014eb910715SAlp Dener . -tao_BNK_gamma1_i -  gamma1 interpolation init factor
1015eb910715SAlp Dener . -tao_BNK_gamma2_i -  gamma2 interpolation init factor
1016eb910715SAlp Dener . -tao_BNK_gamma3_i -  gamma3 interpolation init factor
1017eb910715SAlp Dener . -tao_BNK_gamma4_i -  gamma4 interpolation init factor
1018eb910715SAlp Dener - -tao_BNK_theta_i -  theta interpolation init factor
1019eb910715SAlp Dener 
1020eb910715SAlp Dener   Level: beginner
1021eb910715SAlp Dener M*/
1022eb910715SAlp Dener 
1023eb910715SAlp Dener PetscErrorCode TaoCreate_BNK(Tao tao)
1024eb910715SAlp Dener {
1025eb910715SAlp Dener   TAO_BNK        *bnk;
1026eb910715SAlp Dener   const char     *morethuente_type = TAOLINESEARCHMT;
1027eb910715SAlp Dener   PetscErrorCode ierr;
1028eb910715SAlp Dener 
1029eb910715SAlp Dener   PetscFunctionBegin;
1030eb910715SAlp Dener   ierr = PetscNewLog(tao,&bnk);CHKERRQ(ierr);
1031eb910715SAlp Dener 
1032eb910715SAlp Dener   tao->ops->setup = TaoSetUp_BNK;
1033eb910715SAlp Dener   tao->ops->view = TaoView_BNK;
1034eb910715SAlp Dener   tao->ops->setfromoptions = TaoSetFromOptions_BNK;
1035eb910715SAlp Dener   tao->ops->destroy = TaoDestroy_BNK;
1036eb910715SAlp Dener 
1037eb910715SAlp Dener   /*  Override default settings (unless already changed) */
1038eb910715SAlp Dener   if (!tao->max_it_changed) tao->max_it = 50;
1039eb910715SAlp Dener   if (!tao->trust0_changed) tao->trust0 = 100.0;
1040eb910715SAlp Dener 
1041eb910715SAlp Dener   tao->data = (void*)bnk;
1042eb910715SAlp Dener 
104366ed3702SAlp Dener   /*  Hessian shifting parameters */
1044eb910715SAlp Dener   bnk->sval   = 0.0;
1045eb910715SAlp Dener   bnk->imin   = 1.0e-4;
1046eb910715SAlp Dener   bnk->imax   = 1.0e+2;
1047eb910715SAlp Dener   bnk->imfac  = 1.0e-1;
1048eb910715SAlp Dener 
1049eb910715SAlp Dener   bnk->pmin   = 1.0e-12;
1050eb910715SAlp Dener   bnk->pmax   = 1.0e+2;
1051eb910715SAlp Dener   bnk->pgfac  = 1.0e+1;
1052eb910715SAlp Dener   bnk->psfac  = 4.0e-1;
1053eb910715SAlp Dener   bnk->pmgfac = 1.0e-1;
1054eb910715SAlp Dener   bnk->pmsfac = 1.0e-1;
1055eb910715SAlp Dener 
1056eb910715SAlp Dener   /*  Default values for trust-region radius update based on steplength */
1057eb910715SAlp Dener   bnk->nu1 = 0.25;
1058eb910715SAlp Dener   bnk->nu2 = 0.50;
1059eb910715SAlp Dener   bnk->nu3 = 1.00;
1060eb910715SAlp Dener   bnk->nu4 = 1.25;
1061eb910715SAlp Dener 
1062eb910715SAlp Dener   bnk->omega1 = 0.25;
1063eb910715SAlp Dener   bnk->omega2 = 0.50;
1064eb910715SAlp Dener   bnk->omega3 = 1.00;
1065eb910715SAlp Dener   bnk->omega4 = 2.00;
1066eb910715SAlp Dener   bnk->omega5 = 4.00;
1067eb910715SAlp Dener 
1068eb910715SAlp Dener   /*  Default values for trust-region radius update based on reduction */
1069eb910715SAlp Dener   bnk->eta1 = 1.0e-4;
1070eb910715SAlp Dener   bnk->eta2 = 0.25;
1071eb910715SAlp Dener   bnk->eta3 = 0.50;
1072eb910715SAlp Dener   bnk->eta4 = 0.90;
1073eb910715SAlp Dener 
1074eb910715SAlp Dener   bnk->alpha1 = 0.25;
1075eb910715SAlp Dener   bnk->alpha2 = 0.50;
1076eb910715SAlp Dener   bnk->alpha3 = 1.00;
1077eb910715SAlp Dener   bnk->alpha4 = 2.00;
1078eb910715SAlp Dener   bnk->alpha5 = 4.00;
1079eb910715SAlp Dener 
1080eb910715SAlp Dener   /*  Default values for trust-region radius update based on interpolation */
1081eb910715SAlp Dener   bnk->mu1 = 0.10;
1082eb910715SAlp Dener   bnk->mu2 = 0.50;
1083eb910715SAlp Dener 
1084eb910715SAlp Dener   bnk->gamma1 = 0.25;
1085eb910715SAlp Dener   bnk->gamma2 = 0.50;
1086eb910715SAlp Dener   bnk->gamma3 = 2.00;
1087eb910715SAlp Dener   bnk->gamma4 = 4.00;
1088eb910715SAlp Dener 
1089eb910715SAlp Dener   bnk->theta = 0.05;
1090eb910715SAlp Dener 
1091eb910715SAlp Dener   /*  Default values for trust region initialization based on interpolation */
1092eb910715SAlp Dener   bnk->mu1_i = 0.35;
1093eb910715SAlp Dener   bnk->mu2_i = 0.50;
1094eb910715SAlp Dener 
1095eb910715SAlp Dener   bnk->gamma1_i = 0.0625;
1096eb910715SAlp Dener   bnk->gamma2_i = 0.5;
1097eb910715SAlp Dener   bnk->gamma3_i = 2.0;
1098eb910715SAlp Dener   bnk->gamma4_i = 5.0;
1099eb910715SAlp Dener 
1100eb910715SAlp Dener   bnk->theta_i = 0.25;
1101eb910715SAlp Dener 
1102eb910715SAlp Dener   /*  Remaining parameters */
1103eb910715SAlp Dener   bnk->min_radius = 1.0e-10;
1104eb910715SAlp Dener   bnk->max_radius = 1.0e10;
1105eb910715SAlp Dener   bnk->epsilon = 1.0e-6;
1106eb910715SAlp Dener 
1107eb910715SAlp Dener   bnk->pc_type         = BNK_PC_BFGS;
1108eb910715SAlp Dener   bnk->bfgs_scale_type = BFGS_SCALE_PHESS;
1109eb910715SAlp Dener   bnk->init_type       = BNK_INIT_INTERPOLATION;
1110fed79b8eSAlp Dener   bnk->update_type     = BNK_UPDATE_INTERPOLATION;
1111eb910715SAlp Dener 
1112eb910715SAlp Dener   ierr = TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch);CHKERRQ(ierr);
1113eb910715SAlp Dener   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr);
1114eb910715SAlp Dener   ierr = TaoLineSearchSetType(tao->linesearch,morethuente_type);CHKERRQ(ierr);
1115eb910715SAlp Dener   ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr);
1116eb910715SAlp Dener   ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr);
1117eb910715SAlp Dener 
1118eb910715SAlp Dener   /*  Set linear solver to default for symmetric matrices */
1119eb910715SAlp Dener   ierr = KSPCreate(((PetscObject)tao)->comm,&tao->ksp);CHKERRQ(ierr);
1120eb910715SAlp Dener   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1);CHKERRQ(ierr);
1121eb910715SAlp Dener   ierr = KSPSetOptionsPrefix(tao->ksp,tao->hdr.prefix);CHKERRQ(ierr);
1122eb910715SAlp Dener   ierr = KSPSetType(tao->ksp,KSPCGSTCG);CHKERRQ(ierr);
1123eb910715SAlp Dener   PetscFunctionReturn(0);
1124eb910715SAlp Dener }
1125