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