Actual source code: ssls.h

  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

 10:  Many other reformulations are possible and combinations of
 11:  feasible/infeasible and linesearch/trust region are possible.

 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.

 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.

 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.

 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: */

 37: #pragma once
 38: #include <petsc/private/taoimpl.h>

 40: typedef struct {
 41:   Vec ff;   /* fischer function */
 42:   Vec dpsi; /* gradient of psi */

 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;

 49:   Vec t1; /* work vector */
 50:   Vec t2; /* work vector */

 52:   Vec r1, r2, r3, w; /* work vectors */

 54:   PetscReal merit; /* merit function value (norm(fischer)) */
 55:   PetscReal merit_eqn;
 56:   PetscReal merit_mu;

 58:   PetscReal delta;
 59:   PetscReal rho;

 61:   PetscReal rtol; /* Solution tolerances */
 62:   PetscReal atol;

 64:   PetscReal identifier; /* Active-set identification */

 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 */

 74:   Mat J_sub, Jpre_sub; /* subset of jacobian */
 75:   Vec f;               /* constraint function */

 77:   IS fixed;
 78:   IS free;
 79: } TAO_SSLS;

 81: PetscErrorCode TaoSetFromOptions_SSLS(Tao, PetscOptionItems *);
 82: PetscErrorCode TaoView_SSLS(Tao, PetscViewer);

 84: PetscErrorCode Tao_SSLS_Function(TaoLineSearch, Vec, PetscReal *, void *);
 85: PetscErrorCode Tao_SSLS_FunctionGradient(TaoLineSearch, Vec, PetscReal *, Vec, void *);