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