xref: /petsc/src/tao/complementarity/impls/asls/asfls.c (revision a7e14dcfba0d07adf6226a919460249440ec94c7)
1*a7e14dcfSSatish Balay #include "src/complementarity/impls/ssls/ssls.h"
2*a7e14dcfSSatish Balay /*
3*a7e14dcfSSatish Balay    Context for ASXLS
4*a7e14dcfSSatish Balay      -- active-set	- reduced matrices formed
5*a7e14dcfSSatish Balay   			  - inherit properties of original system
6*a7e14dcfSSatish Balay      -- semismooth (S)  - function not differentiable
7*a7e14dcfSSatish Balay                         - merit function continuously differentiable
8*a7e14dcfSSatish Balay                         - Fischer-Burmeister reformulation of complementarity
9*a7e14dcfSSatish Balay                           - Billups composition for two finite bounds
10*a7e14dcfSSatish Balay      -- infeasible (I)  - iterates not guaranteed to remain within bounds
11*a7e14dcfSSatish Balay      -- feasible (F)    - iterates guaranteed to remain within bounds
12*a7e14dcfSSatish Balay      -- linesearch (LS) - Armijo rule on direction
13*a7e14dcfSSatish Balay 
14*a7e14dcfSSatish Balay    Many other reformulations are possible and combinations of
15*a7e14dcfSSatish Balay    feasible/infeasible and linesearch/trust region are possible.
16*a7e14dcfSSatish Balay 
17*a7e14dcfSSatish Balay    Basic theory
18*a7e14dcfSSatish Balay      Fischer-Burmeister reformulation is semismooth with a continuously
19*a7e14dcfSSatish Balay      differentiable merit function and strongly semismooth if the F has
20*a7e14dcfSSatish Balay      lipschitz continuous derivatives.
21*a7e14dcfSSatish Balay 
22*a7e14dcfSSatish Balay      Every accumulation point generated by the algorithm is a stationary
23*a7e14dcfSSatish Balay      point for the merit function.  Stationary points of the merit function
24*a7e14dcfSSatish Balay      are solutions of the complementarity problem if
25*a7e14dcfSSatish Balay        a.  the stationary point has a BD-regular subdifferential, or
26*a7e14dcfSSatish Balay        b.  the Schur complement F'/F'_ff is a P_0-matrix where ff is the
27*a7e14dcfSSatish Balay            index set corresponding to the free variables.
28*a7e14dcfSSatish Balay 
29*a7e14dcfSSatish Balay      If one of the accumulation points has a BD-regular subdifferential then
30*a7e14dcfSSatish Balay        a.  the entire sequence converges to this accumulation point at
31*a7e14dcfSSatish Balay            a local q-superlinear rate
32*a7e14dcfSSatish Balay        b.  if in addition the reformulation is strongly semismooth near
33*a7e14dcfSSatish Balay            this accumulation point, then the algorithm converges at a
34*a7e14dcfSSatish Balay            local q-quadratic rate.
35*a7e14dcfSSatish Balay 
36*a7e14dcfSSatish Balay    The theory for the feasible version follows from the feasible descent
37*a7e14dcfSSatish Balay    algorithm framework.
38*a7e14dcfSSatish Balay 
39*a7e14dcfSSatish Balay    References:
40*a7e14dcfSSatish Balay      Billups, "Algorithms for Complementarity Problems and Generalized
41*a7e14dcfSSatish Balay        Equations," Ph.D thesis, University of Wisconsin - Madison, 1995.
42*a7e14dcfSSatish Balay      De Luca, Facchinei, Kanzow, "A Semismooth Equation Approach to the
43*a7e14dcfSSatish Balay        Solution of Nonlinear Complementarity Problems," Mathematical
44*a7e14dcfSSatish Balay        Programming, 75, pages 407-439, 1996.
45*a7e14dcfSSatish Balay      Ferris, Kanzow, Munson, "Feasible Descent Algorithms for Mixed
46*a7e14dcfSSatish Balay        Complementarity Problems," Mathematical Programming, 86,
47*a7e14dcfSSatish Balay        pages 475-497, 1999.
48*a7e14dcfSSatish Balay      Fischer, "A Special Newton-type Optimization Method," Optimization,
49*a7e14dcfSSatish Balay        24, pages 269-284, 1992
50*a7e14dcfSSatish Balay      Munson, Facchinei, Ferris, Fischer, Kanzow, "The Semismooth Algorithm
51*a7e14dcfSSatish Balay        for Large Scale Complementarity Problems," Technical Report 99-06,
52*a7e14dcfSSatish Balay        University of Wisconsin - Madison, 1999.
53*a7e14dcfSSatish Balay */
54*a7e14dcfSSatish Balay 
55*a7e14dcfSSatish Balay 
56*a7e14dcfSSatish Balay #undef __FUNCT__
57*a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetUp_ASFLS"
58*a7e14dcfSSatish Balay PetscErrorCode TaoSetUp_ASFLS(TaoSolver tao)
59*a7e14dcfSSatish Balay {
60*a7e14dcfSSatish Balay   TAO_SSLS *asls = (TAO_SSLS *)tao->data;
61*a7e14dcfSSatish Balay   PetscErrorCode ierr;
62*a7e14dcfSSatish Balay 
63*a7e14dcfSSatish Balay   PetscFunctionBegin;
64*a7e14dcfSSatish Balay 
65*a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution,&tao->gradient); CHKERRQ(ierr);
66*a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution,&tao->stepdirection); CHKERRQ(ierr);
67*a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution,&asls->ff); CHKERRQ(ierr);
68*a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution,&asls->dpsi); CHKERRQ(ierr);
69*a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution,&asls->da); CHKERRQ(ierr);
70*a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution,&asls->db); CHKERRQ(ierr);
71*a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution,&asls->t1); CHKERRQ(ierr);
72*a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution,&asls->t2); CHKERRQ(ierr);
73*a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution, &asls->w); CHKERRQ(ierr);
74*a7e14dcfSSatish Balay   asls->fixed = PETSC_NULL;
75*a7e14dcfSSatish Balay   asls->free = PETSC_NULL;
76*a7e14dcfSSatish Balay   asls->J_sub = PETSC_NULL;
77*a7e14dcfSSatish Balay   asls->Jpre_sub = PETSC_NULL;
78*a7e14dcfSSatish Balay   asls->r1 = PETSC_NULL;
79*a7e14dcfSSatish Balay   asls->r2 = PETSC_NULL;
80*a7e14dcfSSatish Balay   asls->r3 = PETSC_NULL;
81*a7e14dcfSSatish Balay   asls->dxfree = PETSC_NULL;
82*a7e14dcfSSatish Balay 
83*a7e14dcfSSatish Balay   PetscFunctionReturn(0);
84*a7e14dcfSSatish Balay }
85*a7e14dcfSSatish Balay 
86*a7e14dcfSSatish Balay #undef __FUNCT__
87*a7e14dcfSSatish Balay #define __FUNCT__ "Tao_ASLS_FunctionGradient"
88*a7e14dcfSSatish Balay static PetscErrorCode Tao_ASLS_FunctionGradient(TaoLineSearch ls, Vec X, PetscReal *fcn,  Vec G, void *ptr)
89*a7e14dcfSSatish Balay {
90*a7e14dcfSSatish Balay   TaoSolver tao = (TaoSolver)ptr;
91*a7e14dcfSSatish Balay   TAO_SSLS *asls = (TAO_SSLS *)tao->data;
92*a7e14dcfSSatish Balay   PetscErrorCode ierr;
93*a7e14dcfSSatish Balay 
94*a7e14dcfSSatish Balay   PetscFunctionBegin;
95*a7e14dcfSSatish Balay 
96*a7e14dcfSSatish Balay   ierr = TaoComputeConstraints(tao, X, tao->constraints); CHKERRQ(ierr);
97*a7e14dcfSSatish Balay   ierr = VecFischer(X,tao->constraints,tao->XL,tao->XU,asls->ff); CHKERRQ(ierr);
98*a7e14dcfSSatish Balay   ierr = VecNorm(asls->ff,NORM_2,&asls->merit); CHKERRQ(ierr);
99*a7e14dcfSSatish Balay   *fcn = 0.5*asls->merit*asls->merit;
100*a7e14dcfSSatish Balay 
101*a7e14dcfSSatish Balay   ierr = TaoComputeJacobian(tao, tao->solution, &tao->jacobian, &tao->jacobian_pre, &asls->matflag); CHKERRQ(ierr);
102*a7e14dcfSSatish Balay 
103*a7e14dcfSSatish Balay   ierr = D_Fischer(tao->jacobian, tao->solution, tao->constraints,
104*a7e14dcfSSatish Balay 		   tao->XL, tao->XU, asls->t1, asls->t2,
105*a7e14dcfSSatish Balay 		   asls->da, asls->db); CHKERRQ(ierr);
106*a7e14dcfSSatish Balay   ierr = VecPointwiseMult(asls->t1, asls->ff, asls->db); CHKERRQ(ierr);
107*a7e14dcfSSatish Balay   ierr = MatMultTranspose(tao->jacobian,asls->t1,G); CHKERRQ(ierr);
108*a7e14dcfSSatish Balay   ierr = VecPointwiseMult(asls->t1, asls->ff, asls->da); CHKERRQ(ierr);
109*a7e14dcfSSatish Balay   ierr = VecAXPY(G,1.0,asls->t1); CHKERRQ(ierr);
110*a7e14dcfSSatish Balay   PetscFunctionReturn(0);
111*a7e14dcfSSatish Balay }
112*a7e14dcfSSatish Balay 
113*a7e14dcfSSatish Balay #undef __FUNCT__
114*a7e14dcfSSatish Balay #define __FUNCT__ "TaoDestroy_ASFLS"
115*a7e14dcfSSatish Balay static PetscErrorCode TaoDestroy_ASFLS(TaoSolver tao)
116*a7e14dcfSSatish Balay {
117*a7e14dcfSSatish Balay   TAO_SSLS *ssls = (TAO_SSLS *)tao->data;
118*a7e14dcfSSatish Balay   PetscErrorCode ierr;
119*a7e14dcfSSatish Balay 
120*a7e14dcfSSatish Balay   PetscFunctionBegin;
121*a7e14dcfSSatish Balay 
122*a7e14dcfSSatish Balay   ierr = VecDestroy(&ssls->ff); CHKERRQ(ierr);
123*a7e14dcfSSatish Balay   ierr = VecDestroy(&ssls->dpsi); CHKERRQ(ierr);
124*a7e14dcfSSatish Balay   ierr = VecDestroy(&ssls->da); CHKERRQ(ierr);
125*a7e14dcfSSatish Balay   ierr = VecDestroy(&ssls->db); CHKERRQ(ierr);
126*a7e14dcfSSatish Balay   ierr = VecDestroy(&ssls->w); CHKERRQ(ierr);
127*a7e14dcfSSatish Balay   ierr = VecDestroy(&ssls->t1); CHKERRQ(ierr);
128*a7e14dcfSSatish Balay   ierr = VecDestroy(&ssls->t2); CHKERRQ(ierr);
129*a7e14dcfSSatish Balay   ierr = VecDestroy(&ssls->r1); CHKERRQ(ierr);
130*a7e14dcfSSatish Balay   ierr = VecDestroy(&ssls->r2); CHKERRQ(ierr);
131*a7e14dcfSSatish Balay   ierr = VecDestroy(&ssls->r3); CHKERRQ(ierr);
132*a7e14dcfSSatish Balay   ierr = VecDestroy(&ssls->dxfree); CHKERRQ(ierr);
133*a7e14dcfSSatish Balay   ierr = MatDestroy(&ssls->J_sub); CHKERRQ(ierr);
134*a7e14dcfSSatish Balay   ierr = MatDestroy(&ssls->Jpre_sub); CHKERRQ(ierr);
135*a7e14dcfSSatish Balay   ierr = ISDestroy(&ssls->fixed); CHKERRQ(ierr);
136*a7e14dcfSSatish Balay   ierr = ISDestroy(&ssls->free); CHKERRQ(ierr);
137*a7e14dcfSSatish Balay   ierr = PetscFree(tao->data); CHKERRQ(ierr);
138*a7e14dcfSSatish Balay 
139*a7e14dcfSSatish Balay   tao->data = PETSC_NULL;
140*a7e14dcfSSatish Balay   PetscFunctionReturn(0);
141*a7e14dcfSSatish Balay 
142*a7e14dcfSSatish Balay }
143*a7e14dcfSSatish Balay #undef __FUNCT__
144*a7e14dcfSSatish Balay #define __FUNCT__ "TaoSolve_ASFLS"
145*a7e14dcfSSatish Balay static PetscErrorCode TaoSolve_ASFLS(TaoSolver tao)
146*a7e14dcfSSatish Balay {
147*a7e14dcfSSatish Balay   TAO_SSLS *asls = (TAO_SSLS *)tao->data;
148*a7e14dcfSSatish Balay   PetscReal psi,ndpsi, normd, innerd, t=0;
149*a7e14dcfSSatish Balay   PetscInt iter=0, nf;
150*a7e14dcfSSatish Balay   PetscErrorCode ierr;
151*a7e14dcfSSatish Balay   TaoSolverTerminationReason reason;
152*a7e14dcfSSatish Balay   TaoLineSearchTerminationReason ls_reason;
153*a7e14dcfSSatish Balay 
154*a7e14dcfSSatish Balay   PetscFunctionBegin;
155*a7e14dcfSSatish Balay 
156*a7e14dcfSSatish Balay   /* Assume that Setup has been called!
157*a7e14dcfSSatish Balay      Set the structure for the Jacobian and create a linear solver. */
158*a7e14dcfSSatish Balay 
159*a7e14dcfSSatish Balay   ierr = TaoComputeVariableBounds(tao); CHKERRQ(ierr);
160*a7e14dcfSSatish Balay   ierr = TaoLineSearchSetObjectiveAndGradientRoutine(tao->linesearch,Tao_ASLS_FunctionGradient,tao); CHKERRQ(ierr);
161*a7e14dcfSSatish Balay   ierr = TaoLineSearchSetObjectiveRoutine(tao->linesearch,Tao_SSLS_Function,tao); CHKERRQ(ierr);
162*a7e14dcfSSatish Balay   ierr = TaoLineSearchSetVariableBounds(tao->linesearch,tao->XL,tao->XU); CHKERRQ(ierr);
163*a7e14dcfSSatish Balay 
164*a7e14dcfSSatish Balay   ierr = VecMedian(tao->XL, tao->solution, tao->XU, tao->solution); CHKERRQ(ierr);
165*a7e14dcfSSatish Balay 
166*a7e14dcfSSatish Balay   /* Calculate the function value and fischer function value at the
167*a7e14dcfSSatish Balay      current iterate */
168*a7e14dcfSSatish Balay   ierr = TaoLineSearchComputeObjectiveAndGradient(tao->linesearch,tao->solution,&psi,asls->dpsi); CHKERRQ(ierr);
169*a7e14dcfSSatish Balay   ierr = VecNorm(asls->dpsi,NORM_2,&ndpsi); CHKERRQ(ierr);
170*a7e14dcfSSatish Balay 
171*a7e14dcfSSatish Balay   while (1) {
172*a7e14dcfSSatish Balay 
173*a7e14dcfSSatish Balay     /* Check the termination criteria */
174*a7e14dcfSSatish Balay     ierr = PetscInfo3(tao,"iter %D, merit: %G, ||dpsi||: %G\n",iter, asls->merit,  ndpsi); CHKERRQ(ierr);
175*a7e14dcfSSatish Balay     ierr = TaoMonitor(tao, iter++, asls->merit, ndpsi, 0.0, t, &reason); CHKERRQ(ierr);
176*a7e14dcfSSatish Balay     if (TAO_CONTINUE_ITERATING != reason) break;
177*a7e14dcfSSatish Balay 
178*a7e14dcfSSatish Balay     /* We are going to solve a linear system of equations.  We need to
179*a7e14dcfSSatish Balay        set the tolerances for the solve so that we maintain an asymptotic
180*a7e14dcfSSatish Balay        rate of convergence that is superlinear.
181*a7e14dcfSSatish Balay        Note: these tolerances are for the reduced system.  We really need
182*a7e14dcfSSatish Balay        to make sure that the full system satisfies the full-space conditions.
183*a7e14dcfSSatish Balay 
184*a7e14dcfSSatish Balay        This rule gives superlinear asymptotic convergence
185*a7e14dcfSSatish Balay        asls->atol = min(0.5, asls->merit*sqrt(asls->merit));
186*a7e14dcfSSatish Balay        asls->rtol = 0.0;
187*a7e14dcfSSatish Balay 
188*a7e14dcfSSatish Balay        This rule gives quadratic asymptotic convergence
189*a7e14dcfSSatish Balay        asls->atol = min(0.5, asls->merit*asls->merit);
190*a7e14dcfSSatish Balay        asls->rtol = 0.0;
191*a7e14dcfSSatish Balay 
192*a7e14dcfSSatish Balay        Calculate a free and fixed set of variables.  The fixed set of
193*a7e14dcfSSatish Balay        variables are those for the d_b is approximately equal to zero.
194*a7e14dcfSSatish Balay        The definition of approximately changes as we approach the solution
195*a7e14dcfSSatish Balay        to the problem.
196*a7e14dcfSSatish Balay 
197*a7e14dcfSSatish Balay        No one rule is guaranteed to work in all cases.  The following
198*a7e14dcfSSatish Balay        definition is based on the norm of the Jacobian matrix.  If the
199*a7e14dcfSSatish Balay        norm is large, the tolerance becomes smaller. */
200*a7e14dcfSSatish Balay     ierr = MatNorm(tao->jacobian,NORM_1,&asls->identifier); CHKERRQ(ierr);
201*a7e14dcfSSatish Balay     asls->identifier = PetscMin(asls->merit, 1e-2) / (1 + asls->identifier);
202*a7e14dcfSSatish Balay 
203*a7e14dcfSSatish Balay     ierr = VecSet(asls->t1,-asls->identifier); CHKERRQ(ierr);
204*a7e14dcfSSatish Balay     ierr = VecSet(asls->t2, asls->identifier); CHKERRQ(ierr);
205*a7e14dcfSSatish Balay 
206*a7e14dcfSSatish Balay     ierr = ISDestroy(&asls->fixed); CHKERRQ(ierr);
207*a7e14dcfSSatish Balay     ierr = ISDestroy(&asls->free); CHKERRQ(ierr);
208*a7e14dcfSSatish Balay     ierr = VecWhichBetweenOrEqual(asls->t1, asls->db, asls->t2, &asls->fixed); CHKERRQ(ierr);
209*a7e14dcfSSatish Balay     ierr = ISCreateComplement(asls->fixed,asls->t1, &asls->free); CHKERRQ(ierr);
210*a7e14dcfSSatish Balay 
211*a7e14dcfSSatish Balay     ierr = ISGetSize(asls->fixed,&nf); CHKERRQ(ierr);
212*a7e14dcfSSatish Balay     ierr = PetscInfo1(tao,"Number of fixed variables: %d\n", nf); CHKERRQ(ierr);
213*a7e14dcfSSatish Balay 
214*a7e14dcfSSatish Balay     /* We now have our partition.  Now calculate the direction in the
215*a7e14dcfSSatish Balay        fixed variable space. */
216*a7e14dcfSSatish Balay     ierr = VecGetSubVec(asls->ff, asls->fixed, tao->subset_type, 0.0, &asls->r1);
217*a7e14dcfSSatish Balay     ierr = VecGetSubVec(asls->da, asls->fixed, tao->subset_type, 1.0, &asls->r2);
218*a7e14dcfSSatish Balay     ierr = VecPointwiseDivide(asls->r1,asls->r1,asls->r2); CHKERRQ(ierr);
219*a7e14dcfSSatish Balay     ierr = VecSet(tao->stepdirection,0.0); CHKERRQ(ierr);
220*a7e14dcfSSatish Balay     ierr = VecReducedXPY(tao->stepdirection,asls->r1, asls->fixed); CHKERRQ(ierr);
221*a7e14dcfSSatish Balay 
222*a7e14dcfSSatish Balay 
223*a7e14dcfSSatish Balay     /* Our direction in the Fixed Variable Set is fixed.  Calculate the
224*a7e14dcfSSatish Balay        information needed for the step in the Free Variable Set.  To
225*a7e14dcfSSatish Balay        do this, we need to know the diagonal perturbation and the
226*a7e14dcfSSatish Balay        right hand side. */
227*a7e14dcfSSatish Balay 
228*a7e14dcfSSatish Balay     ierr = VecGetSubVec(asls->da, asls->free, tao->subset_type, 0.0, &asls->r1); CHKERRQ(ierr);
229*a7e14dcfSSatish Balay     ierr = VecGetSubVec(asls->ff, asls->free, tao->subset_type, 0.0, &asls->r2); CHKERRQ(ierr);
230*a7e14dcfSSatish Balay     ierr = VecGetSubVec(asls->db, asls->free, tao->subset_type, 1.0, &asls->r3); CHKERRQ(ierr);
231*a7e14dcfSSatish Balay     ierr = VecPointwiseDivide(asls->r1,asls->r1, asls->r3); CHKERRQ(ierr);
232*a7e14dcfSSatish Balay     ierr = VecPointwiseDivide(asls->r2,asls->r2, asls->r3); CHKERRQ(ierr);
233*a7e14dcfSSatish Balay 
234*a7e14dcfSSatish Balay     /* r1 is the diagonal perturbation
235*a7e14dcfSSatish Balay        r2 is the right hand side
236*a7e14dcfSSatish Balay        r3 is no longer needed
237*a7e14dcfSSatish Balay 
238*a7e14dcfSSatish Balay        Now need to modify r2 for our direction choice in the fixed
239*a7e14dcfSSatish Balay        variable set:  calculate t1 = J*d, take the reduced vector
240*a7e14dcfSSatish Balay        of t1 and modify r2. */
241*a7e14dcfSSatish Balay 
242*a7e14dcfSSatish Balay     ierr = MatMult(tao->jacobian, tao->stepdirection, asls->t1); CHKERRQ(ierr);
243*a7e14dcfSSatish Balay     ierr = VecGetSubVec(asls->t1,asls->free,tao->subset_type,0.0,&asls->r3); CHKERRQ(ierr);
244*a7e14dcfSSatish Balay     ierr = VecAXPY(asls->r2, -1.0, asls->r3); CHKERRQ(ierr);
245*a7e14dcfSSatish Balay 
246*a7e14dcfSSatish Balay     /* Calculate the reduced problem matrix and the direction */
247*a7e14dcfSSatish Balay     ierr = MatGetSubMat(tao->jacobian, asls->free, asls->w, tao->subset_type,&asls->J_sub); CHKERRQ(ierr);
248*a7e14dcfSSatish Balay     if (tao->jacobian != tao->jacobian_pre) {
249*a7e14dcfSSatish Balay       ierr = MatGetSubMat(tao->jacobian_pre, asls->free, asls->w, tao->subset_type, &asls->Jpre_sub); CHKERRQ(ierr);
250*a7e14dcfSSatish Balay     } else {
251*a7e14dcfSSatish Balay       ierr = MatDestroy(&asls->Jpre_sub); CHKERRQ(ierr);
252*a7e14dcfSSatish Balay       asls->Jpre_sub = asls->J_sub;
253*a7e14dcfSSatish Balay       ierr = PetscObjectReference((PetscObject)(asls->Jpre_sub)); CHKERRQ(ierr);
254*a7e14dcfSSatish Balay     }
255*a7e14dcfSSatish Balay     ierr = MatDiagonalSet(asls->J_sub, asls->r1,ADD_VALUES); CHKERRQ(ierr);
256*a7e14dcfSSatish Balay     ierr = VecGetSubVec(tao->stepdirection, asls->free, tao->subset_type, 0.0, &asls->dxfree); CHKERRQ(ierr);
257*a7e14dcfSSatish Balay     ierr = VecSet(asls->dxfree, 0.0); CHKERRQ(ierr);
258*a7e14dcfSSatish Balay 
259*a7e14dcfSSatish Balay     /* Calculate the reduced direction.  (Really negative of Newton
260*a7e14dcfSSatish Balay        direction.  Therefore, rest of the code uses -d.) */
261*a7e14dcfSSatish Balay     ierr = KSPReset(tao->ksp); CHKERRQ(ierr);
262*a7e14dcfSSatish Balay     ierr = KSPSetOperators(tao->ksp, asls->J_sub, asls->Jpre_sub,  asls->matflag); CHKERRQ(ierr);
263*a7e14dcfSSatish Balay     ierr = KSPSolve(tao->ksp, asls->r2, asls->dxfree); CHKERRQ(ierr);
264*a7e14dcfSSatish Balay 
265*a7e14dcfSSatish Balay     /* Add the direction in the free variables back into the real direction. */
266*a7e14dcfSSatish Balay     ierr = VecReducedXPY(tao->stepdirection, asls->dxfree, asls->free); CHKERRQ(ierr);
267*a7e14dcfSSatish Balay 
268*a7e14dcfSSatish Balay 
269*a7e14dcfSSatish Balay     /* Check the projected real direction for descent and if not, use the negative
270*a7e14dcfSSatish Balay        gradient direction. */
271*a7e14dcfSSatish Balay     ierr = VecCopy(tao->stepdirection, asls->w); CHKERRQ(ierr);
272*a7e14dcfSSatish Balay     ierr = VecScale(asls->w, -1.0); CHKERRQ(ierr);
273*a7e14dcfSSatish Balay     ierr = VecBoundGradientProjection(asls->w, tao->solution, tao->XL, tao->XU, asls->w); CHKERRQ(ierr);
274*a7e14dcfSSatish Balay     ierr = VecNorm(asls->w, NORM_2, &normd); CHKERRQ(ierr);
275*a7e14dcfSSatish Balay     ierr = VecDot(asls->w, asls->dpsi, &innerd); CHKERRQ(ierr);
276*a7e14dcfSSatish Balay 
277*a7e14dcfSSatish Balay     if (innerd >= -asls->delta*pow(normd, asls->rho)) {
278*a7e14dcfSSatish Balay       ierr = PetscInfo1(tao,"Gradient direction: %5.4e.\n", innerd); CHKERRQ(ierr);
279*a7e14dcfSSatish Balay       ierr = PetscInfo1(tao, "Iteration %d: newton direction not descent\n", iter); CHKERRQ(ierr);
280*a7e14dcfSSatish Balay       ierr = VecCopy(asls->dpsi, tao->stepdirection); CHKERRQ(ierr);
281*a7e14dcfSSatish Balay       ierr = VecDot(asls->dpsi, tao->stepdirection, &innerd); CHKERRQ(ierr);
282*a7e14dcfSSatish Balay     }
283*a7e14dcfSSatish Balay 
284*a7e14dcfSSatish Balay     ierr = VecScale(tao->stepdirection, -1.0); CHKERRQ(ierr);
285*a7e14dcfSSatish Balay     innerd = -innerd;
286*a7e14dcfSSatish Balay 
287*a7e14dcfSSatish Balay     /* We now have a correct descent direction.  Apply a linesearch to
288*a7e14dcfSSatish Balay        find the new iterate. */
289*a7e14dcfSSatish Balay     ierr = TaoLineSearchSetInitialStepLength(tao->linesearch, 1.0); CHKERRQ(ierr);
290*a7e14dcfSSatish Balay     ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &psi,
291*a7e14dcfSSatish Balay 		      asls->dpsi, tao->stepdirection, &t, &ls_reason); CHKERRQ(ierr);
292*a7e14dcfSSatish Balay     ierr = VecNorm(asls->dpsi, NORM_2, &ndpsi); CHKERRQ(ierr);
293*a7e14dcfSSatish Balay   }
294*a7e14dcfSSatish Balay 
295*a7e14dcfSSatish Balay   PetscFunctionReturn(0);
296*a7e14dcfSSatish Balay }
297*a7e14dcfSSatish Balay 
298*a7e14dcfSSatish Balay /* ---------------------------------------------------------- */
299*a7e14dcfSSatish Balay EXTERN_C_BEGIN
300*a7e14dcfSSatish Balay #undef __FUNCT__
301*a7e14dcfSSatish Balay #define __FUNCT__ "TaoCreate_ASFLS"
302*a7e14dcfSSatish Balay PetscErrorCode TaoCreate_ASFLS(TaoSolver tao)
303*a7e14dcfSSatish Balay {
304*a7e14dcfSSatish Balay   TAO_SSLS *asls;
305*a7e14dcfSSatish Balay   PetscErrorCode  ierr;
306*a7e14dcfSSatish Balay   const char *armijo_type = TAOLINESEARCH_ARMIJO;
307*a7e14dcfSSatish Balay 
308*a7e14dcfSSatish Balay   PetscFunctionBegin;
309*a7e14dcfSSatish Balay   ierr = PetscNewLog(tao,TAO_SSLS,&asls); CHKERRQ(ierr);
310*a7e14dcfSSatish Balay   tao->data = (void*)asls;
311*a7e14dcfSSatish Balay   tao->ops->solve = TaoSolve_ASFLS;
312*a7e14dcfSSatish Balay   tao->ops->setup = TaoSetUp_ASFLS;
313*a7e14dcfSSatish Balay   tao->ops->view = TaoView_SSLS;
314*a7e14dcfSSatish Balay   tao->ops->setfromoptions = TaoSetFromOptions_SSLS;
315*a7e14dcfSSatish Balay   tao->ops->destroy = TaoDestroy_ASFLS;
316*a7e14dcfSSatish Balay   tao->subset_type = TAO_SUBSET_SUBVEC;
317*a7e14dcfSSatish Balay   asls->delta = 1e-10;
318*a7e14dcfSSatish Balay   asls->rho = 2.1;
319*a7e14dcfSSatish Balay   asls->fixed = PETSC_NULL;
320*a7e14dcfSSatish Balay   asls->free = PETSC_NULL;
321*a7e14dcfSSatish Balay   asls->J_sub = PETSC_NULL;
322*a7e14dcfSSatish Balay   asls->Jpre_sub = PETSC_NULL;
323*a7e14dcfSSatish Balay   asls->w = PETSC_NULL;
324*a7e14dcfSSatish Balay   asls->r1 = PETSC_NULL;
325*a7e14dcfSSatish Balay   asls->r2 = PETSC_NULL;
326*a7e14dcfSSatish Balay   asls->r3 = PETSC_NULL;
327*a7e14dcfSSatish Balay   asls->t1 = PETSC_NULL;
328*a7e14dcfSSatish Balay   asls->t2 = PETSC_NULL;
329*a7e14dcfSSatish Balay   asls->dxfree = PETSC_NULL;
330*a7e14dcfSSatish Balay 
331*a7e14dcfSSatish Balay   asls->identifier = 1e-5;
332*a7e14dcfSSatish Balay 
333*a7e14dcfSSatish Balay   ierr = TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch); CHKERRQ(ierr);
334*a7e14dcfSSatish Balay   ierr = TaoLineSearchSetType(tao->linesearch, armijo_type); CHKERRQ(ierr);
335*a7e14dcfSSatish Balay   ierr = TaoLineSearchSetFromOptions(tao->linesearch); CHKERRQ(ierr);
336*a7e14dcfSSatish Balay 
337*a7e14dcfSSatish Balay   ierr = KSPCreate(((PetscObject)tao)->comm, &tao->ksp); CHKERRQ(ierr);
338*a7e14dcfSSatish Balay   ierr = KSPSetFromOptions(tao->ksp); CHKERRQ(ierr);
339*a7e14dcfSSatish Balay   tao->max_it = 2000;
340*a7e14dcfSSatish Balay   tao->max_funcs = 4000;
341*a7e14dcfSSatish Balay   tao->fatol = 0;
342*a7e14dcfSSatish Balay   tao->frtol = 0;
343*a7e14dcfSSatish Balay   tao->gttol = 0;
344*a7e14dcfSSatish Balay   tao->grtol = 0;
345*a7e14dcfSSatish Balay   tao->gatol = 1.0e-16;
346*a7e14dcfSSatish Balay   tao->fmin = 1.0e-8;
347*a7e14dcfSSatish Balay 
348*a7e14dcfSSatish Balay   PetscFunctionReturn(0);
349*a7e14dcfSSatish Balay }
350*a7e14dcfSSatish Balay EXTERN_C_END
351*a7e14dcfSSatish Balay 
352