xref: /petsc/src/tao/complementarity/impls/asls/asils.c (revision e4cb33bb7dbdbae9285fba102465ca0f1dcb3977)
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
37a7e14dcfSSatish Balay    algorithm framework.
38a7e14dcfSSatish Balay 
39a7e14dcfSSatish Balay    References:
40a7e14dcfSSatish Balay      Billups, "Algorithms for Complementarity Problems and Generalized
41a7e14dcfSSatish Balay        Equations," Ph.D thesis, University of Wisconsin - Madison, 1995.
42a7e14dcfSSatish Balay      De Luca, Facchinei, Kanzow, "A Semismooth Equation Approach to the
43a7e14dcfSSatish Balay        Solution of Nonlinear Complementarity Problems," Mathematical
44a7e14dcfSSatish Balay        Programming, 75, pages 407-439, 1996.
45a7e14dcfSSatish Balay      Ferris, Kanzow, Munson, "Feasible Descent Algorithms for Mixed
46a7e14dcfSSatish Balay        Complementarity Problems," Mathematical Programming, 86,
47a7e14dcfSSatish Balay        pages 475-497, 1999.
48a7e14dcfSSatish Balay      Fischer, "A Special Newton-type Optimization Method," Optimization,
49a7e14dcfSSatish Balay        24, pages 269-284, 1992
50a7e14dcfSSatish Balay      Munson, Facchinei, Ferris, Fischer, Kanzow, "The Semismooth Algorithm
51a7e14dcfSSatish Balay        for Large Scale Complementarity Problems," Technical Report 99-06,
52a7e14dcfSSatish Balay        University of Wisconsin - Madison, 1999.
53a7e14dcfSSatish Balay */
54a7e14dcfSSatish Balay 
55a7e14dcfSSatish Balay 
56a7e14dcfSSatish Balay #undef __FUNCT__
57a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetUp_ASILS"
58441846f8SBarry Smith PetscErrorCode TaoSetUp_ASILS(Tao tao)
59a7e14dcfSSatish Balay {
60a7e14dcfSSatish Balay   TAO_SSLS       *asls = (TAO_SSLS *)tao->data;
61a7e14dcfSSatish Balay   PetscErrorCode ierr;
62a7e14dcfSSatish Balay 
63a7e14dcfSSatish Balay   PetscFunctionBegin;
64a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);
65a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr);
66a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution,&asls->ff);CHKERRQ(ierr);
67a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution,&asls->dpsi);CHKERRQ(ierr);
68a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution,&asls->da);CHKERRQ(ierr);
69a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution,&asls->db);CHKERRQ(ierr);
70a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution,&asls->t1);CHKERRQ(ierr);
71a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution,&asls->t2);CHKERRQ(ierr);
726c23d075SBarry Smith   asls->fixed = NULL;
736c23d075SBarry Smith   asls->free = NULL;
746c23d075SBarry Smith   asls->J_sub = NULL;
756c23d075SBarry Smith   asls->Jpre_sub = NULL;
766c23d075SBarry Smith   asls->w = NULL;
776c23d075SBarry Smith   asls->r1 = NULL;
786c23d075SBarry Smith   asls->r2 = NULL;
796c23d075SBarry Smith   asls->r3 = NULL;
806c23d075SBarry Smith   asls->dxfree = NULL;
81a7e14dcfSSatish Balay   PetscFunctionReturn(0);
82a7e14dcfSSatish Balay }
83a7e14dcfSSatish Balay 
84a7e14dcfSSatish Balay #undef __FUNCT__
85a7e14dcfSSatish Balay #define __FUNCT__ "Tao_ASLS_FunctionGradient"
86a7e14dcfSSatish Balay static PetscErrorCode Tao_ASLS_FunctionGradient(TaoLineSearch ls, Vec X, PetscReal *fcn,  Vec G, void *ptr)
87a7e14dcfSSatish Balay {
88441846f8SBarry Smith   Tao            tao = (Tao)ptr;
89a7e14dcfSSatish Balay   TAO_SSLS       *asls = (TAO_SSLS *)tao->data;
90a7e14dcfSSatish Balay   PetscErrorCode ierr;
91a7e14dcfSSatish Balay 
92a7e14dcfSSatish Balay   PetscFunctionBegin;
93a7e14dcfSSatish Balay   ierr = TaoComputeConstraints(tao, X, tao->constraints);CHKERRQ(ierr);
94a7e14dcfSSatish Balay   ierr = VecFischer(X,tao->constraints,tao->XL,tao->XU,asls->ff);CHKERRQ(ierr);
95a7e14dcfSSatish Balay   ierr = VecNorm(asls->ff,NORM_2,&asls->merit);CHKERRQ(ierr);
96a7e14dcfSSatish Balay   *fcn = 0.5*asls->merit*asls->merit;
97a7e14dcfSSatish Balay 
98a7e14dcfSSatish Balay   ierr = TaoComputeJacobian(tao, tao->solution, &tao->jacobian, &tao->jacobian_pre, &asls->matflag);CHKERRQ(ierr);
9947a47007SBarry Smith   ierr = D_Fischer(tao->jacobian, tao->solution, tao->constraints,tao->XL, tao->XU, asls->t1, asls->t2,asls->da, asls->db);CHKERRQ(ierr);
100a7e14dcfSSatish Balay   ierr = VecPointwiseMult(asls->t1, asls->ff, asls->db);CHKERRQ(ierr);
101a7e14dcfSSatish Balay   ierr = MatMultTranspose(tao->jacobian,asls->t1,G);CHKERRQ(ierr);
102a7e14dcfSSatish Balay   ierr = VecPointwiseMult(asls->t1, asls->ff, asls->da);CHKERRQ(ierr);
103a7e14dcfSSatish Balay   ierr = VecAXPY(G,1.0,asls->t1);CHKERRQ(ierr);
104a7e14dcfSSatish Balay   PetscFunctionReturn(0);
105a7e14dcfSSatish Balay }
106a7e14dcfSSatish Balay 
107a7e14dcfSSatish Balay #undef __FUNCT__
108a7e14dcfSSatish Balay #define __FUNCT__ "TaoDestroy_ASILS"
109441846f8SBarry Smith static PetscErrorCode TaoDestroy_ASILS(Tao tao)
110a7e14dcfSSatish Balay {
111a7e14dcfSSatish Balay   TAO_SSLS       *ssls = (TAO_SSLS *)tao->data;
112a7e14dcfSSatish Balay   PetscErrorCode ierr;
113a7e14dcfSSatish Balay 
114a7e14dcfSSatish Balay   PetscFunctionBegin;
115a7e14dcfSSatish Balay   ierr = VecDestroy(&ssls->ff);CHKERRQ(ierr);
116a7e14dcfSSatish Balay   ierr = VecDestroy(&ssls->dpsi);CHKERRQ(ierr);
117a7e14dcfSSatish Balay   ierr = VecDestroy(&ssls->da);CHKERRQ(ierr);
118a7e14dcfSSatish Balay   ierr = VecDestroy(&ssls->db);CHKERRQ(ierr);
119a7e14dcfSSatish Balay   ierr = VecDestroy(&ssls->w);CHKERRQ(ierr);
120a7e14dcfSSatish Balay   ierr = VecDestroy(&ssls->t1);CHKERRQ(ierr);
121a7e14dcfSSatish Balay   ierr = VecDestroy(&ssls->t2);CHKERRQ(ierr);
122a7e14dcfSSatish Balay   ierr = VecDestroy(&ssls->r1);CHKERRQ(ierr);
123a7e14dcfSSatish Balay   ierr = VecDestroy(&ssls->r2);CHKERRQ(ierr);
124a7e14dcfSSatish Balay   ierr = VecDestroy(&ssls->r3);CHKERRQ(ierr);
125a7e14dcfSSatish Balay   ierr = VecDestroy(&ssls->dxfree);CHKERRQ(ierr);
126a7e14dcfSSatish Balay   ierr = MatDestroy(&ssls->J_sub);CHKERRQ(ierr);
127a7e14dcfSSatish Balay   ierr = MatDestroy(&ssls->Jpre_sub);CHKERRQ(ierr);
128a7e14dcfSSatish Balay   ierr = ISDestroy(&ssls->fixed);CHKERRQ(ierr);
129a7e14dcfSSatish Balay   ierr = ISDestroy(&ssls->free);CHKERRQ(ierr);
130a7e14dcfSSatish Balay   ierr = PetscFree(tao->data);CHKERRQ(ierr);
131a7e14dcfSSatish Balay   PetscFunctionReturn(0);
132a7e14dcfSSatish Balay }
13347a47007SBarry Smith 
134a7e14dcfSSatish Balay #undef __FUNCT__
135a7e14dcfSSatish Balay #define __FUNCT__ "TaoSolve_ASILS"
136441846f8SBarry Smith static PetscErrorCode TaoSolve_ASILS(Tao tao)
137a7e14dcfSSatish Balay {
138a7e14dcfSSatish Balay   TAO_SSLS                     *asls = (TAO_SSLS *)tao->data;
139a7e14dcfSSatish Balay   PetscReal                    psi,ndpsi, normd, innerd, t=0;
140a7e14dcfSSatish Balay   PetscInt                     iter=0, nf;
141a7e14dcfSSatish Balay   PetscErrorCode               ierr;
142*e4cb33bbSBarry Smith   TaoConvergedReason           reason;
143*e4cb33bbSBarry Smith   TaoLineSearchConvergedReason ls_reason;
144a7e14dcfSSatish Balay 
145a7e14dcfSSatish Balay   PetscFunctionBegin;
146a7e14dcfSSatish Balay   /* Assume that Setup has been called!
147a7e14dcfSSatish Balay      Set the structure for the Jacobian and create a linear solver. */
148a7e14dcfSSatish Balay 
149a7e14dcfSSatish Balay   ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr);
150a7e14dcfSSatish Balay   ierr = TaoLineSearchSetObjectiveAndGradientRoutine(tao->linesearch,Tao_ASLS_FunctionGradient,tao);CHKERRQ(ierr);
151a7e14dcfSSatish Balay   ierr = TaoLineSearchSetObjectiveRoutine(tao->linesearch,Tao_SSLS_Function,tao);CHKERRQ(ierr);
152a7e14dcfSSatish Balay 
153a7e14dcfSSatish Balay   /* Calculate the function value and fischer function value at the
154a7e14dcfSSatish Balay      current iterate */
155a7e14dcfSSatish Balay   ierr = TaoLineSearchComputeObjectiveAndGradient(tao->linesearch,tao->solution,&psi,asls->dpsi);CHKERRQ(ierr);
156a7e14dcfSSatish Balay   ierr = VecNorm(asls->dpsi,NORM_2,&ndpsi);CHKERRQ(ierr);
157a7e14dcfSSatish Balay 
158a7e14dcfSSatish Balay   while (1) {
159a7e14dcfSSatish Balay     /* Check the termination criteria */
16047a47007SBarry Smith     ierr = PetscInfo3(tao,"iter %D, merit: %g, ||dpsi||: %g\n",iter, (double)asls->merit,  (double)ndpsi);CHKERRQ(ierr);
161a7e14dcfSSatish Balay     ierr = TaoMonitor(tao, iter++, asls->merit, ndpsi, 0.0, t, &reason);CHKERRQ(ierr);
162a7e14dcfSSatish Balay     if (TAO_CONTINUE_ITERATING != reason) break;
163a7e14dcfSSatish Balay 
164a7e14dcfSSatish Balay     /* We are going to solve a linear system of equations.  We need to
165a7e14dcfSSatish Balay        set the tolerances for the solve so that we maintain an asymptotic
166a7e14dcfSSatish Balay        rate of convergence that is superlinear.
167a7e14dcfSSatish Balay        Note: these tolerances are for the reduced system.  We really need
168a7e14dcfSSatish Balay        to make sure that the full system satisfies the full-space conditions.
169a7e14dcfSSatish Balay 
170a7e14dcfSSatish Balay        This rule gives superlinear asymptotic convergence
171a7e14dcfSSatish Balay        asls->atol = min(0.5, asls->merit*sqrt(asls->merit));
172a7e14dcfSSatish Balay        asls->rtol = 0.0;
173a7e14dcfSSatish Balay 
174a7e14dcfSSatish Balay        This rule gives quadratic asymptotic convergence
175a7e14dcfSSatish Balay        asls->atol = min(0.5, asls->merit*asls->merit);
176a7e14dcfSSatish Balay        asls->rtol = 0.0;
177a7e14dcfSSatish Balay 
178a7e14dcfSSatish Balay        Calculate a free and fixed set of variables.  The fixed set of
179a7e14dcfSSatish Balay        variables are those for the d_b is approximately equal to zero.
180a7e14dcfSSatish Balay        The definition of approximately changes as we approach the solution
181a7e14dcfSSatish Balay        to the problem.
182a7e14dcfSSatish Balay 
183a7e14dcfSSatish Balay        No one rule is guaranteed to work in all cases.  The following
184a7e14dcfSSatish Balay        definition is based on the norm of the Jacobian matrix.  If the
185a7e14dcfSSatish Balay        norm is large, the tolerance becomes smaller. */
186a7e14dcfSSatish Balay     ierr = MatNorm(tao->jacobian,NORM_1,&asls->identifier);CHKERRQ(ierr);
187a7e14dcfSSatish Balay     asls->identifier = PetscMin(asls->merit, 1e-2) / (1 + asls->identifier);
188a7e14dcfSSatish Balay 
189a7e14dcfSSatish Balay     ierr = VecSet(asls->t1,-asls->identifier);CHKERRQ(ierr);
190a7e14dcfSSatish Balay     ierr = VecSet(asls->t2, asls->identifier);CHKERRQ(ierr);
191a7e14dcfSSatish Balay 
192a7e14dcfSSatish Balay     ierr = ISDestroy(&asls->fixed);CHKERRQ(ierr);
193a7e14dcfSSatish Balay     ierr = ISDestroy(&asls->free);CHKERRQ(ierr);
194a7e14dcfSSatish Balay     ierr = VecWhichBetweenOrEqual(asls->t1, asls->db, asls->t2, &asls->fixed);CHKERRQ(ierr);
1954473680cSBarry Smith     ierr = ISComplementVec(asls->fixed,asls->t1, &asls->free);CHKERRQ(ierr);
196a7e14dcfSSatish Balay 
197a7e14dcfSSatish Balay     ierr = ISGetSize(asls->fixed,&nf);CHKERRQ(ierr);
198335036cbSBarry Smith     ierr = PetscInfo1(tao,"Number of fixed variables: %D\n", nf);CHKERRQ(ierr);
199a7e14dcfSSatish Balay 
200a7e14dcfSSatish Balay     /* We now have our partition.  Now calculate the direction in the
201a7e14dcfSSatish Balay        fixed variable space. */
202a7e14dcfSSatish Balay     ierr = VecGetSubVec(asls->ff, asls->fixed, tao->subset_type, 0.0, &asls->r1);
203a7e14dcfSSatish Balay     ierr = VecGetSubVec(asls->da, asls->fixed, tao->subset_type, 1.0, &asls->r2);
204a7e14dcfSSatish Balay     ierr = VecPointwiseDivide(asls->r1,asls->r1,asls->r2);CHKERRQ(ierr);
205a7e14dcfSSatish Balay     ierr = VecSet(tao->stepdirection,0.0);CHKERRQ(ierr);
2064473680cSBarry Smith     ierr = VecISAXPY(tao->stepdirection, asls->fixed,1.0,asls->r1);CHKERRQ(ierr);
207a7e14dcfSSatish Balay 
208a7e14dcfSSatish Balay     /* Our direction in the Fixed Variable Set is fixed.  Calculate the
209a7e14dcfSSatish Balay        information needed for the step in the Free Variable Set.  To
210a7e14dcfSSatish Balay        do this, we need to know the diagonal perturbation and the
211a7e14dcfSSatish Balay        right hand side. */
212a7e14dcfSSatish Balay 
213a7e14dcfSSatish Balay     ierr = VecGetSubVec(asls->da, asls->free, tao->subset_type, 0.0, &asls->r1);CHKERRQ(ierr);
214a7e14dcfSSatish Balay     ierr = VecGetSubVec(asls->ff, asls->free, tao->subset_type, 0.0, &asls->r2);CHKERRQ(ierr);
215a7e14dcfSSatish Balay     ierr = VecGetSubVec(asls->db, asls->free, tao->subset_type, 1.0, &asls->r3);CHKERRQ(ierr);
216a7e14dcfSSatish Balay     ierr = VecPointwiseDivide(asls->r1,asls->r1, asls->r3);CHKERRQ(ierr);
217a7e14dcfSSatish Balay     ierr = VecPointwiseDivide(asls->r2,asls->r2, asls->r3);CHKERRQ(ierr);
218a7e14dcfSSatish Balay 
219a7e14dcfSSatish Balay     /* r1 is the diagonal perturbation
220a7e14dcfSSatish Balay        r2 is the right hand side
221a7e14dcfSSatish Balay        r3 is no longer needed
222a7e14dcfSSatish Balay 
223a7e14dcfSSatish Balay        Now need to modify r2 for our direction choice in the fixed
224a7e14dcfSSatish Balay        variable set:  calculate t1 = J*d, take the reduced vector
225a7e14dcfSSatish Balay        of t1 and modify r2. */
226a7e14dcfSSatish Balay 
227a7e14dcfSSatish Balay     ierr = MatMult(tao->jacobian, tao->stepdirection, asls->t1);CHKERRQ(ierr);
228a7e14dcfSSatish Balay     ierr = VecGetSubVec(asls->t1,asls->free,tao->subset_type,0.0,&asls->r3);CHKERRQ(ierr);
229a7e14dcfSSatish Balay     ierr = VecAXPY(asls->r2, -1.0, asls->r3);CHKERRQ(ierr);
230a7e14dcfSSatish Balay 
231a7e14dcfSSatish Balay     /* Calculate the reduced problem matrix and the direction */
23247a47007SBarry Smith     if (!asls->w && (tao->subset_type == TAO_SUBSET_MASK || tao->subset_type == TAO_SUBSET_MATRIXFREE)) {
233a7e14dcfSSatish Balay       ierr = VecDuplicate(tao->solution, &asls->w);CHKERRQ(ierr);
234a7e14dcfSSatish Balay     }
235a7e14dcfSSatish Balay     ierr = MatGetSubMat(tao->jacobian, asls->free, asls->w, tao->subset_type,&asls->J_sub);CHKERRQ(ierr);
236a7e14dcfSSatish Balay     if (tao->jacobian != tao->jacobian_pre) {
237a7e14dcfSSatish Balay       ierr = MatGetSubMat(tao->jacobian_pre, asls->free, asls->w, tao->subset_type, &asls->Jpre_sub);CHKERRQ(ierr);
238a7e14dcfSSatish Balay     } else {
239a7e14dcfSSatish Balay       ierr = MatDestroy(&asls->Jpre_sub);CHKERRQ(ierr);
240a7e14dcfSSatish Balay       asls->Jpre_sub = asls->J_sub;
241a7e14dcfSSatish Balay       ierr = PetscObjectReference((PetscObject)(asls->Jpre_sub));CHKERRQ(ierr);
242a7e14dcfSSatish Balay     }
243a7e14dcfSSatish Balay     ierr = MatDiagonalSet(asls->J_sub, asls->r1,ADD_VALUES);CHKERRQ(ierr);
244a7e14dcfSSatish Balay     ierr = VecGetSubVec(tao->stepdirection, asls->free, tao->subset_type, 0.0, &asls->dxfree);CHKERRQ(ierr);
245a7e14dcfSSatish Balay     ierr = VecSet(asls->dxfree, 0.0);CHKERRQ(ierr);
246a7e14dcfSSatish Balay 
247a7e14dcfSSatish Balay     /* Calculate the reduced direction.  (Really negative of Newton
248a7e14dcfSSatish Balay        direction.  Therefore, rest of the code uses -d.) */
249a7e14dcfSSatish Balay     ierr = KSPReset(tao->ksp);
250a7e14dcfSSatish Balay     ierr = KSPSetOperators(tao->ksp, asls->J_sub, asls->Jpre_sub,  asls->matflag);CHKERRQ(ierr);
251a7e14dcfSSatish Balay     ierr = KSPSolve(tao->ksp, asls->r2, asls->dxfree);CHKERRQ(ierr);
252a7e14dcfSSatish Balay 
253a7e14dcfSSatish Balay     /* Add the direction in the free variables back into the real direction. */
2544473680cSBarry Smith     ierr = VecISAXPY(tao->stepdirection, asls->free, 1.0,asls->dxfree);CHKERRQ(ierr);
255a7e14dcfSSatish Balay 
256a7e14dcfSSatish Balay     /* Check the real direction for descent and if not, use the negative
257a7e14dcfSSatish Balay        gradient direction. */
258a7e14dcfSSatish Balay     ierr = VecNorm(tao->stepdirection, NORM_2, &normd);CHKERRQ(ierr);
259a7e14dcfSSatish Balay     ierr = VecDot(tao->stepdirection, asls->dpsi, &innerd);CHKERRQ(ierr);
260a7e14dcfSSatish Balay 
261a7e14dcfSSatish Balay     if (innerd <= asls->delta*pow(normd, asls->rho)) {
262335036cbSBarry Smith       ierr = PetscInfo1(tao,"Gradient direction: %5.4e.\n", (double)innerd);CHKERRQ(ierr);
263335036cbSBarry Smith       ierr = PetscInfo1(tao, "Iteration %D: newton direction not descent\n", iter);CHKERRQ(ierr);
264a7e14dcfSSatish Balay       ierr = VecCopy(asls->dpsi, tao->stepdirection);CHKERRQ(ierr);
265a7e14dcfSSatish Balay       ierr = VecDot(asls->dpsi, tao->stepdirection, &innerd);CHKERRQ(ierr);
266a7e14dcfSSatish Balay     }
267a7e14dcfSSatish Balay 
268a7e14dcfSSatish Balay     ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr);
269a7e14dcfSSatish Balay     innerd = -innerd;
270a7e14dcfSSatish Balay 
271a7e14dcfSSatish Balay     /* We now have a correct descent direction.  Apply a linesearch to
272a7e14dcfSSatish Balay        find the new iterate. */
273a7e14dcfSSatish Balay     ierr = TaoLineSearchSetInitialStepLength(tao->linesearch, 1.0);CHKERRQ(ierr);
27447a47007SBarry Smith     ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &psi,asls->dpsi, tao->stepdirection, &t, &ls_reason);CHKERRQ(ierr);
275a7e14dcfSSatish Balay     ierr = VecNorm(asls->dpsi, NORM_2, &ndpsi);CHKERRQ(ierr);
276a7e14dcfSSatish Balay   }
277a7e14dcfSSatish Balay   PetscFunctionReturn(0);
278a7e14dcfSSatish Balay }
279a7e14dcfSSatish Balay 
280a7e14dcfSSatish Balay /* ---------------------------------------------------------- */
281a7e14dcfSSatish Balay EXTERN_C_BEGIN
282a7e14dcfSSatish Balay #undef __FUNCT__
283a7e14dcfSSatish Balay #define __FUNCT__ "TaoCreate_ASILS"
284441846f8SBarry Smith PetscErrorCode TaoCreate_ASILS(Tao tao)
285a7e14dcfSSatish Balay {
286a7e14dcfSSatish Balay   TAO_SSLS       *asls;
287a7e14dcfSSatish Balay   PetscErrorCode ierr;
288a7e14dcfSSatish Balay   const char     *armijo_type = TAOLINESEARCH_ARMIJO;
289a7e14dcfSSatish Balay 
290a7e14dcfSSatish Balay   PetscFunctionBegin;
2913c9e27cfSGeoffrey Irving   ierr = PetscNewLog(tao,&asls);CHKERRQ(ierr);
292a7e14dcfSSatish Balay   tao->data = (void*)asls;
293a7e14dcfSSatish Balay   tao->ops->solve = TaoSolve_ASILS;
294a7e14dcfSSatish Balay   tao->ops->setup = TaoSetUp_ASILS;
295a7e14dcfSSatish Balay   tao->ops->view = TaoView_SSLS;
296a7e14dcfSSatish Balay   tao->ops->setfromoptions = TaoSetFromOptions_SSLS;
297a7e14dcfSSatish Balay   tao->ops->destroy = TaoDestroy_ASILS;
298a7e14dcfSSatish Balay   tao->subset_type = TAO_SUBSET_SUBVEC;
299a7e14dcfSSatish Balay   asls->delta = 1e-10;
300a7e14dcfSSatish Balay   asls->rho = 2.1;
3016c23d075SBarry Smith   asls->fixed = NULL;
3026c23d075SBarry Smith   asls->free = NULL;
3036c23d075SBarry Smith   asls->J_sub = NULL;
3046c23d075SBarry Smith   asls->Jpre_sub = NULL;
3056c23d075SBarry Smith   asls->w = NULL;
3066c23d075SBarry Smith   asls->r1 = NULL;
3076c23d075SBarry Smith   asls->r2 = NULL;
3086c23d075SBarry Smith   asls->r3 = NULL;
3096c23d075SBarry Smith   asls->t1 = NULL;
3106c23d075SBarry Smith   asls->t2 = NULL;
3116c23d075SBarry Smith   asls->dxfree = NULL;
312a7e14dcfSSatish Balay 
313a7e14dcfSSatish Balay   asls->identifier = 1e-5;
314a7e14dcfSSatish Balay 
315a7e14dcfSSatish Balay   ierr = TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch);CHKERRQ(ierr);
316a7e14dcfSSatish Balay   ierr = TaoLineSearchSetType(tao->linesearch, armijo_type);CHKERRQ(ierr);
317a7e14dcfSSatish Balay   ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr);
318a7e14dcfSSatish Balay 
319a7e14dcfSSatish Balay   ierr = KSPCreate(((PetscObject)tao)->comm, &tao->ksp);CHKERRQ(ierr);
320a7e14dcfSSatish Balay   ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
321a7e14dcfSSatish Balay   tao->max_it = 2000;
322a7e14dcfSSatish Balay   tao->max_funcs = 4000;
323a7e14dcfSSatish Balay   tao->fatol = 0;
324a7e14dcfSSatish Balay   tao->frtol = 0;
325a7e14dcfSSatish Balay   tao->gttol = 0;
326a7e14dcfSSatish Balay   tao->grtol = 0;
3276f4723b1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE)
3286f4723b1SBarry Smith   tao->gatol = 1.0e-6;
3296f4723b1SBarry Smith   tao->fmin = 1.0e-4;
3306f4723b1SBarry Smith #else
331a7e14dcfSSatish Balay   tao->gatol = 1.0e-16;
332a7e14dcfSSatish Balay   tao->fmin = 1.0e-8;
3336f4723b1SBarry Smith #endif
334a7e14dcfSSatish Balay   PetscFunctionReturn(0);
335a7e14dcfSSatish Balay }
336a7e14dcfSSatish Balay EXTERN_C_END
337a7e14dcfSSatish Balay 
338