xref: /petsc/src/tao/bound/impls/bnk/bnk.c (revision 66ed3702dc8f7c8e3a7cbd4d5287a373dbf44efa)
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 
22eb910715SAlp Dener PetscErrorCode TaoBNKInitialize(Tao tao)
23eb910715SAlp Dener {
24eb910715SAlp Dener   PetscErrorCode               ierr;
25eb910715SAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
26eb910715SAlp Dener   KSPType                      ksp_type;
27eb910715SAlp Dener   PC                           pc;
28eb910715SAlp Dener 
29eb910715SAlp Dener   PetscReal                    fmin, ftrial, prered, actred, kappa, sigma;
30eb910715SAlp Dener   PetscReal                    tau, tau_1, tau_2, tau_max, tau_min, max_radius;
31eb910715SAlp Dener   PetscReal                    delta, step = 1.0;
32eb910715SAlp Dener 
33eb910715SAlp Dener   PetscInt                     n,N,needH = 1;
34eb910715SAlp Dener 
35eb910715SAlp Dener   PetscInt                     i_max = 5;
36eb910715SAlp Dener   PetscInt                     j_max = 1;
37eb910715SAlp Dener   PetscInt                     i, j;
38eb910715SAlp Dener 
39eb910715SAlp Dener   PetscFunctionBegin;
40eb910715SAlp Dener   /* Number of times ksp stopped because of these reasons */
41eb910715SAlp Dener   bnk->ksp_atol = 0;
42eb910715SAlp Dener   bnk->ksp_rtol = 0;
43eb910715SAlp Dener   bnk->ksp_dtol = 0;
44eb910715SAlp Dener   bnk->ksp_ctol = 0;
45eb910715SAlp Dener   bnk->ksp_negc = 0;
46eb910715SAlp Dener   bnk->ksp_iter = 0;
47eb910715SAlp Dener   bnk->ksp_othr = 0;
48eb910715SAlp Dener 
49fed79b8eSAlp Dener   /* Initialize the Hessian perturbation */
50fed79b8eSAlp Dener   bnk->pert = bnk->sval;
51fed79b8eSAlp Dener 
52eb910715SAlp Dener   /* Initialize trust-region radius when using nash, stcg, or gltr
53eb910715SAlp Dener      Command automatically ignored for other methods
54eb910715SAlp Dener      Will be reset during the first iteration
55eb910715SAlp Dener   */
56eb910715SAlp Dener   ierr = KSPGetType(tao->ksp,&ksp_type);CHKERRQ(ierr);
57eb910715SAlp Dener   ierr = PetscStrcmp(ksp_type,KSPCGNASH,&bnk->is_nash);CHKERRQ(ierr);
58eb910715SAlp Dener   ierr = PetscStrcmp(ksp_type,KSPCGSTCG,&bnk->is_stcg);CHKERRQ(ierr);
59eb910715SAlp Dener   ierr = PetscStrcmp(ksp_type,KSPCGGLTR,&bnk->is_gltr);CHKERRQ(ierr);
60eb910715SAlp Dener 
61eb910715SAlp Dener   ierr = KSPCGSetRadius(tao->ksp,bnk->max_radius);CHKERRQ(ierr);
62eb910715SAlp Dener 
63eb910715SAlp Dener   if (bnk->is_nash || bnk->is_stcg || bnk->is_gltr) {
64eb910715SAlp Dener     if (tao->trust0 < 0.0) SETERRQ(PETSC_COMM_SELF,1,"Initial radius negative");
65eb910715SAlp Dener     tao->trust = tao->trust0;
66eb910715SAlp Dener     tao->trust = PetscMax(tao->trust, bnk->min_radius);
67eb910715SAlp Dener     tao->trust = PetscMin(tao->trust, bnk->max_radius);
68eb910715SAlp Dener   }
69eb910715SAlp Dener 
70eb910715SAlp Dener   /* Get vectors we will need */
71eb910715SAlp Dener   if (BNK_PC_BFGS == bnk->pc_type && !bnk->M) {
72eb910715SAlp Dener     ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr);
73eb910715SAlp Dener     ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr);
74eb910715SAlp Dener     ierr = MatCreateLMVM(((PetscObject)tao)->comm,n,N,&bnk->M);CHKERRQ(ierr);
75eb910715SAlp Dener     ierr = MatLMVMAllocateVectors(bnk->M,tao->solution);CHKERRQ(ierr);
76eb910715SAlp Dener   }
77eb910715SAlp Dener 
78eb910715SAlp Dener   /* create vectors for the limited memory preconditioner */
79eb910715SAlp Dener   if ((BNK_PC_BFGS == bnk->pc_type) && (BFGS_SCALE_BFGS != bnk->bfgs_scale_type)) {
80eb910715SAlp Dener     if (!bnk->Diag) {
81eb910715SAlp Dener       ierr = VecDuplicate(tao->solution,&bnk->Diag);CHKERRQ(ierr);
82eb910715SAlp Dener     }
83eb910715SAlp Dener   }
84eb910715SAlp Dener 
85eb910715SAlp Dener   /* Modify the preconditioner to use the bfgs approximation */
86eb910715SAlp Dener   ierr = KSPGetPC(tao->ksp, &pc);CHKERRQ(ierr);
87eb910715SAlp Dener   switch(bnk->pc_type) {
88eb910715SAlp Dener   case BNK_PC_NONE:
89eb910715SAlp Dener     ierr = PCSetType(pc, PCNONE);CHKERRQ(ierr);
90eb910715SAlp Dener     ierr = PCSetFromOptions(pc);CHKERRQ(ierr);
91eb910715SAlp Dener     break;
92eb910715SAlp Dener 
93eb910715SAlp Dener   case BNK_PC_AHESS:
94eb910715SAlp Dener     ierr = PCSetType(pc, PCJACOBI);CHKERRQ(ierr);
95eb910715SAlp Dener     ierr = PCSetFromOptions(pc);CHKERRQ(ierr);
96eb910715SAlp Dener     ierr = PCJacobiSetUseAbs(pc,PETSC_TRUE);CHKERRQ(ierr);
97eb910715SAlp Dener     break;
98eb910715SAlp Dener 
99eb910715SAlp Dener   case BNK_PC_BFGS:
100eb910715SAlp Dener     ierr = PCSetType(pc, PCSHELL);CHKERRQ(ierr);
101eb910715SAlp Dener     ierr = PCSetFromOptions(pc);CHKERRQ(ierr);
102eb910715SAlp Dener     ierr = PCShellSetName(pc, "bfgs");CHKERRQ(ierr);
103eb910715SAlp Dener     ierr = PCShellSetContext(pc, bnk->M);CHKERRQ(ierr);
104eb910715SAlp Dener     ierr = PCShellSetApply(pc, MatLMVMSolveShell);CHKERRQ(ierr);
105eb910715SAlp Dener     break;
106eb910715SAlp Dener 
107eb910715SAlp Dener   default:
108eb910715SAlp Dener     /* Use the pc method set by pc_type */
109eb910715SAlp Dener     break;
110eb910715SAlp Dener   }
111eb910715SAlp Dener 
112eb910715SAlp Dener   /* Initialize trust-region radius.  The initialization is only performed
113eb910715SAlp Dener      when we are using Nash, Steihaug-Toint or the Generalized Lanczos method. */
114eb910715SAlp Dener   if (bnk->is_nash || bnk->is_stcg || bnk->is_gltr) {
115eb910715SAlp Dener     switch(bnk->init_type) {
116eb910715SAlp Dener     case BNK_INIT_CONSTANT:
117eb910715SAlp Dener       /* Use the initial radius specified */
118eb910715SAlp Dener       break;
119eb910715SAlp Dener 
120eb910715SAlp Dener     case BNK_INIT_INTERPOLATION:
121eb910715SAlp Dener       /* Use the initial radius specified */
122eb910715SAlp Dener       max_radius = 0.0;
123eb910715SAlp Dener 
124eb910715SAlp Dener       for (j = 0; j < j_max; ++j) {
125eb910715SAlp Dener         fmin = bnk->f;
126eb910715SAlp Dener         sigma = 0.0;
127eb910715SAlp Dener 
128eb910715SAlp Dener         if (needH) {
129eb910715SAlp Dener           ierr  = TaoComputeHessian(tao, tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr);
130eb910715SAlp Dener           needH = 0;
131eb910715SAlp Dener         }
132eb910715SAlp Dener 
133eb910715SAlp Dener         for (i = 0; i < i_max; ++i) {
134eb910715SAlp Dener           ierr = VecCopy(tao->solution,bnk->W);CHKERRQ(ierr);
135eb910715SAlp Dener           ierr = VecAXPY(bnk->W,-tao->trust/bnk->gnorm,tao->gradient);CHKERRQ(ierr);
136eb910715SAlp Dener           ierr = TaoComputeObjective(tao, bnk->W, &ftrial);CHKERRQ(ierr);
137eb910715SAlp Dener           if (PetscIsInfOrNanReal(ftrial)) {
138eb910715SAlp Dener             tau = bnk->gamma1_i;
139eb910715SAlp Dener           } else {
140eb910715SAlp Dener             if (ftrial < fmin) {
141eb910715SAlp Dener               fmin = ftrial;
142eb910715SAlp Dener               sigma = -tao->trust / bnk->gnorm;
143eb910715SAlp Dener             }
144eb910715SAlp Dener 
145080d2917SAlp Dener             ierr = MatMult(tao->hessian, tao->gradient, tao->stepdirection);CHKERRQ(ierr);
146080d2917SAlp Dener             ierr = VecDot(tao->gradient, tao->stepdirection, &prered);CHKERRQ(ierr);
147eb910715SAlp Dener 
148eb910715SAlp Dener             prered = tao->trust * (bnk->gnorm - 0.5 * tao->trust * prered / (bnk->gnorm * bnk->gnorm));
149eb910715SAlp Dener             actred = bnk->f - ftrial;
150eb910715SAlp Dener             if ((PetscAbsScalar(actred) <= bnk->epsilon) && (PetscAbsScalar(prered) <= bnk->epsilon)) {
151eb910715SAlp Dener               kappa = 1.0;
152eb910715SAlp Dener             } else {
153eb910715SAlp Dener               kappa = actred / prered;
154eb910715SAlp Dener             }
155eb910715SAlp Dener 
156eb910715SAlp Dener             tau_1 = bnk->theta_i * bnk->gnorm * tao->trust / (bnk->theta_i * bnk->gnorm * tao->trust + (1.0 - bnk->theta_i) * prered - actred);
157eb910715SAlp Dener             tau_2 = bnk->theta_i * bnk->gnorm * tao->trust / (bnk->theta_i * bnk->gnorm * tao->trust - (1.0 + bnk->theta_i) * prered + actred);
158eb910715SAlp Dener             tau_min = PetscMin(tau_1, tau_2);
159eb910715SAlp Dener             tau_max = PetscMax(tau_1, tau_2);
160eb910715SAlp Dener 
161eb910715SAlp Dener             if (PetscAbsScalar(kappa - 1.0) <= bnk->mu1_i) {
162eb910715SAlp Dener               /* Great agreement */
163eb910715SAlp Dener               max_radius = PetscMax(max_radius, tao->trust);
164eb910715SAlp Dener 
165eb910715SAlp Dener               if (tau_max < 1.0) {
166eb910715SAlp Dener                 tau = bnk->gamma3_i;
167eb910715SAlp Dener               } else if (tau_max > bnk->gamma4_i) {
168eb910715SAlp Dener                 tau = bnk->gamma4_i;
169eb910715SAlp Dener               } else if (tau_1 >= 1.0 && tau_1 <= bnk->gamma4_i && tau_2 < 1.0) {
170eb910715SAlp Dener                 tau = tau_1;
171eb910715SAlp Dener               } else if (tau_2 >= 1.0 && tau_2 <= bnk->gamma4_i && tau_1 < 1.0) {
172eb910715SAlp Dener                 tau = tau_2;
173eb910715SAlp Dener               } else {
174eb910715SAlp Dener                 tau = tau_max;
175eb910715SAlp Dener               }
176eb910715SAlp Dener             } else if (PetscAbsScalar(kappa - 1.0) <= bnk->mu2_i) {
177eb910715SAlp Dener               /* Good agreement */
178eb910715SAlp Dener               max_radius = PetscMax(max_radius, tao->trust);
179eb910715SAlp Dener 
180eb910715SAlp Dener               if (tau_max < bnk->gamma2_i) {
181eb910715SAlp Dener                 tau = bnk->gamma2_i;
182eb910715SAlp Dener               } else if (tau_max > bnk->gamma3_i) {
183eb910715SAlp Dener                 tau = bnk->gamma3_i;
184eb910715SAlp Dener               } else {
185eb910715SAlp Dener                 tau = tau_max;
186eb910715SAlp Dener               }
187eb910715SAlp Dener             } else {
188eb910715SAlp Dener               /* Not good agreement */
189eb910715SAlp Dener               if (tau_min > 1.0) {
190eb910715SAlp Dener                 tau = bnk->gamma2_i;
191eb910715SAlp Dener               } else if (tau_max < bnk->gamma1_i) {
192eb910715SAlp Dener                 tau = bnk->gamma1_i;
193eb910715SAlp Dener               } else if ((tau_min < bnk->gamma1_i) && (tau_max >= 1.0)) {
194eb910715SAlp Dener                 tau = bnk->gamma1_i;
195eb910715SAlp Dener               } else if ((tau_1 >= bnk->gamma1_i) && (tau_1 < 1.0) && ((tau_2 < bnk->gamma1_i) || (tau_2 >= 1.0))) {
196eb910715SAlp Dener                 tau = tau_1;
197eb910715SAlp Dener               } else if ((tau_2 >= bnk->gamma1_i) && (tau_2 < 1.0) && ((tau_1 < bnk->gamma1_i) || (tau_2 >= 1.0))) {
198eb910715SAlp Dener                 tau = tau_2;
199eb910715SAlp Dener               } else {
200eb910715SAlp Dener                 tau = tau_max;
201eb910715SAlp Dener               }
202eb910715SAlp Dener             }
203eb910715SAlp Dener           }
204eb910715SAlp Dener           tao->trust = tau * tao->trust;
205eb910715SAlp Dener         }
206eb910715SAlp Dener 
207eb910715SAlp Dener         if (fmin < bnk->f) {
208eb910715SAlp Dener           bnk->f = fmin;
209eb910715SAlp Dener           ierr = VecAXPY(tao->solution,sigma,tao->gradient);CHKERRQ(ierr);
21009164190SAlp Dener           ierr = TaoComputeGradient(tao,tao->solution,bnk->unprojected_gradient);CHKERRQ(ierr);
21109164190SAlp Dener           ierr = VecBoundGradientProjection(bnk->unprojected_gradient,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr);
212eb910715SAlp Dener 
213eb910715SAlp Dener           ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&bnk->gnorm);CHKERRQ(ierr);
214eb910715SAlp Dener           if (PetscIsInfOrNanReal(bnk->gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute gradient generated Inf or NaN");
215eb910715SAlp Dener           needH = 1;
216eb910715SAlp Dener 
217eb910715SAlp Dener           ierr = TaoLogConvergenceHistory(tao,bnk->f,bnk->gnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
218eb910715SAlp Dener           ierr = TaoMonitor(tao,tao->niter,bnk->f,bnk->gnorm,0.0,step);CHKERRQ(ierr);
219eb910715SAlp Dener           ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
220eb910715SAlp Dener           if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
221eb910715SAlp Dener         }
222eb910715SAlp Dener       }
223eb910715SAlp Dener       tao->trust = PetscMax(tao->trust, max_radius);
224eb910715SAlp Dener 
225eb910715SAlp Dener       /* Modify the radius if it is too large or small */
226eb910715SAlp Dener       tao->trust = PetscMax(tao->trust, bnk->min_radius);
227eb910715SAlp Dener       tao->trust = PetscMin(tao->trust, bnk->max_radius);
228eb910715SAlp Dener       break;
229eb910715SAlp Dener 
230eb910715SAlp Dener     default:
231eb910715SAlp Dener       /* Norm of the first direction will initialize radius */
232eb910715SAlp Dener       tao->trust = 0.0;
233eb910715SAlp Dener       break;
234eb910715SAlp Dener     }
235eb910715SAlp Dener   }
236eb910715SAlp Dener 
237eb910715SAlp Dener   /* Set initial scaling for the BFGS preconditioner
238eb910715SAlp Dener      This step is done after computing the initial trust-region radius
239eb910715SAlp Dener      since the function value may have decreased */
240eb910715SAlp Dener   if (BNK_PC_BFGS == bnk->pc_type) {
241eb910715SAlp Dener     if (bnk->f != 0.0) {
242eb910715SAlp Dener       delta = 2.0 * PetscAbsScalar(bnk->f) / (bnk->gnorm*bnk->gnorm);
243eb910715SAlp Dener     } else {
244eb910715SAlp Dener       delta = 2.0 / (bnk->gnorm*bnk->gnorm);
245eb910715SAlp Dener     }
246eb910715SAlp Dener     ierr = MatLMVMSetDelta(bnk->M,delta);CHKERRQ(ierr);
247eb910715SAlp Dener   }
248eb910715SAlp Dener 
249eb910715SAlp Dener   /* Set counter for gradient/reset steps*/
250eb910715SAlp Dener   bnk->newt = 0;
251eb910715SAlp Dener   bnk->bfgs = 0;
252eb910715SAlp Dener   bnk->sgrad = 0;
253eb910715SAlp Dener   bnk->grad = 0;
254eb910715SAlp Dener   PetscFunctionReturn(0);
255eb910715SAlp Dener }
256eb910715SAlp Dener 
257fed79b8eSAlp Dener PetscErrorCode TaoBNKComputeStep(Tao tao, PetscBool safeguard, PetscInt *stepType)
258eb910715SAlp Dener {
259eb910715SAlp Dener   PetscErrorCode               ierr;
260eb910715SAlp Dener   TAO_BNK                      *bnk = (TAO_BNK *)tao->data;
261eb910715SAlp Dener   KSPConvergedReason           ksp_reason;
262eb910715SAlp Dener 
263080d2917SAlp Dener   PetscReal                    gdx, delta, e_min;
264eb910715SAlp Dener 
265eb910715SAlp Dener   PetscInt                     bfgsUpdates = 0;
266eb910715SAlp Dener   PetscInt                     kspits;
267eb910715SAlp Dener 
268eb910715SAlp Dener   PetscFunctionBegin;
269eb910715SAlp Dener   /* Shift the Hessian matrix */
270eb910715SAlp Dener   if (bnk->pert > 0) {
271eb910715SAlp Dener     ierr = MatShift(tao->hessian, bnk->pert);CHKERRQ(ierr);
272eb910715SAlp Dener     if (tao->hessian != tao->hessian_pre) {
273eb910715SAlp Dener       ierr = MatShift(tao->hessian_pre, bnk->pert);CHKERRQ(ierr);
274eb910715SAlp Dener     }
275eb910715SAlp Dener   }
276eb910715SAlp Dener 
27709164190SAlp Dener   /* Determine the inactive set */
27809164190SAlp Dener   ierr = ISDestroy(&bnk->inactive_idx);CHKERRQ(ierr);
27909164190SAlp Dener   ierr = VecWhichInactive(tao->XL,tao->solution,bnk->unprojected_gradient,tao->XU,PETSC_TRUE,&bnk->inactive_idx);CHKERRQ(ierr);
28009164190SAlp Dener 
281eb3ba6a7SAlp Dener   /* Prepare masked matrices for the inactive set */
28209164190SAlp Dener   ierr = MatLMVMSetInactive(bnk->M, bnk->inactive_idx);CHKERRQ(ierr);
28309164190SAlp Dener   ierr = VecSet(tao->stepdirection, 0.0);CHKERRQ(ierr);
284eb3ba6a7SAlp Dener   ierr = TaoMatGetSubMat(tao->hessian, bnk->inactive_idx, bnk->Xwork, TAO_SUBSET_MASK, &bnk->H_inactive);CHKERRQ(ierr);
285eb3ba6a7SAlp Dener   if (tao->hessian == tao->hessian_pre) {
286eb3ba6a7SAlp Dener     bnk->Hpre_inactive = bnk->H_inactive;
287eb3ba6a7SAlp Dener   } else {
288eb3ba6a7SAlp Dener     ierr = TaoMatGetSubMat(tao->hessian_pre, bnk->inactive_idx, bnk->Xwork, TAO_SUBSET_MASK, &bnk->Hpre_inactive);CHKERRQ(ierr);
289eb3ba6a7SAlp Dener   }
29009164190SAlp Dener 
291eb910715SAlp Dener   /* Solve the Newton system of equations */
29209164190SAlp Dener   ierr = KSPSetOperators(tao->ksp,bnk->H_inactive,bnk->Hpre_inactive);CHKERRQ(ierr);
293eb910715SAlp Dener   if (bnk->is_nash || bnk->is_stcg || bnk->is_gltr) {
294fed79b8eSAlp Dener     ierr = KSPCGSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr);
295eb3ba6a7SAlp Dener     ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr);
296eb910715SAlp Dener     ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr);
297eb910715SAlp Dener     tao->ksp_its+=kspits;
298eb910715SAlp Dener     tao->ksp_tot_its+=kspits;
299080d2917SAlp Dener     ierr = KSPCGGetNormD(tao->ksp,&bnk->dnorm);CHKERRQ(ierr);
300eb910715SAlp Dener 
301eb910715SAlp Dener     if (0.0 == tao->trust) {
302eb910715SAlp Dener       /* Radius was uninitialized; use the norm of the direction */
303080d2917SAlp Dener       if (bnk->dnorm > 0.0) {
304080d2917SAlp Dener         tao->trust = bnk->dnorm;
305eb910715SAlp Dener 
306eb910715SAlp Dener         /* Modify the radius if it is too large or small */
307eb910715SAlp Dener         tao->trust = PetscMax(tao->trust, bnk->min_radius);
308eb910715SAlp Dener         tao->trust = PetscMin(tao->trust, bnk->max_radius);
309eb910715SAlp Dener       } else {
310eb910715SAlp Dener         /* The direction was bad; set radius to default value and re-solve
311eb910715SAlp Dener            the trust-region subproblem to get a direction */
312eb910715SAlp Dener         tao->trust = tao->trust0;
313eb910715SAlp Dener 
314eb910715SAlp Dener         /* Modify the radius if it is too large or small */
315eb910715SAlp Dener         tao->trust = PetscMax(tao->trust, bnk->min_radius);
316eb910715SAlp Dener         tao->trust = PetscMin(tao->trust, bnk->max_radius);
317eb910715SAlp Dener 
318fed79b8eSAlp Dener         ierr = KSPCGSetRadius(tao->ksp,tao->trust);CHKERRQ(ierr);
319eb3ba6a7SAlp Dener         ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr);
320eb910715SAlp Dener         ierr = KSPGetIterationNumber(tao->ksp,&kspits);CHKERRQ(ierr);
321eb910715SAlp Dener         tao->ksp_its+=kspits;
322eb910715SAlp Dener         tao->ksp_tot_its+=kspits;
323080d2917SAlp Dener         ierr = KSPCGGetNormD(tao->ksp,&bnk->dnorm);CHKERRQ(ierr);
324eb910715SAlp Dener 
325080d2917SAlp Dener         if (bnk->dnorm == 0.0) SETERRQ(PETSC_COMM_SELF,1, "Initial direction zero");
326eb910715SAlp Dener       }
327eb910715SAlp Dener     }
328eb910715SAlp Dener   } else {
329eb3ba6a7SAlp Dener     ierr = KSPSolve(tao->ksp, tao->gradient, tao->stepdirection);CHKERRQ(ierr);
330eb910715SAlp Dener     ierr = KSPGetIterationNumber(tao->ksp, &kspits);CHKERRQ(ierr);
331eb910715SAlp Dener     tao->ksp_its += kspits;
332eb910715SAlp Dener     tao->ksp_tot_its+=kspits;
333eb910715SAlp Dener   }
334fed79b8eSAlp Dener   ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
335fed79b8eSAlp Dener   ierr = KSPGetConvergedReason(tao->ksp, &ksp_reason);CHKERRQ(ierr);
336fed79b8eSAlp Dener 
337eb3ba6a7SAlp Dener   /* Destroy masked matrices */
338eb3ba6a7SAlp Dener   if (bnk->H_inactive != bnk->Hpre_inactive) {
33909164190SAlp Dener     ierr = MatDestroy(&bnk->Hpre_inactive);CHKERRQ(ierr);
340eb3ba6a7SAlp Dener   }
341eb3ba6a7SAlp Dener   ierr = MatDestroy(&bnk->H_inactive);CHKERRQ(ierr);
34209164190SAlp Dener 
343fed79b8eSAlp Dener   /* Make sure the BFGS preconditioner is healthy */
344fed79b8eSAlp Dener   if (bnk->pc_type == BNK_PC_BFGS) {
345fed79b8eSAlp Dener     ierr = MatLMVMGetUpdates(bnk->M, &bfgsUpdates);CHKERRQ(ierr);
346fed79b8eSAlp Dener     if ((KSP_DIVERGED_INDEFINITE_PC == ksp_reason) && (bfgsUpdates > 1)) {
347fed79b8eSAlp Dener       /* Preconditioner is numerically indefinite; reset the approximation. */
348eb910715SAlp Dener       if (bnk->f != 0.0) {
349eb910715SAlp Dener         delta = 2.0 * PetscAbsScalar(bnk->f) / (bnk->gnorm*bnk->gnorm);
350eb910715SAlp Dener       } else {
351eb910715SAlp Dener         delta = 2.0 / (bnk->gnorm*bnk->gnorm);
352eb910715SAlp Dener       }
353eb910715SAlp Dener       ierr = MatLMVMSetDelta(bnk->M,delta);CHKERRQ(ierr);
354eb910715SAlp Dener       ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr);
35509164190SAlp Dener       ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
356eb910715SAlp Dener       bfgsUpdates = 1;
357eb910715SAlp Dener     }
358fed79b8eSAlp Dener   }
359eb910715SAlp Dener 
360eb910715SAlp Dener   if (KSP_CONVERGED_ATOL == ksp_reason) {
361eb910715SAlp Dener     ++bnk->ksp_atol;
362eb910715SAlp Dener   } else if (KSP_CONVERGED_RTOL == ksp_reason) {
363eb910715SAlp Dener     ++bnk->ksp_rtol;
364eb910715SAlp Dener   } else if (KSP_CONVERGED_CG_CONSTRAINED == ksp_reason) {
365eb910715SAlp Dener     ++bnk->ksp_ctol;
366eb910715SAlp Dener   } else if (KSP_CONVERGED_CG_NEG_CURVE == ksp_reason) {
367eb910715SAlp Dener     ++bnk->ksp_negc;
368eb910715SAlp Dener   } else if (KSP_DIVERGED_DTOL == ksp_reason) {
369eb910715SAlp Dener     ++bnk->ksp_dtol;
370eb910715SAlp Dener   } else if (KSP_DIVERGED_ITS == ksp_reason) {
371eb910715SAlp Dener     ++bnk->ksp_iter;
372eb910715SAlp Dener   } else {
373eb910715SAlp Dener     ++bnk->ksp_othr;
374eb910715SAlp Dener   }
375eb910715SAlp Dener 
376eb910715SAlp Dener   /* Check for success (descent direction) */
377fed79b8eSAlp Dener   if (safeguard) {
378080d2917SAlp Dener     ierr = VecDot(tao->stepdirection, tao->gradient, &gdx);CHKERRQ(ierr);
379eb910715SAlp Dener     if ((gdx >= 0.0) || PetscIsInfOrNanReal(gdx)) {
380eb910715SAlp Dener       /* Newton step is not descent or direction produced Inf or NaN
381eb910715SAlp Dener          Update the perturbation for next time */
382eb910715SAlp Dener       if (bnk->pert <= 0.0) {
383eb910715SAlp Dener         /* Initialize the perturbation */
384eb910715SAlp Dener         bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
385eb910715SAlp Dener         if (bnk->is_gltr) {
386eb910715SAlp Dener           ierr = KSPCGGLTRGetMinEig(tao->ksp,&e_min);CHKERRQ(ierr);
387eb910715SAlp Dener           bnk->pert = PetscMax(bnk->pert, -e_min);
388eb910715SAlp Dener         }
389eb910715SAlp Dener       } else {
390eb910715SAlp Dener         /* Increase the perturbation */
391eb910715SAlp Dener         bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
392eb910715SAlp Dener       }
393eb910715SAlp Dener 
394eb910715SAlp Dener       if (BNK_PC_BFGS != bnk->pc_type) {
395eb910715SAlp Dener         /* We don't have the bfgs matrix around and updated
396eb910715SAlp Dener            Must use gradient direction in this case */
397080d2917SAlp Dener         ierr = VecCopy(tao->gradient, tao->stepdirection);CHKERRQ(ierr);
398080d2917SAlp Dener         ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
399eb910715SAlp Dener         ++bnk->grad;
400eb910715SAlp Dener         *stepType = BNK_GRADIENT;
401eb910715SAlp Dener       } else {
402eb910715SAlp Dener         /* Attempt to use the BFGS direction */
40309164190SAlp Dener         ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
40409164190SAlp Dener         ierr = VecBoundGradientProjection(tao->stepdirection,tao->solution,tao->XL,tao->XU,tao->stepdirection);CHKERRQ(ierr);
405080d2917SAlp Dener         ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
406eb910715SAlp Dener 
407eb910715SAlp Dener         /* Check for success (descent direction) */
408080d2917SAlp Dener         ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr);
409eb910715SAlp Dener         if ((gdx >= 0) || PetscIsInfOrNanReal(gdx)) {
410eb910715SAlp Dener           /* BFGS direction is not descent or direction produced not a number
411eb910715SAlp Dener              We can assert bfgsUpdates > 1 in this case because
412eb910715SAlp Dener              the first solve produces the scaled gradient direction,
413eb910715SAlp Dener              which is guaranteed to be descent */
414eb910715SAlp Dener 
415eb910715SAlp Dener           /* Use steepest descent direction (scaled) */
416eb910715SAlp Dener 
417eb910715SAlp Dener           if (bnk->f != 0.0) {
418eb910715SAlp Dener             delta = 2.0 * PetscAbsScalar(bnk->f) / (bnk->gnorm*bnk->gnorm);
419eb910715SAlp Dener           } else {
420eb910715SAlp Dener             delta = 2.0 / (bnk->gnorm*bnk->gnorm);
421eb910715SAlp Dener           }
422eb910715SAlp Dener           ierr = MatLMVMSetDelta(bnk->M, delta);CHKERRQ(ierr);
423eb910715SAlp Dener           ierr = MatLMVMReset(bnk->M);CHKERRQ(ierr);
42409164190SAlp Dener           ierr = MatLMVMUpdate(bnk->M, tao->solution, bnk->unprojected_gradient);CHKERRQ(ierr);
42509164190SAlp Dener           ierr = MatLMVMSolve(bnk->M, bnk->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
42609164190SAlp Dener           ierr = VecBoundGradientProjection(tao->stepdirection,tao->solution,tao->XL,tao->XU,tao->stepdirection);CHKERRQ(ierr);
427080d2917SAlp Dener           ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
428eb910715SAlp Dener 
429eb910715SAlp Dener           bfgsUpdates = 1;
430eb910715SAlp Dener           ++bnk->sgrad;
431eb910715SAlp Dener           *stepType = BNK_SCALED_GRADIENT;
432eb910715SAlp Dener         } else {
433eb910715SAlp Dener           if (1 == bfgsUpdates) {
434eb910715SAlp Dener             /* The first BFGS direction is always the scaled gradient */
435eb910715SAlp Dener             ++bnk->sgrad;
436eb910715SAlp Dener             *stepType = BNK_SCALED_GRADIENT;
437eb910715SAlp Dener           } else {
438eb910715SAlp Dener             ++bnk->bfgs;
439eb910715SAlp Dener             *stepType = BNK_BFGS;
440eb910715SAlp Dener           }
441eb910715SAlp Dener         }
442eb910715SAlp Dener       }
443eb910715SAlp Dener     } else {
444eb910715SAlp Dener       /* Computed Newton step is descent */
445eb910715SAlp Dener       switch (ksp_reason) {
446eb910715SAlp Dener       case KSP_DIVERGED_NANORINF:
447eb910715SAlp Dener       case KSP_DIVERGED_BREAKDOWN:
448eb910715SAlp Dener       case KSP_DIVERGED_INDEFINITE_MAT:
449eb910715SAlp Dener       case KSP_DIVERGED_INDEFINITE_PC:
450eb910715SAlp Dener       case KSP_CONVERGED_CG_NEG_CURVE:
451eb910715SAlp Dener         /* Matrix or preconditioner is indefinite; increase perturbation */
452eb910715SAlp Dener         if (bnk->pert <= 0.0) {
453eb910715SAlp Dener           /* Initialize the perturbation */
454eb910715SAlp Dener           bnk->pert = PetscMin(bnk->imax, PetscMax(bnk->imin, bnk->imfac * bnk->gnorm));
455eb910715SAlp Dener           if (bnk->is_gltr) {
456eb910715SAlp Dener             ierr = KSPCGGLTRGetMinEig(tao->ksp, &e_min);CHKERRQ(ierr);
457eb910715SAlp Dener             bnk->pert = PetscMax(bnk->pert, -e_min);
458eb910715SAlp Dener           }
459eb910715SAlp Dener         } else {
460eb910715SAlp Dener           /* Increase the perturbation */
461eb910715SAlp Dener           bnk->pert = PetscMin(bnk->pmax, PetscMax(bnk->pgfac * bnk->pert, bnk->pmgfac * bnk->gnorm));
462eb910715SAlp Dener         }
463eb910715SAlp Dener         break;
464eb910715SAlp Dener 
465eb910715SAlp Dener       default:
466eb910715SAlp Dener         /* Newton step computation is good; decrease perturbation */
467eb910715SAlp Dener         bnk->pert = PetscMin(bnk->psfac * bnk->pert, bnk->pmsfac * bnk->gnorm);
468eb910715SAlp Dener         if (bnk->pert < bnk->pmin) {
469eb910715SAlp Dener           bnk->pert = 0.0;
470eb910715SAlp Dener         }
471eb910715SAlp Dener         break;
472eb910715SAlp Dener       }
473eb910715SAlp Dener 
474eb910715SAlp Dener       ++bnk->newt;
475fed79b8eSAlp Dener       *stepType = BNK_NEWTON;
476fed79b8eSAlp Dener     }
477fed79b8eSAlp Dener   } else {
478fed79b8eSAlp Dener     ++bnk->newt;
479fed79b8eSAlp Dener     bnk->pert = 0.0;
480fed79b8eSAlp Dener     *stepType = BNK_NEWTON;
481eb910715SAlp Dener   }
482eb910715SAlp Dener   PetscFunctionReturn(0);
483eb910715SAlp Dener }
484eb910715SAlp Dener 
485080d2917SAlp Dener PetscErrorCode TaoBNKUpdateTrustRadius(Tao tao, PetscReal fold, PetscReal fnew, PetscInt stepType, PetscBool *accept)
486080d2917SAlp Dener {
487080d2917SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
488080d2917SAlp Dener   PetscErrorCode ierr;
489080d2917SAlp Dener 
490080d2917SAlp Dener   PetscReal      step, prered, actred, kappa;
491080d2917SAlp Dener   PetscReal      gdx, tau_1, tau_2, tau_min, tau_max;
492080d2917SAlp Dener 
493080d2917SAlp Dener   PetscFunctionBegin;
494080d2917SAlp Dener   /* Update trust region radius */
495080d2917SAlp Dener   *accept = PETSC_FALSE;
496080d2917SAlp Dener   switch(bnk->update_type) {
497080d2917SAlp Dener   case BNK_UPDATE_STEP:
498080d2917SAlp Dener     if (stepType == BNK_NEWTON) {
499080d2917SAlp Dener       ierr = TaoLineSearchGetStepLength(tao->linesearch, &step);CHKERRQ(ierr);
500080d2917SAlp Dener       if (step < bnk->nu1) {
501080d2917SAlp Dener         /* Very bad step taken; reduce radius */
502080d2917SAlp Dener         tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust);
503080d2917SAlp Dener       } else if (step < bnk->nu2) {
504080d2917SAlp Dener         /* Reasonably bad step taken; reduce radius */
505080d2917SAlp Dener         tao->trust = bnk->omega2 * PetscMin(bnk->dnorm, tao->trust);
506080d2917SAlp Dener       } else if (step < bnk->nu3) {
507080d2917SAlp Dener         /*  Reasonable step was taken; leave radius alone */
508080d2917SAlp Dener         if (bnk->omega3 < 1.0) {
509080d2917SAlp Dener           tao->trust = bnk->omega3 * PetscMin(bnk->dnorm, tao->trust);
510080d2917SAlp Dener         } else if (bnk->omega3 > 1.0) {
511080d2917SAlp Dener           tao->trust = PetscMax(bnk->omega3 * bnk->dnorm, tao->trust);
512080d2917SAlp Dener         }
513080d2917SAlp Dener       } else if (step < bnk->nu4) {
514080d2917SAlp Dener         /*  Full step taken; increase the radius */
515080d2917SAlp Dener         tao->trust = PetscMax(bnk->omega4 * bnk->dnorm, tao->trust);
516080d2917SAlp Dener       } else {
517080d2917SAlp Dener         /*  More than full step taken; increase the radius */
518080d2917SAlp Dener         tao->trust = PetscMax(bnk->omega5 * bnk->dnorm, tao->trust);
519080d2917SAlp Dener       }
520080d2917SAlp Dener     } else {
521080d2917SAlp Dener       /*  Newton step was not good; reduce the radius */
522080d2917SAlp Dener       tao->trust = bnk->omega1 * PetscMin(bnk->dnorm, tao->trust);
523080d2917SAlp Dener     }
524080d2917SAlp Dener     break;
525080d2917SAlp Dener 
526080d2917SAlp Dener   case BNK_UPDATE_REDUCTION:
527080d2917SAlp Dener     if (stepType == BNK_NEWTON) {
528080d2917SAlp Dener       /* Get predicted reduction */
529080d2917SAlp Dener       ierr = KSPCGGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr);
530080d2917SAlp Dener       if (prered >= 0.0) {
531fed79b8eSAlp Dener         /* The predicted reduction has the wrong sign.  This cannot
532fed79b8eSAlp Dener            happen in infinite precision arithmetic.  Step should
533fed79b8eSAlp Dener            be rejected! */
534080d2917SAlp Dener         tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
535fed79b8eSAlp Dener       }
536fed79b8eSAlp Dener       else {
537080d2917SAlp Dener         if (PetscIsInfOrNanReal(fnew)) {
538080d2917SAlp Dener           tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
539080d2917SAlp Dener         } else {
540080d2917SAlp Dener           /* Compute and actual reduction */
541080d2917SAlp Dener           actred = fold - fnew;
542080d2917SAlp Dener           prered = -prered;
543080d2917SAlp Dener           if ((PetscAbsScalar(actred) <= bnk->epsilon) &&
544080d2917SAlp Dener               (PetscAbsScalar(prered) <= bnk->epsilon)) {
545080d2917SAlp Dener             kappa = 1.0;
546fed79b8eSAlp Dener           }
547fed79b8eSAlp Dener           else {
548080d2917SAlp Dener             kappa = actred / prered;
549080d2917SAlp Dener           }
550fed79b8eSAlp Dener 
551fed79b8eSAlp Dener           /* Accept or reject the step and update radius */
552080d2917SAlp Dener           if (kappa < bnk->eta1) {
553fed79b8eSAlp Dener             /* Reject the step */
554080d2917SAlp Dener             tao->trust = bnk->alpha1 * PetscMin(tao->trust, bnk->dnorm);
555fed79b8eSAlp Dener           }
556fed79b8eSAlp Dener           else {
557fed79b8eSAlp Dener             /* Accept the step */
558080d2917SAlp Dener             if (kappa < bnk->eta2) {
559080d2917SAlp Dener               /* Marginal bad step */
560080d2917SAlp Dener               tao->trust = bnk->alpha2 * PetscMin(tao->trust, bnk->dnorm);
561080d2917SAlp Dener             }
562fed79b8eSAlp Dener             else if (kappa < bnk->eta3) {
563fed79b8eSAlp Dener               /* Reasonable step */
564fed79b8eSAlp Dener               tao->trust = bnk->alpha3 * tao->trust;
565fed79b8eSAlp Dener             }
566fed79b8eSAlp Dener             else if (kappa < bnk->eta4) {
567080d2917SAlp Dener               /* Good step */
568080d2917SAlp Dener               tao->trust = PetscMax(bnk->alpha4 * bnk->dnorm, tao->trust);
569fed79b8eSAlp Dener             }
570fed79b8eSAlp Dener             else {
571080d2917SAlp Dener               /* Very good step */
572080d2917SAlp Dener               tao->trust = PetscMax(bnk->alpha5 * bnk->dnorm, tao->trust);
573080d2917SAlp Dener             }
574fed79b8eSAlp Dener             *accept = PETSC_TRUE;
575080d2917SAlp Dener           }
576080d2917SAlp Dener         }
577080d2917SAlp Dener       }
578080d2917SAlp Dener     } else {
579080d2917SAlp Dener       /*  Newton step was not good; reduce the radius */
580080d2917SAlp Dener       tao->trust = bnk->alpha1 * PetscMin(bnk->dnorm, tao->trust);
581080d2917SAlp Dener     }
582080d2917SAlp Dener     break;
583080d2917SAlp Dener 
584080d2917SAlp Dener   default:
585080d2917SAlp Dener     if (stepType == BNK_NEWTON) {
586080d2917SAlp Dener       ierr = KSPCGGetObjFcn(tao->ksp,&prered);CHKERRQ(ierr);
587080d2917SAlp Dener       if (prered >= 0.0) {
588080d2917SAlp Dener         /*  The predicted reduction has the wrong sign.  This cannot */
589080d2917SAlp Dener         /*  happen in infinite precision arithmetic.  Step should */
590080d2917SAlp Dener         /*  be rejected! */
591080d2917SAlp Dener         tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
592080d2917SAlp Dener       } else {
593080d2917SAlp Dener         if (PetscIsInfOrNanReal(fnew)) {
594080d2917SAlp Dener           tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
595080d2917SAlp Dener         } else {
596080d2917SAlp Dener           actred = fold - fnew;
597080d2917SAlp Dener           prered = -prered;
598080d2917SAlp Dener           if ((PetscAbsScalar(actred) <= bnk->epsilon) && (PetscAbsScalar(prered) <= bnk->epsilon)) {
599080d2917SAlp Dener             kappa = 1.0;
600080d2917SAlp Dener           } else {
601080d2917SAlp Dener             kappa = actred / prered;
602080d2917SAlp Dener           }
603080d2917SAlp Dener 
604080d2917SAlp Dener           ierr = VecDot(tao->gradient, tao->stepdirection, &gdx);CHKERRQ(ierr);
605080d2917SAlp Dener           tau_1 = bnk->theta * gdx / (bnk->theta * gdx - (1.0 - bnk->theta) * prered + actred);
606080d2917SAlp Dener           tau_2 = bnk->theta * gdx / (bnk->theta * gdx + (1.0 + bnk->theta) * prered - actred);
607080d2917SAlp Dener           tau_min = PetscMin(tau_1, tau_2);
608080d2917SAlp Dener           tau_max = PetscMax(tau_1, tau_2);
609080d2917SAlp Dener 
610080d2917SAlp Dener           if (kappa >= 1.0 - bnk->mu1) {
611080d2917SAlp Dener             /*  Great agreement */
612080d2917SAlp Dener             *accept = PETSC_TRUE;
613080d2917SAlp Dener             if (tau_max < 1.0) {
614080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm);
615080d2917SAlp Dener             } else if (tau_max > bnk->gamma4) {
616080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, bnk->gamma4 * bnk->dnorm);
617080d2917SAlp Dener             } else {
618080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm);
619080d2917SAlp Dener             }
620080d2917SAlp Dener           } else if (kappa >= 1.0 - bnk->mu2) {
621080d2917SAlp Dener             /*  Good agreement */
622080d2917SAlp Dener             *accept = PETSC_TRUE;
623080d2917SAlp Dener             if (tau_max < bnk->gamma2) {
624080d2917SAlp Dener               tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm);
625080d2917SAlp Dener             } else if (tau_max > bnk->gamma3) {
626080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, bnk->gamma3 * bnk->dnorm);
627080d2917SAlp Dener             } else if (tau_max < 1.0) {
628080d2917SAlp Dener               tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm);
629080d2917SAlp Dener             } else {
630080d2917SAlp Dener               tao->trust = PetscMax(tao->trust, tau_max * bnk->dnorm);
631080d2917SAlp Dener             }
632080d2917SAlp Dener           } else {
633080d2917SAlp Dener             /*  Not good agreement */
634080d2917SAlp Dener             if (tau_min > 1.0) {
635080d2917SAlp Dener               tao->trust = bnk->gamma2 * PetscMin(tao->trust, bnk->dnorm);
636080d2917SAlp Dener             } else if (tau_max < bnk->gamma1) {
637080d2917SAlp Dener               tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
638080d2917SAlp Dener             } else if ((tau_min < bnk->gamma1) && (tau_max >= 1.0)) {
639080d2917SAlp Dener               tao->trust = bnk->gamma1 * PetscMin(tao->trust, bnk->dnorm);
640080d2917SAlp Dener             } else if ((tau_1 >= bnk->gamma1) && (tau_1 < 1.0) && ((tau_2 < bnk->gamma1) || (tau_2 >= 1.0))) {
641080d2917SAlp Dener               tao->trust = tau_1 * PetscMin(tao->trust, bnk->dnorm);
642080d2917SAlp Dener             } else if ((tau_2 >= bnk->gamma1) && (tau_2 < 1.0) && ((tau_1 < bnk->gamma1) || (tau_2 >= 1.0))) {
643080d2917SAlp Dener               tao->trust = tau_2 * PetscMin(tao->trust, bnk->dnorm);
644080d2917SAlp Dener             } else {
645080d2917SAlp Dener               tao->trust = tau_max * PetscMin(tao->trust, bnk->dnorm);
646080d2917SAlp Dener             }
647080d2917SAlp Dener           }
648080d2917SAlp Dener         }
649080d2917SAlp Dener       }
650080d2917SAlp Dener     } else {
651080d2917SAlp Dener       /*  Newton step was not good; reduce the radius */
652080d2917SAlp Dener       tao->trust = bnk->gamma1 * PetscMin(bnk->dnorm, tao->trust);
653080d2917SAlp Dener     }
654080d2917SAlp Dener     /*  The radius may have been increased; modify if it is too large */
655080d2917SAlp Dener     tao->trust = PetscMin(tao->trust, bnk->max_radius);
656080d2917SAlp Dener   }
657fed79b8eSAlp Dener   tao->trust = PetscMax(tao->trust, bnk->min_radius);
658080d2917SAlp Dener   PetscFunctionReturn(0);
659080d2917SAlp Dener }
660080d2917SAlp Dener 
661eb910715SAlp Dener /* ---------------------------------------------------------- */
662eb910715SAlp Dener static PetscErrorCode TaoSetUp_BNK(Tao tao)
663eb910715SAlp Dener {
664eb910715SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
665eb910715SAlp Dener   PetscErrorCode ierr;
666eb910715SAlp Dener 
667eb910715SAlp Dener   PetscFunctionBegin;
668eb910715SAlp Dener   if (!tao->gradient) {ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);}
669eb910715SAlp Dener   if (!tao->stepdirection) {ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr);}
670eb910715SAlp Dener   if (!bnk->W) {ierr = VecDuplicate(tao->solution,&bnk->W);CHKERRQ(ierr);}
671eb910715SAlp Dener   if (!bnk->Xold) {ierr = VecDuplicate(tao->solution,&bnk->Xold);CHKERRQ(ierr);}
672eb910715SAlp Dener   if (!bnk->Gold) {ierr = VecDuplicate(tao->solution,&bnk->Gold);CHKERRQ(ierr);}
67309164190SAlp Dener   if (!bnk->Xwork) {ierr = VecDuplicate(tao->solution,&bnk->Xwork);CHKERRQ(ierr);}
67409164190SAlp Dener   if (!bnk->Gwork) {ierr = VecDuplicate(tao->solution,&bnk->Gwork);CHKERRQ(ierr);}
67509164190SAlp Dener   if (!bnk->unprojected_gradient) {ierr = VecDuplicate(tao->solution,&bnk->unprojected_gradient);CHKERRQ(ierr);}
67609164190SAlp Dener   if (!bnk->unprojected_gradient_old) {ierr = VecDuplicate(tao->solution,&bnk->unprojected_gradient_old);CHKERRQ(ierr);}
677eb910715SAlp Dener   bnk->Diag = 0;
678eb910715SAlp Dener   bnk->M = 0;
67909164190SAlp Dener   bnk->H_inactive = 0;
68009164190SAlp Dener   bnk->Hpre_inactive = 0;
681eb910715SAlp Dener   PetscFunctionReturn(0);
682eb910715SAlp Dener }
683eb910715SAlp Dener 
684eb910715SAlp Dener /*------------------------------------------------------------*/
685eb910715SAlp Dener static PetscErrorCode TaoDestroy_BNK(Tao tao)
686eb910715SAlp Dener {
687eb910715SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
688eb910715SAlp Dener   PetscErrorCode ierr;
689eb910715SAlp Dener 
690eb910715SAlp Dener   PetscFunctionBegin;
691eb910715SAlp Dener   if (tao->setupcalled) {
692eb910715SAlp Dener     ierr = VecDestroy(&bnk->W);CHKERRQ(ierr);
693eb910715SAlp Dener     ierr = VecDestroy(&bnk->Xold);CHKERRQ(ierr);
694eb910715SAlp Dener     ierr = VecDestroy(&bnk->Gold);CHKERRQ(ierr);
69509164190SAlp Dener     ierr = VecDestroy(&bnk->Xwork);CHKERRQ(ierr);
69609164190SAlp Dener     ierr = VecDestroy(&bnk->Gwork);CHKERRQ(ierr);
69709164190SAlp Dener     ierr = VecDestroy(&bnk->unprojected_gradient);CHKERRQ(ierr);
69809164190SAlp Dener     ierr = VecDestroy(&bnk->unprojected_gradient_old);CHKERRQ(ierr);
699eb910715SAlp Dener   }
700eb910715SAlp Dener   ierr = VecDestroy(&bnk->Diag);CHKERRQ(ierr);
701eb910715SAlp Dener   ierr = MatDestroy(&bnk->M);CHKERRQ(ierr);
702eb910715SAlp Dener   ierr = PetscFree(tao->data);CHKERRQ(ierr);
703eb910715SAlp Dener   PetscFunctionReturn(0);
704eb910715SAlp Dener }
705eb910715SAlp Dener 
706eb910715SAlp Dener /*------------------------------------------------------------*/
707eb910715SAlp Dener static PetscErrorCode TaoSetFromOptions_BNK(PetscOptionItems *PetscOptionsObject,Tao tao)
708eb910715SAlp Dener {
709eb910715SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
710eb910715SAlp Dener   PetscErrorCode ierr;
711eb910715SAlp Dener 
712eb910715SAlp Dener   PetscFunctionBegin;
713eb910715SAlp Dener   ierr = PetscOptionsHead(PetscOptionsObject,"Newton line search method for unconstrained optimization");CHKERRQ(ierr);
714eb910715SAlp 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);
715eb910715SAlp 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);
716eb910715SAlp 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);
717eb910715SAlp 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);
718eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_sval", "perturbation starting value", "", bnk->sval, &bnk->sval,NULL);CHKERRQ(ierr);
719eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_imin", "minimum initial perturbation", "", bnk->imin, &bnk->imin,NULL);CHKERRQ(ierr);
720eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_imax", "maximum initial perturbation", "", bnk->imax, &bnk->imax,NULL);CHKERRQ(ierr);
721eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_imfac", "initial merit factor", "", bnk->imfac, &bnk->imfac,NULL);CHKERRQ(ierr);
722eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_pmin", "minimum perturbation", "", bnk->pmin, &bnk->pmin,NULL);CHKERRQ(ierr);
723eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_pmax", "maximum perturbation", "", bnk->pmax, &bnk->pmax,NULL);CHKERRQ(ierr);
724eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_pgfac", "growth factor", "", bnk->pgfac, &bnk->pgfac,NULL);CHKERRQ(ierr);
725eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_psfac", "shrink factor", "", bnk->psfac, &bnk->psfac,NULL);CHKERRQ(ierr);
726eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_pmgfac", "merit growth factor", "", bnk->pmgfac, &bnk->pmgfac,NULL);CHKERRQ(ierr);
727eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_pmsfac", "merit shrink factor", "", bnk->pmsfac, &bnk->pmsfac,NULL);CHKERRQ(ierr);
728eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_eta1", "poor steplength; reduce radius", "", bnk->eta1, &bnk->eta1,NULL);CHKERRQ(ierr);
729eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_eta2", "reasonable steplength; leave radius alone", "", bnk->eta2, &bnk->eta2,NULL);CHKERRQ(ierr);
730eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_eta3", "good steplength; increase radius", "", bnk->eta3, &bnk->eta3,NULL);CHKERRQ(ierr);
731eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_eta4", "excellent steplength; greatly increase radius", "", bnk->eta4, &bnk->eta4,NULL);CHKERRQ(ierr);
732eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_alpha1", "", "", bnk->alpha1, &bnk->alpha1,NULL);CHKERRQ(ierr);
733eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_alpha2", "", "", bnk->alpha2, &bnk->alpha2,NULL);CHKERRQ(ierr);
734eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_alpha3", "", "", bnk->alpha3, &bnk->alpha3,NULL);CHKERRQ(ierr);
735eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_alpha4", "", "", bnk->alpha4, &bnk->alpha4,NULL);CHKERRQ(ierr);
736eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_alpha5", "", "", bnk->alpha5, &bnk->alpha5,NULL);CHKERRQ(ierr);
737eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_nu1", "poor steplength; reduce radius", "", bnk->nu1, &bnk->nu1,NULL);CHKERRQ(ierr);
738eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_nu2", "reasonable steplength; leave radius alone", "", bnk->nu2, &bnk->nu2,NULL);CHKERRQ(ierr);
739eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_nu3", "good steplength; increase radius", "", bnk->nu3, &bnk->nu3,NULL);CHKERRQ(ierr);
740eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_nu4", "excellent steplength; greatly increase radius", "", bnk->nu4, &bnk->nu4,NULL);CHKERRQ(ierr);
741eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_omega1", "", "", bnk->omega1, &bnk->omega1,NULL);CHKERRQ(ierr);
742eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_omega2", "", "", bnk->omega2, &bnk->omega2,NULL);CHKERRQ(ierr);
743eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_omega3", "", "", bnk->omega3, &bnk->omega3,NULL);CHKERRQ(ierr);
744eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_omega4", "", "", bnk->omega4, &bnk->omega4,NULL);CHKERRQ(ierr);
745eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_omega5", "", "", bnk->omega5, &bnk->omega5,NULL);CHKERRQ(ierr);
746eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_mu1_i", "", "", bnk->mu1_i, &bnk->mu1_i,NULL);CHKERRQ(ierr);
747eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_mu2_i", "", "", bnk->mu2_i, &bnk->mu2_i,NULL);CHKERRQ(ierr);
748eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_gamma1_i", "", "", bnk->gamma1_i, &bnk->gamma1_i,NULL);CHKERRQ(ierr);
749eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_gamma2_i", "", "", bnk->gamma2_i, &bnk->gamma2_i,NULL);CHKERRQ(ierr);
750eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_gamma3_i", "", "", bnk->gamma3_i, &bnk->gamma3_i,NULL);CHKERRQ(ierr);
751eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_gamma4_i", "", "", bnk->gamma4_i, &bnk->gamma4_i,NULL);CHKERRQ(ierr);
752eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_theta_i", "", "", bnk->theta_i, &bnk->theta_i,NULL);CHKERRQ(ierr);
753eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_mu1", "", "", bnk->mu1, &bnk->mu1,NULL);CHKERRQ(ierr);
754eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_mu2", "", "", bnk->mu2, &bnk->mu2,NULL);CHKERRQ(ierr);
755eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_gamma1", "", "", bnk->gamma1, &bnk->gamma1,NULL);CHKERRQ(ierr);
756eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_gamma2", "", "", bnk->gamma2, &bnk->gamma2,NULL);CHKERRQ(ierr);
757eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_gamma3", "", "", bnk->gamma3, &bnk->gamma3,NULL);CHKERRQ(ierr);
758eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_gamma4", "", "", bnk->gamma4, &bnk->gamma4,NULL);CHKERRQ(ierr);
759eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_theta", "", "", bnk->theta, &bnk->theta,NULL);CHKERRQ(ierr);
760eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_min_radius", "lower bound on initial radius", "", bnk->min_radius, &bnk->min_radius,NULL);CHKERRQ(ierr);
761eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_max_radius", "upper bound on radius", "", bnk->max_radius, &bnk->max_radius,NULL);CHKERRQ(ierr);
762eb910715SAlp Dener   ierr = PetscOptionsReal("-tao_BNK_epsilon", "tolerance used when computing actual and predicted reduction", "", bnk->epsilon, &bnk->epsilon,NULL);CHKERRQ(ierr);
763eb910715SAlp Dener   ierr = PetscOptionsTail();CHKERRQ(ierr);
764eb910715SAlp Dener   ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr);
765eb910715SAlp Dener   ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
766eb910715SAlp Dener   PetscFunctionReturn(0);
767eb910715SAlp Dener }
768eb910715SAlp Dener 
769eb910715SAlp Dener /*------------------------------------------------------------*/
770eb910715SAlp Dener static PetscErrorCode TaoView_BNK(Tao tao, PetscViewer viewer)
771eb910715SAlp Dener {
772eb910715SAlp Dener   TAO_BNK        *bnk = (TAO_BNK *)tao->data;
773eb910715SAlp Dener   PetscInt       nrejects;
774eb910715SAlp Dener   PetscBool      isascii;
775eb910715SAlp Dener   PetscErrorCode ierr;
776eb910715SAlp Dener 
777eb910715SAlp Dener   PetscFunctionBegin;
778eb910715SAlp Dener   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
779eb910715SAlp Dener   if (isascii) {
780eb910715SAlp Dener     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
781eb910715SAlp Dener     if (BNK_PC_BFGS == bnk->pc_type && bnk->M) {
782eb910715SAlp Dener       ierr = MatLMVMGetRejects(bnk->M,&nrejects);CHKERRQ(ierr);
783eb910715SAlp Dener       ierr = PetscViewerASCIIPrintf(viewer, "Rejected matrix updates: %D\n",nrejects);CHKERRQ(ierr);
784eb910715SAlp Dener     }
785eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "Newton steps: %D\n", bnk->newt);CHKERRQ(ierr);
786eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "BFGS steps: %D\n", bnk->bfgs);CHKERRQ(ierr);
787eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "Scaled gradient steps: %D\n", bnk->sgrad);CHKERRQ(ierr);
788eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", bnk->grad);CHKERRQ(ierr);
789eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "KSP termination reasons:\n");CHKERRQ(ierr);
790eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  atol: %D\n", bnk->ksp_atol);CHKERRQ(ierr);
791eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  rtol: %D\n", bnk->ksp_rtol);CHKERRQ(ierr);
792eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  ctol: %D\n", bnk->ksp_ctol);CHKERRQ(ierr);
793eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  negc: %D\n", bnk->ksp_negc);CHKERRQ(ierr);
794eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  dtol: %D\n", bnk->ksp_dtol);CHKERRQ(ierr);
795eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  iter: %D\n", bnk->ksp_iter);CHKERRQ(ierr);
796eb910715SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "  othr: %D\n", bnk->ksp_othr);CHKERRQ(ierr);
797eb910715SAlp Dener     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
798eb910715SAlp Dener   }
799eb910715SAlp Dener   PetscFunctionReturn(0);
800eb910715SAlp Dener }
801eb910715SAlp Dener 
802eb910715SAlp Dener /* ---------------------------------------------------------- */
803eb910715SAlp Dener /*MC
804eb910715SAlp Dener   TAOBNK - Shared base-type for Bounded Newton-Krylov type algorithms.
805*66ed3702SAlp Dener   At each iteration, the BNK methods solve the symmetric
806eb910715SAlp Dener   system of equations to obtain the step diretion dk:
807eb910715SAlp Dener               Hk dk = -gk
8082b97c8d8SAlp Dener   for free variables only. The step can be globalized either through
8092b97c8d8SAlp Dener   trust-region methods, or a line search, or a heuristic mixture of both.
810eb910715SAlp Dener 
811eb910715SAlp Dener     Options Database Keys:
812eb910715SAlp Dener + -tao_BNK_pc_type - "none","ahess","bfgs","petsc"
813eb910715SAlp Dener . -tao_BNK_bfgs_scale_type - "ahess","phess","bfgs"
814eb910715SAlp Dener . -tao_BNK_init_type - "constant","direction","interpolation"
815eb910715SAlp Dener . -tao_BNK_update_type - "step","direction","interpolation"
816eb910715SAlp Dener . -tao_BNK_sval - perturbation starting value
817eb910715SAlp Dener . -tao_BNK_imin - minimum initial perturbation
818eb910715SAlp Dener . -tao_BNK_imax - maximum initial perturbation
819eb910715SAlp Dener . -tao_BNK_pmin - minimum perturbation
820eb910715SAlp Dener . -tao_BNK_pmax - maximum perturbation
821eb910715SAlp Dener . -tao_BNK_pgfac - growth factor
822eb910715SAlp Dener . -tao_BNK_psfac - shrink factor
823eb910715SAlp Dener . -tao_BNK_imfac - initial merit factor
824eb910715SAlp Dener . -tao_BNK_pmgfac - merit growth factor
825eb910715SAlp Dener . -tao_BNK_pmsfac - merit shrink factor
826eb910715SAlp Dener . -tao_BNK_eta1 - poor steplength; reduce radius
827eb910715SAlp Dener . -tao_BNK_eta2 - reasonable steplength; leave radius
828eb910715SAlp Dener . -tao_BNK_eta3 - good steplength; increase readius
829eb910715SAlp Dener . -tao_BNK_eta4 - excellent steplength; greatly increase radius
830eb910715SAlp Dener . -tao_BNK_alpha1 - alpha1 reduction
831eb910715SAlp Dener . -tao_BNK_alpha2 - alpha2 reduction
832eb910715SAlp Dener . -tao_BNK_alpha3 - alpha3 reduction
833eb910715SAlp Dener . -tao_BNK_alpha4 - alpha4 reduction
834eb910715SAlp Dener . -tao_BNK_alpha - alpha5 reduction
835eb910715SAlp Dener . -tao_BNK_mu1 - mu1 interpolation update
836eb910715SAlp Dener . -tao_BNK_mu2 - mu2 interpolation update
837eb910715SAlp Dener . -tao_BNK_gamma1 - gamma1 interpolation update
838eb910715SAlp Dener . -tao_BNK_gamma2 - gamma2 interpolation update
839eb910715SAlp Dener . -tao_BNK_gamma3 - gamma3 interpolation update
840eb910715SAlp Dener . -tao_BNK_gamma4 - gamma4 interpolation update
841eb910715SAlp Dener . -tao_BNK_theta - theta interpolation update
842eb910715SAlp Dener . -tao_BNK_omega1 - omega1 step update
843eb910715SAlp Dener . -tao_BNK_omega2 - omega2 step update
844eb910715SAlp Dener . -tao_BNK_omega3 - omega3 step update
845eb910715SAlp Dener . -tao_BNK_omega4 - omega4 step update
846eb910715SAlp Dener . -tao_BNK_omega5 - omega5 step update
847eb910715SAlp Dener . -tao_BNK_mu1_i -  mu1 interpolation init factor
848eb910715SAlp Dener . -tao_BNK_mu2_i -  mu2 interpolation init factor
849eb910715SAlp Dener . -tao_BNK_gamma1_i -  gamma1 interpolation init factor
850eb910715SAlp Dener . -tao_BNK_gamma2_i -  gamma2 interpolation init factor
851eb910715SAlp Dener . -tao_BNK_gamma3_i -  gamma3 interpolation init factor
852eb910715SAlp Dener . -tao_BNK_gamma4_i -  gamma4 interpolation init factor
853eb910715SAlp Dener - -tao_BNK_theta_i -  theta interpolation init factor
854eb910715SAlp Dener 
855eb910715SAlp Dener   Level: beginner
856eb910715SAlp Dener M*/
857eb910715SAlp Dener 
858eb910715SAlp Dener PetscErrorCode TaoCreate_BNK(Tao tao)
859eb910715SAlp Dener {
860eb910715SAlp Dener   TAO_BNK        *bnk;
861eb910715SAlp Dener   const char     *morethuente_type = TAOLINESEARCHMT;
862eb910715SAlp Dener   PetscErrorCode ierr;
863eb910715SAlp Dener 
864eb910715SAlp Dener   PetscFunctionBegin;
865eb910715SAlp Dener   ierr = PetscNewLog(tao,&bnk);CHKERRQ(ierr);
866eb910715SAlp Dener 
867eb910715SAlp Dener   tao->ops->setup = TaoSetUp_BNK;
868eb910715SAlp Dener   tao->ops->view = TaoView_BNK;
869eb910715SAlp Dener   tao->ops->setfromoptions = TaoSetFromOptions_BNK;
870eb910715SAlp Dener   tao->ops->destroy = TaoDestroy_BNK;
871eb910715SAlp Dener 
872eb910715SAlp Dener   /*  Override default settings (unless already changed) */
873eb910715SAlp Dener   if (!tao->max_it_changed) tao->max_it = 50;
874eb910715SAlp Dener   if (!tao->trust0_changed) tao->trust0 = 100.0;
875eb910715SAlp Dener 
876eb910715SAlp Dener   tao->data = (void*)bnk;
877eb910715SAlp Dener 
878*66ed3702SAlp Dener   /*  Hessian shifting parameters */
879eb910715SAlp Dener   bnk->sval   = 0.0;
880eb910715SAlp Dener   bnk->imin   = 1.0e-4;
881eb910715SAlp Dener   bnk->imax   = 1.0e+2;
882eb910715SAlp Dener   bnk->imfac  = 1.0e-1;
883eb910715SAlp Dener 
884eb910715SAlp Dener   bnk->pmin   = 1.0e-12;
885eb910715SAlp Dener   bnk->pmax   = 1.0e+2;
886eb910715SAlp Dener   bnk->pgfac  = 1.0e+1;
887eb910715SAlp Dener   bnk->psfac  = 4.0e-1;
888eb910715SAlp Dener   bnk->pmgfac = 1.0e-1;
889eb910715SAlp Dener   bnk->pmsfac = 1.0e-1;
890eb910715SAlp Dener 
891eb910715SAlp Dener   /*  Default values for trust-region radius update based on steplength */
892eb910715SAlp Dener   bnk->nu1 = 0.25;
893eb910715SAlp Dener   bnk->nu2 = 0.50;
894eb910715SAlp Dener   bnk->nu3 = 1.00;
895eb910715SAlp Dener   bnk->nu4 = 1.25;
896eb910715SAlp Dener 
897eb910715SAlp Dener   bnk->omega1 = 0.25;
898eb910715SAlp Dener   bnk->omega2 = 0.50;
899eb910715SAlp Dener   bnk->omega3 = 1.00;
900eb910715SAlp Dener   bnk->omega4 = 2.00;
901eb910715SAlp Dener   bnk->omega5 = 4.00;
902eb910715SAlp Dener 
903eb910715SAlp Dener   /*  Default values for trust-region radius update based on reduction */
904eb910715SAlp Dener   bnk->eta1 = 1.0e-4;
905eb910715SAlp Dener   bnk->eta2 = 0.25;
906eb910715SAlp Dener   bnk->eta3 = 0.50;
907eb910715SAlp Dener   bnk->eta4 = 0.90;
908eb910715SAlp Dener 
909eb910715SAlp Dener   bnk->alpha1 = 0.25;
910eb910715SAlp Dener   bnk->alpha2 = 0.50;
911eb910715SAlp Dener   bnk->alpha3 = 1.00;
912eb910715SAlp Dener   bnk->alpha4 = 2.00;
913eb910715SAlp Dener   bnk->alpha5 = 4.00;
914eb910715SAlp Dener 
915eb910715SAlp Dener   /*  Default values for trust-region radius update based on interpolation */
916eb910715SAlp Dener   bnk->mu1 = 0.10;
917eb910715SAlp Dener   bnk->mu2 = 0.50;
918eb910715SAlp Dener 
919eb910715SAlp Dener   bnk->gamma1 = 0.25;
920eb910715SAlp Dener   bnk->gamma2 = 0.50;
921eb910715SAlp Dener   bnk->gamma3 = 2.00;
922eb910715SAlp Dener   bnk->gamma4 = 4.00;
923eb910715SAlp Dener 
924eb910715SAlp Dener   bnk->theta = 0.05;
925eb910715SAlp Dener 
926eb910715SAlp Dener   /*  Default values for trust region initialization based on interpolation */
927eb910715SAlp Dener   bnk->mu1_i = 0.35;
928eb910715SAlp Dener   bnk->mu2_i = 0.50;
929eb910715SAlp Dener 
930eb910715SAlp Dener   bnk->gamma1_i = 0.0625;
931eb910715SAlp Dener   bnk->gamma2_i = 0.5;
932eb910715SAlp Dener   bnk->gamma3_i = 2.0;
933eb910715SAlp Dener   bnk->gamma4_i = 5.0;
934eb910715SAlp Dener 
935eb910715SAlp Dener   bnk->theta_i = 0.25;
936eb910715SAlp Dener 
937eb910715SAlp Dener   /*  Remaining parameters */
938eb910715SAlp Dener   bnk->min_radius = 1.0e-10;
939eb910715SAlp Dener   bnk->max_radius = 1.0e10;
940eb910715SAlp Dener   bnk->epsilon = 1.0e-6;
941eb910715SAlp Dener 
942eb910715SAlp Dener   bnk->pc_type         = BNK_PC_BFGS;
943eb910715SAlp Dener   bnk->bfgs_scale_type = BFGS_SCALE_PHESS;
944eb910715SAlp Dener   bnk->init_type       = BNK_INIT_INTERPOLATION;
945fed79b8eSAlp Dener   bnk->update_type     = BNK_UPDATE_INTERPOLATION;
946eb910715SAlp Dener 
947eb910715SAlp Dener   ierr = TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch);CHKERRQ(ierr);
948eb910715SAlp Dener   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr);
949eb910715SAlp Dener   ierr = TaoLineSearchSetType(tao->linesearch,morethuente_type);CHKERRQ(ierr);
950eb910715SAlp Dener   ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr);
951eb910715SAlp Dener   ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr);
952eb910715SAlp Dener 
953eb910715SAlp Dener   /*  Set linear solver to default for symmetric matrices */
954eb910715SAlp Dener   ierr = KSPCreate(((PetscObject)tao)->comm,&tao->ksp);CHKERRQ(ierr);
955eb910715SAlp Dener   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1);CHKERRQ(ierr);
956eb910715SAlp Dener   ierr = KSPSetOptionsPrefix(tao->ksp,tao->hdr.prefix);CHKERRQ(ierr);
957eb910715SAlp Dener   ierr = KSPSetType(tao->ksp,KSPCGSTCG);CHKERRQ(ierr);
958eb910715SAlp Dener   PetscFunctionReturn(0);
959eb910715SAlp Dener }
960