1 /* Context for SSXLS
2    -- semismooth (SS) - function not differentiable
3                       - merit function continuously differentiable
4                       - Fischer-Burmeister reformulation of complementarity
5                         - Billups composition for two finite bounds
6    -- infeasible (I)  - iterates not guaranteed to remain within bounds
7    -- feasible (F)    - iterates guaranteed to remain within bounds
8    -- linesearch (LS) - Armijo rule on direction
9 
10  Many other reformulations are possible and combinations of
11  feasible/infeasible and linesearch/trust region are possible.
12 
13  Basic theory
14    Fischer-Burmeister reformulation is semismooth with a continuously
15    differentiable merit function and strongly semismooth if the F has
16    lipschitz continuous derivatives.
17 
18    Every accumulation point generated by the algorithm is a stationary
19    point for the merit function.  Stationary points of the merit function
20    are solutions of the complementarity problem if
21      a.  the stationary point has a BD-regular subdifferential, or
22      b.  the Schur complement F'/F'_ff is a P_0-matrix where ff is the
23          index set corresponding to the free variables.
24 
25    If one of the accumulation points has a BD-regular subdifferential then
26      a.  the entire sequence converges to this accumulation point at
27          a local q-superlinear rate
28      b.  if in addition the reformulation is strongly semismooth near
29          this accumulation point, then the algorithm converges at a
30          local q-quadratic rate.
31 
32  The theory for the feasible version follows from the feasible descent
33  algorithm framework.
34 
35  References:
36    Billups, "Algorithms for Complementarity Problems and Generalized
37      Equations," Ph.D thesis, University of Wisconsin - Madison, 1995.
38    De Luca, Facchinei, Kanzow, "A Semismooth Equation Approach to the
39      Solution of Nonlinear Complementarity Problems," Mathematical
40      Programming, 75, pages 407-439, 1996.
41    Ferris, Kanzow, Munson, "Feasible Descent Algorithms for Mixed
42      Complementarity Problems," Mathematical Programming, 86,
43      pages 475-497, 1999.
44    Fischer, "A Special Newton-type Optimization Method," Optimization,
45      24, pages 269-284, 1992
46    Munson, Facchinei, Ferris, Fischer, Kanzow, "The Semismooth Algorithm
47      for Large Scale Complementarity Problems," Technical Report 99-06,
48      University of Wisconsin - Madison, 1999.
49 */
50 
51 #ifndef __TAO_SSLS_H
52 #define __TAO_SSLS_H
53 #include <petsc/private/taoimpl.h>
54 
55 typedef struct {
56   Vec ff;       /* fischer function */
57   Vec dpsi;     /* gradient of psi */
58 
59   Vec da;       /* work vector for subdifferential calculation (diag pert) */
60   Vec db;       /* work vector for subdifferential calculation (row scale) */
61   Vec dm;   /* work vector for subdifferential calculation (mu vector) */
62   Vec dxfree;
63 
64   Vec t1;       /* work vector */
65   Vec t2;       /* work vector */
66 
67   Vec r1,r2,r3,w; /* work vectors */
68 
69   PetscReal merit; /* merit function value (norm(fischer)) */
70   PetscReal merit_eqn;
71   PetscReal merit_mu;
72 
73   PetscReal delta;
74   PetscReal rho;
75 
76   PetscReal rtol;       /* Solution tolerances */
77   PetscReal atol;
78 
79   PetscReal identifier; /* Active-set identification */
80 
81   /* Interior-point method data */
82   PetscReal mu_init; /* initial smoothing parameter value */
83   PetscReal mu;      /* smoothing parameter */
84   PetscReal dmu;     /* direction in smoothing parameter */
85   PetscReal mucon;   /* smoothing parameter constraint */
86   PetscReal d_mucon; /* derivative of smoothing constraint with respect to mu */
87   PetscReal g_mucon; /* gradient of merit function with respect to mu */
88 
89   Mat J_sub, Jpre_sub; /* subset of jacobian */
90   Vec f;        /* constraint function */
91 
92   IS fixed;
93   IS free;
94 } TAO_SSLS;
95 
96 PetscErrorCode TaoSetFromOptions_SSLS(PetscOptionItems *,Tao);
97 PetscErrorCode TaoView_SSLS(Tao,PetscViewer);
98 
99 PetscErrorCode Tao_SSLS_Function(TaoLineSearch, Vec, PetscReal *, void *);
100 PetscErrorCode Tao_SSLS_FunctionGradient(TaoLineSearch, Vec, PetscReal *, Vec, void *);
101 
102 #endif
103 
104