xref: /petsc/src/tao/bound/impls/blmvm/blmvm.c (revision 5d5277661bc9d5de49c003e79a2b7124bf8a2eb4)
1ba92ff59SBarry Smith #include <petsctaolinesearch.h>
2aaa7dc30SBarry Smith #include <../src/tao/matrix/lmvmmat.h>
3aaa7dc30SBarry Smith #include <../src/tao/bound/impls/blmvm/blmvm.h>
4a7e14dcfSSatish Balay 
5a7e14dcfSSatish Balay /*------------------------------------------------------------*/
6a7e14dcfSSatish Balay #undef __FUNCT__
7a7e14dcfSSatish Balay #define __FUNCT__ "TaoSolve_BLMVM"
8441846f8SBarry Smith static PetscErrorCode TaoSolve_BLMVM(Tao tao)
9a7e14dcfSSatish Balay {
10a7e14dcfSSatish Balay   PetscErrorCode               ierr;
11a7e14dcfSSatish Balay   TAO_BLMVM                    *blmP = (TAO_BLMVM *)tao->data;
12e4cb33bbSBarry Smith   TaoConvergedReason           reason = TAO_CONTINUE_ITERATING;
13e4cb33bbSBarry Smith   TaoLineSearchConvergedReason ls_status = TAOLINESEARCH_CONTINUE_ITERATING;
14a7e14dcfSSatish Balay   PetscReal                    f, fold, gdx, gnorm;
15a7e14dcfSSatish Balay   PetscReal                    stepsize = 1.0,delta;
16a7e14dcfSSatish Balay 
17a7e14dcfSSatish Balay   PetscFunctionBegin;
18a7e14dcfSSatish Balay   /*  Project initial point onto bounds */
19a7e14dcfSSatish Balay   ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr);
20a7e14dcfSSatish Balay   ierr = VecMedian(tao->XL,tao->solution,tao->XU,tao->solution);CHKERRQ(ierr);
21a7e14dcfSSatish Balay   ierr = TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU);CHKERRQ(ierr);
22a7e14dcfSSatish Balay 
23a7e14dcfSSatish Balay   /* Check convergence criteria */
24a7e14dcfSSatish Balay   ierr = TaoComputeObjectiveAndGradient(tao, tao->solution,&f,blmP->unprojected_gradient);CHKERRQ(ierr);
25a7e14dcfSSatish Balay   ierr = VecBoundGradientProjection(blmP->unprojected_gradient,tao->solution, tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr);
26a7e14dcfSSatish Balay 
27a7e14dcfSSatish Balay   ierr = VecNorm(tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr);
2853506e15SBarry Smith   if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf pr NaN");
29a7e14dcfSSatish Balay 
308931d482SJason Sarich   ierr = TaoMonitor(tao, tao->niter, f, gnorm, 0.0, stepsize, &reason);CHKERRQ(ierr);
3153506e15SBarry Smith   if (reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
32a7e14dcfSSatish Balay 
33a7e14dcfSSatish Balay   /* Set initial scaling for the function */
34a7e14dcfSSatish Balay   if (f != 0.0) {
35a7e14dcfSSatish Balay     delta = 2.0*PetscAbsScalar(f) / (gnorm*gnorm);
3653506e15SBarry Smith   } else {
37a7e14dcfSSatish Balay     delta = 2.0 / (gnorm*gnorm);
38a7e14dcfSSatish Balay   }
39a7e14dcfSSatish Balay   ierr = MatLMVMSetDelta(blmP->M,delta);CHKERRQ(ierr);
40a7e14dcfSSatish Balay 
41a7e14dcfSSatish Balay   /* Set counter for gradient/reset steps */
42a7e14dcfSSatish Balay   blmP->grad = 0;
43a7e14dcfSSatish Balay   blmP->reset = 0;
44a7e14dcfSSatish Balay 
45a7e14dcfSSatish Balay   /* Have not converged; continue with Newton method */
46a7e14dcfSSatish Balay   while (reason == TAO_CONTINUE_ITERATING) {
47a7e14dcfSSatish Balay     /* Compute direction */
48a7e14dcfSSatish Balay     ierr = MatLMVMUpdate(blmP->M, tao->solution, tao->gradient);CHKERRQ(ierr);
49a7e14dcfSSatish Balay     ierr = MatLMVMSolve(blmP->M, blmP->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
50a7e14dcfSSatish Balay     ierr = VecBoundGradientProjection(tao->stepdirection,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr);
51a7e14dcfSSatish Balay 
52a7e14dcfSSatish Balay     /* Check for success (descent direction) */
53a7e14dcfSSatish Balay     ierr = VecDot(blmP->unprojected_gradient, tao->gradient, &gdx);CHKERRQ(ierr);
54a7e14dcfSSatish Balay     if (gdx <= 0) {
55a7e14dcfSSatish Balay       /* Step is not descent or solve was not successful
56a7e14dcfSSatish Balay          Use steepest descent direction (scaled) */
57a7e14dcfSSatish Balay       ++blmP->grad;
58a7e14dcfSSatish Balay 
59a7e14dcfSSatish Balay       if (f != 0.0) {
60a7e14dcfSSatish Balay         delta = 2.0*PetscAbsScalar(f) / (gnorm*gnorm);
6153506e15SBarry Smith       } else {
62a7e14dcfSSatish Balay         delta = 2.0 / (gnorm*gnorm);
63a7e14dcfSSatish Balay       }
64a7e14dcfSSatish Balay       ierr = MatLMVMSetDelta(blmP->M,delta);CHKERRQ(ierr);
65a7e14dcfSSatish Balay       ierr = MatLMVMReset(blmP->M);CHKERRQ(ierr);
66a7e14dcfSSatish Balay       ierr = MatLMVMUpdate(blmP->M, tao->solution, blmP->unprojected_gradient);CHKERRQ(ierr);
67a7e14dcfSSatish Balay       ierr = MatLMVMSolve(blmP->M,blmP->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
68a7e14dcfSSatish Balay     }
69a7e14dcfSSatish Balay     ierr = VecScale(tao->stepdirection,-1.0);CHKERRQ(ierr);
70a7e14dcfSSatish Balay 
71a7e14dcfSSatish Balay     /* Perform the linesearch */
72a7e14dcfSSatish Balay     fold = f;
73a7e14dcfSSatish Balay     ierr = VecCopy(tao->solution, blmP->Xold);CHKERRQ(ierr);
74a7e14dcfSSatish Balay     ierr = VecCopy(blmP->unprojected_gradient, blmP->Gold);CHKERRQ(ierr);
75a7e14dcfSSatish Balay     ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,1.0);CHKERRQ(ierr);
76a7e14dcfSSatish Balay     ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, blmP->unprojected_gradient, tao->stepdirection, &stepsize, &ls_status);CHKERRQ(ierr);
77a7e14dcfSSatish Balay     ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
78a7e14dcfSSatish Balay 
79a7e14dcfSSatish Balay     if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) {
80a7e14dcfSSatish Balay       /* Linesearch failed
81a7e14dcfSSatish Balay          Reset factors and use scaled (projected) gradient step */
82a7e14dcfSSatish Balay       ++blmP->reset;
83a7e14dcfSSatish Balay 
84a7e14dcfSSatish Balay       f = fold;
85a7e14dcfSSatish Balay       ierr = VecCopy(blmP->Xold, tao->solution);CHKERRQ(ierr);
86a7e14dcfSSatish Balay       ierr = VecCopy(blmP->Gold, blmP->unprojected_gradient);CHKERRQ(ierr);
87a7e14dcfSSatish Balay 
88a7e14dcfSSatish Balay       if (f != 0.0) {
89a7e14dcfSSatish Balay         delta = 2.0* PetscAbsScalar(f) / (gnorm*gnorm);
9053506e15SBarry Smith       } else {
91a7e14dcfSSatish Balay         delta = 2.0/ (gnorm*gnorm);
92a7e14dcfSSatish Balay       }
93a7e14dcfSSatish Balay       ierr = MatLMVMSetDelta(blmP->M,delta);CHKERRQ(ierr);
94a7e14dcfSSatish Balay       ierr = MatLMVMReset(blmP->M);CHKERRQ(ierr);
95a7e14dcfSSatish Balay       ierr = MatLMVMUpdate(blmP->M, tao->solution, blmP->unprojected_gradient);CHKERRQ(ierr);
96a7e14dcfSSatish Balay       ierr = MatLMVMSolve(blmP->M, blmP->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
97a7e14dcfSSatish Balay       ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
98a7e14dcfSSatish Balay 
99a7e14dcfSSatish Balay       /* This may be incorrect; linesearch has values fo stepmax and stepmin
100a7e14dcfSSatish Balay          that should be reset. */
101302440fdSBarry Smith       ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,1.0);CHKERRQ(ierr);
102a7e14dcfSSatish Balay       ierr = TaoLineSearchApply(tao->linesearch,tao->solution,&f, blmP->unprojected_gradient, tao->stepdirection,  &stepsize, &ls_status);CHKERRQ(ierr);
103a7e14dcfSSatish Balay       ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
104a7e14dcfSSatish Balay 
105a7e14dcfSSatish Balay       if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) {
106a7e14dcfSSatish Balay         tao->reason = TAO_DIVERGED_LS_FAILURE;
107a7e14dcfSSatish Balay         break;
108a7e14dcfSSatish Balay       }
109a7e14dcfSSatish Balay     }
110a7e14dcfSSatish Balay 
111e4cb33bbSBarry Smith     /* Check for converged */
112a7e14dcfSSatish Balay     ierr = VecBoundGradientProjection(blmP->unprojected_gradient, tao->solution, tao->XL, tao->XU, tao->gradient);CHKERRQ(ierr);
113a7e14dcfSSatish Balay     ierr = VecNorm(tao->gradient, NORM_2, &gnorm);CHKERRQ(ierr);
114a7e14dcfSSatish Balay 
115a7e14dcfSSatish Balay 
11653506e15SBarry Smith     if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Not-a-Number");
1178931d482SJason Sarich     tao->niter++;
1188931d482SJason Sarich     ierr = TaoMonitor(tao, tao->niter, f, gnorm, 0.0, stepsize, &reason);CHKERRQ(ierr);
119a7e14dcfSSatish Balay   }
120a7e14dcfSSatish Balay   PetscFunctionReturn(0);
121a7e14dcfSSatish Balay }
122a7e14dcfSSatish Balay 
123a7e14dcfSSatish Balay #undef __FUNCT__
124a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetup_BLMVM"
125441846f8SBarry Smith static PetscErrorCode TaoSetup_BLMVM(Tao tao)
126a7e14dcfSSatish Balay {
127a7e14dcfSSatish Balay   TAO_BLMVM      *blmP = (TAO_BLMVM *)tao->data;
128a7e14dcfSSatish Balay   PetscInt       n,N;
129a7e14dcfSSatish Balay   PetscErrorCode ierr;
130a7e14dcfSSatish Balay 
131a7e14dcfSSatish Balay   PetscFunctionBegin;
132a7e14dcfSSatish Balay   /* Existence of tao->solution checked in TaoSetup() */
133a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution,&blmP->Xold);CHKERRQ(ierr);
134a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution,&blmP->Gold);CHKERRQ(ierr);
135302440fdSBarry Smith   ierr = VecDuplicate(tao->solution, &blmP->unprojected_gradient);CHKERRQ(ierr);
136a7e14dcfSSatish Balay 
137a7e14dcfSSatish Balay   if (!tao->stepdirection) {
13853506e15SBarry Smith     ierr = VecDuplicate(tao->solution, &tao->stepdirection);CHKERRQ(ierr);
139a7e14dcfSSatish Balay   }
140a7e14dcfSSatish Balay   if (!tao->gradient) {
141a7e14dcfSSatish Balay     ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);
142a7e14dcfSSatish Balay   }
143a7e14dcfSSatish Balay   if (!tao->XL) {
144a7e14dcfSSatish Balay     ierr = VecDuplicate(tao->solution,&tao->XL);CHKERRQ(ierr);
145e270355aSBarry Smith     ierr = VecSet(tao->XL,PETSC_NINFINITY);CHKERRQ(ierr);
146a7e14dcfSSatish Balay   }
147a7e14dcfSSatish Balay   if (!tao->XU) {
148a7e14dcfSSatish Balay     ierr = VecDuplicate(tao->solution,&tao->XU);CHKERRQ(ierr);
149e270355aSBarry Smith     ierr = VecSet(tao->XU,PETSC_INFINITY);CHKERRQ(ierr);
150a7e14dcfSSatish Balay   }
151a7e14dcfSSatish Balay   /* Create matrix for the limited memory approximation */
152a7e14dcfSSatish Balay   ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr);
153a7e14dcfSSatish Balay   ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr);
154a7e14dcfSSatish Balay   ierr = MatCreateLMVM(((PetscObject)tao)->comm,n,N,&blmP->M);CHKERRQ(ierr);
155a7e14dcfSSatish Balay   ierr = MatLMVMAllocateVectors(blmP->M,tao->solution);CHKERRQ(ierr);
156a7e14dcfSSatish Balay   PetscFunctionReturn(0);
157a7e14dcfSSatish Balay }
158a7e14dcfSSatish Balay 
159a7e14dcfSSatish Balay /* ---------------------------------------------------------- */
160a7e14dcfSSatish Balay #undef __FUNCT__
161a7e14dcfSSatish Balay #define __FUNCT__ "TaoDestroy_BLMVM"
162441846f8SBarry Smith static PetscErrorCode TaoDestroy_BLMVM(Tao tao)
163a7e14dcfSSatish Balay {
164a7e14dcfSSatish Balay   TAO_BLMVM      *blmP = (TAO_BLMVM *)tao->data;
165a7e14dcfSSatish Balay   PetscErrorCode ierr;
166a7e14dcfSSatish Balay 
167a7e14dcfSSatish Balay   PetscFunctionBegin;
168a7e14dcfSSatish Balay   if (tao->setupcalled) {
169a7e14dcfSSatish Balay     ierr = MatDestroy(&blmP->M);CHKERRQ(ierr);
170a7e14dcfSSatish Balay     ierr = VecDestroy(&blmP->unprojected_gradient);CHKERRQ(ierr);
171a7e14dcfSSatish Balay     ierr = VecDestroy(&blmP->Xold);CHKERRQ(ierr);
172a7e14dcfSSatish Balay     ierr = VecDestroy(&blmP->Gold);CHKERRQ(ierr);
173a7e14dcfSSatish Balay   }
174a7e14dcfSSatish Balay   ierr = PetscFree(tao->data);CHKERRQ(ierr);
175a7e14dcfSSatish Balay   PetscFunctionReturn(0);
176a7e14dcfSSatish Balay }
177a7e14dcfSSatish Balay 
178a7e14dcfSSatish Balay /*------------------------------------------------------------*/
179a7e14dcfSSatish Balay #undef __FUNCT__
180a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetFromOptions_BLMVM"
1818c34d3f5SBarry Smith static PetscErrorCode TaoSetFromOptions_BLMVM(PetscOptions* PetscOptionsObject,Tao tao)
182a7e14dcfSSatish Balay {
183a7e14dcfSSatish Balay   PetscErrorCode ierr;
184a7e14dcfSSatish Balay 
185a7e14dcfSSatish Balay   PetscFunctionBegin;
1861a1499c8SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Limited-memory variable-metric method for bound constrained optimization");CHKERRQ(ierr);
187a7e14dcfSSatish Balay   ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr);
188a7e14dcfSSatish Balay   ierr = PetscOptionsTail();CHKERRQ(ierr);
189a7e14dcfSSatish Balay   PetscFunctionReturn(0);
190a7e14dcfSSatish Balay }
191a7e14dcfSSatish Balay 
192a7e14dcfSSatish Balay 
193a7e14dcfSSatish Balay /*------------------------------------------------------------*/
194a7e14dcfSSatish Balay #undef __FUNCT__
195a7e14dcfSSatish Balay #define __FUNCT__ "TaoView_BLMVM"
196441846f8SBarry Smith static int TaoView_BLMVM(Tao tao, PetscViewer viewer)
197a7e14dcfSSatish Balay {
198a7e14dcfSSatish Balay   TAO_BLMVM      *lmP = (TAO_BLMVM *)tao->data;
199a7e14dcfSSatish Balay   PetscBool      isascii;
200a7e14dcfSSatish Balay   PetscErrorCode ierr;
201a7e14dcfSSatish Balay 
202a7e14dcfSSatish Balay   PetscFunctionBegin;
203a7e14dcfSSatish Balay   ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii);CHKERRQ(ierr);
204a7e14dcfSSatish Balay   if (isascii) {
205a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
206a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", lmP->grad);CHKERRQ(ierr);
207a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
208a7e14dcfSSatish Balay   }
209a7e14dcfSSatish Balay   PetscFunctionReturn(0);
210a7e14dcfSSatish Balay }
211a7e14dcfSSatish Balay 
212a7e14dcfSSatish Balay #undef __FUNCT__
213a7e14dcfSSatish Balay #define __FUNCT__ "TaoComputeDual_BLMVM"
214441846f8SBarry Smith static PetscErrorCode TaoComputeDual_BLMVM(Tao tao, Vec DXL, Vec DXU)
215a7e14dcfSSatish Balay {
216a7e14dcfSSatish Balay   TAO_BLMVM      *blm = (TAO_BLMVM *) tao->data;
217a7e14dcfSSatish Balay   PetscErrorCode ierr;
218a7e14dcfSSatish Balay 
219a7e14dcfSSatish Balay   PetscFunctionBegin;
220441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
221a7e14dcfSSatish Balay   PetscValidHeaderSpecific(DXL,VEC_CLASSID,2);
222a7e14dcfSSatish Balay   PetscValidHeaderSpecific(DXU,VEC_CLASSID,3);
22353506e15SBarry Smith   if (!tao->gradient || !blm->unprojected_gradient) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Dual variables don't exist yet or no longer exist.\n");
224a7e14dcfSSatish Balay 
225a7e14dcfSSatish Balay   ierr = VecCopy(tao->gradient,DXL);CHKERRQ(ierr);
226a7e14dcfSSatish Balay   ierr = VecAXPY(DXL,-1.0,blm->unprojected_gradient);CHKERRQ(ierr);
227a7e14dcfSSatish Balay   ierr = VecSet(DXU,0.0);CHKERRQ(ierr);
228a7e14dcfSSatish Balay   ierr = VecPointwiseMax(DXL,DXL,DXU);CHKERRQ(ierr);
229a7e14dcfSSatish Balay 
230a7e14dcfSSatish Balay   ierr = VecCopy(blm->unprojected_gradient,DXU);CHKERRQ(ierr);
231a7e14dcfSSatish Balay   ierr = VecAXPY(DXU,-1.0,tao->gradient);CHKERRQ(ierr);
232a7e14dcfSSatish Balay   ierr = VecAXPY(DXU,1.0,DXL);CHKERRQ(ierr);
233a7e14dcfSSatish Balay   PetscFunctionReturn(0);
234a7e14dcfSSatish Balay }
235a7e14dcfSSatish Balay 
236a7e14dcfSSatish Balay /* ---------------------------------------------------------- */
2371522df2eSJason Sarich /*MC
2381522df2eSJason Sarich   TAOBLMVM - Bounded limited memory variable metric is a quasi-Newton method
2391522df2eSJason Sarich          for nonlinear minimization with bound constraints. It is an extension
2401522df2eSJason Sarich          of TAOLMVM
2411522df2eSJason Sarich 
2421522df2eSJason Sarich   Options Database Keys:
2431522df2eSJason Sarich +     -tao_lmm_vectors - number of vectors to use for approximation
2441522df2eSJason Sarich .     -tao_lmm_scale_type - "none","scalar","broyden"
2451522df2eSJason Sarich .     -tao_lmm_limit_type - "none","average","relative","absolute"
2461522df2eSJason Sarich .     -tao_lmm_rescale_type - "none","scalar","gl"
2471522df2eSJason Sarich .     -tao_lmm_limit_mu - mu limiting factor
2481522df2eSJason Sarich .     -tao_lmm_limit_nu - nu limiting factor
2491522df2eSJason Sarich .     -tao_lmm_delta_min - minimum delta value
2501522df2eSJason Sarich .     -tao_lmm_delta_max - maximum delta value
2511522df2eSJason Sarich .     -tao_lmm_broyden_phi - phi factor for Broyden scaling
2521522df2eSJason Sarich .     -tao_lmm_scalar_alpha - alpha factor for scalar scaling
2531522df2eSJason Sarich .     -tao_lmm_rescale_alpha - alpha factor for rescaling diagonal
2541522df2eSJason Sarich .     -tao_lmm_rescale_beta - beta factor for rescaling diagonal
2551522df2eSJason Sarich .     -tao_lmm_scalar_history - amount of history for scalar scaling
2561522df2eSJason Sarich .     -tao_lmm_rescale_history - amount of history for rescaling diagonal
2571522df2eSJason Sarich -     -tao_lmm_eps - rejection tolerance
2581522df2eSJason Sarich 
2591eb8069cSJason Sarich   Level: beginner
2601522df2eSJason Sarich M*/
261a7e14dcfSSatish Balay #undef __FUNCT__
262a7e14dcfSSatish Balay #define __FUNCT__ "TaoCreate_BLMVM"
263728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate_BLMVM(Tao tao)
264a7e14dcfSSatish Balay {
265a7e14dcfSSatish Balay   TAO_BLMVM      *blmP;
2668caf6e8cSBarry Smith   const char     *morethuente_type = TAOLINESEARCHMT;
267a7e14dcfSSatish Balay   PetscErrorCode ierr;
268a7e14dcfSSatish Balay 
269a7e14dcfSSatish Balay   PetscFunctionBegin;
270a7e14dcfSSatish Balay   tao->ops->setup = TaoSetup_BLMVM;
271a7e14dcfSSatish Balay   tao->ops->solve = TaoSolve_BLMVM;
272a7e14dcfSSatish Balay   tao->ops->view = TaoView_BLMVM;
273a7e14dcfSSatish Balay   tao->ops->setfromoptions = TaoSetFromOptions_BLMVM;
274a7e14dcfSSatish Balay   tao->ops->destroy = TaoDestroy_BLMVM;
275a7e14dcfSSatish Balay   tao->ops->computedual = TaoComputeDual_BLMVM;
276a7e14dcfSSatish Balay 
2773c9e27cfSGeoffrey Irving   ierr = PetscNewLog(tao,&blmP);CHKERRQ(ierr);
278a7e14dcfSSatish Balay   tao->data = (void*)blmP;
279a7e14dcfSSatish Balay   tao->max_it = 2000;
280a7e14dcfSSatish Balay   tao->max_funcs = 4000;
281a7e14dcfSSatish Balay   tao->fatol = 1e-4;
282a7e14dcfSSatish Balay   tao->frtol = 1e-4;
283a7e14dcfSSatish Balay 
284a7e14dcfSSatish Balay   ierr = TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch);CHKERRQ(ierr);
285a7e14dcfSSatish Balay   ierr = TaoLineSearchSetType(tao->linesearch, morethuente_type);CHKERRQ(ierr);
286441846f8SBarry Smith   ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr);
287*5d527766SPatrick Farrell   ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr);
288a7e14dcfSSatish Balay   PetscFunctionReturn(0);
289a7e14dcfSSatish Balay }
290a7e14dcfSSatish Balay 
291