xref: /petsc/src/tao/bound/impls/blmvm/blmvm.c (revision 691b26d3a7ce3263bd9be9c446af0af2a46feecf)
1f5766c09SAlp Dener #include <petsctaolinesearch.h>
2f5766c09SAlp Dener #include <../src/tao/unconstrained/impls/lmvm/lmvm.h>
3f5766c09SAlp Dener #include <../src/tao/bound/impls/blmvm/blmvm.h>
4a7e14dcfSSatish Balay 
5f5766c09SAlp Dener /*------------------------------------------------------------*/
6f5766c09SAlp Dener static PetscErrorCode TaoSolve_BLMVM(Tao tao)
7f5766c09SAlp Dener {
8f5766c09SAlp Dener   PetscErrorCode               ierr;
9f5766c09SAlp Dener   TAO_BLMVM                    *blmP = (TAO_BLMVM *)tao->data;
10f5766c09SAlp Dener   TaoLineSearchConvergedReason ls_status = TAOLINESEARCH_CONTINUE_ITERATING;
11f5766c09SAlp Dener   PetscReal                    f, fold, gdx, gnorm, gnorm2;
12f5766c09SAlp Dener   PetscReal                    stepsize = 1.0,delta;
13a7e14dcfSSatish Balay 
14f5766c09SAlp Dener   PetscFunctionBegin;
15f5766c09SAlp Dener   /*  Project initial point onto bounds */
16f5766c09SAlp Dener   ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr);
17f5766c09SAlp Dener   ierr = VecMedian(tao->XL,tao->solution,tao->XU,tao->solution);CHKERRQ(ierr);
18f5766c09SAlp Dener   ierr = TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU);CHKERRQ(ierr);
19f5766c09SAlp Dener 
20f5766c09SAlp Dener 
21f5766c09SAlp Dener   /* Check convergence criteria */
22f5766c09SAlp Dener   ierr = TaoComputeObjectiveAndGradient(tao, tao->solution,&f,blmP->unprojected_gradient);CHKERRQ(ierr);
23f5766c09SAlp Dener   ierr = VecBoundGradientProjection(blmP->unprojected_gradient,tao->solution, tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr);
24f5766c09SAlp Dener 
25f5766c09SAlp Dener   ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr);
26*691b26d3SBarry Smith   if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_USER, "User provided compute function generated Inf or NaN");
27f5766c09SAlp Dener 
28f5766c09SAlp Dener   tao->reason = TAO_CONTINUE_ITERATING;
29f5766c09SAlp Dener   ierr = TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
30f5766c09SAlp Dener   ierr = TaoMonitor(tao,tao->niter,f,gnorm,0.0,stepsize);CHKERRQ(ierr);
31f5766c09SAlp Dener   ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
32f5766c09SAlp Dener   if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
33f5766c09SAlp Dener 
34f5766c09SAlp Dener   /* Set counter for gradient/reset steps */
35f5766c09SAlp Dener   if (!blmP->recycle) {
36f5766c09SAlp Dener     blmP->grad = 0;
37f5766c09SAlp Dener     blmP->reset = 0;
38f5766c09SAlp Dener     ierr = MatLMVMReset(blmP->M, PETSC_FALSE);CHKERRQ(ierr);
39f5766c09SAlp Dener   }
40f5766c09SAlp Dener 
41f5766c09SAlp Dener   /* Have not converged; continue with Newton method */
42f5766c09SAlp Dener   while (tao->reason == TAO_CONTINUE_ITERATING) {
43e1e80dc8SAlp Dener     /* Call general purpose update function */
44e1e80dc8SAlp Dener     if (tao->ops->update) {
458fcddce6SStefano Zampini       ierr = (*tao->ops->update)(tao, tao->niter, tao->user_update);CHKERRQ(ierr);
46e1e80dc8SAlp Dener     }
47f5766c09SAlp Dener     /* Compute direction */
48f5766c09SAlp Dener     gnorm2 = gnorm*gnorm;
498cabe928SAlp Dener     if (gnorm2 == 0.0) gnorm2 = PETSC_MACHINE_EPSILON;
508cabe928SAlp Dener     if (f == 0.0) {
518cabe928SAlp Dener       delta = 2.0 / gnorm2;
528cabe928SAlp Dener     } else {
538cabe928SAlp Dener       delta = 2.0 * PetscAbsScalar(f) / gnorm2;
548cabe928SAlp Dener     }
55f5766c09SAlp Dener     ierr = MatSymBrdnSetDelta(blmP->M, delta);CHKERRQ(ierr);
56f5766c09SAlp Dener     ierr = MatLMVMUpdate(blmP->M, tao->solution, tao->gradient);CHKERRQ(ierr);
579515a401SAlp Dener     ierr = MatSolve(blmP->M, blmP->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
58f5766c09SAlp Dener     ierr = VecBoundGradientProjection(tao->stepdirection,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr);
59f5766c09SAlp Dener 
60f5766c09SAlp Dener     /* Check for success (descent direction) */
61f5766c09SAlp Dener     ierr = VecDot(blmP->unprojected_gradient, tao->gradient, &gdx);CHKERRQ(ierr);
62f5766c09SAlp Dener     if (gdx <= 0) {
63f5766c09SAlp Dener       /* Step is not descent or solve was not successful
64f5766c09SAlp Dener          Use steepest descent direction (scaled) */
65f5766c09SAlp Dener       ++blmP->grad;
66f5766c09SAlp Dener 
67f5766c09SAlp Dener       ierr = MatLMVMReset(blmP->M, PETSC_FALSE);CHKERRQ(ierr);
68f5766c09SAlp Dener       ierr = MatLMVMUpdate(blmP->M, tao->solution, blmP->unprojected_gradient);CHKERRQ(ierr);
699515a401SAlp Dener       ierr = MatSolve(blmP->M,blmP->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
70f5766c09SAlp Dener     }
71f5766c09SAlp Dener     ierr = VecScale(tao->stepdirection,-1.0);CHKERRQ(ierr);
72f5766c09SAlp Dener 
73f5766c09SAlp Dener     /* Perform the linesearch */
74f5766c09SAlp Dener     fold = f;
75f5766c09SAlp Dener     ierr = VecCopy(tao->solution, blmP->Xold);CHKERRQ(ierr);
76f5766c09SAlp Dener     ierr = VecCopy(blmP->unprojected_gradient, blmP->Gold);CHKERRQ(ierr);
77f5766c09SAlp Dener     ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,1.0);CHKERRQ(ierr);
78f5766c09SAlp Dener     ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, blmP->unprojected_gradient, tao->stepdirection, &stepsize, &ls_status);CHKERRQ(ierr);
79f5766c09SAlp Dener     ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
80f5766c09SAlp Dener 
81f5766c09SAlp Dener     if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) {
82f5766c09SAlp Dener       /* Linesearch failed
83f5766c09SAlp Dener          Reset factors and use scaled (projected) gradient step */
84f5766c09SAlp Dener       ++blmP->reset;
85f5766c09SAlp Dener 
86f5766c09SAlp Dener       f = fold;
87f5766c09SAlp Dener       ierr = VecCopy(blmP->Xold, tao->solution);CHKERRQ(ierr);
88f5766c09SAlp Dener       ierr = VecCopy(blmP->Gold, blmP->unprojected_gradient);CHKERRQ(ierr);
89f5766c09SAlp Dener 
90f5766c09SAlp Dener       ierr = MatLMVMReset(blmP->M, PETSC_FALSE);CHKERRQ(ierr);
91f5766c09SAlp Dener       ierr = MatLMVMUpdate(blmP->M, tao->solution, blmP->unprojected_gradient);CHKERRQ(ierr);
929515a401SAlp Dener       ierr = MatSolve(blmP->M, blmP->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
93f5766c09SAlp Dener       ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
94f5766c09SAlp Dener 
95f5766c09SAlp Dener       /* This may be incorrect; linesearch has values for stepmax and stepmin
96f5766c09SAlp Dener          that should be reset. */
97f5766c09SAlp Dener       ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,1.0);CHKERRQ(ierr);
98f5766c09SAlp Dener       ierr = TaoLineSearchApply(tao->linesearch,tao->solution,&f, blmP->unprojected_gradient, tao->stepdirection,  &stepsize, &ls_status);CHKERRQ(ierr);
99f5766c09SAlp Dener       ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
100f5766c09SAlp Dener 
101f5766c09SAlp Dener       if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) {
102f5766c09SAlp Dener         tao->reason = TAO_DIVERGED_LS_FAILURE;
103f5766c09SAlp Dener         break;
104f5766c09SAlp Dener       }
105f5766c09SAlp Dener     }
106f5766c09SAlp Dener 
107f5766c09SAlp Dener     /* Check for converged */
108f5766c09SAlp Dener     ierr = VecBoundGradientProjection(blmP->unprojected_gradient, tao->solution, tao->XL, tao->XU, tao->gradient);CHKERRQ(ierr);
109f5766c09SAlp Dener     ierr = TaoGradientNorm(tao, tao->gradient, NORM_2, &gnorm);CHKERRQ(ierr);
110*691b26d3SBarry Smith     if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_USER, "User provided compute function generated Not-a-Number");
111f5766c09SAlp Dener     tao->niter++;
112f5766c09SAlp Dener     ierr = TaoLogConvergenceHistory(tao,f,gnorm,0.0,tao->ksp_its);CHKERRQ(ierr);
113f5766c09SAlp Dener     ierr = TaoMonitor(tao,tao->niter,f,gnorm,0.0,stepsize);CHKERRQ(ierr);
114f5766c09SAlp Dener     ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
115f5766c09SAlp Dener   }
116f5766c09SAlp Dener   PetscFunctionReturn(0);
117f5766c09SAlp Dener }
118f5766c09SAlp Dener 
119f5766c09SAlp Dener static PetscErrorCode TaoSetup_BLMVM(Tao tao)
120f5766c09SAlp Dener {
121f5766c09SAlp Dener   TAO_BLMVM      *blmP = (TAO_BLMVM *)tao->data;
122f5766c09SAlp Dener   PetscErrorCode ierr;
123f5766c09SAlp Dener 
124f5766c09SAlp Dener   PetscFunctionBegin;
125f5766c09SAlp Dener   /* Existence of tao->solution checked in TaoSetup() */
126f5766c09SAlp Dener   ierr = VecDuplicate(tao->solution,&blmP->Xold);CHKERRQ(ierr);
127f5766c09SAlp Dener   ierr = VecDuplicate(tao->solution,&blmP->Gold);CHKERRQ(ierr);
128f5766c09SAlp Dener   ierr = VecDuplicate(tao->solution, &blmP->unprojected_gradient);CHKERRQ(ierr);
129f5766c09SAlp Dener 
130f5766c09SAlp Dener   if (!tao->stepdirection) {
131f5766c09SAlp Dener     ierr = VecDuplicate(tao->solution, &tao->stepdirection);CHKERRQ(ierr);
132f5766c09SAlp Dener   }
133f5766c09SAlp Dener   if (!tao->gradient) {
134f5766c09SAlp Dener     ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);
135f5766c09SAlp Dener   }
136f5766c09SAlp Dener   if (!tao->XL) {
137f5766c09SAlp Dener     ierr = VecDuplicate(tao->solution,&tao->XL);CHKERRQ(ierr);
138f5766c09SAlp Dener     ierr = VecSet(tao->XL,PETSC_NINFINITY);CHKERRQ(ierr);
139f5766c09SAlp Dener   }
140f5766c09SAlp Dener   if (!tao->XU) {
141f5766c09SAlp Dener     ierr = VecDuplicate(tao->solution,&tao->XU);CHKERRQ(ierr);
142f5766c09SAlp Dener     ierr = VecSet(tao->XU,PETSC_INFINITY);CHKERRQ(ierr);
143f5766c09SAlp Dener   }
144f5766c09SAlp Dener   /* Allocate matrix for the limited memory approximation */
145f5766c09SAlp Dener   ierr = MatLMVMAllocate(blmP->M,tao->solution,blmP->unprojected_gradient);CHKERRQ(ierr);
146f5766c09SAlp Dener 
147f5766c09SAlp Dener   /* If the user has set a matrix to solve as the initial H0, set the options prefix here, and set up the KSP */
148f5766c09SAlp Dener   if (blmP->H0) {
149f5766c09SAlp Dener     ierr = MatLMVMSetJ0(blmP->M, blmP->H0);CHKERRQ(ierr);
150f5766c09SAlp Dener   }
151f5766c09SAlp Dener   PetscFunctionReturn(0);
152f5766c09SAlp Dener }
153f5766c09SAlp Dener 
154f5766c09SAlp Dener /* ---------------------------------------------------------- */
155f5766c09SAlp Dener static PetscErrorCode TaoDestroy_BLMVM(Tao tao)
156f5766c09SAlp Dener {
157f5766c09SAlp Dener   TAO_BLMVM      *blmP = (TAO_BLMVM *)tao->data;
158f5766c09SAlp Dener   PetscErrorCode ierr;
159f5766c09SAlp Dener 
160f5766c09SAlp Dener   PetscFunctionBegin;
161f5766c09SAlp Dener   if (tao->setupcalled) {
162f5766c09SAlp Dener     ierr = VecDestroy(&blmP->unprojected_gradient);CHKERRQ(ierr);
163f5766c09SAlp Dener     ierr = VecDestroy(&blmP->Xold);CHKERRQ(ierr);
164f5766c09SAlp Dener     ierr = VecDestroy(&blmP->Gold);CHKERRQ(ierr);
165f5766c09SAlp Dener   }
166f5766c09SAlp Dener   ierr = MatDestroy(&blmP->M);CHKERRQ(ierr);
167f5766c09SAlp Dener   if (blmP->H0) {
168f5766c09SAlp Dener     PetscObjectDereference((PetscObject)blmP->H0);
169f5766c09SAlp Dener   }
170f5766c09SAlp Dener   ierr = PetscFree(tao->data);CHKERRQ(ierr);
171f5766c09SAlp Dener   PetscFunctionReturn(0);
172f5766c09SAlp Dener }
173f5766c09SAlp Dener 
174f5766c09SAlp Dener /*------------------------------------------------------------*/
1752ec5c1acSAlp Dener static PetscErrorCode TaoSetFromOptions_BLMVM(PetscOptionItems* PetscOptionsObject,Tao tao)
176a7e14dcfSSatish Balay {
177f5766c09SAlp Dener   TAO_BLMVM      *blmP = (TAO_BLMVM *)tao->data;
178a7e14dcfSSatish Balay   PetscErrorCode ierr;
179e5fecd4eSAlp Dener   PetscBool      is_spd;
180a7e14dcfSSatish Balay 
181a7e14dcfSSatish Balay   PetscFunctionBegin;
182f5766c09SAlp Dener   ierr = PetscOptionsHead(PetscOptionsObject,"Limited-memory variable-metric method for bound constrained optimization");CHKERRQ(ierr);
183f5766c09SAlp Dener   ierr = PetscOptionsBool("-tao_blmvm_recycle","enable recycling of the BFGS matrix between subsequent TaoSolve() calls","",blmP->recycle,&blmP->recycle,NULL);CHKERRQ(ierr);
184e5fecd4eSAlp Dener   ierr = PetscOptionsTail();CHKERRQ(ierr);
185a7e14dcfSSatish Balay   ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr);
186f5766c09SAlp Dener   ierr = MatSetFromOptions(blmP->M);CHKERRQ(ierr);
187f5766c09SAlp Dener   ierr = MatGetOption(blmP->M, MAT_SPD, &is_spd);CHKERRQ(ierr);
188e5fecd4eSAlp Dener   if (!is_spd) SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "LMVM matrix must be symmetric positive-definite");
189a7e14dcfSSatish Balay   PetscFunctionReturn(0);
190a7e14dcfSSatish Balay }
191a7e14dcfSSatish Balay 
192f5766c09SAlp Dener 
193f5766c09SAlp Dener /*------------------------------------------------------------*/
1945bd1e576SStefano Zampini static PetscErrorCode TaoView_BLMVM(Tao tao, PetscViewer viewer)
195a7e14dcfSSatish Balay {
196f5766c09SAlp Dener   TAO_BLMVM      *lmP = (TAO_BLMVM *)tao->data;
197f5766c09SAlp Dener   PetscBool      isascii;
198a7e14dcfSSatish Balay   PetscErrorCode ierr;
199a7e14dcfSSatish Balay 
200a7e14dcfSSatish Balay   PetscFunctionBegin;
201f5766c09SAlp Dener   ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii);CHKERRQ(ierr);
202f5766c09SAlp Dener   if (isascii) {
203f5766c09SAlp Dener     ierr = PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", lmP->grad);CHKERRQ(ierr);
204f5766c09SAlp Dener     ierr = PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
205f5766c09SAlp Dener     ierr = MatView(lmP->M, viewer);CHKERRQ(ierr);
206f5766c09SAlp Dener     ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
207f5766c09SAlp Dener   }
208f5766c09SAlp Dener   PetscFunctionReturn(0);
209f5766c09SAlp Dener }
210f5766c09SAlp Dener 
211f5766c09SAlp Dener static PetscErrorCode TaoComputeDual_BLMVM(Tao tao, Vec DXL, Vec DXU)
212f5766c09SAlp Dener {
213f5766c09SAlp Dener   TAO_BLMVM      *blm = (TAO_BLMVM *) tao->data;
214f5766c09SAlp Dener   PetscErrorCode ierr;
215f5766c09SAlp Dener 
216f5766c09SAlp Dener   PetscFunctionBegin;
217f5766c09SAlp Dener   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
218f5766c09SAlp Dener   PetscValidHeaderSpecific(DXL,VEC_CLASSID,2);
219f5766c09SAlp Dener   PetscValidHeaderSpecific(DXU,VEC_CLASSID,3);
220f5766c09SAlp Dener   if (!tao->gradient || !blm->unprojected_gradient) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Dual variables don't exist yet or no longer exist.\n");
221f5766c09SAlp Dener 
222f5766c09SAlp Dener   ierr = VecCopy(tao->gradient,DXL);CHKERRQ(ierr);
223f5766c09SAlp Dener   ierr = VecAXPY(DXL,-1.0,blm->unprojected_gradient);CHKERRQ(ierr);
224f5766c09SAlp Dener   ierr = VecSet(DXU,0.0);CHKERRQ(ierr);
225f5766c09SAlp Dener   ierr = VecPointwiseMax(DXL,DXL,DXU);CHKERRQ(ierr);
226f5766c09SAlp Dener 
227f5766c09SAlp Dener   ierr = VecCopy(blm->unprojected_gradient,DXU);CHKERRQ(ierr);
228f5766c09SAlp Dener   ierr = VecAXPY(DXU,-1.0,tao->gradient);CHKERRQ(ierr);
229f5766c09SAlp Dener   ierr = VecAXPY(DXU,1.0,DXL);CHKERRQ(ierr);
230f5766c09SAlp Dener   PetscFunctionReturn(0);
231f5766c09SAlp Dener }
232f5766c09SAlp Dener 
233f5766c09SAlp Dener /* ---------------------------------------------------------- */
234f5766c09SAlp Dener /*MC
235f5766c09SAlp Dener   TAOBLMVM - Bounded limited memory variable metric is a quasi-Newton method
236f5766c09SAlp Dener          for nonlinear minimization with bound constraints. It is an extension
237f5766c09SAlp Dener          of TAOLMVM
238f5766c09SAlp Dener 
239f5766c09SAlp Dener   Options Database Keys:
240f5766c09SAlp Dener .     -tao_lmm_recycle - enable recycling of LMVM information between subsequent TaoSolve calls
241f5766c09SAlp Dener 
242f5766c09SAlp Dener   Level: beginner
243f5766c09SAlp Dener M*/
244f5766c09SAlp Dener PETSC_EXTERN PetscErrorCode TaoCreate_BLMVM(Tao tao)
245f5766c09SAlp Dener {
246f5766c09SAlp Dener   TAO_BLMVM      *blmP;
247f5766c09SAlp Dener   const char     *morethuente_type = TAOLINESEARCHMT;
248f5766c09SAlp Dener   PetscErrorCode ierr;
249f5766c09SAlp Dener 
250f5766c09SAlp Dener   PetscFunctionBegin;
251f5766c09SAlp Dener   tao->ops->setup = TaoSetup_BLMVM;
252f5766c09SAlp Dener   tao->ops->solve = TaoSolve_BLMVM;
253f5766c09SAlp Dener   tao->ops->view = TaoView_BLMVM;
254a7e14dcfSSatish Balay   tao->ops->setfromoptions = TaoSetFromOptions_BLMVM;
255f5766c09SAlp Dener   tao->ops->destroy = TaoDestroy_BLMVM;
256f5766c09SAlp Dener   tao->ops->computedual = TaoComputeDual_BLMVM;
257f5766c09SAlp Dener 
258f5766c09SAlp Dener   ierr = PetscNewLog(tao,&blmP);CHKERRQ(ierr);
259f5766c09SAlp Dener   blmP->H0 = NULL;
260f5766c09SAlp Dener   blmP->recycle = PETSC_FALSE;
261f5766c09SAlp Dener   tao->data = (void*)blmP;
262f5766c09SAlp Dener 
263f5766c09SAlp Dener   /* Override default settings (unless already changed) */
264f5766c09SAlp Dener   if (!tao->max_it_changed) tao->max_it = 2000;
265f5766c09SAlp Dener   if (!tao->max_funcs_changed) tao->max_funcs = 4000;
266f5766c09SAlp Dener 
267f5766c09SAlp Dener   ierr = TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch);CHKERRQ(ierr);
268f5766c09SAlp Dener   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr);
269f5766c09SAlp Dener   ierr = TaoLineSearchSetType(tao->linesearch, morethuente_type);CHKERRQ(ierr);
270f5766c09SAlp Dener   ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr);
271f5766c09SAlp Dener   ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr);
272f5766c09SAlp Dener 
273f5766c09SAlp Dener   ierr = KSPInitializePackage();CHKERRQ(ierr);
274f5766c09SAlp Dener   ierr = MatCreate(((PetscObject)tao)->comm, &blmP->M);CHKERRQ(ierr);
275f5766c09SAlp Dener   ierr = MatSetType(blmP->M, MATLMVMBFGS);CHKERRQ(ierr);
276f5766c09SAlp Dener   ierr = PetscObjectIncrementTabLevel((PetscObject)blmP->M, (PetscObject)tao, 1);CHKERRQ(ierr);
277f5766c09SAlp Dener   ierr = MatSetOptionsPrefix(blmP->M, "tao_blmvm_");CHKERRQ(ierr);
278f5766c09SAlp Dener   PetscFunctionReturn(0);
279f5766c09SAlp Dener }
280f5766c09SAlp Dener 
281b39c12a9SAlp Dener PetscErrorCode TaoLMVMRecycle(Tao tao, PetscBool flg)
282b39c12a9SAlp Dener {
283b39c12a9SAlp Dener   TAO_LMVM       *lmP;
284b39c12a9SAlp Dener   TAO_BLMVM      *blmP;
285b39c12a9SAlp Dener   TaoType        type;
286b39c12a9SAlp Dener   PetscBool      is_lmvm, is_blmvm;
287b39c12a9SAlp Dener   PetscErrorCode ierr;
288b39c12a9SAlp Dener 
289b39c12a9SAlp Dener   PetscFunctionBegin;
290b39c12a9SAlp Dener   ierr = TaoGetType(tao, &type);CHKERRQ(ierr);
291b39c12a9SAlp Dener   ierr = PetscStrcmp(type, TAOLMVM,  &is_lmvm);CHKERRQ(ierr);
292b39c12a9SAlp Dener   ierr = PetscStrcmp(type, TAOBLMVM, &is_blmvm);CHKERRQ(ierr);
293b39c12a9SAlp Dener 
294b39c12a9SAlp Dener   if (is_lmvm) {
295b39c12a9SAlp Dener     lmP = (TAO_LMVM *)tao->data;
296b39c12a9SAlp Dener     lmP->recycle = flg;
297b39c12a9SAlp Dener   } else if (is_blmvm) {
298b39c12a9SAlp Dener     blmP = (TAO_BLMVM *)tao->data;
299b39c12a9SAlp Dener     blmP->recycle = flg;
300b39c12a9SAlp Dener   } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "This routine applies to TAO_LMVM and TAO_BLMVM.");
301b39c12a9SAlp Dener   PetscFunctionReturn(0);
302b39c12a9SAlp Dener }
303b39c12a9SAlp Dener 
304f5766c09SAlp Dener PetscErrorCode TaoLMVMSetH0(Tao tao, Mat H0)
305f5766c09SAlp Dener {
306f5766c09SAlp Dener   TAO_LMVM       *lmP;
307f5766c09SAlp Dener   TAO_BLMVM      *blmP;
308f5766c09SAlp Dener   TaoType        type;
309f5766c09SAlp Dener   PetscBool      is_lmvm, is_blmvm;
310f5766c09SAlp Dener   PetscErrorCode ierr;
311f5766c09SAlp Dener 
312b39c12a9SAlp Dener   PetscFunctionBegin;
313f5766c09SAlp Dener   ierr = TaoGetType(tao, &type);CHKERRQ(ierr);
314f5766c09SAlp Dener   ierr = PetscStrcmp(type, TAOLMVM,  &is_lmvm);CHKERRQ(ierr);
315f5766c09SAlp Dener   ierr = PetscStrcmp(type, TAOBLMVM, &is_blmvm);CHKERRQ(ierr);
316f5766c09SAlp Dener 
317f5766c09SAlp Dener   if (is_lmvm) {
318f5766c09SAlp Dener     lmP = (TAO_LMVM *)tao->data;
319f5766c09SAlp Dener     ierr = PetscObjectReference((PetscObject)H0);CHKERRQ(ierr);
320f5766c09SAlp Dener     lmP->H0 = H0;
321f5766c09SAlp Dener   } else if (is_blmvm) {
322f5766c09SAlp Dener     blmP = (TAO_BLMVM *)tao->data;
323f5766c09SAlp Dener     ierr = PetscObjectReference((PetscObject)H0);CHKERRQ(ierr);
324f5766c09SAlp Dener     blmP->H0 = H0;
325f5766c09SAlp Dener   } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "This routine applies to TAO_LMVM and TAO_BLMVM.");
326f5766c09SAlp Dener   PetscFunctionReturn(0);
327f5766c09SAlp Dener }
328f5766c09SAlp Dener 
329f5766c09SAlp Dener PetscErrorCode TaoLMVMGetH0(Tao tao, Mat *H0)
330f5766c09SAlp Dener {
331f5766c09SAlp Dener   TAO_LMVM       *lmP;
332f5766c09SAlp Dener   TAO_BLMVM      *blmP;
333f5766c09SAlp Dener   TaoType        type;
334f5766c09SAlp Dener   PetscBool      is_lmvm, is_blmvm;
335f5766c09SAlp Dener   Mat            M;
336f5766c09SAlp Dener 
337f5766c09SAlp Dener   PetscErrorCode ierr;
338f5766c09SAlp Dener 
339b39c12a9SAlp Dener   PetscFunctionBegin;
340f5766c09SAlp Dener   ierr = TaoGetType(tao, &type);CHKERRQ(ierr);
341f5766c09SAlp Dener   ierr = PetscStrcmp(type, TAOLMVM,  &is_lmvm);CHKERRQ(ierr);
342f5766c09SAlp Dener   ierr = PetscStrcmp(type, TAOBLMVM, &is_blmvm);CHKERRQ(ierr);
343f5766c09SAlp Dener 
344f5766c09SAlp Dener   if (is_lmvm) {
345f5766c09SAlp Dener     lmP = (TAO_LMVM *)tao->data;
346f5766c09SAlp Dener     M = lmP->M;
347f5766c09SAlp Dener   } else if (is_blmvm) {
348f5766c09SAlp Dener     blmP = (TAO_BLMVM *)tao->data;
349f5766c09SAlp Dener     M = blmP->M;
350f5766c09SAlp Dener   } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "This routine applies to TAO_LMVM and TAO_BLMVM.");
351f5766c09SAlp Dener   ierr = MatLMVMGetJ0(M, H0);CHKERRQ(ierr);
352f5766c09SAlp Dener   PetscFunctionReturn(0);
353f5766c09SAlp Dener }
354f5766c09SAlp Dener 
355f5766c09SAlp Dener PetscErrorCode TaoLMVMGetH0KSP(Tao tao, KSP *ksp)
356f5766c09SAlp Dener {
357f5766c09SAlp Dener   TAO_LMVM       *lmP;
358f5766c09SAlp Dener   TAO_BLMVM      *blmP;
359f5766c09SAlp Dener   TaoType        type;
360f5766c09SAlp Dener   PetscBool      is_lmvm, is_blmvm;
361f5766c09SAlp Dener   Mat            M;
362f5766c09SAlp Dener   PetscErrorCode ierr;
363f5766c09SAlp Dener 
364f5766c09SAlp Dener   ierr = TaoGetType(tao, &type);CHKERRQ(ierr);
365f5766c09SAlp Dener   ierr = PetscStrcmp(type, TAOLMVM,  &is_lmvm);CHKERRQ(ierr);
366f5766c09SAlp Dener   ierr = PetscStrcmp(type, TAOBLMVM, &is_blmvm);CHKERRQ(ierr);
367f5766c09SAlp Dener 
368f5766c09SAlp Dener   if (is_lmvm) {
369f5766c09SAlp Dener     lmP = (TAO_LMVM *)tao->data;
370f5766c09SAlp Dener     M = lmP->M;
371f5766c09SAlp Dener   } else if (is_blmvm) {
372f5766c09SAlp Dener     blmP = (TAO_BLMVM *)tao->data;
373f5766c09SAlp Dener     M = blmP->M;
374f5766c09SAlp Dener   } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "This routine applies to TAO_LMVM and TAO_BLMVM.");
375f5766c09SAlp Dener   ierr = MatLMVMGetJ0KSP(M, ksp);CHKERRQ(ierr);
376a9603a14SPatrick Farrell   PetscFunctionReturn(0);
377a9603a14SPatrick Farrell }
378