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