xref: /petsc/src/tao/unconstrained/impls/owlqn/owlqn.c (revision 5a884c48ab0c46bab83cd9bb8710f380fa6d8bcf)
1ba92ff59SBarry Smith #include <petsctaolinesearch.h>
2aaa7dc30SBarry Smith #include <../src/tao/unconstrained/impls/owlqn/owlqn.h>
3a7e14dcfSSatish Balay 
4a7e14dcfSSatish Balay #define OWLQN_BFGS            0
5a7e14dcfSSatish Balay #define OWLQN_SCALED_GRADIENT 1
6a7e14dcfSSatish Balay #define OWLQN_GRADIENT        2
7a7e14dcfSSatish Balay 
8d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProjDirect_OWLQN(Vec d, Vec g)
9d71ae5a4SJacob Faibussowitsch {
105e081366SBarry Smith   const PetscReal *gptr;
115e081366SBarry Smith   PetscReal       *dptr;
12a7e14dcfSSatish Balay   PetscInt         low, high, low1, high1, i;
13a7e14dcfSSatish Balay 
14a7e14dcfSSatish Balay   PetscFunctionBegin;
159566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(d, &low, &high));
169566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(g, &low1, &high1));
17a7e14dcfSSatish Balay 
189566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(g, &gptr));
199566063dSJacob Faibussowitsch   PetscCall(VecGetArray(d, &dptr));
20a7e14dcfSSatish Balay   for (i = 0; i < high - low; i++) {
21ad540459SPierre Jolivet     if (dptr[i] * gptr[i] <= 0.0) dptr[i] = 0.0;
22a7e14dcfSSatish Balay   }
239566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(d, &dptr));
249566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(g, &gptr));
253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
26a7e14dcfSSatish Balay }
27a7e14dcfSSatish Balay 
28d71ae5a4SJacob Faibussowitsch static PetscErrorCode ComputePseudoGrad_OWLQN(Vec x, Vec gv, PetscReal lambda)
29d71ae5a4SJacob Faibussowitsch {
305e081366SBarry Smith   const PetscReal *xptr;
315e081366SBarry Smith   PetscReal       *gptr;
32a7e14dcfSSatish Balay   PetscInt         low, high, low1, high1, i;
33a7e14dcfSSatish Balay 
34a7e14dcfSSatish Balay   PetscFunctionBegin;
359566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(x, &low, &high));
369566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(gv, &low1, &high1));
37a7e14dcfSSatish Balay 
389566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(x, &xptr));
399566063dSJacob Faibussowitsch   PetscCall(VecGetArray(gv, &gptr));
40a7e14dcfSSatish Balay   for (i = 0; i < high - low; i++) {
4153506e15SBarry Smith     if (xptr[i] < 0.0) gptr[i] = gptr[i] - lambda;
4253506e15SBarry Smith     else if (xptr[i] > 0.0) gptr[i] = gptr[i] + lambda;
4353506e15SBarry Smith     else if (gptr[i] + lambda < 0.0) gptr[i] = gptr[i] + lambda;
4453506e15SBarry Smith     else if (gptr[i] - lambda > 0.0) gptr[i] = gptr[i] - lambda;
4553506e15SBarry Smith     else gptr[i] = 0.0;
46a7e14dcfSSatish Balay   }
479566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(gv, &gptr));
489566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(x, &xptr));
493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
50a7e14dcfSSatish Balay }
51a7e14dcfSSatish Balay 
52d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSolve_OWLQN(Tao tao)
53d71ae5a4SJacob Faibussowitsch {
54a7e14dcfSSatish Balay   TAO_OWLQN                   *lmP = (TAO_OWLQN *)tao->data;
55a7e14dcfSSatish Balay   PetscReal                    f, fold, gdx, gnorm;
56a7e14dcfSSatish Balay   PetscReal                    step = 1.0;
57a7e14dcfSSatish Balay   PetscReal                    delta;
58a7e14dcfSSatish Balay   PetscInt                     stepType;
59e4cb33bbSBarry Smith   TaoLineSearchConvergedReason ls_status = TAOLINESEARCH_CONTINUE_ITERATING;
60a7e14dcfSSatish Balay 
61a7e14dcfSSatish Balay   PetscFunctionBegin;
6248a46eb9SPierre Jolivet   if (tao->XL || tao->XU || tao->ops->computebounds) PetscCall(PetscInfo(tao, "WARNING: Variable bounds have been set but will be ignored by owlqn algorithm\n"));
63a7e14dcfSSatish Balay 
64a7e14dcfSSatish Balay   /* Check convergence criteria */
659566063dSJacob Faibussowitsch   PetscCall(TaoComputeObjectiveAndGradient(tao, tao->solution, &f, tao->gradient));
669566063dSJacob Faibussowitsch   PetscCall(VecCopy(tao->gradient, lmP->GV));
679566063dSJacob Faibussowitsch   PetscCall(ComputePseudoGrad_OWLQN(tao->solution, lmP->GV, lmP->lambda));
689566063dSJacob Faibussowitsch   PetscCall(VecNorm(lmP->GV, NORM_2, &gnorm));
6976c63389SBarry Smith   PetscCheck(!PetscIsInfOrNanReal(f) && !PetscIsInfOrNanReal(gnorm), PetscObjectComm((PetscObject)tao), PETSC_ERR_USER, "User provided compute function generated infinity or NaN");
70a7e14dcfSSatish Balay 
713ecd9318SAlp Dener   tao->reason = TAO_CONTINUE_ITERATING;
729566063dSJacob Faibussowitsch   PetscCall(TaoLogConvergenceHistory(tao, f, gnorm, 0.0, tao->ksp_its));
73*efab7034SHan Liu   PetscCall(TaoMonitor(tao, tao->niter, f, gnorm, 0.0, step));
74dbbe0bcdSBarry Smith   PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
753ba16761SJacob Faibussowitsch   if (tao->reason != TAO_CONTINUE_ITERATING) PetscFunctionReturn(PETSC_SUCCESS);
76a7e14dcfSSatish Balay 
77a7e14dcfSSatish Balay   /* Set initial scaling for the function */
78cd929ea3SAlp Dener   delta = 2.0 * PetscMax(1.0, PetscAbsScalar(f)) / (gnorm * gnorm);
799566063dSJacob Faibussowitsch   PetscCall(MatLMVMSetJ0Scale(lmP->M, delta));
80a7e14dcfSSatish Balay 
81a7e14dcfSSatish Balay   /* Set counter for gradient/reset steps */
82a7e14dcfSSatish Balay   lmP->bfgs  = 0;
83a7e14dcfSSatish Balay   lmP->sgrad = 0;
84a7e14dcfSSatish Balay   lmP->grad  = 0;
85a7e14dcfSSatish Balay 
86a7e14dcfSSatish Balay   /* Have not converged; continue with Newton method */
873ecd9318SAlp Dener   while (tao->reason == TAO_CONTINUE_ITERATING) {
88e1e80dc8SAlp Dener     /* Call general purpose update function */
89dbbe0bcdSBarry Smith     PetscTryTypeMethod(tao, update, tao->niter, tao->user_update);
90e1e80dc8SAlp Dener 
91a7e14dcfSSatish Balay     /* Compute direction */
929566063dSJacob Faibussowitsch     PetscCall(MatLMVMUpdate(lmP->M, tao->solution, tao->gradient));
939566063dSJacob Faibussowitsch     PetscCall(MatSolve(lmP->M, lmP->GV, lmP->D));
94a7e14dcfSSatish Balay 
959566063dSJacob Faibussowitsch     PetscCall(ProjDirect_OWLQN(lmP->D, lmP->GV));
96a7e14dcfSSatish Balay 
97a7e14dcfSSatish Balay     ++lmP->bfgs;
98a7e14dcfSSatish Balay 
99a7e14dcfSSatish Balay     /* Check for success (descent direction) */
1009566063dSJacob Faibussowitsch     PetscCall(VecDot(lmP->D, lmP->GV, &gdx));
101a7e14dcfSSatish Balay     if ((gdx <= 0.0) || PetscIsInfOrNanReal(gdx)) {
102a7e14dcfSSatish Balay       /* Step is not descent or direction produced not a number
103a7e14dcfSSatish Balay          We can assert bfgsUpdates > 1 in this case because
104a7e14dcfSSatish Balay          the first solve produces the scaled gradient direction,
105a7e14dcfSSatish Balay          which is guaranteed to be descent
106a7e14dcfSSatish Balay 
107a7e14dcfSSatish Balay          Use steepest descent direction (scaled) */
108a7e14dcfSSatish Balay       ++lmP->grad;
109a7e14dcfSSatish Balay 
110cd929ea3SAlp Dener       delta = 2.0 * PetscMax(1.0, PetscAbsScalar(f)) / (gnorm * gnorm);
1119566063dSJacob Faibussowitsch       PetscCall(MatLMVMSetJ0Scale(lmP->M, delta));
1129566063dSJacob Faibussowitsch       PetscCall(MatLMVMReset(lmP->M, PETSC_FALSE));
1139566063dSJacob Faibussowitsch       PetscCall(MatLMVMUpdate(lmP->M, tao->solution, tao->gradient));
1149566063dSJacob Faibussowitsch       PetscCall(MatSolve(lmP->M, lmP->GV, lmP->D));
115a7e14dcfSSatish Balay 
1169566063dSJacob Faibussowitsch       PetscCall(ProjDirect_OWLQN(lmP->D, lmP->GV));
117a7e14dcfSSatish Balay 
118a7e14dcfSSatish Balay       lmP->bfgs = 1;
119a7e14dcfSSatish Balay       ++lmP->sgrad;
120a7e14dcfSSatish Balay       stepType = OWLQN_SCALED_GRADIENT;
12153506e15SBarry Smith     } else {
122a7e14dcfSSatish Balay       if (1 == lmP->bfgs) {
123a7e14dcfSSatish Balay         /* The first BFGS direction is always the scaled gradient */
124a7e14dcfSSatish Balay         ++lmP->sgrad;
125a7e14dcfSSatish Balay         stepType = OWLQN_SCALED_GRADIENT;
12653506e15SBarry Smith       } else {
127a7e14dcfSSatish Balay         ++lmP->bfgs;
128a7e14dcfSSatish Balay         stepType = OWLQN_BFGS;
129a7e14dcfSSatish Balay       }
130a7e14dcfSSatish Balay     }
131a7e14dcfSSatish Balay 
1329566063dSJacob Faibussowitsch     PetscCall(VecScale(lmP->D, -1.0));
133a7e14dcfSSatish Balay 
134a7e14dcfSSatish Balay     /* Perform the linesearch */
135a7e14dcfSSatish Balay     fold = f;
1369566063dSJacob Faibussowitsch     PetscCall(VecCopy(tao->solution, lmP->Xold));
1379566063dSJacob Faibussowitsch     PetscCall(VecCopy(tao->gradient, lmP->Gold));
138a7e14dcfSSatish Balay 
1399566063dSJacob Faibussowitsch     PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &f, lmP->GV, lmP->D, &step, &ls_status));
1409566063dSJacob Faibussowitsch     PetscCall(TaoAddLineSearchCounts(tao));
141a7e14dcfSSatish Balay 
142a7e14dcfSSatish Balay     while (((int)ls_status < 0) && (stepType != OWLQN_GRADIENT)) {
143a7e14dcfSSatish Balay       /* Reset factors and use scaled gradient step */
144a7e14dcfSSatish Balay       f = fold;
1459566063dSJacob Faibussowitsch       PetscCall(VecCopy(lmP->Xold, tao->solution));
1469566063dSJacob Faibussowitsch       PetscCall(VecCopy(lmP->Gold, tao->gradient));
1479566063dSJacob Faibussowitsch       PetscCall(VecCopy(tao->gradient, lmP->GV));
148a7e14dcfSSatish Balay 
1499566063dSJacob Faibussowitsch       PetscCall(ComputePseudoGrad_OWLQN(tao->solution, lmP->GV, lmP->lambda));
150a7e14dcfSSatish Balay 
151a7e14dcfSSatish Balay       switch (stepType) {
152a7e14dcfSSatish Balay       case OWLQN_BFGS:
153a7e14dcfSSatish Balay         /* Failed to obtain acceptable iterate with BFGS step
154a7e14dcfSSatish Balay            Attempt to use the scaled gradient direction */
155a7e14dcfSSatish Balay 
156cd929ea3SAlp Dener         delta = 2.0 * PetscMax(1.0, PetscAbsScalar(f)) / (gnorm * gnorm);
1579566063dSJacob Faibussowitsch         PetscCall(MatLMVMSetJ0Scale(lmP->M, delta));
1589566063dSJacob Faibussowitsch         PetscCall(MatLMVMReset(lmP->M, PETSC_FALSE));
1599566063dSJacob Faibussowitsch         PetscCall(MatLMVMUpdate(lmP->M, tao->solution, tao->gradient));
1609566063dSJacob Faibussowitsch         PetscCall(MatSolve(lmP->M, lmP->GV, lmP->D));
161a7e14dcfSSatish Balay 
1629566063dSJacob Faibussowitsch         PetscCall(ProjDirect_OWLQN(lmP->D, lmP->GV));
163a7e14dcfSSatish Balay 
164a7e14dcfSSatish Balay         lmP->bfgs = 1;
165a7e14dcfSSatish Balay         ++lmP->sgrad;
166a7e14dcfSSatish Balay         stepType = OWLQN_SCALED_GRADIENT;
167a7e14dcfSSatish Balay         break;
168a7e14dcfSSatish Balay 
169a7e14dcfSSatish Balay       case OWLQN_SCALED_GRADIENT:
170a7e14dcfSSatish Balay         /* The scaled gradient step did not produce a new iterate;
171a7e14dcfSSatish Balay            attempt to use the gradient direction.
172a7e14dcfSSatish Balay            Need to make sure we are not using a different diagonal scaling */
1739566063dSJacob Faibussowitsch         PetscCall(MatLMVMSetJ0Scale(lmP->M, 1.0));
1749566063dSJacob Faibussowitsch         PetscCall(MatLMVMReset(lmP->M, PETSC_FALSE));
1759566063dSJacob Faibussowitsch         PetscCall(MatLMVMUpdate(lmP->M, tao->solution, tao->gradient));
1769566063dSJacob Faibussowitsch         PetscCall(MatSolve(lmP->M, lmP->GV, lmP->D));
177a7e14dcfSSatish Balay 
1789566063dSJacob Faibussowitsch         PetscCall(ProjDirect_OWLQN(lmP->D, lmP->GV));
179a7e14dcfSSatish Balay 
180a7e14dcfSSatish Balay         lmP->bfgs = 1;
181a7e14dcfSSatish Balay         ++lmP->grad;
182a7e14dcfSSatish Balay         stepType = OWLQN_GRADIENT;
183a7e14dcfSSatish Balay         break;
184a7e14dcfSSatish Balay       }
1859566063dSJacob Faibussowitsch       PetscCall(VecScale(lmP->D, -1.0));
186a7e14dcfSSatish Balay 
187a7e14dcfSSatish Balay       /* Perform the linesearch */
1889566063dSJacob Faibussowitsch       PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &f, lmP->GV, lmP->D, &step, &ls_status));
1899566063dSJacob Faibussowitsch       PetscCall(TaoAddLineSearchCounts(tao));
190a7e14dcfSSatish Balay     }
191a7e14dcfSSatish Balay 
192a7e14dcfSSatish Balay     if ((int)ls_status < 0) {
193a7e14dcfSSatish Balay       /* Failed to find an improving point*/
194a7e14dcfSSatish Balay       f = fold;
1959566063dSJacob Faibussowitsch       PetscCall(VecCopy(lmP->Xold, tao->solution));
1969566063dSJacob Faibussowitsch       PetscCall(VecCopy(lmP->Gold, tao->gradient));
1979566063dSJacob Faibussowitsch       PetscCall(VecCopy(tao->gradient, lmP->GV));
198a7e14dcfSSatish Balay       step = 0.0;
19953506e15SBarry Smith     } else {
200a7e14dcfSSatish Balay       /* a little hack here, because that gv is used to store g */
2019566063dSJacob Faibussowitsch       PetscCall(VecCopy(lmP->GV, tao->gradient));
202a7e14dcfSSatish Balay     }
203a7e14dcfSSatish Balay 
2049566063dSJacob Faibussowitsch     PetscCall(ComputePseudoGrad_OWLQN(tao->solution, lmP->GV, lmP->lambda));
205a7e14dcfSSatish Balay 
206a7e14dcfSSatish Balay     /* Check for termination */
207a7e14dcfSSatish Balay 
2089566063dSJacob Faibussowitsch     PetscCall(VecNorm(lmP->GV, NORM_2, &gnorm));
209a7e14dcfSSatish Balay 
210*efab7034SHan Liu     ++tao->niter;
2119566063dSJacob Faibussowitsch     PetscCall(TaoLogConvergenceHistory(tao, f, gnorm, 0.0, tao->ksp_its));
212*efab7034SHan Liu     PetscCall(TaoMonitor(tao, tao->niter, f, gnorm, 0.0, step));
213dbbe0bcdSBarry Smith     PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
214a7e14dcfSSatish Balay 
21553506e15SBarry Smith     if ((int)ls_status < 0) break;
216a7e14dcfSSatish Balay   }
2173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
218a7e14dcfSSatish Balay }
219a7e14dcfSSatish Balay 
220d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSetUp_OWLQN(Tao tao)
221d71ae5a4SJacob Faibussowitsch {
222a7e14dcfSSatish Balay   TAO_OWLQN *lmP = (TAO_OWLQN *)tao->data;
223a7e14dcfSSatish Balay   PetscInt   n, N;
224a7e14dcfSSatish Balay 
225a7e14dcfSSatish Balay   PetscFunctionBegin;
226441846f8SBarry Smith   /* Existence of tao->solution checked in TaoSetUp() */
2279566063dSJacob Faibussowitsch   if (!tao->gradient) PetscCall(VecDuplicate(tao->solution, &tao->gradient));
2289566063dSJacob Faibussowitsch   if (!tao->stepdirection) PetscCall(VecDuplicate(tao->solution, &tao->stepdirection));
2299566063dSJacob Faibussowitsch   if (!lmP->D) PetscCall(VecDuplicate(tao->solution, &lmP->D));
2309566063dSJacob Faibussowitsch   if (!lmP->GV) PetscCall(VecDuplicate(tao->solution, &lmP->GV));
2319566063dSJacob Faibussowitsch   if (!lmP->Xold) PetscCall(VecDuplicate(tao->solution, &lmP->Xold));
2329566063dSJacob Faibussowitsch   if (!lmP->Gold) PetscCall(VecDuplicate(tao->solution, &lmP->Gold));
233a7e14dcfSSatish Balay 
234a7e14dcfSSatish Balay   /* Create matrix for the limited memory approximation */
2359566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(tao->solution, &n));
2369566063dSJacob Faibussowitsch   PetscCall(VecGetSize(tao->solution, &N));
2379566063dSJacob Faibussowitsch   PetscCall(MatCreateLMVMBFGS(((PetscObject)tao)->comm, n, N, &lmP->M));
2389566063dSJacob Faibussowitsch   PetscCall(MatLMVMAllocate(lmP->M, tao->solution, tao->gradient));
2393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
240a7e14dcfSSatish Balay }
241a7e14dcfSSatish Balay 
242d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoDestroy_OWLQN(Tao tao)
243d71ae5a4SJacob Faibussowitsch {
244a7e14dcfSSatish Balay   TAO_OWLQN *lmP = (TAO_OWLQN *)tao->data;
245a7e14dcfSSatish Balay 
246a7e14dcfSSatish Balay   PetscFunctionBegin;
247a7e14dcfSSatish Balay   if (tao->setupcalled) {
2489566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&lmP->Xold));
2499566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&lmP->Gold));
2509566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&lmP->D));
2519566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&lmP->M));
2529566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&lmP->GV));
253a7e14dcfSSatish Balay   }
2549566063dSJacob Faibussowitsch   PetscCall(PetscFree(tao->data));
2553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
256a7e14dcfSSatish Balay }
257a7e14dcfSSatish Balay 
258ce78bad3SBarry Smith static PetscErrorCode TaoSetFromOptions_OWLQN(Tao tao, PetscOptionItems PetscOptionsObject)
259d71ae5a4SJacob Faibussowitsch {
260a7e14dcfSSatish Balay   TAO_OWLQN *lmP = (TAO_OWLQN *)tao->data;
261a7e14dcfSSatish Balay 
262a7e14dcfSSatish Balay   PetscFunctionBegin;
263d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Orthant-Wise Limited-memory method for Quasi-Newton unconstrained optimization");
2649566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_owlqn_lambda", "regulariser weight", "", 100, &lmP->lambda, NULL));
265d0609cedSBarry Smith   PetscOptionsHeadEnd();
2669566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetFromOptions(tao->linesearch));
2673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
268a7e14dcfSSatish Balay }
269a7e14dcfSSatish Balay 
270d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoView_OWLQN(Tao tao, PetscViewer viewer)
271d71ae5a4SJacob Faibussowitsch {
272a7e14dcfSSatish Balay   TAO_OWLQN *lm = (TAO_OWLQN *)tao->data;
273a7e14dcfSSatish Balay   PetscBool  isascii;
274a7e14dcfSSatish Balay 
275a7e14dcfSSatish Balay   PetscFunctionBegin;
2769566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
277a7e14dcfSSatish Balay   if (isascii) {
2789566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
27963a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "BFGS steps: %" PetscInt_FMT "\n", lm->bfgs));
28063a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Scaled gradient steps: %" PetscInt_FMT "\n", lm->sgrad));
28163a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Gradient steps: %" PetscInt_FMT "\n", lm->grad));
2829566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
283a7e14dcfSSatish Balay   }
2843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
285a7e14dcfSSatish Balay }
286a7e14dcfSSatish Balay 
2871522df2eSJason Sarich /*MC
2881522df2eSJason Sarich   TAOOWLQN - orthant-wise limited memory quasi-newton algorithm
2891522df2eSJason Sarich 
2901522df2eSJason Sarich . - tao_owlqn_lambda - regulariser weight
2911522df2eSJason Sarich 
2921eb8069cSJason Sarich   Level: beginner
2931522df2eSJason Sarich M*/
2941522df2eSJason Sarich 
295d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TaoCreate_OWLQN(Tao tao)
296d71ae5a4SJacob Faibussowitsch {
297a7e14dcfSSatish Balay   TAO_OWLQN  *lmP;
2988caf6e8cSBarry Smith   const char *owarmijo_type = TAOLINESEARCHOWARMIJO;
299a7e14dcfSSatish Balay 
300a7e14dcfSSatish Balay   PetscFunctionBegin;
301a7e14dcfSSatish Balay   tao->ops->setup          = TaoSetUp_OWLQN;
302a7e14dcfSSatish Balay   tao->ops->solve          = TaoSolve_OWLQN;
303a7e14dcfSSatish Balay   tao->ops->view           = TaoView_OWLQN;
304a7e14dcfSSatish Balay   tao->ops->setfromoptions = TaoSetFromOptions_OWLQN;
305a7e14dcfSSatish Balay   tao->ops->destroy        = TaoDestroy_OWLQN;
306a7e14dcfSSatish Balay 
3074dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&lmP));
30883c8fe1dSLisandro Dalcin   lmP->D      = NULL;
30983c8fe1dSLisandro Dalcin   lmP->M      = NULL;
31083c8fe1dSLisandro Dalcin   lmP->GV     = NULL;
31183c8fe1dSLisandro Dalcin   lmP->Xold   = NULL;
31283c8fe1dSLisandro Dalcin   lmP->Gold   = NULL;
313a7e14dcfSSatish Balay   lmP->lambda = 1.0;
314a7e14dcfSSatish Balay 
315a7e14dcfSSatish Balay   tao->data = (void *)lmP;
3166552cf8aSJason Sarich   /* Override default settings (unless already changed) */
317606f75f6SBarry Smith   PetscCall(TaoParametersInitialize(tao));
318606f75f6SBarry Smith   PetscObjectParameterSetDefault(tao, max_it, 2000);
319606f75f6SBarry Smith   PetscObjectParameterSetDefault(tao, max_funcs, 4000);
320a7e14dcfSSatish Balay 
3219566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch));
3229566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1));
3239566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetType(tao->linesearch, owarmijo_type));
3249566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchUseTaoRoutines(tao->linesearch, tao));
3259566063dSJacob Faibussowitsch   PetscCall(TaoLineSearchSetOptionsPrefix(tao->linesearch, tao->hdr.prefix));
3263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
327a7e14dcfSSatish Balay }
328