1aaa7dc30SBarry Smith #include <../src/tao/complementarity/impls/ssls/ssls.h>
2a7e14dcfSSatish Balay /*
3a7e14dcfSSatish Balay Context for ASXLS
4a7e14dcfSSatish Balay -- active-set - reduced matrices formed
5a7e14dcfSSatish Balay - inherit properties of original system
6a7e14dcfSSatish Balay -- semismooth (S) - function not differentiable
7a7e14dcfSSatish Balay - merit function continuously differentiable
8a7e14dcfSSatish Balay - Fischer-Burmeister reformulation of complementarity
9a7e14dcfSSatish Balay - Billups composition for two finite bounds
10a7e14dcfSSatish Balay -- infeasible (I) - iterates not guaranteed to remain within bounds
11a7e14dcfSSatish Balay -- feasible (F) - iterates guaranteed to remain within bounds
12a7e14dcfSSatish Balay -- linesearch (LS) - Armijo rule on direction
13a7e14dcfSSatish Balay
14a7e14dcfSSatish Balay Many other reformulations are possible and combinations of
15a7e14dcfSSatish Balay feasible/infeasible and linesearch/trust region are possible.
16a7e14dcfSSatish Balay
17a7e14dcfSSatish Balay Basic theory
18a7e14dcfSSatish Balay Fischer-Burmeister reformulation is semismooth with a continuously
19a7e14dcfSSatish Balay differentiable merit function and strongly semismooth if the F has
20a7e14dcfSSatish Balay lipschitz continuous derivatives.
21a7e14dcfSSatish Balay
22a7e14dcfSSatish Balay Every accumulation point generated by the algorithm is a stationary
23a7e14dcfSSatish Balay point for the merit function. Stationary points of the merit function
24a7e14dcfSSatish Balay are solutions of the complementarity problem if
25a7e14dcfSSatish Balay a. the stationary point has a BD-regular subdifferential, or
26a7e14dcfSSatish Balay b. the Schur complement F'/F'_ff is a P_0-matrix where ff is the
27a7e14dcfSSatish Balay index set corresponding to the free variables.
28a7e14dcfSSatish Balay
29a7e14dcfSSatish Balay If one of the accumulation points has a BD-regular subdifferential then
30a7e14dcfSSatish Balay a. the entire sequence converges to this accumulation point at
31a7e14dcfSSatish Balay a local q-superlinear rate
32a7e14dcfSSatish Balay b. if in addition the reformulation is strongly semismooth near
33a7e14dcfSSatish Balay this accumulation point, then the algorithm converges at a
34a7e14dcfSSatish Balay local q-quadratic rate.
35a7e14dcfSSatish Balay
36a7e14dcfSSatish Balay The theory for the feasible version follows from the feasible descent
371d27aa22SBarry Smith algorithm framework. See {cite}`billups:algorithms`, {cite}`deluca.facchinei.ea:semismooth`,
381d27aa22SBarry Smith {cite}`ferris.kanzow.ea:feasible`, {cite}`fischer:special`, and {cite}`munson.facchinei.ea:semismooth`.
39a7e14dcfSSatish Balay */
40a7e14dcfSSatish Balay
TaoSetUp_ASILS(Tao tao)41d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSetUp_ASILS(Tao tao)
42d71ae5a4SJacob Faibussowitsch {
43a7e14dcfSSatish Balay TAO_SSLS *asls = (TAO_SSLS *)tao->data;
44a7e14dcfSSatish Balay
45a7e14dcfSSatish Balay PetscFunctionBegin;
469566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &tao->gradient));
479566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &tao->stepdirection));
489566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &asls->ff));
499566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &asls->dpsi));
509566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &asls->da));
519566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &asls->db));
529566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &asls->t1));
539566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &asls->t2));
546c23d075SBarry Smith asls->fixed = NULL;
556c23d075SBarry Smith asls->free = NULL;
566c23d075SBarry Smith asls->J_sub = NULL;
576c23d075SBarry Smith asls->Jpre_sub = NULL;
586c23d075SBarry Smith asls->w = NULL;
596c23d075SBarry Smith asls->r1 = NULL;
606c23d075SBarry Smith asls->r2 = NULL;
616c23d075SBarry Smith asls->r3 = NULL;
626c23d075SBarry Smith asls->dxfree = NULL;
633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
64a7e14dcfSSatish Balay }
65a7e14dcfSSatish Balay
Tao_ASLS_FunctionGradient(TaoLineSearch ls,Vec X,PetscReal * fcn,Vec G,void * ptr)66d71ae5a4SJacob Faibussowitsch static PetscErrorCode Tao_ASLS_FunctionGradient(TaoLineSearch ls, Vec X, PetscReal *fcn, Vec G, void *ptr)
67d71ae5a4SJacob Faibussowitsch {
68441846f8SBarry Smith Tao tao = (Tao)ptr;
69a7e14dcfSSatish Balay TAO_SSLS *asls = (TAO_SSLS *)tao->data;
70a7e14dcfSSatish Balay
71a7e14dcfSSatish Balay PetscFunctionBegin;
729566063dSJacob Faibussowitsch PetscCall(TaoComputeConstraints(tao, X, tao->constraints));
739566063dSJacob Faibussowitsch PetscCall(VecFischer(X, tao->constraints, tao->XL, tao->XU, asls->ff));
749566063dSJacob Faibussowitsch PetscCall(VecNorm(asls->ff, NORM_2, &asls->merit));
75a7e14dcfSSatish Balay *fcn = 0.5 * asls->merit * asls->merit;
76a7e14dcfSSatish Balay
779566063dSJacob Faibussowitsch PetscCall(TaoComputeJacobian(tao, tao->solution, tao->jacobian, tao->jacobian_pre));
789566063dSJacob Faibussowitsch PetscCall(MatDFischer(tao->jacobian, tao->solution, tao->constraints, tao->XL, tao->XU, asls->t1, asls->t2, asls->da, asls->db));
799566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(asls->t1, asls->ff, asls->db));
809566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(tao->jacobian, asls->t1, G));
819566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(asls->t1, asls->ff, asls->da));
829566063dSJacob Faibussowitsch PetscCall(VecAXPY(G, 1.0, asls->t1));
833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
84a7e14dcfSSatish Balay }
85a7e14dcfSSatish Balay
TaoDestroy_ASILS(Tao tao)86d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoDestroy_ASILS(Tao tao)
87d71ae5a4SJacob Faibussowitsch {
88a7e14dcfSSatish Balay TAO_SSLS *ssls = (TAO_SSLS *)tao->data;
89a7e14dcfSSatish Balay
90a7e14dcfSSatish Balay PetscFunctionBegin;
919566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ssls->ff));
929566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ssls->dpsi));
939566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ssls->da));
949566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ssls->db));
959566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ssls->w));
969566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ssls->t1));
979566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ssls->t2));
989566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ssls->r1));
999566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ssls->r2));
1009566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ssls->r3));
1019566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ssls->dxfree));
1029566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ssls->J_sub));
1039566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ssls->Jpre_sub));
1049566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ssls->fixed));
1059566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ssls->free));
106a958fbfcSStefano Zampini PetscCall(KSPDestroy(&tao->ksp));
1079566063dSJacob Faibussowitsch PetscCall(PetscFree(tao->data));
1083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
109a7e14dcfSSatish Balay }
11047a47007SBarry Smith
TaoSolve_ASILS(Tao tao)111d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSolve_ASILS(Tao tao)
112d71ae5a4SJacob Faibussowitsch {
113a7e14dcfSSatish Balay TAO_SSLS *asls = (TAO_SSLS *)tao->data;
114a7e14dcfSSatish Balay PetscReal psi, ndpsi, normd, innerd, t = 0;
1158931d482SJason Sarich PetscInt nf;
116e4cb33bbSBarry Smith TaoLineSearchConvergedReason ls_reason;
117a7e14dcfSSatish Balay
118a7e14dcfSSatish Balay PetscFunctionBegin;
119a7e14dcfSSatish Balay /* Assume that Setup has been called!
120a7e14dcfSSatish Balay Set the structure for the Jacobian and create a linear solver. */
121a7e14dcfSSatish Balay
1229566063dSJacob Faibussowitsch PetscCall(TaoComputeVariableBounds(tao));
1239566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetObjectiveAndGradientRoutine(tao->linesearch, Tao_ASLS_FunctionGradient, tao));
1249566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetObjectiveRoutine(tao->linesearch, Tao_SSLS_Function, tao));
125a7e14dcfSSatish Balay
126a7e14dcfSSatish Balay /* Calculate the function value and fischer function value at the
127a7e14dcfSSatish Balay current iterate */
1289566063dSJacob Faibussowitsch PetscCall(TaoLineSearchComputeObjectiveAndGradient(tao->linesearch, tao->solution, &psi, asls->dpsi));
1299566063dSJacob Faibussowitsch PetscCall(VecNorm(asls->dpsi, NORM_2, &ndpsi));
130a7e14dcfSSatish Balay
131763847b4SAlp Dener tao->reason = TAO_CONTINUE_ITERATING;
132a7e14dcfSSatish Balay while (1) {
133a7e14dcfSSatish Balay /* Check the termination criteria */
13463a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(tao, "iter %" PetscInt_FMT ", merit: %g, ||dpsi||: %g\n", tao->niter, (double)asls->merit, (double)ndpsi));
1359566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao, asls->merit, ndpsi, 0.0, tao->ksp_its));
1369566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao, tao->niter, asls->merit, ndpsi, 0.0, t));
137dbbe0bcdSBarry Smith PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
138763847b4SAlp Dener if (TAO_CONTINUE_ITERATING != tao->reason) break;
139e1e80dc8SAlp Dener
140e1e80dc8SAlp Dener /* Call general purpose update function */
141dbbe0bcdSBarry Smith PetscTryTypeMethod(tao, update, tao->niter, tao->user_update);
142e6d4cb7fSJason Sarich tao->niter++;
143a7e14dcfSSatish Balay
144a7e14dcfSSatish Balay /* We are going to solve a linear system of equations. We need to
145a7e14dcfSSatish Balay set the tolerances for the solve so that we maintain an asymptotic
146a7e14dcfSSatish Balay rate of convergence that is superlinear.
147a7e14dcfSSatish Balay Note: these tolerances are for the reduced system. We really need
148a7e14dcfSSatish Balay to make sure that the full system satisfies the full-space conditions.
149a7e14dcfSSatish Balay
150a7e14dcfSSatish Balay This rule gives superlinear asymptotic convergence
151a7e14dcfSSatish Balay asls->atol = min(0.5, asls->merit*sqrt(asls->merit));
152a7e14dcfSSatish Balay asls->rtol = 0.0;
153a7e14dcfSSatish Balay
154a7e14dcfSSatish Balay This rule gives quadratic asymptotic convergence
155a7e14dcfSSatish Balay asls->atol = min(0.5, asls->merit*asls->merit);
156a7e14dcfSSatish Balay asls->rtol = 0.0;
157a7e14dcfSSatish Balay
158a7e14dcfSSatish Balay Calculate a free and fixed set of variables. The fixed set of
159a7e14dcfSSatish Balay variables are those for the d_b is approximately equal to zero.
160a7e14dcfSSatish Balay The definition of approximately changes as we approach the solution
161a7e14dcfSSatish Balay to the problem.
162a7e14dcfSSatish Balay
163a7e14dcfSSatish Balay No one rule is guaranteed to work in all cases. The following
164a7e14dcfSSatish Balay definition is based on the norm of the Jacobian matrix. If the
165a7e14dcfSSatish Balay norm is large, the tolerance becomes smaller. */
1669566063dSJacob Faibussowitsch PetscCall(MatNorm(tao->jacobian, NORM_1, &asls->identifier));
167a7e14dcfSSatish Balay asls->identifier = PetscMin(asls->merit, 1e-2) / (1 + asls->identifier);
168a7e14dcfSSatish Balay
1699566063dSJacob Faibussowitsch PetscCall(VecSet(asls->t1, -asls->identifier));
1709566063dSJacob Faibussowitsch PetscCall(VecSet(asls->t2, asls->identifier));
171a7e14dcfSSatish Balay
1729566063dSJacob Faibussowitsch PetscCall(ISDestroy(&asls->fixed));
1739566063dSJacob Faibussowitsch PetscCall(ISDestroy(&asls->free));
1749566063dSJacob Faibussowitsch PetscCall(VecWhichBetweenOrEqual(asls->t1, asls->db, asls->t2, &asls->fixed));
1759566063dSJacob Faibussowitsch PetscCall(ISComplementVec(asls->fixed, asls->t1, &asls->free));
176a7e14dcfSSatish Balay
1779566063dSJacob Faibussowitsch PetscCall(ISGetSize(asls->fixed, &nf));
17863a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(tao, "Number of fixed variables: %" PetscInt_FMT "\n", nf));
179a7e14dcfSSatish Balay
180a7e14dcfSSatish Balay /* We now have our partition. Now calculate the direction in the
181a7e14dcfSSatish Balay fixed variable space. */
1829566063dSJacob Faibussowitsch PetscCall(TaoVecGetSubVec(asls->ff, asls->fixed, tao->subset_type, 0.0, &asls->r1));
1839566063dSJacob Faibussowitsch PetscCall(TaoVecGetSubVec(asls->da, asls->fixed, tao->subset_type, 1.0, &asls->r2));
1849566063dSJacob Faibussowitsch PetscCall(VecPointwiseDivide(asls->r1, asls->r1, asls->r2));
1859566063dSJacob Faibussowitsch PetscCall(VecSet(tao->stepdirection, 0.0));
1869566063dSJacob Faibussowitsch PetscCall(VecISAXPY(tao->stepdirection, asls->fixed, 1.0, asls->r1));
187a7e14dcfSSatish Balay
188a7e14dcfSSatish Balay /* Our direction in the Fixed Variable Set is fixed. Calculate the
189a7e14dcfSSatish Balay information needed for the step in the Free Variable Set. To
190a7e14dcfSSatish Balay do this, we need to know the diagonal perturbation and the
191dd8e379bSPierre Jolivet right-hand side. */
192a7e14dcfSSatish Balay
1939566063dSJacob Faibussowitsch PetscCall(TaoVecGetSubVec(asls->da, asls->free, tao->subset_type, 0.0, &asls->r1));
1949566063dSJacob Faibussowitsch PetscCall(TaoVecGetSubVec(asls->ff, asls->free, tao->subset_type, 0.0, &asls->r2));
1959566063dSJacob Faibussowitsch PetscCall(TaoVecGetSubVec(asls->db, asls->free, tao->subset_type, 1.0, &asls->r3));
1969566063dSJacob Faibussowitsch PetscCall(VecPointwiseDivide(asls->r1, asls->r1, asls->r3));
1979566063dSJacob Faibussowitsch PetscCall(VecPointwiseDivide(asls->r2, asls->r2, asls->r3));
198a7e14dcfSSatish Balay
199a7e14dcfSSatish Balay /* r1 is the diagonal perturbation
200dd8e379bSPierre Jolivet r2 is the right-hand side
201a7e14dcfSSatish Balay r3 is no longer needed
202a7e14dcfSSatish Balay
203a7e14dcfSSatish Balay Now need to modify r2 for our direction choice in the fixed
204a7e14dcfSSatish Balay variable set: calculate t1 = J*d, take the reduced vector
205a7e14dcfSSatish Balay of t1 and modify r2. */
206a7e14dcfSSatish Balay
2079566063dSJacob Faibussowitsch PetscCall(MatMult(tao->jacobian, tao->stepdirection, asls->t1));
2089566063dSJacob Faibussowitsch PetscCall(TaoVecGetSubVec(asls->t1, asls->free, tao->subset_type, 0.0, &asls->r3));
2099566063dSJacob Faibussowitsch PetscCall(VecAXPY(asls->r2, -1.0, asls->r3));
210a7e14dcfSSatish Balay
211a7e14dcfSSatish Balay /* Calculate the reduced problem matrix and the direction */
21248a46eb9SPierre Jolivet if (!asls->w && (tao->subset_type == TAO_SUBSET_MASK || tao->subset_type == TAO_SUBSET_MATRIXFREE)) PetscCall(VecDuplicate(tao->solution, &asls->w));
2139566063dSJacob Faibussowitsch PetscCall(TaoMatGetSubMat(tao->jacobian, asls->free, asls->w, tao->subset_type, &asls->J_sub));
214a7e14dcfSSatish Balay if (tao->jacobian != tao->jacobian_pre) {
2159566063dSJacob Faibussowitsch PetscCall(TaoMatGetSubMat(tao->jacobian_pre, asls->free, asls->w, tao->subset_type, &asls->Jpre_sub));
216a7e14dcfSSatish Balay } else {
2179566063dSJacob Faibussowitsch PetscCall(MatDestroy(&asls->Jpre_sub));
218a7e14dcfSSatish Balay asls->Jpre_sub = asls->J_sub;
219f4f49eeaSPierre Jolivet PetscCall(PetscObjectReference((PetscObject)asls->Jpre_sub));
220a7e14dcfSSatish Balay }
2219566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet(asls->J_sub, asls->r1, ADD_VALUES));
2229566063dSJacob Faibussowitsch PetscCall(TaoVecGetSubVec(tao->stepdirection, asls->free, tao->subset_type, 0.0, &asls->dxfree));
2239566063dSJacob Faibussowitsch PetscCall(VecSet(asls->dxfree, 0.0));
224a7e14dcfSSatish Balay
225a7e14dcfSSatish Balay /* Calculate the reduced direction. (Really negative of Newton
226a7e14dcfSSatish Balay direction. Therefore, rest of the code uses -d.) */
2279566063dSJacob Faibussowitsch PetscCall(KSPReset(tao->ksp));
2289566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(tao->ksp, asls->J_sub, asls->Jpre_sub));
2299566063dSJacob Faibussowitsch PetscCall(KSPSolve(tao->ksp, asls->r2, asls->dxfree));
2309566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(tao->ksp, &tao->ksp_its));
231b0026674SJason Sarich tao->ksp_tot_its += tao->ksp_its;
232a7e14dcfSSatish Balay
233a7e14dcfSSatish Balay /* Add the direction in the free variables back into the real direction. */
2349566063dSJacob Faibussowitsch PetscCall(VecISAXPY(tao->stepdirection, asls->free, 1.0, asls->dxfree));
235a7e14dcfSSatish Balay
236a7e14dcfSSatish Balay /* Check the real direction for descent and if not, use the negative
237a7e14dcfSSatish Balay gradient direction. */
2389566063dSJacob Faibussowitsch PetscCall(VecNorm(tao->stepdirection, NORM_2, &normd));
2399566063dSJacob Faibussowitsch PetscCall(VecDot(tao->stepdirection, asls->dpsi, &innerd));
240a7e14dcfSSatish Balay
2411118d4bcSLisandro Dalcin if (innerd <= asls->delta * PetscPowReal(normd, asls->rho)) {
2429566063dSJacob Faibussowitsch PetscCall(PetscInfo(tao, "Gradient direction: %5.4e.\n", (double)innerd));
24363a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(tao, "Iteration %" PetscInt_FMT ": newton direction not descent\n", tao->niter));
2449566063dSJacob Faibussowitsch PetscCall(VecCopy(asls->dpsi, tao->stepdirection));
2459566063dSJacob Faibussowitsch PetscCall(VecDot(asls->dpsi, tao->stepdirection, &innerd));
246a7e14dcfSSatish Balay }
247a7e14dcfSSatish Balay
2489566063dSJacob Faibussowitsch PetscCall(VecScale(tao->stepdirection, -1.0));
249a7e14dcfSSatish Balay innerd = -innerd;
250a7e14dcfSSatish Balay
251a7e14dcfSSatish Balay /* We now have a correct descent direction. Apply a linesearch to
252a7e14dcfSSatish Balay find the new iterate. */
2539566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetInitialStepLength(tao->linesearch, 1.0));
2549566063dSJacob Faibussowitsch PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &psi, asls->dpsi, tao->stepdirection, &t, &ls_reason));
2559566063dSJacob Faibussowitsch PetscCall(VecNorm(asls->dpsi, NORM_2, &ndpsi));
256a7e14dcfSSatish Balay }
2573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
258a7e14dcfSSatish Balay }
259a7e14dcfSSatish Balay
2601522df2eSJason Sarich /*MC
2611d27aa22SBarry Smith TAOASILS - Active-set infeasible linesearch algorithm for solving complementarity constraints
2621522df2eSJason Sarich
2631522df2eSJason Sarich Options Database Keys:
2641522df2eSJason Sarich + -tao_ssls_delta - descent test fraction
2651522df2eSJason Sarich - -tao_ssls_rho - descent test power
2661522df2eSJason Sarich
2671eb8069cSJason Sarich Level: beginner
2681d27aa22SBarry Smith
2691d27aa22SBarry Smith Note:
2701d27aa22SBarry Smith See {cite}`billups:algorithms`, {cite}`deluca.facchinei.ea:semismooth`,
2711d27aa22SBarry Smith {cite}`ferris.kanzow.ea:feasible`, {cite}`fischer:special`, and {cite}`munson.facchinei.ea:semismooth`.
2721d27aa22SBarry Smith
2731d27aa22SBarry Smith .seealso: `Tao`, `TaoType`, `TAOASFLS`
2741522df2eSJason Sarich M*/
TaoCreate_ASILS(Tao tao)275d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TaoCreate_ASILS(Tao tao)
276d71ae5a4SJacob Faibussowitsch {
277a7e14dcfSSatish Balay TAO_SSLS *asls;
2788caf6e8cSBarry Smith const char *armijo_type = TAOLINESEARCHARMIJO;
279a7e14dcfSSatish Balay
280a7e14dcfSSatish Balay PetscFunctionBegin;
2814dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&asls));
282a7e14dcfSSatish Balay tao->data = (void *)asls;
283a7e14dcfSSatish Balay tao->ops->solve = TaoSolve_ASILS;
284a7e14dcfSSatish Balay tao->ops->setup = TaoSetUp_ASILS;
285a7e14dcfSSatish Balay tao->ops->view = TaoView_SSLS;
286a7e14dcfSSatish Balay tao->ops->setfromoptions = TaoSetFromOptions_SSLS;
287a7e14dcfSSatish Balay tao->ops->destroy = TaoDestroy_ASILS;
288a7e14dcfSSatish Balay tao->subset_type = TAO_SUBSET_SUBVEC;
289a7e14dcfSSatish Balay asls->delta = 1e-10;
290a7e14dcfSSatish Balay asls->rho = 2.1;
2916c23d075SBarry Smith asls->fixed = NULL;
2926c23d075SBarry Smith asls->free = NULL;
2936c23d075SBarry Smith asls->J_sub = NULL;
2946c23d075SBarry Smith asls->Jpre_sub = NULL;
2956c23d075SBarry Smith asls->w = NULL;
2966c23d075SBarry Smith asls->r1 = NULL;
2976c23d075SBarry Smith asls->r2 = NULL;
2986c23d075SBarry Smith asls->r3 = NULL;
2996c23d075SBarry Smith asls->t1 = NULL;
3006c23d075SBarry Smith asls->t2 = NULL;
3016c23d075SBarry Smith asls->dxfree = NULL;
302a7e14dcfSSatish Balay
303a7e14dcfSSatish Balay asls->identifier = 1e-5;
304a7e14dcfSSatish Balay
305*606f75f6SBarry Smith PetscCall(TaoParametersInitialize(tao));
306*606f75f6SBarry Smith
3079566063dSJacob Faibussowitsch PetscCall(TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch));
3089566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1));
3099566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetType(tao->linesearch, armijo_type));
3109566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetOptionsPrefix(tao->linesearch, tao->hdr.prefix));
3119566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetFromOptions(tao->linesearch));
312a7e14dcfSSatish Balay
3139566063dSJacob Faibussowitsch PetscCall(KSPCreate(((PetscObject)tao)->comm, &tao->ksp));
3149566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1));
3159566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(tao->ksp, tao->hdr.prefix));
3169566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(tao->ksp));
3176552cf8aSJason Sarich
3186552cf8aSJason Sarich /* Override default settings (unless already changed) */
319*606f75f6SBarry Smith PetscObjectParameterSetDefault(tao, max_it, 2000);
320*606f75f6SBarry Smith PetscObjectParameterSetDefault(tao, max_funcs, 4000);
321*606f75f6SBarry Smith PetscObjectParameterSetDefault(tao, gttol, 0);
322*606f75f6SBarry Smith PetscObjectParameterSetDefault(tao, grtol, 0);
323*606f75f6SBarry Smith PetscObjectParameterSetDefault(tao, gatol, PetscDefined(USE_REAL_SINGLE) ? 1.0e-6 : 1.0e-16);
324*606f75f6SBarry Smith PetscObjectParameterSetDefault(tao, fmin, PetscDefined(USE_REAL_SINGLE) ? 1.0e-4 : 1.0e-8);
3253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
326a7e14dcfSSatish Balay }
327