xref: /petsc/src/tao/bound/impls/blmvm/blmvm.c (revision 4f4fdda464867e01e2f7b514a4f52c421e99e383)
1 #include <petsctaolinesearch.h>
2 #include <../src/tao/unconstrained/impls/lmvm/lmvm.h>
3 #include <../src/tao/bound/impls/blmvm/blmvm.h>
4 
5 #define BLMVM_STEP_BFGS     0
6 #define BLMVM_STEP_GRAD     1
7 
8 #define BLMVM_AS_NONE       0
9 #define BLMVM_AS_BERTSEKAS  1
10 #define BLMVM_AS_SIZE       2
11 
12 static const char *BLMVM_AS_TYPE[64] = {"none", "bertsekas"};
13 
14 PETSC_INTERN PetscErrorCode TaoBLMVMEstimateActiveSet(Tao tao, PetscInt asType)
15 {
16   PetscErrorCode               ierr;
17   TAO_BLMVM                     *blmP = (TAO_BLMVM *)tao->data;
18 
19   PetscFunctionBegin;
20   if (!tao->bounded) PetscFunctionReturn(0);
21   switch (asType) {
22   case BLMVM_AS_NONE:
23     ierr = ISDestroy(&blmP->inactive_idx);CHKERRQ(ierr);
24     ierr = VecWhichInactive(tao->XL, tao->solution, blmP->unprojected_gradient, tao->XU, PETSC_TRUE, &blmP->inactive_idx);CHKERRQ(ierr);
25     ierr = ISDestroy(&blmP->active_idx);CHKERRQ(ierr);
26     ierr = ISComplementVec(blmP->inactive_idx, tao->solution, &blmP->active_idx);CHKERRQ(ierr);
27     break;
28 
29   case BLMVM_AS_BERTSEKAS:
30     /* Use gradient descent to estimate the active set */
31     ierr = VecCopy(blmP->unprojected_gradient, blmP->W);CHKERRQ(ierr);
32     ierr = VecScale(blmP->W, -1.0);CHKERRQ(ierr);
33     ierr = TaoEstimateActiveBounds(tao->solution, tao->XL, tao->XU, blmP->unprojected_gradient, blmP->W, blmP->work, blmP->as_step, &blmP->as_tol,
34                                    &blmP->active_lower, &blmP->active_upper, &blmP->active_fixed, &blmP->active_idx, &blmP->inactive_idx);CHKERRQ(ierr);
35     break;
36 
37   default:
38     break;
39   }
40   PetscFunctionReturn(0);
41 }
42 
43 PETSC_INTERN PetscErrorCode TaoBLMVMBoundStep(Tao tao, PetscInt asType, Vec step)
44 {
45   PetscErrorCode               ierr;
46   TAO_BLMVM                     *blmP = (TAO_BLMVM *)tao->data;
47 
48   PetscFunctionBegin;
49   if (!tao->bounded) PetscFunctionReturn(0);
50   switch (asType) {
51   case BLMVM_AS_NONE:
52     ierr = VecISSet(step, blmP->active_idx, 0.0);CHKERRQ(ierr);
53     break;
54 
55   case BLMVM_AS_BERTSEKAS:
56     ierr = TaoBoundStep(tao->solution, tao->XL, tao->XU, blmP->active_lower, blmP->active_upper, blmP->active_fixed, 1.0, step);CHKERRQ(ierr);
57     break;
58 
59   default:
60     break;
61   }
62   PetscFunctionReturn(0);
63 }
64 
65 /*------------------------------------------------------------*/
66 PETSC_INTERN PetscErrorCode TaoSolve_BLMVM(Tao tao)
67 {
68   PetscErrorCode               ierr;
69   TAO_BLMVM                    *blmP = (TAO_BLMVM *)tao->data;
70   TaoLineSearchConvergedReason ls_status = TAOLINESEARCH_CONTINUE_ITERATING;
71   PetscReal                    f, fold, gdx, gnorm, resnorm;
72   PetscReal                    stepsize = 1.0;
73   PetscInt                     nDiff, nupdates, stepType = BLMVM_STEP_GRAD;
74 
75   PetscFunctionBegin;
76   /*  Project initial point onto bounds */
77   ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr);
78   ierr = TaoBoundSolution(tao->solution, tao->XL,tao->XU, 0.0, &nDiff, tao->solution);CHKERRQ(ierr);
79   if (tao->bounded) {
80     ierr = TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU);CHKERRQ(ierr);
81   }
82 
83   /* Check convergence criteria */
84   ierr = TaoComputeObjectiveAndGradient(tao, tao->solution,&f,blmP->unprojected_gradient);CHKERRQ(ierr);
85   ierr = VecBoundGradientProjection(blmP->unprojected_gradient,tao->solution, tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr);
86 
87   ierr = TaoGradientNorm(tao, blmP->unprojected_gradient,NORM_2,&gnorm);CHKERRQ(ierr);
88   if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
89 
90   tao->reason = TAO_CONTINUE_ITERATING;
91   ierr = VecFischer(tao->solution, blmP->unprojected_gradient, tao->XL, tao->XU, blmP->W);CHKERRQ(ierr);
92   ierr = VecNorm(blmP->W, NORM_2, &resnorm);CHKERRQ(ierr);
93   ierr = TaoLogConvergenceHistory(tao,f,resnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
94   ierr = TaoMonitor(tao,tao->niter,f,resnorm,0.0,stepsize);CHKERRQ(ierr);
95   ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
96   if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
97 
98   /* Set counter for gradient/reset steps */
99   if (!blmP->recycle) {
100     blmP->bfgs = 0;
101     blmP->grad = 0;
102     ierr = MatLMVMReset(blmP->M, PETSC_FALSE);CHKERRQ(ierr);
103   }
104 
105   /* Have not converged; continue with Newton method */
106   while (tao->reason == TAO_CONTINUE_ITERATING) {
107     /* Estimate active set at the current iterate */
108     ierr = TaoBLMVMEstimateActiveSet(tao, blmP->as_type);CHKERRQ(ierr);
109 
110     /* Compute direction */
111     if (blmP->H0) {
112       ierr = MatLMVMSetJ0(blmP->M, blmP->H0);CHKERRQ(ierr);
113       stepType = BLMVM_STEP_BFGS;
114     } else if (!blmP->no_scale) {
115       ierr = MatLMVMSetJ0(blmP->M, blmP->Mscale);CHKERRQ(ierr);
116     }
117     ierr = MatLMVMUpdate(blmP->M, tao->solution, blmP->unprojected_gradient);CHKERRQ(ierr);
118     ierr = VecCopy(blmP->unprojected_gradient, blmP->red_grad);CHKERRQ(ierr);
119     ierr = VecISSet(blmP->red_grad, blmP->active_idx, 0.0);CHKERRQ(ierr);
120     ierr = MatSolve(blmP->M, blmP->red_grad, tao->stepdirection);CHKERRQ(ierr);
121     ierr = VecScale(tao->stepdirection,-1.0);CHKERRQ(ierr);
122     ierr = TaoBLMVMBoundStep(tao, blmP->as_type, tao->stepdirection);CHKERRQ(ierr);
123     ierr = MatLMVMGetUpdateCount(blmP->M, &nupdates);CHKERRQ(ierr);
124     if (nupdates > 0) stepType = BLMVM_STEP_BFGS;
125 
126     /* Check for success (descent direction) */
127     ierr = VecDot(tao->stepdirection, blmP->red_grad, &gdx);CHKERRQ(ierr);
128     if (gdx >= 0 || PetscIsInfOrNanReal(gdx)) {
129       /* Step is not descent or solve was not successful
130          Use steepest descent direction (scaled) */
131       ierr = MatLMVMReset(blmP->M, PETSC_FALSE);CHKERRQ(ierr);
132       ierr = MatLMVMClearJ0(blmP->M);CHKERRQ(ierr);
133       ierr = MatLMVMUpdate(blmP->M, tao->solution, blmP->unprojected_gradient);CHKERRQ(ierr);
134       ierr = MatSolve(blmP->M, blmP->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
135       ierr = VecScale(tao->stepdirection,-1.0);CHKERRQ(ierr);
136       ierr = TaoBLMVMBoundStep(tao, blmP->as_type, tao->stepdirection);CHKERRQ(ierr);
137       stepType = BLMVM_STEP_GRAD;
138     }
139 
140     /* Perform the linesearch */
141     fold = f;
142     ierr = VecCopy(tao->solution, blmP->Xold);CHKERRQ(ierr);
143     ierr = VecCopy(blmP->unprojected_gradient, blmP->Gold);CHKERRQ(ierr);
144     ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, blmP->unprojected_gradient, tao->stepdirection, &stepsize, &ls_status);CHKERRQ(ierr);
145     ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
146 
147     if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER && stepType != BLMVM_STEP_GRAD) {
148       /* Linesearch failed
149          Reset factors and use scaled (projected) gradient step */
150       f = fold;
151       ierr = VecCopy(blmP->Xold, tao->solution);CHKERRQ(ierr);
152       ierr = VecCopy(blmP->Gold, blmP->unprojected_gradient);CHKERRQ(ierr);
153 
154       ierr = MatLMVMReset(blmP->M, PETSC_FALSE);CHKERRQ(ierr);
155       ierr = MatLMVMClearJ0(blmP->M);CHKERRQ(ierr);
156       ierr = MatLMVMUpdate(blmP->M, tao->solution, blmP->unprojected_gradient);CHKERRQ(ierr);
157       ierr = MatSolve(blmP->M, blmP->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
158       ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
159       ierr = TaoBLMVMBoundStep(tao, blmP->as_type, tao->stepdirection);CHKERRQ(ierr);
160       stepType = BLMVM_STEP_GRAD;
161 
162       ierr = TaoLineSearchApply(tao->linesearch,tao->solution,&f, blmP->unprojected_gradient, tao->stepdirection,  &stepsize, &ls_status);CHKERRQ(ierr);
163       ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
164     }
165 
166     if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) {
167       /* Line search failed on a gradient step, so just mark reason for divergence */
168       f = fold;
169       ierr = VecCopy(blmP->Xold, tao->solution);CHKERRQ(ierr);
170       ierr = VecCopy(blmP->Gold, blmP->unprojected_gradient);CHKERRQ(ierr);
171       tao->reason = TAO_DIVERGED_LS_FAILURE;
172     } else {
173       /* LS found valid step, so tally the step type and compute projected gradient */
174       switch (stepType) {
175       case BLMVM_STEP_BFGS:
176         ++blmP->bfgs;
177         break;
178       case BLMVM_STEP_GRAD:
179         ++blmP->grad;
180         break;
181       default:
182         break;
183       }
184       ierr = VecBoundGradientProjection(blmP->unprojected_gradient, tao->solution, tao->XL, tao->XU, tao->gradient);CHKERRQ(ierr);
185       ierr = TaoGradientNorm(tao, blmP->unprojected_gradient,NORM_2,&gnorm);CHKERRQ(ierr);
186       if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Not-a-Number");
187     }
188 
189     /* Check for converged */
190     tao->niter++;
191     ierr = VecFischer(tao->solution, blmP->unprojected_gradient, tao->XL, tao->XU, blmP->W);CHKERRQ(ierr);
192     ierr = VecNorm(blmP->W, NORM_2, &resnorm);CHKERRQ(ierr);
193     ierr = TaoLogConvergenceHistory(tao,f,resnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
194     ierr = TaoMonitor(tao,tao->niter,f,resnorm,0.0,stepsize);CHKERRQ(ierr);
195     ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
196   }
197   PetscFunctionReturn(0);
198 }
199 
200 PETSC_INTERN PetscErrorCode TaoSetup_BLMVM(Tao tao)
201 {
202   TAO_BLMVM      *blmP = (TAO_BLMVM *)tao->data;
203   PetscInt       n,N;
204   PetscErrorCode ierr;
205   PetscBool      is_spd;
206 
207   PetscFunctionBegin;
208   /* Existence of tao->solution checked in TaoSetup() */
209   ierr = VecDuplicate(tao->solution,&blmP->Xold);CHKERRQ(ierr);
210   ierr = VecDuplicate(tao->solution,&blmP->Gold);CHKERRQ(ierr);
211   ierr = VecDuplicate(tao->solution, &blmP->unprojected_gradient);CHKERRQ(ierr);
212   ierr = VecDuplicate(tao->solution, &blmP->W);CHKERRQ(ierr);
213   ierr = VecDuplicate(tao->solution, &blmP->work);CHKERRQ(ierr);
214   ierr = VecDuplicate(tao->solution, &blmP->red_grad);CHKERRQ(ierr);
215 
216   if (!tao->stepdirection) {
217     ierr = VecDuplicate(tao->solution, &tao->stepdirection);CHKERRQ(ierr);
218   }
219   if (!tao->gradient) {
220     ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);
221   }
222   /* Create matrix for the limited memory approximation */
223   ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr);
224   ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr);
225   ierr = MatSetSizes(blmP->M, n, n, N, N);CHKERRQ(ierr);
226   ierr = MatLMVMAllocate(blmP->M,tao->solution,blmP->unprojected_gradient);CHKERRQ(ierr);
227   ierr = MatGetOption(blmP->M, MAT_SPD, &is_spd);CHKERRQ(ierr);
228   if (!is_spd) SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM matrix is not symmetric positive-definite.");
229 
230   /* If the user has set a matrix to solve as the initial H0, set the options prefix here, and set up the KSP */
231   if (blmP->H0) {
232     ierr = MatLMVMSetJ0(blmP->M, blmP->H0);CHKERRQ(ierr);
233   } else if (!blmP->no_scale) {
234     if (!blmP->Mscale) {
235       ierr = MatCreateLMVMDiagBrdn(PetscObjectComm((PetscObject)blmP->M), n, N, &blmP->Mscale);CHKERRQ(ierr);
236       ierr = MatSetOptionsPrefix(blmP->Mscale, "tao_blmvm_scale_");CHKERRQ(ierr);
237       ierr = MatLMVMAllocate(blmP->Mscale, tao->solution, blmP->unprojected_gradient);CHKERRQ(ierr);
238     }
239     ierr = MatLMVMSetJ0(blmP->M, blmP->Mscale);CHKERRQ(ierr);
240   }
241   PetscFunctionReturn(0);
242 }
243 
244 /* ---------------------------------------------------------- */
245 PETSC_INTERN PetscErrorCode TaoDestroy_BLMVM(Tao tao)
246 {
247   TAO_BLMVM      *blmP = (TAO_BLMVM *)tao->data;
248   PetscErrorCode ierr;
249 
250   PetscFunctionBegin;
251   if (tao->setupcalled) {
252     ierr = VecDestroy(&blmP->unprojected_gradient);CHKERRQ(ierr);
253     ierr = VecDestroy(&blmP->Xold);CHKERRQ(ierr);
254     ierr = VecDestroy(&blmP->Gold);CHKERRQ(ierr);
255     ierr = VecDestroy(&blmP->W);CHKERRQ(ierr);
256     ierr = VecDestroy(&blmP->work);CHKERRQ(ierr);
257     ierr = VecDestroy(&blmP->red_grad);CHKERRQ(ierr);
258   }
259   ierr = ISDestroy(&blmP->active_lower);CHKERRQ(ierr);
260   ierr = ISDestroy(&blmP->active_upper);CHKERRQ(ierr);
261   ierr = ISDestroy(&blmP->active_fixed);CHKERRQ(ierr);
262   ierr = ISDestroy(&blmP->active_idx);CHKERRQ(ierr);
263   ierr = ISDestroy(&blmP->inactive_idx);CHKERRQ(ierr);
264   ierr = MatDestroy(&blmP->M);CHKERRQ(ierr);
265   if (blmP->H0) {
266     PetscObjectDereference((PetscObject)blmP->H0);
267   }
268   if (blmP->Mscale) {
269     ierr = MatDestroy(&blmP->Mscale);CHKERRQ(ierr);
270   }
271   ierr = PetscFree(tao->data);CHKERRQ(ierr);
272   PetscFunctionReturn(0);
273 }
274 
275 /*------------------------------------------------------------*/
276 PETSC_INTERN PetscErrorCode TaoSetFromOptions_BLMVM(PetscOptionItems* PetscOptionsObject,Tao tao)
277 {
278   TAO_BLMVM      *blmP = (TAO_BLMVM *)tao->data;
279   PetscErrorCode ierr;
280   PetscBool      is_lmvm, is_spd;
281 
282   PetscFunctionBegin;
283   ierr = PetscOptionsHead(PetscOptionsObject,"Limited-memory variable-metric method for bound constrained optimization");CHKERRQ(ierr);
284   ierr = PetscOptionsBool("-tao_blmvm_recycle","enable recycling of the BFGS matrix between subsequent TaoSolve() calls","",blmP->recycle,&blmP->recycle,NULL);CHKERRQ(ierr);
285   ierr = PetscOptionsBool("-tao_blmvm_no_scale","(developer) disable the diagonal Broyden scaling of the BFGS approximation","",blmP->no_scale,&blmP->no_scale,NULL);CHKERRQ(ierr);
286   ierr = PetscOptionsEList("-tao_blmvm_as_type","active set estimation method", "", BLMVM_AS_TYPE, BLMVM_AS_SIZE, BLMVM_AS_TYPE[blmP->as_type], &blmP->as_type,NULL);CHKERRQ(ierr);
287   ierr = PetscOptionsReal("-tao_blmvm_as_tol", "initial tolerance used when estimating actively bounded variables","",blmP->as_tol,&blmP->as_tol,NULL);CHKERRQ(ierr);
288   ierr = PetscOptionsReal("-tao_blmvm_as_step", "step length used when estimating actively bounded variables","",blmP->as_step,&blmP->as_step,NULL);CHKERRQ(ierr);
289   ierr = PetscOptionsTail();CHKERRQ(ierr);
290   ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr);
291   ierr = MatSetFromOptions(blmP->M);CHKERRQ(ierr);
292   ierr = PetscObjectBaseTypeCompare((PetscObject)blmP->M, MATLMVM, &is_lmvm);CHKERRQ(ierr);
293   if (!is_lmvm) SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "Matrix must be an LMVM-type");
294   ierr = MatGetOption(blmP->M, MAT_SPD, &is_spd);CHKERRQ(ierr);
295   if (!is_spd) SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM matrix must be symmetric positive-definite");
296   PetscFunctionReturn(0);
297 }
298 
299 /*------------------------------------------------------------*/
300 PETSC_INTERN PetscErrorCode TaoView_BLMVM(Tao tao, PetscViewer viewer)
301 {
302   TAO_BLMVM      *lmP = (TAO_BLMVM *)tao->data;
303   PetscBool      isascii;
304   PetscErrorCode ierr;
305   PetscInt       recycled_its;
306 
307   PetscFunctionBegin;
308   ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii);CHKERRQ(ierr);
309   if (isascii) {
310     ierr = PetscViewerASCIIPrintf(viewer, "  Gradient steps: %D\n", lmP->grad);CHKERRQ(ierr);
311     if (lmP->recycle) {
312       ierr = PetscViewerASCIIPrintf(viewer, "  Recycle: on\n");CHKERRQ(ierr);
313       recycled_its = lmP->bfgs + lmP->grad;
314       ierr = PetscViewerASCIIPrintf(viewer, "  Total recycled iterations: %D\n", recycled_its);CHKERRQ(ierr);
315     }
316     ierr = PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
317     ierr = MatView(lmP->M, viewer);CHKERRQ(ierr);
318     ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
319   }
320   PetscFunctionReturn(0);
321 }
322 
323 PETSC_INTERN PetscErrorCode TaoComputeDual_BLMVM(Tao tao, Vec DXL, Vec DXU)
324 {
325   TAO_BLMVM      *blm = (TAO_BLMVM *) tao->data;
326   PetscErrorCode ierr;
327 
328   PetscFunctionBegin;
329   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
330   PetscValidHeaderSpecific(DXL,VEC_CLASSID,2);
331   PetscValidHeaderSpecific(DXU,VEC_CLASSID,3);
332   if (!tao->gradient || !blm->unprojected_gradient) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Dual variables don't exist yet or no longer exist.\n");
333 
334   ierr = VecCopy(tao->gradient,DXL);CHKERRQ(ierr);
335   ierr = VecAXPY(DXL,-1.0,blm->unprojected_gradient);CHKERRQ(ierr);
336   ierr = VecSet(DXU,0.0);CHKERRQ(ierr);
337   ierr = VecPointwiseMax(DXL,DXL,DXU);CHKERRQ(ierr);
338 
339   ierr = VecCopy(blm->unprojected_gradient,DXU);CHKERRQ(ierr);
340   ierr = VecAXPY(DXU,-1.0,tao->gradient);CHKERRQ(ierr);
341   ierr = VecAXPY(DXU,1.0,DXL);CHKERRQ(ierr);
342   PetscFunctionReturn(0);
343 }
344 
345 /* ---------------------------------------------------------- */
346 /*MC
347   TAOBLMVM - Bounded limited memory variable metric is a quasi-Newton method
348          for nonlinear minimization with bound constraints. It is an extension
349          of TAOLMVM
350 
351   Options Database Keys:
352 .   -tao_blmvm_recycle - enable recycling LMVM updates between TaoSolve() calls
353 .   -tao_blmvm_no_scale - (developer) disables diagonal Broyden scaling on the LMVM approximation
354 .   -tao_blmvm_as_type - (developer) active set estimation method <none,bertsekas>
355 .   -tao_blmvm_as_tol - (developer) initial distance tolerance used when estimating the active set
356 .   -tao_blmvm_as_step - (developer) trial step length used when estimating the active set
357 
358   Level: beginner
359 M*/
360 PETSC_EXTERN PetscErrorCode TaoCreate_BLMVM(Tao tao)
361 {
362   TAO_BLMVM      *blmP;
363   const char     *morethuente_type = TAOLINESEARCHMT;
364   PetscErrorCode ierr;
365 
366   PetscFunctionBegin;
367   tao->ops->setup = TaoSetup_BLMVM;
368   tao->ops->solve = TaoSolve_BLMVM;
369   tao->ops->view = TaoView_BLMVM;
370   tao->ops->setfromoptions = TaoSetFromOptions_BLMVM;
371   tao->ops->destroy = TaoDestroy_BLMVM;
372   tao->ops->computedual = TaoComputeDual_BLMVM;
373 
374   ierr = PetscNewLog(tao,&blmP);CHKERRQ(ierr);
375   blmP->H0 = NULL;
376   blmP->as_step = 0.001;
377   blmP->as_tol = 0.001;
378   blmP->as_type = BLMVM_AS_BERTSEKAS;
379   blmP->recycle = PETSC_FALSE;
380   blmP->no_scale = PETSC_FALSE;
381 
382   tao->data = (void*)blmP;
383 
384   /* Override default settings (unless already changed) */
385   if (!tao->max_it_changed) tao->max_it = 2000;
386   if (!tao->max_funcs_changed) tao->max_funcs = 4000;
387 
388   ierr = TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch);CHKERRQ(ierr);
389   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr);
390   ierr = TaoLineSearchSetType(tao->linesearch, morethuente_type);CHKERRQ(ierr);
391   ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr);
392   ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr);
393 
394   ierr = KSPInitializePackage();CHKERRQ(ierr);
395   ierr = MatCreate(((PetscObject)tao)->comm, &blmP->M);CHKERRQ(ierr);
396   ierr = PetscObjectIncrementTabLevel((PetscObject)blmP->M, (PetscObject)tao, 1);CHKERRQ(ierr);
397   ierr = MatSetType(blmP->M, MATLMVMBFGS);CHKERRQ(ierr);
398   ierr = MatSetOptionsPrefix(blmP->M, "tao_blmvm_");CHKERRQ(ierr);
399   PetscFunctionReturn(0);
400 }
401 
402 PETSC_EXTERN PetscErrorCode TaoBLMVMSetH0(Tao tao, Mat H0)
403 {
404   TAO_LMVM       *lmP;
405   TAO_BLMVM      *blmP;
406   TaoType        type;
407   PetscBool      is_lmvm, is_blmvm;
408   PetscErrorCode ierr;
409 
410   ierr = TaoGetType(tao, &type);CHKERRQ(ierr);
411   ierr = PetscStrcmp(type, TAOLMVM,  &is_lmvm);CHKERRQ(ierr);
412   ierr = PetscStrcmp(type, TAOBLMVM, &is_blmvm);CHKERRQ(ierr);
413 
414   if (is_lmvm) {
415     lmP = (TAO_LMVM *)tao->data;
416     ierr = PetscObjectReference((PetscObject)H0);CHKERRQ(ierr);
417     lmP->H0 = H0;
418   } else if (is_blmvm) {
419     blmP = (TAO_BLMVM *)tao->data;
420     ierr = PetscObjectReference((PetscObject)H0);CHKERRQ(ierr);
421     blmP->H0 = H0;
422   } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "This routine applies to TAO_LMVM and TAO_BLMVM.");
423   PetscFunctionReturn(0);
424 }
425 
426 PETSC_EXTERN PetscErrorCode TaoBLMVMGetH0(Tao tao, Mat *H0)
427 {
428   TAO_LMVM       *lmP;
429   TAO_BLMVM      *blmP;
430   TaoType        type;
431   PetscBool      is_lmvm, is_blmvm;
432   Mat            M;
433 
434   PetscErrorCode ierr;
435 
436   ierr = TaoGetType(tao, &type);CHKERRQ(ierr);
437   ierr = PetscStrcmp(type, TAOLMVM,  &is_lmvm);CHKERRQ(ierr);
438   ierr = PetscStrcmp(type, TAOBLMVM, &is_blmvm);CHKERRQ(ierr);
439 
440   if (is_lmvm) {
441     lmP = (TAO_LMVM *)tao->data;
442     M = lmP->M;
443   } else if (is_blmvm) {
444     blmP = (TAO_BLMVM *)tao->data;
445     M = blmP->M;
446   } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "This routine applies to TAO_LMVM and TAO_BLMVM.");
447   ierr = MatLMVMGetJ0(M, H0);CHKERRQ(ierr);
448   PetscFunctionReturn(0);
449 }
450 
451 PETSC_EXTERN PetscErrorCode TaoBLMVMGetH0KSP(Tao tao, KSP *ksp)
452 {
453   TAO_LMVM       *lmP;
454   TAO_BLMVM      *blmP;
455   TaoType        type;
456   PetscBool      is_lmvm, is_blmvm;
457   Mat            M;
458   PetscErrorCode ierr;
459 
460   ierr = TaoGetType(tao, &type);CHKERRQ(ierr);
461   ierr = PetscStrcmp(type, TAOLMVM,  &is_lmvm);CHKERRQ(ierr);
462   ierr = PetscStrcmp(type, TAOBLMVM, &is_blmvm);CHKERRQ(ierr);
463 
464   if (is_lmvm) {
465     lmP = (TAO_LMVM *)tao->data;
466     M = lmP->M;
467   } else if (is_blmvm) {
468     blmP = (TAO_BLMVM *)tao->data;
469     M = blmP->M;
470   } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "This routine applies to TAO_LMVM and TAO_BLMVM.");
471   ierr = MatLMVMGetJ0KSP(M, ksp);CHKERRQ(ierr);
472   PetscFunctionReturn(0);
473 }
474