xref: /petsc/src/tao/bound/impls/blmvm/blmvm.c (revision cd929ea3f739fd9f7b6394f772cb40b9d4e6d97c)
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,delta;
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 initial scaling for the function */
99   if (!blmP->H0) {
100     delta = 2.0 * PetscMax(1.0, PetscAbsScalar(f)) / (gnorm*gnorm);
101     ierr = MatLMVMSetJ0Scale(blmP->M, delta);CHKERRQ(ierr);
102   }
103 
104 
105   /* Set counter for gradient/reset steps */
106   if (!blmP->recycle) {
107     blmP->bfgs = 0;
108     blmP->grad = 0;
109     ierr = MatLMVMReset(blmP->M, PETSC_FALSE);CHKERRQ(ierr);
110   }
111 
112   /* Have not converged; continue with Newton method */
113   while (tao->reason == TAO_CONTINUE_ITERATING) {
114     /* Estimate active set at the current iterate */
115     ierr = TaoBLMVMEstimateActiveSet(tao, blmP->as_type);CHKERRQ(ierr);
116 
117     /* Compute direction */
118     if (blmP->H0) {
119       ierr = MatLMVMSetJ0(blmP->M, blmP->H0);CHKERRQ(ierr);
120       stepType = BLMVM_STEP_BFGS;
121     }
122     ierr = MatLMVMUpdate(blmP->M, tao->solution, tao->gradient);CHKERRQ(ierr);
123     ierr = MatSolve(blmP->M, blmP->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
124     ierr = MatLMVMGetUpdateCount(blmP->M, &nupdates);CHKERRQ(ierr);
125     if (nupdates > 0) stepType = BLMVM_STEP_BFGS;
126 
127     /* Check for success (descent direction) */
128     ierr = VecDot(tao->stepdirection, tao->gradient, &gdx);CHKERRQ(ierr);
129     if (gdx <= 0 || PetscIsInfOrNanReal(gdx)) {
130       /* Step is not descent or solve was not successful
131          Use steepest descent direction (scaled) */
132 
133       delta = 2.0 * PetscMax(1.0, PetscAbsScalar(f)) / (gnorm*gnorm);
134       ierr = MatLMVMSetJ0Scale(blmP->M, delta);CHKERRQ(ierr);
135       ierr = MatLMVMReset(blmP->M, PETSC_FALSE);CHKERRQ(ierr);
136       ierr = MatLMVMUpdate(blmP->M, tao->solution, blmP->unprojected_gradient);CHKERRQ(ierr);
137       ierr = MatSolve(blmP->M, blmP->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
138       stepType = BLMVM_STEP_GRAD;
139     }
140     ierr = VecScale(tao->stepdirection,-1.0);CHKERRQ(ierr);
141     ierr = TaoBLMVMBoundStep(tao, blmP->as_type, tao->stepdirection);CHKERRQ(ierr);
142 
143     /* Perform the linesearch */
144     fold = f;
145     ierr = VecCopy(tao->solution, blmP->Xold);CHKERRQ(ierr);
146     ierr = VecCopy(blmP->unprojected_gradient, blmP->Gold);CHKERRQ(ierr);
147     ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, blmP->unprojected_gradient, tao->stepdirection, &stepsize, &ls_status);CHKERRQ(ierr);
148     ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
149 
150     if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER && stepType != BLMVM_STEP_GRAD) {
151       /* Linesearch failed
152          Reset factors and use scaled (projected) gradient step */
153       f = fold;
154       ierr = VecCopy(blmP->Xold, tao->solution);CHKERRQ(ierr);
155       ierr = VecCopy(blmP->Gold, blmP->unprojected_gradient);CHKERRQ(ierr);
156 
157       delta = 2.0 * PetscMax(1.0, PetscAbsScalar(f)) / (gnorm*gnorm);
158       ierr = MatLMVMSetJ0Scale(blmP->M, delta);CHKERRQ(ierr);
159       ierr = MatLMVMReset(blmP->M, PETSC_FALSE);CHKERRQ(ierr);
160       ierr = MatLMVMUpdate(blmP->M, tao->solution, blmP->unprojected_gradient);CHKERRQ(ierr);
161       ierr = MatSolve(blmP->M, blmP->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
162       ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
163       ierr = TaoBLMVMBoundStep(tao, blmP->as_type, tao->stepdirection);CHKERRQ(ierr);
164       stepType = BLMVM_STEP_GRAD;
165 
166       ierr = TaoLineSearchApply(tao->linesearch,tao->solution,&f, blmP->unprojected_gradient, tao->stepdirection,  &stepsize, &ls_status);CHKERRQ(ierr);
167       ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
168     }
169 
170     if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) {
171       /* Line search failed on a gradient step, so just mark reason for divergence */
172       f = fold;
173       ierr = VecCopy(blmP->Xold, tao->solution);CHKERRQ(ierr);
174       ierr = VecCopy(blmP->Gold, blmP->unprojected_gradient);CHKERRQ(ierr);
175       tao->reason = TAO_DIVERGED_LS_FAILURE;
176     } else {
177       /* LS found valid step, so tally the step type and compute projected gradient */
178       switch (stepType) {
179       case BLMVM_STEP_BFGS:
180         ++blmP->bfgs;
181         break;
182       case BLMVM_STEP_GRAD:
183         ++blmP->grad;
184         break;
185       default:
186         break;
187       }
188       ierr = VecBoundGradientProjection(blmP->unprojected_gradient, tao->solution, tao->XL, tao->XU, tao->gradient);CHKERRQ(ierr);
189       ierr = TaoGradientNorm(tao, blmP->unprojected_gradient,NORM_2,&gnorm);CHKERRQ(ierr);
190       if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Not-a-Number");
191     }
192 
193     /* Check for converged */
194     tao->niter++;
195     ierr = VecFischer(tao->solution, blmP->unprojected_gradient, tao->XL, tao->XU, blmP->W);CHKERRQ(ierr);
196     ierr = VecNorm(blmP->W, NORM_2, &resnorm);CHKERRQ(ierr);
197     ierr = TaoLogConvergenceHistory(tao,f,resnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
198     ierr = TaoMonitor(tao,tao->niter,f,resnorm,0.0,stepsize);CHKERRQ(ierr);
199     ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
200   }
201   PetscFunctionReturn(0);
202 }
203 
204 PETSC_INTERN PetscErrorCode TaoSetup_BLMVM(Tao tao)
205 {
206   TAO_BLMVM      *blmP = (TAO_BLMVM *)tao->data;
207   PetscInt       n,N;
208   PetscErrorCode ierr;
209 
210   PetscFunctionBegin;
211   /* Existence of tao->solution checked in TaoSetup() */
212   ierr = VecDuplicate(tao->solution,&blmP->Xold);CHKERRQ(ierr);
213   ierr = VecDuplicate(tao->solution,&blmP->Gold);CHKERRQ(ierr);
214   ierr = VecDuplicate(tao->solution, &blmP->unprojected_gradient);CHKERRQ(ierr);
215   ierr = VecDuplicate(tao->solution, &blmP->W);CHKERRQ(ierr);
216   ierr = VecDuplicate(tao->solution, &blmP->work);CHKERRQ(ierr);
217 
218   if (!tao->stepdirection) {
219     ierr = VecDuplicate(tao->solution, &tao->stepdirection);CHKERRQ(ierr);
220   }
221   if (!tao->gradient) {
222     ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);
223   }
224   /* Create matrix for the limited memory approximation */
225   ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr);
226   ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr);
227   ierr = MatSetSizes(blmP->M, n, n, N, N);CHKERRQ(ierr);
228   ierr = MatLMVMAllocate(blmP->M,tao->solution,blmP->unprojected_gradient);CHKERRQ(ierr);
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   }
234   PetscFunctionReturn(0);
235 }
236 
237 /* ---------------------------------------------------------- */
238 PETSC_INTERN PetscErrorCode TaoDestroy_BLMVM(Tao tao)
239 {
240   TAO_BLMVM      *blmP = (TAO_BLMVM *)tao->data;
241   PetscErrorCode ierr;
242 
243   PetscFunctionBegin;
244   if (tao->setupcalled) {
245     ierr = MatDestroy(&blmP->M);CHKERRQ(ierr);
246     ierr = VecDestroy(&blmP->unprojected_gradient);CHKERRQ(ierr);
247     ierr = VecDestroy(&blmP->Xold);CHKERRQ(ierr);
248     ierr = VecDestroy(&blmP->Gold);CHKERRQ(ierr);
249     ierr = VecDestroy(&blmP->W);CHKERRQ(ierr);
250     ierr = VecDestroy(&blmP->work);CHKERRQ(ierr);
251   }
252   ierr = ISDestroy(&blmP->active_lower);CHKERRQ(ierr);
253   ierr = ISDestroy(&blmP->active_upper);CHKERRQ(ierr);
254   ierr = ISDestroy(&blmP->active_fixed);CHKERRQ(ierr);
255   ierr = ISDestroy(&blmP->active_idx);CHKERRQ(ierr);
256   ierr = ISDestroy(&blmP->inactive_idx);CHKERRQ(ierr);
257   if (blmP->H0) {
258     PetscObjectDereference((PetscObject)blmP->H0);
259   }
260 
261   ierr = PetscFree(tao->data);CHKERRQ(ierr);
262   PetscFunctionReturn(0);
263 }
264 
265 /*------------------------------------------------------------*/
266 PETSC_INTERN PetscErrorCode TaoSetFromOptions_BLMVM(PetscOptionItems* PetscOptionsObject,Tao tao)
267 {
268   TAO_BLMVM      *blmP = (TAO_BLMVM *)tao->data;
269   PetscErrorCode ierr;
270 
271   PetscFunctionBegin;
272   ierr = PetscOptionsHead(PetscOptionsObject,"Limited-memory variable-metric method for bound constrained optimization");CHKERRQ(ierr);
273   ierr = PetscOptionsBool("-tao_blmvm_recycle","enable recycling of the BFGS matrix between subsequent TaoSolve() calls","",blmP->recycle,&blmP->recycle,NULL);CHKERRQ(ierr);
274   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);
275   ierr = PetscOptionsReal("-tao_blmvm_as_tol", "initial tolerance used when estimating actively bounded variables","",blmP->as_tol,&blmP->as_tol,NULL);CHKERRQ(ierr);
276   ierr = PetscOptionsReal("-tao_blmvm_as_step", "step length used when estimating actively bounded variables","",blmP->as_step,&blmP->as_step,NULL);CHKERRQ(ierr);
277   ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr);
278   ierr = PetscOptionsTail();CHKERRQ(ierr);
279   PetscFunctionReturn(0);
280 }
281 
282 
283 /*------------------------------------------------------------*/
284 PETSC_INTERN PetscErrorCode TaoView_BLMVM(Tao tao, PetscViewer viewer)
285 {
286   TAO_BLMVM      *lmP = (TAO_BLMVM *)tao->data;
287   PetscBool      isascii;
288   PetscErrorCode ierr;
289   PetscInt       recycled_its;
290 
291   PetscFunctionBegin;
292   ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii);CHKERRQ(ierr);
293   if (isascii) {
294     ierr = PetscViewerASCIIPrintf(viewer, "  Gradient steps: %D\n", lmP->grad);CHKERRQ(ierr);
295     if (lmP->recycle) {
296       ierr = PetscViewerASCIIPrintf(viewer, "  Recycle: on\n");CHKERRQ(ierr);
297       recycled_its = lmP->bfgs + lmP->grad;
298       ierr = PetscViewerASCIIPrintf(viewer, "  Total recycled iterations: %D\n", recycled_its);CHKERRQ(ierr);
299     }
300   }
301   PetscFunctionReturn(0);
302 }
303 
304 PETSC_INTERN PetscErrorCode TaoComputeDual_BLMVM(Tao tao, Vec DXL, Vec DXU)
305 {
306   TAO_BLMVM      *blm = (TAO_BLMVM *) tao->data;
307   PetscErrorCode ierr;
308 
309   PetscFunctionBegin;
310   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
311   PetscValidHeaderSpecific(DXL,VEC_CLASSID,2);
312   PetscValidHeaderSpecific(DXU,VEC_CLASSID,3);
313   if (!tao->gradient || !blm->unprojected_gradient) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Dual variables don't exist yet or no longer exist.\n");
314 
315   ierr = VecCopy(tao->gradient,DXL);CHKERRQ(ierr);
316   ierr = VecAXPY(DXL,-1.0,blm->unprojected_gradient);CHKERRQ(ierr);
317   ierr = VecSet(DXU,0.0);CHKERRQ(ierr);
318   ierr = VecPointwiseMax(DXL,DXL,DXU);CHKERRQ(ierr);
319 
320   ierr = VecCopy(blm->unprojected_gradient,DXU);CHKERRQ(ierr);
321   ierr = VecAXPY(DXU,-1.0,tao->gradient);CHKERRQ(ierr);
322   ierr = VecAXPY(DXU,1.0,DXL);CHKERRQ(ierr);
323   PetscFunctionReturn(0);
324 }
325 
326 /* ---------------------------------------------------------- */
327 /*MC
328   TAOBLMVM - Bounded limited memory variable metric is a quasi-Newton method
329          for nonlinear minimization with bound constraints. It is an extension
330          of TAOLMVM
331 
332   Level: beginner
333 M*/
334 PETSC_EXTERN PetscErrorCode TaoCreate_BLMVM(Tao tao)
335 {
336   TAO_BLMVM      *blmP;
337   const char     *morethuente_type = TAOLINESEARCHMT;
338   PetscErrorCode ierr;
339 
340   PetscFunctionBegin;
341   tao->ops->setup = TaoSetup_BLMVM;
342   tao->ops->solve = TaoSolve_BLMVM;
343   tao->ops->view = TaoView_BLMVM;
344   tao->ops->setfromoptions = TaoSetFromOptions_BLMVM;
345   tao->ops->destroy = TaoDestroy_BLMVM;
346   tao->ops->computedual = TaoComputeDual_BLMVM;
347 
348   ierr = PetscNewLog(tao,&blmP);CHKERRQ(ierr);
349   blmP->H0 = NULL;
350   blmP->as_step = 0.001;
351   blmP->as_tol = 0.001;
352   blmP->as_type = BLMVM_AS_BERTSEKAS;
353   blmP->recycle = PETSC_FALSE;
354 
355   tao->data = (void*)blmP;
356 
357   /* Override default settings (unless already changed) */
358   if (!tao->max_it_changed) tao->max_it = 2000;
359   if (!tao->max_funcs_changed) tao->max_funcs = 4000;
360 
361   ierr = TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch);CHKERRQ(ierr);
362   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr);
363   ierr = TaoLineSearchSetType(tao->linesearch, morethuente_type);CHKERRQ(ierr);
364   ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr);
365   ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr);
366 
367   ierr = KSPInitializePackage();CHKERRQ(ierr);
368   ierr = MatCreate(((PetscObject)tao)->comm, &blmP->M);CHKERRQ(ierr);
369   ierr = PetscObjectIncrementTabLevel((PetscObject)blmP->M, (PetscObject)tao, 1);CHKERRQ(ierr);
370   ierr = MatSetType(blmP->M, MATLBFGS);CHKERRQ(ierr);
371   ierr = MatSetOptionsPrefix(blmP->M, "tao_blmvm_");CHKERRQ(ierr);
372   PetscFunctionReturn(0);
373 }
374 
375 PETSC_EXTERN PetscErrorCode TaoBLMVMSetH0(Tao tao, Mat H0)
376 {
377   TAO_LMVM       *lmP;
378   TAO_BLMVM      *blmP;
379   TaoType        type;
380   PetscBool      is_lmvm, is_blmvm;
381   PetscErrorCode ierr;
382 
383   ierr = TaoGetType(tao, &type);CHKERRQ(ierr);
384   ierr = PetscStrcmp(type, TAOLMVM,  &is_lmvm);CHKERRQ(ierr);
385   ierr = PetscStrcmp(type, TAOBLMVM, &is_blmvm);CHKERRQ(ierr);
386 
387   if (is_lmvm) {
388     lmP = (TAO_LMVM *)tao->data;
389     ierr = PetscObjectReference((PetscObject)H0);CHKERRQ(ierr);
390     lmP->H0 = H0;
391   } else if (is_blmvm) {
392     blmP = (TAO_BLMVM *)tao->data;
393     ierr = PetscObjectReference((PetscObject)H0);CHKERRQ(ierr);
394     blmP->H0 = H0;
395   } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "This routine applies to TAO_LMVM and TAO_BLMVM.");
396   PetscFunctionReturn(0);
397 }
398 
399 PETSC_EXTERN PetscErrorCode TaoBLMVMGetH0(Tao tao, Mat *H0)
400 {
401   TAO_LMVM       *lmP;
402   TAO_BLMVM      *blmP;
403   TaoType        type;
404   PetscBool      is_lmvm, is_blmvm;
405   Mat            M;
406 
407   PetscErrorCode ierr;
408 
409   ierr = TaoGetType(tao, &type);CHKERRQ(ierr);
410   ierr = PetscStrcmp(type, TAOLMVM,  &is_lmvm);CHKERRQ(ierr);
411   ierr = PetscStrcmp(type, TAOBLMVM, &is_blmvm);CHKERRQ(ierr);
412 
413   if (is_lmvm) {
414     lmP = (TAO_LMVM *)tao->data;
415     M = lmP->M;
416   } else if (is_blmvm) {
417     blmP = (TAO_BLMVM *)tao->data;
418     M = blmP->M;
419   } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "This routine applies to TAO_LMVM and TAO_BLMVM.");
420   ierr = MatLMVMGetJ0(M, H0);CHKERRQ(ierr);
421   PetscFunctionReturn(0);
422 }
423 
424 PETSC_EXTERN PetscErrorCode TaoBLMVMGetH0KSP(Tao tao, KSP *ksp)
425 {
426   TAO_LMVM       *lmP;
427   TAO_BLMVM      *blmP;
428   TaoType        type;
429   PetscBool      is_lmvm, is_blmvm;
430   Mat            M;
431   PetscErrorCode ierr;
432 
433   ierr = TaoGetType(tao, &type);CHKERRQ(ierr);
434   ierr = PetscStrcmp(type, TAOLMVM,  &is_lmvm);CHKERRQ(ierr);
435   ierr = PetscStrcmp(type, TAOBLMVM, &is_blmvm);CHKERRQ(ierr);
436 
437   if (is_lmvm) {
438     lmP = (TAO_LMVM *)tao->data;
439     M = lmP->M;
440   } else if (is_blmvm) {
441     blmP = (TAO_BLMVM *)tao->data;
442     M = blmP->M;
443   } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "This routine applies to TAO_LMVM and TAO_BLMVM.");
444   ierr = MatLMVMGetJ0KSP(M, ksp);CHKERRQ(ierr);
445   PetscFunctionReturn(0);
446 }
447