xref: /petsc/src/tao/unconstrained/impls/lmvm/lmvm.c (revision 4d6623b4ac38e61a86b5c0cf094187b5a5ae29d2)
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 
9441846f8SBarry Smith static PetscErrorCode TaoSolve_LMVM(Tao tao)
10a7e14dcfSSatish Balay {
11a7e14dcfSSatish Balay   TAO_LMVM                     *lmP = (TAO_LMVM *)tao->data;
12a7e14dcfSSatish Balay   PetscReal                    f, fold, gdx, gnorm;
13a7e14dcfSSatish Balay   PetscReal                    step = 1.0;
14a7e14dcfSSatish Balay   PetscReal                    delta;
15a7e14dcfSSatish Balay   PetscErrorCode               ierr;
16*4d6623b4SAlp Dener   PetscInt                     stepType, nupdates;
17de6ffafeSAlp Dener   PetscBool                    recycle;
18e4cb33bbSBarry Smith   TaoConvergedReason           reason = TAO_CONTINUE_ITERATING;
19e4cb33bbSBarry Smith   TaoLineSearchConvergedReason ls_status = TAOLINESEARCH_CONTINUE_ITERATING;
20a7e14dcfSSatish Balay 
21a7e14dcfSSatish Balay   PetscFunctionBegin;
22a7e14dcfSSatish Balay 
23a7e14dcfSSatish Balay   if (tao->XL || tao->XU || tao->ops->computebounds) {
24a7e14dcfSSatish Balay     ierr = PetscPrintf(((PetscObject)tao)->comm,"WARNING: Variable bounds have been set but will be ignored by lmvm algorithm\n");CHKERRQ(ierr);
25a7e14dcfSSatish Balay   }
26a7e14dcfSSatish Balay 
27a7e14dcfSSatish Balay   /*  Check convergence criteria */
28a7e14dcfSSatish Balay   ierr = TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient);CHKERRQ(ierr);
29a9603a14SPatrick Farrell   ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr);
30a9603a14SPatrick Farrell 
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 */
45de6ffafeSAlp Dener   ierr = MatLMVMGetRecycleFlag(lmP->M, &recycle);CHKERRQ(ierr);
46de6ffafeSAlp Dener   if (!recycle) {
47a7e14dcfSSatish Balay     lmP->bfgs = 0;
48a7e14dcfSSatish Balay     lmP->sgrad = 0;
49a7e14dcfSSatish Balay     lmP->grad = 0;
50de6ffafeSAlp Dener   }
51a7e14dcfSSatish Balay 
52a7e14dcfSSatish Balay   /*  Have not converged; continue with Newton method */
53a7e14dcfSSatish Balay   while (reason == TAO_CONTINUE_ITERATING) {
54a7e14dcfSSatish Balay     /*  Compute direction */
55a7e14dcfSSatish Balay     ierr = MatLMVMUpdate(lmP->M,tao->solution,tao->gradient);CHKERRQ(ierr);
56a7e14dcfSSatish Balay     ierr = MatLMVMSolve(lmP->M, tao->gradient, lmP->D);CHKERRQ(ierr);
57a7e14dcfSSatish Balay 
58a7e14dcfSSatish Balay     /*  Check for success (descent direction) */
59a7e14dcfSSatish Balay     ierr = VecDot(lmP->D, tao->gradient, &gdx);CHKERRQ(ierr);
60a7e14dcfSSatish Balay     if ((gdx <= 0.0) || PetscIsInfOrNanReal(gdx)) {
61a7e14dcfSSatish Balay       /* Step is not descent or direction produced not a number
62a7e14dcfSSatish Balay          We can assert bfgsUpdates > 1 in this case because
63a7e14dcfSSatish Balay          the first solve produces the scaled gradient direction,
64a7e14dcfSSatish Balay          which is guaranteed to be descent
65a7e14dcfSSatish Balay 
66a7e14dcfSSatish Balay          Use steepest descent direction (scaled)
67a7e14dcfSSatish Balay       */
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       ++lmP->sgrad;
82a7e14dcfSSatish Balay       stepType = LMVM_SCALED_GRADIENT;
8387f595a5SBarry Smith     } else {
84*4d6623b4SAlp Dener       ierr = MatLMVMGetUpdates(lmP->M, &nupdates); CHKERRQ(ierr);
85*4d6623b4SAlp Dener       if (1 == nupdates) {
86a7e14dcfSSatish Balay         /*  The first BFGS direction is always the scaled gradient */
87a7e14dcfSSatish Balay         ++lmP->sgrad;
88a7e14dcfSSatish Balay         stepType = LMVM_SCALED_GRADIENT;
8987f595a5SBarry Smith       } else {
90a7e14dcfSSatish Balay         ++lmP->bfgs;
91a7e14dcfSSatish Balay         stepType = LMVM_BFGS;
92a7e14dcfSSatish Balay       }
93a7e14dcfSSatish Balay     }
94a7e14dcfSSatish Balay     ierr = VecScale(lmP->D, -1.0);CHKERRQ(ierr);
95a7e14dcfSSatish Balay 
96a7e14dcfSSatish Balay     /*  Perform the linesearch */
97a7e14dcfSSatish Balay     fold = f;
98a7e14dcfSSatish Balay     ierr = VecCopy(tao->solution, lmP->Xold);CHKERRQ(ierr);
99a7e14dcfSSatish Balay     ierr = VecCopy(tao->gradient, lmP->Gold);CHKERRQ(ierr);
100a7e14dcfSSatish Balay 
101a7e14dcfSSatish Balay     ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, lmP->D, &step,&ls_status);CHKERRQ(ierr);
102a7e14dcfSSatish Balay     ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
103a7e14dcfSSatish Balay 
10487f595a5SBarry Smith     while (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER && (stepType != LMVM_GRADIENT)) {
105a7e14dcfSSatish Balay       /*  Linesearch failed */
106a7e14dcfSSatish Balay       /*  Reset factors and use scaled gradient step */
107a7e14dcfSSatish Balay       f = fold;
108a7e14dcfSSatish Balay       ierr = VecCopy(lmP->Xold, tao->solution);CHKERRQ(ierr);
109a7e14dcfSSatish Balay       ierr = VecCopy(lmP->Gold, tao->gradient);CHKERRQ(ierr);
110a7e14dcfSSatish Balay 
111a7e14dcfSSatish Balay       switch(stepType) {
112a7e14dcfSSatish Balay       case LMVM_BFGS:
113a7e14dcfSSatish Balay         /*  Failed to obtain acceptable iterate with BFGS step */
114a7e14dcfSSatish Balay         /*  Attempt to use the scaled gradient direction */
115a7e14dcfSSatish Balay 
116a7e14dcfSSatish Balay         if (f != 0.0) {
117a7e14dcfSSatish Balay           delta = 2.0 * PetscAbsScalar(f) / (gnorm*gnorm);
11887f595a5SBarry Smith         } else {
119a7e14dcfSSatish Balay           delta = 2.0 / (gnorm*gnorm);
120a7e14dcfSSatish Balay         }
121a7e14dcfSSatish Balay         ierr = MatLMVMSetDelta(lmP->M, delta);CHKERRQ(ierr);
122a7e14dcfSSatish Balay         ierr = MatLMVMReset(lmP->M);CHKERRQ(ierr);
123a7e14dcfSSatish Balay         ierr = MatLMVMUpdate(lmP->M, tao->solution, tao->gradient);CHKERRQ(ierr);
124a7e14dcfSSatish Balay         ierr = MatLMVMSolve(lmP->M, tao->gradient, lmP->D);CHKERRQ(ierr);
125a7e14dcfSSatish Balay 
126a7e14dcfSSatish Balay         /* On a reset, the direction cannot be not a number; it is a
127a7e14dcfSSatish Balay            scaled gradient step.  No need to check for this condition. */
128*4d6623b4SAlp Dener         --lmP->bfgs;
129a7e14dcfSSatish Balay         ++lmP->sgrad;
130a7e14dcfSSatish Balay         stepType = LMVM_SCALED_GRADIENT;
131a7e14dcfSSatish Balay         break;
132a7e14dcfSSatish Balay 
133a7e14dcfSSatish Balay       case LMVM_SCALED_GRADIENT:
134a7e14dcfSSatish Balay         /* The scaled gradient step did not produce a new iterate;
135a7e14dcfSSatish Balay            attempt to use the gradient direction.
136a7e14dcfSSatish Balay            Need to make sure we are not using a different diagonal scaling */
137a7e14dcfSSatish Balay         ierr = MatLMVMSetDelta(lmP->M, 1.0);CHKERRQ(ierr);
138a7e14dcfSSatish Balay         ierr = MatLMVMReset(lmP->M);CHKERRQ(ierr);
139a7e14dcfSSatish Balay         ierr = MatLMVMUpdate(lmP->M, tao->solution, tao->gradient);CHKERRQ(ierr);
140a7e14dcfSSatish Balay         ierr = MatLMVMSolve(lmP->M, tao->gradient, lmP->D);CHKERRQ(ierr);
141a7e14dcfSSatish Balay 
142*4d6623b4SAlp Dener         --lmP->sgrad;
143a7e14dcfSSatish Balay         ++lmP->grad;
144a7e14dcfSSatish Balay         stepType = LMVM_GRADIENT;
145a7e14dcfSSatish Balay         break;
146a7e14dcfSSatish Balay       }
147a7e14dcfSSatish Balay       ierr = VecScale(lmP->D, -1.0);CHKERRQ(ierr);
148a7e14dcfSSatish Balay 
149a7e14dcfSSatish Balay       /*  Perform the linesearch */
150a7e14dcfSSatish Balay       ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &f, tao->gradient, lmP->D, &step, &ls_status);CHKERRQ(ierr);
151a7e14dcfSSatish Balay       ierr = TaoAddLineSearchCounts(tao);CHKERRQ(ierr);
152a7e14dcfSSatish Balay     }
153a7e14dcfSSatish Balay 
15487f595a5SBarry Smith     if (ls_status != TAOLINESEARCH_SUCCESS && ls_status != TAOLINESEARCH_SUCCESS_USER) {
155a7e14dcfSSatish Balay       /*  Failed to find an improving point */
156a7e14dcfSSatish Balay       f = fold;
157a7e14dcfSSatish Balay       ierr = VecCopy(lmP->Xold, tao->solution);CHKERRQ(ierr);
158a7e14dcfSSatish Balay       ierr = VecCopy(lmP->Gold, tao->gradient);CHKERRQ(ierr);
159a7e14dcfSSatish Balay       step = 0.0;
160a7e14dcfSSatish Balay       reason = TAO_DIVERGED_LS_FAILURE;
161a7e14dcfSSatish Balay       tao->reason = TAO_DIVERGED_LS_FAILURE;
162a7e14dcfSSatish Balay     }
163a9603a14SPatrick Farrell 
164a7e14dcfSSatish Balay     /*  Check for termination */
165a9603a14SPatrick Farrell     ierr = TaoGradientNorm(tao, tao->gradient,NORM_2,&gnorm);CHKERRQ(ierr);
166a9603a14SPatrick Farrell 
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 
173441846f8SBarry Smith static PetscErrorCode TaoSetUp_LMVM(Tao tao)
174a7e14dcfSSatish Balay {
175a7e14dcfSSatish Balay   TAO_LMVM       *lmP = (TAO_LMVM *)tao->data;
176a7e14dcfSSatish Balay   PetscInt       n,N;
177a7e14dcfSSatish Balay   PetscErrorCode ierr;
178a9603a14SPatrick Farrell   KSP            H0ksp;
179a7e14dcfSSatish Balay 
180a7e14dcfSSatish Balay   PetscFunctionBegin;
181a7e14dcfSSatish Balay   /* Existence of tao->solution checked in TaoSetUp() */
182a7e14dcfSSatish Balay   if (!tao->gradient) {ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);  }
183a7e14dcfSSatish Balay   if (!tao->stepdirection) {ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr);  }
184a7e14dcfSSatish Balay   if (!lmP->D) {ierr = VecDuplicate(tao->solution,&lmP->D);CHKERRQ(ierr);  }
185a7e14dcfSSatish Balay   if (!lmP->Xold) {ierr = VecDuplicate(tao->solution,&lmP->Xold);CHKERRQ(ierr);  }
186a7e14dcfSSatish Balay   if (!lmP->Gold) {ierr = VecDuplicate(tao->solution,&lmP->Gold);CHKERRQ(ierr);  }
187a7e14dcfSSatish Balay 
188a7e14dcfSSatish Balay   /*  Create matrix for the limited memory approximation */
189a7e14dcfSSatish Balay   ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr);
190a7e14dcfSSatish Balay   ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr);
191a7e14dcfSSatish Balay   ierr = MatCreateLMVM(((PetscObject)tao)->comm,n,N,&lmP->M);CHKERRQ(ierr);
192a7e14dcfSSatish Balay   ierr = MatLMVMAllocateVectors(lmP->M,tao->solution);CHKERRQ(ierr);
193a9603a14SPatrick Farrell 
194a9603a14SPatrick Farrell   /* If the user has set a matrix to solve as the initial H0, set the options prefix here, and set up the KSP */
195a9603a14SPatrick Farrell   if (lmP->H0) {
196a9603a14SPatrick Farrell     const char *prefix;
197a9603a14SPatrick Farrell     PC H0pc;
198a9603a14SPatrick Farrell 
199a9603a14SPatrick Farrell     ierr = MatLMVMSetH0(lmP->M, lmP->H0);CHKERRQ(ierr);
200a9603a14SPatrick Farrell     ierr = MatLMVMGetH0KSP(lmP->M, &H0ksp);CHKERRQ(ierr);
201a9603a14SPatrick Farrell 
202a9603a14SPatrick Farrell     ierr = TaoGetOptionsPrefix(tao, &prefix);CHKERRQ(ierr);
203a9603a14SPatrick Farrell     ierr = KSPSetOptionsPrefix(H0ksp, prefix);CHKERRQ(ierr);
204a9603a14SPatrick Farrell     ierr = PetscObjectAppendOptionsPrefix((PetscObject)H0ksp, "tao_h0_");CHKERRQ(ierr);
205a9603a14SPatrick Farrell     ierr = KSPGetPC(H0ksp, &H0pc);CHKERRQ(ierr);
206a9603a14SPatrick Farrell     ierr = PetscObjectAppendOptionsPrefix((PetscObject)H0pc,  "tao_h0_");CHKERRQ(ierr);
207a9603a14SPatrick Farrell 
208a9603a14SPatrick Farrell     ierr = KSPSetFromOptions(H0ksp);CHKERRQ(ierr);
209a9603a14SPatrick Farrell     ierr = KSPSetUp(H0ksp);CHKERRQ(ierr);
210a9603a14SPatrick Farrell   }
211a9603a14SPatrick Farrell 
212a7e14dcfSSatish Balay   PetscFunctionReturn(0);
213a7e14dcfSSatish Balay }
214a7e14dcfSSatish Balay 
215a7e14dcfSSatish Balay /* ---------------------------------------------------------- */
216441846f8SBarry Smith static PetscErrorCode TaoDestroy_LMVM(Tao tao)
217a7e14dcfSSatish Balay {
218a7e14dcfSSatish Balay   TAO_LMVM       *lmP = (TAO_LMVM *)tao->data;
219a7e14dcfSSatish Balay   PetscErrorCode ierr;
220a7e14dcfSSatish Balay 
221a7e14dcfSSatish Balay   PetscFunctionBegin;
222a7e14dcfSSatish Balay   if (tao->setupcalled) {
223a7e14dcfSSatish Balay     ierr = VecDestroy(&lmP->Xold);CHKERRQ(ierr);
224a7e14dcfSSatish Balay     ierr = VecDestroy(&lmP->Gold);CHKERRQ(ierr);
225a7e14dcfSSatish Balay     ierr = VecDestroy(&lmP->D);CHKERRQ(ierr);
226a7e14dcfSSatish Balay     ierr = MatDestroy(&lmP->M);CHKERRQ(ierr);
227a7e14dcfSSatish Balay   }
228a9603a14SPatrick Farrell 
229a9603a14SPatrick Farrell   if (lmP->H0) {
230a9603a14SPatrick Farrell     ierr = PetscObjectDereference((PetscObject)lmP->H0);CHKERRQ(ierr);
231a9603a14SPatrick Farrell   }
232a9603a14SPatrick Farrell 
233a7e14dcfSSatish Balay   ierr = PetscFree(tao->data);CHKERRQ(ierr);
234a9603a14SPatrick Farrell 
235a7e14dcfSSatish Balay   PetscFunctionReturn(0);
236a7e14dcfSSatish Balay }
237a7e14dcfSSatish Balay 
238a7e14dcfSSatish Balay /*------------------------------------------------------------*/
2394416b707SBarry Smith static PetscErrorCode TaoSetFromOptions_LMVM(PetscOptionItems *PetscOptionsObject,Tao tao)
240a7e14dcfSSatish Balay {
241a7e14dcfSSatish Balay   PetscErrorCode ierr;
242a7e14dcfSSatish Balay 
243a7e14dcfSSatish Balay   PetscFunctionBegin;
2441a1499c8SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Limited-memory variable-metric method for unconstrained optimization");CHKERRQ(ierr);
245114d2d62SAlp Dener   ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr);
246288b7216SAlp Dener   ierr = PetscOptionsTail();CHKERRQ(ierr);
247a7e14dcfSSatish Balay   PetscFunctionReturn(0);
248a7e14dcfSSatish Balay }
249a7e14dcfSSatish Balay 
250a7e14dcfSSatish Balay /*------------------------------------------------------------*/
251441846f8SBarry Smith static PetscErrorCode TaoView_LMVM(Tao tao, PetscViewer viewer)
252a7e14dcfSSatish Balay {
253a7e14dcfSSatish Balay   TAO_LMVM       *lm = (TAO_LMVM *)tao->data;
254de6ffafeSAlp Dener   PetscBool      isascii, recycle;
255*4d6623b4SAlp Dener   PetscInt       recycled_its;
256a7e14dcfSSatish Balay   PetscErrorCode ierr;
257a7e14dcfSSatish Balay 
258a7e14dcfSSatish Balay   PetscFunctionBegin;
259a7e14dcfSSatish Balay   ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii);CHKERRQ(ierr);
260a7e14dcfSSatish Balay   if (isascii) {
261a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
262a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "BFGS steps: %D\n", lm->bfgs);CHKERRQ(ierr);
263a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "Scaled gradient steps: %D\n", lm->sgrad);CHKERRQ(ierr);
264a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPrintf(viewer, "Gradient steps: %D\n", lm->grad);CHKERRQ(ierr);
265de6ffafeSAlp Dener     ierr = MatLMVMGetRecycleFlag(lm->M, &recycle);CHKERRQ(ierr);
266de6ffafeSAlp Dener     if (recycle) {
267288b7216SAlp Dener       ierr = PetscViewerASCIIPrintf(viewer, "Recycle: on\n");CHKERRQ(ierr);
268*4d6623b4SAlp Dener       recycled_its = lm->bfgs + lm->sgrad + lm->grad;
269*4d6623b4SAlp Dener       ierr = PetscViewerASCIIPrintf(viewer, "Total recycled iterations: %D\n", recycled_its);CHKERRQ(ierr);
270a0bfee83SAlp Dener     }
271a7e14dcfSSatish Balay     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
272a7e14dcfSSatish Balay   }
273a7e14dcfSSatish Balay   PetscFunctionReturn(0);
274a7e14dcfSSatish Balay }
275a7e14dcfSSatish Balay 
276a7e14dcfSSatish Balay /* ---------------------------------------------------------- */
277a7e14dcfSSatish Balay 
2784aa34175SJason Sarich /*MC
2794aa34175SJason Sarich      TAOLMVM - Limited Memory Variable Metric method is a quasi-Newton
2804aa34175SJason Sarich      optimization solver for unconstrained minimization. It solves
2814aa34175SJason Sarich      the Newton step
2824aa34175SJason Sarich               Hkdk = - gk
2834aa34175SJason Sarich 
2844aa34175SJason Sarich      using an approximation Bk in place of Hk, where Bk is composed using
2854aa34175SJason Sarich      the BFGS update formula. A More-Thuente line search is then used
2864aa34175SJason Sarich      to computed the steplength in the dk direction
2874aa34175SJason Sarich   Options Database Keys:
2884aa34175SJason Sarich +     -tao_lmm_vectors - number of vectors to use for approximation
2894aa34175SJason Sarich .     -tao_lmm_scale_type - "none","scalar","broyden"
2904aa34175SJason Sarich .     -tao_lmm_limit_type - "none","average","relative","absolute"
2914aa34175SJason Sarich .     -tao_lmm_rescale_type - "none","scalar","gl"
2924aa34175SJason Sarich .     -tao_lmm_limit_mu - mu limiting factor
2934aa34175SJason Sarich .     -tao_lmm_limit_nu - nu limiting factor
2944aa34175SJason Sarich .     -tao_lmm_delta_min - minimum delta value
2954aa34175SJason Sarich .     -tao_lmm_delta_max - maximum delta value
2964aa34175SJason Sarich .     -tao_lmm_broyden_phi - phi factor for Broyden scaling
2974aa34175SJason Sarich .     -tao_lmm_scalar_alpha - alpha factor for scalar scaling
2984aa34175SJason Sarich .     -tao_lmm_rescale_alpha - alpha factor for rescaling diagonal
2994aa34175SJason Sarich .     -tao_lmm_rescale_beta - beta factor for rescaling diagonal
3004aa34175SJason Sarich .     -tao_lmm_scalar_history - amount of history for scalar scaling
3014aa34175SJason Sarich .     -tao_lmm_rescale_history - amount of history for rescaling diagonal
3024aa34175SJason Sarich -     -tao_lmm_eps - rejection tolerance
3034aa34175SJason Sarich 
3041eb8069cSJason Sarich   Level: beginner
3054aa34175SJason Sarich M*/
3064aa34175SJason Sarich 
307728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate_LMVM(Tao tao)
308a7e14dcfSSatish Balay {
309a7e14dcfSSatish Balay   TAO_LMVM       *lmP;
3108caf6e8cSBarry Smith   const char     *morethuente_type = TAOLINESEARCHMT;
311a7e14dcfSSatish Balay   PetscErrorCode ierr;
312a7e14dcfSSatish Balay 
313a7e14dcfSSatish Balay   PetscFunctionBegin;
314a7e14dcfSSatish Balay   tao->ops->setup = TaoSetUp_LMVM;
315a7e14dcfSSatish Balay   tao->ops->solve = TaoSolve_LMVM;
316a7e14dcfSSatish Balay   tao->ops->view = TaoView_LMVM;
317a7e14dcfSSatish Balay   tao->ops->setfromoptions = TaoSetFromOptions_LMVM;
318a7e14dcfSSatish Balay   tao->ops->destroy = TaoDestroy_LMVM;
319a7e14dcfSSatish Balay 
3203c9e27cfSGeoffrey Irving   ierr = PetscNewLog(tao,&lmP);CHKERRQ(ierr);
321a7e14dcfSSatish Balay   lmP->D = 0;
322a7e14dcfSSatish Balay   lmP->M = 0;
323a7e14dcfSSatish Balay   lmP->Xold = 0;
324a7e14dcfSSatish Balay   lmP->Gold = 0;
325a9603a14SPatrick Farrell   lmP->H0   = NULL;
326a7e14dcfSSatish Balay 
327a7e14dcfSSatish Balay   tao->data = (void*)lmP;
3286552cf8aSJason Sarich   /* Override default settings (unless already changed) */
3296552cf8aSJason Sarich   if (!tao->max_it_changed) tao->max_it = 2000;
3306552cf8aSJason Sarich   if (!tao->max_funcs_changed) tao->max_funcs = 4000;
331a7e14dcfSSatish Balay 
332a7e14dcfSSatish Balay   ierr = TaoLineSearchCreate(((PetscObject)tao)->comm,&tao->linesearch);CHKERRQ(ierr);
33363b15415SAlp Dener   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr);
334a7e14dcfSSatish Balay   ierr = TaoLineSearchSetType(tao->linesearch,morethuente_type);CHKERRQ(ierr);
335441846f8SBarry Smith   ierr = TaoLineSearchUseTaoRoutines(tao->linesearch,tao);CHKERRQ(ierr);
3365d527766SPatrick Farrell   ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr);
337a7e14dcfSSatish Balay   PetscFunctionReturn(0);
338a7e14dcfSSatish Balay }
339