xref: /petsc/src/tao/complementarity/impls/ssls/ssls.h (revision 6d8694c4fbab79f9439f1ad13c0386ba7ee1ca4b)
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. See {cite}`billups:algorithms`, {cite}`deluca.facchinei.ea:semismooth`,
34   {cite}`ferris.kanzow.ea:feasible`, {cite}`fischer:special`, and {cite}`munson.facchinei.ea:semismooth`.
35 */
36 
37 #pragma once
38 #include <petsc/private/taoimpl.h>
39 
40 typedef struct {
41   Vec ff;   /* fischer function */
42   Vec dpsi; /* gradient of psi */
43 
44   Vec da; /* work vector for subdifferential calculation (diag pert) */
45   Vec db; /* work vector for subdifferential calculation (row scale) */
46   Vec dm; /* work vector for subdifferential calculation (mu vector) */
47   Vec dxfree;
48 
49   Vec t1; /* work vector */
50   Vec t2; /* work vector */
51 
52   Vec r1, r2, r3, w; /* work vectors */
53 
54   PetscReal merit; /* merit function value (norm(fischer)) */
55   PetscReal merit_eqn;
56   PetscReal merit_mu;
57 
58   PetscReal delta;
59   PetscReal rho;
60 
61   PetscReal rtol; /* Solution tolerances */
62   PetscReal atol;
63 
64   PetscReal identifier; /* Active-set identification */
65 
66   /* Interior-point method data */
67   PetscReal mu_init; /* initial smoothing parameter value */
68   PetscReal mu;      /* smoothing parameter */
69   PetscReal dmu;     /* direction in smoothing parameter */
70   PetscReal mucon;   /* smoothing parameter constraint */
71   PetscReal d_mucon; /* derivative of smoothing constraint with respect to mu */
72   PetscReal g_mucon; /* gradient of merit function with respect to mu */
73 
74   Mat J_sub, Jpre_sub; /* subset of jacobian */
75   Vec f;               /* constraint function */
76 
77   IS fixed;
78   IS free;
79 } TAO_SSLS;
80 
81 PetscErrorCode TaoSetFromOptions_SSLS(Tao, PetscOptionItems);
82 PetscErrorCode TaoView_SSLS(Tao, PetscViewer);
83 
84 PetscErrorCode Tao_SSLS_Function(TaoLineSearch, Vec, PetscReal *, void *);
85 PetscErrorCode Tao_SSLS_FunctionGradient(TaoLineSearch, Vec, PetscReal *, Vec, void *);
86