xref: /petsc/src/tao/bound/impls/blmvm/blmvm.c (revision 63b1541575edd3f8a34e081896782a8d7a30de53)
1ba92ff59SBarry Smith #include <petsctaolinesearch.h>
2aaa7dc30SBarry Smith #include <../src/tao/matrix/lmvmmat.h>
3a9603a14SPatrick Farrell #include <../src/tao/unconstrained/impls/lmvm/lmvm.h>
4aaa7dc30SBarry Smith #include <../src/tao/bound/impls/blmvm/blmvm.h>
5a7e14dcfSSatish Balay 
6a7e14dcfSSatish Balay /*------------------------------------------------------------*/
7441846f8SBarry Smith static PetscErrorCode TaoSolve_BLMVM(Tao tao)
8a7e14dcfSSatish Balay {
9a7e14dcfSSatish Balay   PetscErrorCode               ierr;
10a7e14dcfSSatish Balay   TAO_BLMVM                    *blmP = (TAO_BLMVM *)tao->data;
11e4cb33bbSBarry Smith   TaoConvergedReason           reason = TAO_CONTINUE_ITERATING;
12e4cb33bbSBarry Smith   TaoLineSearchConvergedReason ls_status = TAOLINESEARCH_CONTINUE_ITERATING;
13a7e14dcfSSatish Balay   PetscReal                    f, fold, gdx, gnorm;
14a7e14dcfSSatish Balay   PetscReal                    stepsize = 1.0,delta;
15a7e14dcfSSatish Balay 
16a7e14dcfSSatish Balay   PetscFunctionBegin;
17a7e14dcfSSatish Balay   /*  Project initial point onto bounds */
18a7e14dcfSSatish Balay   ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr);
19a7e14dcfSSatish Balay   ierr = VecMedian(tao->XL,tao->solution,tao->XU,tao->solution);CHKERRQ(ierr);
20a7e14dcfSSatish Balay   ierr = TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU);CHKERRQ(ierr);
21a7e14dcfSSatish Balay 
22a9603a14SPatrick Farrell 
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 
27a9603a14SPatrick Farrell   ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr);
283fcc4fd3STodd Munson   if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or 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);
401feb5a2fSTodd Munson   ierr = MatLMVMReset(blmP->M);CHKERRQ(ierr);
41a7e14dcfSSatish Balay 
42a7e14dcfSSatish Balay   /* Set counter for gradient/reset steps */
43a7e14dcfSSatish Balay   blmP->grad = 0;
44a7e14dcfSSatish Balay   blmP->reset = 0;
45a7e14dcfSSatish Balay 
46a7e14dcfSSatish Balay   /* Have not converged; continue with Newton method */
47a7e14dcfSSatish Balay   while (reason == TAO_CONTINUE_ITERATING) {
48a7e14dcfSSatish Balay     /* Compute direction */
49a7e14dcfSSatish Balay     ierr = MatLMVMUpdate(blmP->M, tao->solution, tao->gradient);CHKERRQ(ierr);
50a7e14dcfSSatish Balay     ierr = MatLMVMSolve(blmP->M, blmP->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
51a7e14dcfSSatish Balay     ierr = VecBoundGradientProjection(tao->stepdirection,tao->solution,tao->XL,tao->XU,tao->gradient);CHKERRQ(ierr);
52a7e14dcfSSatish Balay 
53a7e14dcfSSatish Balay     /* Check for success (descent direction) */
54a7e14dcfSSatish Balay     ierr = VecDot(blmP->unprojected_gradient, tao->gradient, &gdx);CHKERRQ(ierr);
55a7e14dcfSSatish Balay     if (gdx <= 0) {
56a7e14dcfSSatish Balay       /* Step is not descent or solve was not successful
57a7e14dcfSSatish Balay          Use steepest descent direction (scaled) */
58a7e14dcfSSatish Balay       ++blmP->grad;
59a7e14dcfSSatish Balay 
60a7e14dcfSSatish Balay       if (f != 0.0) {
61a7e14dcfSSatish Balay         delta = 2.0*PetscAbsScalar(f) / (gnorm*gnorm);
6253506e15SBarry Smith       } else {
63a7e14dcfSSatish Balay         delta = 2.0 / (gnorm*gnorm);
64a7e14dcfSSatish Balay       }
65a7e14dcfSSatish Balay       ierr = MatLMVMSetDelta(blmP->M,delta);CHKERRQ(ierr);
66a7e14dcfSSatish Balay       ierr = MatLMVMReset(blmP->M);CHKERRQ(ierr);
67a7e14dcfSSatish Balay       ierr = MatLMVMUpdate(blmP->M, tao->solution, blmP->unprojected_gradient);CHKERRQ(ierr);
68a7e14dcfSSatish Balay       ierr = MatLMVMSolve(blmP->M,blmP->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
69a7e14dcfSSatish Balay     }
70a7e14dcfSSatish Balay     ierr = VecScale(tao->stepdirection,-1.0);CHKERRQ(ierr);
71a7e14dcfSSatish Balay 
72a7e14dcfSSatish Balay     /* Perform the linesearch */
73a7e14dcfSSatish Balay     fold = f;
74a7e14dcfSSatish Balay     ierr = VecCopy(tao->solution, blmP->Xold);CHKERRQ(ierr);
75a7e14dcfSSatish Balay     ierr = VecCopy(blmP->unprojected_gradient, blmP->Gold);CHKERRQ(ierr);
76a7e14dcfSSatish Balay     ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,1.0);CHKERRQ(ierr);
77a7e14dcfSSatish Balay     ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, blmP->unprojected_gradient, tao->stepdirection, &stepsize, &ls_status);CHKERRQ(ierr);
78a7e14dcfSSatish Balay     ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
79a7e14dcfSSatish Balay 
80a7e14dcfSSatish Balay     if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) {
81a7e14dcfSSatish Balay       /* Linesearch failed
82a7e14dcfSSatish Balay          Reset factors and use scaled (projected) gradient step */
83a7e14dcfSSatish Balay       ++blmP->reset;
84a7e14dcfSSatish Balay 
85a7e14dcfSSatish Balay       f = fold;
86a7e14dcfSSatish Balay       ierr = VecCopy(blmP->Xold, tao->solution);CHKERRQ(ierr);
87a7e14dcfSSatish Balay       ierr = VecCopy(blmP->Gold, blmP->unprojected_gradient);CHKERRQ(ierr);
88a7e14dcfSSatish Balay 
89a7e14dcfSSatish Balay       if (f != 0.0) {
90a7e14dcfSSatish Balay         delta = 2.0* PetscAbsScalar(f) / (gnorm*gnorm);
9153506e15SBarry Smith       } else {
92a7e14dcfSSatish Balay         delta = 2.0/ (gnorm*gnorm);
93a7e14dcfSSatish Balay       }
94a7e14dcfSSatish Balay       ierr = MatLMVMSetDelta(blmP->M,delta);CHKERRQ(ierr);
95a7e14dcfSSatish Balay       ierr = MatLMVMReset(blmP->M);CHKERRQ(ierr);
96a7e14dcfSSatish Balay       ierr = MatLMVMUpdate(blmP->M, tao->solution, blmP->unprojected_gradient);CHKERRQ(ierr);
97a7e14dcfSSatish Balay       ierr = MatLMVMSolve(blmP->M, blmP->unprojected_gradient, tao->stepdirection);CHKERRQ(ierr);
98a7e14dcfSSatish Balay       ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
99a7e14dcfSSatish Balay 
1003fcc4fd3STodd Munson       /* This may be incorrect; linesearch has values for stepmax and stepmin
101a7e14dcfSSatish Balay          that should be reset. */
102302440fdSBarry Smith       ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,1.0);CHKERRQ(ierr);
103a7e14dcfSSatish Balay       ierr = TaoLineSearchApply(tao->linesearch,tao->solution,&f, blmP->unprojected_gradient, tao->stepdirection,  &stepsize, &ls_status);CHKERRQ(ierr);
104a7e14dcfSSatish Balay       ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
105a7e14dcfSSatish Balay 
106a7e14dcfSSatish Balay       if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) {
107a7e14dcfSSatish Balay         tao->reason = TAO_DIVERGED_LS_FAILURE;
108a7e14dcfSSatish Balay         break;
109a7e14dcfSSatish Balay       }
110a7e14dcfSSatish Balay     }
111a7e14dcfSSatish Balay 
112e4cb33bbSBarry Smith     /* Check for converged */
113a7e14dcfSSatish Balay     ierr = VecBoundGradientProjection(blmP->unprojected_gradient, tao->solution, tao->XL, tao->XU, tao->gradient);CHKERRQ(ierr);
114a9603a14SPatrick Farrell     ierr = TaoGradientNorm(tao, tao->gradient, NORM_2, &gnorm);CHKERRQ(ierr);
115a7e14dcfSSatish Balay 
116a7e14dcfSSatish Balay 
11753506e15SBarry Smith     if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Not-a-Number");
1188931d482SJason Sarich     tao->niter++;
1198931d482SJason Sarich     ierr = TaoMonitor(tao, tao->niter, f, gnorm, 0.0, stepsize, &reason);CHKERRQ(ierr);
120a7e14dcfSSatish Balay   }
121a7e14dcfSSatish Balay   PetscFunctionReturn(0);
122a7e14dcfSSatish Balay }
123a7e14dcfSSatish Balay 
124441846f8SBarry Smith static PetscErrorCode TaoSetup_BLMVM(Tao tao)
125a7e14dcfSSatish Balay {
126a7e14dcfSSatish Balay   TAO_BLMVM      *blmP = (TAO_BLMVM *)tao->data;
127a7e14dcfSSatish Balay   PetscInt       n,N;
128a7e14dcfSSatish Balay   PetscErrorCode ierr;
129a9603a14SPatrick Farrell   KSP            H0ksp;
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);
156a9603a14SPatrick Farrell 
157a9603a14SPatrick Farrell   /* If the user has set a matrix to solve as the initial H0, set the options prefix here, and set up the KSP */
158a9603a14SPatrick Farrell   if (blmP->H0) {
159a9603a14SPatrick Farrell     const char *prefix;
160a9603a14SPatrick Farrell     ierr = MatLMVMSetH0(blmP->M, blmP->H0);CHKERRQ(ierr);
161a9603a14SPatrick Farrell     ierr = MatLMVMGetH0KSP(blmP->M, &H0ksp);CHKERRQ(ierr);
162a9603a14SPatrick Farrell 
163a9603a14SPatrick Farrell     ierr = TaoGetOptionsPrefix(tao, &prefix);CHKERRQ(ierr);
164a9603a14SPatrick Farrell     ierr = KSPSetOptionsPrefix(H0ksp, prefix);CHKERRQ(ierr);
16517f6384fSBarry Smith     ierr = KSPAppendOptionsPrefix(H0ksp, "tao_h0_");CHKERRQ(ierr);
166a9603a14SPatrick Farrell     ierr = KSPSetFromOptions(H0ksp);CHKERRQ(ierr);
167a9603a14SPatrick Farrell     ierr = KSPSetUp(H0ksp);CHKERRQ(ierr);
168a9603a14SPatrick Farrell   }
169a7e14dcfSSatish Balay   PetscFunctionReturn(0);
170a7e14dcfSSatish Balay }
171a7e14dcfSSatish Balay 
172a7e14dcfSSatish Balay /* ---------------------------------------------------------- */
173441846f8SBarry Smith static PetscErrorCode TaoDestroy_BLMVM(Tao tao)
174a7e14dcfSSatish Balay {
175a7e14dcfSSatish Balay   TAO_BLMVM      *blmP = (TAO_BLMVM *)tao->data;
176a7e14dcfSSatish Balay   PetscErrorCode ierr;
177a7e14dcfSSatish Balay 
178a7e14dcfSSatish Balay   PetscFunctionBegin;
179a7e14dcfSSatish Balay   if (tao->setupcalled) {
180a7e14dcfSSatish Balay     ierr = MatDestroy(&blmP->M);CHKERRQ(ierr);
181a7e14dcfSSatish Balay     ierr = VecDestroy(&blmP->unprojected_gradient);CHKERRQ(ierr);
182a7e14dcfSSatish Balay     ierr = VecDestroy(&blmP->Xold);CHKERRQ(ierr);
183a7e14dcfSSatish Balay     ierr = VecDestroy(&blmP->Gold);CHKERRQ(ierr);
184a7e14dcfSSatish Balay   }
185a9603a14SPatrick Farrell 
186a9603a14SPatrick Farrell   if (blmP->H0) {
187a9603a14SPatrick Farrell     PetscObjectDereference((PetscObject)blmP->H0);
188a9603a14SPatrick Farrell   }
189a9603a14SPatrick Farrell 
190a7e14dcfSSatish Balay   ierr = PetscFree(tao->data);CHKERRQ(ierr);
191a7e14dcfSSatish Balay   PetscFunctionReturn(0);
192a7e14dcfSSatish Balay }
193a7e14dcfSSatish Balay 
194a7e14dcfSSatish Balay /*------------------------------------------------------------*/
1954416b707SBarry Smith static PetscErrorCode TaoSetFromOptions_BLMVM(PetscOptionItems* PetscOptionsObject,Tao tao)
196a7e14dcfSSatish Balay {
197a7e14dcfSSatish Balay   PetscErrorCode ierr;
198a7e14dcfSSatish Balay 
199a7e14dcfSSatish Balay   PetscFunctionBegin;
2001a1499c8SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Limited-memory variable-metric method for bound constrained optimization");CHKERRQ(ierr);
201a7e14dcfSSatish Balay   ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr);
202a7e14dcfSSatish Balay   ierr = PetscOptionsTail();CHKERRQ(ierr);
203a7e14dcfSSatish Balay   PetscFunctionReturn(0);
204a7e14dcfSSatish Balay }
205a7e14dcfSSatish Balay 
206a7e14dcfSSatish Balay 
207a7e14dcfSSatish Balay /*------------------------------------------------------------*/
208441846f8SBarry Smith static int TaoView_BLMVM(Tao tao, PetscViewer viewer)
209a7e14dcfSSatish Balay {
210a7e14dcfSSatish Balay   TAO_BLMVM      *lmP = (TAO_BLMVM *)tao->data;
211a7e14dcfSSatish Balay   PetscBool      isascii;
212a7e14dcfSSatish Balay   PetscErrorCode ierr;
213a7e14dcfSSatish Balay 
214a7e14dcfSSatish Balay   PetscFunctionBegin;
215a7e14dcfSSatish Balay   ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii);CHKERRQ(ierr);
216a7e14dcfSSatish Balay   if (isascii) {
217a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", lmP->grad);CHKERRQ(ierr);
218a7e14dcfSSatish Balay   }
219a7e14dcfSSatish Balay   PetscFunctionReturn(0);
220a7e14dcfSSatish Balay }
221a7e14dcfSSatish Balay 
222441846f8SBarry Smith static PetscErrorCode TaoComputeDual_BLMVM(Tao tao, Vec DXL, Vec DXU)
223a7e14dcfSSatish Balay {
224a7e14dcfSSatish Balay   TAO_BLMVM      *blm = (TAO_BLMVM *) tao->data;
225a7e14dcfSSatish Balay   PetscErrorCode ierr;
226a7e14dcfSSatish Balay 
227a7e14dcfSSatish Balay   PetscFunctionBegin;
228441846f8SBarry Smith   PetscValidHeaderSpecific(tao,TAO_CLASSID,1);
229a7e14dcfSSatish Balay   PetscValidHeaderSpecific(DXL,VEC_CLASSID,2);
230a7e14dcfSSatish Balay   PetscValidHeaderSpecific(DXU,VEC_CLASSID,3);
23153506e15SBarry 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");
232a7e14dcfSSatish Balay 
233a7e14dcfSSatish Balay   ierr = VecCopy(tao->gradient,DXL);CHKERRQ(ierr);
234a7e14dcfSSatish Balay   ierr = VecAXPY(DXL,-1.0,blm->unprojected_gradient);CHKERRQ(ierr);
235a7e14dcfSSatish Balay   ierr = VecSet(DXU,0.0);CHKERRQ(ierr);
236a7e14dcfSSatish Balay   ierr = VecPointwiseMax(DXL,DXL,DXU);CHKERRQ(ierr);
237a7e14dcfSSatish Balay 
238a7e14dcfSSatish Balay   ierr = VecCopy(blm->unprojected_gradient,DXU);CHKERRQ(ierr);
239a7e14dcfSSatish Balay   ierr = VecAXPY(DXU,-1.0,tao->gradient);CHKERRQ(ierr);
240a7e14dcfSSatish Balay   ierr = VecAXPY(DXU,1.0,DXL);CHKERRQ(ierr);
241a7e14dcfSSatish Balay   PetscFunctionReturn(0);
242a7e14dcfSSatish Balay }
243a7e14dcfSSatish Balay 
244a7e14dcfSSatish Balay /* ---------------------------------------------------------- */
2451522df2eSJason Sarich /*MC
2461522df2eSJason Sarich   TAOBLMVM - Bounded limited memory variable metric is a quasi-Newton method
2471522df2eSJason Sarich          for nonlinear minimization with bound constraints. It is an extension
2481522df2eSJason Sarich          of TAOLMVM
2491522df2eSJason Sarich 
2501522df2eSJason Sarich   Options Database Keys:
2511522df2eSJason Sarich +     -tao_lmm_vectors - number of vectors to use for approximation
2521522df2eSJason Sarich .     -tao_lmm_scale_type - "none","scalar","broyden"
2531522df2eSJason Sarich .     -tao_lmm_limit_type - "none","average","relative","absolute"
2541522df2eSJason Sarich .     -tao_lmm_rescale_type - "none","scalar","gl"
2551522df2eSJason Sarich .     -tao_lmm_limit_mu - mu limiting factor
2561522df2eSJason Sarich .     -tao_lmm_limit_nu - nu limiting factor
2571522df2eSJason Sarich .     -tao_lmm_delta_min - minimum delta value
2581522df2eSJason Sarich .     -tao_lmm_delta_max - maximum delta value
2591522df2eSJason Sarich .     -tao_lmm_broyden_phi - phi factor for Broyden scaling
2601522df2eSJason Sarich .     -tao_lmm_scalar_alpha - alpha factor for scalar scaling
2611522df2eSJason Sarich .     -tao_lmm_rescale_alpha - alpha factor for rescaling diagonal
2621522df2eSJason Sarich .     -tao_lmm_rescale_beta - beta factor for rescaling diagonal
2631522df2eSJason Sarich .     -tao_lmm_scalar_history - amount of history for scalar scaling
2641522df2eSJason Sarich .     -tao_lmm_rescale_history - amount of history for rescaling diagonal
2651522df2eSJason Sarich -     -tao_lmm_eps - rejection tolerance
2661522df2eSJason Sarich 
2671eb8069cSJason Sarich   Level: beginner
2681522df2eSJason Sarich M*/
269728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate_BLMVM(Tao tao)
270a7e14dcfSSatish Balay {
271a7e14dcfSSatish Balay   TAO_BLMVM      *blmP;
2728caf6e8cSBarry Smith   const char     *morethuente_type = TAOLINESEARCHMT;
273a7e14dcfSSatish Balay   PetscErrorCode ierr;
274a7e14dcfSSatish Balay 
275a7e14dcfSSatish Balay   PetscFunctionBegin;
276a7e14dcfSSatish Balay   tao->ops->setup = TaoSetup_BLMVM;
277a7e14dcfSSatish Balay   tao->ops->solve = TaoSolve_BLMVM;
278a7e14dcfSSatish Balay   tao->ops->view = TaoView_BLMVM;
279a7e14dcfSSatish Balay   tao->ops->setfromoptions = TaoSetFromOptions_BLMVM;
280a7e14dcfSSatish Balay   tao->ops->destroy = TaoDestroy_BLMVM;
281a7e14dcfSSatish Balay   tao->ops->computedual = TaoComputeDual_BLMVM;
282a7e14dcfSSatish Balay 
2833c9e27cfSGeoffrey Irving   ierr = PetscNewLog(tao,&blmP);CHKERRQ(ierr);
284a9603a14SPatrick Farrell   blmP->H0 = NULL;
285a7e14dcfSSatish Balay   tao->data = (void*)blmP;
2866552cf8aSJason Sarich 
2876552cf8aSJason Sarich   /* Override default settings (unless already changed) */
2886552cf8aSJason Sarich   if (!tao->max_it_changed) tao->max_it = 2000;
2896552cf8aSJason Sarich   if (!tao->max_funcs_changed) tao->max_funcs = 4000;
290a7e14dcfSSatish Balay 
291a7e14dcfSSatish Balay   ierr = TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch);CHKERRQ(ierr);
292*63b15415SAlp Dener   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr);
293a7e14dcfSSatish Balay   ierr = TaoLineSearchSetType(tao->linesearch, morethuente_type);CHKERRQ(ierr);
294441846f8SBarry Smith   ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr);
2955d527766SPatrick Farrell   ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr);
296a7e14dcfSSatish Balay   PetscFunctionReturn(0);
297a7e14dcfSSatish Balay }
298a7e14dcfSSatish Balay 
299a9603a14SPatrick Farrell PETSC_EXTERN PetscErrorCode TaoLMVMSetH0(Tao tao, Mat H0)
300a9603a14SPatrick Farrell {
301a9603a14SPatrick Farrell   TAO_LMVM       *lmP;
302a9603a14SPatrick Farrell   TAO_BLMVM      *blmP;
303a9603a14SPatrick Farrell   const TaoType  type;
304a9603a14SPatrick Farrell   PetscBool      is_lmvm, is_blmvm;
305a9603a14SPatrick Farrell   PetscErrorCode ierr;
306a9603a14SPatrick Farrell 
307a9603a14SPatrick Farrell   ierr = TaoGetType(tao, &type);CHKERRQ(ierr);
308a9603a14SPatrick Farrell   ierr = PetscStrcmp(type, TAOLMVM,  &is_lmvm);CHKERRQ(ierr);
309a9603a14SPatrick Farrell   ierr = PetscStrcmp(type, TAOBLMVM, &is_blmvm);CHKERRQ(ierr);
310a9603a14SPatrick Farrell 
311a9603a14SPatrick Farrell   if (is_lmvm) {
312a9603a14SPatrick Farrell     lmP = (TAO_LMVM *)tao->data;
313a9603a14SPatrick Farrell     ierr = PetscObjectReference((PetscObject)H0);CHKERRQ(ierr);
314a9603a14SPatrick Farrell     lmP->H0 = H0;
315a9603a14SPatrick Farrell   } else if (is_blmvm) {
316a9603a14SPatrick Farrell     blmP = (TAO_BLMVM *)tao->data;
317a9603a14SPatrick Farrell     ierr = PetscObjectReference((PetscObject)H0);CHKERRQ(ierr);
318a9603a14SPatrick Farrell     blmP->H0 = H0;
31974c66251SPatrick Farrell   } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "This routine applies to TAO_LMVM and TAO_BLMVM.");
320a9603a14SPatrick Farrell   PetscFunctionReturn(0);
321a9603a14SPatrick Farrell }
322a9603a14SPatrick Farrell 
323a9603a14SPatrick Farrell PETSC_EXTERN PetscErrorCode TaoLMVMGetH0(Tao tao, Mat *H0)
324a9603a14SPatrick Farrell {
325a9603a14SPatrick Farrell   TAO_LMVM       *lmP;
326a9603a14SPatrick Farrell   TAO_BLMVM      *blmP;
327a9603a14SPatrick Farrell   const TaoType  type;
328a9603a14SPatrick Farrell   PetscBool      is_lmvm, is_blmvm;
329a9603a14SPatrick Farrell   Mat            M;
330a9603a14SPatrick Farrell 
331a9603a14SPatrick Farrell   PetscErrorCode ierr;
332a9603a14SPatrick Farrell 
333a9603a14SPatrick Farrell   ierr = TaoGetType(tao, &type);CHKERRQ(ierr);
334a9603a14SPatrick Farrell   ierr = PetscStrcmp(type, TAOLMVM,  &is_lmvm);CHKERRQ(ierr);
335a9603a14SPatrick Farrell   ierr = PetscStrcmp(type, TAOBLMVM, &is_blmvm);CHKERRQ(ierr);
336a9603a14SPatrick Farrell 
337a9603a14SPatrick Farrell   if (is_lmvm) {
338a9603a14SPatrick Farrell     lmP = (TAO_LMVM *)tao->data;
339a9603a14SPatrick Farrell     M = lmP->M;
340a9603a14SPatrick Farrell   } else if (is_blmvm) {
341a9603a14SPatrick Farrell     blmP = (TAO_BLMVM *)tao->data;
342a9603a14SPatrick Farrell     M = blmP->M;
34374c66251SPatrick Farrell   } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "This routine applies to TAO_LMVM and TAO_BLMVM.");
344a9603a14SPatrick Farrell   ierr = MatLMVMGetH0(M, H0);CHKERRQ(ierr);
345a9603a14SPatrick Farrell   PetscFunctionReturn(0);
346a9603a14SPatrick Farrell }
347a9603a14SPatrick Farrell 
348a9603a14SPatrick Farrell PETSC_EXTERN PetscErrorCode TaoLMVMGetH0KSP(Tao tao, KSP *ksp)
349a9603a14SPatrick Farrell {
350a9603a14SPatrick Farrell   TAO_LMVM       *lmP;
351a9603a14SPatrick Farrell   TAO_BLMVM      *blmP;
352a9603a14SPatrick Farrell   const TaoType  type;
353a9603a14SPatrick Farrell   PetscBool      is_lmvm, is_blmvm;
354a9603a14SPatrick Farrell   Mat            M;
355a9603a14SPatrick Farrell   PetscErrorCode ierr;
356a9603a14SPatrick Farrell 
357a9603a14SPatrick Farrell   ierr = TaoGetType(tao, &type);CHKERRQ(ierr);
358a9603a14SPatrick Farrell   ierr = PetscStrcmp(type, TAOLMVM,  &is_lmvm);CHKERRQ(ierr);
359a9603a14SPatrick Farrell   ierr = PetscStrcmp(type, TAOBLMVM, &is_blmvm);CHKERRQ(ierr);
360a9603a14SPatrick Farrell 
361a9603a14SPatrick Farrell   if (is_lmvm) {
362a9603a14SPatrick Farrell     lmP = (TAO_LMVM *)tao->data;
363a9603a14SPatrick Farrell     M = lmP->M;
364a9603a14SPatrick Farrell   } else if (is_blmvm) {
365a9603a14SPatrick Farrell     blmP = (TAO_BLMVM *)tao->data;
366a9603a14SPatrick Farrell     M = blmP->M;
36774c66251SPatrick Farrell   } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "This routine applies to TAO_LMVM and TAO_BLMVM.");
368a9603a14SPatrick Farrell   ierr = MatLMVMGetH0KSP(M, ksp);CHKERRQ(ierr);
369a9603a14SPatrick Farrell   PetscFunctionReturn(0);
370a9603a14SPatrick Farrell }
371