xref: /petsc/src/tao/bound/impls/bncg/bncg.c (revision 8ebe3e4e9e00d86ece2e9fcd0cc84910b0ad437c)
1 #include <petsctaolinesearch.h>
2 #include <../src/tao/bound/impls/bncg/bncg.h> /*I "petsctao.h" I*/
3 #include <petscksp.h>
4 
5 #define CG_GradientDescent      0
6 #define CG_HestenesStiefel      1
7 #define CG_FletcherReeves       2
8 #define CG_PolakRibierePolyak   3
9 #define CG_PolakRibierePlus     4
10 #define CG_DaiYuan              5
11 #define CG_HagerZhang           6
12 #define CG_DaiKou               7
13 #define CG_KouDai               8
14 #define CG_SSML_BFGS            9
15 #define CG_SSML_DFP             10
16 #define CG_SSML_BROYDEN         11
17 #define CG_PCGradientDescent    12
18 #define CGTypes                 13
19 
20 static const char *CG_Table[64] = {"gd", "hs", "fr", "pr", "prp", "dy", "hz", "dk", "kd", "ssml_bfgs", "ssml_dfp", "ssml_brdn", "pcgd"};
21 
22 #define CG_AS_NONE       0
23 #define CG_AS_BERTSEKAS  1
24 #define CG_AS_SIZE       2
25 
26 static const char *CG_AS_TYPE[64] = {"none", "bertsekas"};
27 
28 PetscErrorCode TaoBNCGEstimateActiveSet(Tao tao, PetscInt asType)
29 {
30   PetscErrorCode               ierr;
31   TAO_BNCG                     *cg = (TAO_BNCG *)tao->data;
32 
33   PetscFunctionBegin;
34   ierr = ISDestroy(&cg->inactive_old);CHKERRQ(ierr);
35   if (cg->inactive_idx) {
36     ierr = ISDuplicate(cg->inactive_idx, &cg->inactive_old);CHKERRQ(ierr);
37     ierr = ISCopy(cg->inactive_idx, cg->inactive_old);CHKERRQ(ierr);
38   }
39   switch (asType) {
40   case CG_AS_NONE:
41     ierr = ISDestroy(&cg->inactive_idx);CHKERRQ(ierr);
42     ierr = VecWhichInactive(tao->XL, tao->solution, cg->unprojected_gradient, tao->XU, PETSC_TRUE, &cg->inactive_idx);CHKERRQ(ierr);
43     ierr = ISDestroy(&cg->active_idx);CHKERRQ(ierr);
44     ierr = ISComplementVec(cg->inactive_idx, tao->solution, &cg->active_idx);CHKERRQ(ierr);
45     break;
46 
47   case CG_AS_BERTSEKAS:
48     /* Use gradient descent to estimate the active set */
49     ierr = VecCopy(cg->unprojected_gradient, cg->W);CHKERRQ(ierr);
50     ierr = VecScale(cg->W, -1.0);CHKERRQ(ierr);
51     ierr = TaoEstimateActiveBounds(tao->solution, tao->XL, tao->XU, cg->unprojected_gradient, cg->W, cg->work, cg->as_step, &cg->as_tol,
52                                    &cg->active_lower, &cg->active_upper, &cg->active_fixed, &cg->active_idx, &cg->inactive_idx);CHKERRQ(ierr);
53     break;
54 
55   default:
56     break;
57   }
58   PetscFunctionReturn(0);
59 }
60 
61 PetscErrorCode TaoBNCGBoundStep(Tao tao, PetscInt asType, Vec step)
62 {
63   PetscErrorCode               ierr;
64   TAO_BNCG                     *cg = (TAO_BNCG *)tao->data;
65 
66   PetscFunctionBegin;
67   switch (asType) {
68   case CG_AS_NONE:
69     ierr = VecISSet(step, cg->active_idx, 0.0);CHKERRQ(ierr);
70     break;
71 
72   case CG_AS_BERTSEKAS:
73     ierr = TaoBoundStep(tao->solution, tao->XL, tao->XU, cg->active_lower, cg->active_upper, cg->active_fixed, 1.0, step);CHKERRQ(ierr);
74     break;
75 
76   default:
77     break;
78   }
79   PetscFunctionReturn(0);
80 }
81 
82 static PetscErrorCode TaoSolve_BNCG(Tao tao)
83 {
84   TAO_BNCG                     *cg = (TAO_BNCG*)tao->data;
85   PetscErrorCode               ierr;
86   PetscReal                    step=1.0,gnorm,gnorm2, resnorm;
87   PetscInt                     nDiff;
88 
89   PetscFunctionBegin;
90   /*   Project the current point onto the feasible set */
91   ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr);
92   ierr = TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU);CHKERRQ(ierr);
93 
94   /* Project the initial point onto the feasible region */
95   ierr = TaoBoundSolution(tao->solution, tao->XL,tao->XU, 0.0, &nDiff, tao->solution);CHKERRQ(ierr);
96 
97   if (nDiff > 0 || !tao->recycle) {
98     ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &cg->f, cg->unprojected_gradient);CHKERRQ(ierr);
99   }
100   ierr = VecNorm(cg->unprojected_gradient,NORM_2,&gnorm);CHKERRQ(ierr);
101   if (PetscIsInfOrNanReal(cg->f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_USER, "User provided compute function generated Inf or NaN");
102 
103   /* Estimate the active set and compute the projected gradient */
104   ierr = TaoBNCGEstimateActiveSet(tao, cg->as_type);CHKERRQ(ierr);
105 
106   /* Project the gradient and calculate the norm */
107   ierr = VecCopy(cg->unprojected_gradient, tao->gradient);CHKERRQ(ierr);
108   ierr = VecISSet(tao->gradient, cg->active_idx, 0.0);CHKERRQ(ierr);
109   ierr = VecNorm(tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr);
110   gnorm2 = gnorm*gnorm;
111 
112   /* Initialize counters */
113   tao->niter = 0;
114   cg->ls_fails = cg->descent_error = 0;
115   cg->resets = -1;
116   cg->skipped_updates = cg->pure_gd_steps = 0;
117   cg->iter_quad = 0;
118 
119   /* Convergence test at the starting point. */
120   tao->reason = TAO_CONTINUE_ITERATING;
121 
122   ierr = VecFischer(tao->solution, cg->unprojected_gradient, tao->XL, tao->XU, cg->W);CHKERRQ(ierr);
123   ierr = VecNorm(cg->W, NORM_2, &resnorm);CHKERRQ(ierr);
124   if (PetscIsInfOrNanReal(resnorm)) SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_USER, "User provided compute function generated Inf or NaN");
125   ierr = TaoLogConvergenceHistory(tao, cg->f, resnorm, 0.0, tao->ksp_its);CHKERRQ(ierr);
126   ierr = TaoMonitor(tao, tao->niter, cg->f, resnorm, 0.0, step);CHKERRQ(ierr);
127   ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
128   if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
129   /* Calculate initial direction. */
130   if (!tao->recycle) {
131     /* We are not recycling a solution/history from a past TaoSolve */
132     ierr = TaoBNCGResetUpdate(tao, gnorm2);CHKERRQ(ierr);
133   }
134   /* Initial gradient descent step. Scaling by 1.0 also does a decent job for some problems. */
135   while (1) {
136     /* Call general purpose update function */
137     if (tao->ops->update) {
138       ierr = (*tao->ops->update)(tao, tao->niter, tao->user_update);CHKERRQ(ierr);
139     }
140     ierr = TaoBNCGConductIteration(tao, gnorm);CHKERRQ(ierr);
141     if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
142   }
143 }
144 
145 static PetscErrorCode TaoSetUp_BNCG(Tao tao)
146 {
147   TAO_BNCG         *cg = (TAO_BNCG*)tao->data;
148   PetscErrorCode ierr;
149 
150   PetscFunctionBegin;
151   if (!tao->gradient) {
152     ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);
153   }
154   if (!tao->stepdirection) {
155     ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr);
156   }
157   if (!cg->W) {
158     ierr = VecDuplicate(tao->solution,&cg->W);CHKERRQ(ierr);
159   }
160   if (!cg->work) {
161     ierr = VecDuplicate(tao->solution,&cg->work);CHKERRQ(ierr);
162   }
163   if (!cg->sk) {
164     ierr = VecDuplicate(tao->solution,&cg->sk);CHKERRQ(ierr);
165   }
166   if (!cg->yk) {
167     ierr = VecDuplicate(tao->gradient,&cg->yk);CHKERRQ(ierr);
168   }
169   if (!cg->X_old) {
170     ierr = VecDuplicate(tao->solution,&cg->X_old);CHKERRQ(ierr);
171   }
172   if (!cg->G_old) {
173     ierr = VecDuplicate(tao->gradient,&cg->G_old);CHKERRQ(ierr);
174   }
175   if (cg->diag_scaling) {
176     ierr = VecDuplicate(tao->solution,&cg->d_work);CHKERRQ(ierr);
177     ierr = VecDuplicate(tao->solution,&cg->y_work);CHKERRQ(ierr);
178     ierr = VecDuplicate(tao->solution,&cg->g_work);CHKERRQ(ierr);
179   }
180   if (!cg->unprojected_gradient) {
181     ierr = VecDuplicate(tao->gradient,&cg->unprojected_gradient);CHKERRQ(ierr);
182   }
183   if (!cg->unprojected_gradient_old) {
184     ierr = VecDuplicate(tao->gradient,&cg->unprojected_gradient_old);CHKERRQ(ierr);
185   }
186   ierr = MatLMVMAllocate(cg->B, cg->sk, cg->yk);CHKERRQ(ierr);
187   if (cg->pc) {
188     ierr = MatLMVMSetJ0(cg->B, cg->pc);CHKERRQ(ierr);
189   }
190   PetscFunctionReturn(0);
191 }
192 
193 static PetscErrorCode TaoDestroy_BNCG(Tao tao)
194 {
195   TAO_BNCG       *cg = (TAO_BNCG*) tao->data;
196   PetscErrorCode ierr;
197 
198   PetscFunctionBegin;
199   if (tao->setupcalled) {
200     ierr = VecDestroy(&cg->W);CHKERRQ(ierr);
201     ierr = VecDestroy(&cg->work);CHKERRQ(ierr);
202     ierr = VecDestroy(&cg->X_old);CHKERRQ(ierr);
203     ierr = VecDestroy(&cg->G_old);CHKERRQ(ierr);
204     ierr = VecDestroy(&cg->unprojected_gradient);CHKERRQ(ierr);
205     ierr = VecDestroy(&cg->unprojected_gradient_old);CHKERRQ(ierr);
206     ierr = VecDestroy(&cg->g_work);CHKERRQ(ierr);
207     ierr = VecDestroy(&cg->d_work);CHKERRQ(ierr);
208     ierr = VecDestroy(&cg->y_work);CHKERRQ(ierr);
209     ierr = VecDestroy(&cg->sk);CHKERRQ(ierr);
210     ierr = VecDestroy(&cg->yk);CHKERRQ(ierr);
211   }
212   ierr = ISDestroy(&cg->active_lower);CHKERRQ(ierr);
213   ierr = ISDestroy(&cg->active_upper);CHKERRQ(ierr);
214   ierr = ISDestroy(&cg->active_fixed);CHKERRQ(ierr);
215   ierr = ISDestroy(&cg->active_idx);CHKERRQ(ierr);
216   ierr = ISDestroy(&cg->inactive_idx);CHKERRQ(ierr);
217   ierr = ISDestroy(&cg->inactive_old);CHKERRQ(ierr);
218   ierr = ISDestroy(&cg->new_inactives);CHKERRQ(ierr);
219   ierr = MatDestroy(&cg->B);CHKERRQ(ierr);
220   if (cg->pc) {
221     ierr = MatDestroy(&cg->pc);CHKERRQ(ierr);
222   }
223   ierr = PetscFree(tao->data);CHKERRQ(ierr);
224   PetscFunctionReturn(0);
225 }
226 
227 static PetscErrorCode TaoSetFromOptions_BNCG(PetscOptionItems *PetscOptionsObject,Tao tao)
228 {
229   TAO_BNCG       *cg = (TAO_BNCG*)tao->data;
230   PetscErrorCode ierr;
231 
232   PetscFunctionBegin;
233   ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr);
234   ierr = PetscOptionsHead(PetscOptionsObject,"Nonlinear Conjugate Gradient method for unconstrained optimization");CHKERRQ(ierr);
235   ierr = PetscOptionsEList("-tao_bncg_type","cg formula", "", CG_Table, CGTypes, CG_Table[cg->cg_type], &cg->cg_type,NULL);CHKERRQ(ierr);
236   if (cg->cg_type != CG_SSML_BFGS) cg->alpha = -1.0; /* Setting defaults for non-BFGS methods. User can change it below. */
237   if (CG_GradientDescent == cg->cg_type) {
238     cg->cg_type = CG_PCGradientDescent;
239     /* Set scaling equal to none or, at best, scalar scaling. */
240     cg->unscaled_restart = PETSC_TRUE;
241     cg->diag_scaling = PETSC_FALSE;
242   }
243   ierr = PetscOptionsEList("-tao_bncg_as_type","active set estimation method", "", CG_AS_TYPE, CG_AS_SIZE, CG_AS_TYPE[cg->cg_type], &cg->cg_type,NULL);CHKERRQ(ierr);
244   ierr = PetscOptionsReal("-tao_bncg_hz_eta","(developer) cutoff tolerance for HZ", "", cg->hz_eta,&cg->hz_eta,NULL);CHKERRQ(ierr);
245   ierr = PetscOptionsReal("-tao_bncg_eps","(developer) cutoff value for restarts", "", cg->epsilon,&cg->epsilon,NULL);CHKERRQ(ierr);
246   ierr = PetscOptionsReal("-tao_bncg_dk_eta","(developer) cutoff tolerance for DK", "", cg->dk_eta,&cg->dk_eta,NULL);CHKERRQ(ierr);
247   ierr = PetscOptionsReal("-tao_bncg_xi","(developer) Parameter in the KD method", "", cg->xi,&cg->xi,NULL);CHKERRQ(ierr);
248   ierr = PetscOptionsReal("-tao_bncg_theta", "(developer) update parameter for the Broyden method", "", cg->theta, &cg->theta, NULL);CHKERRQ(ierr);
249   ierr = PetscOptionsReal("-tao_bncg_hz_theta", "(developer) parameter for the HZ (2006) method", "", cg->hz_theta, &cg->hz_theta, NULL);CHKERRQ(ierr);
250   ierr = PetscOptionsReal("-tao_bncg_alpha","(developer) parameter for the scalar scaling","",cg->alpha,&cg->alpha,NULL);CHKERRQ(ierr);
251   ierr = PetscOptionsReal("-tao_bncg_bfgs_scale", "(developer) update parameter for bfgs/brdn CG methods", "", cg->bfgs_scale, &cg->bfgs_scale, NULL);CHKERRQ(ierr);
252   ierr = PetscOptionsReal("-tao_bncg_dfp_scale", "(developer) update parameter for bfgs/brdn CG methods", "", cg->dfp_scale, &cg->dfp_scale, NULL);CHKERRQ(ierr);
253   ierr = PetscOptionsBool("-tao_bncg_diag_scaling","Enable diagonal Broyden-like preconditioning","",cg->diag_scaling,&cg->diag_scaling,NULL);CHKERRQ(ierr);
254   ierr = PetscOptionsBool("-tao_bncg_dynamic_restart","(developer) use dynamic restarts as in HZ, DK, KD","",cg->use_dynamic_restart,&cg->use_dynamic_restart,NULL);CHKERRQ(ierr);
255   ierr = PetscOptionsBool("-tao_bncg_unscaled_restart","(developer) use unscaled gradient restarts","",cg->unscaled_restart,&cg->unscaled_restart,NULL);CHKERRQ(ierr);
256   ierr = PetscOptionsReal("-tao_bncg_zeta", "(developer) Free parameter for the Kou-Dai method", "", cg->zeta, &cg->zeta, NULL);CHKERRQ(ierr);
257   ierr = PetscOptionsInt("-tao_bncg_min_quad", "(developer) Number of iterations with approximate quadratic behavior needed for restart", "", cg->min_quad, &cg->min_quad, NULL);CHKERRQ(ierr);
258   ierr = PetscOptionsInt("-tao_bncg_min_restart_num", "(developer) Number of iterations between restarts (times dimension)", "", cg->min_restart_num, &cg->min_restart_num, NULL);CHKERRQ(ierr);
259   ierr = PetscOptionsBool("-tao_bncg_spaced_restart","(developer) Enable regular steepest descent restarting every fixed number of iterations","",cg->spaced_restart,&cg->spaced_restart,NULL);CHKERRQ(ierr);
260   ierr = PetscOptionsBool("-tao_bncg_no_scaling","Disable all scaling except in restarts","",cg->no_scaling,&cg->no_scaling,NULL);CHKERRQ(ierr);
261   if (cg->no_scaling) {
262     cg->diag_scaling = PETSC_FALSE;
263     cg->alpha = -1.0;
264   }
265   if (cg->alpha == -1.0 && cg->cg_type == CG_KouDai && !cg->diag_scaling) { /* Some more default options that appear to be good. */
266     cg->neg_xi = PETSC_TRUE;
267   }
268   ierr = PetscOptionsBool("-tao_bncg_neg_xi","(developer) Use negative xi when it might be a smaller descent direction than necessary","",cg->neg_xi,&cg->neg_xi,NULL);CHKERRQ(ierr);
269   ierr = PetscOptionsReal("-tao_bncg_as_tol", "(developer) initial tolerance used when estimating actively bounded variables","",cg->as_tol,&cg->as_tol,NULL);CHKERRQ(ierr);
270   ierr = PetscOptionsReal("-tao_bncg_as_step", "(developer) step length used when estimating actively bounded variables","",cg->as_step,&cg->as_step,NULL);CHKERRQ(ierr);
271   ierr = PetscOptionsReal("-tao_bncg_delta_min", "(developer) minimum scaling factor used for scaled gradient restarts","",cg->delta_min,&cg->delta_min,NULL);CHKERRQ(ierr);
272   ierr = PetscOptionsReal("-tao_bncg_delta_max", "(developer) maximum scaling factor used for scaled gradient restarts","",cg->delta_max,&cg->delta_max,NULL);CHKERRQ(ierr);
273 
274   ierr = PetscOptionsTail();CHKERRQ(ierr);
275   ierr = MatSetOptionsPrefix(cg->B, ((PetscObject)tao)->prefix);CHKERRQ(ierr);
276   ierr = MatAppendOptionsPrefix(cg->B, "tao_bncg_");CHKERRQ(ierr);
277   ierr = MatSetFromOptions(cg->B);CHKERRQ(ierr);
278   PetscFunctionReturn(0);
279 }
280 
281 static PetscErrorCode TaoView_BNCG(Tao tao, PetscViewer viewer)
282 {
283   PetscBool      isascii;
284   TAO_BNCG       *cg = (TAO_BNCG*)tao->data;
285   PetscErrorCode ierr;
286 
287   PetscFunctionBegin;
288   ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii);CHKERRQ(ierr);
289   if (isascii) {
290     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
291     ierr = PetscViewerASCIIPrintf(viewer, "CG Type: %s\n", CG_Table[cg->cg_type]);CHKERRQ(ierr);
292     ierr = PetscViewerASCIIPrintf(viewer, "Skipped Stepdirection Updates: %i\n", cg->skipped_updates);CHKERRQ(ierr);
293     ierr = PetscViewerASCIIPrintf(viewer, "Scaled gradient steps: %i\n", cg->resets);CHKERRQ(ierr);
294     ierr = PetscViewerASCIIPrintf(viewer, "Pure gradient steps: %i\n", cg->pure_gd_steps);CHKERRQ(ierr);
295     ierr = PetscViewerASCIIPrintf(viewer, "Not a descent direction: %i\n", cg->descent_error);CHKERRQ(ierr);
296     ierr = PetscViewerASCIIPrintf(viewer, "Line search fails: %i\n", cg->ls_fails);CHKERRQ(ierr);
297     if (cg->diag_scaling) {
298       ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
299       if (isascii) {
300         ierr = PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
301         ierr = MatView(cg->B, viewer);CHKERRQ(ierr);
302         ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
303       }
304     }
305     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
306   }
307   PetscFunctionReturn(0);
308 }
309 
310 PetscErrorCode TaoBNCGComputeScalarScaling(PetscReal yty, PetscReal yts, PetscReal sts, PetscReal *scale, PetscReal alpha)
311 {
312   PetscReal a, b, c, sig1, sig2;
313 
314   PetscFunctionBegin;
315   *scale = 0.0;
316   if (1.0 == alpha) *scale = yts/yty;
317   else if (0.0 == alpha) *scale = sts/yts;
318   else if (-1.0 == alpha) *scale = 1.0;
319   else {
320     a = yty;
321     b = yts;
322     c = sts;
323     a *= alpha;
324     b *= -(2.0*alpha - 1.0);
325     c *= alpha - 1.0;
326     sig1 = (-b + PetscSqrtReal(b*b - 4.0*a*c))/(2.0*a);
327     sig2 = (-b - PetscSqrtReal(b*b - 4.0*a*c))/(2.0*a);
328     /* accept the positive root as the scalar */
329     if (sig1 > 0.0) *scale = sig1;
330     else if (sig2 > 0.0) *scale = sig2;
331     else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_CONV_FAILED, "Cannot find positive scalar");
332   }
333   PetscFunctionReturn(0);
334 }
335 
336 /*MC
337   TAOBNCG - Bound-constrained Nonlinear Conjugate Gradient method.
338 
339   Options Database Keys:
340 + -tao_bncg_recycle - enable recycling the latest calculated gradient vector in subsequent TaoSolve() calls (currently disabled)
341 . -tao_bncg_eta <r> - restart tolerance
342 . -tao_bncg_type <taocg_type> - cg formula
343 . -tao_bncg_as_type <none,bertsekas> - active set estimation method
344 . -tao_bncg_as_tol <r> - tolerance used in Bertsekas active-set estimation
345 . -tao_bncg_as_step <r> - trial step length used in Bertsekas active-set estimation
346 . -tao_bncg_eps <r> - cutoff used for determining whether or not we restart based on steplength each iteration, as well as determining whether or not we continue using the last stepdirection. Defaults to machine precision.
347 . -tao_bncg_theta <r> - convex combination parameter for the Broyden method
348 . -tao_bncg_hz_eta <r> - cutoff tolerance for the beta term in the HZ, DK methods
349 . -tao_bncg_dk_eta <r> - cutoff tolerance for the beta term in the HZ, DK methods
350 . -tao_bncg_xi <r> - Multiplicative constant of the gamma term in the KD method
351 . -tao_bncg_hz_theta <r> - Multiplicative constant of the theta term for the HZ method
352 . -tao_bncg_bfgs_scale <r> - Scaling parameter of the bfgs contribution to the scalar Broyden method
353 . -tao_bncg_dfp_scale <r> - Scaling parameter of the dfp contribution to the scalar Broyden method
354 . -tao_bncg_diag_scaling <b> - Whether or not to use diagonal initialization/preconditioning for the CG methods. Default True.
355 . -tao_bncg_dynamic_restart <b> - use dynamic restart strategy in the HZ, DK, KD methods
356 . -tao_bncg_unscaled_restart <b> - whether or not to scale the gradient when doing gradient descent restarts
357 . -tao_bncg_zeta <r> - Scaling parameter in the KD method
358 . -tao_bncg_delta_min <r> - Minimum bound for rescaling during restarted gradient descent steps
359 . -tao_bncg_delta_max <r> - Maximum bound for rescaling during restarted gradient descent steps
360 . -tao_bncg_min_quad <i> - Number of quadratic-like steps in a row necessary to do a dynamic restart
361 . -tao_bncg_min_restart_num <i> - This number, x, makes sure there is a gradient descent step every x*n iterations, where n is the dimension of the problem
362 . -tao_bncg_spaced_restart <b> - whether or not to do gradient descent steps every x*n iterations
363 . -tao_bncg_no_scaling <b> - If true, eliminates all scaling, including defaults.
364 - -tao_bncg_neg_xi <b> - Whether or not to use negative xi in the KD method under certain conditions
365 
366   Notes:
367     CG formulas are:
368 + "gd" - Gradient Descent
369 . "fr" - Fletcher-Reeves
370 . "pr" - Polak-Ribiere-Polyak
371 . "prp" - Polak-Ribiere-Plus
372 . "hs" - Hestenes-Steifel
373 . "dy" - Dai-Yuan
374 . "ssml_bfgs" - Self-Scaling Memoryless BFGS
375 . "ssml_dfp"  - Self-Scaling Memoryless DFP
376 . "ssml_brdn" - Self-Scaling Memoryless Broyden
377 . "hz" - Hager-Zhang (CG_DESCENT 5.3)
378 . "dk" - Dai-Kou (2013)
379 - "kd" - Kou-Dai (2015)
380 
381   Level: beginner
382 M*/
383 
384 PETSC_EXTERN PetscErrorCode TaoCreate_BNCG(Tao tao)
385 {
386   TAO_BNCG       *cg;
387   const char     *morethuente_type = TAOLINESEARCHMT;
388   PetscErrorCode ierr;
389 
390   PetscFunctionBegin;
391   tao->ops->setup = TaoSetUp_BNCG;
392   tao->ops->solve = TaoSolve_BNCG;
393   tao->ops->view = TaoView_BNCG;
394   tao->ops->setfromoptions = TaoSetFromOptions_BNCG;
395   tao->ops->destroy = TaoDestroy_BNCG;
396 
397   /* Override default settings (unless already changed) */
398   if (!tao->max_it_changed) tao->max_it = 2000;
399   if (!tao->max_funcs_changed) tao->max_funcs = 4000;
400 
401   /*  Note: nondefault values should be used for nonlinear conjugate gradient  */
402   /*  method.  In particular, gtol should be less that 0.5; the value used in  */
403   /*  Nocedal and Wright is 0.10.  We use the default values for the  */
404   /*  linesearch because it seems to work better. */
405   ierr = TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch);CHKERRQ(ierr);
406   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr);
407   ierr = TaoLineSearchSetType(tao->linesearch, morethuente_type);CHKERRQ(ierr);
408   ierr = TaoLineSearchUseTaoRoutines(tao->linesearch, tao);CHKERRQ(ierr);
409 
410   ierr = PetscNewLog(tao,&cg);CHKERRQ(ierr);
411   tao->data = (void*)cg;
412   ierr = KSPInitializePackage();CHKERRQ(ierr);
413   ierr = MatCreate(PetscObjectComm((PetscObject)tao), &cg->B);CHKERRQ(ierr);
414   ierr = PetscObjectIncrementTabLevel((PetscObject)cg->B, (PetscObject)tao, 1);CHKERRQ(ierr);
415   ierr = MatSetType(cg->B, MATLMVMDIAGBROYDEN);CHKERRQ(ierr);
416 
417   cg->pc = NULL;
418 
419   cg->dk_eta = 0.5;
420   cg->hz_eta = 0.4;
421   cg->dynamic_restart = PETSC_FALSE;
422   cg->unscaled_restart = PETSC_FALSE;
423   cg->no_scaling = PETSC_FALSE;
424   cg->delta_min = 1e-7;
425   cg->delta_max = 100;
426   cg->theta = 1.0;
427   cg->hz_theta = 1.0;
428   cg->dfp_scale = 1.0;
429   cg->bfgs_scale = 1.0;
430   cg->zeta = 0.1;
431   cg->min_quad = 6;
432   cg->min_restart_num = 6; /* As in CG_DESCENT and KD2015*/
433   cg->xi = 1.0;
434   cg->neg_xi = PETSC_TRUE;
435   cg->spaced_restart = PETSC_FALSE;
436   cg->tol_quad = 1e-8;
437   cg->as_step = 0.001;
438   cg->as_tol = 0.001;
439   cg->eps_23 = PetscPowReal(PETSC_MACHINE_EPSILON, 2.0/3.0); /* Just a little tighter*/
440   cg->as_type = CG_AS_BERTSEKAS;
441   cg->cg_type = CG_SSML_BFGS;
442   cg->alpha = 1.0;
443   cg->diag_scaling = PETSC_TRUE;
444   PetscFunctionReturn(0);
445 }
446 
447 PetscErrorCode TaoBNCGResetUpdate(Tao tao, PetscReal gnormsq)
448 {
449    TAO_BNCG          *cg = (TAO_BNCG*)tao->data;
450    PetscErrorCode    ierr;
451    PetscReal         scaling;
452 
453    PetscFunctionBegin;
454    ++cg->resets;
455    scaling = 2.0 * PetscMax(1.0, PetscAbsScalar(cg->f)) / PetscMax(gnormsq, cg->eps_23);
456    scaling = PetscMin(cg->delta_max, PetscMax(cg->delta_min, scaling));
457    if (cg->unscaled_restart) {
458      scaling = 1.0;
459      ++cg->pure_gd_steps;
460    }
461    ierr = VecAXPBY(tao->stepdirection, -scaling, 0.0, tao->gradient);CHKERRQ(ierr);
462    /* Also want to reset our diagonal scaling with each restart */
463    if (cg->diag_scaling) {
464      ierr = MatLMVMReset(cg->B, PETSC_FALSE);CHKERRQ(ierr);
465    }
466    PetscFunctionReturn(0);
467  }
468 
469 PetscErrorCode TaoBNCGCheckDynamicRestart(Tao tao, PetscReal stepsize, PetscReal gd, PetscReal gd_old, PetscBool *dynrestart, PetscReal fold)
470 {
471    TAO_BNCG          *cg = (TAO_BNCG*)tao->data;
472    PetscReal         quadinterp;
473 
474    PetscFunctionBegin;
475    if (cg->f < cg->min_quad/10) {
476      *dynrestart = PETSC_FALSE;
477      PetscFunctionReturn(0);
478    } /* just skip this since this strategy doesn't work well for functions near zero */
479    quadinterp = 2.0*(cg->f - fold)/(stepsize*(gd + gd_old));
480    if (PetscAbs(quadinterp - 1.0) < cg->tol_quad) ++cg->iter_quad;
481    else {
482      cg->iter_quad = 0;
483      *dynrestart = PETSC_FALSE;
484    }
485    if (cg->iter_quad >= cg->min_quad) {
486      cg->iter_quad = 0;
487      *dynrestart = PETSC_TRUE;
488    }
489    PetscFunctionReturn(0);
490  }
491 
492 PETSC_INTERN PetscErrorCode TaoBNCGStepDirectionUpdate(Tao tao, PetscReal gnorm2, PetscReal step, PetscReal fold, PetscReal gnorm2_old, PetscReal dnorm, PetscBool pcgd_fallback)
493 {
494   TAO_BNCG          *cg = (TAO_BNCG*)tao->data;
495   PetscErrorCode    ierr;
496   PetscReal         gamma = 1.0, tau_k, beta;
497   PetscReal         tmp = 1.0, ynorm, ynorm2 = 1.0, snorm = 1.0, dk_yk=1.0, gd;
498   PetscReal         gkp1_yk, gd_old, tau_bfgs, tau_dfp, gkp1D_yk, gtDg;
499   PetscInt          dim;
500   PetscBool         cg_restart = PETSC_FALSE;
501   PetscFunctionBegin;
502 
503   /* Local curvature check to see if we need to restart */
504   if (tao->niter >= 1 || tao->recycle) {
505     ierr = VecWAXPY(cg->yk, -1.0, cg->G_old, tao->gradient);CHKERRQ(ierr);
506     ierr = VecNorm(cg->yk, NORM_2, &ynorm);CHKERRQ(ierr);
507     ynorm2 = ynorm*ynorm;
508     ierr = VecDot(cg->yk, tao->stepdirection, &dk_yk);CHKERRQ(ierr);
509     if (step*dnorm < PETSC_MACHINE_EPSILON || step*dk_yk < PETSC_MACHINE_EPSILON) {
510       cg_restart = PETSC_TRUE;
511       ++cg->skipped_updates;
512     }
513     if (cg->spaced_restart) {
514       ierr = VecGetSize(tao->gradient, &dim);CHKERRQ(ierr);
515       if (tao->niter % (dim*cg->min_restart_num)) cg_restart = PETSC_TRUE;
516     }
517   }
518   /* If the user wants regular restarts, do it every 6n iterations, where n=dimension */
519   if (cg->spaced_restart) {
520     ierr = VecGetSize(tao->gradient, &dim);CHKERRQ(ierr);
521     if (0 == tao->niter % (6*dim)) cg_restart = PETSC_TRUE;
522   }
523   /* Compute the diagonal scaling vector if applicable */
524   if (cg->diag_scaling) {
525     ierr = MatLMVMUpdate(cg->B, tao->solution, tao->gradient);CHKERRQ(ierr);
526   }
527 
528   /* A note on diagonal scaling (to be added to paper):
529    For the FR, PR, PRP, and DY methods, the diagonally scaled versions
530    must be derived as a preconditioned CG method rather than as
531    a Hessian initialization like in the Broyden methods. */
532 
533   /* In that case, one writes the objective function as
534    f(x) \equiv f(Ay). Gradient evaluations yield g(x_k) = A g(Ay_k) = A g(x_k).
535    Furthermore, the direction d_k \equiv (x_k - x_{k-1})/step according to
536    HZ (2006) becomes A^{-1} d_k, such that d_k^T g_k remains the
537    same under preconditioning. Note that A is diagonal, such that A^T = A. */
538 
539   /* This yields questions like what the dot product d_k^T y_k
540    should look like. HZ mistakenly treats that as the same under
541    preconditioning, but that is not necessarily true. */
542 
543   /* Observe y_k \equiv g_k - g_{k-1}, and under the P.C. transformation,
544    we get d_k^T y_k = (d_k^T A_k^{-T} A_k g_k - d_k^T A_k^{-T} A_{k-1} g_{k-1}),
545    yielding d_k^T y_k = d_k^T g_k - d_k^T (A_k^{-T} A_{k-1} g_{k-1}), which is
546    NOT the same if our preconditioning matrix is updated between iterations.
547    This same issue is found when considering dot products of the form g_{k+1}^T y_k. */
548 
549   /* Compute CG step direction */
550   if (cg_restart) {
551     ierr = TaoBNCGResetUpdate(tao, gnorm2);CHKERRQ(ierr);
552   } else if (pcgd_fallback) {
553     /* Just like preconditioned CG */
554     ierr = MatSolve(cg->B, tao->gradient, cg->g_work);CHKERRQ(ierr);
555     ierr = VecAXPBY(tao->stepdirection, -1.0, 0.0, cg->g_work);CHKERRQ(ierr);
556   } else if (ynorm2 > PETSC_MACHINE_EPSILON) {
557     switch (cg->cg_type) {
558     case CG_PCGradientDescent:
559       if (!cg->diag_scaling) {
560         if (!cg->no_scaling) {
561         cg->sts = step*step*dnorm*dnorm;
562         ierr = TaoBNCGComputeScalarScaling(ynorm2, step*dk_yk, cg->sts, &tau_k, cg->alpha);CHKERRQ(ierr);
563         } else {
564           tau_k = 1.0;
565           ++cg->pure_gd_steps;
566         }
567         ierr = VecAXPBY(tao->stepdirection, -tau_k, 0.0, tao->gradient);CHKERRQ(ierr);
568       } else {
569         ierr = MatSolve(cg->B, tao->gradient, cg->g_work);CHKERRQ(ierr);
570         ierr = VecAXPBY(tao->stepdirection, -1.0, 0.0, cg->g_work);CHKERRQ(ierr);
571       }
572       break;
573 
574     case CG_HestenesStiefel:
575       /* Classic Hestenes-Stiefel method, modified with scalar and diagonal preconditioning. */
576       if (!cg->diag_scaling) {
577         cg->sts = step*step*dnorm*dnorm;
578         ierr = VecDot(cg->yk, tao->gradient, &gkp1_yk);CHKERRQ(ierr);
579         ierr = TaoBNCGComputeScalarScaling(ynorm2, step*dk_yk, cg->sts, &tau_k, cg->alpha);CHKERRQ(ierr);
580         beta = tau_k*gkp1_yk/dk_yk;
581         ierr = VecAXPBY(tao->stepdirection, -tau_k, beta, tao->gradient);CHKERRQ(ierr);
582       } else {
583         ierr = MatSolve(cg->B, tao->gradient, cg->g_work);CHKERRQ(ierr);
584         ierr = VecDot(cg->yk, cg->g_work, &gkp1_yk);CHKERRQ(ierr);
585         beta = gkp1_yk/dk_yk;
586         ierr = VecAXPBY(tao->stepdirection, -1.0, beta, cg->g_work);CHKERRQ(ierr);
587       }
588       break;
589 
590     case CG_FletcherReeves:
591       ierr = VecDot(cg->G_old, cg->G_old, &gnorm2_old);CHKERRQ(ierr);
592       ierr = VecWAXPY(cg->yk, -1.0, cg->G_old, tao->gradient);CHKERRQ(ierr);
593       ierr = VecNorm(cg->yk, NORM_2, &ynorm);CHKERRQ(ierr);
594       ynorm2 = ynorm*ynorm;
595       ierr = VecDot(cg->yk, tao->stepdirection, &dk_yk);CHKERRQ(ierr);
596       if (!cg->diag_scaling) {
597         ierr = TaoBNCGComputeScalarScaling(ynorm2, step*dk_yk, step*step*dnorm*dnorm, &tau_k, cg->alpha);CHKERRQ(ierr);
598         beta = tau_k*gnorm2/gnorm2_old;
599         ierr = VecAXPBY(tao->stepdirection, -tau_k, beta, tao->gradient);CHKERRQ(ierr);
600       } else {
601         ierr = VecDot(cg->G_old, cg->g_work, &gnorm2_old);CHKERRQ(ierr); /* Before it's updated */
602         ierr = MatSolve(cg->B, tao->gradient, cg->g_work);CHKERRQ(ierr);
603         ierr = VecDot(tao->gradient, cg->g_work, &tmp);CHKERRQ(ierr);
604         beta = tmp/gnorm2_old;
605         ierr = VecAXPBY(tao->stepdirection, -1.0, beta, cg->g_work);CHKERRQ(ierr);
606       }
607       break;
608 
609     case CG_PolakRibierePolyak:
610       snorm = step*dnorm;
611       if (!cg->diag_scaling) {
612         ierr = VecDot(cg->G_old, cg->G_old, &gnorm2_old);CHKERRQ(ierr);
613         ierr = VecDot(cg->yk, tao->gradient, &gkp1_yk);CHKERRQ(ierr);
614         ierr = TaoBNCGComputeScalarScaling(ynorm2, step*dk_yk, snorm*snorm, &tau_k, cg->alpha);CHKERRQ(ierr);
615         beta = tau_k*gkp1_yk/gnorm2_old;
616         ierr = VecAXPBY(tao->stepdirection, -tau_k, beta, tao->gradient);CHKERRQ(ierr);
617       } else {
618         ierr = VecDot(cg->G_old, cg->g_work, &gnorm2_old);CHKERRQ(ierr);
619         ierr = MatSolve(cg->B, tao->gradient, cg->g_work);CHKERRQ(ierr);
620         ierr = VecDot(cg->g_work, cg->yk, &gkp1_yk);CHKERRQ(ierr);
621         beta = gkp1_yk/gnorm2_old;
622         ierr = VecAXPBY(tao->stepdirection, -1.0, beta, cg->g_work);CHKERRQ(ierr);
623       }
624       break;
625 
626     case CG_PolakRibierePlus:
627       ierr = VecWAXPY(cg->yk, -1.0, cg->G_old, tao->gradient);CHKERRQ(ierr);
628       ierr = VecNorm(cg->yk, NORM_2, &ynorm);CHKERRQ(ierr);
629       ynorm2 = ynorm*ynorm;
630       if (!cg->diag_scaling) {
631         ierr = VecDot(cg->G_old, cg->G_old, &gnorm2_old);CHKERRQ(ierr);
632         ierr = VecDot(cg->yk, tao->gradient, &gkp1_yk);CHKERRQ(ierr);
633         ierr = TaoBNCGComputeScalarScaling(ynorm2, step*dk_yk, snorm*snorm, &tau_k, cg->alpha);CHKERRQ(ierr);
634         beta = tau_k*gkp1_yk/gnorm2_old;
635         beta = PetscMax(beta, 0.0);
636         ierr = VecAXPBY(tao->stepdirection, -tau_k, beta, tao->gradient);CHKERRQ(ierr);
637       } else {
638         ierr = VecDot(cg->G_old, cg->g_work, &gnorm2_old);CHKERRQ(ierr); /* Old gtDg */
639         ierr = MatSolve(cg->B, tao->gradient, cg->g_work);CHKERRQ(ierr);
640         ierr = VecDot(cg->g_work, cg->yk, &gkp1_yk);CHKERRQ(ierr);
641         beta = gkp1_yk/gnorm2_old;
642         beta = PetscMax(beta, 0.0);
643         ierr = VecAXPBY(tao->stepdirection, -1.0, beta, cg->g_work);CHKERRQ(ierr);
644       }
645       break;
646 
647     case CG_DaiYuan:
648       /* Dai, Yu-Hong, and Yaxiang Yuan. "A nonlinear conjugate gradient method with a strong global convergence property."
649          SIAM Journal on optimization 10, no. 1 (1999): 177-182. */
650       if (!cg->diag_scaling) {
651         ierr = VecDot(tao->stepdirection, tao->gradient, &gd);CHKERRQ(ierr);
652         ierr = VecDot(cg->G_old, tao->stepdirection, &gd_old);CHKERRQ(ierr);
653         ierr = TaoBNCGComputeScalarScaling(ynorm2, step*dk_yk, cg->yts, &tau_k, cg->alpha);CHKERRQ(ierr);
654         beta = tau_k*gnorm2/(gd - gd_old);
655         ierr = VecAXPBY(tao->stepdirection, -tau_k, beta, tao->gradient);CHKERRQ(ierr);
656       } else {
657         ierr = MatMult(cg->B, tao->stepdirection, cg->d_work);CHKERRQ(ierr);
658         ierr = MatSolve(cg->B, tao->gradient, cg->g_work);CHKERRQ(ierr);
659         ierr = VecDot(cg->g_work, tao->gradient, &gtDg);CHKERRQ(ierr);
660         ierr = VecDot(tao->stepdirection, cg->G_old, &gd_old);CHKERRQ(ierr);
661         ierr = VecDot(cg->d_work, cg->g_work, &dk_yk);CHKERRQ(ierr);
662         dk_yk = dk_yk - gd_old;
663         beta = gtDg/dk_yk;
664         ierr = VecScale(cg->d_work, beta);CHKERRQ(ierr);
665         ierr = VecWAXPY(tao->stepdirection, -1.0, cg->g_work, cg->d_work);CHKERRQ(ierr);
666       }
667       break;
668 
669     case CG_HagerZhang:
670       /* Hager, William W., and Hongchao Zhang. "Algorithm 851: CG_DESCENT, a conjugate gradient method with guaranteed descent."
671          ACM Transactions on Mathematical Software (TOMS) 32, no. 1 (2006): 113-137. */
672       ierr = VecDot(tao->gradient, tao->stepdirection, &gd);CHKERRQ(ierr);
673       ierr = VecDot(cg->G_old, tao->stepdirection, &gd_old);CHKERRQ(ierr);
674       ierr = VecWAXPY(cg->sk, -1.0, cg->X_old, tao->solution);CHKERRQ(ierr);
675       snorm = dnorm*step;
676       cg->yts = step*dk_yk;
677       if (cg->use_dynamic_restart) {
678         ierr = TaoBNCGCheckDynamicRestart(tao, step, gd, gd_old, &cg->dynamic_restart, fold);CHKERRQ(ierr);
679       }
680       if (cg->dynamic_restart) {
681         ierr = TaoBNCGResetUpdate(tao, gnorm2);CHKERRQ(ierr);
682       } else {
683         if (!cg->diag_scaling) {
684           ierr = VecDot(cg->yk, tao->gradient, &gkp1_yk);CHKERRQ(ierr);
685           ierr = TaoBNCGComputeScalarScaling(ynorm2, cg->yts, snorm*snorm, &tau_k, cg->alpha);CHKERRQ(ierr);
686           /* Supplying cg->alpha = -1.0 will give the CG_DESCENT 5.3 special case of tau_k = 1.0 */
687           tmp = gd/dk_yk;
688           beta = tau_k*(gkp1_yk/dk_yk - ynorm2*gd/(dk_yk*dk_yk));
689           /* Bound beta as in CG_DESCENT 5.3, as implemented, with the third comparison from DK 2013 */
690           beta = PetscMax(PetscMax(beta, cg->hz_eta*tau_k*gd_old/(dnorm*dnorm)), cg->dk_eta*tau_k*gd/(dnorm*dnorm));
691           /* d <- -t*g + beta*t*d */
692           ierr = VecAXPBY(tao->stepdirection, -tau_k, beta, tao->gradient);CHKERRQ(ierr);
693         } else {
694           /* We have diagonal scaling enabled and are taking a diagonally-scaled memoryless BFGS step */
695           cg->yty = ynorm2;
696           cg->sts = snorm*snorm;
697           /* Apply the diagonal scaling to all my vectors */
698           ierr = MatSolve(cg->B, tao->gradient, cg->g_work);CHKERRQ(ierr);
699           ierr = MatSolve(cg->B, cg->yk, cg->y_work);CHKERRQ(ierr);
700           ierr = MatSolve(cg->B, tao->stepdirection, cg->d_work);CHKERRQ(ierr);
701           /* Construct the constant ytDgkp1 */
702           ierr = VecDot(cg->yk, cg->g_work, &gkp1_yk);CHKERRQ(ierr);
703           /* Construct the constant for scaling Dkyk in the update */
704           tmp = gd/dk_yk;
705           ierr = VecDot(cg->yk, cg->y_work, &tau_k);CHKERRQ(ierr);
706           tau_k = -tau_k*gd/(dk_yk*dk_yk);
707           /* beta is the constant which adds the dk contribution */
708           beta = gkp1_yk/dk_yk + cg->hz_theta*tau_k; /* HZ; (1.15) from DK 2013 */
709           /* From HZ2013, modified to account for diagonal scaling*/
710           ierr = VecDot(cg->G_old, cg->d_work, &gd_old);CHKERRQ(ierr);
711           ierr = VecDot(tao->stepdirection, cg->g_work, &gd);CHKERRQ(ierr);
712           beta = PetscMax(PetscMax(beta, cg->hz_eta*gd_old/(dnorm*dnorm)), cg->dk_eta*gd/(dnorm*dnorm));
713           /* Do the update */
714           ierr = VecAXPBY(tao->stepdirection, -1.0, beta, cg->g_work);CHKERRQ(ierr);
715         }
716       }
717       break;
718 
719     case CG_DaiKou:
720       /* Dai, Yu-Hong, and Cai-Xia Kou. "A nonlinear conjugate gradient algorithm with an optimal property and an improved Wolfe line search."
721          SIAM Journal on Optimization 23, no. 1 (2013): 296-320. */
722       ierr = VecDot(tao->gradient, tao->stepdirection, &gd);CHKERRQ(ierr);
723       ierr = VecDot(cg->G_old, tao->stepdirection, &gd_old);CHKERRQ(ierr);
724       ierr = VecWAXPY(cg->sk, -1.0, cg->X_old, tao->solution);CHKERRQ(ierr);
725       snorm = step*dnorm;
726       cg->yts = dk_yk*step;
727       if (!cg->diag_scaling) {
728         ierr = VecDot(cg->yk, tao->gradient, &gkp1_yk);CHKERRQ(ierr);
729         ierr = TaoBNCGComputeScalarScaling(ynorm2, cg->yts, snorm*snorm, &tau_k, cg->alpha);CHKERRQ(ierr);
730         /* Use cg->alpha = -1.0 to get tau_k = 1.0 as in CG_DESCENT 5.3 */
731         tmp = gd/dk_yk;
732         beta = tau_k*(gkp1_yk/dk_yk - ynorm2*gd/(dk_yk*dk_yk) + gd/(dnorm*dnorm)) - step*gd/dk_yk;
733         beta = PetscMax(PetscMax(beta, cg->hz_eta*tau_k*gd_old/(dnorm*dnorm)), cg->dk_eta*tau_k*gd/(dnorm*dnorm));
734         /* d <- -t*g + beta*t*d */
735         ierr = VecAXPBYPCZ(tao->stepdirection, -tau_k, 0.0, beta, tao->gradient, cg->yk);CHKERRQ(ierr);
736       } else {
737         /* We have diagonal scaling enabled and are taking a diagonally-scaled memoryless BFGS step */
738         cg->yty = ynorm2;
739         cg->sts = snorm*snorm;
740         ierr = MatSolve(cg->B, tao->gradient, cg->g_work);CHKERRQ(ierr);
741         ierr = MatSolve(cg->B, cg->yk, cg->y_work);CHKERRQ(ierr);
742         ierr = MatSolve(cg->B, tao->stepdirection, cg->d_work);CHKERRQ(ierr);
743         /* Construct the constant ytDgkp1 */
744         ierr = VecDot(cg->yk, cg->g_work, &gkp1_yk);CHKERRQ(ierr);
745         ierr = VecDot(cg->yk, cg->y_work, &tau_k);CHKERRQ(ierr);
746         tau_k = tau_k*gd/(dk_yk*dk_yk);
747         tmp = gd/dk_yk;
748         /* beta is the constant which adds the dk contribution */
749         beta = gkp1_yk/dk_yk - step*tmp - tau_k;
750         /* Update this for the last term in beta */
751         ierr = VecDot(cg->y_work, tao->stepdirection, &dk_yk);CHKERRQ(ierr);
752         beta += tmp*dk_yk/(dnorm*dnorm); /* projection of y_work onto dk */
753         ierr = VecDot(tao->stepdirection, cg->g_work, &gd);CHKERRQ(ierr);
754         ierr = VecDot(cg->G_old, cg->d_work, &gd_old);CHKERRQ(ierr);
755         beta = PetscMax(PetscMax(beta, cg->hz_eta*gd_old/(dnorm*dnorm)), cg->dk_eta*gd/(dnorm*dnorm));
756         /* Do the update */
757         ierr = VecAXPBY(tao->stepdirection, -1.0, beta, cg->g_work);CHKERRQ(ierr);
758       }
759       break;
760 
761     case CG_KouDai:
762       /* Kou, Cai-Xia, and Yu-Hong Dai. "A modified self-scaling memoryless Broyden-Fletcher-Goldfarb-Shanno method for unconstrained optimization."
763          Journal of Optimization Theory and Applications 165, no. 1 (2015): 209-224. */
764       ierr = VecDot(tao->gradient, tao->stepdirection, &gd);CHKERRQ(ierr);
765       ierr = VecDot(cg->G_old, tao->stepdirection, &gd_old);CHKERRQ(ierr);
766       ierr = VecWAXPY(cg->sk, -1.0, cg->X_old, tao->solution);CHKERRQ(ierr);
767       snorm = step*dnorm;
768       cg->yts = dk_yk*step;
769       if (cg->use_dynamic_restart) {
770         ierr = TaoBNCGCheckDynamicRestart(tao, step, gd, gd_old, &cg->dynamic_restart, fold);CHKERRQ(ierr);
771       }
772       if (cg->dynamic_restart) {
773         ierr = TaoBNCGResetUpdate(tao, gnorm2);CHKERRQ(ierr);
774       } else {
775         if (!cg->diag_scaling) {
776           ierr = VecDot(cg->yk, tao->gradient, &gkp1_yk);CHKERRQ(ierr);
777           ierr = TaoBNCGComputeScalarScaling(ynorm2, cg->yts, snorm*snorm, &tau_k, cg->alpha);CHKERRQ(ierr);
778           beta = tau_k*(gkp1_yk/dk_yk - ynorm2*gd/(dk_yk*dk_yk)) - step*gd/dk_yk;
779           if (beta < cg->zeta*tau_k*gd/(dnorm*dnorm)) /* 0.1 is KD's zeta parameter */
780           {
781             beta = cg->zeta*tau_k*gd/(dnorm*dnorm);
782             gamma = 0.0;
783           } else {
784             if (gkp1_yk < 0 && cg->neg_xi) gamma = -1.0*gd/dk_yk;
785             /* This seems to be very effective when there's no tau_k scaling.
786                This guarantees a large descent step every iteration, going through DK 2015 Lemma 3.1's proof but allowing for negative xi */
787             else {
788               gamma = cg->xi*gd/dk_yk;
789             }
790           }
791           /* d <- -t*g + beta*t*d + t*tmp*yk */
792           ierr = VecAXPBYPCZ(tao->stepdirection, -tau_k, gamma*tau_k, beta, tao->gradient, cg->yk);CHKERRQ(ierr);
793         } else {
794           /* We have diagonal scaling enabled and are taking a diagonally-scaled memoryless BFGS step */
795           cg->yty = ynorm2;
796           cg->sts = snorm*snorm;
797           ierr = MatSolve(cg->B, tao->gradient, cg->g_work);CHKERRQ(ierr);
798           ierr = MatSolve(cg->B, cg->yk, cg->y_work);CHKERRQ(ierr);
799           /* Construct the constant ytDgkp1 */
800           ierr = VecDot(cg->yk, cg->g_work, &gkp1D_yk);CHKERRQ(ierr);
801           /* Construct the constant for scaling Dkyk in the update */
802           gamma = gd/dk_yk;
803           /* tau_k = -ytDy/(ytd)^2 * gd */
804           ierr = VecDot(cg->yk, cg->y_work, &tau_k);CHKERRQ(ierr);
805           tau_k = tau_k*gd/(dk_yk*dk_yk);
806           /* beta is the constant which adds the d_k contribution */
807           beta = gkp1D_yk/dk_yk - step*gamma - tau_k;
808           /* Here is the requisite check */
809           ierr = VecDot(tao->stepdirection, cg->g_work, &tmp);CHKERRQ(ierr);
810           if (cg->neg_xi) {
811             /* modified KD implementation */
812             if (gkp1D_yk/dk_yk < 0) gamma = -1.0*gd/dk_yk;
813             else {
814               gamma = cg->xi*gd/dk_yk;
815             }
816             if (beta < cg->zeta*tmp/(dnorm*dnorm)) {
817               beta = cg->zeta*tmp/(dnorm*dnorm);
818               gamma = 0.0;
819             }
820           } else { /* original KD 2015 implementation */
821             if (beta < cg->zeta*tmp/(dnorm*dnorm)) {
822               beta = cg->zeta*tmp/(dnorm*dnorm);
823               gamma = 0.0;
824             } else {
825               gamma = cg->xi*gd/dk_yk;
826             }
827           }
828           /* Do the update in two steps */
829           ierr = VecAXPBY(tao->stepdirection, -1.0, beta, cg->g_work);CHKERRQ(ierr);
830           ierr = VecAXPY(tao->stepdirection, gamma, cg->y_work);CHKERRQ(ierr);
831         }
832       }
833       break;
834 
835     case CG_SSML_BFGS:
836       /* Perry, J. M. "A class of conjugate gradient algorithms with a two-step variable-metric memory."
837          Discussion Papers 269 (1977). */
838       ierr = VecDot(tao->gradient, tao->stepdirection, &gd);CHKERRQ(ierr);
839       ierr = VecWAXPY(cg->sk, -1.0, cg->X_old, tao->solution);CHKERRQ(ierr);
840       snorm = step*dnorm;
841       cg->yts = dk_yk*step;
842       cg->yty = ynorm2;
843       cg->sts = snorm*snorm;
844       if (!cg->diag_scaling) {
845         ierr = VecDot(cg->yk, tao->gradient, &gkp1_yk);CHKERRQ(ierr);
846         ierr = TaoBNCGComputeScalarScaling(cg->yty, cg->yts, cg->sts, &tau_k, cg->alpha);CHKERRQ(ierr);
847         tmp = gd/dk_yk;
848         beta = tau_k*(gkp1_yk/dk_yk - cg->yty*gd/(dk_yk*dk_yk)) - step*tmp;
849         /* d <- -t*g + beta*t*d + t*tmp*yk */
850         ierr = VecAXPBYPCZ(tao->stepdirection, -tau_k, tmp*tau_k, beta, tao->gradient, cg->yk);CHKERRQ(ierr);
851       } else {
852         /* We have diagonal scaling enabled and are taking a diagonally-scaled memoryless BFGS step */
853         ierr = MatSolve(cg->B, tao->gradient, cg->g_work);CHKERRQ(ierr);
854         ierr = MatSolve(cg->B, cg->yk, cg->y_work);CHKERRQ(ierr);
855         /* compute scalar gamma */
856         ierr = VecDot(cg->g_work, cg->yk, &gkp1_yk);CHKERRQ(ierr);
857         ierr = VecDot(cg->y_work, cg->yk, &tmp);CHKERRQ(ierr);
858         gamma = gd/dk_yk;
859         /* Compute scalar beta */
860         beta = (gkp1_yk/dk_yk - gd*tmp/(dk_yk*dk_yk)) - step*gd/dk_yk;
861         /* Compute stepdirection d_kp1 = gamma*Dkyk + beta*dk - Dkgkp1 */
862         ierr = VecAXPBYPCZ(tao->stepdirection, -1.0, gamma, beta, cg->g_work, cg->y_work);CHKERRQ(ierr);
863       }
864       break;
865 
866     case CG_SSML_DFP:
867       ierr = VecDot(tao->gradient, tao->stepdirection, &gd);CHKERRQ(ierr);
868       ierr = VecWAXPY(cg->sk, -1.0, cg->X_old, tao->solution);CHKERRQ(ierr);
869       snorm = step*dnorm;
870       cg->yts = dk_yk*step;
871       cg->yty = ynorm2;
872       cg->sts = snorm*snorm;
873       if (!cg->diag_scaling) {
874         /* Instead of a regular convex combination, we will solve a quadratic formula. */
875         ierr = TaoBNCGComputeScalarScaling(cg->yty, cg->yts, cg->sts, &tau_k, cg->alpha);CHKERRQ(ierr);
876         ierr = VecDot(cg->yk, tao->gradient, &gkp1_yk);CHKERRQ(ierr);
877         tau_k = cg->dfp_scale*tau_k;
878         tmp = tau_k*gkp1_yk/cg->yty;
879         beta = -step*gd/dk_yk;
880         /* d <- -t*g + beta*d + tmp*yk */
881         ierr = VecAXPBYPCZ(tao->stepdirection, -tau_k, tmp, beta, tao->gradient, cg->yk);CHKERRQ(ierr);
882       } else {
883         /* We have diagonal scaling enabled and are taking a diagonally-scaled memoryless DFP step */
884         ierr = MatSolve(cg->B, tao->gradient, cg->g_work);CHKERRQ(ierr);
885         ierr = MatSolve(cg->B, cg->yk, cg->y_work);CHKERRQ(ierr);
886         /* compute scalar gamma */
887         ierr = VecDot(cg->g_work, cg->yk, &gkp1_yk);CHKERRQ(ierr);
888         ierr = VecDot(cg->y_work, cg->yk, &tmp);CHKERRQ(ierr);
889         gamma = (gkp1_yk/tmp);
890         /* Compute scalar beta */
891         beta = -step*gd/dk_yk;
892         /* Compute stepdirection d_kp1 = gamma*Dkyk + beta*dk - Dkgkp1 */
893         ierr = VecAXPBYPCZ(tao->stepdirection, -1.0, gamma, beta, cg->g_work, cg->y_work);CHKERRQ(ierr);
894       }
895       break;
896 
897     case CG_SSML_BROYDEN:
898       ierr = VecDot(tao->gradient, tao->stepdirection, &gd);CHKERRQ(ierr);
899       ierr = VecWAXPY(cg->sk, -1.0, cg->X_old, tao->solution);CHKERRQ(ierr);
900       snorm = step*dnorm;
901       cg->yts = step*dk_yk;
902       cg->yty = ynorm2;
903       cg->sts = snorm*snorm;
904       if (!cg->diag_scaling) {
905         /* Instead of a regular convex combination, we will solve a quadratic formula. */
906         ierr = TaoBNCGComputeScalarScaling(cg->yty, step*dk_yk, snorm*snorm, &tau_bfgs, cg->bfgs_scale);CHKERRQ(ierr);
907         ierr = TaoBNCGComputeScalarScaling(cg->yty, step*dk_yk, snorm*snorm, &tau_dfp, cg->dfp_scale);CHKERRQ(ierr);
908         ierr = VecDot(cg->yk, tao->gradient, &gkp1_yk);CHKERRQ(ierr);
909         tau_k = cg->theta*tau_bfgs + (1.0-cg->theta)*tau_dfp;
910         /* If bfgs_scale = 1.0, it should reproduce the bfgs tau_bfgs. If bfgs_scale = 0.0,
911            it should reproduce the tau_dfp scaling. Same with dfp_scale.   */
912         tmp = cg->theta*tau_bfgs*gd/dk_yk + (1-cg->theta)*tau_dfp*gkp1_yk/cg->yty;
913         beta = cg->theta*tau_bfgs*(gkp1_yk/dk_yk - cg->yty*gd/(dk_yk*dk_yk)) - step*gd/dk_yk;
914         /* d <- -t*g + beta*d + tmp*yk */
915         ierr = VecAXPBYPCZ(tao->stepdirection, -tau_k, tmp, beta, tao->gradient, cg->yk);CHKERRQ(ierr);
916       } else {
917         /* We have diagonal scaling enabled */
918         ierr = MatSolve(cg->B, tao->gradient, cg->g_work);CHKERRQ(ierr);
919         ierr = MatSolve(cg->B, cg->yk, cg->y_work);CHKERRQ(ierr);
920         /* compute scalar gamma */
921         ierr = VecDot(cg->g_work, cg->yk, &gkp1_yk);CHKERRQ(ierr);
922         ierr = VecDot(cg->y_work, cg->yk, &tmp);CHKERRQ(ierr);
923         gamma = cg->theta*gd/dk_yk + (1-cg->theta)*(gkp1_yk/tmp);
924         /* Compute scalar beta */
925         beta = cg->theta*(gkp1_yk/dk_yk - gd*tmp/(dk_yk*dk_yk)) - step*gd/dk_yk;
926         /* Compute stepdirection dkp1 = gamma*Dkyk + beta*dk - Dkgkp1 */
927         ierr = VecAXPBYPCZ(tao->stepdirection, -1.0, gamma, beta, cg->g_work, cg->y_work);CHKERRQ(ierr);
928       }
929       break;
930 
931     default:
932       break;
933 
934     }
935   }
936   PetscFunctionReturn(0);
937 }
938 
939 PETSC_INTERN PetscErrorCode TaoBNCGConductIteration(Tao tao, PetscReal gnorm)
940 {
941   TAO_BNCG                     *cg = (TAO_BNCG*)tao->data;
942   PetscErrorCode               ierr;
943   TaoLineSearchConvergedReason ls_status = TAOLINESEARCH_CONTINUE_ITERATING;
944   PetscReal                    step=1.0,gnorm2,gd,dnorm=0.0;
945   PetscReal                    gnorm2_old,f_old,resnorm, gnorm_old;
946   PetscBool                    pcgd_fallback = PETSC_FALSE;
947 
948   PetscFunctionBegin;
949   /* We are now going to perform a line search along the direction. */
950   /* Store solution and gradient info before it changes */
951   ierr = VecCopy(tao->solution, cg->X_old);CHKERRQ(ierr);
952   ierr = VecCopy(tao->gradient, cg->G_old);CHKERRQ(ierr);
953   ierr = VecCopy(cg->unprojected_gradient, cg->unprojected_gradient_old);CHKERRQ(ierr);
954 
955   gnorm_old = gnorm;
956   gnorm2_old = gnorm_old*gnorm_old;
957   f_old = cg->f;
958   /* Perform bounded line search. If we are recycling a solution from a previous */
959   /* TaoSolve, then we want to immediately skip to calculating a new direction rather than performing a linesearch */
960   if (!(tao->recycle && 0 == tao->niter)) {
961     /* Above logic: the below code happens every iteration, except for the first iteration of a recycled TaoSolve */
962     ierr = TaoLineSearchSetInitialStepLength(tao->linesearch, 1.0);CHKERRQ(ierr);
963     ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &cg->f, cg->unprojected_gradient, tao->stepdirection, &step, &ls_status);CHKERRQ(ierr);
964     ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
965 
966     /*  Check linesearch failure */
967     if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) {
968       ++cg->ls_fails;
969       if (cg->cg_type == CG_GradientDescent) {
970         /* Nothing left to do but fail out of the optimization */
971         step = 0.0;
972         tao->reason = TAO_DIVERGED_LS_FAILURE;
973       } else {
974         /* Restore previous point, perform preconditioned GD and regular GD steps at the last good point */
975         ierr = VecCopy(cg->X_old, tao->solution);CHKERRQ(ierr);
976         ierr = VecCopy(cg->G_old, tao->gradient);CHKERRQ(ierr);
977         ierr = VecCopy(cg->unprojected_gradient_old, cg->unprojected_gradient);CHKERRQ(ierr);
978         gnorm = gnorm_old;
979         gnorm2 = gnorm2_old;
980         cg->f = f_old;
981 
982         /* Fall back on preconditioned CG (so long as you're not already using it) */
983         if (cg->cg_type != CG_PCGradientDescent && cg->diag_scaling) {
984           pcgd_fallback = PETSC_TRUE;
985           ierr = TaoBNCGStepDirectionUpdate(tao, gnorm2, step, f_old, gnorm2_old, dnorm, pcgd_fallback);CHKERRQ(ierr);
986 
987           ierr = TaoBNCGResetUpdate(tao, gnorm2);CHKERRQ(ierr);
988           ierr = TaoBNCGBoundStep(tao, cg->as_type, tao->stepdirection);CHKERRQ(ierr);
989 
990           ierr = TaoLineSearchSetInitialStepLength(tao->linesearch, 1.0);CHKERRQ(ierr);
991           ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &cg->f, cg->unprojected_gradient, tao->stepdirection, &step, &ls_status);CHKERRQ(ierr);
992           ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
993 
994           pcgd_fallback = PETSC_FALSE;
995           if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) {
996             /* Going to perform a regular gradient descent step. */
997             ++cg->ls_fails;
998             step = 0.0;
999           }
1000         }
1001         /* Fall back on the scaled gradient step */
1002         if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) {
1003           ++cg->ls_fails;
1004           ierr = TaoBNCGResetUpdate(tao, gnorm2);CHKERRQ(ierr);
1005           ierr = TaoBNCGBoundStep(tao, cg->as_type, tao->stepdirection);CHKERRQ(ierr);
1006           ierr = TaoLineSearchSetInitialStepLength(tao->linesearch, 1.0);CHKERRQ(ierr);
1007           ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &cg->f, cg->unprojected_gradient, tao->stepdirection, &step, &ls_status);CHKERRQ(ierr);
1008           ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
1009         }
1010 
1011         if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) {
1012           /* Nothing left to do but fail out of the optimization */
1013           ++cg->ls_fails;
1014           step = 0.0;
1015           tao->reason = TAO_DIVERGED_LS_FAILURE;
1016         } else {
1017           /* One of the fallbacks worked. Set them both back equal to false. */
1018           pcgd_fallback = PETSC_FALSE;
1019         }
1020       }
1021     }
1022     /* Convergence test for line search failure */
1023     if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
1024 
1025     /* Standard convergence test */
1026     ierr = VecFischer(tao->solution, cg->unprojected_gradient, tao->XL, tao->XU, cg->W);CHKERRQ(ierr);
1027     ierr = VecNorm(cg->W, NORM_2, &resnorm);CHKERRQ(ierr);
1028     if (PetscIsInfOrNanReal(resnorm)) SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_USER, "User provided compute function generated Inf or NaN");
1029     ierr = TaoLogConvergenceHistory(tao, cg->f, resnorm, 0.0, tao->ksp_its);CHKERRQ(ierr);
1030     ierr = TaoMonitor(tao, tao->niter, cg->f, resnorm, 0.0, step);CHKERRQ(ierr);
1031     ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
1032     if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
1033   }
1034   /* Assert we have an updated step and we need at least one more iteration. */
1035   /* Calculate the next direction */
1036   /* Estimate the active set at the new solution */
1037   ierr = TaoBNCGEstimateActiveSet(tao, cg->as_type);CHKERRQ(ierr);
1038   /* Compute the projected gradient and its norm */
1039   ierr = VecCopy(cg->unprojected_gradient, tao->gradient);CHKERRQ(ierr);
1040   ierr = VecISSet(tao->gradient, cg->active_idx, 0.0);CHKERRQ(ierr);
1041   ierr = VecNorm(tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr);
1042   gnorm2 = gnorm*gnorm;
1043 
1044   /* Calculate some quantities used in the StepDirectionUpdate. */
1045   ierr = VecNorm(tao->stepdirection, NORM_2, &dnorm);CHKERRQ(ierr);
1046   /* Update the step direction. */
1047   ierr = TaoBNCGStepDirectionUpdate(tao, gnorm2, step, f_old, gnorm2_old, dnorm, pcgd_fallback);CHKERRQ(ierr);
1048   ++tao->niter;
1049   ierr = TaoBNCGBoundStep(tao, cg->as_type, tao->stepdirection);CHKERRQ(ierr);
1050 
1051   if (cg->cg_type != CG_GradientDescent) {
1052     /* Figure out which previously active variables became inactive this iteration */
1053     ierr = ISDestroy(&cg->new_inactives);CHKERRQ(ierr);
1054     if (cg->inactive_idx && cg->inactive_old) {
1055       ierr = ISDifference(cg->inactive_idx, cg->inactive_old, &cg->new_inactives);CHKERRQ(ierr);
1056     }
1057     /* Selectively reset the CG step those freshly inactive variables */
1058     if (cg->new_inactives) {
1059       ierr = VecGetSubVector(tao->stepdirection, cg->new_inactives, &cg->inactive_step);CHKERRQ(ierr);
1060       ierr = VecGetSubVector(cg->unprojected_gradient, cg->new_inactives, &cg->inactive_grad);CHKERRQ(ierr);
1061       ierr = VecCopy(cg->inactive_grad, cg->inactive_step);CHKERRQ(ierr);
1062       ierr = VecScale(cg->inactive_step, -1.0);CHKERRQ(ierr);
1063       ierr = VecRestoreSubVector(tao->stepdirection, cg->new_inactives, &cg->inactive_step);CHKERRQ(ierr);
1064       ierr = VecRestoreSubVector(cg->unprojected_gradient, cg->new_inactives, &cg->inactive_grad);CHKERRQ(ierr);
1065     }
1066     /* Verify that this is a descent direction */
1067     ierr = VecDot(tao->gradient, tao->stepdirection, &gd);CHKERRQ(ierr);
1068     ierr = VecNorm(tao->stepdirection, NORM_2, &dnorm);CHKERRQ(ierr);
1069     if (PetscIsInfOrNanReal(gd) || (gd/(dnorm*dnorm) <= -1e10 || gd/(dnorm*dnorm) >= -1e-10)) {
1070       /* Not a descent direction, so we reset back to projected gradient descent */
1071       ierr = TaoBNCGResetUpdate(tao, gnorm2);CHKERRQ(ierr);
1072       ierr = TaoBNCGBoundStep(tao, cg->as_type, tao->stepdirection);CHKERRQ(ierr);
1073       ++cg->descent_error;
1074     } else {
1075     }
1076   }
1077   PetscFunctionReturn(0);
1078 }
1079 
1080 PetscErrorCode TaoBNCGSetH0(Tao tao, Mat H0)
1081 {
1082   TAO_BNCG                     *cg = (TAO_BNCG*)tao->data;
1083   PetscErrorCode               ierr;
1084 
1085   PetscFunctionBegin;
1086   ierr = PetscObjectReference((PetscObject)H0);CHKERRQ(ierr);
1087   cg->pc = H0;
1088   PetscFunctionReturn(0);
1089 }
1090