xref: /petsc/src/tao/unconstrained/impls/lmvm/lmvm.c (revision 5d5277661bc9d5de49c003e79a2b7124bf8a2eb4)
1ba92ff59SBarry Smith #include <petsctaolinesearch.h>
2aaa7dc30SBarry Smith #include <../src/tao/matrix/lmvmmat.h>
3aaa7dc30SBarry Smith #include <../src/tao/unconstrained/impls/lmvm/lmvm.h>
4a7e14dcfSSatish Balay 
5a7e14dcfSSatish Balay #define LMVM_BFGS                0
6a7e14dcfSSatish Balay #define LMVM_SCALED_GRADIENT     1
7a7e14dcfSSatish Balay #define LMVM_GRADIENT            2
8a7e14dcfSSatish Balay 
9a7e14dcfSSatish Balay #undef __FUNCT__
10a7e14dcfSSatish Balay #define __FUNCT__ "TaoSolve_LMVM"
11441846f8SBarry Smith static PetscErrorCode TaoSolve_LMVM(Tao tao)
12a7e14dcfSSatish Balay {
13a7e14dcfSSatish Balay   TAO_LMVM                     *lmP = (TAO_LMVM *)tao->data;
14a7e14dcfSSatish Balay   PetscReal                    f, fold, gdx, gnorm;
15a7e14dcfSSatish Balay   PetscReal                    step = 1.0;
16a7e14dcfSSatish Balay   PetscReal                    delta;
17a7e14dcfSSatish Balay   PetscErrorCode               ierr;
18a7e14dcfSSatish Balay   PetscInt                     stepType;
19e4cb33bbSBarry Smith   TaoConvergedReason           reason = TAO_CONTINUE_ITERATING;
20e4cb33bbSBarry Smith   TaoLineSearchConvergedReason ls_status = TAOLINESEARCH_CONTINUE_ITERATING;
21a7e14dcfSSatish Balay 
22a7e14dcfSSatish Balay   PetscFunctionBegin;
23a7e14dcfSSatish Balay 
24a7e14dcfSSatish Balay   if (tao->XL || tao->XU || tao->ops->computebounds) {
25a7e14dcfSSatish Balay     ierr = PetscPrintf(((PetscObject)tao)->comm,"WARNING: Variable bounds have been set but will be ignored by lmvm algorithm\n");CHKERRQ(ierr);
26a7e14dcfSSatish Balay   }
27a7e14dcfSSatish Balay 
28a7e14dcfSSatish Balay   /*  Check convergence criteria */
29a7e14dcfSSatish Balay   ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient);CHKERRQ(ierr);
30a7e14dcfSSatish Balay   ierr = VecNorm(tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr);
3187f595a5SBarry Smith   if (PetscIsInfOrNanReal(f) || PetscIsInfOrNanReal(gnorm)) SETERRQ(PETSC_COMM_SELF,1, "User provided compute function generated Inf or NaN");
32a7e14dcfSSatish Balay 
338931d482SJason Sarich   ierr = TaoMonitor(tao, tao->niter, f, gnorm, 0.0, step, &reason);CHKERRQ(ierr);
3487f595a5SBarry Smith   if (reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(0);
35a7e14dcfSSatish Balay 
36a7e14dcfSSatish Balay   /*  Set initial scaling for the function */
37a7e14dcfSSatish Balay   if (f != 0.0) {
38a7e14dcfSSatish Balay     delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm);
3987f595a5SBarry Smith   } else {
40a7e14dcfSSatish Balay     delta = 2.0 / (gnorm*gnorm);
41a7e14dcfSSatish Balay   }
42a7e14dcfSSatish Balay   ierr = MatLMVMSetDelta(lmP->M,delta);CHKERRQ(ierr);
43a7e14dcfSSatish Balay 
44a7e14dcfSSatish Balay   /*  Set counter for gradient/reset steps */
45a7e14dcfSSatish Balay   lmP->bfgs = 0;
46a7e14dcfSSatish Balay   lmP->sgrad = 0;
47a7e14dcfSSatish Balay   lmP->grad = 0;
48a7e14dcfSSatish Balay 
49a7e14dcfSSatish Balay   /*  Have not converged; continue with Newton method */
50a7e14dcfSSatish Balay   while (reason == TAO_CONTINUE_ITERATING) {
51a7e14dcfSSatish Balay     /*  Compute direction */
52a7e14dcfSSatish Balay     ierr = MatLMVMUpdate(lmP->M,tao->solution,tao->gradient);CHKERRQ(ierr);
53a7e14dcfSSatish Balay     ierr = MatLMVMSolve(lmP->M, tao->gradient, lmP->D);CHKERRQ(ierr);
54a7e14dcfSSatish Balay     ++lmP->bfgs;
55a7e14dcfSSatish Balay 
56a7e14dcfSSatish Balay     /*  Check for success (descent direction) */
57a7e14dcfSSatish Balay     ierr = VecDot(lmP->D, tao->gradient, &gdx);CHKERRQ(ierr);
58a7e14dcfSSatish Balay     if ((gdx <= 0.0) || PetscIsInfOrNanReal(gdx)) {
59a7e14dcfSSatish Balay       /* Step is not descent or direction produced not a number
60a7e14dcfSSatish Balay          We can assert bfgsUpdates > 1 in this case because
61a7e14dcfSSatish Balay          the first solve produces the scaled gradient direction,
62a7e14dcfSSatish Balay          which is guaranteed to be descent
63a7e14dcfSSatish Balay 
64a7e14dcfSSatish Balay          Use steepest descent direction (scaled)
65a7e14dcfSSatish Balay       */
66a7e14dcfSSatish Balay 
67a7e14dcfSSatish Balay       ++lmP->grad;
68a7e14dcfSSatish Balay 
69a7e14dcfSSatish Balay       if (f != 0.0) {
70a7e14dcfSSatish Balay         delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm);
7187f595a5SBarry Smith       } else {
72a7e14dcfSSatish Balay         delta = 2.0 / (gnorm*gnorm);
73a7e14dcfSSatish Balay       }
74a7e14dcfSSatish Balay       ierr = MatLMVMSetDelta(lmP->M, delta);CHKERRQ(ierr);
75a7e14dcfSSatish Balay       ierr = MatLMVMReset(lmP->M);CHKERRQ(ierr);
76a7e14dcfSSatish Balay       ierr = MatLMVMUpdate(lmP->M, tao->solution, tao->gradient);CHKERRQ(ierr);
77a7e14dcfSSatish Balay       ierr = MatLMVMSolve(lmP->M,tao->gradient, lmP->D);CHKERRQ(ierr);
78a7e14dcfSSatish Balay 
79a7e14dcfSSatish Balay       /* On a reset, the direction cannot be not a number; it is a
80a7e14dcfSSatish Balay          scaled gradient step.  No need to check for this condition. */
81a7e14dcfSSatish Balay 
82a7e14dcfSSatish Balay       lmP->bfgs = 1;
83a7e14dcfSSatish Balay       ++lmP->sgrad;
84a7e14dcfSSatish Balay       stepType = LMVM_SCALED_GRADIENT;
8587f595a5SBarry Smith     } else {
86a7e14dcfSSatish Balay       if (1 == lmP->bfgs) {
87a7e14dcfSSatish Balay         /*  The first BFGS direction is always the scaled gradient */
88a7e14dcfSSatish Balay         ++lmP->sgrad;
89a7e14dcfSSatish Balay         stepType = LMVM_SCALED_GRADIENT;
9087f595a5SBarry Smith       } else {
91a7e14dcfSSatish Balay         ++lmP->bfgs;
92a7e14dcfSSatish Balay         stepType = LMVM_BFGS;
93a7e14dcfSSatish Balay       }
94a7e14dcfSSatish Balay     }
95a7e14dcfSSatish Balay     ierr = VecScale(lmP->D, -1.0);CHKERRQ(ierr);
96a7e14dcfSSatish Balay 
97a7e14dcfSSatish Balay     /*  Perform the linesearch */
98a7e14dcfSSatish Balay     fold = f;
99a7e14dcfSSatish Balay     ierr = VecCopy(tao->solution, lmP->Xold);CHKERRQ(ierr);
100a7e14dcfSSatish Balay     ierr = VecCopy(tao->gradient, lmP->Gold);CHKERRQ(ierr);
101a7e14dcfSSatish Balay 
102a7e14dcfSSatish Balay     ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, lmP->D, &step,&ls_status);CHKERRQ(ierr);
103a7e14dcfSSatish Balay     ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
104a7e14dcfSSatish Balay 
10587f595a5SBarry Smith     while (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER && (stepType != LMVM_GRADIENT)) {
106a7e14dcfSSatish Balay       /*  Linesearch failed */
107a7e14dcfSSatish Balay       /*  Reset factors and use scaled gradient step */
108a7e14dcfSSatish Balay       f = fold;
109a7e14dcfSSatish Balay       ierr = VecCopy(lmP->Xold, tao->solution);CHKERRQ(ierr);
110a7e14dcfSSatish Balay       ierr = VecCopy(lmP->Gold, tao->gradient);CHKERRQ(ierr);
111a7e14dcfSSatish Balay 
112a7e14dcfSSatish Balay       switch(stepType) {
113a7e14dcfSSatish Balay       case LMVM_BFGS:
114a7e14dcfSSatish Balay         /*  Failed to obtain acceptable iterate with BFGS step */
115a7e14dcfSSatish Balay         /*  Attempt to use the scaled gradient direction */
116a7e14dcfSSatish Balay 
117a7e14dcfSSatish Balay         if (f != 0.0) {
118a7e14dcfSSatish Balay           delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm);
11987f595a5SBarry Smith         } else {
120a7e14dcfSSatish Balay           delta = 2.0 / (gnorm*gnorm);
121a7e14dcfSSatish Balay         }
122a7e14dcfSSatish Balay         ierr = MatLMVMSetDelta(lmP->M, delta);CHKERRQ(ierr);
123a7e14dcfSSatish Balay         ierr = MatLMVMReset(lmP->M);CHKERRQ(ierr);
124a7e14dcfSSatish Balay         ierr = MatLMVMUpdate(lmP->M, tao->solution, tao->gradient);CHKERRQ(ierr);
125a7e14dcfSSatish Balay         ierr = MatLMVMSolve(lmP->M, tao->gradient, lmP->D);CHKERRQ(ierr);
126a7e14dcfSSatish Balay 
127a7e14dcfSSatish Balay         /* On a reset, the direction cannot be not a number; it is a
128a7e14dcfSSatish Balay            scaled gradient step.  No need to check for this condition. */
129a7e14dcfSSatish Balay 
130a7e14dcfSSatish Balay         lmP->bfgs = 1;
131a7e14dcfSSatish Balay         ++lmP->sgrad;
132a7e14dcfSSatish Balay         stepType = LMVM_SCALED_GRADIENT;
133a7e14dcfSSatish Balay         break;
134a7e14dcfSSatish Balay 
135a7e14dcfSSatish Balay       case LMVM_SCALED_GRADIENT:
136a7e14dcfSSatish Balay         /* The scaled gradient step did not produce a new iterate;
137a7e14dcfSSatish Balay            attempt to use the gradient direction.
138a7e14dcfSSatish Balay            Need to make sure we are not using a different diagonal scaling */
139a7e14dcfSSatish Balay         ierr = MatLMVMSetDelta(lmP->M, 1.0);CHKERRQ(ierr);
140a7e14dcfSSatish Balay         ierr = MatLMVMReset(lmP->M);CHKERRQ(ierr);
141a7e14dcfSSatish Balay         ierr = MatLMVMUpdate(lmP->M, tao->solution, tao->gradient);CHKERRQ(ierr);
142a7e14dcfSSatish Balay         ierr = MatLMVMSolve(lmP->M, tao->gradient, lmP->D);CHKERRQ(ierr);
143a7e14dcfSSatish Balay 
144a7e14dcfSSatish Balay         lmP->bfgs = 1;
145a7e14dcfSSatish Balay         ++lmP->grad;
146a7e14dcfSSatish Balay         stepType = LMVM_GRADIENT;
147a7e14dcfSSatish Balay         break;
148a7e14dcfSSatish Balay       }
149a7e14dcfSSatish Balay       ierr = VecScale(lmP->D, -1.0);CHKERRQ(ierr);
150a7e14dcfSSatish Balay 
151a7e14dcfSSatish Balay       /*  Perform the linesearch */
152a7e14dcfSSatish Balay       ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, lmP->D, &step, &ls_status);CHKERRQ(ierr);
153a7e14dcfSSatish Balay       ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
154a7e14dcfSSatish Balay     }
155a7e14dcfSSatish Balay 
15687f595a5SBarry Smith     if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) {
157a7e14dcfSSatish Balay       /*  Failed to find an improving point */
158a7e14dcfSSatish Balay       f = fold;
159a7e14dcfSSatish Balay       ierr = VecCopy(lmP->Xold, tao->solution);CHKERRQ(ierr);
160a7e14dcfSSatish Balay       ierr = VecCopy(lmP->Gold, tao->gradient);CHKERRQ(ierr);
161a7e14dcfSSatish Balay       step = 0.0;
162a7e14dcfSSatish Balay       reason = TAO_DIVERGED_LS_FAILURE;
163a7e14dcfSSatish Balay       tao->reason = TAO_DIVERGED_LS_FAILURE;
164a7e14dcfSSatish Balay     }
165a7e14dcfSSatish Balay     /*  Check for termination */
166a7e14dcfSSatish Balay     ierr = VecNorm(tao->gradient, NORM_2, &gnorm);CHKERRQ(ierr);
1678931d482SJason Sarich     tao->niter++;
1688931d482SJason Sarich     ierr = TaoMonitor(tao,tao->niter,f,gnorm,0.0,step,&reason);CHKERRQ(ierr);
169a7e14dcfSSatish Balay   }
170a7e14dcfSSatish Balay   PetscFunctionReturn(0);
171a7e14dcfSSatish Balay }
17287f595a5SBarry Smith 
173a7e14dcfSSatish Balay #undef __FUNCT__
174a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetUp_LMVM"
175441846f8SBarry Smith static PetscErrorCode TaoSetUp_LMVM(Tao tao)
176a7e14dcfSSatish Balay {
177a7e14dcfSSatish Balay   TAO_LMVM       *lmP = (TAO_LMVM *)tao->data;
178a7e14dcfSSatish Balay   PetscInt       n,N;
179a7e14dcfSSatish Balay   PetscErrorCode ierr;
180a7e14dcfSSatish Balay 
181a7e14dcfSSatish Balay   PetscFunctionBegin;
182a7e14dcfSSatish Balay   /* Existence of tao->solution checked in TaoSetUp() */
183a7e14dcfSSatish Balay   if (!tao->gradient) {ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);  }
184a7e14dcfSSatish Balay   if (!tao->stepdirection) {ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr);  }
185a7e14dcfSSatish Balay   if (!lmP->D) {ierr = VecDuplicate(tao->solution,&lmP->D);CHKERRQ(ierr);  }
186a7e14dcfSSatish Balay   if (!lmP->Xold) {ierr = VecDuplicate(tao->solution,&lmP->Xold);CHKERRQ(ierr);  }
187a7e14dcfSSatish Balay   if (!lmP->Gold) {ierr = VecDuplicate(tao->solution,&lmP->Gold);CHKERRQ(ierr);  }
188a7e14dcfSSatish Balay 
189a7e14dcfSSatish Balay   /*  Create matrix for the limited memory approximation */
190a7e14dcfSSatish Balay   ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr);
191a7e14dcfSSatish Balay   ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr);
192a7e14dcfSSatish Balay   ierr = MatCreateLMVM(((PetscObject)tao)->comm,n,N,&lmP->M);CHKERRQ(ierr);
193a7e14dcfSSatish Balay   ierr = MatLMVMAllocateVectors(lmP->M,tao->solution);CHKERRQ(ierr);
194a7e14dcfSSatish Balay   PetscFunctionReturn(0);
195a7e14dcfSSatish Balay }
196a7e14dcfSSatish Balay 
197a7e14dcfSSatish Balay /* ---------------------------------------------------------- */
198a7e14dcfSSatish Balay #undef __FUNCT__
199a7e14dcfSSatish Balay #define __FUNCT__ "TaoDestroy_LMVM"
200441846f8SBarry Smith static PetscErrorCode TaoDestroy_LMVM(Tao tao)
201a7e14dcfSSatish Balay {
202a7e14dcfSSatish Balay   TAO_LMVM       *lmP = (TAO_LMVM *)tao->data;
203a7e14dcfSSatish Balay   PetscErrorCode ierr;
204a7e14dcfSSatish Balay 
205a7e14dcfSSatish Balay   PetscFunctionBegin;
206a7e14dcfSSatish Balay   if (tao->setupcalled) {
207a7e14dcfSSatish Balay     ierr = VecDestroy(&lmP->Xold);CHKERRQ(ierr);
208a7e14dcfSSatish Balay     ierr = VecDestroy(&lmP->Gold);CHKERRQ(ierr);
209a7e14dcfSSatish Balay     ierr = VecDestroy(&lmP->D);CHKERRQ(ierr);
210a7e14dcfSSatish Balay     ierr = MatDestroy(&lmP->M);CHKERRQ(ierr);
211a7e14dcfSSatish Balay   }
212a7e14dcfSSatish Balay   ierr = PetscFree(tao->data);CHKERRQ(ierr);
213a7e14dcfSSatish Balay   PetscFunctionReturn(0);
214a7e14dcfSSatish Balay }
215a7e14dcfSSatish Balay 
216a7e14dcfSSatish Balay /*------------------------------------------------------------*/
217a7e14dcfSSatish Balay #undef __FUNCT__
218a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetFromOptions_LMVM"
2198c34d3f5SBarry Smith static PetscErrorCode TaoSetFromOptions_LMVM(PetscOptions *PetscOptionsObject,Tao tao)
220a7e14dcfSSatish Balay {
221a7e14dcfSSatish Balay   PetscErrorCode ierr;
222a7e14dcfSSatish Balay 
223a7e14dcfSSatish Balay   PetscFunctionBegin;
2241a1499c8SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Limited-memory variable-metric method for unconstrained optimization");CHKERRQ(ierr);
225a7e14dcfSSatish Balay   ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr);
226a7e14dcfSSatish Balay   ierr = PetscOptionsTail();CHKERRQ(ierr);
227a7e14dcfSSatish Balay   PetscFunctionReturn(0);
228a7e14dcfSSatish Balay }
229a7e14dcfSSatish Balay 
230a7e14dcfSSatish Balay /*------------------------------------------------------------*/
231a7e14dcfSSatish Balay #undef __FUNCT__
232a7e14dcfSSatish Balay #define __FUNCT__ "TaoView_LMVM"
233441846f8SBarry Smith static PetscErrorCode TaoView_LMVM(Tao tao, PetscViewer viewer)
234a7e14dcfSSatish Balay {
235a7e14dcfSSatish Balay   TAO_LMVM       *lm = (TAO_LMVM *)tao->data;
236a7e14dcfSSatish Balay   PetscBool      isascii;
237a7e14dcfSSatish Balay   PetscErrorCode ierr;
238a7e14dcfSSatish Balay 
239a7e14dcfSSatish Balay   PetscFunctionBegin;
240a7e14dcfSSatish Balay   ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii);CHKERRQ(ierr);
241a7e14dcfSSatish Balay   if (isascii) {
242a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
243a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "BFGS steps: %D\n", lm->bfgs);CHKERRQ(ierr);
244a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "Scaled gradient steps: %D\n", lm->sgrad);CHKERRQ(ierr);
245a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", lm->grad);CHKERRQ(ierr);
246a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
247a7e14dcfSSatish Balay   }
248a7e14dcfSSatish Balay   PetscFunctionReturn(0);
249a7e14dcfSSatish Balay }
250a7e14dcfSSatish Balay 
251a7e14dcfSSatish Balay /* ---------------------------------------------------------- */
252a7e14dcfSSatish Balay 
2534aa34175SJason Sarich /*MC
2544aa34175SJason Sarich      TAOLMVM - Limited Memory Variable Metric method is a quasi-Newton
2554aa34175SJason Sarich      optimization solver for unconstrained minimization. It solves
2564aa34175SJason Sarich      the Newton step
2574aa34175SJason Sarich               Hkdk = - gk
2584aa34175SJason Sarich 
2594aa34175SJason Sarich      using an approximation Bk in place of Hk, where Bk is composed using
2604aa34175SJason Sarich      the BFGS update formula. A More-Thuente line search is then used
2614aa34175SJason Sarich      to computed the steplength in the dk direction
2624aa34175SJason Sarich   Options Database Keys:
2634aa34175SJason Sarich +     -tao_lmm_vectors - number of vectors to use for approximation
2644aa34175SJason Sarich .     -tao_lmm_scale_type - "none","scalar","broyden"
2654aa34175SJason Sarich .     -tao_lmm_limit_type - "none","average","relative","absolute"
2664aa34175SJason Sarich .     -tao_lmm_rescale_type - "none","scalar","gl"
2674aa34175SJason Sarich .     -tao_lmm_limit_mu - mu limiting factor
2684aa34175SJason Sarich .     -tao_lmm_limit_nu - nu limiting factor
2694aa34175SJason Sarich .     -tao_lmm_delta_min - minimum delta value
2704aa34175SJason Sarich .     -tao_lmm_delta_max - maximum delta value
2714aa34175SJason Sarich .     -tao_lmm_broyden_phi - phi factor for Broyden scaling
2724aa34175SJason Sarich .     -tao_lmm_scalar_alpha - alpha factor for scalar scaling
2734aa34175SJason Sarich .     -tao_lmm_rescale_alpha - alpha factor for rescaling diagonal
2744aa34175SJason Sarich .     -tao_lmm_rescale_beta - beta factor for rescaling diagonal
2754aa34175SJason Sarich .     -tao_lmm_scalar_history - amount of history for scalar scaling
2764aa34175SJason Sarich .     -tao_lmm_rescale_history - amount of history for rescaling diagonal
2774aa34175SJason Sarich -     -tao_lmm_eps - rejection tolerance
2784aa34175SJason Sarich 
2791eb8069cSJason Sarich   Level: beginner
2804aa34175SJason Sarich M*/
2814aa34175SJason Sarich 
282a7e14dcfSSatish Balay #undef __FUNCT__
283a7e14dcfSSatish Balay #define __FUNCT__ "TaoCreate_LMVM"
284728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate_LMVM(Tao tao)
285a7e14dcfSSatish Balay {
286a7e14dcfSSatish Balay   TAO_LMVM       *lmP;
2878caf6e8cSBarry Smith   const char     *morethuente_type = TAOLINESEARCHMT;
288a7e14dcfSSatish Balay   PetscErrorCode ierr;
289a7e14dcfSSatish Balay 
290a7e14dcfSSatish Balay   PetscFunctionBegin;
291a7e14dcfSSatish Balay   tao->ops->setup = TaoSetUp_LMVM;
292a7e14dcfSSatish Balay   tao->ops->solve = TaoSolve_LMVM;
293a7e14dcfSSatish Balay   tao->ops->view = TaoView_LMVM;
294a7e14dcfSSatish Balay   tao->ops->setfromoptions = TaoSetFromOptions_LMVM;
295a7e14dcfSSatish Balay   tao->ops->destroy = TaoDestroy_LMVM;
296a7e14dcfSSatish Balay 
2973c9e27cfSGeoffrey Irving   ierr = PetscNewLog(tao,&lmP);CHKERRQ(ierr);
298a7e14dcfSSatish Balay   lmP->D = 0;
299a7e14dcfSSatish Balay   lmP->M = 0;
300a7e14dcfSSatish Balay   lmP->Xold = 0;
301a7e14dcfSSatish Balay   lmP->Gold = 0;
302a7e14dcfSSatish Balay 
303a7e14dcfSSatish Balay   tao->data = (void*)lmP;
304a7e14dcfSSatish Balay   tao->max_it = 2000;
305a7e14dcfSSatish Balay   tao->max_funcs = 4000;
306a7e14dcfSSatish Balay   tao->fatol = 1e-4;
307a7e14dcfSSatish Balay   tao->frtol = 1e-4;
308a7e14dcfSSatish Balay 
309a7e14dcfSSatish Balay   ierr = TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch);CHKERRQ(ierr);
310a7e14dcfSSatish Balay   ierr = TaoLineSearchSetType(tao->linesearch,morethuente_type);CHKERRQ(ierr);
311441846f8SBarry Smith   ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr);
312*5d527766SPatrick Farrell   ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr);
313a7e14dcfSSatish Balay   PetscFunctionReturn(0);
314a7e14dcfSSatish Balay }
315728e0ed0SBarry Smith 
316