xref: /petsc/src/tao/bound/impls/blmvm/blmvm.c (revision 3c9e27cfca911a7d7e3219758be42726e83c4ab2)
1a7e14dcfSSatish Balay #include "taolinesearch.h"
2f89ca46fSSatish Balay #include "../src/tao/matrix/lmvmmat.h"
3a7e14dcfSSatish Balay #include "blmvm.h"
4a7e14dcfSSatish Balay 
5a7e14dcfSSatish Balay /*------------------------------------------------------------*/
6a7e14dcfSSatish Balay #undef __FUNCT__
7a7e14dcfSSatish Balay #define __FUNCT__ "TaoSolve_BLMVM"
8a7e14dcfSSatish Balay static PetscErrorCode TaoSolve_BLMVM(TaoSolver tao)
9a7e14dcfSSatish Balay {
10a7e14dcfSSatish Balay   PetscErrorCode ierr;
11a7e14dcfSSatish Balay   TAO_BLMVM *blmP = (TAO_BLMVM *)tao->data;
12a7e14dcfSSatish Balay 
13a7e14dcfSSatish Balay   TaoSolverTerminationReason reason = TAO_CONTINUE_ITERATING;
14a7e14dcfSSatish Balay   TaoLineSearchTerminationReason ls_status = TAOLINESEARCH_CONTINUE_ITERATING;
15a7e14dcfSSatish Balay 
16a7e14dcfSSatish Balay   PetscReal f, fold, gdx, gnorm;
17a7e14dcfSSatish Balay   PetscReal stepsize = 1.0,delta;
18a7e14dcfSSatish Balay 
19a7e14dcfSSatish Balay   PetscInt iter = 0;
20a7e14dcfSSatish Balay 
21a7e14dcfSSatish Balay 
22a7e14dcfSSatish Balay   PetscFunctionBegin;
23a7e14dcfSSatish Balay 
24a7e14dcfSSatish Balay   /*  Project initial point onto bounds */
25a7e14dcfSSatish Balay   ierr = TaoComputeVariableBounds(tao); CHKERRQ(ierr);
26a7e14dcfSSatish Balay   ierr = VecMedian(tao->XL,tao->solution,tao->XU,tao->solution); CHKERRQ(ierr);
27a7e14dcfSSatish Balay   ierr = TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU); CHKERRQ(ierr);
28a7e14dcfSSatish Balay 
29a7e14dcfSSatish Balay   /* Check convergence criteria */
30a7e14dcfSSatish Balay   ierr = TaoComputeObjectiveAndGradient(tao, tao->solution,&f,blmP->unprojected_gradient); CHKERRQ(ierr);
31a7e14dcfSSatish Balay   ierr = VecBoundGradientProjection(blmP->unprojected_gradient,tao->solution, tao->XL,tao->XU,tao->gradient); CHKERRQ(ierr);
32a7e14dcfSSatish Balay 
33a7e14dcfSSatish Balay 
34a7e14dcfSSatish Balay   ierr = VecNorm(tao->gradient,NORM_2,&gnorm); CHKERRQ(ierr);
35a7e14dcfSSatish Balay   if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) {
36a7e14dcfSSatish Balay     SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf pr NaN");
37a7e14dcfSSatish Balay   }
38a7e14dcfSSatish Balay 
39a7e14dcfSSatish Balay   ierr = TaoMonitor(tao, iter, f, gnorm, 0.0, stepsize, &reason); CHKERRQ(ierr);
40a7e14dcfSSatish Balay   if (reason != TAO_CONTINUE_ITERATING) {
41a7e14dcfSSatish Balay     PetscFunctionReturn(0);
42a7e14dcfSSatish Balay   }
43a7e14dcfSSatish Balay 
44a7e14dcfSSatish Balay   /* Set initial scaling for the function */
45a7e14dcfSSatish Balay   if (f != 0.0) {
46a7e14dcfSSatish Balay     delta = 2.0*PetscAbsScalar(f) / (gnorm*gnorm);
47a7e14dcfSSatish Balay   }
48a7e14dcfSSatish Balay   else {
49a7e14dcfSSatish Balay     delta = 2.0 / (gnorm*gnorm);
50a7e14dcfSSatish Balay   }
51a7e14dcfSSatish Balay   ierr = MatLMVMSetDelta(blmP->M,delta); CHKERRQ(ierr);
52a7e14dcfSSatish Balay 
53a7e14dcfSSatish Balay   /* Set counter for gradient/reset steps */
54a7e14dcfSSatish Balay   blmP->grad = 0;
55a7e14dcfSSatish Balay   blmP->reset = 0;
56a7e14dcfSSatish Balay 
57a7e14dcfSSatish Balay   /* Have not converged; continue with Newton method */
58a7e14dcfSSatish Balay   while (reason == TAO_CONTINUE_ITERATING) {
59a7e14dcfSSatish Balay 
60a7e14dcfSSatish Balay     /* Compute direction */
61a7e14dcfSSatish Balay     ierr = MatLMVMUpdate(blmP->M, tao->solution, tao->gradient); CHKERRQ(ierr);
62a7e14dcfSSatish Balay     ierr = MatLMVMSolve(blmP->M, blmP->unprojected_gradient, tao->stepdirection); CHKERRQ(ierr);
63a7e14dcfSSatish Balay     ierr = VecBoundGradientProjection(tao->stepdirection,tao->solution,tao->XL,tao->XU,tao->gradient); CHKERRQ(ierr);
64a7e14dcfSSatish Balay 
65a7e14dcfSSatish Balay     /* Check for success (descent direction) */
66a7e14dcfSSatish Balay     ierr = VecDot(blmP->unprojected_gradient, tao->gradient, &gdx); CHKERRQ(ierr);
67a7e14dcfSSatish Balay     if (gdx <= 0) {
68a7e14dcfSSatish Balay       /* Step is not descent or solve was not successful
69a7e14dcfSSatish Balay 	 Use steepest descent direction (scaled) */
70a7e14dcfSSatish Balay       ++blmP->grad;
71a7e14dcfSSatish Balay 
72a7e14dcfSSatish Balay       if (f != 0.0) {
73a7e14dcfSSatish Balay 	delta = 2.0*PetscAbsScalar(f) / (gnorm*gnorm);
74a7e14dcfSSatish Balay       }
75a7e14dcfSSatish Balay       else {
76a7e14dcfSSatish Balay 	delta = 2.0 / (gnorm*gnorm);
77a7e14dcfSSatish Balay       }
78a7e14dcfSSatish Balay       ierr = MatLMVMSetDelta(blmP->M,delta); CHKERRQ(ierr);
79a7e14dcfSSatish Balay       ierr = MatLMVMReset(blmP->M); CHKERRQ(ierr);
80a7e14dcfSSatish Balay       ierr = MatLMVMUpdate(blmP->M, tao->solution, blmP->unprojected_gradient); CHKERRQ(ierr);
81a7e14dcfSSatish Balay       ierr = MatLMVMSolve(blmP->M,blmP->unprojected_gradient, tao->stepdirection); CHKERRQ(ierr);
82a7e14dcfSSatish Balay     }
83a7e14dcfSSatish Balay     ierr = VecScale(tao->stepdirection,-1.0); CHKERRQ(ierr);
84a7e14dcfSSatish Balay 
85a7e14dcfSSatish Balay     /* Perform the linesearch */
86a7e14dcfSSatish Balay     fold = f;
87a7e14dcfSSatish Balay     ierr = VecCopy(tao->solution, blmP->Xold); CHKERRQ(ierr);
88a7e14dcfSSatish Balay     ierr = VecCopy(blmP->unprojected_gradient, blmP->Gold); CHKERRQ(ierr);
89a7e14dcfSSatish Balay     ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,1.0); CHKERRQ(ierr);
90a7e14dcfSSatish Balay     ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, blmP->unprojected_gradient, tao->stepdirection, &stepsize, &ls_status); CHKERRQ(ierr);
91a7e14dcfSSatish Balay     ierr = TaoAddLineSearchCounts(tao); CHKERRQ(ierr);
92a7e14dcfSSatish Balay 
93a7e14dcfSSatish Balay     if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) {
94a7e14dcfSSatish Balay       /* Linesearch failed
95a7e14dcfSSatish Balay 	 Reset factors and use scaled (projected) gradient step */
96a7e14dcfSSatish Balay       ++blmP->reset;
97a7e14dcfSSatish Balay 
98a7e14dcfSSatish Balay       f = fold;
99a7e14dcfSSatish Balay       ierr = VecCopy(blmP->Xold, tao->solution); CHKERRQ(ierr);
100a7e14dcfSSatish Balay       ierr = VecCopy(blmP->Gold, blmP->unprojected_gradient); CHKERRQ(ierr);
101a7e14dcfSSatish Balay 
102a7e14dcfSSatish Balay       if (f != 0.0) {
103a7e14dcfSSatish Balay 	delta = 2.0* PetscAbsScalar(f) / (gnorm*gnorm);
104a7e14dcfSSatish Balay       }
105a7e14dcfSSatish Balay       else {
106a7e14dcfSSatish Balay 	delta = 2.0/ (gnorm*gnorm);
107a7e14dcfSSatish Balay       }
108a7e14dcfSSatish Balay       ierr = MatLMVMSetDelta(blmP->M,delta); CHKERRQ(ierr);
109a7e14dcfSSatish Balay       ierr = MatLMVMReset(blmP->M); CHKERRQ(ierr);
110a7e14dcfSSatish Balay       ierr = MatLMVMUpdate(blmP->M, tao->solution, blmP->unprojected_gradient); CHKERRQ(ierr);
111a7e14dcfSSatish Balay       ierr = MatLMVMSolve(blmP->M, blmP->unprojected_gradient, tao->stepdirection); CHKERRQ(ierr);
112a7e14dcfSSatish Balay       ierr = VecScale(tao->stepdirection, -1.0); CHKERRQ(ierr);
113a7e14dcfSSatish Balay 
114a7e14dcfSSatish Balay       /* This may be incorrect; linesearch has values fo stepmax and stepmin
115a7e14dcfSSatish Balay 	 that should be reset. */
116a7e14dcfSSatish Balay       ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,1.0);
117a7e14dcfSSatish Balay       ierr = TaoLineSearchApply(tao->linesearch,tao->solution,&f, blmP->unprojected_gradient, tao->stepdirection,  &stepsize, &ls_status); CHKERRQ(ierr);
118a7e14dcfSSatish Balay       ierr = TaoAddLineSearchCounts(tao); CHKERRQ(ierr);
119a7e14dcfSSatish Balay 
120a7e14dcfSSatish Balay       if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) {
121a7e14dcfSSatish Balay 	tao->reason = TAO_DIVERGED_LS_FAILURE;
122a7e14dcfSSatish Balay 	break;
123a7e14dcfSSatish Balay       }
124a7e14dcfSSatish Balay     }
125a7e14dcfSSatish Balay 
126a7e14dcfSSatish Balay     /* Check for termination */
127a7e14dcfSSatish Balay     ierr = VecBoundGradientProjection(blmP->unprojected_gradient, tao->solution, tao->XL, tao->XU, tao->gradient); CHKERRQ(ierr);
128a7e14dcfSSatish Balay     ierr = VecNorm(tao->gradient, NORM_2, &gnorm); CHKERRQ(ierr);
129a7e14dcfSSatish Balay 
130a7e14dcfSSatish Balay 
131a7e14dcfSSatish Balay     if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) {
132a7e14dcfSSatish Balay       SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Not-a-Number");
133a7e14dcfSSatish Balay     }
134a7e14dcfSSatish Balay     iter++;
135a7e14dcfSSatish Balay     ierr = TaoMonitor(tao, iter, f, gnorm, 0.0, stepsize, &reason); CHKERRQ(ierr);
136a7e14dcfSSatish Balay   }
137a7e14dcfSSatish Balay 
138a7e14dcfSSatish Balay 
139a7e14dcfSSatish Balay 
140a7e14dcfSSatish Balay   PetscFunctionReturn(0);
141a7e14dcfSSatish Balay }
142a7e14dcfSSatish Balay 
143a7e14dcfSSatish Balay #undef __FUNCT__
144a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetup_BLMVM"
145a7e14dcfSSatish Balay static PetscErrorCode TaoSetup_BLMVM(TaoSolver tao)
146a7e14dcfSSatish Balay {
147a7e14dcfSSatish Balay   TAO_BLMVM *blmP = (TAO_BLMVM *)tao->data;
148a7e14dcfSSatish Balay   PetscInt n,N;
149a7e14dcfSSatish Balay   PetscErrorCode ierr;
150a7e14dcfSSatish Balay 
151a7e14dcfSSatish Balay   PetscFunctionBegin;
152a7e14dcfSSatish Balay   /* Existence of tao->solution checked in TaoSetup() */
153a7e14dcfSSatish Balay 
154a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution,&blmP->Xold); CHKERRQ(ierr);
155a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution,&blmP->Gold); CHKERRQ(ierr);
156a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution, &blmP->unprojected_gradient);
157a7e14dcfSSatish Balay 
158a7e14dcfSSatish Balay   if (!tao->stepdirection) {
159a7e14dcfSSatish Balay       ierr = VecDuplicate(tao->solution, &tao->stepdirection);
160a7e14dcfSSatish Balay       CHKERRQ(ierr);
161a7e14dcfSSatish Balay   }
162a7e14dcfSSatish Balay   if (!tao->gradient) {
163a7e14dcfSSatish Balay       ierr = VecDuplicate(tao->solution,&tao->gradient); CHKERRQ(ierr);
164a7e14dcfSSatish Balay   }
165a7e14dcfSSatish Balay   if (!tao->XL) {
166a7e14dcfSSatish Balay       ierr = VecDuplicate(tao->solution,&tao->XL); CHKERRQ(ierr);
167a7e14dcfSSatish Balay       ierr = VecSet(tao->XL,TAO_NINFINITY); CHKERRQ(ierr);
168a7e14dcfSSatish Balay   }
169a7e14dcfSSatish Balay   if (!tao->XU) {
170a7e14dcfSSatish Balay       ierr = VecDuplicate(tao->solution,&tao->XU); CHKERRQ(ierr);
171a7e14dcfSSatish Balay       ierr = VecSet(tao->XU,TAO_INFINITY); CHKERRQ(ierr);
172a7e14dcfSSatish Balay   }
173a7e14dcfSSatish Balay   /* Create matrix for the limited memory approximation */
174a7e14dcfSSatish Balay   ierr = VecGetLocalSize(tao->solution,&n); CHKERRQ(ierr);
175a7e14dcfSSatish Balay   ierr = VecGetSize(tao->solution,&N); CHKERRQ(ierr);
176a7e14dcfSSatish Balay   ierr = MatCreateLMVM(((PetscObject)tao)->comm,n,N,&blmP->M); CHKERRQ(ierr);
177a7e14dcfSSatish Balay   ierr = MatLMVMAllocateVectors(blmP->M,tao->solution); CHKERRQ(ierr);
178a7e14dcfSSatish Balay 
179a7e14dcfSSatish Balay    PetscFunctionReturn(0);
180a7e14dcfSSatish Balay }
181a7e14dcfSSatish Balay 
182a7e14dcfSSatish Balay /* ---------------------------------------------------------- */
183a7e14dcfSSatish Balay #undef __FUNCT__
184a7e14dcfSSatish Balay #define __FUNCT__ "TaoDestroy_BLMVM"
185a7e14dcfSSatish Balay static PetscErrorCode TaoDestroy_BLMVM(TaoSolver tao)
186a7e14dcfSSatish Balay {
187a7e14dcfSSatish Balay   TAO_BLMVM *blmP = (TAO_BLMVM *)tao->data;
188a7e14dcfSSatish Balay   PetscErrorCode ierr;
189a7e14dcfSSatish Balay 
190a7e14dcfSSatish Balay   PetscFunctionBegin;
191a7e14dcfSSatish Balay   if (tao->setupcalled) {
192a7e14dcfSSatish Balay     ierr = MatDestroy(&blmP->M); CHKERRQ(ierr);
193a7e14dcfSSatish Balay     ierr = VecDestroy(&blmP->unprojected_gradient); CHKERRQ(ierr);
194a7e14dcfSSatish Balay     ierr = VecDestroy(&blmP->Xold); CHKERRQ(ierr);
195a7e14dcfSSatish Balay     ierr = VecDestroy(&blmP->Gold); CHKERRQ(ierr);
196a7e14dcfSSatish Balay   }
197a7e14dcfSSatish Balay   ierr = PetscFree(tao->data); CHKERRQ(ierr);
198a7e14dcfSSatish Balay   tao->data = PETSC_NULL;
199a7e14dcfSSatish Balay 
200a7e14dcfSSatish Balay   PetscFunctionReturn(0);
201a7e14dcfSSatish Balay }
202a7e14dcfSSatish Balay 
203a7e14dcfSSatish Balay /*------------------------------------------------------------*/
204a7e14dcfSSatish Balay #undef __FUNCT__
205a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetFromOptions_BLMVM"
206a7e14dcfSSatish Balay static PetscErrorCode TaoSetFromOptions_BLMVM(TaoSolver tao)
207a7e14dcfSSatish Balay {
208a7e14dcfSSatish Balay 
209a7e14dcfSSatish Balay   PetscErrorCode ierr;
210a7e14dcfSSatish Balay 
211a7e14dcfSSatish Balay   PetscFunctionBegin;
212a7e14dcfSSatish Balay   ierr = PetscOptionsHead("Limited-memory variable-metric method for bound constrained optimization"); CHKERRQ(ierr);
213a7e14dcfSSatish Balay   ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr);
214a7e14dcfSSatish Balay   ierr = PetscOptionsTail();CHKERRQ(ierr);
215a7e14dcfSSatish Balay   PetscFunctionReturn(0);
216a7e14dcfSSatish Balay }
217a7e14dcfSSatish Balay 
218a7e14dcfSSatish Balay 
219a7e14dcfSSatish Balay /*------------------------------------------------------------*/
220a7e14dcfSSatish Balay #undef __FUNCT__
221a7e14dcfSSatish Balay #define __FUNCT__ "TaoView_BLMVM"
222a7e14dcfSSatish Balay static int TaoView_BLMVM(TaoSolver tao, PetscViewer viewer)
223a7e14dcfSSatish Balay {
224a7e14dcfSSatish Balay 
225a7e14dcfSSatish Balay 
226a7e14dcfSSatish Balay     TAO_BLMVM *lmP = (TAO_BLMVM *)tao->data;
227a7e14dcfSSatish Balay     PetscBool isascii;
228a7e14dcfSSatish Balay     PetscErrorCode ierr;
229a7e14dcfSSatish Balay 
230a7e14dcfSSatish Balay 
231a7e14dcfSSatish Balay     PetscFunctionBegin;
232a7e14dcfSSatish Balay     ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii); CHKERRQ(ierr);
233a7e14dcfSSatish Balay     if (isascii) {
234a7e14dcfSSatish Balay       ierr = PetscViewerASCIIPushTab(viewer); CHKERRQ(ierr);
235a7e14dcfSSatish Balay       ierr = PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", lmP->grad); CHKERRQ(ierr);
236a7e14dcfSSatish Balay       ierr = PetscViewerASCIIPopTab(viewer); CHKERRQ(ierr);
237a7e14dcfSSatish Balay     } else {
238a7e14dcfSSatish Balay       SETERRQ1(((PetscObject)tao)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for TAO BLMVM",((PetscObject)viewer)->type_name);
239a7e14dcfSSatish Balay     }
240a7e14dcfSSatish Balay     PetscFunctionReturn(0);
241a7e14dcfSSatish Balay }
242a7e14dcfSSatish Balay 
243a7e14dcfSSatish Balay #undef __FUNCT__
244a7e14dcfSSatish Balay #define __FUNCT__ "TaoComputeDual_BLMVM"
245a7e14dcfSSatish Balay static PetscErrorCode TaoComputeDual_BLMVM(TaoSolver tao, Vec DXL, Vec DXU)
246a7e14dcfSSatish Balay {
247a7e14dcfSSatish Balay   TAO_BLMVM *blm = (TAO_BLMVM *) tao->data;
248a7e14dcfSSatish Balay   PetscErrorCode ierr;
249a7e14dcfSSatish Balay 
250a7e14dcfSSatish Balay   PetscFunctionBegin;
251a7e14dcfSSatish Balay   PetscValidHeaderSpecific(tao,TAOSOLVER_CLASSID,1);
252a7e14dcfSSatish Balay   PetscValidHeaderSpecific(DXL,VEC_CLASSID,2);
253a7e14dcfSSatish Balay   PetscValidHeaderSpecific(DXU,VEC_CLASSID,3);
254a7e14dcfSSatish Balay 
255a7e14dcfSSatish Balay   if (!tao->gradient || !blm->unprojected_gradient) {
256a7e14dcfSSatish Balay       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Dual variables don't exist yet or no longer exist.\n");
257a7e14dcfSSatish Balay   }
258a7e14dcfSSatish Balay 
259a7e14dcfSSatish Balay   ierr = VecCopy(tao->gradient,DXL); CHKERRQ(ierr);
260a7e14dcfSSatish Balay   ierr = VecAXPY(DXL,-1.0,blm->unprojected_gradient); CHKERRQ(ierr);
261a7e14dcfSSatish Balay   ierr = VecSet(DXU,0.0); CHKERRQ(ierr);
262a7e14dcfSSatish Balay   ierr = VecPointwiseMax(DXL,DXL,DXU); CHKERRQ(ierr);
263a7e14dcfSSatish Balay 
264a7e14dcfSSatish Balay   ierr = VecCopy(blm->unprojected_gradient,DXU); CHKERRQ(ierr);
265a7e14dcfSSatish Balay   ierr = VecAXPY(DXU,-1.0,tao->gradient); CHKERRQ(ierr);
266a7e14dcfSSatish Balay   ierr = VecAXPY(DXU,1.0,DXL); CHKERRQ(ierr);
267a7e14dcfSSatish Balay 
268a7e14dcfSSatish Balay   PetscFunctionReturn(0);
269a7e14dcfSSatish Balay }
270a7e14dcfSSatish Balay 
271a7e14dcfSSatish Balay /* ---------------------------------------------------------- */
272a7e14dcfSSatish Balay EXTERN_C_BEGIN
273a7e14dcfSSatish Balay #undef __FUNCT__
274a7e14dcfSSatish Balay #define __FUNCT__ "TaoCreate_BLMVM"
275a7e14dcfSSatish Balay PetscErrorCode TaoCreate_BLMVM(TaoSolver tao)
276a7e14dcfSSatish Balay {
277a7e14dcfSSatish Balay   TAO_BLMVM *blmP;
278a7e14dcfSSatish Balay   const char *morethuente_type = TAOLINESEARCH_MT;
279a7e14dcfSSatish Balay   PetscErrorCode ierr;
280a7e14dcfSSatish Balay 
281a7e14dcfSSatish Balay   PetscFunctionBegin;
282a7e14dcfSSatish Balay   tao->ops->setup = TaoSetup_BLMVM;
283a7e14dcfSSatish Balay   tao->ops->solve = TaoSolve_BLMVM;
284a7e14dcfSSatish Balay   tao->ops->view = TaoView_BLMVM;
285a7e14dcfSSatish Balay   tao->ops->setfromoptions = TaoSetFromOptions_BLMVM;
286a7e14dcfSSatish Balay   tao->ops->destroy = TaoDestroy_BLMVM;
287a7e14dcfSSatish Balay   tao->ops->computedual = TaoComputeDual_BLMVM;
288a7e14dcfSSatish Balay 
289*3c9e27cfSGeoffrey Irving   ierr = PetscNewLog(tao,&blmP); CHKERRQ(ierr);
290a7e14dcfSSatish Balay   tao->data = (void*)blmP;
291a7e14dcfSSatish Balay   tao->max_it = 2000;
292a7e14dcfSSatish Balay   tao->max_funcs = 4000;
293a7e14dcfSSatish Balay   tao->fatol = 1e-4;
294a7e14dcfSSatish Balay   tao->frtol = 1e-4;
295a7e14dcfSSatish Balay 
296a7e14dcfSSatish Balay   ierr = TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch); CHKERRQ(ierr);
297a7e14dcfSSatish Balay   ierr = TaoLineSearchSetType(tao->linesearch, morethuente_type); CHKERRQ(ierr);
298a7e14dcfSSatish Balay   ierr = TaoLineSearchUseTaoSolverRoutines(tao->linesearch,tao); CHKERRQ(ierr);
299a7e14dcfSSatish Balay   PetscFunctionReturn(0);
300a7e14dcfSSatish Balay 
301a7e14dcfSSatish Balay }
302a7e14dcfSSatish Balay EXTERN_C_END
303a7e14dcfSSatish Balay 
304