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_ASFLS(Tao tao)41d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSetUp_ASFLS(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));
549566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &asls->w));
556c23d075SBarry Smith asls->fixed = NULL;
566c23d075SBarry Smith asls->free = NULL;
576c23d075SBarry Smith asls->J_sub = NULL;
586c23d075SBarry Smith asls->Jpre_sub = 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;
769566063dSJacob Faibussowitsch PetscCall(TaoComputeJacobian(tao, tao->solution, tao->jacobian, tao->jacobian_pre));
77a7e14dcfSSatish Balay
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_ASFLS(Tao tao)86d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoDestroy_ASFLS(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_ASFLS(Tao tao)111d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSolve_ASFLS(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));
1259566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetVariableBounds(tao->linesearch, tao->XL, tao->XU));
126a7e14dcfSSatish Balay
1279566063dSJacob Faibussowitsch PetscCall(VecMedian(tao->XL, tao->solution, tao->XU, tao->solution));
128a7e14dcfSSatish Balay
129a7e14dcfSSatish Balay /* Calculate the function value and fischer function value at the
130a7e14dcfSSatish Balay current iterate */
1319566063dSJacob Faibussowitsch PetscCall(TaoLineSearchComputeObjectiveAndGradient(tao->linesearch, tao->solution, &psi, asls->dpsi));
1329566063dSJacob Faibussowitsch PetscCall(VecNorm(asls->dpsi, NORM_2, &ndpsi));
133a7e14dcfSSatish Balay
134763847b4SAlp Dener tao->reason = TAO_CONTINUE_ITERATING;
135a7e14dcfSSatish Balay while (1) {
136e4cb33bbSBarry Smith /* Check the converged criteria */
13763a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(tao, "iter %" PetscInt_FMT ", merit: %g, ||dpsi||: %g\n", tao->niter, (double)asls->merit, (double)ndpsi));
1389566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao, asls->merit, ndpsi, 0.0, tao->ksp_its));
1399566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao, tao->niter, asls->merit, ndpsi, 0.0, t));
140dbbe0bcdSBarry Smith PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
141763847b4SAlp Dener if (TAO_CONTINUE_ITERATING != tao->reason) break;
142e1e80dc8SAlp Dener
143e1e80dc8SAlp Dener /* Call general purpose update function */
144dbbe0bcdSBarry Smith PetscTryTypeMethod(tao, update, tao->niter, tao->user_update);
145e6d4cb7fSJason Sarich tao->niter++;
146a7e14dcfSSatish Balay
147a7e14dcfSSatish Balay /* We are going to solve a linear system of equations. We need to
148a7e14dcfSSatish Balay set the tolerances for the solve so that we maintain an asymptotic
149a7e14dcfSSatish Balay rate of convergence that is superlinear.
150a7e14dcfSSatish Balay Note: these tolerances are for the reduced system. We really need
151a7e14dcfSSatish Balay to make sure that the full system satisfies the full-space conditions.
152a7e14dcfSSatish Balay
153a7e14dcfSSatish Balay This rule gives superlinear asymptotic convergence
154a7e14dcfSSatish Balay asls->atol = min(0.5, asls->merit*sqrt(asls->merit));
155a7e14dcfSSatish Balay asls->rtol = 0.0;
156a7e14dcfSSatish Balay
157a7e14dcfSSatish Balay This rule gives quadratic asymptotic convergence
158a7e14dcfSSatish Balay asls->atol = min(0.5, asls->merit*asls->merit);
159a7e14dcfSSatish Balay asls->rtol = 0.0;
160a7e14dcfSSatish Balay
161a7e14dcfSSatish Balay Calculate a free and fixed set of variables. The fixed set of
162a7e14dcfSSatish Balay variables are those for the d_b is approximately equal to zero.
163a7e14dcfSSatish Balay The definition of approximately changes as we approach the solution
164a7e14dcfSSatish Balay to the problem.
165a7e14dcfSSatish Balay
166a7e14dcfSSatish Balay No one rule is guaranteed to work in all cases. The following
167a7e14dcfSSatish Balay definition is based on the norm of the Jacobian matrix. If the
168a7e14dcfSSatish Balay norm is large, the tolerance becomes smaller. */
1699566063dSJacob Faibussowitsch PetscCall(MatNorm(tao->jacobian, NORM_1, &asls->identifier));
170a7e14dcfSSatish Balay asls->identifier = PetscMin(asls->merit, 1e-2) / (1 + asls->identifier);
171a7e14dcfSSatish Balay
1729566063dSJacob Faibussowitsch PetscCall(VecSet(asls->t1, -asls->identifier));
1739566063dSJacob Faibussowitsch PetscCall(VecSet(asls->t2, asls->identifier));
174a7e14dcfSSatish Balay
1759566063dSJacob Faibussowitsch PetscCall(ISDestroy(&asls->fixed));
1769566063dSJacob Faibussowitsch PetscCall(ISDestroy(&asls->free));
1779566063dSJacob Faibussowitsch PetscCall(VecWhichBetweenOrEqual(asls->t1, asls->db, asls->t2, &asls->fixed));
1789566063dSJacob Faibussowitsch PetscCall(ISComplementVec(asls->fixed, asls->t1, &asls->free));
179a7e14dcfSSatish Balay
1809566063dSJacob Faibussowitsch PetscCall(ISGetSize(asls->fixed, &nf));
18163a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(tao, "Number of fixed variables: %" PetscInt_FMT "\n", nf));
182a7e14dcfSSatish Balay
183a7e14dcfSSatish Balay /* We now have our partition. Now calculate the direction in the
184a7e14dcfSSatish Balay fixed variable space. */
1859566063dSJacob Faibussowitsch PetscCall(TaoVecGetSubVec(asls->ff, asls->fixed, tao->subset_type, 0.0, &asls->r1));
1869566063dSJacob Faibussowitsch PetscCall(TaoVecGetSubVec(asls->da, asls->fixed, tao->subset_type, 1.0, &asls->r2));
1879566063dSJacob Faibussowitsch PetscCall(VecPointwiseDivide(asls->r1, asls->r1, asls->r2));
1889566063dSJacob Faibussowitsch PetscCall(VecSet(tao->stepdirection, 0.0));
1899566063dSJacob Faibussowitsch PetscCall(VecISAXPY(tao->stepdirection, asls->fixed, 1.0, asls->r1));
190a7e14dcfSSatish Balay
191a7e14dcfSSatish Balay /* Our direction in the Fixed Variable Set is fixed. Calculate the
192a7e14dcfSSatish Balay information needed for the step in the Free Variable Set. To
193a7e14dcfSSatish Balay do this, we need to know the diagonal perturbation and the
194dd8e379bSPierre Jolivet right-hand side. */
195a7e14dcfSSatish Balay
1969566063dSJacob Faibussowitsch PetscCall(TaoVecGetSubVec(asls->da, asls->free, tao->subset_type, 0.0, &asls->r1));
1979566063dSJacob Faibussowitsch PetscCall(TaoVecGetSubVec(asls->ff, asls->free, tao->subset_type, 0.0, &asls->r2));
1989566063dSJacob Faibussowitsch PetscCall(TaoVecGetSubVec(asls->db, asls->free, tao->subset_type, 1.0, &asls->r3));
1999566063dSJacob Faibussowitsch PetscCall(VecPointwiseDivide(asls->r1, asls->r1, asls->r3));
2009566063dSJacob Faibussowitsch PetscCall(VecPointwiseDivide(asls->r2, asls->r2, asls->r3));
201a7e14dcfSSatish Balay
202a7e14dcfSSatish Balay /* r1 is the diagonal perturbation
203dd8e379bSPierre Jolivet r2 is the right-hand side
204a7e14dcfSSatish Balay r3 is no longer needed
205a7e14dcfSSatish Balay
206a7e14dcfSSatish Balay Now need to modify r2 for our direction choice in the fixed
207a7e14dcfSSatish Balay variable set: calculate t1 = J*d, take the reduced vector
208a7e14dcfSSatish Balay of t1 and modify r2. */
209a7e14dcfSSatish Balay
2109566063dSJacob Faibussowitsch PetscCall(MatMult(tao->jacobian, tao->stepdirection, asls->t1));
2119566063dSJacob Faibussowitsch PetscCall(TaoVecGetSubVec(asls->t1, asls->free, tao->subset_type, 0.0, &asls->r3));
2129566063dSJacob Faibussowitsch PetscCall(VecAXPY(asls->r2, -1.0, asls->r3));
213a7e14dcfSSatish Balay
214a7e14dcfSSatish Balay /* Calculate the reduced problem matrix and the direction */
2159566063dSJacob Faibussowitsch PetscCall(TaoMatGetSubMat(tao->jacobian, asls->free, asls->w, tao->subset_type, &asls->J_sub));
216a7e14dcfSSatish Balay if (tao->jacobian != tao->jacobian_pre) {
2179566063dSJacob Faibussowitsch PetscCall(TaoMatGetSubMat(tao->jacobian_pre, asls->free, asls->w, tao->subset_type, &asls->Jpre_sub));
218a7e14dcfSSatish Balay } else {
2199566063dSJacob Faibussowitsch PetscCall(MatDestroy(&asls->Jpre_sub));
220a7e14dcfSSatish Balay asls->Jpre_sub = asls->J_sub;
221f4f49eeaSPierre Jolivet PetscCall(PetscObjectReference((PetscObject)asls->Jpre_sub));
222a7e14dcfSSatish Balay }
2239566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet(asls->J_sub, asls->r1, ADD_VALUES));
2249566063dSJacob Faibussowitsch PetscCall(TaoVecGetSubVec(tao->stepdirection, asls->free, tao->subset_type, 0.0, &asls->dxfree));
2259566063dSJacob Faibussowitsch PetscCall(VecSet(asls->dxfree, 0.0));
226a7e14dcfSSatish Balay
227a7e14dcfSSatish Balay /* Calculate the reduced direction. (Really negative of Newton
228a7e14dcfSSatish Balay direction. Therefore, rest of the code uses -d.) */
2299566063dSJacob Faibussowitsch PetscCall(KSPReset(tao->ksp));
2309566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(tao->ksp, asls->J_sub, asls->Jpre_sub));
2319566063dSJacob Faibussowitsch PetscCall(KSPSolve(tao->ksp, asls->r2, asls->dxfree));
2329566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(tao->ksp, &tao->ksp_its));
233b0026674SJason Sarich tao->ksp_tot_its += tao->ksp_its;
234a7e14dcfSSatish Balay
235a7e14dcfSSatish Balay /* Add the direction in the free variables back into the real direction. */
2369566063dSJacob Faibussowitsch PetscCall(VecISAXPY(tao->stepdirection, asls->free, 1.0, asls->dxfree));
237a7e14dcfSSatish Balay
238a7e14dcfSSatish Balay /* Check the projected real direction for descent and if not, use the negative
239a7e14dcfSSatish Balay gradient direction. */
2409566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->stepdirection, asls->w));
2419566063dSJacob Faibussowitsch PetscCall(VecScale(asls->w, -1.0));
2429566063dSJacob Faibussowitsch PetscCall(VecBoundGradientProjection(asls->w, tao->solution, tao->XL, tao->XU, asls->w));
2439566063dSJacob Faibussowitsch PetscCall(VecNorm(asls->w, NORM_2, &normd));
2449566063dSJacob Faibussowitsch PetscCall(VecDot(asls->w, asls->dpsi, &innerd));
245a7e14dcfSSatish Balay
246d90ca52eSBarry Smith if (innerd >= -asls->delta * PetscPowReal(normd, asls->rho)) {
2479566063dSJacob Faibussowitsch PetscCall(PetscInfo(tao, "Gradient direction: %5.4e.\n", (double)innerd));
24863a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(tao, "Iteration %" PetscInt_FMT ": newton direction not descent\n", tao->niter));
2499566063dSJacob Faibussowitsch PetscCall(VecCopy(asls->dpsi, tao->stepdirection));
2509566063dSJacob Faibussowitsch PetscCall(VecDot(asls->dpsi, tao->stepdirection, &innerd));
251a7e14dcfSSatish Balay }
252a7e14dcfSSatish Balay
2539566063dSJacob Faibussowitsch PetscCall(VecScale(tao->stepdirection, -1.0));
254a7e14dcfSSatish Balay innerd = -innerd;
255a7e14dcfSSatish Balay
256a7e14dcfSSatish Balay /* We now have a correct descent direction. Apply a linesearch to
257a7e14dcfSSatish Balay find the new iterate. */
2589566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetInitialStepLength(tao->linesearch, 1.0));
2599566063dSJacob Faibussowitsch PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &psi, asls->dpsi, tao->stepdirection, &t, &ls_reason));
2609566063dSJacob Faibussowitsch PetscCall(VecNorm(asls->dpsi, NORM_2, &ndpsi));
261a7e14dcfSSatish Balay }
2623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
263a7e14dcfSSatish Balay }
264a7e14dcfSSatish Balay
2651522df2eSJason Sarich /*MC
2661d27aa22SBarry Smith TAOASFLS - Active-set feasible linesearch algorithm for solving complementarity constraints
2671522df2eSJason Sarich
2681522df2eSJason Sarich Options Database Keys:
2691522df2eSJason Sarich + -tao_ssls_delta - descent test fraction
2701522df2eSJason Sarich - -tao_ssls_rho - descent test power
2711522df2eSJason Sarich
2721eb8069cSJason Sarich Level: beginner
2731d27aa22SBarry Smith
2741d27aa22SBarry Smith Note:
2751d27aa22SBarry Smith See {cite}`billups:algorithms`, {cite}`deluca.facchinei.ea:semismooth`,
2761d27aa22SBarry Smith {cite}`ferris.kanzow.ea:feasible`, {cite}`fischer:special`, and {cite}`munson.facchinei.ea:semismooth`.
2771d27aa22SBarry Smith
2781d27aa22SBarry Smith .seealso: `Tao`, `TaoType`, `TAOASILS`
2791522df2eSJason Sarich M*/
TaoCreate_ASFLS(Tao tao)280d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TaoCreate_ASFLS(Tao tao)
281d71ae5a4SJacob Faibussowitsch {
282a7e14dcfSSatish Balay TAO_SSLS *asls;
2838caf6e8cSBarry Smith const char *armijo_type = TAOLINESEARCHARMIJO;
284a7e14dcfSSatish Balay
285a7e14dcfSSatish Balay PetscFunctionBegin;
2864dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&asls));
287a7e14dcfSSatish Balay tao->data = (void *)asls;
288a7e14dcfSSatish Balay tao->ops->solve = TaoSolve_ASFLS;
289a7e14dcfSSatish Balay tao->ops->setup = TaoSetUp_ASFLS;
290a7e14dcfSSatish Balay tao->ops->view = TaoView_SSLS;
291a7e14dcfSSatish Balay tao->ops->setfromoptions = TaoSetFromOptions_SSLS;
292a7e14dcfSSatish Balay tao->ops->destroy = TaoDestroy_ASFLS;
293a7e14dcfSSatish Balay tao->subset_type = TAO_SUBSET_SUBVEC;
294a7e14dcfSSatish Balay asls->delta = 1e-10;
295a7e14dcfSSatish Balay asls->rho = 2.1;
2966c23d075SBarry Smith asls->fixed = NULL;
2976c23d075SBarry Smith asls->free = NULL;
2986c23d075SBarry Smith asls->J_sub = NULL;
2996c23d075SBarry Smith asls->Jpre_sub = NULL;
3006c23d075SBarry Smith asls->w = NULL;
3016c23d075SBarry Smith asls->r1 = NULL;
3026c23d075SBarry Smith asls->r2 = NULL;
3036c23d075SBarry Smith asls->r3 = NULL;
3046c23d075SBarry Smith asls->t1 = NULL;
3056c23d075SBarry Smith asls->t2 = NULL;
3066c23d075SBarry Smith asls->dxfree = NULL;
307a7e14dcfSSatish Balay asls->identifier = 1e-5;
308a7e14dcfSSatish Balay
309*606f75f6SBarry Smith PetscCall(TaoParametersInitialize(tao));
310*606f75f6SBarry Smith
3119566063dSJacob Faibussowitsch PetscCall(TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch));
3129566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1));
3139566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetType(tao->linesearch, armijo_type));
3149566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetOptionsPrefix(tao->linesearch, tao->hdr.prefix));
3159566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetFromOptions(tao->linesearch));
316a7e14dcfSSatish Balay
3179566063dSJacob Faibussowitsch PetscCall(KSPCreate(((PetscObject)tao)->comm, &tao->ksp));
3189566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1));
3199566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(tao->ksp, tao->hdr.prefix));
3209566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(tao->ksp));
3216552cf8aSJason Sarich
3226552cf8aSJason Sarich /* Override default settings (unless already changed) */
323*606f75f6SBarry Smith PetscObjectParameterSetDefault(tao, max_it, 2000);
324*606f75f6SBarry Smith PetscObjectParameterSetDefault(tao, max_funcs, 4000);
325*606f75f6SBarry Smith PetscObjectParameterSetDefault(tao, gttol, 0);
326*606f75f6SBarry Smith PetscObjectParameterSetDefault(tao, grtol, 0);
327*606f75f6SBarry Smith PetscObjectParameterSetDefault(tao, gatol, PetscDefined(USE_REAL_SINGLE) ? 1.0e-6 : 1.0e-16);
328*606f75f6SBarry Smith PetscObjectParameterSetDefault(tao, fmin, PetscDefined(USE_REAL_SINGLE) ? 1.0e-4 : 1.0e-8);
3293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
330a7e14dcfSSatish Balay }
331